30 #ifndef KOKKOS_EXPERIMENTAL_VIEW_SACADO_FAD_CONTIGUOUS_HPP
31 #define KOKKOS_EXPERIMENTAL_VIEW_SACADO_FAD_CONTIGUOUS_HPP
34 #if defined(HAVE_SACADO_KOKKOSCORE)
42 template <
typename ViewType,
typename Enabled =
void>
43 struct ThreadLocalScalarType {
44 typedef typename ViewType::non_const_value_type type;
47 template <
typename ViewType>
48 struct ViewScalarStride {
49 static const unsigned stride =
50 Impl::LayoutScalarStride< typename ViewType::array_layout>::stride;
51 static const bool is_unit_stride =
52 Impl::LayoutScalarStride< typename ViewType::array_layout>::is_unit_stride;
62 template <
unsigned Size = 0>
64 static const unsigned PartitionSize = Size;
68 template<
typename iType0 ,
typename iType1 >
70 Partition(
const iType0 & i0 ,
const iType1 & i1 ) :
71 offset(i0), stride(i1) {
76 struct is_fad_partition {
77 static const bool value =
false;
80 template <
unsigned Str
ide>
81 struct is_fad_partition< Partition<Stride> > {
82 static const bool value =
true;
88 template <
typename T,
unsigned Str
ide = 0>
89 struct LocalScalarType {
92 template <
typename T,
unsigned Str
ide>
93 struct LocalScalarType<const
T, Stride> {
94 typedef typename LocalScalarType<T,Stride>::type lst;
95 typedef const lst type;
106 template <
typename T,
int N>
class StaticStorage;
107 template <
typename S>
class GeneralFad;
110 template <
typename T,
int N,
unsigned Str
ide>
111 struct LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >,
113 static const int Ns = (N+Stride-1) / Stride;
114 typedef Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,Ns> > type;
116 #ifndef SACADO_NEW_FAD_DESIGN_IS_DEFAULT
118 template <
typename T,
int N>
class SLFad;
120 template <
typename T,
int N,
unsigned Str
ide>
121 struct LocalScalarType< Fad::
SLFad<T,N>, Stride > {
122 static const int Ns = (N+Stride-1) / Stride;
123 typedef Fad::SLFad<T,Ns> type;
133 template <
typename T,
int N>
class StaticFixedStorage;
134 template <
typename T,
int N>
class StaticStorage;
135 template <
typename S>
class GeneralFad;
138 template <
typename T,
int N,
unsigned Str
ide>
139 struct LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >,
141 static const int Ns = (N+Stride-1) / Stride;
142 typedef typename std::conditional<
144 Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,Ns> > ,
145 Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,Ns> >
149 #ifndef SACADO_NEW_FAD_DESIGN_IS_DEFAULT
151 template <
typename T,
int N>
class SFad;
153 template <
typename T,
int N,
unsigned Str
ide>
154 struct LocalScalarType< Fad::
SFad<T,N>, Stride > {
155 static const int Ns = (N+Stride-1) / Stride;
156 typedef typename std::conditional< Ns == N/Stride , Fad::SFad<T,Ns> , Fad::SLFad<T,Ns> >::type type;
160 template <
unsigned Str
ide,
typename T>
162 const T& partition_scalar(
const T& x) {
return x; }
167 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
170 #include "Kokkos_Core.hpp"
171 #include "impl/Kokkos_ViewMapping.hpp"
177 #if defined(__CUDA_ARCH__)
180 template <
typename T,
typename U>
class DynamicStorage;
181 template <
typename T,
int N>
class StaticFixedStorage;
182 template <
typename T,
int N>
class StaticStorage;
183 template <
typename S>
class GeneralFad;
186 #ifndef SACADO_VIEW_CUDA_HIERARCHICAL_DFAD
187 template <
unsigned Str
ide,
typename T,
typename U>
189 typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::DynamicStorage<T,U> >, Stride >::type
190 partition_scalar(
const Fad::Exp::GeneralFad< Fad::Exp::DynamicStorage<T,U> >& x) {
191 typedef typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::DynamicStorage<T,U> >, Stride >::type ret_type;
192 const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
193 const int offset = threadIdx.x;
194 ret_type xp(size, x.val());
200 const T*
dx = x.dx();
201 for (
int i=0; i<size; ++i)
202 xp.fastAccessDx(i) = dx[offset+i*Stride];
207 template <
unsigned Str
ide,
typename T,
int N>
209 typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >, Stride >::type
210 partition_scalar(
const Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >& x) {
211 typedef typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >, Stride >::type ret_type;
212 const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
213 const int offset = threadIdx.x;
214 ret_type xp(size, x.val());
215 for (
int i=0; i<size; ++i)
216 xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
219 template <
unsigned Str
ide,
typename T,
int N>
221 typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >, Stride >::type
222 partition_scalar(
const Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >& x) {
223 typedef typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >, Stride >::type ret_type;
224 const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
225 const int offset = threadIdx.x;
226 ret_type xp(size, x.val());
227 for (
int i=0; i<size; ++i)
228 xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
232 #ifndef SACADO_NEW_FAD_DESIGN_IS_DEFAULT
234 template <
typename T>
class DFad;
235 template <
typename T,
int N>
class SLFad;
236 template <
typename T,
int N>
class SFad;
238 #ifndef SACADO_VIEW_CUDA_HIERARCHICAL_DFAD
239 template <
unsigned Str
ide,
typename T>
241 typename LocalScalarType< Fad::DFad<T>, Stride >::type
242 partition_scalar(
const Fad::DFad<T>& x) {
243 typedef typename LocalScalarType< Fad::DFad<T>, Stride >::type ret_type;
244 const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
245 const int offset = threadIdx.x;
246 ret_type xp(size, x.val());
252 const T* dx = x.dx();
253 for (
int i=0; i<size; ++i)
254 xp.fastAccessDx(i) = dx[offset+i*Stride];
259 template <
unsigned Str
ide,
typename T,
int N>
261 typename LocalScalarType< Fad::SLFad<T,N>, Stride >::type
262 partition_scalar(
const Fad::SLFad<T,N>& x) {
263 typedef typename LocalScalarType< Fad::SLFad<T,N>, Stride >::type ret_type;
264 const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
265 const int offset = threadIdx.x;
266 ret_type xp(size, x.val());
267 for (
int i=0; i<size; ++i)
268 xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
271 template <
unsigned Str
ide,
typename T,
int N>
273 typename LocalScalarType< Fad::SFad<T,N>, Stride >::type
274 partition_scalar(
const Fad::SFad<T,N>& x) {
275 typedef typename LocalScalarType< Fad::SFad<T,N>, Stride >::type ret_type;
276 const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
277 const int offset = threadIdx.x;
278 ret_type xp(size, x.val());
279 for (
int i=0; i<size; ++i)
280 xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
283 #endif // SACADO_NEW_FAD_DESIGN_IS_DEFAULT
285 #endif // defined(__CUDA_ARCH__)
293 template<
unsigned Stride,
typename D,
typename ... P >
295 typename Kokkos::Impl::ViewMapping<
void,
typename Kokkos::ViewTraits<
D,P...>, Sacado::Fad::Partition<Stride> >::type
296 partition(
const Kokkos::View<D,P...> & src ,
297 const unsigned offset ,
298 const unsigned stride )
300 typedef Kokkos::ViewTraits<
D,P...> traits;
301 typedef typename Kokkos::Impl::ViewMapping< void, traits, Sacado::Fad::Partition<Stride> >::type DstViewType;
302 const Sacado::Fad::Partition<Stride> part( offset , stride );
303 return DstViewType(src, part);
306 template <
typename ViewType>
307 struct ThreadLocalScalarType<
309 typename std::enable_if< is_view_fad_contiguous<ViewType>::value >::type > {
310 typedef typename ViewType::traits TraitsType;
311 typedef Impl::ViewMapping<TraitsType, typename TraitsType::specialize> MappingType;
312 typedef typename MappingType::thread_local_scalar_type type;
317 #if defined (KOKKOS_ENABLE_CUDA) && defined(SACADO_VIEW_CUDA_HIERARCHICAL)
318 template<
class OutputView >
319 struct SacadoViewFill<
321 typename std::enable_if<
322 ( Kokkos::is_view_fad_contiguous<OutputView>::value &&
323 std::is_same<typename OutputView::execution_space, Kokkos::Cuda>::value &&
324 !Kokkos::ViewScalarStride<OutputView>::is_unit_stride )
328 typedef typename OutputView::const_value_type const_value_type ;
329 typedef typename OutputView::execution_space execution_space ;
330 typedef Kokkos::TeamPolicy< execution_space> team_policy;
331 typedef typename team_policy::member_type team_impl_handle;
332 typedef typename Kokkos::ThreadLocalScalarType<OutputView>::type local_scalar_type;
333 static const unsigned stride = Kokkos::ViewScalarStride<OutputView>::stride;
335 const OutputView output ;
336 const_value_type input ;
339 void operator()(
const size_t i0 )
const
341 local_scalar_type input_stride = Sacado::partition_scalar<stride>(input);
343 const size_t n1 = output.extent(1);
344 const size_t n2 = output.extent(2);
345 const size_t n3 = output.extent(3);
346 const size_t n4 = output.extent(4);
347 const size_t n5 = output.extent(5);
348 const size_t n6 = output.extent(6);
349 const size_t n7 = output.extent(7);
351 for (
size_t i1 = 0 ; i1 < n1 ; ++i1 ) {
352 for (
size_t i2 = 0 ; i2 < n2 ; ++i2 ) {
353 for (
size_t i3 = 0 ; i3 < n3 ; ++i3 ) {
354 for (
size_t i4 = 0 ; i4 < n4 ; ++i4 ) {
355 for (
size_t i5 = 0 ; i5 < n5 ; ++i5 ) {
356 for (
size_t i6 = 0 ; i6 < n6 ; ++i6 ) {
357 for (
size_t i7 = 0 ; i7 < n7 ; ++i7 ) {
358 output.access(i0,i1,i2,i3,i4,i5,i6,i7) = input_stride ;
363 void operator()(
const team_impl_handle& team )
const
365 const size_t i0 = team.league_rank()*team.team_size() + team.team_rank();
366 if (i0 < output.extent(0))
370 SacadoViewFill(
const OutputView & arg_out , const_value_type & arg_in )
371 : output( arg_out ), input( arg_in )
373 const size_t team_size = 256 / stride;
374 team_policy policy( (output.extent(0)+team_size-1)/team_size ,
375 team_size , stride );
376 Kokkos::parallel_for( policy, *
this );
377 execution_space::fence();
390 template<
class ... Args >
391 struct is_ViewSpecializeSacadoFadContiguous {
enum { value =
false }; };
393 template<
class D ,
class ... P ,
class ... Args >
394 struct is_ViewSpecializeSacadoFadContiguous< Kokkos::
View<D,P...> , Args... > {
396 std::is_same<
typename Kokkos::ViewTraits<
D,P...>::specialize
397 , ViewSpecializeSacadoFadContiguous >::value
399 ( (
sizeof...(Args) == 0 ) ||
400 is_ViewSpecializeSacadoFadContiguous< Args... >::value ) };
414 constexpr
unsigned computeFadPartitionSize(
unsigned size,
unsigned stride)
417 ((size+stride-1)/stride) == (size/stride) ? ((size+stride-1)/stride) : 0;
423 template <
unsigned rank,
unsigned static_dim,
typename Layout>
425 typename std::enable_if< !std::is_same<Layout, LayoutLeft>::value &&
426 !std::is_same<Layout, LayoutStride>::value,
428 create_fad_array_layout(
const Layout& layout)
431 for (
int i=0; i<8; ++i)
432 dims[i] = layout.dimension[i];
434 dims[rank] = static_dim+1;
435 return Layout( dims[0], dims[1], dims[2], dims[3],
436 dims[4], dims[5], dims[6], dims[7] );
442 template <
unsigned rank,
unsigned static_dim,
typename Layout>
444 typename std::enable_if< std::is_same<Layout, LayoutStride>::value, Layout>::type
445 create_fad_array_layout(
const Layout& layout)
447 size_t dims[8], strides[8];
448 for (
int i=0; i<8; ++i) {
449 dims[i] = layout.dimension[i];
450 strides[i] = layout.stride[i];
452 if (static_dim > 0) {
453 dims[rank] = static_dim+1;
456 return Layout( dims[0], strides[0],
463 dims[7], strides[7] );
469 template <
unsigned rank,
unsigned static_dim,
typename Layout>
471 typename std::enable_if< std::is_same<Layout, LayoutLeft>::value, Layout >::type
472 create_fad_array_layout(
const Layout& layout)
475 for (
int i=0; i<8; ++i)
476 dims[i] = layout.dimension[i];
477 size_t fad_dim = static_dim == 0 ? dims[rank] : static_dim+1;
478 for (
int i=rank; i>=1; --i)
481 return Layout( dims[0], dims[1], dims[2], dims[3],
482 dims[4], dims[5], dims[6], dims[7] );
485 template <
unsigned Rank,
typename Dimension,
typename Layout>
487 typename std::enable_if< !std::is_same<Layout, LayoutLeft>::value,
size_t>::type
488 getFadDimension(
const ViewOffset<Dimension,Layout,void>& offset)
491 ( Rank == 0 ? offset.dimension_0() :
492 ( Rank == 1 ? offset.dimension_1() :
493 ( Rank == 2 ? offset.dimension_2() :
494 ( Rank == 3 ? offset.dimension_3() :
495 ( Rank == 4 ? offset.dimension_4() :
496 ( Rank == 5 ? offset.dimension_5() :
497 ( Rank == 6 ? offset.dimension_6() :
498 offset.dimension_7() )))))));
501 template <
unsigned Rank,
typename Dimension,
typename Layout>
503 typename std::enable_if< std::is_same<Layout, LayoutLeft>::value,
size_t >::type
504 getFadDimension(
const ViewOffset<Dimension,Layout,void>& offset)
506 return offset.dimension_0();
509 template<
class Traits >
510 class ViewMapping< Traits ,
511 typename std::enable_if<
512 ( std::is_same< typename Traits::specialize
513 , ViewSpecializeSacadoFadContiguous >::value
515 ( std::is_same< typename Traits::array_layout
516 , Kokkos::LayoutLeft >::value
518 std::is_same< typename Traits::array_layout
519 , Kokkos::LayoutRight >::value
521 std::is_same< typename Traits::array_layout
522 , Kokkos::LayoutStride >::value
525 , typename Traits::specialize
530 template< class ,
class ... >
friend class ViewMapping ;
531 template< class ,
class ... >
friend class Kokkos::View ;
533 typedef typename Traits::value_type fad_type ;
536 std::add_const< fad_value_type >::type const_fad_value_type ;
541 enum { PartitionedFadStride = Traits::array_layout::scalar_stride };
545 enum { PartitionedFadStaticDimension =
546 computeFadPartitionSize(FadStaticDimension,PartitionedFadStride) };
548 #ifdef KOKKOS_ENABLE_CUDA
549 typedef typename Sacado::LocalScalarType< fad_type, unsigned(PartitionedFadStride) >::type strided_scalar_type;
550 typedef typename std::conditional< std::is_same<typename Traits::execution_space, Kokkos::Cuda>::value, strided_scalar_type, fad_type >::type thread_local_scalar_type;
552 typedef fad_type thread_local_scalar_type;
558 typedef fad_value_type * handle_type ;
560 typedef ViewArrayAnalysis< typename Traits::data_type > array_analysis ;
562 typedef ViewOffset<
typename Traits::dimension
563 ,
typename Traits::array_layout
568 static constexpr
bool is_layout_left =
569 std::is_same< typename Traits::array_layout, Kokkos::LayoutLeft>::value;
571 typename std::conditional<
573 typename array_analysis::dimension::
574 template prepend<0>::type,
575 typename array_analysis::dimension::
576 template append<0>::type >::type,
577 typename Traits::array_layout,
581 handle_type m_impl_handle ;
582 offset_type m_impl_offset ;
583 array_offset_type m_array_offset ;
584 sacado_size_type m_fad_size ;
587 unsigned m_original_fad_size ;
588 unsigned m_fad_stride ;
589 unsigned m_fad_index ;
596 enum { Rank = Traits::dimension::rank };
599 enum { Sacado_Rank = std::is_same< typename Traits::array_layout, Kokkos::LayoutLeft >::value ? 0 : Rank+1 };
602 template<
typename iType >
604 {
return m_impl_offset.m_dim.extent(r); }
607 typename Traits::array_layout layout()
const
608 {
return m_impl_offset.layout(); }
611 {
return m_impl_offset.dimension_0(); }
613 {
return m_impl_offset.dimension_1(); }
615 {
return m_impl_offset.dimension_2(); }
617 {
return m_impl_offset.dimension_3(); }
619 {
return m_impl_offset.dimension_4(); }
621 {
return m_impl_offset.dimension_5(); }
623 {
return m_impl_offset.dimension_6(); }
625 {
return m_impl_offset.dimension_7(); }
630 using is_regular = std::false_type ;
634 {
return m_impl_offset.stride_0(); }
636 {
return m_impl_offset.stride_1(); }
638 {
return m_impl_offset.stride_2(); }
640 {
return m_impl_offset.stride_3(); }
642 {
return m_impl_offset.stride_4(); }
644 {
return m_impl_offset.stride_5(); }
646 {
return m_impl_offset.stride_6(); }
648 {
return m_impl_offset.stride_7(); }
650 template<
typename iType >
652 { m_impl_offset.stride(s); }
656 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
657 {
return PartitionedFadStaticDimension ? PartitionedFadStaticDimension+1 : (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x + 1; }
659 {
return m_fad_size.value+1; }
664 {
return m_fad_stride; }
669 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
683 typedef fad_value_type * pointer_type ;
687 {
return m_array_offset.span(); }
691 {
return m_array_offset.span_is_contiguous() && (m_fad_stride == 1); }
695 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
696 {
return m_impl_handle + threadIdx.x; }
698 {
return m_impl_handle + m_fad_index; }
707 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
708 const unsigned index = threadIdx.x;
709 const unsigned strd = blockDim.x;
710 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
712 const unsigned index = m_fad_index;
713 const unsigned strd = m_fad_stride;
714 const unsigned size = m_fad_size.value;
716 return reference_type( m_impl_handle + index
717 , m_impl_handle + m_original_fad_size
721 template<
typename I0 >
723 typename std::enable_if< Kokkos::Impl::are_integral<I0>::value &&
724 is_layout_left, reference_type>::type
725 reference(
const I0 & i0 )
const
726 { pointer_type beg = m_impl_handle + m_array_offset(0,i0);
727 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
728 const unsigned index = threadIdx.x;
729 const unsigned strd = blockDim.x;
730 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
732 const unsigned index = m_fad_index;
733 const unsigned strd = m_fad_stride;
734 const unsigned size = m_fad_size.value;
736 return reference_type( beg + index
737 , beg + m_original_fad_size
741 template<
typename I0 >
743 typename std::enable_if< Kokkos::Impl::are_integral<I0>::value &&
744 !is_layout_left, reference_type>::type
745 reference(
const I0 & i0 )
const
746 { pointer_type beg = m_impl_handle + m_array_offset(i0,0);
747 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
748 const unsigned index = threadIdx.x;
749 const unsigned strd = blockDim.x;
750 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
752 const unsigned index = m_fad_index;
753 const unsigned strd = m_fad_stride;
754 const unsigned size = m_fad_size.value;
756 return reference_type( beg + index
757 , beg + m_original_fad_size
761 template<
typename I0 ,
typename I1 >
763 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1>::value &&
764 is_layout_left, reference_type>::type
765 reference(
const I0 & i0 ,
const I1 & i1 )
const
766 { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1);
767 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
768 const unsigned index = threadIdx.x;
769 const unsigned strd = blockDim.x;
770 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
772 const unsigned index = m_fad_index;
773 const unsigned strd = m_fad_stride;
774 const unsigned size = m_fad_size.value;
776 return reference_type( beg + index
777 , beg + m_original_fad_size
781 template<
typename I0 ,
typename I1 >
783 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1>::value &&
784 !is_layout_left, reference_type>::type
785 reference(
const I0 & i0 ,
const I1 & i1 )
const
786 { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,0);
787 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
788 const unsigned index = threadIdx.x;
789 const unsigned strd = blockDim.x;
790 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
792 const unsigned index = m_fad_index;
793 const unsigned strd = m_fad_stride;
794 const unsigned size = m_fad_size.value;
796 return reference_type( beg + index
797 , beg + m_original_fad_size
802 template<
typename I0 ,
typename I1 ,
typename I2 >
804 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2>::value &&
805 is_layout_left, reference_type>::type
806 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 )
const
807 { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2);
808 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
809 const unsigned index = threadIdx.x;
810 const unsigned strd = blockDim.x;
811 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
813 const unsigned index = m_fad_index;
814 const unsigned strd = m_fad_stride;
815 const unsigned size = m_fad_size.value;
817 return reference_type( beg + index
818 , beg + m_original_fad_size
822 template<
typename I0 ,
typename I1 ,
typename I2 >
824 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2>::value &&
825 !is_layout_left, reference_type>::type
826 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 )
const
827 { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,0);
828 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
829 const unsigned index = threadIdx.x;
830 const unsigned strd = blockDim.x;
831 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
833 const unsigned index = m_fad_index;
834 const unsigned strd = m_fad_stride;
835 const unsigned size = m_fad_size.value;
837 return reference_type( beg + index
838 , beg + m_original_fad_size
842 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3 >
844 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3>::value &&
845 is_layout_left, reference_type>::type
846 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3 )
const
847 { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3);
848 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
849 const unsigned index = threadIdx.x;
850 const unsigned strd = blockDim.x;
851 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
853 const unsigned index = m_fad_index;
854 const unsigned strd = m_fad_stride;
855 const unsigned size = m_fad_size.value;
857 return reference_type( beg + index
858 , beg + m_original_fad_size
862 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3 >
864 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3>::value &&
865 !is_layout_left, reference_type>::type
866 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3 )
const
867 { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,0);
868 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
869 const unsigned index = threadIdx.x;
870 const unsigned strd = blockDim.x;
871 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
873 const unsigned index = m_fad_index;
874 const unsigned strd = m_fad_stride;
875 const unsigned size = m_fad_size.value;
877 return reference_type( beg + index
878 , beg + m_original_fad_size
882 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3
885 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4>::value &&
886 is_layout_left, reference_type>::type
887 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3
888 ,
const I4 & i4 )
const
889 { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4);
890 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
891 const unsigned index = threadIdx.x;
892 const unsigned strd = blockDim.x;
893 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
895 const unsigned index = m_fad_index;
896 const unsigned strd = m_fad_stride;
897 const unsigned size = m_fad_size.value;
899 return reference_type( beg + index
900 , beg + m_original_fad_size
904 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3
907 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4>::value &&
908 !is_layout_left, reference_type>::type
909 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3
910 ,
const I4 & i4 )
const
911 { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,0);
912 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
913 const unsigned index = threadIdx.x;
914 const unsigned strd = blockDim.x;
915 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
917 const unsigned index = m_fad_index;
918 const unsigned strd = m_fad_stride;
919 const unsigned size = m_fad_size.value;
921 return reference_type( beg + index
922 , beg + m_original_fad_size
926 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3
927 ,
typename I4 ,
typename I5 >
929 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5>::value &&
930 is_layout_left, reference_type>::type
931 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3
932 ,
const I4 & i4 ,
const I5 & i5 )
const
933 { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5);
934 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
935 const unsigned index = threadIdx.x;
936 const unsigned strd = blockDim.x;
937 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
939 const unsigned index = m_fad_index;
940 const unsigned strd = m_fad_stride;
941 const unsigned size = m_fad_size.value;
943 return reference_type( beg + index
944 , beg + m_original_fad_size
948 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3
949 ,
typename I4 ,
typename I5 >
951 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5>::value &&
952 !is_layout_left, reference_type>::type
953 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3
954 ,
const I4 & i4 ,
const I5 & i5 )
const
955 { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,0);
956 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
957 const unsigned index = threadIdx.x;
958 const unsigned strd = blockDim.x;
959 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
961 const unsigned index = m_fad_index;
962 const unsigned strd = m_fad_stride;
963 const unsigned size = m_fad_size.value;
965 return reference_type( beg + index
966 , beg + m_original_fad_size
970 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3
971 ,
typename I4 ,
typename I5 ,
typename I6 >
973 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5,I6>::value &&
974 is_layout_left, reference_type>::type
975 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3
976 ,
const I4 & i4 ,
const I5 & i5 ,
const I6 & i6 )
const
977 { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5,i6);
978 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
979 const unsigned index = threadIdx.x;
980 const unsigned strd = blockDim.x;
981 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
983 const unsigned index = m_fad_index;
984 const unsigned strd = m_fad_stride;
985 const unsigned size = m_fad_size.value;
987 return reference_type( beg + index
988 , beg + m_original_fad_size
992 template<
typename I0 ,
typename I1 ,
typename I2 ,
typename I3
993 ,
typename I4 ,
typename I5 ,
typename I6 >
995 typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5,I6>::value &&
996 !is_layout_left, reference_type>::type
997 reference(
const I0 & i0 ,
const I1 & i1 ,
const I2 & i2 ,
const I3 & i3
998 ,
const I4 & i4 ,
const I5 & i5 ,
const I6 & i6 )
const
999 { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,i6,0);
1000 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
1001 const unsigned index = threadIdx.x;
1002 const unsigned strd = blockDim.x;
1003 const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1005 const unsigned index = m_fad_index;
1006 const unsigned strd = m_fad_stride;
1007 const unsigned size = m_fad_size.value;
1009 return reference_type( beg + index
1010 , beg + m_original_fad_size
1018 static constexpr
size_t memory_span(
typename Traits::array_layout
const & layout )
1021 typedef std::integral_constant< unsigned , 0 > padding ;
1022 return array_offset_type(
1024 create_fad_array_layout<
unsigned(Rank),
unsigned(FadStaticDimension)>( layout ) ).span() *
sizeof(fad_value_type);
1030 KOKKOS_INLINE_FUNCTION ViewMapping() : m_impl_handle(0) , m_impl_offset() , m_array_offset() , m_fad_size(0) , m_original_fad_size(0) , m_fad_stride(1) , m_fad_index(0) {}
1038 template<
class ... P >
1041 ( ViewCtorProp< P ... >
const & prop
1042 ,
typename Traits::array_layout
const & local_layout
1044 : m_impl_handle( ( (ViewCtorProp<
void,pointer_type> const &) prop ).value )
1045 , m_impl_offset( std::integral_constant< unsigned , 0 >()
1048 std::integral_constant< unsigned , 0 >()
1049 , create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( local_layout ) )
1050 , m_fad_size( getFadDimension<unsigned(Rank)>( m_array_offset ) - 1 )
1051 , m_original_fad_size( m_fad_size.value )
1055 const unsigned fad_dim =
1056 getFadDimension<unsigned(Rank)>( m_array_offset );
1057 if (
unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1058 Kokkos::abort(
"invalid fad dimension (0) supplied!");
1066 template<
class ... P >
1067 SharedAllocationRecord<> *
1068 allocate_shared( ViewCtorProp< P... >
const & prop
1069 ,
typename Traits::array_layout
const & local_layout )
1071 typedef ViewCtorProp< P... > ctor_prop ;
1073 typedef typename ctor_prop::execution_space execution_space ;
1074 typedef typename Traits::memory_space memory_space ;
1075 typedef ViewValueFunctor< execution_space , fad_value_type > functor_type ;
1076 typedef SharedAllocationRecord< memory_space , functor_type > record_type ;
1079 typedef std::integral_constant< unsigned , 0 > padding ;
1082 enum { test_traits_check = Kokkos::Impl::check_has_common_view_alloc_prop< P... >::value };
1084 typename Traits::array_layout internal_layout =
1085 (test_traits_check ==
true)
1086 ? Kokkos::Impl::appendFadToLayoutViewAllocHelper< Traits, P... >::returnNewLayoutPlusFad(prop, local_layout)
1089 m_impl_offset = offset_type( padding(), internal_layout );
1092 array_offset_type( padding() ,
1093 create_fad_array_layout<
unsigned(Rank),
unsigned(FadStaticDimension)>( internal_layout ) );
1094 const unsigned fad_dim =
1095 getFadDimension<unsigned(Rank)>( m_array_offset );
1096 if (
unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1097 Kokkos::abort(
"invalid fad dimension (0) supplied!");
1098 m_fad_size = fad_dim - 1 ;
1099 m_original_fad_size = m_fad_size.value ;
1103 const size_t alloc_size = m_array_offset.span() *
sizeof(fad_value_type);
1106 record_type *
const record =
1107 record_type::allocate( ( (ViewCtorProp<void,memory_space>
const &) prop ).value
1108 , ( (ViewCtorProp<void,std::string>
const &) prop ).value
1115 m_impl_handle = handle_type( reinterpret_cast< pointer_type >( record->data() ) );
1117 if ( ctor_prop::initialize ) {
1120 record->m_destroy = functor_type( ( (ViewCtorProp<void,execution_space>
const &) prop).value
1121 , (fad_value_type *) m_impl_handle
1122 , m_array_offset.span()
1126 record->m_destroy.construct_shared_allocation();
1147 template<
class DstTraits ,
class SrcTraits >
1148 class ViewMapping< DstTraits , SrcTraits ,
1149 typename std::enable_if<(
1150 Kokkos::Impl::MemorySpaceAccess
1151 < typename DstTraits::memory_space
1152 , typename SrcTraits::memory_space >::assignable
1155 std::is_same< typename DstTraits::specialize
1156 , ViewSpecializeSacadoFadContiguous >::value
1159 std::is_same< typename SrcTraits::specialize
1160 , ViewSpecializeSacadoFadContiguous >::value
1162 , typename DstTraits::specialize
1167 enum { is_assignable =
true };
1169 typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1170 typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1171 typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1173 template<
class DstType >
1175 void assign( DstType & dst
1176 ,
const SrcFadType & src
1177 ,
const TrackType & )
1181 std::is_same<
typename DstTraits::array_layout
1182 , Kokkos::LayoutLeft >::value ||
1183 std::is_same<
typename DstTraits::array_layout
1184 , Kokkos::LayoutRight >::value ||
1185 std::is_same<
typename DstTraits::array_layout
1186 , Kokkos::LayoutStride >::value
1190 std::is_same<
typename SrcTraits::array_layout
1191 , Kokkos::LayoutLeft >::value ||
1192 std::is_same<
typename SrcTraits::array_layout
1193 , Kokkos::LayoutRight >::value ||
1194 std::is_same<
typename SrcTraits::array_layout
1195 , Kokkos::LayoutStride >::value
1197 ,
"View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1200 std::is_same<
typename DstTraits::array_layout
1201 ,
typename SrcTraits::array_layout >::value ||
1202 std::is_same<
typename DstTraits::array_layout
1203 , Kokkos::LayoutStride >::value ,
1204 "View assignment must have compatible layout" );
1207 std::is_same<
typename DstTraits::scalar_array_type
1208 ,
typename SrcTraits::scalar_array_type >::value ||
1209 std::is_same<
typename DstTraits::scalar_array_type
1210 ,
typename SrcTraits::const_scalar_array_type >::value ,
1211 "View assignment must have same value type or const = non-const" );
1214 ViewDimensionAssignable
1215 <
typename DstType::offset_type::dimension_type
1216 ,
typename SrcFadType::offset_type::dimension_type >::value ,
1217 "View assignment must have compatible dimensions" );
1220 ViewDimensionAssignable
1221 <
typename DstType::array_offset_type::dimension_type
1222 ,
typename SrcFadType::array_offset_type::dimension_type >::value ,
1223 "View assignment must have compatible dimensions" );
1225 typedef typename DstType::offset_type dst_offset_type ;
1226 typedef typename DstType::array_offset_type dst_array_offset_type ;
1228 dst.m_impl_handle = src.m_impl_handle ;
1229 dst.m_impl_offset = dst_offset_type( src.m_impl_offset );
1230 dst.m_array_offset = dst_array_offset_type( src.m_array_offset );
1231 dst.m_fad_size = src.m_fad_size.value ;
1232 dst.m_original_fad_size = src.m_original_fad_size ;
1233 dst.m_fad_stride = src.m_fad_stride ;
1234 dst.m_fad_index = src.m_fad_index ;
1242 template<
class DstTraits ,
class SrcTraits >
1243 class ViewMapping< DstTraits , SrcTraits ,
1244 typename std::enable_if<(
1245 std::is_same< typename DstTraits::memory_space
1246 , typename SrcTraits::memory_space >::value
1249 std::is_same< typename DstTraits::specialize
1250 , ViewSpecializeSacadoFad >::value
1253 std::is_same< typename SrcTraits::specialize
1254 , ViewSpecializeSacadoFadContiguous >::value
1257 std::is_same< typename DstTraits::array_layout
1258 , Kokkos::LayoutStride >::value
1260 , typename DstTraits::specialize
1265 enum { is_assignable =
true };
1267 typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1268 typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1269 typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1271 template<
class DstType >
1273 void assign( DstType & dst
1274 ,
const SrcFadType & src
1275 ,
const TrackType & )
1278 std::is_same<
typename SrcTraits::array_layout
1279 , Kokkos::LayoutLeft >::value ||
1280 std::is_same<
typename SrcTraits::array_layout
1281 , Kokkos::LayoutRight >::value ||
1282 std::is_same<
typename SrcTraits::array_layout
1283 , Kokkos::LayoutStride >::value ,
1284 "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1287 std::is_same<
typename DstTraits::value_type
1288 ,
typename SrcTraits::value_type >::value ||
1289 std::is_same<
typename DstTraits::value_type
1290 ,
typename SrcTraits::const_value_type >::value ,
1291 "View assignment must have same value type or const = non-const" );
1294 DstTraits::dimension::rank == SrcTraits::dimension::rank,
1295 "View assignment must have same rank" );
1297 typedef typename DstType::array_offset_type dst_offset_type ;
1299 dst.m_impl_handle = src.m_impl_handle ;
1300 dst.m_fad_size = src.m_fad_size.value ;
1301 dst.m_fad_stride = src.m_fad_stride ;
1302 dst.m_impl_offset = src.m_impl_offset;
1305 N[0] = src.m_array_offset.dimension_0();
1306 N[1] = src.m_array_offset.dimension_1();
1307 N[2] = src.m_array_offset.dimension_2();
1308 N[3] = src.m_array_offset.dimension_3();
1309 N[4] = src.m_array_offset.dimension_4();
1310 N[5] = src.m_array_offset.dimension_5();
1311 N[6] = src.m_array_offset.dimension_6();
1312 N[7] = src.m_array_offset.dimension_7();
1313 S[0] = src.m_array_offset.stride_0();
1314 S[1] = src.m_array_offset.stride_1();
1315 S[2] = src.m_array_offset.stride_2();
1316 S[3] = src.m_array_offset.stride_3();
1317 S[4] = src.m_array_offset.stride_4();
1318 S[5] = src.m_array_offset.stride_5();
1319 S[6] = src.m_array_offset.stride_6();
1320 S[7] = src.m_array_offset.stride_7();
1324 if (std::is_same<
typename SrcTraits::array_layout
1325 , Kokkos::LayoutLeft >::value)
1327 const size_t N_fad = N[0];
1328 const size_t S_fad = S[0];
1329 for (
int i=0; i<7; ++i) {
1333 N[DstTraits::dimension::rank] = N_fad;
1334 S[DstTraits::dimension::rank] = S_fad;
1336 Kokkos::LayoutStride ls( N[0], S[0],
1344 dst.m_array_offset = dst_offset_type(std::integral_constant<unsigned,0>(), ls);
1352 template<
class DstTraits ,
class SrcTraits >
1353 class ViewMapping< DstTraits , SrcTraits ,
1354 typename std::enable_if<(
1355 std::is_same< typename DstTraits::memory_space
1356 , typename SrcTraits::memory_space >::value
1359 std::is_same< typename DstTraits::specialize , void >::value
1362 std::is_same< typename SrcTraits::specialize
1363 , ViewSpecializeSacadoFadContiguous >::value
1365 , typename DstTraits::specialize
1370 enum { is_assignable =
true };
1372 typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1373 typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1374 typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1379 template <
class DstType,
class SrcFadType,
class Enable =
void >
1380 struct AssignOffset;
1382 template <
class DstType,
class SrcFadType >
1383 struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank != (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1387 static void assign( DstType & dst,
const SrcFadType & src )
1389 typedef typename SrcTraits::value_type TraitsValueType;
1396 typedef typename DstType::offset_type::array_layout DstLayoutType;
1398 typedef typename SrcFadType::array_offset_type::dimension_type SrcViewDimension;
1403 static constexpr
bool is_layout_left =
1404 std::is_same< DstLayoutType, Kokkos::LayoutLeft>::value;
1406 typedef typename std::conditional< is_layout_left,
1407 typename SrcViewDimension:: template prepend< InnerStaticDim+1 >::type,
1408 typename SrcViewDimension:: template append < InnerStaticDim+1 >::type
1409 >::type SrcViewDimensionAppended;
1411 typedef std::integral_constant< unsigned , 0 > padding ;
1413 typedef ViewOffset< SrcViewDimensionAppended, DstLayoutType > TmpOffsetType;
1415 auto src_layout = src.m_array_offset.layout();
1417 if ( is_layout_left ) {
1418 auto prepend_layout = Kokkos::Impl::prependFadToLayout< DstLayoutType >::returnNewLayoutPlusFad(src_layout, InnerStaticDim+1);
1419 TmpOffsetType offset_tmp( padding(), prepend_layout );
1420 dst.m_impl_offset = offset_tmp;
1423 TmpOffsetType offset_tmp( padding(), src_layout );
1424 dst.m_impl_offset = offset_tmp;
1427 Kokkos::abort(
"Sacado error: Applying AssignOffset for case with nested Fads, but without nested Fads - something went wrong");
1432 template <
class DstType,
class SrcFadType >
1433 struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank == (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1436 static void assign( DstType & dst,
const SrcFadType & src )
1438 typedef typename DstType::offset_type dst_offset_type ;
1439 dst.m_impl_offset = dst_offset_type( src.m_array_offset );
1443 template<
class DstType >
1445 void assign( DstType & dst
1446 ,
const SrcFadType & src
1452 std::is_same<
typename DstTraits::array_layout
1453 , Kokkos::LayoutLeft >::value ||
1454 std::is_same<
typename DstTraits::array_layout
1455 , Kokkos::LayoutRight >::value ||
1456 std::is_same<
typename DstTraits::array_layout
1457 , Kokkos::LayoutStride >::value
1461 std::is_same<
typename SrcTraits::array_layout
1462 , Kokkos::LayoutLeft >::value ||
1463 std::is_same<
typename SrcTraits::array_layout
1464 , Kokkos::LayoutRight >::value ||
1465 std::is_same<
typename SrcTraits::array_layout
1466 , Kokkos::LayoutStride >::value
1468 ,
"View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1471 std::is_same<
typename DstTraits::array_layout
1472 ,
typename SrcTraits::array_layout >::value ||
1473 std::is_same<
typename DstTraits::array_layout
1474 , Kokkos::LayoutStride >::value ,
1475 "View assignment must have compatible layout" );
1477 if ( src.m_fad_index != 0 || src.m_fad_stride != 1 ) {
1478 Kokkos::abort(
"\n\n ****** Kokkos::View< Sacado::Fad ... > Cannot assign to array with partitioned view ******\n\n");
1481 AssignOffset< DstType, SrcFadType >::assign( dst, src );
1482 dst.m_impl_handle =
reinterpret_cast< typename DstType::handle_type
>(src.m_impl_handle) ;
1496 template<
class LayoutDest,
class LayoutSrc,
int RankDest,
int RankSrc,
int CurrentArg,
class ... SubViewArgs>
1497 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1498 enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1501 template<
class LayoutDest,
class LayoutSrc,
int RankDest,
int RankSrc,
int CurrentArg,
class ... SubViewArgs>
1502 struct SubviewLegalArgsCompileTime<LayoutDest,Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1503 enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1506 template<
class LayoutDest,
class LayoutSrc,
int RankDest,
int RankSrc,
int CurrentArg,
class ... SubViewArgs>
1507 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,
Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1508 enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1513 template<
class SrcTraits ,
class Arg0 ,
class ... Args >
1515 < typename std::enable_if<(
1517 std::is_same< typename SrcTraits::specialize
1518 , ViewSpecializeSacadoFadContiguous >::value
1521 std::is_same< typename SrcTraits::array_layout
1522 , Kokkos::LayoutLeft >::value ||
1523 std::is_same< typename SrcTraits::array_layout
1524 , Kokkos::LayoutRight >::value ||
1525 std::is_same< typename SrcTraits::array_layout
1526 , Kokkos::LayoutStride >::value
1528 && !Sacado::Fad::is_fad_partition<Arg0>::value
1536 static_assert( SrcTraits::rank ==
sizeof...(Args)+1 ,
"" );
1540 , R0 = bool(is_integral_extent<0,Arg0,Args...>::value)
1541 , R1 = bool(is_integral_extent<1,Arg0,Args...>::value)
1542 , R2 = bool(is_integral_extent<2,Arg0,Args...>::value)
1543 , R3 = bool(is_integral_extent<3,Arg0,Args...>::value)
1544 , R4 = bool(is_integral_extent<4,Arg0,Args...>::value)
1545 , R5 = bool(is_integral_extent<5,Arg0,Args...>::value)
1546 , R6 = bool(is_integral_extent<6,Arg0,Args...>::value)
1550 enum { rank = unsigned(R0) + unsigned(R1) + unsigned(R2) + unsigned(R3)
1551 + unsigned(R4) + unsigned(R5) + unsigned(R6) };
1554 enum { R0_rev = ( 0 == SrcTraits::rank ? RZ : (
1555 1 == SrcTraits::rank ? R0 : (
1556 2 == SrcTraits::rank ? R1 : (
1557 3 == SrcTraits::rank ? R2 : (
1558 4 == SrcTraits::rank ? R3 : (
1559 5 == SrcTraits::rank ? R4 : (
1560 6 == SrcTraits::rank ? R5 : R6 ))))))) };
1563 typedef typename std::conditional<
1569 ( rank <= 2 && R0 && std::is_same< typename SrcTraits::array_layout , Kokkos::LayoutLeft >::value )
1573 ( rank <= 2 && R0_rev && std::is_same< typename SrcTraits::array_layout , Kokkos::LayoutRight >::value )
1575 >::type array_layout ;
1577 typedef typename SrcTraits::value_type fad_type ;
1579 typedef typename std::conditional< rank == 0 , fad_type ,
1580 typename std::conditional< rank == 1 , fad_type * ,
1581 typename std::conditional< rank == 2 , fad_type ** ,
1582 typename std::conditional< rank == 3 , fad_type *** ,
1583 typename std::conditional< rank == 4 , fad_type **** ,
1584 typename std::conditional< rank == 5 , fad_type ***** ,
1585 typename std::conditional< rank == 6 , fad_type ****** ,
1587 >::type >::type >::type >::type >::type >::type >::type
1592 typedef Kokkos::ViewTraits
1595 ,
typename SrcTraits::device_type
1596 ,
typename SrcTraits::memory_traits > traits_type ;
1598 typedef Kokkos::View
1601 ,
typename SrcTraits::device_type
1602 ,
typename SrcTraits::memory_traits > type ;
1606 static void assign( ViewMapping< traits_type , typename traits_type::specialize > & dst
1607 , ViewMapping< SrcTraits , typename SrcTraits::specialize >
const & src
1608 , Arg0 arg0 , Args ... args )
1610 typedef ViewMapping< traits_type , typename traits_type::specialize > DstType ;
1611 typedef typename DstType::offset_type dst_offset_type ;
1612 typedef typename DstType::array_offset_type dst_array_offset_type ;
1613 typedef typename DstType::handle_type dst_handle_type ;
1616 if (std::is_same< typename SrcTraits::array_layout, LayoutLeft >::value) {
1617 const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1618 array_extents( src.m_array_offset.m_dim , Kokkos::ALL() , arg0 , args... );
1619 offset = src.m_array_offset( array_extents.domain_offset(0)
1620 , array_extents.domain_offset(1)
1621 , array_extents.domain_offset(2)
1622 , array_extents.domain_offset(3)
1623 , array_extents.domain_offset(4)
1624 , array_extents.domain_offset(5)
1625 , array_extents.domain_offset(6)
1626 , array_extents.domain_offset(7) );
1627 dst.m_array_offset = dst_array_offset_type( src.m_array_offset ,
1631 const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1632 array_extents( src.m_array_offset.m_dim , arg0 , args... , Kokkos::ALL() );
1633 offset = src.m_array_offset( array_extents.domain_offset(0)
1634 , array_extents.domain_offset(1)
1635 , array_extents.domain_offset(2)
1636 , array_extents.domain_offset(3)
1637 , array_extents.domain_offset(4)
1638 , array_extents.domain_offset(5)
1639 , array_extents.domain_offset(6)
1640 , array_extents.domain_offset(7) );
1641 dst.m_array_offset = dst_array_offset_type( src.m_array_offset ,
1645 const SubviewExtents< SrcTraits::rank , rank >
1646 extents( src.m_impl_offset.m_dim , arg0 , args... );
1648 dst.m_impl_offset = dst_offset_type( src.m_impl_offset , extents );
1649 dst.m_impl_handle = dst_handle_type( src.m_impl_handle + offset );
1650 dst.m_fad_size = src.m_fad_size;
1651 dst.m_original_fad_size = src.m_original_fad_size;
1652 dst.m_fad_stride = src.m_fad_stride;
1653 dst.m_fad_index = src.m_fad_index;
1668 template<
class DataType,
class ...P,
unsigned Stride >
1671 ViewTraits<DataType,P...> ,
1672 Sacado::Fad::Partition<Stride>
1677 enum { is_assignable =
true };
1679 typedef ViewTraits<DataType,P...> src_traits;
1680 typedef ViewMapping< src_traits , typename src_traits::specialize > src_type ;
1682 typedef typename src_type::offset_type::dimension_type src_dimension;
1683 typedef typename src_traits::value_type fad_type;
1684 typedef typename Sacado::LocalScalarType<fad_type,Stride>::type strided_fad_type;
1686 ViewDataType< strided_fad_type , src_dimension >::type strided_data_type;
1687 typedef ViewTraits<strided_data_type,P...> dst_traits;
1688 typedef View<strided_data_type,P...> type;
1689 typedef ViewMapping< dst_traits , typename dst_traits::specialize > dst_type ;
1692 void assign( dst_type & dst
1693 ,
const src_type & src
1694 ,
const Sacado::Fad::Partition<Stride> & part )
1696 if ( Stride != part.stride && Stride != 0 ) {
1697 Kokkos::abort(
"\n\n ****** Kokkos::View< Sacado::Fad ... > Invalid size in partitioned view assignment ******\n\n");
1699 if ( src.m_fad_stride != 1 ) {
1700 Kokkos::abort(
"\n\n ****** Kokkos::View< Sacado::Fad ... > Can't partition already partitioned view ******\n\n");
1703 dst.m_impl_handle = src.m_impl_handle ;
1704 dst.m_impl_offset = src.m_impl_offset ;
1705 dst.m_array_offset = src.m_array_offset ;
1710 (src.m_fad_size.value + part.stride-part.offset-1) / part.stride ;
1712 dst.m_original_fad_size = src.m_original_fad_size ;
1713 dst.m_fad_stride = part.stride ;
1714 dst.m_fad_index = part.offset ;
1721 #endif // defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
1723 #endif // defined(HAVE_SACADO_KOKKOSCORE)
Base template specification for whether a type is a Fad type.
GeneralFad< StaticStorage< T, Num > > SLFad
Base template specification for static size.
#define KOKKOS_INLINE_FUNCTION
GeneralFad< DynamicStorage< T > > DFad
#define KOKKOS_FORCEINLINE_FUNCTION
GeneralFad< StaticFixedStorage< T, Num > > SFad
Base template specification for testing whether type is statically sized.
Get view type for any Fad type.