Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
KokkosExp_View_Fad_Contiguous.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Sacado Package
4 //
5 // Copyright 2006 NTESS and the Sacado contributors.
6 // SPDX-License-Identifier: LGPL-2.1-or-later
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef KOKKOS_EXPERIMENTAL_VIEW_SACADO_FAD_CONTIGUOUS_HPP
11 #define KOKKOS_EXPERIMENTAL_VIEW_SACADO_FAD_CONTIGUOUS_HPP
12 
13 #include "Sacado_ConfigDefs.h"
14 #if defined(HAVE_SACADO_KOKKOS)
15 
17 
18 // Some default traits definitions that need to be defined even if the view
19 // specialization is disabled
20 namespace Kokkos {
21 
22 template <typename ViewType, typename Enabled = void>
23 struct ThreadLocalScalarType {
24  typedef typename ViewType::non_const_value_type type;
25 };
26 
27 template <typename ViewType>
28 struct ViewScalarStride {
29  static const unsigned stride =
30  Impl::LayoutScalarStride< typename ViewType::array_layout>::stride;
31  static const bool is_unit_stride =
32  Impl::LayoutScalarStride< typename ViewType::array_layout>::is_unit_stride;
33 };
34 
35 } // namespace Kokkos
36 
37 namespace Sacado {
38 
39  namespace Fad {
40 
41  /* Define a partition of a View of Sacado Fad type */
42  template <unsigned Size = 0>
43  struct Partition {
44  static const unsigned PartitionSize = Size;
45  unsigned offset ;
46  unsigned stride ;
47 
48  template< typename iType0 , typename iType1 >
49  KOKKOS_INLINE_FUNCTION
50  Partition( const iType0 & i0 , const iType1 & i1 ) :
51  offset(i0), stride(i1) {
52  }
53  };
54 
55  template <typename T>
56  struct is_fad_partition {
57  static const bool value = false;
58  };
59 
60  template <unsigned Stride>
61  struct is_fad_partition< Partition<Stride> > {
62  static const bool value = true;
63  };
64 
65  }
66 
67  // Type of local scalar type when partitioning a view
68  template <typename T, unsigned Stride = 0>
69  struct LocalScalarType {
70  typedef T type;
71  };
72  template <typename T, unsigned Stride>
73  struct LocalScalarType<const T, Stride> {
74  typedef typename LocalScalarType<T,Stride>::type lst;
75  typedef const lst type;
76  };
77 
78  // For DFad, the size is not part of the type, so the default implementation
79  // is sufficient
80 
81  // Type of local scalar type when partitioning a view
82  //
83  // For SLFad, divde the array size by the given stride
84  namespace Fad {
85  namespace Exp {
86  template <typename T, int N> class StaticStorage;
87  template <typename S> class GeneralFad;
88  }
89  }
90  template <typename T, int N, unsigned Stride>
91  struct LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >,
92  Stride > {
93  static const int Ns = (N+Stride-1) / Stride;
94  typedef Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,Ns> > type;
95  };
96 #ifndef SACADO_NEW_FAD_DESIGN_IS_DEFAULT
97  namespace Fad {
98  template <typename T, int N> class SLFad;
99  }
100  template <typename T, int N, unsigned Stride>
101  struct LocalScalarType< Fad::SLFad<T,N>, Stride > {
102  static const int Ns = (N+Stride-1) / Stride;
103  typedef Fad::SLFad<T,Ns> type;
104  };
105 #endif
106 
107  // Type of local scalar type when partitioning a view
108  //
109  // For SFad, divde the array size by the given stride. If it divides evenly,
110  // use SFad, otherwise use SLFad
111  namespace Fad {
112  namespace Exp {
113  template <typename T, int N> class StaticFixedStorage;
114  template <typename T, int N> class StaticStorage;
115  template <typename S> class GeneralFad;
116  }
117  }
118  template <typename T, int N, unsigned Stride>
119  struct LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >,
120  Stride > {
121  static const int Ns = (N+Stride-1) / Stride;
122  typedef typename std::conditional<
123  Ns == N/Stride ,
124  Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,Ns> > ,
125  Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,Ns> >
126  >::type type;
127  };
128 
129 #ifndef SACADO_NEW_FAD_DESIGN_IS_DEFAULT
130  namespace Fad {
131  template <typename T, int N> class SFad;
132  }
133  template <typename T, int N, unsigned Stride>
134  struct LocalScalarType< Fad::SFad<T,N>, Stride > {
135  static const int Ns = (N+Stride-1) / Stride;
136  typedef typename std::conditional< Ns == N/Stride , Fad::SFad<T,Ns> , Fad::SLFad<T,Ns> >::type type;
137  };
138 #endif
139 
140  template <unsigned Stride, typename T>
141  KOKKOS_INLINE_FUNCTION
142  const T& partition_scalar(const T& x) { return x; }
143 
144 } // namespace Sacado
145 
146 // Make sure the user really wants these View specializations
147 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
148 
149 #include "Sacado_Traits.hpp"
150 #include "Kokkos_Core.hpp"
151 #if KOKKOS_VERSION >= 40499
152 #include "View/Kokkos_ViewMapping.hpp"
153 #else
154 #include "impl/Kokkos_ViewMapping.hpp"
155 #endif
156 
157 //----------------------------------------------------------------------------
158 
159 namespace Sacado {
160 
161 #if defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__)
162  namespace Fad {
163  namespace Exp {
164  template <typename T, typename U> class DynamicStorage;
165  template <typename T, int N> class StaticFixedStorage;
166  template <typename T, int N> class StaticStorage;
167  template <typename S> class GeneralFad;
168  }
169  }
170 #ifndef SACADO_VIEW_CUDA_HIERARCHICAL_DFAD
171  template <unsigned Stride, typename T, typename U>
172  KOKKOS_INLINE_FUNCTION
173  typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::DynamicStorage<T,U> >, Stride >::type
174  partition_scalar(const Fad::Exp::GeneralFad< Fad::Exp::DynamicStorage<T,U> >& x) {
175  typedef typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::DynamicStorage<T,U> >, Stride >::type ret_type;
176  const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
177  const int offset = threadIdx.x;
178  ret_type xp(size, x.val());
179 
180  // Note: we can't use x.dx(offset+i*Stride) if
181  // SACADO_VIEW_CUDA_HIERARCHICAL_DFAD_STRIDED is defined because it already
182  // uses blockDim.x in its index calculation. This approach should work
183  // regardless
184  const T* dx = x.dx();
185  for (int i=0; i<size; ++i)
186  xp.fastAccessDx(i) = dx[offset+i*Stride];
187 
188  return xp;
189  }
190 #endif
191  template <unsigned Stride, typename T, int N>
192  KOKKOS_INLINE_FUNCTION
193  typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >, Stride >::type
194  partition_scalar(const Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >& x) {
195  typedef typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >, Stride >::type ret_type;
196  const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
197  const int offset = threadIdx.x;
198  ret_type xp(size, x.val());
199  for (int i=0; i<size; ++i)
200  xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
201  return xp;
202  }
203  template <unsigned Stride, typename T, int N>
204  KOKKOS_INLINE_FUNCTION
205  typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >, Stride >::type
206  partition_scalar(const Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >& x) {
207  typedef typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticFixedStorage<T,N> >, Stride >::type ret_type;
208  const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
209  const int offset = threadIdx.x;
210  ret_type xp(size, x.val());
211  for (int i=0; i<size; ++i)
212  xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
213  return xp;
214  }
215 
216 #ifndef SACADO_NEW_FAD_DESIGN_IS_DEFAULT
217  namespace Fad {
218  template <typename T> class DFad;
219  template <typename T, int N> class SLFad;
220  template <typename T, int N> class SFad;
221  }
222 #ifndef SACADO_VIEW_CUDA_HIERARCHICAL_DFAD
223  template <unsigned Stride, typename T>
224  KOKKOS_INLINE_FUNCTION
225  typename LocalScalarType< Fad::DFad<T>, Stride >::type
226  partition_scalar(const Fad::DFad<T>& x) {
227  typedef typename LocalScalarType< Fad::DFad<T>, Stride >::type ret_type;
228  const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
229  const int offset = threadIdx.x;
230  ret_type xp(size, x.val());
231 
232  // Note: we can't use x.dx(offset+i*Stride) if
233  // SACADO_VIEW_CUDA_HIERARCHICAL_DFAD_STRIDED is defined because it already
234  // uses blockDim.x in its index calculation. This approach should work
235  // regardless
236  const T* dx = x.dx();
237  for (int i=0; i<size; ++i)
238  xp.fastAccessDx(i) = dx[offset+i*Stride];
239 
240  return xp;
241  }
242 #endif
243  template <unsigned Stride, typename T, int N>
244  KOKKOS_INLINE_FUNCTION
245  typename LocalScalarType< Fad::SLFad<T,N>, Stride >::type
246  partition_scalar(const Fad::SLFad<T,N>& x) {
247  typedef typename LocalScalarType< Fad::SLFad<T,N>, Stride >::type ret_type;
248  const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
249  const int offset = threadIdx.x;
250  ret_type xp(size, x.val());
251  for (int i=0; i<size; ++i)
252  xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
253  return xp;
254  }
255  template <unsigned Stride, typename T, int N>
256  KOKKOS_INLINE_FUNCTION
257  typename LocalScalarType< Fad::SFad<T,N>, Stride >::type
258  partition_scalar(const Fad::SFad<T,N>& x) {
259  typedef typename LocalScalarType< Fad::SFad<T,N>, Stride >::type ret_type;
260  const int size = (x.size()+blockDim.x-threadIdx.x-1) / blockDim.x;
261  const int offset = threadIdx.x;
262  ret_type xp(size, x.val());
263  for (int i=0; i<size; ++i)
264  xp.fastAccessDx(i) = x.fastAccessDx(offset+i*Stride);
265  return xp;
266  }
267 #endif // SACADO_NEW_FAD_DESIGN_IS_DEFAULT
268 
269 #endif // defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__)
270 
271 } // namespace Sacado
272 
273 //----------------------------------------------------------------------------
274 
275 namespace Kokkos {
276 
277 template< unsigned Stride, typename D, typename ... P >
278 KOKKOS_INLINE_FUNCTION
279 typename Kokkos::Impl::ViewMapping< void, typename Kokkos::ViewTraits<D,P...>, Sacado::Fad::Partition<Stride> >::type
280 partition( const Kokkos::View<D,P...> & src ,
281  const unsigned offset ,
282  const unsigned stride )
283 {
284  typedef Kokkos::ViewTraits<D,P...> traits;
285  typedef typename Kokkos::Impl::ViewMapping< void, traits, Sacado::Fad::Partition<Stride> >::type DstViewType;
286  const Sacado::Fad::Partition<Stride> part( offset , stride );
287  return DstViewType(src, part);
288 }
289 
290 template <typename ViewType>
291 struct ThreadLocalScalarType<
292  ViewType,
293  typename std::enable_if< is_view_fad_contiguous<ViewType>::value >::type > {
294  typedef typename ViewType::traits TraitsType;
295  typedef Impl::ViewMapping<TraitsType, typename TraitsType::specialize> MappingType;
296  typedef typename MappingType::thread_local_scalar_type type;
297 };
298 
299 namespace Impl {
300 
301 #if defined (KOKKOS_ENABLE_CUDA) && defined(SACADO_VIEW_CUDA_HIERARCHICAL)
302 template< class OutputView >
303 struct SacadoViewFill<
304  OutputView,
305  typename std::enable_if<
306  ( Kokkos::is_view_fad_contiguous<OutputView>::value &&
307  std::is_same<typename OutputView::execution_space, Kokkos::Cuda>::value &&
308  !Kokkos::ViewScalarStride<OutputView>::is_unit_stride )
309  >::type
310  >
311 {
312  typedef typename OutputView::const_value_type const_value_type ;
313  typedef typename OutputView::execution_space execution_space ;
314  typedef Kokkos::TeamPolicy< execution_space> team_policy;
315  typedef typename team_policy::member_type team_impl_handle;
316  typedef typename Kokkos::ThreadLocalScalarType<OutputView>::type local_scalar_type;
317  static const unsigned stride = Kokkos::ViewScalarStride<OutputView>::stride;
318 
319  const OutputView output ;
320  const_value_type input ;
321 
322  KOKKOS_INLINE_FUNCTION
323  void operator()( const size_t i0 ) const
324  {
325  local_scalar_type input_stride = Sacado::partition_scalar<stride>(input);
326 
327  const size_t n1 = output.extent(1);
328  const size_t n2 = output.extent(2);
329  const size_t n3 = output.extent(3);
330  const size_t n4 = output.extent(4);
331  const size_t n5 = output.extent(5);
332  const size_t n6 = output.extent(6);
333 
334  for ( size_t i1 = 0 ; i1 < n1 ; ++i1 ) {
335  for ( size_t i2 = 0 ; i2 < n2 ; ++i2 ) {
336  for ( size_t i3 = 0 ; i3 < n3 ; ++i3 ) {
337  for ( size_t i4 = 0 ; i4 < n4 ; ++i4 ) {
338  for ( size_t i5 = 0 ; i5 < n5 ; ++i5 ) {
339  for ( size_t i6 = 0 ; i6 < n6 ; ++i6 ) {
340  output.access(i0,i1,i2,i3,i4,i5,i6) = input_stride ;
341  }}}}}}
342  }
343 
344  KOKKOS_INLINE_FUNCTION
345  void operator()( const team_impl_handle& team ) const
346  {
347  const size_t i0 = team.league_rank()*team.team_size() + team.team_rank();
348  if (i0 < output.extent(0))
349  (*this)(i0);
350  }
351 
352  SacadoViewFill( const OutputView & arg_out , const_value_type & arg_in )
353  : output( arg_out ), input( arg_in )
354  {
355  const size_t team_size = 256 / stride;
356  team_policy policy( (output.extent(0)+team_size-1)/team_size ,
357  team_size , stride );
358  Kokkos::parallel_for( policy, *this );
359  }
360 };
361 #endif
362 
363 #if defined (KOKKOS_ENABLE_HIP) && defined(SACADO_VIEW_CUDA_HIERARCHICAL)
364 template< class OutputView >
365 struct SacadoViewFill<
366  OutputView,
367  typename std::enable_if<
368  ( Kokkos::is_view_fad_contiguous<OutputView>::value &&
369  std::is_same<typename OutputView::execution_space, Kokkos::HIP>::value &&
370  !Kokkos::ViewScalarStride<OutputView>::is_unit_stride )
371  >::type
372  >
373 {
374  typedef typename OutputView::const_value_type const_value_type ;
375  typedef typename OutputView::execution_space execution_space ;
376  typedef Kokkos::TeamPolicy< execution_space> team_policy;
377  typedef typename team_policy::member_type team_impl_handle;
378  typedef typename Kokkos::ThreadLocalScalarType<OutputView>::type local_scalar_type;
379  static const unsigned stride = Kokkos::ViewScalarStride<OutputView>::stride;
380 
381  const OutputView output ;
382  const_value_type input ;
383 
384  KOKKOS_INLINE_FUNCTION
385  void operator()( const size_t i0 ) const
386  {
387  local_scalar_type input_stride = Sacado::partition_scalar<stride>(input);
388 
389  const size_t n1 = output.extent(1);
390  const size_t n2 = output.extent(2);
391  const size_t n3 = output.extent(3);
392  const size_t n4 = output.extent(4);
393  const size_t n5 = output.extent(5);
394  const size_t n6 = output.extent(6);
395  const size_t n7 = output.extent(7);
396 
397  for ( size_t i1 = 0 ; i1 < n1 ; ++i1 ) {
398  for ( size_t i2 = 0 ; i2 < n2 ; ++i2 ) {
399  for ( size_t i3 = 0 ; i3 < n3 ; ++i3 ) {
400  for ( size_t i4 = 0 ; i4 < n4 ; ++i4 ) {
401  for ( size_t i5 = 0 ; i5 < n5 ; ++i5 ) {
402  for ( size_t i6 = 0 ; i6 < n6 ; ++i6 ) {
403  output.access(i0,i1,i2,i3,i4,i5,i6) = input_stride ;
404  }}}}}}
405  }
406 
407  KOKKOS_INLINE_FUNCTION
408  void operator()( const team_impl_handle& team ) const
409  {
410  const size_t i0 = team.league_rank()*team.team_size() + team.team_rank();
411  if (i0 < output.extent(0))
412  (*this)(i0);
413  }
414 
415  SacadoViewFill( const OutputView & arg_out , const_value_type & arg_in )
416  : output( arg_out ), input( arg_in )
417  {
418  const size_t team_size = 256 / stride;
419  team_policy policy( (output.extent(0)+team_size-1)/team_size ,
420  team_size , stride );
421  Kokkos::parallel_for( policy, *this );
422  }
423 };
424 #endif
425 
426 } // namespace Impl
427 } // namespace Kokkos
428 
429 //----------------------------------------------------------------------------
430 
431 namespace Kokkos {
432 namespace Impl {
433 
434 template< class ... Args >
435 struct is_ViewSpecializeSacadoFadContiguous { enum { value = false }; };
436 
437 template< class D , class ... P , class ... Args >
438 struct is_ViewSpecializeSacadoFadContiguous< Kokkos::View<D,P...> , Args... > {
439  enum { value =
440  std::is_same< typename Kokkos::ViewTraits<D,P...>::specialize
441  , ViewSpecializeSacadoFadContiguous >::value
442  &&
443  ( ( sizeof...(Args) == 0 ) ||
444  is_ViewSpecializeSacadoFadContiguous< Args... >::value ) };
445 };
446 
447 } // namespace Impl
448 } // namespace Kokkos
449 
450 //----------------------------------------------------------------------------
451 
452 namespace Kokkos {
453 namespace Impl {
454 
455 // Compute a partitioned fad size given a stride. Return 0 if the stride
456 // does not evenly divide the size
457 KOKKOS_INLINE_FUNCTION
458 constexpr unsigned computeFadPartitionSize(unsigned size, unsigned stride)
459 {
460  return
461  ((size+stride-1)/stride) == (size/stride) ? ((size+stride-1)/stride) : 0;
462 }
463 
464 // Create new Layout with Fad dimension set to the last
465 // Note: we use enable_if<> here to handle both LayoutLeft and
466 // LayoutContiguous<LayoutLeft>
467 template <unsigned rank, unsigned static_dim, typename Layout>
468 KOKKOS_INLINE_FUNCTION
471  Layout>::type
472 create_fad_array_layout(const Layout& layout)
473 {
474  size_t dims[8];
475  for (int i=0; i<8; ++i)
476  dims[i] = layout.dimension[i];
477  if (static_dim > 0)
478  dims[rank] = static_dim+1;
479  return Layout( dims[0], dims[1], dims[2], dims[3],
480  dims[4], dims[5], dims[6], dims[7] );
481 }
482 
483 // Create new Layout with Fad dimension set to the last
484 // Note: we use enable_if<> here to handle both LayoutStride and
485 // LayoutContiguous<LayoutStride>
486 template <unsigned rank, unsigned static_dim, typename Layout>
487 KOKKOS_INLINE_FUNCTION
489 create_fad_array_layout(const Layout& layout)
490 {
491  size_t dims[8], strides[8];
492  for (int i=0; i<8; ++i) {
493  dims[i] = layout.dimension[i];
494  strides[i] = layout.stride[i];
495  }
496  if (static_dim > 0) {
497  dims[rank] = static_dim+1;
498  strides[rank] = 1;
499  }
500  return Layout( dims[0], strides[0],
501  dims[1], strides[1],
502  dims[2], strides[2],
503  dims[3], strides[3],
504  dims[4], strides[4],
505  dims[5], strides[5],
506  dims[6], strides[6],
507  dims[7], strides[7] );
508 }
509 
510 // Create new LayoutLeft with Fad dimension shuffled to the first
511 // Note: we use enable_if<> here to handle both LayoutLeft and
512 // LayoutContiguous<LayoutLeft>
513  template <unsigned rank, unsigned static_dim, typename Layout>
514 KOKKOS_INLINE_FUNCTION
516 create_fad_array_layout(const Layout& layout)
517 {
518  size_t dims[8];
519  for (int i=0; i<8; ++i)
520  dims[i] = layout.dimension[i];
521  size_t fad_dim = static_dim == 0 ? dims[rank] : static_dim+1;
522  for (int i=rank; i>=1; --i)
523  dims[i] = dims[i-1];
524  dims[0] = fad_dim;
525  return Layout( dims[0], dims[1], dims[2], dims[3],
526  dims[4], dims[5], dims[6], dims[7] );
527 }
528 
529 template <unsigned Rank, typename Dimension, typename Layout>
530 KOKKOS_INLINE_FUNCTION
532 getFadDimension(const ViewOffset<Dimension,Layout,void>& offset)
533 {
534  return
535  ( Rank == 0 ? offset.dimension_0() :
536  ( Rank == 1 ? offset.dimension_1() :
537  ( Rank == 2 ? offset.dimension_2() :
538  ( Rank == 3 ? offset.dimension_3() :
539  ( Rank == 4 ? offset.dimension_4() :
540  ( Rank == 5 ? offset.dimension_5() :
541  ( Rank == 6 ? offset.dimension_6() :
542  offset.dimension_7() )))))));
543 }
544 
545 template <unsigned Rank, typename Dimension, typename Layout>
546 KOKKOS_INLINE_FUNCTION
548 getFadDimension(const ViewOffset<Dimension,Layout,void>& offset)
549 {
550  return offset.dimension_0();
551 }
552 
553 template< class Traits >
554 class ViewMapping< Traits , /* View internal mapping */
555  typename std::enable_if<
556  ( std::is_same< typename Traits::specialize
557  , ViewSpecializeSacadoFadContiguous >::value
558  &&
559  ( std::is_same< typename Traits::array_layout
560  , Kokkos::LayoutLeft >::value
561  ||
562  std::is_same< typename Traits::array_layout
563  , Kokkos::LayoutRight >::value
564  ||
565  std::is_same< typename Traits::array_layout
566  , Kokkos::LayoutStride >::value
567  )
568  )
569  , typename Traits::specialize
570  >::type >
571 {
572 private:
573 
574  template< class , class ... > friend class ViewMapping ;
575  template< class , class ... > friend class Kokkos::View ;
576 
577  typedef typename Traits::value_type fad_type ;
578  typedef typename Sacado::ValueType< fad_type >::type fad_value_type ;
579  typedef typename
580  std::add_const< fad_value_type >::type const_fad_value_type ;
581 
582 public:
583  enum { is_assignable_data_type = true };
584 
585  enum { FadStaticDimension = Sacado::StaticSize< fad_type >::value };
586  enum { PartitionedFadStride = Traits::array_layout::scalar_stride };
587 
588  // The partitioned static size -- this will be 0 if ParitionedFadStride
589  // does not evenly divide FadStaticDimension
590  enum { PartitionedFadStaticDimension =
591  computeFadPartitionSize(FadStaticDimension,PartitionedFadStride) };
592 
593 #ifdef KOKKOS_ENABLE_CUDA
594  typedef typename Sacado::LocalScalarType< fad_type, unsigned(PartitionedFadStride) >::type strided_scalar_type;
595  typedef typename std::conditional< std::is_same<typename Traits::execution_space, Kokkos::Cuda>::value, strided_scalar_type, fad_type >::type thread_local_scalar_type;
596 #elif defined(KOKKOS_ENABLE_HIP)
597  typedef typename Sacado::LocalScalarType< fad_type, unsigned(PartitionedFadStride) >::type strided_scalar_type;
598  typedef typename std::conditional< std::is_same<typename Traits::execution_space, Kokkos::HIP>::value, strided_scalar_type, fad_type >::type thread_local_scalar_type;
599 #else
600  typedef fad_type thread_local_scalar_type;
601 #endif
602 
603 private:
605 
606  typedef fad_value_type * handle_type ;
607 
608  typedef ViewArrayAnalysis< typename Traits::data_type > array_analysis ;
609 
610  typedef ViewOffset< typename Traits::dimension
611  , typename Traits::array_layout
612  , void
613  > offset_type ;
614 
615  // Prepend/append the fad dimension for the internal offset mapping.
616  static constexpr bool is_layout_left =
618  typedef ViewOffset<
619  typename std::conditional<
620  is_layout_left,
621  typename array_analysis::dimension::
622  template prepend<0>::type,
623  typename array_analysis::dimension::
624  template append<0>::type >::type,
625  typename Traits::array_layout,
626  void >
627  array_offset_type ;
628 
629  handle_type m_impl_handle ;
630  offset_type m_impl_offset ;
631  array_offset_type m_array_offset ;
632  sacado_size_type m_fad_size ;
633 
634  // These are for manual partitioning, and will likely be removed
635  unsigned m_original_fad_size ;
636  unsigned m_fad_stride ;
637  unsigned m_fad_index ;
638 
639 public:
640 
641  //----------------------------------------
642  // Domain dimensions
643 
644  enum { Rank = Traits::dimension::rank };
645 
646  // Rank corresponding to the sacado dimension
648 
649  // Using the internal offset mapping so limit to public rank:
650  template< typename iType >
651  KOKKOS_INLINE_FUNCTION constexpr size_t extent( const iType & r ) const
652  { return m_impl_offset.m_dim.extent(r); }
653 
654  KOKKOS_INLINE_FUNCTION constexpr
655  typename Traits::array_layout layout() const
656  { return m_impl_offset.layout(); }
657 
658  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_0() const
659  { return m_impl_offset.dimension_0(); }
660  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_1() const
661  { return m_impl_offset.dimension_1(); }
662  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_2() const
663  { return m_impl_offset.dimension_2(); }
664  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_3() const
665  { return m_impl_offset.dimension_3(); }
666  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_4() const
667  { return m_impl_offset.dimension_4(); }
668  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_5() const
669  { return m_impl_offset.dimension_5(); }
670  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_6() const
671  { return m_impl_offset.dimension_6(); }
672  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_7() const
673  { return m_impl_offset.dimension_7(); }
674 
675  // Is a regular layout with uniform striding for each index.
676  // Since we allow for striding within the data type, we can't guarantee
677  // regular striding
678  using is_regular = std::false_type ;
679 
680  // FIXME: Adjust these for m_stride
681  KOKKOS_INLINE_FUNCTION constexpr size_t stride_0() const
682  { return m_impl_offset.stride_0(); }
683  KOKKOS_INLINE_FUNCTION constexpr size_t stride_1() const
684  { return m_impl_offset.stride_1(); }
685  KOKKOS_INLINE_FUNCTION constexpr size_t stride_2() const
686  { return m_impl_offset.stride_2(); }
687  KOKKOS_INLINE_FUNCTION constexpr size_t stride_3() const
688  { return m_impl_offset.stride_3(); }
689  KOKKOS_INLINE_FUNCTION constexpr size_t stride_4() const
690  { return m_impl_offset.stride_4(); }
691  KOKKOS_INLINE_FUNCTION constexpr size_t stride_5() const
692  { return m_impl_offset.stride_5(); }
693  KOKKOS_INLINE_FUNCTION constexpr size_t stride_6() const
694  { return m_impl_offset.stride_6(); }
695  KOKKOS_INLINE_FUNCTION constexpr size_t stride_7() const
696  { return m_impl_offset.stride_7(); }
697 
698  template< typename iType >
699  KOKKOS_INLINE_FUNCTION constexpr std::enable_if_t<std::is_integral_v<iType>,
700  size_t>
701  stride( iType const s ) const
702  { return m_impl_offset.stride(s) ; }
703 
704  template< typename iType >
705  KOKKOS_INLINE_FUNCTION void stride( iType * const s ) const
706  { m_impl_offset.stride(s); }
707 
708  // Size of sacado scalar dimension
709  KOKKOS_FORCEINLINE_FUNCTION constexpr unsigned dimension_scalar() const
710 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
711  { return PartitionedFadStaticDimension ? PartitionedFadStaticDimension+1 : (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x + 1; }
712 #else
713  { return m_fad_size.value+1; }
714 #endif
715 
716  // trode of sacado scalar dimension
717  KOKKOS_FORCEINLINE_FUNCTION constexpr unsigned stride_scalar() const
718  { return m_fad_stride; }
719 
720  //----------------------------------------
721  // Range of mapping
722 
723 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
724  // Return type of reference operators
725  // this only works if you are using a team-parallel operation on Cuda or HIP!
726  // typedef typename
727  // Sacado::ViewFadType< thread_local_scalar_type , PartitionedFadStaticDimension , (unsigned(ParitionedFadStride) > 1 ? PartitionedFadStride : 0) >::type reference_type ;
728  typedef typename
730 #else
731  // Return type of reference operators
732  typedef typename
734 #endif
735 
737  typedef fad_value_type * pointer_type ;
738 
740  KOKKOS_INLINE_FUNCTION constexpr size_t span() const
741  { return m_array_offset.span(); }
742 
744  KOKKOS_INLINE_FUNCTION constexpr bool span_is_contiguous() const
745  { return m_array_offset.span_is_contiguous() && (m_fad_stride == 1); }
746 
748  KOKKOS_INLINE_FUNCTION constexpr pointer_type data() const
749 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && (defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
750  { return m_impl_handle + threadIdx.x; }
751 #else
752  { return m_impl_handle + m_fad_index; }
753 #endif
754 
755  //----------------------------------------
756 
757  KOKKOS_FORCEINLINE_FUNCTION
758  reference_type
759  reference() const
760  {
761 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
762  const unsigned index = threadIdx.x;
763  const unsigned strd = blockDim.x;
764  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
765 #else
766  const unsigned index = m_fad_index;
767  const unsigned strd = m_fad_stride;
768  const unsigned size = m_fad_size.value;
769 #endif
770  return reference_type( m_impl_handle + index
771  , m_impl_handle + m_original_fad_size
772  , size
773  , strd ); }
774 
775  template< typename I0 >
776  KOKKOS_FORCEINLINE_FUNCTION
778  is_layout_left, reference_type>::type
779  reference( const I0 & i0 ) const
780  { pointer_type beg = m_impl_handle + m_array_offset(0,i0);
781 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
782  const unsigned index = threadIdx.x;
783  const unsigned strd = blockDim.x;
784  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
785 #else
786  const unsigned index = m_fad_index;
787  const unsigned strd = m_fad_stride;
788  const unsigned size = m_fad_size.value;
789 #endif
790  return reference_type( beg + index
791  , beg + m_original_fad_size
792  , size
793  , strd ); }
794 
795  template< typename I0 >
796  KOKKOS_FORCEINLINE_FUNCTION
798  !is_layout_left, reference_type>::type
799  reference( const I0 & i0 ) const
800  { pointer_type beg = m_impl_handle + m_array_offset(i0,0);
801 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
802  const unsigned index = threadIdx.x;
803  const unsigned strd = blockDim.x;
804  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
805 #else
806  const unsigned index = m_fad_index;
807  const unsigned strd = m_fad_stride;
808  const unsigned size = m_fad_size.value;
809 #endif
810  return reference_type( beg + index
811  , beg + m_original_fad_size
812  , size
813  , strd ); }
814 
815  template< typename I0 , typename I1 >
816  KOKKOS_FORCEINLINE_FUNCTION
818  is_layout_left, reference_type>::type
819  reference( const I0 & i0 , const I1 & i1 ) const
820  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1);
821 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && (defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
822  const unsigned index = threadIdx.x;
823  const unsigned strd = blockDim.x;
824  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
825 #else
826  const unsigned index = m_fad_index;
827  const unsigned strd = m_fad_stride;
828  const unsigned size = m_fad_size.value;
829 #endif
830  return reference_type( beg + index
831  , beg + m_original_fad_size
832  , size
833  , strd ); }
834 
835  template< typename I0 , typename I1 >
836  KOKKOS_FORCEINLINE_FUNCTION
838  !is_layout_left, reference_type>::type
839  reference( const I0 & i0 , const I1 & i1 ) const
840  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,0);
841 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
842  const unsigned index = threadIdx.x;
843  const unsigned strd = blockDim.x;
844  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
845 #else
846  const unsigned index = m_fad_index;
847  const unsigned strd = m_fad_stride;
848  const unsigned size = m_fad_size.value;
849 #endif
850  return reference_type( beg + index
851  , beg + m_original_fad_size
852  , size
853  , strd ); }
854 
855 
856  template< typename I0 , typename I1 , typename I2 >
857  KOKKOS_FORCEINLINE_FUNCTION
859  is_layout_left, reference_type>::type
860  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 ) const
861  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2);
862 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
863  const unsigned index = threadIdx.x;
864  const unsigned strd = blockDim.x;
865  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
866 #else
867  const unsigned index = m_fad_index;
868  const unsigned strd = m_fad_stride;
869  const unsigned size = m_fad_size.value;
870 #endif
871  return reference_type( beg + index
872  , beg + m_original_fad_size
873  , size
874  , strd ); }
875 
876  template< typename I0 , typename I1 , typename I2 >
877  KOKKOS_FORCEINLINE_FUNCTION
879  !is_layout_left, reference_type>::type
880  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 ) const
881  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,0);
882 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
883  const unsigned index = threadIdx.x;
884  const unsigned strd = blockDim.x;
885  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
886 #else
887  const unsigned index = m_fad_index;
888  const unsigned strd = m_fad_stride;
889  const unsigned size = m_fad_size.value;
890 #endif
891  return reference_type( beg + index
892  , beg + m_original_fad_size
893  , size
894  , strd ); }
895 
896  template< typename I0 , typename I1 , typename I2 , typename I3 >
897  KOKKOS_FORCEINLINE_FUNCTION
899  is_layout_left, reference_type>::type
900  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3 ) const
901  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3);
902 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
903  const unsigned index = threadIdx.x;
904  const unsigned strd = blockDim.x;
905  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
906 #else
907  const unsigned index = m_fad_index;
908  const unsigned strd = m_fad_stride;
909  const unsigned size = m_fad_size.value;
910 #endif
911  return reference_type( beg + index
912  , beg + m_original_fad_size
913  , size
914  , strd ); }
915 
916  template< typename I0 , typename I1 , typename I2 , typename I3 >
917  KOKKOS_FORCEINLINE_FUNCTION
919  !is_layout_left, reference_type>::type
920  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3 ) const
921  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,0);
922 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
923  const unsigned index = threadIdx.x;
924  const unsigned strd = blockDim.x;
925  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
926 #else
927  const unsigned index = m_fad_index;
928  const unsigned strd = m_fad_stride;
929  const unsigned size = m_fad_size.value;
930 #endif
931  return reference_type( beg + index
932  , beg + m_original_fad_size
933  , size
934  , strd ); }
935 
936  template< typename I0 , typename I1 , typename I2 , typename I3
937  , typename I4 >
938  KOKKOS_FORCEINLINE_FUNCTION
940  is_layout_left, reference_type>::type
941  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
942  , const I4 & i4 ) const
943  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4);
944 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
945  const unsigned index = threadIdx.x;
946  const unsigned strd = blockDim.x;
947  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
948 #else
949  const unsigned index = m_fad_index;
950  const unsigned strd = m_fad_stride;
951  const unsigned size = m_fad_size.value;
952 #endif
953  return reference_type( beg + index
954  , beg + m_original_fad_size
955  , size
956  , strd ); }
957 
958  template< typename I0 , typename I1 , typename I2 , typename I3
959  , typename I4 >
960  KOKKOS_FORCEINLINE_FUNCTION
962  !is_layout_left, reference_type>::type
963  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
964  , const I4 & i4 ) const
965  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,0);
966 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
967  const unsigned index = threadIdx.x;
968  const unsigned strd = blockDim.x;
969  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
970 #else
971  const unsigned index = m_fad_index;
972  const unsigned strd = m_fad_stride;
973  const unsigned size = m_fad_size.value;
974 #endif
975  return reference_type( beg + index
976  , beg + m_original_fad_size
977  , size
978  , strd ); }
979 
980  template< typename I0 , typename I1 , typename I2 , typename I3
981  , typename I4 , typename I5 >
982  KOKKOS_FORCEINLINE_FUNCTION
984  is_layout_left, reference_type>::type
985  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
986  , const I4 & i4 , const I5 & i5 ) const
987  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5);
988 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
989  const unsigned index = threadIdx.x;
990  const unsigned strd = blockDim.x;
991  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
992 #else
993  const unsigned index = m_fad_index;
994  const unsigned strd = m_fad_stride;
995  const unsigned size = m_fad_size.value;
996 #endif
997  return reference_type( beg + index
998  , beg + m_original_fad_size
999  , size
1000  , strd ); }
1001 
1002  template< typename I0 , typename I1 , typename I2 , typename I3
1003  , typename I4 , typename I5 >
1004  KOKKOS_FORCEINLINE_FUNCTION
1006  !is_layout_left, reference_type>::type
1007  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
1008  , const I4 & i4 , const I5 & i5 ) const
1009  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,0);
1010 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
1011  const unsigned index = threadIdx.x;
1012  const unsigned strd = blockDim.x;
1013  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1014 #else
1015  const unsigned index = m_fad_index;
1016  const unsigned strd = m_fad_stride;
1017  const unsigned size = m_fad_size.value;
1018 #endif
1019  return reference_type( beg + index
1020  , beg + m_original_fad_size
1021  , size
1022  , strd ); }
1023 
1024  template< typename I0 , typename I1 , typename I2 , typename I3
1025  , typename I4 , typename I5 , typename I6 >
1026  KOKKOS_FORCEINLINE_FUNCTION
1028  is_layout_left, reference_type>::type
1029  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
1030  , const I4 & i4 , const I5 & i5 , const I6 & i6 ) const
1031  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5,i6);
1032 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
1033  const unsigned index = threadIdx.x;
1034  const unsigned strd = blockDim.x;
1035  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1036 #else
1037  const unsigned index = m_fad_index;
1038  const unsigned strd = m_fad_stride;
1039  const unsigned size = m_fad_size.value;
1040 #endif
1041  return reference_type( beg + index
1042  , beg + m_original_fad_size
1043  , size
1044  , strd ); }
1045 
1046  template< typename I0 , typename I1 , typename I2 , typename I3
1047  , typename I4 , typename I5 , typename I6 >
1048  KOKKOS_FORCEINLINE_FUNCTION
1050  !is_layout_left, reference_type>::type
1051  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
1052  , const I4 & i4 , const I5 & i5 , const I6 & i6 ) const
1053  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,i6,0);
1054 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
1055  const unsigned index = threadIdx.x;
1056  const unsigned strd = blockDim.x;
1057  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1058 #else
1059  const unsigned index = m_fad_index;
1060  const unsigned strd = m_fad_stride;
1061  const unsigned size = m_fad_size.value;
1062 #endif
1063  return reference_type( beg + index
1064  , beg + m_original_fad_size
1065  , size
1066  , strd ); }
1067 
1068  //----------------------------------------
1069 
1071  KOKKOS_INLINE_FUNCTION
1072  static constexpr size_t memory_span( typename Traits::array_layout const & layout )
1073  {
1074  // Do not introduce padding...
1075  typedef std::integral_constant< unsigned , 0 > padding ;
1076  return array_offset_type(
1077  padding() ,
1078  create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( layout ) ).span() * sizeof(fad_value_type);
1079  }
1080 
1081  //----------------------------------------
1082 
1083  KOKKOS_DEFAULTED_FUNCTION ~ViewMapping() = default ;
1084  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) {}
1085 
1086  KOKKOS_DEFAULTED_FUNCTION ViewMapping( const ViewMapping & ) = default ;
1087  KOKKOS_DEFAULTED_FUNCTION ViewMapping & operator = ( const ViewMapping & ) = default ;
1088 
1089  KOKKOS_DEFAULTED_FUNCTION ViewMapping( ViewMapping && ) = default ;
1090  KOKKOS_DEFAULTED_FUNCTION ViewMapping & operator = ( ViewMapping && ) = default ;
1091 
1092  template< class ... P >
1093  KOKKOS_INLINE_FUNCTION
1094  ViewMapping
1095  ( ViewCtorProp< P ... > const & prop
1096  , typename Traits::array_layout const & local_layout
1097  )
1098  : m_impl_handle( ( (ViewCtorProp<void,pointer_type> const &) prop ).value )
1099  , m_impl_offset( std::integral_constant< unsigned , 0 >()
1100  , local_layout )
1101  , m_array_offset(
1102  std::integral_constant< unsigned , 0 >()
1103  , create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( local_layout ) )
1104  , m_fad_size( getFadDimension<unsigned(Rank)>( m_array_offset ) - 1 )
1105  , m_original_fad_size( m_fad_size.value )
1106  , m_fad_stride( 1 )
1107  , m_fad_index( 0 )
1108  {
1109  const unsigned fad_dim =
1110  getFadDimension<unsigned(Rank)>( m_array_offset );
1111  if (unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1112  Kokkos::abort("invalid fad dimension (0) supplied!");
1113  }
1114 
1115  //----------------------------------------
1116  /* Allocate and construct mapped array.
1117  * Allocate via shared allocation record and
1118  * return that record for allocation tracking.
1119  */
1120  template< class ... P >
1121  SharedAllocationRecord<> *
1122  allocate_shared( ViewCtorProp< P... > const & prop
1123  , typename Traits::array_layout const & local_layout
1124  , bool execution_space_specified)
1125  {
1126  typedef ViewCtorProp< P... > ctor_prop ;
1127 
1128  typedef typename ctor_prop::execution_space execution_space ;
1129  typedef typename Traits::memory_space memory_space ;
1130  typedef ViewValueFunctor< execution_space , fad_value_type > functor_type ;
1131  typedef SharedAllocationRecord< memory_space , functor_type > record_type ;
1132 
1133  // Disallow padding
1134  typedef std::integral_constant< unsigned , 0 > padding ;
1135 
1136  // Check if ViewCtorProp has CommonViewAllocProp - if so, retrieve the fad_size and append to layout
1137  enum { test_traits_check = Kokkos::Impl::check_has_common_view_alloc_prop< P... >::value };
1138 
1139  typename Traits::array_layout internal_layout =
1140  (test_traits_check == true)
1141  ? Kokkos::Impl::appendFadToLayoutViewAllocHelper< Traits, P... >::returnNewLayoutPlusFad(prop, local_layout)
1142  : local_layout;
1143 
1144  m_impl_offset = offset_type( padding(), internal_layout );
1145 
1146  m_array_offset =
1147  array_offset_type( padding() ,
1148  create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( internal_layout ) );
1149  const unsigned fad_dim =
1150  getFadDimension<unsigned(Rank)>( m_array_offset );
1151  if (unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1152  Kokkos::abort("invalid fad dimension (0) supplied!");
1153  m_fad_size = fad_dim - 1 ;
1154  m_original_fad_size = m_fad_size.value ;
1155  m_fad_stride = 1;
1156  m_fad_index = 0;
1157 
1158  const size_t alloc_size = m_array_offset.span() * sizeof(fad_value_type);
1159 
1160  // Create shared memory tracking record with allocate memory from the memory space
1161  record_type * const record =
1162  record_type::allocate( ( (ViewCtorProp<void,memory_space> const &) prop ).value
1163  , ( (ViewCtorProp<void,std::string> const &) prop ).value
1164  , alloc_size );
1165 
1166  // Only set the the pointer and initialize if the allocation is non-zero.
1167  // May be zero if one of the dimensions is zero.
1168  if ( alloc_size ) {
1169 
1170  m_impl_handle = handle_type( reinterpret_cast< pointer_type >( record->data() ) );
1171 
1172  if ( ctor_prop::initialize ) {
1173  // Assume destruction is only required when construction is requested.
1174  // The ViewValueFunctor has both value construction and destruction operators.
1175  if (execution_space_specified)
1176  record->m_destroy = functor_type( ( (ViewCtorProp<void,execution_space> const &) prop).value
1177  , (fad_value_type *) m_impl_handle
1178  , m_array_offset.span()
1179  , record->get_label()
1180  );
1181  else
1182  record->m_destroy = functor_type((fad_value_type *) m_impl_handle
1183  , m_array_offset.span()
1184  , record->get_label()
1185  );
1186 
1187  // Construct values
1188  record->m_destroy.construct_shared_allocation();
1189  }
1190  }
1191 
1192  return record ;
1193  }
1194 
1195 };
1196 
1197 } // namespace Impl
1198 } // namespace Kokkos
1199 
1200 //----------------------------------------------------------------------------
1201 
1202 namespace Kokkos {
1203 namespace Impl {
1204 
1209 template< class DstTraits , class SrcTraits >
1210 class ViewMapping< DstTraits , SrcTraits ,
1211  typename std::enable_if<(
1212  Kokkos::Impl::MemorySpaceAccess
1213  < typename DstTraits::memory_space
1214  , typename SrcTraits::memory_space >::assignable
1215  &&
1216  // Destination view has FAD
1217  std::is_same< typename DstTraits::specialize
1218  , ViewSpecializeSacadoFadContiguous >::value
1219  &&
1220  // Source view has FAD
1221  std::is_same< typename SrcTraits::specialize
1222  , ViewSpecializeSacadoFadContiguous >::value
1223  )
1224  , typename DstTraits::specialize
1225  >::type >
1226 {
1227 public:
1228 
1229  enum { is_assignable = true };
1230  enum { is_assignable_data_type = true };
1231 
1232  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1233  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1234  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1235 
1236  template< class DstType >
1237  KOKKOS_INLINE_FUNCTION static
1238  void assign( DstType & dst
1239  , const SrcFadType & src
1240  , const TrackType & )
1241  {
1242  static_assert(
1243  (
1244  std::is_same< typename DstTraits::array_layout
1245  , Kokkos::LayoutLeft >::value ||
1246  std::is_same< typename DstTraits::array_layout
1247  , Kokkos::LayoutRight >::value ||
1248  std::is_same< typename DstTraits::array_layout
1249  , Kokkos::LayoutStride >::value
1250  )
1251  &&
1252  (
1253  std::is_same< typename SrcTraits::array_layout
1254  , Kokkos::LayoutLeft >::value ||
1255  std::is_same< typename SrcTraits::array_layout
1256  , Kokkos::LayoutRight >::value ||
1257  std::is_same< typename SrcTraits::array_layout
1258  , Kokkos::LayoutStride >::value
1259  )
1260  , "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1261 
1262  static_assert(
1263  std::is_same< typename DstTraits::array_layout
1264  , typename SrcTraits::array_layout >::value ||
1265  std::is_same< typename DstTraits::array_layout
1266  , Kokkos::LayoutStride >::value ,
1267  "View assignment must have compatible layout" );
1268 
1269  static_assert(
1270  std::is_same< typename DstTraits::value_type
1271  , typename SrcTraits::value_type >::value ||
1272  std::is_same< typename DstTraits::value_type
1273  , typename SrcTraits::const_value_type >::value ,
1274  "View assignment must have same value type or const = non-const" );
1275 
1276  static_assert(
1277  ViewDimensionAssignable
1278  < typename DstType::offset_type::dimension_type
1279  , typename SrcFadType::offset_type::dimension_type >::value ,
1280  "View assignment must have compatible dimensions" );
1281 
1282  static_assert(
1283  ViewDimensionAssignable
1284  < typename DstType::array_offset_type::dimension_type
1285  , typename SrcFadType::array_offset_type::dimension_type >::value ,
1286  "View assignment must have compatible dimensions" );
1287 
1288  typedef typename DstType::offset_type dst_offset_type ;
1289  typedef typename DstType::array_offset_type dst_array_offset_type ;
1290 
1291  dst.m_impl_handle = src.m_impl_handle ;
1292  dst.m_impl_offset = dst_offset_type( src.m_impl_offset );
1293  dst.m_array_offset = dst_array_offset_type( src.m_array_offset );
1294  dst.m_fad_size = src.m_fad_size.value ;
1295  dst.m_original_fad_size = src.m_original_fad_size ;
1296  dst.m_fad_stride = src.m_fad_stride ;
1297  dst.m_fad_index = src.m_fad_index ;
1298  }
1299 };
1300 
1305 template< class DstTraits , class SrcTraits >
1306 class ViewMapping< DstTraits , SrcTraits ,
1307  typename std::enable_if<(
1308  std::is_same< typename DstTraits::memory_space
1309  , typename SrcTraits::memory_space >::value
1310  &&
1311  // Destination view has FAD
1312  std::is_same< typename DstTraits::specialize
1313  , ViewSpecializeSacadoFad >::value
1314  &&
1315  // Source view has FAD contiguous
1316  std::is_same< typename SrcTraits::specialize
1317  , ViewSpecializeSacadoFadContiguous >::value
1318  &&
1319  // Destination view is LayoutStride
1320  std::is_same< typename DstTraits::array_layout
1321  , Kokkos::LayoutStride >::value
1322  )
1323  , typename DstTraits::specialize
1324  >::type >
1325 {
1326 public:
1327 
1328  enum { is_assignable = true };
1329  enum { is_assignable_data_type = true };
1330 
1331  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1332  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1333  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1334 
1335  template< class DstType >
1336  KOKKOS_INLINE_FUNCTION static
1337  void assign( DstType & dst
1338  , const SrcFadType & src
1339  , const TrackType & )
1340  {
1341  static_assert(
1342  std::is_same< typename SrcTraits::array_layout
1343  , Kokkos::LayoutLeft >::value ||
1344  std::is_same< typename SrcTraits::array_layout
1345  , Kokkos::LayoutRight >::value ||
1346  std::is_same< typename SrcTraits::array_layout
1347  , Kokkos::LayoutStride >::value ,
1348  "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1349 
1350  static_assert(
1351  std::is_same< typename DstTraits::value_type
1352  , typename SrcTraits::value_type >::value ||
1353  std::is_same< typename DstTraits::value_type
1354  , typename SrcTraits::const_value_type >::value ,
1355  "View assignment must have same value type or const = non-const" );
1356 
1357  static_assert(
1358  DstTraits::dimension::rank == SrcTraits::dimension::rank,
1359  "View assignment must have same rank" );
1360 
1361  typedef typename DstType::array_offset_type dst_offset_type ;
1362 
1363  dst.m_impl_handle = src.m_impl_handle ;
1364  dst.m_fad_size = src.m_fad_size.value ;
1365  dst.m_fad_stride = src.m_fad_stride ;
1366  dst.m_impl_offset = src.m_impl_offset;
1367 
1368  size_t N[8], S[8];
1369  N[0] = src.m_array_offset.dimension_0();
1370  N[1] = src.m_array_offset.dimension_1();
1371  N[2] = src.m_array_offset.dimension_2();
1372  N[3] = src.m_array_offset.dimension_3();
1373  N[4] = src.m_array_offset.dimension_4();
1374  N[5] = src.m_array_offset.dimension_5();
1375  N[6] = src.m_array_offset.dimension_6();
1376  N[7] = src.m_array_offset.dimension_7();
1377  S[0] = src.m_array_offset.stride_0();
1378  S[1] = src.m_array_offset.stride_1();
1379  S[2] = src.m_array_offset.stride_2();
1380  S[3] = src.m_array_offset.stride_3();
1381  S[4] = src.m_array_offset.stride_4();
1382  S[5] = src.m_array_offset.stride_5();
1383  S[6] = src.m_array_offset.stride_6();
1384  S[7] = src.m_array_offset.stride_7();
1385 
1386  // For LayoutLeft, we have to move the Sacado dimension from the first
1387  // to the last
1388  if (std::is_same< typename SrcTraits::array_layout
1389  , Kokkos::LayoutLeft >::value)
1390  {
1391  const size_t N_fad = N[0];
1392  const size_t S_fad = S[0];
1393  for (int i=0; i<7; ++i) {
1394  N[i] = N[i+1];
1395  S[i] = S[i+1];
1396  }
1397  N[DstTraits::dimension::rank] = N_fad;
1398  S[DstTraits::dimension::rank] = S_fad;
1399  }
1400  Kokkos::LayoutStride ls( N[0], S[0],
1401  N[1], S[1],
1402  N[2], S[2],
1403  N[3], S[3],
1404  N[4], S[4],
1405  N[5], S[5],
1406  N[6], S[6],
1407  N[7], S[7] );
1408  dst.m_array_offset = dst_offset_type(std::integral_constant<unsigned,0>(), ls);
1409  }
1410 };
1411 
1416 template< class DstTraits , class SrcTraits >
1417 class ViewMapping< DstTraits , SrcTraits ,
1418  typename std::enable_if<(
1419  std::is_same< typename DstTraits::memory_space
1420  , typename SrcTraits::memory_space >::value
1421  &&
1422  // Destination view has ordinary
1423  std::is_same< typename DstTraits::specialize , void >::value
1424  &&
1425  // Source view has FAD only
1426  std::is_same< typename SrcTraits::specialize
1427  , ViewSpecializeSacadoFadContiguous >::value
1428  )
1429  , typename DstTraits::specialize
1430  >::type >
1431 {
1432 public:
1433 
1434  enum { is_assignable = true };
1435  enum { is_assignable_data_type = true };
1436 
1437  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1438  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1439  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1440 
1441 
1442  // Helpers to assign, and generate if necessary, ViewOffset to the dst map
1443  // These are necessary to use Kokkos' deep_copy with nested fads
1444  template < class DstType, class SrcFadType, class Enable = void >
1445  struct AssignOffset;
1446 
1447  template < class DstType, class SrcFadType >
1448  struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank != (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1449  {
1450  // ViewOffset's Dimensions Ranks do not match
1451  KOKKOS_INLINE_FUNCTION
1452  static void assign( DstType & dst, const SrcFadType & src )
1453  {
1454  typedef typename SrcTraits::value_type TraitsValueType;
1455 
1458  )
1459  {
1460 
1461  typedef typename DstType::offset_type::array_layout DstLayoutType;
1462  //typedef typename ViewArrayLayoutSelector<typename DstType::offset_type::array_layout>::type DstLayoutType;
1463  typedef typename SrcFadType::array_offset_type::dimension_type SrcViewDimension;
1464 
1465  // This is the static dimension of the inner fad, missing from ViewDimension
1467 
1468  static constexpr bool is_layout_left =
1470 
1471  typedef typename std::conditional< is_layout_left,
1472  typename SrcViewDimension:: template prepend< InnerStaticDim+1 >::type,
1473  typename SrcViewDimension:: template append < InnerStaticDim+1 >::type
1474  >::type SrcViewDimensionAppended;
1475 
1476  typedef std::integral_constant< unsigned , 0 > padding ;
1477 
1478  typedef ViewOffset< SrcViewDimensionAppended, DstLayoutType > TmpOffsetType;
1479 
1480  auto src_layout = src.m_array_offset.layout();
1481 
1482  if ( is_layout_left ) {
1483  auto prepend_layout = Kokkos::Impl::prependFadToLayout< DstLayoutType >::returnNewLayoutPlusFad(src_layout, InnerStaticDim+1);
1484  TmpOffsetType offset_tmp( padding(), prepend_layout );
1485  dst.m_impl_offset = offset_tmp;
1486  }
1487  else {
1488  TmpOffsetType offset_tmp( padding(), src_layout );
1489  dst.m_impl_offset = offset_tmp;
1490  }
1491  } else {
1492  Kokkos::abort("Sacado error: Applying AssignOffset for case with nested Fads, but without nested Fads - something went wrong");
1493  }
1494  }
1495  };
1496 
1497  template < class DstType, class SrcFadType >
1498  struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank == (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1499  {
1500  KOKKOS_INLINE_FUNCTION
1501  static void assign( DstType & dst, const SrcFadType & src )
1502  {
1503  typedef typename DstType::offset_type dst_offset_type ;
1504  dst.m_impl_offset = dst_offset_type( src.m_array_offset );
1505  }
1506  };
1507 
1508  template< class DstType >
1509  KOKKOS_INLINE_FUNCTION static
1510  void assign( DstType & dst
1511  , const SrcFadType & src
1512  , const TrackType &
1513  )
1514  {
1515  static_assert(
1516  (
1517  std::is_same< typename DstTraits::array_layout
1518  , Kokkos::LayoutLeft >::value ||
1519  std::is_same< typename DstTraits::array_layout
1520  , Kokkos::LayoutRight >::value ||
1521  std::is_same< typename DstTraits::array_layout
1522  , Kokkos::LayoutStride >::value
1523  )
1524  &&
1525  (
1526  std::is_same< typename SrcTraits::array_layout
1527  , Kokkos::LayoutLeft >::value ||
1528  std::is_same< typename SrcTraits::array_layout
1529  , Kokkos::LayoutRight >::value ||
1530  std::is_same< typename SrcTraits::array_layout
1531  , Kokkos::LayoutStride >::value
1532  )
1533  , "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1534 
1535  static_assert(
1536  std::is_same< typename DstTraits::array_layout
1537  , typename SrcTraits::array_layout >::value ||
1538  std::is_same< typename DstTraits::array_layout
1539  , Kokkos::LayoutStride >::value ,
1540  "View assignment must have compatible layout" );
1541 
1542  if ( src.m_fad_index != 0 || src.m_fad_stride != 1 ) {
1543  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Cannot assign to array with partitioned view ******\n\n");
1544  }
1545 
1546  AssignOffset< DstType, SrcFadType >::assign( dst, src );
1547  dst.m_impl_handle = reinterpret_cast< typename DstType::handle_type >(src.m_impl_handle) ;
1548  }
1549 };
1550 
1551 } // namespace Impl
1552 } // namespace Kokkos
1553 
1554 //----------------------------------------------------------------------------
1555 
1556 namespace Kokkos {
1557 namespace Impl {
1558 
1559 // Rules for subview arguments and layouts matching
1560 
1561 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1562 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1563  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1564 };
1565 
1566 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1567 struct SubviewLegalArgsCompileTime<LayoutDest,Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1568  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1569 };
1570 
1571 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1572 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1573  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1574 };
1575 
1576 // Subview mapping
1577 
1578 template< class SrcTraits , class Arg0 , class ... Args >
1579 struct ViewMapping
1580  < typename std::enable_if<(
1581  // Source view has FAD only
1582  std::is_same< typename SrcTraits::specialize
1583  , ViewSpecializeSacadoFadContiguous >::value
1584  &&
1585  (
1586  std::is_same< typename SrcTraits::array_layout
1587  , Kokkos::LayoutLeft >::value ||
1588  std::is_same< typename SrcTraits::array_layout
1589  , Kokkos::LayoutRight >::value ||
1590  std::is_same< typename SrcTraits::array_layout
1591  , Kokkos::LayoutStride >::value
1592  )
1593  && !Sacado::Fad::is_fad_partition<Arg0>::value
1594  )
1595  >::type
1596  , SrcTraits
1597  , Arg0, Args ... >
1598 {
1599 private:
1600 
1601  static_assert( SrcTraits::rank == sizeof...(Args)+1 , "" );
1602 
1603  enum
1604  { RZ = false
1612  };
1613 
1614  // Public rank
1615  enum { rank = unsigned(R0) + unsigned(R1) + unsigned(R2) + unsigned(R3)
1616  + unsigned(R4) + unsigned(R5) + unsigned(R6) };
1617 
1618  // Whether right-most non-FAD rank is a range.
1619  enum { R0_rev = ( 0 == SrcTraits::rank ? RZ : (
1620  1 == SrcTraits::rank ? R0 : (
1621  2 == SrcTraits::rank ? R1 : (
1622  3 == SrcTraits::rank ? R2 : (
1623  4 == SrcTraits::rank ? R3 : (
1624  5 == SrcTraits::rank ? R4 : (
1625  6 == SrcTraits::rank ? R5 : R6 ))))))) };
1626 
1627  // Subview's layout
1628  typedef typename std::conditional<
1629  ( /* Same array layout IF */
1630  ( rank == 0 ) /* output rank zero */
1631  ||
1632  // OutputRank 1, InputLayout Left, Interval 0
1633  // because single stride one
1635  ||
1636  // OutputRank 1, InputLayout Right, Interval [InputRank-1]
1637  // because single stride one
1640  >::type array_layout ;
1641 
1642  typedef typename SrcTraits::value_type fad_type ;
1643 
1644  typedef typename std::conditional< rank == 0 , fad_type ,
1645  typename std::conditional< rank == 1 , fad_type * ,
1646  typename std::conditional< rank == 2 , fad_type ** ,
1647  typename std::conditional< rank == 3 , fad_type *** ,
1648  typename std::conditional< rank == 4 , fad_type **** ,
1649  typename std::conditional< rank == 5 , fad_type ***** ,
1650  typename std::conditional< rank == 6 , fad_type ****** ,
1651  fad_type *******
1652  >::type >::type >::type >::type >::type >::type >::type
1653  data_type ;
1654 
1655 public:
1656 
1657  typedef Kokkos::ViewTraits
1658  < data_type
1659  , array_layout
1660  , typename SrcTraits::device_type
1661  , typename SrcTraits::memory_traits > traits_type ;
1662 
1663  typedef Kokkos::View
1664  < data_type
1665  , array_layout
1666  , typename SrcTraits::device_type
1667  , typename SrcTraits::memory_traits > type ;
1668 
1669 
1670  KOKKOS_INLINE_FUNCTION
1671  static void assign( ViewMapping< traits_type , typename traits_type::specialize > & dst
1672  , ViewMapping< SrcTraits , typename SrcTraits::specialize > const & src
1673  , Arg0 arg0 , Args ... args )
1674  {
1675  typedef ViewMapping< traits_type , typename traits_type::specialize > DstType ;
1676  typedef typename DstType::offset_type dst_offset_type ;
1677  typedef typename DstType::array_offset_type dst_array_offset_type ;
1678  typedef typename DstType::handle_type dst_handle_type ;
1679 
1680  size_t offset;
1682  const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1683  array_extents( src.m_array_offset.m_dim , Kokkos::ALL() , arg0 , args... );
1684  offset = src.m_array_offset( array_extents.domain_offset(0)
1685  , array_extents.domain_offset(1)
1686  , array_extents.domain_offset(2)
1687  , array_extents.domain_offset(3)
1688  , array_extents.domain_offset(4)
1689  , array_extents.domain_offset(5)
1690  , array_extents.domain_offset(6)
1691  , array_extents.domain_offset(7) );
1692  dst_array_offset_type dst_array_offset( src.m_array_offset ,
1693  array_extents );
1694  // For LayoutStride, we always use LayoutRight indexing (because we
1695  // don't know whether the original array was Left or Right), so we
1696  // need to swap the Fad dimension to the last and shift all of the
1697  // other dimensions left by 1
1699  {
1700  Kokkos::LayoutStride ls(
1701  dst_array_offset.m_dim.N0, dst_array_offset.m_stride.S0,
1702  dst_array_offset.m_dim.N1, dst_array_offset.m_stride.S1,
1703  dst_array_offset.m_dim.N2, dst_array_offset.m_stride.S2,
1704  dst_array_offset.m_dim.N3, dst_array_offset.m_stride.S3,
1705  dst_array_offset.m_dim.N4, dst_array_offset.m_stride.S4,
1706  dst_array_offset.m_dim.N5, dst_array_offset.m_stride.S5,
1707  dst_array_offset.m_dim.N6, dst_array_offset.m_stride.S6,
1708  dst_array_offset.m_dim.N7, dst_array_offset.m_stride.S7);
1709  auto t1 = ls.dimension[0];
1710  for (unsigned i=0; i<rank; ++i)
1711  ls.dimension[i] = ls.dimension[i+1];
1712  ls.dimension[rank] = t1;
1713  auto t2 = ls.stride[0];
1714  for (unsigned i=0; i<rank; ++i)
1715  ls.stride[i] = ls.stride[i+1];
1716  ls.stride[rank] = t2;
1717  dst.m_array_offset = dst_array_offset_type(std::integral_constant<unsigned, 0>(), ls);
1718  }
1719  else
1720  dst.m_array_offset = dst_array_offset;
1721  }
1722  else {
1723  const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1724  array_extents( src.m_array_offset.m_dim , arg0 , args... , Kokkos::ALL() );
1725  offset = src.m_array_offset( array_extents.domain_offset(0)
1726  , array_extents.domain_offset(1)
1727  , array_extents.domain_offset(2)
1728  , array_extents.domain_offset(3)
1729  , array_extents.domain_offset(4)
1730  , array_extents.domain_offset(5)
1731  , array_extents.domain_offset(6)
1732  , array_extents.domain_offset(7) );
1733  dst.m_array_offset = dst_array_offset_type( src.m_array_offset ,
1734  array_extents );
1735  }
1736 
1737  const SubviewExtents< SrcTraits::rank , rank >
1738  extents( src.m_impl_offset.m_dim , arg0 , args... );
1739 
1740  dst.m_impl_offset = dst_offset_type( src.m_impl_offset , extents );
1741  dst.m_impl_handle = dst_handle_type( src.m_impl_handle + offset );
1742  dst.m_fad_size = src.m_fad_size;
1743  dst.m_original_fad_size = src.m_original_fad_size;
1744  dst.m_fad_stride = src.m_fad_stride;
1745  dst.m_fad_index = src.m_fad_index;
1746  }
1747 
1748 };
1749 
1750 } // namespace Impl
1751 } // namespace Kokkos
1752 
1753 //---------------------------------------------------------------------------
1754 
1755 namespace Kokkos {
1756 namespace Impl {
1757 
1758 // Partition mapping
1759 
1760 template< class DataType, class ...P, unsigned Stride >
1761 class ViewMapping<
1762  void,
1763  ViewTraits<DataType,P...> ,
1764  Sacado::Fad::Partition<Stride>
1765  >
1766 {
1767 public:
1768 
1769  enum { is_assignable = true };
1770  enum { is_assignable_data_type = true };
1771 
1772  typedef ViewTraits<DataType,P...> src_traits;
1773  typedef ViewMapping< src_traits , typename src_traits::specialize > src_type ;
1774 
1775  typedef typename src_type::offset_type::dimension_type src_dimension;
1776  typedef typename src_traits::value_type fad_type;
1777  typedef typename Sacado::LocalScalarType<fad_type,Stride>::type strided_fad_type;
1778  typedef typename
1779  ViewDataType< strided_fad_type , src_dimension >::type strided_data_type;
1780  typedef ViewTraits<strided_data_type,P...> dst_traits;
1781  typedef View<strided_data_type,P...> type;
1782  typedef ViewMapping< dst_traits , typename dst_traits::specialize > dst_type ;
1783 
1784  KOKKOS_INLINE_FUNCTION static
1785  void assign( dst_type & dst
1786  , const src_type & src
1787  , const Sacado::Fad::Partition<Stride> & part )
1788  {
1789  if ( Stride != part.stride && Stride != 0 ) {
1790  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Invalid size in partitioned view assignment ******\n\n");
1791  }
1792  if ( src.m_fad_stride != 1 ) {
1793  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Can't partition already partitioned view ******\n\n");
1794  }
1795 
1796  dst.m_impl_handle = src.m_impl_handle ;
1797  dst.m_impl_offset = src.m_impl_offset ;
1798  dst.m_array_offset = src.m_array_offset ;
1799 
1800  // Assuming the alignment was choosen correctly for the partitioning,
1801  // each partition should get the same size. This allows the use of SFad.
1802  dst.m_fad_size =
1803  (src.m_fad_size.value + part.stride-part.offset-1) / part.stride ;
1804 
1805  dst.m_original_fad_size = src.m_original_fad_size ;
1806  dst.m_fad_stride = part.stride ;
1807  dst.m_fad_index = part.offset ;
1808  }
1809 };
1810 
1811 } // namespace Impl
1812 } // namespace Kokkos
1813 
1814 #endif // defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
1815 
1816 #endif // defined(HAVE_SACADO_KOKKOS)
1817 
1818 #endif /* #ifndef KOKKOS_EXPERIMENTAL_VIEW_SACADO_FAD_HPP */
expr expr dx(i)
Base template specification for whether a type is a Fad type.
GeneralFad< StaticStorage< T, Num > > SLFad
Base template specification for static size.
#define T
Definition: Sacado_rad.hpp:553
#define D
Definition: Sacado_rad.hpp:557
GeneralFad< DynamicStorage< T > > DFad
int value
const int N
void
Definition: uninit.c:105
const int fad_dim
expr expr expr bar false
GeneralFad< StaticFixedStorage< T, Num > > SFad
if(first)
Definition: uninit.c:110
Base template specification for testing whether type is statically sized.
Get view type for any Fad type.