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  enum { is_assignable_data_type = true };
540 
541  enum { FadStaticDimension = Sacado::StaticSize< fad_type >::value };
542  enum { PartitionedFadStride = Traits::array_layout::scalar_stride };
543 
544  // The partitioned static size -- this will be 0 if ParitionedFadStride
545  // does not evenly divide FadStaticDimension
546  enum { PartitionedFadStaticDimension =
547  computeFadPartitionSize(FadStaticDimension,PartitionedFadStride) };
548 
549 #ifdef KOKKOS_ENABLE_CUDA
550  typedef typename Sacado::LocalScalarType< fad_type, unsigned(PartitionedFadStride) >::type strided_scalar_type;
551  typedef typename std::conditional< std::is_same<typename Traits::execution_space, Kokkos::Cuda>::value, strided_scalar_type, fad_type >::type thread_local_scalar_type;
552 #else
553  typedef fad_type thread_local_scalar_type;
554 #endif
555 
556 private:
558 
559  typedef fad_value_type * handle_type ;
560 
561  typedef ViewArrayAnalysis< typename Traits::data_type > array_analysis ;
562 
563  typedef ViewOffset< typename Traits::dimension
564  , typename Traits::array_layout
565  , void
566  > offset_type ;
567 
568  // Prepend/append the fad dimension for the internal offset mapping.
569  static constexpr bool is_layout_left =
570  std::is_same< typename Traits::array_layout, Kokkos::LayoutLeft>::value;
571  typedef ViewOffset<
572  typename std::conditional<
573  is_layout_left,
574  typename array_analysis::dimension::
575  template prepend<0>::type,
576  typename array_analysis::dimension::
577  template append<0>::type >::type,
578  typename Traits::array_layout,
579  void >
580  array_offset_type ;
581 
582  handle_type m_impl_handle ;
583  offset_type m_impl_offset ;
584  array_offset_type m_array_offset ;
585  sacado_size_type m_fad_size ;
586 
587  // These are for manual partitioning, and will likely be removed
588  unsigned m_original_fad_size ;
589  unsigned m_fad_stride ;
590  unsigned m_fad_index ;
591 
592 public:
593 
594  //----------------------------------------
595  // Domain dimensions
596 
597  enum { Rank = Traits::dimension::rank };
598 
599  // Rank corresponding to the sacado dimension
600  enum { Sacado_Rank = std::is_same< typename Traits::array_layout, Kokkos::LayoutLeft >::value ? 0 : Rank+1 };
601 
602  // Using the internal offset mapping so limit to public rank:
603  template< typename iType >
604  KOKKOS_INLINE_FUNCTION constexpr size_t extent( const iType & r ) const
605  { return m_impl_offset.m_dim.extent(r); }
606 
607  KOKKOS_INLINE_FUNCTION constexpr
608  typename Traits::array_layout layout() const
609  { return m_impl_offset.layout(); }
610 
611  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_0() const
612  { return m_impl_offset.dimension_0(); }
613  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_1() const
614  { return m_impl_offset.dimension_1(); }
615  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_2() const
616  { return m_impl_offset.dimension_2(); }
617  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_3() const
618  { return m_impl_offset.dimension_3(); }
619  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_4() const
620  { return m_impl_offset.dimension_4(); }
621  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_5() const
622  { return m_impl_offset.dimension_5(); }
623  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_6() const
624  { return m_impl_offset.dimension_6(); }
625  KOKKOS_INLINE_FUNCTION constexpr size_t dimension_7() const
626  { return m_impl_offset.dimension_7(); }
627 
628  // Is a regular layout with uniform striding for each index.
629  // Since we allow for striding within the data type, we can't guarantee
630  // regular striding
631  using is_regular = std::false_type ;
632 
633  // FIXME: Adjust these for m_stride
634  KOKKOS_INLINE_FUNCTION constexpr size_t stride_0() const
635  { return m_impl_offset.stride_0(); }
636  KOKKOS_INLINE_FUNCTION constexpr size_t stride_1() const
637  { return m_impl_offset.stride_1(); }
638  KOKKOS_INLINE_FUNCTION constexpr size_t stride_2() const
639  { return m_impl_offset.stride_2(); }
640  KOKKOS_INLINE_FUNCTION constexpr size_t stride_3() const
641  { return m_impl_offset.stride_3(); }
642  KOKKOS_INLINE_FUNCTION constexpr size_t stride_4() const
643  { return m_impl_offset.stride_4(); }
644  KOKKOS_INLINE_FUNCTION constexpr size_t stride_5() const
645  { return m_impl_offset.stride_5(); }
646  KOKKOS_INLINE_FUNCTION constexpr size_t stride_6() const
647  { return m_impl_offset.stride_6(); }
648  KOKKOS_INLINE_FUNCTION constexpr size_t stride_7() const
649  { return m_impl_offset.stride_7(); }
650 
651  template< typename iType >
652  KOKKOS_INLINE_FUNCTION void stride( iType * const s ) const
653  { m_impl_offset.stride(s); }
654 
655  // Size of sacado scalar dimension
656  KOKKOS_FORCEINLINE_FUNCTION constexpr unsigned dimension_scalar() const
657 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
658  { return PartitionedFadStaticDimension ? PartitionedFadStaticDimension+1 : (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x + 1; }
659 #else
660  { return m_fad_size.value+1; }
661 #endif
662 
663  // trode of sacado scalar dimension
664  KOKKOS_FORCEINLINE_FUNCTION constexpr unsigned stride_scalar() const
665  { return m_fad_stride; }
666 
667  //----------------------------------------
668  // Range of mapping
669 
670 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
671  // Return type of reference operators
672  // this only works if you are using a team-parallel operation on Cuda!
673  // typedef typename
674  // Sacado::ViewFadType< thread_local_scalar_type , PartitionedFadStaticDimension , (unsigned(ParitionedFadStride) > 1 ? PartitionedFadStride : 0) >::type reference_type ;
675  typedef typename
677 #else
678  // Return type of reference operators
679  typedef typename
681 #endif
682 
684  typedef fad_value_type * pointer_type ;
685 
687  KOKKOS_INLINE_FUNCTION constexpr size_t span() const
688  { return m_array_offset.span(); }
689 
691  KOKKOS_INLINE_FUNCTION constexpr bool span_is_contiguous() const
692  { return m_array_offset.span_is_contiguous() && (m_fad_stride == 1); }
693 
695  KOKKOS_INLINE_FUNCTION constexpr pointer_type data() const
696 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
697  { return m_impl_handle + threadIdx.x; }
698 #else
699  { return m_impl_handle + m_fad_index; }
700 #endif
701 
702  //----------------------------------------
703 
705  reference_type
706  reference() const
707  {
708 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
709  const unsigned index = threadIdx.x;
710  const unsigned strd = blockDim.x;
711  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
712 #else
713  const unsigned index = m_fad_index;
714  const unsigned strd = m_fad_stride;
715  const unsigned size = m_fad_size.value;
716 #endif
717  return reference_type( m_impl_handle + index
718  , m_impl_handle + m_original_fad_size
719  , size
720  , strd ); }
721 
722  template< typename I0 >
724  typename std::enable_if< Kokkos::Impl::are_integral<I0>::value &&
725  is_layout_left, reference_type>::type
726  reference( const I0 & i0 ) const
727  { pointer_type beg = m_impl_handle + m_array_offset(0,i0);
728 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
729  const unsigned index = threadIdx.x;
730  const unsigned strd = blockDim.x;
731  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
732 #else
733  const unsigned index = m_fad_index;
734  const unsigned strd = m_fad_stride;
735  const unsigned size = m_fad_size.value;
736 #endif
737  return reference_type( beg + index
738  , beg + m_original_fad_size
739  , size
740  , strd ); }
741 
742  template< typename I0 >
744  typename std::enable_if< Kokkos::Impl::are_integral<I0>::value &&
745  !is_layout_left, reference_type>::type
746  reference( const I0 & i0 ) const
747  { pointer_type beg = m_impl_handle + m_array_offset(i0,0);
748 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
749  const unsigned index = threadIdx.x;
750  const unsigned strd = blockDim.x;
751  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
752 #else
753  const unsigned index = m_fad_index;
754  const unsigned strd = m_fad_stride;
755  const unsigned size = m_fad_size.value;
756 #endif
757  return reference_type( beg + index
758  , beg + m_original_fad_size
759  , size
760  , strd ); }
761 
762  template< typename I0 , typename I1 >
764  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1>::value &&
765  is_layout_left, reference_type>::type
766  reference( const I0 & i0 , const I1 & i1 ) const
767  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1);
768 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
769  const unsigned index = threadIdx.x;
770  const unsigned strd = blockDim.x;
771  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
772 #else
773  const unsigned index = m_fad_index;
774  const unsigned strd = m_fad_stride;
775  const unsigned size = m_fad_size.value;
776 #endif
777  return reference_type( beg + index
778  , beg + m_original_fad_size
779  , size
780  , strd ); }
781 
782  template< typename I0 , typename I1 >
784  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1>::value &&
785  !is_layout_left, reference_type>::type
786  reference( const I0 & i0 , const I1 & i1 ) const
787  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,0);
788 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
789  const unsigned index = threadIdx.x;
790  const unsigned strd = blockDim.x;
791  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
792 #else
793  const unsigned index = m_fad_index;
794  const unsigned strd = m_fad_stride;
795  const unsigned size = m_fad_size.value;
796 #endif
797  return reference_type( beg + index
798  , beg + m_original_fad_size
799  , size
800  , strd ); }
801 
802 
803  template< typename I0 , typename I1 , typename I2 >
805  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2>::value &&
806  is_layout_left, reference_type>::type
807  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 ) const
808  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2);
809 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
810  const unsigned index = threadIdx.x;
811  const unsigned strd = blockDim.x;
812  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
813 #else
814  const unsigned index = m_fad_index;
815  const unsigned strd = m_fad_stride;
816  const unsigned size = m_fad_size.value;
817 #endif
818  return reference_type( beg + index
819  , beg + m_original_fad_size
820  , size
821  , strd ); }
822 
823  template< typename I0 , typename I1 , typename I2 >
825  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2>::value &&
826  !is_layout_left, reference_type>::type
827  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 ) const
828  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,0);
829 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
830  const unsigned index = threadIdx.x;
831  const unsigned strd = blockDim.x;
832  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
833 #else
834  const unsigned index = m_fad_index;
835  const unsigned strd = m_fad_stride;
836  const unsigned size = m_fad_size.value;
837 #endif
838  return reference_type( beg + index
839  , beg + m_original_fad_size
840  , size
841  , strd ); }
842 
843  template< typename I0 , typename I1 , typename I2 , typename I3 >
845  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3>::value &&
846  is_layout_left, reference_type>::type
847  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3 ) const
848  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3);
849 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
850  const unsigned index = threadIdx.x;
851  const unsigned strd = blockDim.x;
852  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
853 #else
854  const unsigned index = m_fad_index;
855  const unsigned strd = m_fad_stride;
856  const unsigned size = m_fad_size.value;
857 #endif
858  return reference_type( beg + index
859  , beg + m_original_fad_size
860  , size
861  , strd ); }
862 
863  template< typename I0 , typename I1 , typename I2 , typename I3 >
865  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3>::value &&
866  !is_layout_left, reference_type>::type
867  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3 ) const
868  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,0);
869 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
870  const unsigned index = threadIdx.x;
871  const unsigned strd = blockDim.x;
872  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
873 #else
874  const unsigned index = m_fad_index;
875  const unsigned strd = m_fad_stride;
876  const unsigned size = m_fad_size.value;
877 #endif
878  return reference_type( beg + index
879  , beg + m_original_fad_size
880  , size
881  , strd ); }
882 
883  template< typename I0 , typename I1 , typename I2 , typename I3
884  , typename I4 >
886  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4>::value &&
887  is_layout_left, reference_type>::type
888  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
889  , const I4 & i4 ) const
890  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4);
891 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
892  const unsigned index = threadIdx.x;
893  const unsigned strd = blockDim.x;
894  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
895 #else
896  const unsigned index = m_fad_index;
897  const unsigned strd = m_fad_stride;
898  const unsigned size = m_fad_size.value;
899 #endif
900  return reference_type( beg + index
901  , beg + m_original_fad_size
902  , size
903  , strd ); }
904 
905  template< typename I0 , typename I1 , typename I2 , typename I3
906  , typename I4 >
908  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4>::value &&
909  !is_layout_left, reference_type>::type
910  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
911  , const I4 & i4 ) const
912  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,0);
913 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
914  const unsigned index = threadIdx.x;
915  const unsigned strd = blockDim.x;
916  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
917 #else
918  const unsigned index = m_fad_index;
919  const unsigned strd = m_fad_stride;
920  const unsigned size = m_fad_size.value;
921 #endif
922  return reference_type( beg + index
923  , beg + m_original_fad_size
924  , size
925  , strd ); }
926 
927  template< typename I0 , typename I1 , typename I2 , typename I3
928  , typename I4 , typename I5 >
930  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5>::value &&
931  is_layout_left, reference_type>::type
932  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
933  , const I4 & i4 , const I5 & i5 ) const
934  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5);
935 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
936  const unsigned index = threadIdx.x;
937  const unsigned strd = blockDim.x;
938  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
939 #else
940  const unsigned index = m_fad_index;
941  const unsigned strd = m_fad_stride;
942  const unsigned size = m_fad_size.value;
943 #endif
944  return reference_type( beg + index
945  , beg + m_original_fad_size
946  , size
947  , strd ); }
948 
949  template< typename I0 , typename I1 , typename I2 , typename I3
950  , typename I4 , typename I5 >
952  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5>::value &&
953  !is_layout_left, reference_type>::type
954  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
955  , const I4 & i4 , const I5 & i5 ) const
956  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,0);
957 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
958  const unsigned index = threadIdx.x;
959  const unsigned strd = blockDim.x;
960  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
961 #else
962  const unsigned index = m_fad_index;
963  const unsigned strd = m_fad_stride;
964  const unsigned size = m_fad_size.value;
965 #endif
966  return reference_type( beg + index
967  , beg + m_original_fad_size
968  , size
969  , strd ); }
970 
971  template< typename I0 , typename I1 , typename I2 , typename I3
972  , typename I4 , typename I5 , typename I6 >
974  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5,I6>::value &&
975  is_layout_left, reference_type>::type
976  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
977  , const I4 & i4 , const I5 & i5 , const I6 & i6 ) const
978  { pointer_type beg = m_impl_handle + m_array_offset(0,i0,i1,i2,i3,i4,i5,i6);
979 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
980  const unsigned index = threadIdx.x;
981  const unsigned strd = blockDim.x;
982  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
983 #else
984  const unsigned index = m_fad_index;
985  const unsigned strd = m_fad_stride;
986  const unsigned size = m_fad_size.value;
987 #endif
988  return reference_type( beg + index
989  , beg + m_original_fad_size
990  , size
991  , strd ); }
992 
993  template< typename I0 , typename I1 , typename I2 , typename I3
994  , typename I4 , typename I5 , typename I6 >
996  typename std::enable_if< Kokkos::Impl::are_integral<I0,I1,I2,I3,I4,I5,I6>::value &&
997  !is_layout_left, reference_type>::type
998  reference( const I0 & i0 , const I1 & i1 , const I2 & i2 , const I3 & i3
999  , const I4 & i4 , const I5 & i5 , const I6 & i6 ) const
1000  { pointer_type beg = m_impl_handle + m_array_offset(i0,i1,i2,i3,i4,i5,i6,0);
1001 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) && defined(__CUDA_ARCH__)
1002  const unsigned index = threadIdx.x;
1003  const unsigned strd = blockDim.x;
1004  const unsigned size = (m_fad_size.value+blockDim.x-threadIdx.x-1) / blockDim.x;
1005 #else
1006  const unsigned index = m_fad_index;
1007  const unsigned strd = m_fad_stride;
1008  const unsigned size = m_fad_size.value;
1009 #endif
1010  return reference_type( beg + index
1011  , beg + m_original_fad_size
1012  , size
1013  , strd ); }
1014 
1015  //----------------------------------------
1016 
1019  static constexpr size_t memory_span( typename Traits::array_layout const & layout )
1020  {
1021  // Do not introduce padding...
1022  typedef std::integral_constant< unsigned , 0 > padding ;
1023  return array_offset_type(
1024  padding() ,
1025  create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( layout ) ).span() * sizeof(fad_value_type);
1026  }
1027 
1028  //----------------------------------------
1029 
1030  KOKKOS_DEFAULTED_FUNCTION ~ViewMapping() = default ;
1031  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) {}
1032 
1033  KOKKOS_DEFAULTED_FUNCTION ViewMapping( const ViewMapping & ) = default ;
1034  KOKKOS_DEFAULTED_FUNCTION ViewMapping & operator = ( const ViewMapping & ) = default ;
1035 
1036  KOKKOS_DEFAULTED_FUNCTION ViewMapping( ViewMapping && ) = default ;
1037  KOKKOS_DEFAULTED_FUNCTION ViewMapping & operator = ( ViewMapping && ) = default ;
1038 
1039  template< class ... P >
1041  ViewMapping
1042  ( ViewCtorProp< P ... > const & prop
1043  , typename Traits::array_layout const & local_layout
1044  )
1045  : m_impl_handle( ( (ViewCtorProp<void,pointer_type> const &) prop ).value )
1046  , m_impl_offset( std::integral_constant< unsigned , 0 >()
1047  , local_layout )
1048  , m_array_offset(
1049  std::integral_constant< unsigned , 0 >()
1050  , create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( local_layout ) )
1051  , m_fad_size( getFadDimension<unsigned(Rank)>( m_array_offset ) - 1 )
1052  , m_original_fad_size( m_fad_size.value )
1053  , m_fad_stride( 1 )
1054  , m_fad_index( 0 )
1055  {
1056  const unsigned fad_dim =
1057  getFadDimension<unsigned(Rank)>( m_array_offset );
1058  if (unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1059  Kokkos::abort("invalid fad dimension (0) supplied!");
1060  }
1061 
1062  //----------------------------------------
1063  /* Allocate and construct mapped array.
1064  * Allocate via shared allocation record and
1065  * return that record for allocation tracking.
1066  */
1067  template< class ... P >
1068  SharedAllocationRecord<> *
1069  allocate_shared( ViewCtorProp< P... > const & prop
1070  , typename Traits::array_layout const & local_layout )
1071  {
1072  typedef ViewCtorProp< P... > ctor_prop ;
1073 
1074  typedef typename ctor_prop::execution_space execution_space ;
1075  typedef typename Traits::memory_space memory_space ;
1076  typedef ViewValueFunctor< execution_space , fad_value_type > functor_type ;
1077  typedef SharedAllocationRecord< memory_space , functor_type > record_type ;
1078 
1079  // Disallow padding
1080  typedef std::integral_constant< unsigned , 0 > padding ;
1081 
1082  // Check if ViewCtorProp has CommonViewAllocProp - if so, retrieve the fad_size and append to layout
1083  enum { test_traits_check = Kokkos::Impl::check_has_common_view_alloc_prop< P... >::value };
1084 
1085  typename Traits::array_layout internal_layout =
1086  (test_traits_check == true)
1087  ? Kokkos::Impl::appendFadToLayoutViewAllocHelper< Traits, P... >::returnNewLayoutPlusFad(prop, local_layout)
1088  : local_layout;
1089 
1090  m_impl_offset = offset_type( padding(), internal_layout );
1091 
1092  m_array_offset =
1093  array_offset_type( padding() ,
1094  create_fad_array_layout<unsigned(Rank), unsigned(FadStaticDimension)>( internal_layout ) );
1095  const unsigned fad_dim =
1096  getFadDimension<unsigned(Rank)>( m_array_offset );
1097  if (unsigned(FadStaticDimension) == 0 && fad_dim == 0)
1098  Kokkos::abort("invalid fad dimension (0) supplied!");
1099  m_fad_size = fad_dim - 1 ;
1100  m_original_fad_size = m_fad_size.value ;
1101  m_fad_stride = 1;
1102  m_fad_index = 0;
1103 
1104  const size_t alloc_size = m_array_offset.span() * sizeof(fad_value_type);
1105 
1106  // Create shared memory tracking record with allocate memory from the memory space
1107  record_type * const record =
1108  record_type::allocate( ( (ViewCtorProp<void,memory_space> const &) prop ).value
1109  , ( (ViewCtorProp<void,std::string> const &) prop ).value
1110  , alloc_size );
1111 
1112  // Only set the the pointer and initialize if the allocation is non-zero.
1113  // May be zero if one of the dimensions is zero.
1114  if ( alloc_size ) {
1115 
1116  m_impl_handle = handle_type( reinterpret_cast< pointer_type >( record->data() ) );
1117 
1118  if ( ctor_prop::initialize ) {
1119  // Assume destruction is only required when construction is requested.
1120  // The ViewValueFunctor has both value construction and destruction operators.
1121  record->m_destroy = functor_type( ( (ViewCtorProp<void,execution_space> const &) prop).value
1122  , (fad_value_type *) m_impl_handle
1123  , m_array_offset.span()
1124  );
1125 
1126  // Construct values
1127  record->m_destroy.construct_shared_allocation();
1128  }
1129  }
1130 
1131  return record ;
1132  }
1133 
1134 };
1135 
1136 } // namespace Impl
1137 } // namespace Kokkos
1138 
1139 //----------------------------------------------------------------------------
1140 
1141 namespace Kokkos {
1142 namespace Impl {
1143 
1148 template< class DstTraits , class SrcTraits >
1149 class ViewMapping< DstTraits , SrcTraits ,
1150  typename std::enable_if<(
1151  Kokkos::Impl::MemorySpaceAccess
1152  < typename DstTraits::memory_space
1153  , typename SrcTraits::memory_space >::assignable
1154  &&
1155  // Destination view has FAD
1156  std::is_same< typename DstTraits::specialize
1157  , ViewSpecializeSacadoFadContiguous >::value
1158  &&
1159  // Source view has FAD
1160  std::is_same< typename SrcTraits::specialize
1161  , ViewSpecializeSacadoFadContiguous >::value
1162  )
1163  , typename DstTraits::specialize
1164  >::type >
1165 {
1166 public:
1167 
1168  enum { is_assignable = true };
1169  enum { is_assignable_data_type = true };
1170 
1171  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1172  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1173  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1174 
1175  template< class DstType >
1176  KOKKOS_INLINE_FUNCTION static
1177  void assign( DstType & dst
1178  , const SrcFadType & src
1179  , const TrackType & )
1180  {
1181  static_assert(
1182  (
1183  std::is_same< typename DstTraits::array_layout
1184  , Kokkos::LayoutLeft >::value ||
1185  std::is_same< typename DstTraits::array_layout
1186  , Kokkos::LayoutRight >::value ||
1187  std::is_same< typename DstTraits::array_layout
1188  , Kokkos::LayoutStride >::value
1189  )
1190  &&
1191  (
1192  std::is_same< typename SrcTraits::array_layout
1193  , Kokkos::LayoutLeft >::value ||
1194  std::is_same< typename SrcTraits::array_layout
1195  , Kokkos::LayoutRight >::value ||
1196  std::is_same< typename SrcTraits::array_layout
1197  , Kokkos::LayoutStride >::value
1198  )
1199  , "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1200 
1201  static_assert(
1202  std::is_same< typename DstTraits::array_layout
1203  , typename SrcTraits::array_layout >::value ||
1204  std::is_same< typename DstTraits::array_layout
1205  , Kokkos::LayoutStride >::value ,
1206  "View assignment must have compatible layout" );
1207 
1208  static_assert(
1209  std::is_same< typename DstTraits::scalar_array_type
1210  , typename SrcTraits::scalar_array_type >::value ||
1211  std::is_same< typename DstTraits::scalar_array_type
1212  , typename SrcTraits::const_scalar_array_type >::value ,
1213  "View assignment must have same value type or const = non-const" );
1214 
1215  static_assert(
1216  ViewDimensionAssignable
1217  < typename DstType::offset_type::dimension_type
1218  , typename SrcFadType::offset_type::dimension_type >::value ,
1219  "View assignment must have compatible dimensions" );
1220 
1221  static_assert(
1222  ViewDimensionAssignable
1223  < typename DstType::array_offset_type::dimension_type
1224  , typename SrcFadType::array_offset_type::dimension_type >::value ,
1225  "View assignment must have compatible dimensions" );
1226 
1227  typedef typename DstType::offset_type dst_offset_type ;
1228  typedef typename DstType::array_offset_type dst_array_offset_type ;
1229 
1230  dst.m_impl_handle = src.m_impl_handle ;
1231  dst.m_impl_offset = dst_offset_type( src.m_impl_offset );
1232  dst.m_array_offset = dst_array_offset_type( src.m_array_offset );
1233  dst.m_fad_size = src.m_fad_size.value ;
1234  dst.m_original_fad_size = src.m_original_fad_size ;
1235  dst.m_fad_stride = src.m_fad_stride ;
1236  dst.m_fad_index = src.m_fad_index ;
1237  }
1238 };
1239 
1244 template< class DstTraits , class SrcTraits >
1245 class ViewMapping< DstTraits , SrcTraits ,
1246  typename std::enable_if<(
1247  std::is_same< typename DstTraits::memory_space
1248  , typename SrcTraits::memory_space >::value
1249  &&
1250  // Destination view has FAD
1251  std::is_same< typename DstTraits::specialize
1252  , ViewSpecializeSacadoFad >::value
1253  &&
1254  // Source view has FAD contiguous
1255  std::is_same< typename SrcTraits::specialize
1256  , ViewSpecializeSacadoFadContiguous >::value
1257  &&
1258  // Destination view is LayoutStride
1259  std::is_same< typename DstTraits::array_layout
1260  , Kokkos::LayoutStride >::value
1261  )
1262  , typename DstTraits::specialize
1263  >::type >
1264 {
1265 public:
1266 
1267  enum { is_assignable = true };
1268  enum { is_assignable_data_type = true };
1269 
1270  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1271  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1272  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1273 
1274  template< class DstType >
1275  KOKKOS_INLINE_FUNCTION static
1276  void assign( DstType & dst
1277  , const SrcFadType & src
1278  , const TrackType & )
1279  {
1280  static_assert(
1281  std::is_same< typename SrcTraits::array_layout
1282  , Kokkos::LayoutLeft >::value ||
1283  std::is_same< typename SrcTraits::array_layout
1284  , Kokkos::LayoutRight >::value ||
1285  std::is_same< typename SrcTraits::array_layout
1286  , Kokkos::LayoutStride >::value ,
1287  "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1288 
1289  static_assert(
1290  std::is_same< typename DstTraits::value_type
1291  , typename SrcTraits::value_type >::value ||
1292  std::is_same< typename DstTraits::value_type
1293  , typename SrcTraits::const_value_type >::value ,
1294  "View assignment must have same value type or const = non-const" );
1295 
1296  static_assert(
1297  DstTraits::dimension::rank == SrcTraits::dimension::rank,
1298  "View assignment must have same rank" );
1299 
1300  typedef typename DstType::array_offset_type dst_offset_type ;
1301 
1302  dst.m_impl_handle = src.m_impl_handle ;
1303  dst.m_fad_size = src.m_fad_size.value ;
1304  dst.m_fad_stride = src.m_fad_stride ;
1305  dst.m_impl_offset = src.m_impl_offset;
1306 
1307  size_t N[8], S[8];
1308  N[0] = src.m_array_offset.dimension_0();
1309  N[1] = src.m_array_offset.dimension_1();
1310  N[2] = src.m_array_offset.dimension_2();
1311  N[3] = src.m_array_offset.dimension_3();
1312  N[4] = src.m_array_offset.dimension_4();
1313  N[5] = src.m_array_offset.dimension_5();
1314  N[6] = src.m_array_offset.dimension_6();
1315  N[7] = src.m_array_offset.dimension_7();
1316  S[0] = src.m_array_offset.stride_0();
1317  S[1] = src.m_array_offset.stride_1();
1318  S[2] = src.m_array_offset.stride_2();
1319  S[3] = src.m_array_offset.stride_3();
1320  S[4] = src.m_array_offset.stride_4();
1321  S[5] = src.m_array_offset.stride_5();
1322  S[6] = src.m_array_offset.stride_6();
1323  S[7] = src.m_array_offset.stride_7();
1324 
1325  // For LayoutLeft, we have to move the Sacado dimension from the first
1326  // to the last
1327  if (std::is_same< typename SrcTraits::array_layout
1328  , Kokkos::LayoutLeft >::value)
1329  {
1330  const size_t N_fad = N[0];
1331  const size_t S_fad = S[0];
1332  for (int i=0; i<7; ++i) {
1333  N[i] = N[i+1];
1334  S[i] = S[i+1];
1335  }
1336  N[DstTraits::dimension::rank] = N_fad;
1337  S[DstTraits::dimension::rank] = S_fad;
1338  }
1339  Kokkos::LayoutStride ls( N[0], S[0],
1340  N[1], S[1],
1341  N[2], S[2],
1342  N[3], S[3],
1343  N[4], S[4],
1344  N[5], S[5],
1345  N[6], S[6],
1346  N[7], S[7] );
1347  dst.m_array_offset = dst_offset_type(std::integral_constant<unsigned,0>(), ls);
1348  }
1349 };
1350 
1355 template< class DstTraits , class SrcTraits >
1356 class ViewMapping< DstTraits , SrcTraits ,
1357  typename std::enable_if<(
1358  std::is_same< typename DstTraits::memory_space
1359  , typename SrcTraits::memory_space >::value
1360  &&
1361  // Destination view has ordinary
1362  std::is_same< typename DstTraits::specialize , void >::value
1363  &&
1364  // Source view has FAD only
1365  std::is_same< typename SrcTraits::specialize
1366  , ViewSpecializeSacadoFadContiguous >::value
1367  )
1368  , typename DstTraits::specialize
1369  >::type >
1370 {
1371 public:
1372 
1373  enum { is_assignable = true };
1374  enum { is_assignable_data_type = true };
1375 
1376  typedef Kokkos::Impl::SharedAllocationTracker TrackType ;
1377  typedef ViewMapping< DstTraits , typename DstTraits::specialize > DstType ;
1378  typedef ViewMapping< SrcTraits , typename SrcTraits::specialize > SrcFadType ;
1379 
1380 
1381  // Helpers to assign, and generate if necessary, ViewOffset to the dst map
1382  // These are necessary to use Kokkos' deep_copy with nested fads
1383  template < class DstType, class SrcFadType, class Enable = void >
1384  struct AssignOffset;
1385 
1386  template < class DstType, class SrcFadType >
1387  struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank != (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1388  {
1389  // ViewOffset's Dimensions Ranks do not match
1391  static void assign( DstType & dst, const SrcFadType & src )
1392  {
1393  typedef typename SrcTraits::value_type TraitsValueType;
1394 
1397  )
1398  {
1399 
1400  typedef typename DstType::offset_type::array_layout DstLayoutType;
1401  //typedef typename ViewArrayLayoutSelector<typename DstType::offset_type::array_layout>::type DstLayoutType;
1402  typedef typename SrcFadType::array_offset_type::dimension_type SrcViewDimension;
1403 
1404  // This is the static dimension of the inner fad, missing from ViewDimension
1405  const size_t InnerStaticDim = Sacado::StaticSize< typename Sacado::ValueType< TraitsValueType >::type >::value;
1406 
1407  static constexpr bool is_layout_left =
1408  std::is_same< DstLayoutType, Kokkos::LayoutLeft>::value;
1409 
1410  typedef typename std::conditional< is_layout_left,
1411  typename SrcViewDimension:: template prepend< InnerStaticDim+1 >::type,
1412  typename SrcViewDimension:: template append < InnerStaticDim+1 >::type
1413  >::type SrcViewDimensionAppended;
1414 
1415  typedef std::integral_constant< unsigned , 0 > padding ;
1416 
1417  typedef ViewOffset< SrcViewDimensionAppended, DstLayoutType > TmpOffsetType;
1418 
1419  auto src_layout = src.m_array_offset.layout();
1420 
1421  if ( is_layout_left ) {
1422  auto prepend_layout = Kokkos::Impl::prependFadToLayout< DstLayoutType >::returnNewLayoutPlusFad(src_layout, InnerStaticDim+1);
1423  TmpOffsetType offset_tmp( padding(), prepend_layout );
1424  dst.m_impl_offset = offset_tmp;
1425  }
1426  else {
1427  TmpOffsetType offset_tmp( padding(), src_layout );
1428  dst.m_impl_offset = offset_tmp;
1429  }
1430  } else {
1431  Kokkos::abort("Sacado error: Applying AssignOffset for case with nested Fads, but without nested Fads - something went wrong");
1432  }
1433  }
1434  };
1435 
1436  template < class DstType, class SrcFadType >
1437  struct AssignOffset< DstType, SrcFadType, typename std::enable_if< ((int)DstType::offset_type::dimension_type::rank == (int)SrcFadType::array_offset_type::dimension_type::rank) >::type >
1438  {
1440  static void assign( DstType & dst, const SrcFadType & src )
1441  {
1442  typedef typename DstType::offset_type dst_offset_type ;
1443  dst.m_impl_offset = dst_offset_type( src.m_array_offset );
1444  }
1445  };
1446 
1447  template< class DstType >
1448  KOKKOS_INLINE_FUNCTION static
1449  void assign( DstType & dst
1450  , const SrcFadType & src
1451  , const TrackType &
1452  )
1453  {
1454  static_assert(
1455  (
1456  std::is_same< typename DstTraits::array_layout
1457  , Kokkos::LayoutLeft >::value ||
1458  std::is_same< typename DstTraits::array_layout
1459  , Kokkos::LayoutRight >::value ||
1460  std::is_same< typename DstTraits::array_layout
1461  , Kokkos::LayoutStride >::value
1462  )
1463  &&
1464  (
1465  std::is_same< typename SrcTraits::array_layout
1466  , Kokkos::LayoutLeft >::value ||
1467  std::is_same< typename SrcTraits::array_layout
1468  , Kokkos::LayoutRight >::value ||
1469  std::is_same< typename SrcTraits::array_layout
1470  , Kokkos::LayoutStride >::value
1471  )
1472  , "View of FAD requires LayoutLeft, LayoutRight, or LayoutStride" );
1473 
1474  static_assert(
1475  std::is_same< typename DstTraits::array_layout
1476  , typename SrcTraits::array_layout >::value ||
1477  std::is_same< typename DstTraits::array_layout
1478  , Kokkos::LayoutStride >::value ,
1479  "View assignment must have compatible layout" );
1480 
1481  if ( src.m_fad_index != 0 || src.m_fad_stride != 1 ) {
1482  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Cannot assign to array with partitioned view ******\n\n");
1483  }
1484 
1485  AssignOffset< DstType, SrcFadType >::assign( dst, src );
1486  dst.m_impl_handle = reinterpret_cast< typename DstType::handle_type >(src.m_impl_handle) ;
1487  }
1488 };
1489 
1490 } // namespace Impl
1491 } // namespace Kokkos
1492 
1493 //----------------------------------------------------------------------------
1494 
1495 namespace Kokkos {
1496 namespace Impl {
1497 
1498 // Rules for subview arguments and layouts matching
1499 
1500 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1501 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1502  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1503 };
1504 
1505 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1506 struct SubviewLegalArgsCompileTime<LayoutDest,Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1507  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1508 };
1509 
1510 template<class LayoutDest, class LayoutSrc, int RankDest, int RankSrc, int CurrentArg, class ... SubViewArgs>
1511 struct SubviewLegalArgsCompileTime<Kokkos::LayoutContiguous<LayoutDest>,Kokkos::LayoutContiguous<LayoutSrc>,RankDest,RankSrc,CurrentArg,SubViewArgs...> {
1512  enum { value = SubviewLegalArgsCompileTime<LayoutDest,LayoutSrc,RankDest,RankSrc,CurrentArg,SubViewArgs...>::value };
1513 };
1514 
1515 // Subview mapping
1516 
1517 template< class SrcTraits , class Arg0 , class ... Args >
1518 struct ViewMapping
1519  < typename std::enable_if<(
1520  // Source view has FAD only
1521  std::is_same< typename SrcTraits::specialize
1522  , ViewSpecializeSacadoFadContiguous >::value
1523  &&
1524  (
1525  std::is_same< typename SrcTraits::array_layout
1526  , Kokkos::LayoutLeft >::value ||
1527  std::is_same< typename SrcTraits::array_layout
1528  , Kokkos::LayoutRight >::value ||
1529  std::is_same< typename SrcTraits::array_layout
1530  , Kokkos::LayoutStride >::value
1531  )
1532  && !Sacado::Fad::is_fad_partition<Arg0>::value
1533  )
1534  >::type
1535  , SrcTraits
1536  , Arg0, Args ... >
1537 {
1538 private:
1539 
1540  static_assert( SrcTraits::rank == sizeof...(Args)+1 , "" );
1541 
1542  enum
1543  { RZ = false
1544  , R0 = bool(is_integral_extent<0,Arg0,Args...>::value)
1545  , R1 = bool(is_integral_extent<1,Arg0,Args...>::value)
1546  , R2 = bool(is_integral_extent<2,Arg0,Args...>::value)
1547  , R3 = bool(is_integral_extent<3,Arg0,Args...>::value)
1548  , R4 = bool(is_integral_extent<4,Arg0,Args...>::value)
1549  , R5 = bool(is_integral_extent<5,Arg0,Args...>::value)
1550  , R6 = bool(is_integral_extent<6,Arg0,Args...>::value)
1551  };
1552 
1553  // Public rank
1554  enum { rank = unsigned(R0) + unsigned(R1) + unsigned(R2) + unsigned(R3)
1555  + unsigned(R4) + unsigned(R5) + unsigned(R6) };
1556 
1557  // Whether right-most non-FAD rank is a range.
1558  enum { R0_rev = ( 0 == SrcTraits::rank ? RZ : (
1559  1 == SrcTraits::rank ? R0 : (
1560  2 == SrcTraits::rank ? R1 : (
1561  3 == SrcTraits::rank ? R2 : (
1562  4 == SrcTraits::rank ? R3 : (
1563  5 == SrcTraits::rank ? R4 : (
1564  6 == SrcTraits::rank ? R5 : R6 ))))))) };
1565 
1566  // Subview's layout
1567  typedef typename std::conditional<
1568  ( /* Same array layout IF */
1569  ( rank == 0 ) /* output rank zero */
1570  ||
1571  // OutputRank 1 or 2, InputLayout Left, Interval 0
1572  // because single stride one or second index has a stride.
1573  ( rank <= 2 && R0 && std::is_same< typename SrcTraits::array_layout , Kokkos::LayoutLeft >::value )
1574  ||
1575  // OutputRank 1 or 2, InputLayout Right, Interval [InputRank-1]
1576  // because single stride one or second index has a stride.
1577  ( rank <= 2 && R0_rev && std::is_same< typename SrcTraits::array_layout , Kokkos::LayoutRight >::value )
1579  >::type array_layout ;
1580 
1581  typedef typename SrcTraits::value_type fad_type ;
1582 
1583  typedef typename std::conditional< rank == 0 , fad_type ,
1584  typename std::conditional< rank == 1 , fad_type * ,
1585  typename std::conditional< rank == 2 , fad_type ** ,
1586  typename std::conditional< rank == 3 , fad_type *** ,
1587  typename std::conditional< rank == 4 , fad_type **** ,
1588  typename std::conditional< rank == 5 , fad_type ***** ,
1589  typename std::conditional< rank == 6 , fad_type ****** ,
1590  fad_type *******
1591  >::type >::type >::type >::type >::type >::type >::type
1592  data_type ;
1593 
1594 public:
1595 
1596  typedef Kokkos::ViewTraits
1597  < data_type
1598  , array_layout
1599  , typename SrcTraits::device_type
1600  , typename SrcTraits::memory_traits > traits_type ;
1601 
1602  typedef Kokkos::View
1603  < data_type
1604  , array_layout
1605  , typename SrcTraits::device_type
1606  , typename SrcTraits::memory_traits > type ;
1607 
1608 
1610  static void assign( ViewMapping< traits_type , typename traits_type::specialize > & dst
1611  , ViewMapping< SrcTraits , typename SrcTraits::specialize > const & src
1612  , Arg0 arg0 , Args ... args )
1613  {
1614  typedef ViewMapping< traits_type , typename traits_type::specialize > DstType ;
1615  typedef typename DstType::offset_type dst_offset_type ;
1616  typedef typename DstType::array_offset_type dst_array_offset_type ;
1617  typedef typename DstType::handle_type dst_handle_type ;
1618 
1619  size_t offset;
1620  if (std::is_same< typename SrcTraits::array_layout, LayoutLeft >::value) {
1621  const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1622  array_extents( src.m_array_offset.m_dim , Kokkos::ALL() , arg0 , args... );
1623  offset = src.m_array_offset( array_extents.domain_offset(0)
1624  , array_extents.domain_offset(1)
1625  , array_extents.domain_offset(2)
1626  , array_extents.domain_offset(3)
1627  , array_extents.domain_offset(4)
1628  , array_extents.domain_offset(5)
1629  , array_extents.domain_offset(6)
1630  , array_extents.domain_offset(7) );
1631  dst.m_array_offset = dst_array_offset_type( src.m_array_offset ,
1632  array_extents );
1633  }
1634  else {
1635  const SubviewExtents< SrcTraits::rank + 1 , rank + 1 >
1636  array_extents( src.m_array_offset.m_dim , arg0 , args... , Kokkos::ALL() );
1637  offset = src.m_array_offset( array_extents.domain_offset(0)
1638  , array_extents.domain_offset(1)
1639  , array_extents.domain_offset(2)
1640  , array_extents.domain_offset(3)
1641  , array_extents.domain_offset(4)
1642  , array_extents.domain_offset(5)
1643  , array_extents.domain_offset(6)
1644  , array_extents.domain_offset(7) );
1645  dst.m_array_offset = dst_array_offset_type( src.m_array_offset ,
1646  array_extents );
1647  }
1648 
1649  const SubviewExtents< SrcTraits::rank , rank >
1650  extents( src.m_impl_offset.m_dim , arg0 , args... );
1651 
1652  dst.m_impl_offset = dst_offset_type( src.m_impl_offset , extents );
1653  dst.m_impl_handle = dst_handle_type( src.m_impl_handle + offset );
1654  dst.m_fad_size = src.m_fad_size;
1655  dst.m_original_fad_size = src.m_original_fad_size;
1656  dst.m_fad_stride = src.m_fad_stride;
1657  dst.m_fad_index = src.m_fad_index;
1658  }
1659 
1660 };
1661 
1662 } // namespace Impl
1663 } // namespace Kokkos
1664 
1665 //---------------------------------------------------------------------------
1666 
1667 namespace Kokkos {
1668 namespace Impl {
1669 
1670 // Partition mapping
1671 
1672 template< class DataType, class ...P, unsigned Stride >
1673 class ViewMapping<
1674  void,
1675  ViewTraits<DataType,P...> ,
1676  Sacado::Fad::Partition<Stride>
1677  >
1678 {
1679 public:
1680 
1681  enum { is_assignable = true };
1682  enum { is_assignable_data_type = true };
1683 
1684  typedef ViewTraits<DataType,P...> src_traits;
1685  typedef ViewMapping< src_traits , typename src_traits::specialize > src_type ;
1686 
1687  typedef typename src_type::offset_type::dimension_type src_dimension;
1688  typedef typename src_traits::value_type fad_type;
1689  typedef typename Sacado::LocalScalarType<fad_type,Stride>::type strided_fad_type;
1690  typedef typename
1691  ViewDataType< strided_fad_type , src_dimension >::type strided_data_type;
1692  typedef ViewTraits<strided_data_type,P...> dst_traits;
1693  typedef View<strided_data_type,P...> type;
1694  typedef ViewMapping< dst_traits , typename dst_traits::specialize > dst_type ;
1695 
1696  KOKKOS_INLINE_FUNCTION static
1697  void assign( dst_type & dst
1698  , const src_type & src
1699  , const Sacado::Fad::Partition<Stride> & part )
1700  {
1701  if ( Stride != part.stride && Stride != 0 ) {
1702  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Invalid size in partitioned view assignment ******\n\n");
1703  }
1704  if ( src.m_fad_stride != 1 ) {
1705  Kokkos::abort("\n\n ****** Kokkos::View< Sacado::Fad ... > Can't partition already partitioned view ******\n\n");
1706  }
1707 
1708  dst.m_impl_handle = src.m_impl_handle ;
1709  dst.m_impl_offset = src.m_impl_offset ;
1710  dst.m_array_offset = src.m_array_offset ;
1711 
1712  // Assuming the alignment was choosen correctly for the partitioning,
1713  // each partition should get the same size. This allows the use of SFad.
1714  dst.m_fad_size =
1715  (src.m_fad_size.value + part.stride-part.offset-1) / part.stride ;
1716 
1717  dst.m_original_fad_size = src.m_original_fad_size ;
1718  dst.m_fad_stride = part.stride ;
1719  dst.m_fad_index = part.offset ;
1720  }
1721 };
1722 
1723 } // namespace Impl
1724 } // namespace Kokkos
1725 
1726 #endif // defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
1727 
1728 #endif // defined(HAVE_SACADO_KOKKOSCORE)
1729 
1730 #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 KOKKOS_DEFAULTED_FUNCTION
#define D
Definition: Sacado_rad.hpp:577
GeneralFad< DynamicStorage< T > > DFad
const int N
void
Definition: uninit.c:96
const int fad_dim
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.