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 void stride( iType * const s ) const
700  { m_impl_offset.stride(s); }
701 
702  // Size of sacado scalar dimension
703  KOKKOS_FORCEINLINE_FUNCTION constexpr unsigned dimension_scalar() const
704 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
705  { return PartitionedFadStaticDimension ? PartitionedFadStaticDimension+1 : (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x + 1; }
706 #else
707  { return m_fad_size.value+1; }
708 #endif
709 
710  // trode of sacado scalar dimension
711  KOKKOS_FORCEINLINE_FUNCTION constexpr unsigned stride_scalar() const
712  { return m_fad_stride; }
713 
714  //----------------------------------------
715  // Range of mapping
716 
717 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
718  // Return type of reference operators
719  // this only works if you are using a team-parallel operation on Cuda or HIP!
720  // typedef typename
721  // Sacado::ViewFadType< thread_local_scalar_type , PartitionedFadStaticDimension , (unsigned(ParitionedFadStride) > 1 ? PartitionedFadStride : 0) >::type reference_type ;
722  typedef typename
724 #else
725  // Return type of reference operators
726  typedef typename
728 #endif
729 
731  typedef fad_value_type * pointer_type ;
732 
734  KOKKOS_INLINE_FUNCTION constexpr size_t span() const
735  { return m_array_offset.span(); }
736 
738  KOKKOS_INLINE_FUNCTION constexpr bool span_is_contiguous() const
739  { return m_array_offset.span_is_contiguous() && (m_fad_stride == 1); }
740 
742  KOKKOS_INLINE_FUNCTION constexpr pointer_type data() const
743 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && (defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
744  { return m_impl_handle + threadIdx.x; }
745 #else
746  { return m_impl_handle + m_fad_index; }
747 #endif
748 
749  //----------------------------------------
750 
751  KOKKOS_FORCEINLINE_FUNCTION
752  reference_type
753  reference() const
754  {
755 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
756  const unsigned index = threadIdx.x;
757  const unsigned strd = blockDim.x;
758  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
759 #else
760  const unsigned index = m_fad_index;
761  const unsigned strd = m_fad_stride;
762  const unsigned size = m_fad_size.value;
763 #endif
764  return reference_type( m_impl_handle + index
765  , m_impl_handle + m_original_fad_size
766  , size
767  , strd ); }
768 
769  template< typename I0 >
770  KOKKOS_FORCEINLINE_FUNCTION
772  is_layout_left, reference_type>::type
773  reference( const I0 & i0 ) const
774  { pointer_type beg = m_impl_handle + m_array_offset(0,i0);
775 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
776  const unsigned index = threadIdx.x;
777  const unsigned strd = blockDim.x;
778  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
779 #else
780  const unsigned index = m_fad_index;
781  const unsigned strd = m_fad_stride;
782  const unsigned size = m_fad_size.value;
783 #endif
784  return reference_type( beg + index
785  , beg + m_original_fad_size
786  , size
787  , strd ); }
788 
789  template< typename I0 >
790  KOKKOS_FORCEINLINE_FUNCTION
792  !is_layout_left, reference_type>::type
793  reference( const I0 & i0 ) const
794  { pointer_type beg = m_impl_handle + m_array_offset(i0,0);
795 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
796  const unsigned index = threadIdx.x;
797  const unsigned strd = blockDim.x;
798  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
799 #else
800  const unsigned index = m_fad_index;
801  const unsigned strd = m_fad_stride;
802  const unsigned size = m_fad_size.value;
803 #endif
804  return reference_type( beg + index
805  , beg + m_original_fad_size
806  , size
807  , strd ); }
808 
809  template< typename I0 , typename I1 >
810  KOKKOS_FORCEINLINE_FUNCTION
812  is_layout_left, reference_type>::type
813  reference( const I0 & i0 , const I1 & i1 ) const
814  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1);
815 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && (defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
816  const unsigned index = threadIdx.x;
817  const unsigned strd = blockDim.x;
818  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
819 #else
820  const unsigned index = m_fad_index;
821  const unsigned strd = m_fad_stride;
822  const unsigned size = m_fad_size.value;
823 #endif
824  return reference_type( beg + index
825  , beg + m_original_fad_size
826  , size
827  , strd ); }
828 
829  template< typename I0 , typename I1 >
830  KOKKOS_FORCEINLINE_FUNCTION
832  !is_layout_left, reference_type>::type
833  reference( const I0 & i0 , const I1 & i1 ) const
834  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,0);
835 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
836  const unsigned index = threadIdx.x;
837  const unsigned strd = blockDim.x;
838  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
839 #else
840  const unsigned index = m_fad_index;
841  const unsigned strd = m_fad_stride;
842  const unsigned size = m_fad_size.value;
843 #endif
844  return reference_type( beg + index
845  , beg + m_original_fad_size
846  , size
847  , strd ); }
848 
849 
850  template< typename I0 , typename I1 , typename I2 >
851  KOKKOS_FORCEINLINE_FUNCTION
853  is_layout_left, reference_type>::type
854  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 ) const
855  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2);
856 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
857  const unsigned index = threadIdx.x;
858  const unsigned strd = blockDim.x;
859  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
860 #else
861  const unsigned index = m_fad_index;
862  const unsigned strd = m_fad_stride;
863  const unsigned size = m_fad_size.value;
864 #endif
865  return reference_type( beg + index
866  , beg + m_original_fad_size
867  , size
868  , strd ); }
869 
870  template< typename I0 , typename I1 , typename I2 >
871  KOKKOS_FORCEINLINE_FUNCTION
873  !is_layout_left, reference_type>::type
874  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 ) const
875  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,0);
876 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
877  const unsigned index = threadIdx.x;
878  const unsigned strd = blockDim.x;
879  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
880 #else
881  const unsigned index = m_fad_index;
882  const unsigned strd = m_fad_stride;
883  const unsigned size = m_fad_size.value;
884 #endif
885  return reference_type( beg + index
886  , beg + m_original_fad_size
887  , size
888  , strd ); }
889 
890  template< typename I0 , typename I1 , typename I2 , typename I3 >
891  KOKKOS_FORCEINLINE_FUNCTION
893  is_layout_left, reference_type>::type
894  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3 ) const
895  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3);
896 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
897  const unsigned index = threadIdx.x;
898  const unsigned strd = blockDim.x;
899  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
900 #else
901  const unsigned index = m_fad_index;
902  const unsigned strd = m_fad_stride;
903  const unsigned size = m_fad_size.value;
904 #endif
905  return reference_type( beg + index
906  , beg + m_original_fad_size
907  , size
908  , strd ); }
909 
910  template< typename I0 , typename I1 , typename I2 , typename I3 >
911  KOKKOS_FORCEINLINE_FUNCTION
913  !is_layout_left, reference_type>::type
914  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3 ) const
915  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,0);
916 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
917  const unsigned index = threadIdx.x;
918  const unsigned strd = blockDim.x;
919  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
920 #else
921  const unsigned index = m_fad_index;
922  const unsigned strd = m_fad_stride;
923  const unsigned size = m_fad_size.value;
924 #endif
925  return reference_type( beg + index
926  , beg + m_original_fad_size
927  , size
928  , strd ); }
929 
930  template< typename I0 , typename I1 , typename I2 , typename I3
931  , typename I4 >
932  KOKKOS_FORCEINLINE_FUNCTION
934  is_layout_left, reference_type>::type
935  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
936  , const I4 & i4 ) const
937  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4);
938 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
939  const unsigned index = threadIdx.x;
940  const unsigned strd = blockDim.x;
941  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
942 #else
943  const unsigned index = m_fad_index;
944  const unsigned strd = m_fad_stride;
945  const unsigned size = m_fad_size.value;
946 #endif
947  return reference_type( beg + index
948  , beg + m_original_fad_size
949  , size
950  , strd ); }
951 
952  template< typename I0 , typename I1 , typename I2 , typename I3
953  , typename I4 >
954  KOKKOS_FORCEINLINE_FUNCTION
956  !is_layout_left, reference_type>::type
957  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
958  , const I4 & i4 ) const
959  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,0);
960 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
961  const unsigned index = threadIdx.x;
962  const unsigned strd = blockDim.x;
963  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
964 #else
965  const unsigned index = m_fad_index;
966  const unsigned strd = m_fad_stride;
967  const unsigned size = m_fad_size.value;
968 #endif
969  return reference_type( beg + index
970  , beg + m_original_fad_size
971  , size
972  , strd ); }
973 
974  template< typename I0 , typename I1 , typename I2 , typename I3
975  , typename I4 , typename I5 >
976  KOKKOS_FORCEINLINE_FUNCTION
978  is_layout_left, reference_type>::type
979  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
980  , const I4 & i4 , const I5 & i5 ) const
981  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5);
982 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
983  const unsigned index = threadIdx.x;
984  const unsigned strd = blockDim.x;
985  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
986 #else
987  const unsigned index = m_fad_index;
988  const unsigned strd = m_fad_stride;
989  const unsigned size = m_fad_size.value;
990 #endif
991  return reference_type( beg + index
992  , beg + m_original_fad_size
993  , size
994  , strd ); }
995 
996  template< typename I0 , typename I1 , typename I2 , typename I3
997  , typename I4 , typename I5 >
998  KOKKOS_FORCEINLINE_FUNCTION
1000  !is_layout_left, reference_type>::type
1001  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
1002  , const I4 & i4 , const I5 & i5 ) const
1003  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,0);
1004 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
1005  const unsigned index = threadIdx.x;
1006  const unsigned strd = blockDim.x;
1007  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1008 #else
1009  const unsigned index = m_fad_index;
1010  const unsigned strd = m_fad_stride;
1011  const unsigned size = m_fad_size.value;
1012 #endif
1013  return reference_type( beg + index
1014  , beg + m_original_fad_size
1015  , size
1016  , strd ); }
1017 
1018  template< typename I0 , typename I1 , typename I2 , typename I3
1019  , typename I4 , typename I5 , typename I6 >
1020  KOKKOS_FORCEINLINE_FUNCTION
1022  is_layout_left, reference_type>::type
1023  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
1024  , const I4 & i4 , const I5 & i5 , const I6 & i6 ) const
1025  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5,i6);
1026 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
1027  const unsigned index = threadIdx.x;
1028  const unsigned strd = blockDim.x;
1029  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1030 #else
1031  const unsigned index = m_fad_index;
1032  const unsigned strd = m_fad_stride;
1033  const unsigned size = m_fad_size.value;
1034 #endif
1035  return reference_type( beg + index
1036  , beg + m_original_fad_size
1037  , size
1038  , strd ); }
1039 
1040  template< typename I0 , typename I1 , typename I2 , typename I3
1041  , typename I4 , typename I5 , typename I6 >
1042  KOKKOS_FORCEINLINE_FUNCTION
1044  !is_layout_left, reference_type>::type
1045  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
1046  , const I4 & i4 , const I5 & i5 , const I6 & i6 ) const
1047  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,i6,0);
1048 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
1049  const unsigned index = threadIdx.x;
1050  const unsigned strd = blockDim.x;
1051  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1052 #else
1053  const unsigned index = m_fad_index;
1054  const unsigned strd = m_fad_stride;
1055  const unsigned size = m_fad_size.value;
1056 #endif
1057  return reference_type( beg + index
1058  , beg + m_original_fad_size
1059  , size
1060  , strd ); }
1061 
1062  //----------------------------------------
1063 
1065  KOKKOS_INLINE_FUNCTION
1066  static constexpr size_t memory_span( typename Traits::array_layout const & layout )
1067  {
1068  // Do not introduce padding...
1069  typedef std::integral_constant< unsigned , 0 > padding ;
1070  return array_offset_type(
1071  padding() ,
1072  create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( layout ) ).span() * sizeof(fad_value_type);
1073  }
1074 
1075  //----------------------------------------
1076 
1077  KOKKOS_DEFAULTED_FUNCTION ~ViewMapping() = default ;
1078  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) {}
1079 
1080  KOKKOS_DEFAULTED_FUNCTION ViewMapping( const ViewMapping & ) = default ;
1081  KOKKOS_DEFAULTED_FUNCTION ViewMapping & operator = ( const ViewMapping & ) = default ;
1082 
1083  KOKKOS_DEFAULTED_FUNCTION ViewMapping( ViewMapping && ) = default ;
1084  KOKKOS_DEFAULTED_FUNCTION ViewMapping & operator = ( ViewMapping && ) = default ;
1085 
1086  template< class ... P >
1087  KOKKOS_INLINE_FUNCTION
1088  ViewMapping
1089  ( ViewCtorProp< P ... > const & prop
1090  , typename Traits::array_layout const & local_layout
1091  )
1092  : m_impl_handle( ( (ViewCtorProp<void,pointer_type> const &) prop ).value )
1093  , m_impl_offset( std::integral_constant< unsigned , 0 >()
1094  , local_layout )
1095  , m_array_offset(
1096  std::integral_constant< unsigned , 0 >()
1097  , create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( local_layout ) )
1098  , m_fad_size( getFadDimension<unsigned(Rank)>( m_array_offset ) - 1 )
1099  , m_original_fad_size( m_fad_size.value )
1100  , m_fad_stride( 1 )
1101  , m_fad_index( 0 )
1102  {
1103  const unsigned fad_dim =
1104  getFadDimension<unsigned(Rank)>( m_array_offset );
1105  if (unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1106  Kokkos::abort("invalid fad dimension (0) supplied!");
1107  }
1108 
1109  //----------------------------------------
1110  /* Allocate and construct mapped array.
1111  * Allocate via shared allocation record and
1112  * return that record for allocation tracking.
1113  */
1114  template< class ... P >
1115  SharedAllocationRecord<> *
1116  allocate_shared( ViewCtorProp< P... > const & prop
1117  , typename Traits::array_layout const & local_layout
1118  , bool execution_space_specified)
1119  {
1120  typedef ViewCtorProp< P... > ctor_prop ;
1121 
1122  typedef typename ctor_prop::execution_space execution_space ;
1123  typedef typename Traits::memory_space memory_space ;
1124  typedef ViewValueFunctor< execution_space , fad_value_type > functor_type ;
1125  typedef SharedAllocationRecord< memory_space , functor_type > record_type ;
1126 
1127  // Disallow padding
1128  typedef std::integral_constant< unsigned , 0 > padding ;
1129 
1130  // Check if ViewCtorProp has CommonViewAllocProp - if so, retrieve the fad_size and append to layout
1131  enum { test_traits_check = Kokkos::Impl::check_has_common_view_alloc_prop< P... >::value };
1132 
1133  typename Traits::array_layout internal_layout =
1134  (test_traits_check == true)
1135  ? Kokkos::Impl::appendFadToLayoutViewAllocHelper< Traits, P... >::returnNewLayoutPlusFad(prop, local_layout)
1136  : local_layout;
1137 
1138  m_impl_offset = offset_type( padding(), internal_layout );
1139 
1140  m_array_offset =
1141  array_offset_type( padding() ,
1142  create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( internal_layout ) );
1143  const unsigned fad_dim =
1144  getFadDimension<unsigned(Rank)>( m_array_offset );
1145  if (unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1146  Kokkos::abort("invalid fad dimension (0) supplied!");
1147  m_fad_size = fad_dim - 1 ;
1148  m_original_fad_size = m_fad_size.value ;
1149  m_fad_stride = 1;
1150  m_fad_index = 0;
1151 
1152  const size_t alloc_size = m_array_offset.span() * sizeof(fad_value_type);
1153 
1154  // Create shared memory tracking record with allocate memory from the memory space
1155  record_type * const record =
1156  record_type::allocate( ( (ViewCtorProp<void,memory_space> const &) prop ).value
1157  , ( (ViewCtorProp<void,std::string> const &) prop ).value
1158  , alloc_size );
1159 
1160  // Only set the the pointer and initialize if the allocation is non-zero.
1161  // May be zero if one of the dimensions is zero.
1162  if ( alloc_size ) {
1163 
1164  m_impl_handle = handle_type( reinterpret_cast< pointer_type >( record->data() ) );
1165 
1166  if ( ctor_prop::initialize ) {
1167  // Assume destruction is only required when construction is requested.
1168  // The ViewValueFunctor has both value construction and destruction operators.
1169  if (execution_space_specified)
1170  record->m_destroy = functor_type( ( (ViewCtorProp<void,execution_space> const &) prop).value
1171  , (fad_value_type *) m_impl_handle
1172  , m_array_offset.span()
1173  , record->get_label()
1174  );
1175  else
1176  record->m_destroy = functor_type((fad_value_type *) m_impl_handle
1177  , m_array_offset.span()
1178  , record->get_label()
1179  );
1180 
1181  // Construct values
1182  record->m_destroy.construct_shared_allocation();
1183  }
1184  }
1185 
1186  return record ;
1187  }
1188 
1189 };
1190 
1191 } // namespace Impl
1192 } // namespace Kokkos
1193 
1194 //----------------------------------------------------------------------------
1195 
1196 namespace Kokkos {
1197 namespace Impl {
1198 
1203 template< class DstTraits , class SrcTraits >
1204 class ViewMapping< DstTraits , SrcTraits ,
1205  typename std::enable_if<(
1206  Kokkos::Impl::MemorySpaceAccess
1207  < typename DstTraits::memory_space
1208  , typename SrcTraits::memory_space >::assignable
1209  &&
1210  // Destination view has FAD
1211  std::is_same< typename DstTraits::specialize
1212  , ViewSpecializeSacadoFadContiguous >::value
1213  &&
1214  // Source view has FAD
1215  std::is_same< typename SrcTraits::specialize
1216  , ViewSpecializeSacadoFadContiguous >::value
1217  )
1218  , typename DstTraits::specialize
1219  >::type >
1220 {
1221 public:
1222 
1223  enum { is_assignable = true };
1224  enum { is_assignable_data_type = true };
1225 
1226  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1227  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1228  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1229 
1230  template< class DstType >
1231  KOKKOS_INLINE_FUNCTION static
1232  void assign( DstType & dst
1233  , const SrcFadType & src
1234  , const TrackType & )
1235  {
1236  static_assert(
1237  (
1238  std::is_same< typename DstTraits::array_layout
1239  , Kokkos::LayoutLeft >::value ||
1240  std::is_same< typename DstTraits::array_layout
1241  , Kokkos::LayoutRight >::value ||
1242  std::is_same< typename DstTraits::array_layout
1243  , Kokkos::LayoutStride >::value
1244  )
1245  &&
1246  (
1247  std::is_same< typename SrcTraits::array_layout
1248  , Kokkos::LayoutLeft >::value ||
1249  std::is_same< typename SrcTraits::array_layout
1250  , Kokkos::LayoutRight >::value ||
1251  std::is_same< typename SrcTraits::array_layout
1252  , Kokkos::LayoutStride >::value
1253  )
1254  , "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1255 
1256  static_assert(
1257  std::is_same< typename DstTraits::array_layout
1258  , typename SrcTraits::array_layout >::value ||
1259  std::is_same< typename DstTraits::array_layout
1260  , Kokkos::LayoutStride >::value ,
1261  "View assignment must have compatible layout" );
1262 
1263  static_assert(
1264  std::is_same< typename DstTraits::value_type
1265  , typename SrcTraits::value_type >::value ||
1266  std::is_same< typename DstTraits::value_type
1267  , typename SrcTraits::const_value_type >::value ,
1268  "View assignment must have same value type or const = non-const" );
1269 
1270  static_assert(
1271  ViewDimensionAssignable
1272  < typename DstType::offset_type::dimension_type
1273  , typename SrcFadType::offset_type::dimension_type >::value ,
1274  "View assignment must have compatible dimensions" );
1275 
1276  static_assert(
1277  ViewDimensionAssignable
1278  < typename DstType::array_offset_type::dimension_type
1279  , typename SrcFadType::array_offset_type::dimension_type >::value ,
1280  "View assignment must have compatible dimensions" );
1281 
1282  typedef typename DstType::offset_type dst_offset_type ;
1283  typedef typename DstType::array_offset_type dst_array_offset_type ;
1284 
1285  dst.m_impl_handle = src.m_impl_handle ;
1286  dst.m_impl_offset = dst_offset_type( src.m_impl_offset );
1287  dst.m_array_offset = dst_array_offset_type( src.m_array_offset );
1288  dst.m_fad_size = src.m_fad_size.value ;
1289  dst.m_original_fad_size = src.m_original_fad_size ;
1290  dst.m_fad_stride = src.m_fad_stride ;
1291  dst.m_fad_index = src.m_fad_index ;
1292  }
1293 };
1294 
1299 template< class DstTraits , class SrcTraits >
1300 class ViewMapping< DstTraits , SrcTraits ,
1301  typename std::enable_if<(
1302  std::is_same< typename DstTraits::memory_space
1303  , typename SrcTraits::memory_space >::value
1304  &&
1305  // Destination view has FAD
1306  std::is_same< typename DstTraits::specialize
1307  , ViewSpecializeSacadoFad >::value
1308  &&
1309  // Source view has FAD contiguous
1310  std::is_same< typename SrcTraits::specialize
1311  , ViewSpecializeSacadoFadContiguous >::value
1312  &&
1313  // Destination view is LayoutStride
1314  std::is_same< typename DstTraits::array_layout
1315  , Kokkos::LayoutStride >::value
1316  )
1317  , typename DstTraits::specialize
1318  >::type >
1319 {
1320 public:
1321 
1322  enum { is_assignable = true };
1323  enum { is_assignable_data_type = true };
1324 
1325  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1326  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1327  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1328 
1329  template< class DstType >
1330  KOKKOS_INLINE_FUNCTION static
1331  void assign( DstType & dst
1332  , const SrcFadType & src
1333  , const TrackType & )
1334  {
1335  static_assert(
1336  std::is_same< typename SrcTraits::array_layout
1337  , Kokkos::LayoutLeft >::value ||
1338  std::is_same< typename SrcTraits::array_layout
1339  , Kokkos::LayoutRight >::value ||
1340  std::is_same< typename SrcTraits::array_layout
1341  , Kokkos::LayoutStride >::value ,
1342  "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1343 
1344  static_assert(
1345  std::is_same< typename DstTraits::value_type
1346  , typename SrcTraits::value_type >::value ||
1347  std::is_same< typename DstTraits::value_type
1348  , typename SrcTraits::const_value_type >::value ,
1349  "View assignment must have same value type or const = non-const" );
1350 
1351  static_assert(
1352  DstTraits::dimension::rank == SrcTraits::dimension::rank,
1353  "View assignment must have same rank" );
1354 
1355  typedef typename DstType::array_offset_type dst_offset_type ;
1356 
1357  dst.m_impl_handle = src.m_impl_handle ;
1358  dst.m_fad_size = src.m_fad_size.value ;
1359  dst.m_fad_stride = src.m_fad_stride ;
1360  dst.m_impl_offset = src.m_impl_offset;
1361 
1362  size_t N[8], S[8];
1363  N[0] = src.m_array_offset.dimension_0();
1364  N[1] = src.m_array_offset.dimension_1();
1365  N[2] = src.m_array_offset.dimension_2();
1366  N[3] = src.m_array_offset.dimension_3();
1367  N[4] = src.m_array_offset.dimension_4();
1368  N[5] = src.m_array_offset.dimension_5();
1369  N[6] = src.m_array_offset.dimension_6();
1370  N[7] = src.m_array_offset.dimension_7();
1371  S[0] = src.m_array_offset.stride_0();
1372  S[1] = src.m_array_offset.stride_1();
1373  S[2] = src.m_array_offset.stride_2();
1374  S[3] = src.m_array_offset.stride_3();
1375  S[4] = src.m_array_offset.stride_4();
1376  S[5] = src.m_array_offset.stride_5();
1377  S[6] = src.m_array_offset.stride_6();
1378  S[7] = src.m_array_offset.stride_7();
1379 
1380  // For LayoutLeft, we have to move the Sacado dimension from the first
1381  // to the last
1382  if (std::is_same< typename SrcTraits::array_layout
1383  , Kokkos::LayoutLeft >::value)
1384  {
1385  const size_t N_fad = N[0];
1386  const size_t S_fad = S[0];
1387  for (int i=0; i<7; ++i) {
1388  N[i] = N[i+1];
1389  S[i] = S[i+1];
1390  }
1391  N[DstTraits::dimension::rank] = N_fad;
1392  S[DstTraits::dimension::rank] = S_fad;
1393  }
1394  Kokkos::LayoutStride ls( N[0], S[0],
1395  N[1], S[1],
1396  N[2], S[2],
1397  N[3], S[3],
1398  N[4], S[4],
1399  N[5], S[5],
1400  N[6], S[6],
1401  N[7], S[7] );
1402  dst.m_array_offset = dst_offset_type(std::integral_constant<unsigned,0>(), ls);
1403  }
1404 };
1405 
1410 template< class DstTraits , class SrcTraits >
1411 class ViewMapping< DstTraits , SrcTraits ,
1412  typename std::enable_if<(
1413  std::is_same< typename DstTraits::memory_space
1414  , typename SrcTraits::memory_space >::value
1415  &&
1416  // Destination view has ordinary
1417  std::is_same< typename DstTraits::specialize , void >::value
1418  &&
1419  // Source view has FAD only
1420  std::is_same< typename SrcTraits::specialize
1421  , ViewSpecializeSacadoFadContiguous >::value
1422  )
1423  , typename DstTraits::specialize
1424  >::type >
1425 {
1426 public:
1427 
1428  enum { is_assignable = true };
1429  enum { is_assignable_data_type = true };
1430 
1431  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1432  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1433  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1434 
1435 
1436  // Helpers to assign, and generate if necessary, ViewOffset to the dst map
1437  // These are necessary to use Kokkos' deep_copy with nested fads
1438  template < class DstType, class SrcFadType, class Enable = void >
1439  struct AssignOffset;
1440 
1441  template < class DstType, class SrcFadType >
1442  struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank != (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1443  {
1444  // ViewOffset's Dimensions Ranks do not match
1445  KOKKOS_INLINE_FUNCTION
1446  static void assign( DstType & dst, const SrcFadType & src )
1447  {
1448  typedef typename SrcTraits::value_type TraitsValueType;
1449 
1452  )
1453  {
1454 
1455  typedef typename DstType::offset_type::array_layout DstLayoutType;
1456  //typedef typename ViewArrayLayoutSelector<typename DstType::offset_type::array_layout>::type DstLayoutType;
1457  typedef typename SrcFadType::array_offset_type::dimension_type SrcViewDimension;
1458 
1459  // This is the static dimension of the inner fad, missing from ViewDimension
1461 
1462  static constexpr bool is_layout_left =
1464 
1465  typedef typename std::conditional< is_layout_left,
1466  typename SrcViewDimension:: template prepend< InnerStaticDim+1 >::type,
1467  typename SrcViewDimension:: template append < InnerStaticDim+1 >::type
1468  >::type SrcViewDimensionAppended;
1469 
1470  typedef std::integral_constant< unsigned , 0 > padding ;
1471 
1472  typedef ViewOffset< SrcViewDimensionAppended, DstLayoutType > TmpOffsetType;
1473 
1474  auto src_layout = src.m_array_offset.layout();
1475 
1476  if ( is_layout_left ) {
1477  auto prepend_layout = Kokkos::Impl::prependFadToLayout< DstLayoutType >::returnNewLayoutPlusFad(src_layout, InnerStaticDim+1);
1478  TmpOffsetType offset_tmp( padding(), prepend_layout );
1479  dst.m_impl_offset = offset_tmp;
1480  }
1481  else {
1482  TmpOffsetType offset_tmp( padding(), src_layout );
1483  dst.m_impl_offset = offset_tmp;
1484  }
1485  } else {
1486  Kokkos::abort("Sacado error: Applying AssignOffset for case with nested Fads, but without nested Fads - something went wrong");
1487  }
1488  }
1489  };
1490 
1491  template < class DstType, class SrcFadType >
1492  struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank == (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1493  {
1494  KOKKOS_INLINE_FUNCTION
1495  static void assign( DstType & dst, const SrcFadType & src )
1496  {
1497  typedef typename DstType::offset_type dst_offset_type ;
1498  dst.m_impl_offset = dst_offset_type( src.m_array_offset );
1499  }
1500  };
1501 
1502  template< class DstType >
1503  KOKKOS_INLINE_FUNCTION static
1504  void assign( DstType & dst
1505  , const SrcFadType & src
1506  , const TrackType &
1507  )
1508  {
1509  static_assert(
1510  (
1511  std::is_same< typename DstTraits::array_layout
1512  , Kokkos::LayoutLeft >::value ||
1513  std::is_same< typename DstTraits::array_layout
1514  , Kokkos::LayoutRight >::value ||
1515  std::is_same< typename DstTraits::array_layout
1516  , Kokkos::LayoutStride >::value
1517  )
1518  &&
1519  (
1520  std::is_same< typename SrcTraits::array_layout
1521  , Kokkos::LayoutLeft >::value ||
1522  std::is_same< typename SrcTraits::array_layout
1523  , Kokkos::LayoutRight >::value ||
1524  std::is_same< typename SrcTraits::array_layout
1525  , Kokkos::LayoutStride >::value
1526  )
1527  , "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1528 
1529  static_assert(
1530  std::is_same< typename DstTraits::array_layout
1531  , typename SrcTraits::array_layout >::value ||
1532  std::is_same< typename DstTraits::array_layout
1533  , Kokkos::LayoutStride >::value ,
1534  "View assignment must have compatible layout" );
1535 
1536  if ( src.m_fad_index != 0 || src.m_fad_stride != 1 ) {
1537  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Cannot assign to array with partitioned view ******\n\n");
1538  }
1539 
1540  AssignOffset< DstType, SrcFadType >::assign( dst, src );
1541  dst.m_impl_handle = reinterpret_cast< typename DstType::handle_type >(src.m_impl_handle) ;
1542  }
1543 };
1544 
1545 } // namespace Impl
1546 } // namespace Kokkos
1547 
1548 //----------------------------------------------------------------------------
1549 
1550 namespace Kokkos {
1551 namespace Impl {
1552 
1553 // Rules for subview arguments and layouts matching
1554 
1555 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1556 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1557  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1558 };
1559 
1560 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1561 struct SubviewLegalArgsCompileTime<LayoutDest,Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1562  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1563 };
1564 
1565 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1566 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1567  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1568 };
1569 
1570 // Subview mapping
1571 
1572 template< class SrcTraits , class Arg0 , class ... Args >
1573 struct ViewMapping
1574  < typename std::enable_if<(
1575  // Source view has FAD only
1576  std::is_same< typename SrcTraits::specialize
1577  , ViewSpecializeSacadoFadContiguous >::value
1578  &&
1579  (
1580  std::is_same< typename SrcTraits::array_layout
1581  , Kokkos::LayoutLeft >::value ||
1582  std::is_same< typename SrcTraits::array_layout
1583  , Kokkos::LayoutRight >::value ||
1584  std::is_same< typename SrcTraits::array_layout
1585  , Kokkos::LayoutStride >::value
1586  )
1587  && !Sacado::Fad::is_fad_partition<Arg0>::value
1588  )
1589  >::type
1590  , SrcTraits
1591  , Arg0, Args ... >
1592 {
1593 private:
1594 
1595  static_assert( SrcTraits::rank == sizeof...(Args)+1 , "" );
1596 
1597  enum
1598  { RZ = false
1606  };
1607 
1608  // Public rank
1609  enum { rank = unsigned(R0) + unsigned(R1) + unsigned(R2) + unsigned(R3)
1610  + unsigned(R4) + unsigned(R5) + unsigned(R6) };
1611 
1612  // Whether right-most non-FAD rank is a range.
1613  enum { R0_rev = ( 0 == SrcTraits::rank ? RZ : (
1614  1 == SrcTraits::rank ? R0 : (
1615  2 == SrcTraits::rank ? R1 : (
1616  3 == SrcTraits::rank ? R2 : (
1617  4 == SrcTraits::rank ? R3 : (
1618  5 == SrcTraits::rank ? R4 : (
1619  6 == SrcTraits::rank ? R5 : R6 ))))))) };
1620 
1621  // Subview's layout
1622  typedef typename std::conditional<
1623  ( /* Same array layout IF */
1624  ( rank == 0 ) /* output rank zero */
1625  ||
1626  // OutputRank 1, InputLayout Left, Interval 0
1627  // because single stride one
1629  ||
1630  // OutputRank 1, InputLayout Right, Interval [InputRank-1]
1631  // because single stride one
1634  >::type array_layout ;
1635 
1636  typedef typename SrcTraits::value_type fad_type ;
1637 
1638  typedef typename std::conditional< rank == 0 , fad_type ,
1639  typename std::conditional< rank == 1 , fad_type * ,
1640  typename std::conditional< rank == 2 , fad_type ** ,
1641  typename std::conditional< rank == 3 , fad_type *** ,
1642  typename std::conditional< rank == 4 , fad_type **** ,
1643  typename std::conditional< rank == 5 , fad_type ***** ,
1644  typename std::conditional< rank == 6 , fad_type ****** ,
1645  fad_type *******
1646  >::type >::type >::type >::type >::type >::type >::type
1647  data_type ;
1648 
1649 public:
1650 
1651  typedef Kokkos::ViewTraits
1652  < data_type
1653  , array_layout
1654  , typename SrcTraits::device_type
1655  , typename SrcTraits::memory_traits > traits_type ;
1656 
1657  typedef Kokkos::View
1658  < data_type
1659  , array_layout
1660  , typename SrcTraits::device_type
1661  , typename SrcTraits::memory_traits > type ;
1662 
1663 
1664  KOKKOS_INLINE_FUNCTION
1665  static void assign( ViewMapping< traits_type , typename traits_type::specialize > & dst
1666  , ViewMapping< SrcTraits , typename SrcTraits::specialize > const & src
1667  , Arg0 arg0 , Args ... args )
1668  {
1669  typedef ViewMapping< traits_type , typename traits_type::specialize > DstType ;
1670  typedef typename DstType::offset_type dst_offset_type ;
1671  typedef typename DstType::array_offset_type dst_array_offset_type ;
1672  typedef typename DstType::handle_type dst_handle_type ;
1673 
1674  size_t offset;
1676  const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1677  array_extents( src.m_array_offset.m_dim , Kokkos::ALL() , arg0 , args... );
1678  offset = src.m_array_offset( array_extents.domain_offset(0)
1679  , array_extents.domain_offset(1)
1680  , array_extents.domain_offset(2)
1681  , array_extents.domain_offset(3)
1682  , array_extents.domain_offset(4)
1683  , array_extents.domain_offset(5)
1684  , array_extents.domain_offset(6)
1685  , array_extents.domain_offset(7) );
1686  dst_array_offset_type dst_array_offset( src.m_array_offset ,
1687  array_extents );
1688  // For LayoutStride, we always use LayoutRight indexing (because we
1689  // don't know whether the original array was Left or Right), so we
1690  // need to swap the Fad dimension to the last and shift all of the
1691  // other dimensions left by 1
1693  {
1694  Kokkos::LayoutStride ls(
1695  dst_array_offset.m_dim.N0, dst_array_offset.m_stride.S0,
1696  dst_array_offset.m_dim.N1, dst_array_offset.m_stride.S1,
1697  dst_array_offset.m_dim.N2, dst_array_offset.m_stride.S2,
1698  dst_array_offset.m_dim.N3, dst_array_offset.m_stride.S3,
1699  dst_array_offset.m_dim.N4, dst_array_offset.m_stride.S4,
1700  dst_array_offset.m_dim.N5, dst_array_offset.m_stride.S5,
1701  dst_array_offset.m_dim.N6, dst_array_offset.m_stride.S6,
1702  dst_array_offset.m_dim.N7, dst_array_offset.m_stride.S7);
1703  auto t1 = ls.dimension[0];
1704  for (unsigned i=0; i<rank; ++i)
1705  ls.dimension[i] = ls.dimension[i+1];
1706  ls.dimension[rank] = t1;
1707  auto t2 = ls.stride[0];
1708  for (unsigned i=0; i<rank; ++i)
1709  ls.stride[i] = ls.stride[i+1];
1710  ls.stride[rank] = t2;
1711  dst.m_array_offset = dst_array_offset_type(std::integral_constant<unsigned, 0>(), ls);
1712  }
1713  else
1714  dst.m_array_offset = dst_array_offset;
1715  }
1716  else {
1717  const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1718  array_extents( src.m_array_offset.m_dim , arg0 , args... , Kokkos::ALL() );
1719  offset = src.m_array_offset( array_extents.domain_offset(0)
1720  , array_extents.domain_offset(1)
1721  , array_extents.domain_offset(2)
1722  , array_extents.domain_offset(3)
1723  , array_extents.domain_offset(4)
1724  , array_extents.domain_offset(5)
1725  , array_extents.domain_offset(6)
1726  , array_extents.domain_offset(7) );
1727  dst.m_array_offset = dst_array_offset_type( src.m_array_offset ,
1728  array_extents );
1729  }
1730 
1731  const SubviewExtents< SrcTraits::rank , rank >
1732  extents( src.m_impl_offset.m_dim , arg0 , args... );
1733 
1734  dst.m_impl_offset = dst_offset_type( src.m_impl_offset , extents );
1735  dst.m_impl_handle = dst_handle_type( src.m_impl_handle + offset );
1736  dst.m_fad_size = src.m_fad_size;
1737  dst.m_original_fad_size = src.m_original_fad_size;
1738  dst.m_fad_stride = src.m_fad_stride;
1739  dst.m_fad_index = src.m_fad_index;
1740  }
1741 
1742 };
1743 
1744 } // namespace Impl
1745 } // namespace Kokkos
1746 
1747 //---------------------------------------------------------------------------
1748 
1749 namespace Kokkos {
1750 namespace Impl {
1751 
1752 // Partition mapping
1753 
1754 template< class DataType, class ...P, unsigned Stride >
1755 class ViewMapping<
1756  void,
1757  ViewTraits<DataType,P...> ,
1758  Sacado::Fad::Partition<Stride>
1759  >
1760 {
1761 public:
1762 
1763  enum { is_assignable = true };
1764  enum { is_assignable_data_type = true };
1765 
1766  typedef ViewTraits<DataType,P...> src_traits;
1767  typedef ViewMapping< src_traits , typename src_traits::specialize > src_type ;
1768 
1769  typedef typename src_type::offset_type::dimension_type src_dimension;
1770  typedef typename src_traits::value_type fad_type;
1771  typedef typename Sacado::LocalScalarType<fad_type,Stride>::type strided_fad_type;
1772  typedef typename
1773  ViewDataType< strided_fad_type , src_dimension >::type strided_data_type;
1774  typedef ViewTraits<strided_data_type,P...> dst_traits;
1775  typedef View<strided_data_type,P...> type;
1776  typedef ViewMapping< dst_traits , typename dst_traits::specialize > dst_type ;
1777 
1778  KOKKOS_INLINE_FUNCTION static
1779  void assign( dst_type & dst
1780  , const src_type & src
1781  , const Sacado::Fad::Partition<Stride> & part )
1782  {
1783  if ( Stride != part.stride && Stride != 0 ) {
1784  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Invalid size in partitioned view assignment ******\n\n");
1785  }
1786  if ( src.m_fad_stride != 1 ) {
1787  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Can't partition already partitioned view ******\n\n");
1788  }
1789 
1790  dst.m_impl_handle = src.m_impl_handle ;
1791  dst.m_impl_offset = src.m_impl_offset ;
1792  dst.m_array_offset = src.m_array_offset ;
1793 
1794  // Assuming the alignment was choosen correctly for the partitioning,
1795  // each partition should get the same size. This allows the use of SFad.
1796  dst.m_fad_size =
1797  (src.m_fad_size.value + part.stride-part.offset-1) / part.stride ;
1798 
1799  dst.m_original_fad_size = src.m_original_fad_size ;
1800  dst.m_fad_stride = part.stride ;
1801  dst.m_fad_index = part.offset ;
1802  }
1803 };
1804 
1805 } // namespace Impl
1806 } // namespace Kokkos
1807 
1808 #endif // defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
1809 
1810 #endif // defined(HAVE_SACADO_KOKKOS)
1811 
1812 #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
const int N
void
Definition: uninit.c:105
const int fad_dim
int value
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.