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