Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Fad_KokkosTests.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
30 
31 #include "Sacado.hpp"
33 
34 template <typename T>
35 struct is_sfad {
36  static const bool value = false;
37 };
38 
39 template <typename T, int N>
40 struct is_sfad< Sacado::Fad::SFad<T,N> > {
41  static const bool value = true;
42 };
43 
44 template <typename T>
45 struct is_dfad {
46  static const bool value = false;
47 };
48 
49 template <typename T>
50 struct is_dfad< Sacado::Fad::DFad<T> > {
51  static const bool value = true;
52 };
53 
54 template <typename FadType1, typename FadType2>
55 bool checkFads(const FadType1& x, const FadType2& x2,
56  Teuchos::FancyOStream& out, double tol = 1.0e-15)
57 {
58  bool success = true;
59 
60  // Check sizes match
61  TEUCHOS_TEST_EQUALITY(x.size(), x2.size(), out, success);
62 
63  // Check values match
64  TEUCHOS_TEST_FLOATING_EQUALITY(x.val(), x2.val(), tol, out, success);
65 
66  // Check derivatives match
67  for (int i=0; i<x.size(); ++i)
68  TEUCHOS_TEST_FLOATING_EQUALITY(x.dx(i), x2.dx(i), tol, out, success);
69 
70  return success;
71 }
72 
73 template <typename fadtype, typename ordinal>
74 inline
75 fadtype generate_fad( const ordinal num_rows,
76  const ordinal num_cols,
77  const ordinal fad_size,
78  const ordinal row,
79  const ordinal col )
80 {
81  typedef typename fadtype::value_type scalar;
82  fadtype x(fad_size, scalar(0.0));
83 
84  const scalar x_row = 100.0 + scalar(num_rows) / scalar(row+1);
85  const scalar x_col = 10.0 + scalar(num_cols) / scalar(col+1);
86  x.val() = x_row + x_col;
87  for (ordinal i=0; i<fad_size; ++i) {
88  const scalar x_fad = 1.0 + scalar(fad_size) / scalar(i+1);
89  x.fastAccessDx(i) = x_row + x_col + x_fad;
90  }
91  return x;
92 }
93 
94 #ifndef GLOBAL_FAD_SIZE
95 #define GLOBAL_FAD_SIZE 5
96 #endif
97 const int global_num_rows = 11;
98 const int global_num_cols = 7;
100 
101 // Kernel to multiply two views
102 template <typename InputViewType1,
103  typename InputViewType2 = InputViewType1,
104  typename OutputViewType = InputViewType1>
106  typedef typename InputViewType1::execution_space execution_space;
107  typedef typename InputViewType1::size_type size_type;
108  typedef Kokkos::RangePolicy< execution_space> range_policy_type;
109  typedef Kokkos::TeamPolicy< execution_space> team_policy_type;
110  typedef typename team_policy_type::member_type team_handle;
111 
112  const InputViewType1 m_v1;
113  const InputViewType2 m_v2;
114  const OutputViewType m_v3;
115  const bool m_update;
116 
117  MultiplyKernel(const InputViewType1 v1,
118  const InputViewType2 v2,
119  const OutputViewType v3,
120  const bool update) :
121  m_v1(v1), m_v2(v2), m_v3(v3), m_update(update) {};
122 
123  // Multiply entries for row 'i' with a value
124  KOKKOS_INLINE_FUNCTION
125  void operator() (const size_type i) const {
126  if (m_update)
127  m_v3(i) += m_v1(i)*m_v2(i);
128  else
129  m_v3(i) = m_v1(i)*m_v2(i);
130  }
131 
132  KOKKOS_INLINE_FUNCTION
133  void operator()( const team_handle& team ) const
134  {
135  const size_type i = team.league_rank()*team.team_size() + team.team_rank();
136  if (i < m_v1.extent(0))
137  (*this)(i);
138  }
139 
140  // Kernel launch
141  static void apply(const InputViewType1 v1,
142  const InputViewType2 v2,
143  const OutputViewType v3,
144  const bool update = false) {
145  const size_type nrow = v1.extent(0);
146 
147 #if defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
148  const size_type stride = Kokkos::ViewScalarStride<InputViewType1>::stride;
149  const bool use_team =
153  ( stride > 1 );
154 #elif defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
155  const size_type stride = team_policy_type::vector_length_max(); // 32
156  const bool use_team =
161 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
162  const size_type stride = Kokkos::ViewScalarStride<InputViewType1>::stride;
163  const bool use_team =
167  ( stride > 1 );
168 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
169  const size_type stride = team_policy_type::vector_length_max(); // 64
170  const bool use_team =
175 #else
176  const size_type stride = 1;
177  const bool use_team = false;
178 #endif
179 
180  if (use_team) {
181  const size_type team_size = 256 / stride;
182  team_policy_type policy( (nrow+team_size-1)/team_size, team_size, stride );
183  Kokkos::parallel_for( policy, MultiplyKernel(v1,v2,v3,update) );
184  }
185  else {
186  range_policy_type policy( 0, nrow );
187  Kokkos::parallel_for( policy, MultiplyKernel(v1,v2,v3,update) );
188  }
189  }
190 };
191 
192 // Kernel to assign a constant to a view
193 template <typename ViewType>
195  typedef typename ViewType::execution_space execution_space;
196  typedef typename ViewType::size_type size_type;
197  typedef typename ViewType::value_type::value_type ScalarType;
198  typedef Kokkos::TeamPolicy< execution_space> team_policy_type;
199  typedef Kokkos::RangePolicy< execution_space> range_policy_type;
200  typedef typename team_policy_type::member_type team_handle;
201  static const size_type stride = Kokkos::ViewScalarStride<ViewType>::stride;
202 
203  const ViewType m_v;
205 
206  ScalarAssignKernel(const ViewType& v, const ScalarType& s) :
207  m_v(v), m_s(s) {};
208 
209  // Multiply entries for row 'i' with a value
210  KOKKOS_INLINE_FUNCTION
211  void operator() (const size_type i) const {
212  m_v(i) = m_s;
213  }
214 
215  KOKKOS_INLINE_FUNCTION
216  void operator()( const team_handle& team ) const
217  {
218  const size_type i = team.league_rank()*team.team_size() + team.team_rank();
219  if (i < m_v.extent(0))
220  (*this)(i);
221  }
222 
223  // Kernel launch
224  static void apply(const ViewType& v, const ScalarType& s) {
225  const size_type nrow = v.extent(0);
226 
227 #if defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
228  const bool use_team =
232  ( stride > 1 );
233 #elif defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
234  const bool use_team =
239 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
240  const bool use_team =
244  ( stride > 1 );
245 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
246  const bool use_team =
251 #else
252  const bool use_team = false;
253 #endif
254 
255  if (use_team) {
256  const size_type team_size = 256 / stride;
257  team_policy_type policy( (nrow+team_size-1)/team_size, team_size, stride );
258  Kokkos::parallel_for( policy, ScalarAssignKernel(v,s) );
259  }
260  else {
261  range_policy_type policy( 0, nrow );
262  Kokkos::parallel_for( policy, ScalarAssignKernel(v,s) );
263  }
264  }
265 };
266 
267 // Kernel to assign a constant to a view
268 template <typename ViewType, typename ScalarViewType>
270  typedef typename ViewType::execution_space execution_space;
271  typedef typename ViewType::size_type size_type;
272  typedef typename ViewType::value_type ValueType;
273  typedef Kokkos::TeamPolicy< execution_space> team_policy_type;
274  typedef Kokkos::RangePolicy< execution_space> range_policy_type;
275  typedef typename team_policy_type::member_type team_handle;
276  typedef typename Kokkos::ThreadLocalScalarType<ViewType>::type local_scalar_type;
277  static const size_type stride = Kokkos::ViewScalarStride<ViewType>::stride;
278 
279  const ViewType m_v;
280  const ScalarViewType m_s;
281 
282  ValueAssignKernel(const ViewType& v, const ScalarViewType& s) :
283  m_v(v), m_s(s) {};
284 
285  // Multiply entries for row 'i' with a value
286  KOKKOS_INLINE_FUNCTION
287  void operator() (const size_type i) const {
288  local_scalar_type s = Sacado::partition_scalar<stride>(m_s());
289  m_v(i) = s;
290  }
291 
292  KOKKOS_INLINE_FUNCTION
293  void operator()( const team_handle& team ) const
294  {
295  const size_type i = team.league_rank()*team.team_size() + team.team_rank();
296  if (i < m_v.extent(0))
297  (*this)(i);
298  }
299 
300  // Kernel launch
301  static void apply(const ViewType& v, const ScalarViewType& s) {
302  const size_type nrow = v.extent(0);
303 
304 #if defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
305  const bool use_team =
309  ( stride > 1 );
310 #elif defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
311  const bool use_team =
316 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
317  const bool use_team =
321  ( stride > 1 );
322 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
323  const bool use_team =
328 #else
329  const bool use_team = false;
330 #endif
331 
332  if (use_team) {
333  const size_type team_size = 256 / stride;
334  team_policy_type policy( (nrow+team_size-1)/team_size, team_size, stride );
335  Kokkos::parallel_for( policy, ValueAssignKernel(v,s) );
336  }
337  else {
338  range_policy_type policy( 0, nrow );
339  Kokkos::parallel_for( policy, ValueAssignKernel(v,s) );
340  }
341  }
342 };
343 
344 // Kernel to assign a column of a rank-2 to a rank-1
345 template <typename InputViewType,
346  typename OutputViewType,
347  typename Enabled = void>
349  typedef typename InputViewType::execution_space execution_space;
350  typedef typename InputViewType::size_type size_type;
351  typedef Kokkos::TeamPolicy< execution_space> team_policy_type;
352  typedef Kokkos::RangePolicy< execution_space> range_policy_type;
353  typedef typename team_policy_type::member_type team_handle;
354  static const size_type stride = Kokkos::ViewScalarStride<InputViewType>::stride;
355 
356  const InputViewType m_v1;
357  const OutputViewType m_v2;
359 
360  AssignRank2Rank1Kernel(const InputViewType v1,
361  const OutputViewType v2,
362  const size_type col) :
363  m_v1(v1), m_v2(v2), m_col(col) {
364  static_assert( unsigned(InputViewType::rank) == 2 ,
365  "Require rank-2 input view" );
366  static_assert( unsigned(OutputViewType::rank) == 1 ,
367  "Require rank-1 output view" );
368  };
369 
370  // Multiply entries for row 'i' with a value
371  KOKKOS_INLINE_FUNCTION
372  void operator() (const size_type i) const {
373  m_v2(i) = m_v1(i,m_col);
374  }
375 
376  KOKKOS_INLINE_FUNCTION
377  void operator()( const team_handle& team ) const
378  {
379  const size_type i = team.league_rank()*team.team_size() + team.team_rank();
380  if (i < m_v1.extent(0))
381  (*this)(i);
382  }
383 
384  // Kernel launch
385  static void apply(const InputViewType v1,
386  const OutputViewType v2,
387  const size_type col) {
388  const size_type nrow = v1.extent(0);
389 
390 #if defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
391  const bool use_team =
395  ( stride > 1 );
396 #elif defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
397  const bool use_team =
402 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
403  const bool use_team =
407  ( stride > 1 );
408 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
409  const bool use_team =
414 #else
415  const bool use_team = false;
416 #endif
417 
418  if (use_team) {
419  const size_type team_size = 256 / stride;
420  team_policy_type policy( (nrow+team_size-1)/team_size, team_size, stride );
421  Kokkos::parallel_for( policy, AssignRank2Rank1Kernel(v1,v2,col) );
422  }
423  else {
424  range_policy_type policy( 0, nrow );
425  Kokkos::parallel_for( policy, AssignRank2Rank1Kernel(v1,v2,col) );
426  }
427  }
428 };
429 
430 // Kernel to test atomic_add
431 template <typename ViewType, typename ScalarViewType>
433  typedef typename ViewType::execution_space execution_space;
434  typedef typename ViewType::size_type size_type;
435  typedef Kokkos::TeamPolicy< execution_space> team_policy_type;
436  typedef Kokkos::RangePolicy< execution_space> range_policy_type;
437  typedef typename team_policy_type::member_type team_handle;
438  typedef typename Kokkos::ThreadLocalScalarType<ViewType>::type local_scalar_type;
439  static const size_type stride = Kokkos::ViewScalarStride<ViewType>::stride;
440 
441  const ViewType m_v;
442  const ScalarViewType m_s;
443 
444  AtomicAddKernel(const ViewType& v, const ScalarViewType& s) :
445  m_v(v), m_s(s) {};
446 
447  // Multiply entries for row 'i' with a value
448  KOKKOS_INLINE_FUNCTION
449  void operator() (const size_type i) const {
450  local_scalar_type x = m_v(i);
451  Kokkos::atomic_add(&(m_s()), x);
452  }
453 
454  KOKKOS_INLINE_FUNCTION
455  void operator()( const team_handle& team ) const
456  {
457  const size_type i = team.league_rank()*team.team_size() + team.team_rank();
458  if (i < m_v.extent(0))
459  (*this)(i);
460  }
461 
462  // Kernel launch
463  static void apply(const ViewType& v, const ScalarViewType& s) {
464  const size_type nrow = v.extent(0);
465 
466 #if defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
467  const bool use_team =
471  ( stride > 1 );
472 #elif defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
473  const bool use_team =
478 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
479  const bool use_team =
483  ( stride > 1 );
484 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
485  const bool use_team =
490 #else
491  const bool use_team = false;
492 #endif
493 
494  if (use_team) {
495  const size_type team_size = 256 / stride;
496  team_policy_type policy( (nrow+team_size-1)/team_size, team_size, stride );
497  Kokkos::parallel_for( policy, AtomicAddKernel(v,s) );
498  }
499  else {
500  range_policy_type policy( 0, nrow );
501  Kokkos::parallel_for( policy, AtomicAddKernel(v,s) );
502  }
503  }
504 };
505 
507  Kokkos_View_Fad, Size, FadType, Layout, Device )
508 {
509  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
510  typedef typename ViewType::size_type size_type;
511 
512  const size_type num_rows = global_num_rows;
513 
514  // Create and fill view
515  ViewType v;
516 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
517  v = ViewType("view", num_rows);
518 #else
519  const size_type fad_size = global_fad_size;
520  v = ViewType("view", num_rows, fad_size+1);
521 #endif
522  TEUCHOS_TEST_EQUALITY(v.size(), num_rows, out, success);
523 }
524 
526  Kokkos_View_Fad, DeepCopy, FadType, Layout, Device )
527 {
528  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
529  typedef typename ViewType::size_type size_type;
530  typedef typename ViewType::HostMirror host_view_type;
531 
532  const size_type num_rows = global_num_rows;
533  const size_type num_cols = global_num_cols;
534  const size_type fad_size = global_fad_size;
535 
536  // Create and fill view
537  ViewType v;
538 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
539  v = ViewType ("view", num_rows, num_cols);
540 #else
541  v = ViewType ("view", num_rows, num_cols, fad_size+1);
542 #endif
543  host_view_type h_v = Kokkos::create_mirror_view(v);
544  for (size_type i=0; i<num_rows; ++i)
545  for (size_type j=0; j<num_cols; ++j)
546  h_v(i,j) = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
547  Kokkos::deep_copy(v, h_v);
548 
549  // Copy back
550  host_view_type h_v2 = Kokkos::create_mirror_view(v);
551  Kokkos::deep_copy(h_v2, v);
552 
553  // Check
554  success = true;
555  for (size_type i=0; i<num_rows; ++i) {
556  for (size_type j=0; j<num_cols; ++j) {
557  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
558  success = success && checkFads(f, h_v2(i,j), out);
559  }
560  }
561 }
562 
564  Kokkos_View_Fad, DeepCopy_ConstantScalar, FadType, Layout, Device )
565 {
566  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
567  typedef typename ViewType::size_type size_type;
568  typedef typename ViewType::HostMirror host_view_type;
569  typedef typename FadType::value_type value_type;
570 
571  const size_type num_rows = global_num_rows;
572  const size_type num_cols = global_num_cols;
573 
574  // Create and fill view
575  ViewType v;
576 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
577  v = ViewType ("view", num_rows, num_cols);
578 #else
579  const size_type fad_size = global_fad_size;
580  v = ViewType ("view", num_rows, num_cols, fad_size+1);
581 #endif
582  typename ViewType::array_type va = v;
583  Kokkos::deep_copy( va, 1.0 );
584 
585  // Deep copy a constant scalar
586  value_type a = 2.3456;
587  Kokkos::deep_copy( v, a );
588 
589  // Copy to host
590  host_view_type hv = Kokkos::create_mirror_view(v);
591  Kokkos::deep_copy(hv, v);
592 
593  // Check
594  success = true;
595  for (size_type i=0; i<num_rows; ++i) {
596  for (size_type j=0; j<num_cols; ++j) {
597 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
598  FadType f = FadType(fad_size, a);
599 #else
600  FadType f = a;
601 #endif
602  success = success && checkFads(f, hv(i,j), out);
603  }
604  }
605 }
606 
608  Kokkos_View_Fad, DeepCopy_ConstantZero, FadType, Layout, Device )
609 {
610  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
611  typedef typename ViewType::size_type size_type;
612  typedef typename ViewType::HostMirror host_view_type;
613  typedef typename FadType::value_type value_type;
614 
615  const size_type num_rows = global_num_rows;
616  const size_type num_cols = global_num_cols;
617 
618  // Create and fill view
619  ViewType v;
620 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
621  v = ViewType ("view", num_rows, num_cols);
622 #else
623  const size_type fad_size = global_fad_size;
624  v = ViewType ("view", num_rows, num_cols, fad_size+1);
625 #endif
626  typename ViewType::array_type va = v;
627  Kokkos::deep_copy( va, 1.0 );
628 
629  // Deep copy a constant scalar
630  value_type a = 0.0;
631  Kokkos::deep_copy( v, a );
632 
633  // Copy to host
634  host_view_type hv = Kokkos::create_mirror_view(v);
635  Kokkos::deep_copy(hv, v);
636 
637  // Check
638  success = true;
639  for (size_type i=0; i<num_rows; ++i) {
640  for (size_type j=0; j<num_cols; ++j) {
641 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
642  FadType f = FadType(fad_size, a);
643 #else
644  FadType f = a;
645 #endif
646  success = success && checkFads(f, hv(i,j), out);
647  }
648  }
649 }
650 
652  Kokkos_View_Fad, DeepCopy_ConstantFad, FadType, Layout, Device )
653 {
654  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
655  typedef typename ViewType::size_type size_type;
656  typedef typename ViewType::HostMirror host_view_type;
657 
658  const size_type num_rows = global_num_rows;
659  const size_type num_cols = global_num_cols;
660 
661  // Create and fill view
662  ViewType v;
663 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
664  v = ViewType ("view", num_rows, num_cols);
665 #else
666  const size_type fad_size = global_fad_size;
667  v = ViewType ("view", num_rows, num_cols, fad_size+1);
668 #endif
669  typename ViewType::array_type va = v;
670  Kokkos::deep_copy( va, 1.0 );
671 
672  // Deep copy a constant scalar
673  FadType a = 2.3456;
674  Kokkos::deep_copy( v, a );
675 
676  // Copy to host
677  host_view_type hv = Kokkos::create_mirror_view(v);
678  Kokkos::deep_copy(hv, v);
679 
680  // Check
681  success = true;
682  for (size_type i=0; i<num_rows; ++i) {
683  for (size_type j=0; j<num_cols; ++j) {
684 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
685  FadType f = FadType(fad_size, a.val());
686 #else
687  FadType f = a;
688 #endif
689  success = success && checkFads(f, hv(i,j), out);
690  }
691  }
692 }
693 
695  Kokkos_View_Fad, DeepCopy_ConstantFadFull, FadType, Layout, Device )
696 {
697  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
698  typedef typename ViewType::size_type size_type;
699  typedef typename ViewType::HostMirror host_view_type;
700 
701  const size_type num_rows = global_num_rows;
702  const size_type num_cols = global_num_cols;
703  const size_type fad_size = global_fad_size;
704 
705  // Create and fill view
706  ViewType v;
707 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
708  v = ViewType ("view", num_rows, num_cols);
709 #else
710  v = ViewType ("view", num_rows, num_cols, fad_size+1);
711 #endif
712  typename ViewType::array_type va = v;
713  Kokkos::deep_copy( va, 1.0 );
714 
715  // Deep copy a constant Fad
716  FadType a(fad_size, 2.3456);
717  for (size_type i=0; i<fad_size; ++i)
718  a.fastAccessDx(i) = 7.89 + (i+1);
719 
720  // Copy to host
721  host_view_type hv = Kokkos::create_mirror_view(v);
722  Kokkos::deep_copy(hv, a);
723 
724  // Check
725  success = true;
726  for (size_type i=0; i<num_rows; ++i) {
727  for (size_type j=0; j<num_cols; ++j) {
728  success = success && checkFads(a, hv(i,j), out);
729  }
730  }
731 }
732 
734  Kokkos_View_Fad, LocalDeepCopy, FadType, Layout, Device )
735 {
736  typedef Kokkos::View<FadType***,Layout,Device> ViewType;
737  typedef Kokkos::View<FadType,Layout,Device> ScalarViewType;
738  typedef typename ViewType::size_type size_type;
739  typedef typename ViewType::HostMirror host_view_type;
740  typedef typename ScalarViewType::HostMirror host_scalar_view_type;
741 
742  const size_type num_rows = global_num_rows;
743  const size_type num_cols = global_num_cols;
744  const size_type num_slices = 10;
745  const size_type fad_size = global_fad_size;
746 
747  // Create and fill view
748  ViewType v;
749 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
750  v = ViewType ("view", num_rows, num_cols, num_slices);
751 #else
752  v = ViewType ("view", num_rows, num_cols, num_slices, fad_size+1);
753 #endif
754  typename ViewType::array_type va = v;
755  Kokkos::deep_copy( va, 1.0 );
756 
757  // Deep copy a constant Fad to the device
758  // Can't deep_copy directly because that doesn't work with DFad
759  FadType a(fad_size, 2.3456);
760  for (size_type i=0; i<fad_size; ++i)
761  a.fastAccessDx(i) = 7.89 + (i+1);
762 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
763  ScalarViewType a_view("a");
764 #else
765  ScalarViewType a_view("a", fad_size+1);
766 #endif
767  host_scalar_view_type ha_view = Kokkos::create_mirror_view(a_view);
768  ha_view() = a;
769  Kokkos::deep_copy( a_view, ha_view );
770 
771  // Excersize local_deep_copy by setting each row of s to a
772  Kokkos::parallel_for(Kokkos::RangePolicy<Device>(0,num_rows),
773  KOKKOS_LAMBDA(const int i)
774  {
775  auto s = Kokkos::subview(v,i,Kokkos::ALL,Kokkos::ALL);
776  Kokkos::Experimental::local_deep_copy(s,a_view());
777  });
778 
779  // Copy back to host
780  host_view_type hv = Kokkos::create_mirror_view(v);
781  Kokkos::deep_copy(hv, a);
782 
783  // Check
784  success = true;
785  for (size_type i=0; i<num_rows; ++i) {
786  for (size_type j=0; j<num_cols; ++j) {
787  for (size_type k=0; k<num_slices; ++k) {
788  success = success && checkFads(a, hv(i,j,k), out);
789  }
790  }
791  }
792 }
793 
795  Kokkos_View_Fad, LocalDeepCopyTeam, FadType, Layout, Device )
796 {
797  typedef Kokkos::View<FadType***,Layout,Device> ViewType;
798  typedef Kokkos::View<FadType,Layout,Device> ScalarViewType;
799  typedef typename ViewType::size_type size_type;
800  typedef typename ViewType::HostMirror host_view_type;
801  typedef typename ScalarViewType::HostMirror host_scalar_view_type;
802 
803  const size_type num_rows = global_num_rows;
804  const size_type num_cols = global_num_cols;
805  const size_type num_slices = 10;
806  const size_type fad_size = global_fad_size;
807 
808  // Create and fill view
809  ViewType v;
810 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
811  v = ViewType ("view", num_rows, num_cols, num_slices);
812 #else
813  v = ViewType ("view", num_rows, num_cols, num_slices, fad_size+1);
814 #endif
815  typename ViewType::array_type va = v;
816  Kokkos::deep_copy( va, 1.0 );
817 
818  // Deep copy a constant Fad to the device
819  // Can't deep_copy directly because that doesn't work with DFad
820  FadType a(fad_size, 2.3456);
821  for (size_type i=0; i<fad_size; ++i)
822  a.fastAccessDx(i) = 7.89 + (i+1);
823 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
824  ScalarViewType a_view("a");
825 #else
826  ScalarViewType a_view("a", fad_size+1);
827 #endif
828  host_scalar_view_type ha_view = Kokkos::create_mirror_view(a_view);
829  ha_view() = a;
830  Kokkos::deep_copy( a_view, ha_view );
831 
832  // Excersize local_deep_copy by setting each row of s to a
833  typedef Kokkos::TeamPolicy<Device> Policy;
834  static const size_type stride = Kokkos::ViewScalarStride<ViewType>::stride;
835  Kokkos::parallel_for(Policy(num_rows,Kokkos::AUTO,stride),
836  KOKKOS_LAMBDA(const typename Policy::member_type& team)
837  {
838  int i = team.league_rank();
839  auto s = Kokkos::subview(v,i,Kokkos::ALL,Kokkos::ALL);
840  Kokkos::Experimental::local_deep_copy(team,s,a_view());
841  });
842 
843  // Copy back to host
844  host_view_type hv = Kokkos::create_mirror_view(v);
845  Kokkos::deep_copy(hv, a);
846 
847  // Check
848  success = true;
849  for (size_type i=0; i<num_rows; ++i) {
850  for (size_type j=0; j<num_cols; ++j) {
851  for (size_type k=0; k<num_slices; ++k) {
852  success = success && checkFads(a, hv(i,j,k), out);
853  }
854  }
855  }
856 }
857 
859  Kokkos_View_Fad, ScalarAssign, FadType, Layout, Device )
860 {
861  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
862  typedef typename ViewType::size_type size_type;
863  typedef typename ViewType::HostMirror host_view_type;
864  typedef typename FadType::value_type value_type;
865 
866  const size_type num_rows = global_num_rows;
867 
868  // Create and fill view
869  ViewType v;
870 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
871  v = ViewType ("view", num_rows);
872 #else
873  const size_type fad_size = global_fad_size;
874  v = ViewType ("view", num_rows, fad_size+1);
875 #endif
876  typename ViewType::array_type va = v;
877  Kokkos::deep_copy( va, 1.0 );
878 
879  // Deep copy a constant scalar
880  value_type a = 2.3456;
882 
883  // Copy to host
884  host_view_type hv = Kokkos::create_mirror_view(v);
885  Kokkos::deep_copy(hv, v);
886 
887  // Check
888  success = true;
889  for (size_type i=0; i<num_rows; ++i) {
890 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
891  FadType f = FadType(fad_size, a);
892 #else
893  FadType f = a;
894 #endif
895  success = success && checkFads(f, hv(i), out);
896  }
897 }
898 
900  Kokkos_View_Fad, ValueAssign, FadType, Layout, Device )
901 {
902  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
903  typedef Kokkos::View<FadType,Layout,Device> ScalarViewType;
904  typedef typename ViewType::size_type size_type;
905  typedef typename ViewType::HostMirror host_view_type;
906  typedef typename ScalarViewType::HostMirror host_scalar_view_type;
907 
908  const size_type num_rows = global_num_rows;
909  const size_type fad_size = global_fad_size;
910 
911  // Create and fill view
912  ViewType v;
913  ScalarViewType a;
914 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
915  v = ViewType ("view", num_rows);
916  a = ScalarViewType ("fad");
917 #else
918  v = ViewType ("view", num_rows, fad_size+1);
919  a = ScalarViewType ("fad", fad_size+1);
920 #endif
921  typename ViewType::array_type va = v;
922  Kokkos::deep_copy( va, 1.0 );
923 
924  // Deep copy a constant scalar
925  Kokkos::deep_copy(a, 2.3456);
926 
927  Kokkos::parallel_for(Kokkos::RangePolicy< Device>(0, fad_size), KOKKOS_LAMBDA(const int i) {
928  a().fastAccessDx(i) = 7.89 + i;
929  });
930  Kokkos::fence();
931 
933  Kokkos::fence();
934 
935  // Copy to host
936  host_view_type hv = Kokkos::create_mirror_view(v);
937  Kokkos::deep_copy(hv, v);
938 
939  host_scalar_view_type ha = Kokkos::create_mirror_view(a);
940  Kokkos::deep_copy(ha, a);
941 
942  // Check
943  success = true;
944  for (size_type i=0; i<num_rows; ++i) {
945  success = success && checkFads(ha(), hv(i), out);
946  }
947 }
948 
950  Kokkos_View_Fad, Resize, FadType, Layout, Device )
951 {
952  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
953  typedef typename ViewType::size_type size_type;
954  typedef typename ViewType::HostMirror host_view_type;
955 
956  const size_type num_rows = global_num_rows;
957  const size_type num_cols = global_num_cols;
958  const size_type fad_size = global_fad_size;
959 
960  // Create and fill view
961  ViewType v;
962 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
963  v = ViewType ("view", num_rows, num_cols);
964 #else
965  v = ViewType ("view", num_rows, num_cols, fad_size+1);
966 #endif
967  host_view_type h_v = Kokkos::create_mirror_view(v);
968  for (size_type i=0; i<num_rows; ++i)
969  for (size_type j=0; j<num_cols; ++j)
970  h_v(i,j) = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
971  Kokkos::deep_copy(v, h_v);
972 
973  // Resize
974 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
975  Kokkos::resize(v, num_rows, num_cols+1);
976 #else
977  Kokkos::resize(v, num_rows, num_cols+1, fad_size+1);
978 #endif
979 
980  // Copy back
981  host_view_type h_v2 = Kokkos::create_mirror_view(v);
982  Kokkos::deep_copy(h_v2, v);
983 
984  // Check
985  success = true;
986  for (size_type i=0; i<num_rows; ++i) {
987  for (size_type j=0; j<num_cols; ++j) {
988  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
989  success = success && checkFads(f, h_v2(i,j), out);
990  }
991 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
992  FadType f = 0.0;
993 #else
994  FadType f(fad_size, 0.0);
995 #endif
996  success = success && checkFads(f, h_v2(i,num_cols), out);
997  }
998 }
999 
1001  Kokkos_View_Fad, Multiply, FadType, Layout, Device )
1002 {
1003  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
1004  typedef typename ViewType::size_type size_type;
1005  typedef typename ViewType::HostMirror host_view_type;
1006 
1007  const size_type num_rows = global_num_rows;
1008  const size_type fad_size = global_fad_size;
1009 
1010  // Create and fill views
1011  ViewType v1, v2;
1012 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1013  v1 = ViewType ("view1", num_rows);
1014  v2 = ViewType ("view2", num_rows);
1015 #else
1016  v1 = ViewType ("view1", num_rows, fad_size+1);
1017  v2 = ViewType ("view2", num_rows, fad_size+1);
1018 #endif
1019  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
1020  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
1021  for (size_type i=0; i<num_rows; ++i) {
1022  h_v1(i) = generate_fad<FadType>(
1023  num_rows, size_type(2), fad_size, i, size_type(0));
1024  h_v2(i) = generate_fad<FadType>(
1025  num_rows, size_type(2), fad_size, i, size_type(1));
1026  }
1027  Kokkos::deep_copy(v1, h_v1);
1028  Kokkos::deep_copy(v2, h_v2);
1029 
1030  // Launch kernel
1031  ViewType v3;
1032 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1033  v3 = ViewType ("view3", num_rows);
1034 #else
1035  v3 = ViewType ("view3", num_rows, fad_size+1);
1036 #endif
1038 
1039  // Copy back
1040  host_view_type h_v3 = Kokkos::create_mirror_view(v3);
1041  Kokkos::deep_copy(h_v3, v3);
1042 
1043  // Check
1044  success = true;
1045  for (size_type i=0; i<num_rows; ++i) {
1046  FadType f1 =
1047  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(0));
1048  FadType f2 =
1049  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(1));
1050  FadType f3 = f1*f2;
1051  success = success && checkFads(f3, h_v3(i), out);
1052  }
1053 }
1054 
1056  Kokkos_View_Fad, MultiplyUpdate, FadType, Layout, Device )
1057 {
1058  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
1059  typedef typename ViewType::size_type size_type;
1060  typedef typename ViewType::HostMirror host_view_type;
1061 
1062  const size_type num_rows = global_num_rows;
1063  const size_type fad_size = global_fad_size;
1064 
1065  // Create and fill views
1066  ViewType v1, v2;
1067 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1068  v1 = ViewType ("view1", num_rows);
1069  v2 = ViewType ("view2", num_rows);
1070 #else
1071  v1 = ViewType ("view1", num_rows, fad_size+1);
1072  v2 = ViewType ("view2", num_rows, fad_size+1);
1073 #endif
1074  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
1075  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
1076  for (size_type i=0; i<num_rows; ++i) {
1077  h_v1(i) = generate_fad<FadType>(
1078  num_rows, size_type(2), fad_size, i, size_type(0));
1079  h_v2(i) = generate_fad<FadType>(
1080  num_rows, size_type(2), fad_size, i, size_type(1));
1081  }
1082  Kokkos::deep_copy(v1, h_v1);
1083  Kokkos::deep_copy(v2, h_v2);
1084 
1085  // Launch kernel
1086  ViewType v3;
1087 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1088  v3 = ViewType ("view3", num_rows);
1089 #else
1090  v3 = ViewType ("view3", num_rows, fad_size+1);
1091 #endif
1092  Kokkos::deep_copy(v3, 1.0);
1093  MultiplyKernel<ViewType>::apply(v1,v2,v3,true);
1094 
1095  // Copy back
1096  host_view_type h_v3 = Kokkos::create_mirror_view(v3);
1097  Kokkos::deep_copy(h_v3, v3);
1098 
1099  // Check
1100  success = true;
1101  for (size_type i=0; i<num_rows; ++i) {
1102  FadType f1 =
1103  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(0));
1104  FadType f2 =
1105  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(1));
1106  FadType f3 = 1.0 + f1*f2;
1107  success = success && checkFads(f3, h_v3(i), out);
1108  }
1109 }
1110 
1112  Kokkos_View_Fad, MultiplyConst, FadType, Layout, Device )
1113 {
1114  typedef Kokkos::View<const FadType*,Layout,Device,Kokkos::MemoryUnmanaged> ConstViewType;
1115  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
1116  typedef typename ViewType::size_type size_type;
1117  typedef typename ViewType::HostMirror host_view_type;
1118 
1119  const size_type num_rows = global_num_rows;
1120  const size_type fad_size = global_fad_size;
1121 
1122  // Create and fill views
1123  ViewType v1, v2;
1124 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1125  v1 = ViewType ("view1", num_rows);
1126  v2 = ViewType ("view2", num_rows);
1127 #else
1128  v1 = ViewType ("view1", num_rows, fad_size+1);
1129  v2 = ViewType ("view2", num_rows, fad_size+1);
1130 #endif
1131  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
1132  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
1133  for (size_type i=0; i<num_rows; ++i) {
1134  h_v1(i) = generate_fad<FadType>(
1135  num_rows, size_type(2), fad_size, i, size_type(0));
1136  h_v2(i) = generate_fad<FadType>(
1137  num_rows, size_type(2), fad_size, i, size_type(1));
1138  }
1139  Kokkos::deep_copy(v1, h_v1);
1140  Kokkos::deep_copy(v2, h_v2);
1141 
1142  ConstViewType cv1 = v1;
1143 
1144  // Launch kernel
1145  ViewType v3;
1146 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1147  v3 = ViewType ("view3", num_rows);
1148 #else
1149  v3 = ViewType ("view3", num_rows, fad_size+1);
1150 #endif
1152 
1153  // Copy back
1154  host_view_type h_v3 = Kokkos::create_mirror_view(v3);
1155  Kokkos::deep_copy(h_v3, v3);
1156 
1157  // Check
1158  success = true;
1159  for (size_type i=0; i<num_rows; ++i) {
1160  FadType f1 =
1161  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(0));
1162  FadType f2 =
1163  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(1));
1164  FadType f3 = f1*f2;
1165  success = success && checkFads(f3, h_v3(i), out);
1166  }
1167 }
1168 
1170  Kokkos_View_Fad, MultiplyMixed, FadType, Layout, Device )
1171 {
1172  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
1173  typedef typename ViewType::size_type size_type;
1174  typedef typename ViewType::HostMirror host_view_type;
1175 
1176  const size_type num_rows = 2;
1177  const size_type fad_size = global_fad_size;
1178 
1179  // Create and fill views -- do everything on the host for this test
1180  FadType f0 = generate_fad<FadType>(
1181  num_rows, size_type(2), fad_size, size_type(0), size_type(0));
1182  FadType f1 = generate_fad<FadType>(
1183  num_rows, size_type(2), fad_size, size_type(1), size_type(0));
1184  host_view_type h_v;
1185 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1186  h_v = host_view_type ("view1", num_rows);
1187 #else
1188  h_v = host_view_type ("view1", num_rows, fad_size+1);
1189 #endif
1190  h_v(0) = f0;
1191  h_v(1) = f1;
1192 
1193  FadType f2 = f0 * h_v(1);
1194 
1195  // Check
1196  FadType f3 = f0 * f1;
1197  success = checkFads(f3, f2, out);
1198 }
1199 
1201  Kokkos_View_Fad, AtomicAdd, FadType, Layout, Device )
1202 {
1203  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
1204  typedef Kokkos::View<FadType,Layout,Device> ScalarViewType;
1205  typedef typename ViewType::size_type size_type;
1206  typedef typename ScalarViewType::HostMirror host_scalar_view_type;
1207 
1208  const size_type num_rows = global_num_rows;
1209  const size_type fad_size = global_fad_size;
1210 
1211  // Create and fill view
1212  ViewType v;
1213 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1214  v = ViewType ("view", num_rows);
1215 #else
1216  v = ViewType ("view", num_rows, fad_size+1);
1217 #endif
1218  Kokkos::deep_copy(v, 2.3456);
1219 
1220  Kokkos::parallel_for(Kokkos::RangePolicy<Device>(0, num_rows), KOKKOS_LAMBDA(const size_type i) {
1221  for (size_type j = 0; j < fad_size; ++j)
1222  v(i).fastAccessDx(j) = 7.89 + j;
1223  });
1224 
1225  // Create scalar view
1226  ScalarViewType s;
1227 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1228  s = ScalarViewType ("scalar view");
1229 #else
1230  s = ScalarViewType ("scalar view", fad_size+1);
1231 #endif
1232 
1233  // Call atomic_add kernel, which adds up entries in v
1235 
1236  // Copy to host
1237  host_scalar_view_type hs = Kokkos::create_mirror_view(s);
1238  Kokkos::deep_copy(hs, s);
1239 
1240  // Check
1241  auto hv = Kokkos::create_mirror_view(v);
1242  Kokkos::deep_copy(hv, v);
1243 
1244  FadType b = num_rows*hv(0);
1245  success = checkFads(b, hs(), out);
1246 }
1247 
1249  Kokkos_View_Fad, Rank8, FadType, Layout, Device )
1250 {
1251  typedef Kokkos::View<FadType*******,Layout,Device> ViewType;
1252  typedef typename ViewType::size_type size_type;
1253  typedef typename ViewType::HostMirror host_view_type;
1254 
1255  const size_type fad_size = global_fad_size;
1256 
1257  // Create and fill view
1258  ViewType v;
1259 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1260  v = ViewType ("view", 100, 1, 2, 3, 4, 5, 6);
1261 #else
1262  v = ViewType ("view", 100, 1, 2, 3, 4, 5, 6, fad_size+1);
1263 #endif
1264  host_view_type h_v = Kokkos::create_mirror_view(v);
1265  typename host_view_type::array_type h_a = h_v;
1266  Kokkos::deep_copy(h_a, 1.0);
1267 
1268  FadType f1 = FadType(fad_size, 2.0);
1269  h_v(99,0,1,2,3,4,5) = f1;
1270  FadType f2 = h_v(99,0,1,2,3,4,5);
1271 
1272  // Check
1273  success = checkFads(f1, f2, out);
1274 }
1275 
1277  Kokkos_View_Fad, Roger, FadType, Layout, Device )
1278 {
1279  Kokkos::View<FadType*,Layout,Device> a;
1280  Kokkos::View<FadType**,Layout,Device> b;
1281  Kokkos::View<FadType***,Layout,Device> c;
1282  Kokkos::View<FadType****,Layout,Device> d;
1283  Kokkos::View<FadType*****,Layout,Device> e;
1284  Kokkos::View<FadType******,Layout,Device> f;
1285  Kokkos::View<FadType*******,Layout,Device> g;
1286 
1287 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1288  a = Kokkos::View<FadType*,Layout,Device>("a",4);
1289  b = Kokkos::View<FadType**,Layout,Device> ("b",4,4);
1290  c = Kokkos::View<FadType***,Layout,Device> ("c",4,4,4);
1291  d = Kokkos::View<FadType****,Layout,Device> ("d",4,4,4,4);
1292  e = Kokkos::View<FadType*****,Layout,Device> ("e",4,4,4,4,4);
1293  f = Kokkos::View<FadType******,Layout,Device> ("f",4,4,4,4,4,4);
1294  g = Kokkos::View<FadType*******,Layout,Device> ("g",4,4,4,4,4,4,4);
1295 #else
1296  const unsigned fad_size = global_fad_size;
1297  a = Kokkos::View<FadType*,Layout,Device>("a",4,fad_size+1);
1298  b = Kokkos::View<FadType**,Layout,Device> ("b",4,4,fad_size+1);
1299  c = Kokkos::View<FadType***,Layout,Device> ("c",4,4,4,fad_size+1);
1300  d = Kokkos::View<FadType****,Layout,Device> ("d",4,4,4,4,fad_size+1);
1301  e = Kokkos::View<FadType*****,Layout,Device> ("e",4,4,4,4,4,fad_size+1);
1302  f = Kokkos::View<FadType******,Layout,Device> ("f",4,4,4,4,4,4,fad_size+1);
1303  g = Kokkos::View<FadType*******,Layout,Device> ("g",4,4,4,4,4,4,4,fad_size+1);
1304 #endif
1305 
1306  typedef typename Device::memory_space memory_space;
1307  const bool is_accessible =
1308  Kokkos::Impl::MemorySpaceAccess<Kokkos::HostSpace,
1309  memory_space>::accessible;
1310  if (is_accessible) {
1311  a(0) = FadType(1.0);
1312  f(0,0,0,0,0,0) = FadType(1.0);
1313  g(0,0,0,0,0,0,0) = FadType(1.0);
1314  }
1315 
1316  // Check
1317  success = true;
1318 }
1319 
1321  Kokkos_View_Fad, AssignDifferentStrides, FadType, Layout, Device )
1322 {
1323  typedef Kokkos::View<FadType**,Layout,Device> ViewType1;
1324  typedef Kokkos::View<FadType*,Layout,Device> ViewType2;
1325  typedef typename ViewType1::size_type size_type;
1326  typedef typename ViewType1::HostMirror host_view_type1;
1327  typedef typename ViewType2::HostMirror host_view_type2;
1328 
1329  const size_type num_rows = global_num_rows;
1330  const size_type num_cols = global_num_cols;
1331  const size_type fad_size = global_fad_size;
1332 
1333  // Create and fill views
1334  ViewType1 v1;
1335 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1336  v1 = ViewType1 ("view1", num_rows, num_cols);
1337 #else
1338  v1 = ViewType1 ("view1", num_rows, num_cols, fad_size+1);
1339 #endif
1340  host_view_type1 h_v1 = Kokkos::create_mirror_view(v1);
1341  for (size_type i=0; i<num_rows; ++i) {
1342  for (size_type j=0; j<num_cols; ++j) {
1343  h_v1(i,j) = generate_fad<FadType>(
1344  num_rows, num_cols, fad_size, i, j);
1345  }
1346  }
1347  Kokkos::deep_copy(v1, h_v1);
1348 
1349  // Launch kernel
1350  ViewType2 v2;
1351 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1352  v2 = ViewType2 ("view2", num_rows);
1353 #else
1354  v2 = ViewType2 ("view2", num_rows, fad_size+1);
1355 #endif
1357 
1358  // Copy back
1359  host_view_type2 h_v2 = Kokkos::create_mirror_view(v2);
1360  Kokkos::deep_copy(h_v2, v2);
1361 
1362  // Check
1363  success = true;
1364  for (size_type i=0; i<num_rows; ++i) {
1365  FadType f =
1366  generate_fad<FadType>(num_rows, num_cols, fad_size, i, size_type(1));
1367  success = success && checkFads(f, h_v2(i), out);
1368  }
1369 }
1370 
1372  Kokkos_View_Fad, ScalarValue, FadType, Layout, Device )
1373 {
1374  typedef typename Sacado::ScalarType<FadType>::type ScalarType;
1375  typedef Kokkos::View<FadType,Layout,Device> ViewType1;
1376  typedef Kokkos::View<ScalarType,Layout,Device> ViewType2;
1377  typedef typename ViewType1::size_type size_type;
1378  typedef typename ViewType1::HostMirror host_view_type1;
1379  typedef typename ViewType2::HostMirror host_view_type2;
1380 
1381  const int fad_size = global_fad_size;
1382 
1383  // Create and fill views
1384  ViewType1 v1;
1385 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1386  v1 = ViewType1 ("view1");
1387 #else
1388  v1 = ViewType1 ("view1", fad_size+1);
1389 #endif
1390  host_view_type1 h_v1 = Kokkos::create_mirror_view(v1);
1391  h_v1() = generate_fad<FadType>(1, 1, fad_size, 0, 0);
1392  Kokkos::deep_copy(v1, h_v1);
1393 
1394  // Launch kernel
1395  ViewType2 v2 = ViewType2 ("view2");
1396  Kokkos::parallel_for(Kokkos::RangePolicy<Device>(0,1),
1397  KOKKOS_LAMBDA(const size_type i)
1398  {
1399  v2() = Sacado::scalarValue(v1());
1400  });
1401 
1402  // Copy back
1403  host_view_type2 h_v2 = Kokkos::create_mirror_view(v2);
1404  Kokkos::deep_copy(h_v2, v2);
1405 
1406  // Check
1407  success = true;
1408  TEUCHOS_TEST_EQUALITY(h_v1().val(), h_v2(), out, success);
1409 }
1410 
1411 #if defined(HAVE_SACADO_KOKKOS) && defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
1412 
1414  Kokkos_View_Fad, DynRankDimensionScalar, FadType, Layout, Device )
1415 {
1416  typedef Kokkos::DynRankView<double,Layout,Device> DoubleViewType;
1417  typedef Kokkos::DynRankView<FadType,Layout,Device> FadViewType;
1418  typedef typename FadViewType::size_type size_type;
1419 
1420  const size_type num_rows = global_num_rows;
1421  const size_type fad_size = global_fad_size;
1422 
1423  // Create views
1424  DoubleViewType v1("view1", num_rows);
1425  FadViewType v2 ("view2", num_rows, fad_size+1);
1426 
1427  // Check dimension scalar works
1428  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(v1), 0, out, success);
1429  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(v2), fad_size+1, out, success);
1430 }
1431 
1433  Kokkos_View_Fad, DynRankAssignStatic0, FadType, Layout, Device )
1434 {
1435  typedef Kokkos::View<FadType,Layout,Device> StaticViewType;
1436  typedef Kokkos::DynRankView<FadType,Layout,Device> DynamicViewType;
1437  typedef typename StaticViewType::size_type size_type;
1438 
1439  const size_type num_rows = global_num_rows;
1440  const size_type num_cols = global_num_cols;
1441  const size_type fad_size = global_fad_size;
1442 
1443  // Create and fill views
1444  StaticViewType v1("view", fad_size+1);
1445  auto h_v1 = Kokkos::create_mirror_view(v1);
1446  h_v1() = generate_fad<FadType>(num_rows, num_cols, fad_size, size_type(0), size_type(0));
1447  Kokkos::deep_copy(v1, h_v1);
1448 
1449  // Assign static to dynamic
1450  DynamicViewType v2 = v1;
1451 
1452  // Copy back
1453  auto h_v2 = Kokkos::create_mirror_view(v2);
1454  Kokkos::deep_copy(h_v2, v2);
1455 
1456  // Check dimensions are correct
1457  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(v2), fad_size+1, out, success);
1458  TEUCHOS_TEST_EQUALITY(v2.stride_0(), v1.stride_0(), out, success);
1459 
1460  // Check values
1461  FadType f =
1462  generate_fad<FadType>(num_rows, num_cols, fad_size, size_type(0), size_type(0));
1463  success = success && checkFads(f, h_v2(), out);
1464 }
1465 
1467  Kokkos_View_Fad, DynRankAssignStatic1, FadType, Layout, Device )
1468 {
1469  typedef Kokkos::View<FadType*,Layout,Device> StaticViewType;
1470  typedef Kokkos::DynRankView<FadType,Layout,Device> DynamicViewType;
1471  typedef typename StaticViewType::size_type size_type;
1472 
1473  const size_type num_rows = global_num_rows;
1474  const size_type num_cols = global_num_cols;
1475  const size_type fad_size = global_fad_size;
1476 
1477  // Create and fill views
1478  StaticViewType v1("view", num_rows, fad_size+1);
1479  auto h_v1 = Kokkos::create_mirror_view(v1);
1480  for (size_type i=0; i<num_rows; ++i)
1481  h_v1(i) =
1482  generate_fad<FadType>(num_rows, num_cols, fad_size, i, size_type(0));
1483  Kokkos::deep_copy(v1, h_v1);
1484 
1485  // Assign static to dynamic
1486  DynamicViewType v2 = v1;
1487 
1488  // Copy back
1489  auto h_v2 = Kokkos::create_mirror_view(v2);
1490  Kokkos::deep_copy(h_v2, v2);
1491 
1492  // Check dimensions are correct
1493  TEUCHOS_TEST_EQUALITY(v2.extent(0), num_rows, out, success);
1494  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(v2), fad_size+1, out, success);
1495  TEUCHOS_TEST_EQUALITY(v2.stride_0(), v1.stride_0(), out, success);
1496  TEUCHOS_TEST_EQUALITY(v2.stride_1(), v1.stride_1(), out, success);
1497 
1498  // Check values
1499  for (size_type i=0; i<num_rows; ++i) {
1500  FadType f =
1501  generate_fad<FadType>(num_rows, num_cols, fad_size, i, size_type(0));
1502  success = success && checkFads(f, h_v2(i), out);
1503  }
1504 }
1505 
1507  Kokkos_View_Fad, DynRankAssignStatic2, FadType, Layout, Device )
1508 {
1509  typedef Kokkos::View<FadType**,Layout,Device> StaticViewType;
1510  typedef Kokkos::DynRankView<FadType,Layout,Device> DynamicViewType;
1511  typedef typename StaticViewType::size_type size_type;
1512 
1513  const size_type num_rows = global_num_rows;
1514  const size_type num_cols = global_num_cols;
1515  const size_type fad_size = global_fad_size;
1516 
1517  // Create and fill views
1518  StaticViewType v1("view", num_rows, num_cols, fad_size+1);
1519  auto h_v1 = Kokkos::create_mirror_view(v1);
1520  for (size_type i=0; i<num_rows; ++i)
1521  for (size_type j=0; j<num_cols; ++j)
1522  h_v1(i,j) = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
1523  Kokkos::deep_copy(v1, h_v1);
1524 
1525  // Assign static to dynamic
1526  DynamicViewType v2 = v1;
1527 
1528  // Copy back
1529  auto h_v2 = Kokkos::create_mirror_view(v2);
1530  Kokkos::deep_copy(h_v2, v2);
1531 
1532  // Check dimensions are correct
1533  TEUCHOS_TEST_EQUALITY(v2.extent(0), num_rows, out, success);
1534  TEUCHOS_TEST_EQUALITY(v2.extent(1), num_cols, out, success);
1535  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(v2), fad_size+1, out, success);
1536  TEUCHOS_TEST_EQUALITY(v2.stride_0(), v1.stride_0(), out, success);
1537  TEUCHOS_TEST_EQUALITY(v2.stride_1(), v1.stride_1(), out, success);
1538  TEUCHOS_TEST_EQUALITY(v2.stride_2(), v1.stride_2(), out, success);
1539 
1540  // Check values
1541  for (size_type i=0; i<num_rows; ++i) {
1542  for (size_type j=0; j<num_cols; ++j) {
1543  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
1544  success = success && checkFads(f, h_v2(i,j), out);
1545  }
1546  }
1547 }
1548 
1550  Kokkos_View_Fad, DynRankMultiply, FadType, Layout, Device )
1551 {
1552  typedef Kokkos::DynRankView<FadType,Layout,Device> ViewType;
1553  typedef typename ViewType::size_type size_type;
1554  typedef typename ViewType::HostMirror host_view_type;
1555 
1556  const size_type num_rows = global_num_rows;
1557  const size_type fad_size = global_fad_size;
1558 
1559  // Create and fill views
1560  ViewType v1("view1", num_rows, fad_size+1);
1561  ViewType v2("view2", num_rows, fad_size+1);
1562  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
1563  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
1564  for (size_type i=0; i<num_rows; ++i) {
1565  h_v1(i) = generate_fad<FadType>(
1566  num_rows, size_type(2), fad_size, i, size_type(0));
1567  h_v2(i) = generate_fad<FadType>(
1568  num_rows, size_type(2), fad_size, i, size_type(1));
1569  }
1570  Kokkos::deep_copy(v1, h_v1);
1571  Kokkos::deep_copy(v2, h_v2);
1572 
1573  // Launch kernel
1574  ViewType v3("view3", num_rows, fad_size+1);
1576 
1577  // Copy back
1578  host_view_type h_v3 = Kokkos::create_mirror_view(v3);
1579  Kokkos::deep_copy(h_v3, v3);
1580 
1581  // Check
1582  success = true;
1583  TEUCHOS_TEST_EQUALITY(v3.rank(), 1, out, success);
1584  for (size_type i=0; i<num_rows; ++i) {
1585  FadType f1 =
1586  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(0));
1587  FadType f2 =
1588  generate_fad<FadType>(num_rows, size_type(2), fad_size, i, size_type(1));
1589  FadType f3 = f1*f2;
1590  success = success && checkFads(f3, h_v3(i), out);
1591  }
1592 }
1593 
1595  Kokkos_View_Fad, SubdynrankviewCol, FadType, Layout, Device )
1596 {
1597  typedef Kokkos::DynRankView<FadType,Layout,Device> ViewType;
1598  typedef typename ViewType::size_type size_type;
1599  typedef typename ViewType::HostMirror host_view_type;
1600 
1601  const size_type num_rows = global_num_rows;
1602  const size_type num_cols = global_num_cols;
1603  const size_type fad_size = global_fad_size;
1604 
1605  // Create and fill view
1606  ViewType v("view", num_rows, num_cols, fad_size+1);
1607  host_view_type h_v = Kokkos::create_mirror_view(v);
1608  for (size_type i=0; i<num_rows; ++i) {
1609  for (size_type j=0; j<num_cols; ++j) {
1610  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
1611  h_v(i,j) = f;
1612  }
1613  }
1614  Kokkos::deep_copy(v, h_v);
1615 
1616  // Create subview of first column
1617  size_type col = 1;
1618  auto s = Kokkos::subdynrankview(v, Kokkos::ALL(), col);
1619 
1620  // Copy back
1621  typedef decltype(s) SubviewType;
1622  typedef typename SubviewType::HostMirror HostSubviewType;
1623 
1624  // Note: don't create h_s through create_mirror_view and deep_copy
1625  // since Kokkos doesn't support deep_copy of non-contiguous views
1626  //HostSubviewType h_s = Kokkos::create_mirror_view(s);
1627  //Kokkos::deep_copy(h_s, s);
1628  HostSubviewType h_s = Kokkos::subdynrankview(h_v, Kokkos::ALL(), col);
1629 
1630  // Check
1631  success = true;
1632  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(s), fad_size+1, out, success);
1633  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_s), fad_size+1, out, success);
1634  TEUCHOS_TEST_EQUALITY(h_s.extent(0), num_rows, out, success);
1635  TEUCHOS_TEST_EQUALITY(h_s.extent(1), 1, out, success);
1636  TEUCHOS_TEST_EQUALITY(h_s.extent(7), 1, out, success);
1637 
1638  for (size_type i=0; i<num_rows; ++i) {
1639  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, col);
1640  success = success && checkFads(f, h_s(i), out);
1641  }
1642 }
1643 
1645  Kokkos_View_Fad, SubdynrankviewRow, FadType, Layout, Device )
1646 {
1647  typedef Kokkos::DynRankView<FadType,Layout,Device> ViewType;
1648  typedef typename ViewType::size_type size_type;
1649  typedef typename ViewType::HostMirror host_view_type;
1650 
1651  const size_type num_rows = global_num_rows;
1652  const size_type num_cols = global_num_cols;
1653  const size_type num_planes = 9;
1654  const size_type fad_size = global_fad_size;
1655 
1656  // Create and fill view
1657  ViewType v("view", num_rows, num_cols, num_planes, fad_size+1);
1658  host_view_type h_v = Kokkos::create_mirror_view(v);
1659  for (size_type i=0; i<num_rows; ++i) {
1660  for (size_type j=0; j<num_cols; ++j) {
1661  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
1662  for (size_type k=0; k<num_planes; ++k) {
1663  h_v(i,j,k) = (k+1)*f;
1664  }
1665  }
1666  }
1667  Kokkos::deep_copy(v, h_v);
1668 
1669  // Create subview of first column
1670  size_type row = 2;
1671  auto s = Kokkos::subdynrankview(v, row, Kokkos::ALL(), Kokkos::ALL());
1672 
1673  // Copy back
1674  typedef decltype(s) SubviewType;
1675  typedef typename SubviewType::HostMirror HostSubviewType;
1676 
1677  // Note: don't create h_s through create_mirror_view and deep_copy
1678  // since Kokkos doesn't support deep_copy of non-contiguous views
1679  //HostSubviewType h_s = Kokkos::create_mirror_view(s);
1680  //Kokkos::deep_copy(h_s, s);
1681  HostSubviewType h_s =
1682  Kokkos::subdynrankview(h_v, row, Kokkos::ALL(), Kokkos::ALL());
1683 
1684  // Check
1685  success = true;
1686  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(s), fad_size+1, out, success);
1687  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_s), fad_size+1, out, success);
1688  TEUCHOS_TEST_EQUALITY(h_s.extent(0), num_cols, out, success);
1689  TEUCHOS_TEST_EQUALITY(h_s.extent(1), num_planes, out, success);
1690  TEUCHOS_TEST_EQUALITY(h_s.extent(2), 1, out, success);
1691  TEUCHOS_TEST_EQUALITY(h_s.extent(7), 1, out, success);
1692 
1693  for (size_type j=0; j<num_cols; ++j) {
1694  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, row, j);
1695  for (size_type k=0; k<num_planes; ++k) {
1696  FadType g = (k+1)*f;
1697  success = success && checkFads(g, h_s(j,k), out);
1698  }
1699  }
1700 }
1701 
1703  Kokkos_View_Fad, SubdynrankviewScalar, FadType, Layout, Device )
1704 {
1705  typedef Kokkos::DynRankView<FadType,Layout,Device> ViewType;
1706  typedef typename ViewType::size_type size_type;
1707  typedef typename ViewType::HostMirror host_view_type;
1708 
1709  const size_type num_rows = global_num_rows;
1710  const size_type num_cols = global_num_cols;
1711  const size_type fad_size = global_fad_size;
1712 
1713  // Create and fill view
1714  ViewType v("view", num_rows, num_cols, fad_size+1);
1715  host_view_type h_v = Kokkos::create_mirror_view(v);
1716  for (size_type i=0; i<num_rows; ++i) {
1717  for (size_type j=0; j<num_cols; ++j) {
1718  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
1719  h_v(i,j) = f;
1720  }
1721  }
1722  Kokkos::deep_copy(v, h_v);
1723 
1724  // Create subview of first column
1725  size_type row = 3;
1726  size_type col = 1;
1727  auto s = Kokkos::subdynrankview(v, row, col);
1728 
1729  // Copy back
1730  typedef decltype(s) SubviewType;
1731  typedef typename SubviewType::HostMirror HostSubviewType;
1732 
1733  // Note: don't create h_s through create_mirror_view and deep_copy
1734  // since Kokkos doesn't support deep_copy of non-contiguous views
1735  //HostSubviewType h_s = Kokkos::create_mirror_view(s);
1736  //Kokkos::deep_copy(h_s, s);
1737  HostSubviewType h_s = Kokkos::subdynrankview(h_v, row, col);
1738 
1739  // Check
1740  success = true;
1741  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(s), fad_size+1, out, success);
1742  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_s), fad_size+1, out, success);
1743  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, row, col);
1744  success = success && checkFads(f, h_s(), out);
1745 }
1746 
1747 #else
1748 
1750  Kokkos_View_Fad, DynRankDimensionScalar, FadType, Layout, Device ) {}
1752  Kokkos_View_Fad, DynRankAssignStatic0, FadType, Layout, Device ) {}
1754  Kokkos_View_Fad, DynRankAssignStatic1, FadType, Layout, Device ) {}
1756  Kokkos_View_Fad, DynRankAssignStatic2, FadType, Layout, Device ) {}
1758  Kokkos_View_Fad, DynRankMultiply, FadType, Layout, Device ) {}
1760  Kokkos_View_Fad, SubdynrankviewCol, FadType, Layout, Device ) {}
1762  Kokkos_View_Fad, SubdynrankviewRow, FadType, Layout, Device ) {}
1764  Kokkos_View_Fad, SubdynrankviewScalar, FadType, Layout, Device ) {}
1765 
1766 #endif
1767 
1769  Kokkos_View_Fad, Subview, FadType, Layout, Device )
1770 {
1771  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
1772  typedef typename ViewType::size_type size_type;
1773  typedef typename ViewType::HostMirror host_view_type;
1774 
1775  const size_type num_rows = global_num_rows;
1776  const size_type num_cols = global_num_cols;
1777  const size_type fad_size = global_fad_size;
1778 
1779  // Create and fill view
1780  ViewType v;
1781 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1782  v = ViewType ("view", num_rows, num_cols);
1783 #else
1784  v = ViewType ("view", num_rows, num_cols, fad_size+1);
1785 #endif
1786  host_view_type h_v = Kokkos::create_mirror_view(v);
1787  for (size_type i=0; i<num_rows; ++i) {
1788  for (size_type j=0; j<num_cols; ++j) {
1789  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
1790  h_v(i,j) = f;
1791  }
1792  }
1793  Kokkos::deep_copy(v, h_v);
1794 
1795  // Create subview of first column
1796  size_type col = 1;
1797  auto s = Kokkos::subview(v, Kokkos::ALL(), col);
1798 
1799  // Copy back
1800  typedef decltype(s) SubviewType;
1801  typedef typename SubviewType::HostMirror HostSubviewType;
1802 
1803  // Note: don't create h_s through create_mirror_view and deep_copy
1804  // since Kokkos doesn't support deep_copy of non-contiguous views
1805  //HostSubviewType h_s = Kokkos::create_mirror_view(s);
1806  //Kokkos::deep_copy(h_s, s);
1807  HostSubviewType h_s = Kokkos::subview(h_v, Kokkos::ALL(), col);
1808 
1809  // Check
1810  success = true;
1811 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
1812  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(s), fad_size+1, out, success);
1813  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_s), fad_size+1, out, success);
1814 #endif
1815  for (size_type i=0; i<num_rows; ++i) {
1816  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, col);
1817  success = success && checkFads(f, h_s(i), out);
1818  }
1819 }
1820 
1822  Kokkos_View_Fad, Subview2, FadType, Layout, Device )
1823 {
1824  typedef Kokkos::View<FadType***,Layout,Device> ViewType;
1825  typedef typename ViewType::HostMirror host_view_type;
1826 
1827  // Test various subview operations to check the resulting indexing is correct.
1828  // We only need to run these tests on the host because the indexing does
1829  // not depend on the device.
1830 
1831  const int num_cell = 5;
1832  const int num_qp = 4;
1833  const int num_dim = 3;
1834  const int num_deriv = 2;
1835 
1836  // Create and fill view
1837  host_view_type v;
1838 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1839  v = host_view_type ("view", num_cell, num_qp, num_dim);
1840 #else
1841  v = host_view_type ("view", num_cell, num_qp, num_dim, num_deriv+1);
1842 #endif
1843  for (int cell=0; cell < num_cell; ++cell) {
1844  for (int qp=0; qp < num_qp; ++qp) {
1845  for (int dim = 0; dim < num_dim; ++dim) {
1846  v(cell,qp,dim).val() = 100.*cell + 10.*qp + 1.*dim;
1847  for (int deriv = 0; deriv < num_deriv; ++deriv) {
1848  v(cell,qp,dim).fastAccessDx(deriv) = v(cell,qp,dim).val() + (1.0*deriv)/10.;
1849  }
1850  }
1851  }
1852  }
1853 
1854  success = true;
1855 
1856  out << "checking subview(v,ALL,*,*)..." << std::endl;
1857  for (int qp=0; qp < num_qp; ++qp) {
1858  for (int dim=0; dim < num_dim; ++dim) {
1859  auto v_tmp = subview(v,Kokkos::ALL(),qp,dim);
1860  for (int cell=0; cell < num_cell; ++cell) {
1861  out << "\tChecking (" << cell << "," << qp << "," << dim << ")" << std::endl;
1862  success = success && checkFads(v(cell,qp,dim), v_tmp(cell), out);
1863  }
1864  }
1865  }
1866 
1867  out << "checking subview(v,*,ALL,*)..." << std::endl;
1868  for (int cell=0; cell < num_cell; ++cell) {
1869  for (int dim=0; dim < num_dim; ++dim) {
1870  auto v_tmp = subview(v,cell,Kokkos::ALL(),dim);
1871  for (int qp=0; qp < num_qp; ++qp) {
1872  out << "\tChecking (" << cell << "," << qp << "," << dim << ")" << std::endl;
1873  success = success && checkFads(v(cell,qp,dim), v_tmp(qp), out);
1874  }
1875  }
1876  }
1877 
1878  out << "checking subview(v,*,*,ALL)..." << std::endl;
1879  for (int cell=0; cell < num_cell; ++cell) {
1880  for (int qp=0; qp < num_qp; ++qp) {
1881  auto v_tmp = subview(v,cell,qp,Kokkos::ALL());
1882  for (int dim=0; dim < num_dim; ++dim) {
1883  out << "\tChecking (" << cell << "," << qp << "," << dim << ")" << std::endl;
1884  success = success && checkFads(v(cell,qp,dim), v_tmp(dim), out);
1885  }
1886  }
1887  }
1888 
1889  out << "checking subview(v,ALL,ALL,*)..." << std::endl;
1890  for (int dim=0; dim < num_dim; ++dim) {
1891  auto v_tmp = subview(v,Kokkos::ALL(),Kokkos::ALL(),dim);
1892  for (int cell=0; cell < num_cell; ++cell) {
1893  for (int qp=0; qp < num_qp; ++qp) {
1894  out << "\tChecking (" << cell << "," << qp << "," << dim << ")" << std::endl;
1895  success = success && checkFads(v(cell,qp,dim), v_tmp(cell,qp), out);
1896  }
1897  }
1898  }
1899 
1900  out << "checking subview(v,*,ALL,ALL)..." << std::endl;
1901  for (int cell=0; cell < num_cell; ++cell) {
1902  auto v_tmp = subview(v,cell,Kokkos::ALL(),Kokkos::ALL());
1903  for (int qp=0; qp < num_qp; ++qp) {
1904  for (int dim=0; dim < num_dim; ++dim) {
1905  out << "\tChecking (" << cell << "," << qp << "," << dim << ")" << std::endl;
1906  success = success && checkFads(v(cell,qp,dim), v_tmp(qp,dim), out);
1907  }
1908  }
1909  }
1910 
1911  out << "checking subview(v,ALL,*,ALL)..." << std::endl;
1912  for (int qp=0; qp < num_qp; ++qp) {
1913  auto v_tmp = subview(v,Kokkos::ALL(),qp,Kokkos::ALL());
1914  for (int cell=0; cell < num_cell; ++cell) {
1915  for (int dim=0; dim < num_dim; ++dim) {
1916  out << "\tChecking (" << cell << "," << qp << "," << dim << ")" << std::endl;
1917  success = success && checkFads(v(cell,qp,dim), v_tmp(cell,dim), out);
1918  }
1919  }
1920  }
1921 
1922  out << "checking subview(v,range,range,range)..." << std::endl;
1923  auto v_tmp = subview(v,std::make_pair(1,5),std::make_pair(1,4),std::make_pair(1,3));
1924  for (int cell=1; cell < num_cell; ++cell) {
1925  for (int qp=1; qp < num_qp; ++qp) {
1926  for (int dim=1; dim < num_dim; ++dim) {
1927  out << "\tChecking (" << cell << "," << qp << "," << dim << ")" << std::endl;
1928  success = success && checkFads(v(cell,qp,dim), v_tmp(cell-1,qp-1,dim-1), out);
1929  }
1930  }
1931  }
1932 }
1933 
1934 #ifdef SACADO_NEW_FAD_DESIGN_IS_DEFAULT
1936  Kokkos_View_Fad, ConstViewAssign, FadType, Layout, Device )
1937 {
1938  typedef Kokkos::View<FadType*,Layout,Device> ViewType;
1939  typedef Kokkos::View<const FadType,Layout,Device> ConstViewType;
1940  typedef typename ViewType::size_type size_type;
1941  typedef typename ViewType::HostMirror host_view_type;
1942  typedef typename ViewType::execution_space exec_space;
1943 
1944  const size_type num_rows = global_num_rows;
1945  const size_type fad_size = global_fad_size;
1946 
1947  // Create and fill view
1948 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1949  ViewType v1("view1", num_rows);
1950 #else
1951  ViewType v1("view1", num_rows, fad_size+1);
1952 #endif
1953  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
1954  for (size_type i=0; i<num_rows; ++i) {
1955  FadType f = generate_fad<FadType>(num_rows, size_type(1), fad_size, i,
1956  size_type(0));
1957  h_v1(i) = f;
1958  }
1959  Kokkos::deep_copy(v1, h_v1);
1960 
1961 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
1962  ViewType v2("view2", num_rows);
1963 #else
1964  ViewType v2("view2", num_rows, fad_size+1);
1965 #endif
1966 
1967  static const size_type stride = Kokkos::ViewScalarStride<ViewType>::stride;
1968 #if defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
1969  const bool use_team =
1973  ( stride > 1 );
1974 #elif defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
1975  const bool use_team =
1980 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
1981  const bool use_team =
1985  ( stride > 1 );
1986 #elif defined (KOKKOS_ENABLE_HIP) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
1987  const bool use_team =
1992 #else
1993  const bool use_team = false;
1994 #endif
1995 
1996  if (use_team) {
1997  typedef Kokkos::TeamPolicy<exec_space> team_policy;
1998  Kokkos::parallel_for(team_policy(num_rows, 1, stride),
1999  KOKKOS_LAMBDA(typename team_policy::member_type team)
2000  {
2001  const int i = team.league_rank();
2002  typename ConstViewType::reference_type x = v1(i);
2003  v2(i) = x;
2004  });
2005  }
2006  else {
2007  Kokkos::parallel_for(Kokkos::RangePolicy<exec_space>(0,num_rows),
2008  KOKKOS_LAMBDA(const int i)
2009  {
2010  typename ConstViewType::reference_type x = v1(i);
2011  v2(i) = x;
2012  });
2013  }
2014 
2015  // Copy back
2016  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
2017  Kokkos::deep_copy(h_v2, v2);
2018 
2019  // Check
2020  success = true;
2021  for (size_type i=0; i<num_rows; ++i) {
2022  FadType f = generate_fad<FadType>(num_rows, size_type(1), fad_size, i,
2023  size_type(0));
2024  success = success && checkFads(f, h_v2(i), out);
2025  }
2026 }
2027 #else
2029  Kokkos_View_Fad, ConstViewAssign, FadType, Layout, Device ) {}
2030 #endif
2031 
2032 // Tests that require view spec
2033 
2034 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
2036  Kokkos_View_Fad, ShmemSize, FadType, Layout, Device )
2037 {
2038  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
2039  typedef typename FadType::value_type value_type;
2040  typedef typename ViewType::size_type size_type;
2041 
2042  const size_type num_rows = global_num_rows;
2043  const size_type num_cols = global_num_cols;
2044  const size_type fad_size = global_fad_size;
2045 
2046  // Compute shared memory size for View
2047  const size_type shmem_size =
2048  ViewType::shmem_size(num_rows, num_cols, fad_size+1);
2049 
2050  // Check
2051  const size_type align = 8;
2052  const size_type mask = align - 1;
2053  ViewType v;
2054 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2055  v = ViewType ("view", num_rows, num_cols);
2056 #else
2057  v = ViewType ("view", num_rows, num_cols, fad_size+1);
2058 #endif
2059  const size_type shmem_size_expected =
2060  (( sizeof(value_type) * global_num_rows * global_num_cols * (fad_size+1) + mask ) & ~mask) + sizeof(typename ViewType::traits::value_type);
2061  TEUCHOS_TEST_EQUALITY(shmem_size, shmem_size_expected, out, success);
2062 }
2063 
2065  Kokkos_View_Fad, Unmanaged, FadType, Layout, Device )
2066 {
2067  // For LayoutContiguous or LayoutNatural, strip out the layout they are templated on
2068  typedef typename Kokkos::inner_layout<Layout>::type TestLayout;
2069 
2070  typedef typename FadType::value_type scalar_type;
2071  typedef Kokkos::View<scalar_type***,TestLayout,Device> ViewType;
2072  typedef Kokkos::View<FadType**,TestLayout,Device,Kokkos::MemoryUnmanaged> FadViewType;
2073  typedef typename ViewType::size_type size_type;
2074  typedef typename ViewType::HostMirror host_view_type;
2075  typedef typename FadViewType::HostMirror fad_host_view_type;
2076 
2077  const size_type num_rows = global_num_rows;
2078  const size_type num_cols = global_num_cols;
2079  const size_type fad_size = global_fad_size;
2080 
2081  // Create and fill view
2082  ViewType v;
2083  host_view_type h_v;
2086  v = ViewType ("view", fad_size+1, num_rows, num_cols);
2087  h_v = Kokkos::create_mirror_view(v);
2088  for (size_type i=0; i<num_rows; ++i) {
2089  for (size_type j=0; j<num_cols; ++j) {
2090  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2091  for (size_type k=0; k<fad_size; k++)
2092  h_v(k,i,j) = f.dx(k);
2093  h_v(fad_size,i,j) = f.val();
2094  }
2095  }
2096  }
2097  else {
2098  v = ViewType ("view", num_rows, num_cols, fad_size+1);
2099  h_v = Kokkos::create_mirror_view(v);
2100  for (size_type i=0; i<num_rows; ++i) {
2101  for (size_type j=0; j<num_cols; ++j) {
2102  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2103  for (size_type k=0; k<fad_size; k++)
2104  h_v(i,j,k) = f.dx(k);
2105  h_v(i,j,fad_size) = f.val();
2106  }
2107  }
2108  }
2109  Kokkos::deep_copy(v, h_v);
2110 
2111  // Create unmanaged view
2112  FadViewType v_fad;//( v.data(), num_rows, num_cols, fad_size+1);
2113  fad_host_view_type h_v_fad;
2114 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2115  v_fad = FadViewType ( v.data(), num_rows, num_cols);
2116  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols);
2117 #else
2118  v_fad = FadViewType ( v.data(), num_rows, num_cols, fad_size+1);
2119  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols, fad_size+1);
2120 #endif
2121 
2122  // Copy back -- can't use create_mirror_view() because v_fad is unmanaged
2123  Kokkos::deep_copy(h_v_fad, v_fad);
2124 
2125  // Check
2126  success = true;
2127  for (size_type i=0; i<num_rows; ++i) {
2128  for (size_type j=0; j<num_cols; ++j) {
2129  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2130  success = success && checkFads(f, h_v_fad(i,j), out);
2131  }
2132  }
2133 }
2134 
2136  Kokkos_View_Fad, Unmanaged2, FadType, Layout, Device )
2137 {
2138  // For LayoutContiguous or LayoutNatural, strip out the layout they are templated on
2139  typedef typename Kokkos::inner_layout<Layout>::type TestLayout;
2140 
2141  typedef typename FadType::value_type scalar_type;
2142  typedef Kokkos::View<scalar_type***,TestLayout,Device> ViewType;
2143  typedef Kokkos::View<FadType**,TestLayout,Device> FadViewType;
2144  typedef typename ViewType::size_type size_type;
2145  typedef typename ViewType::HostMirror host_view_type;
2146  typedef typename FadViewType::HostMirror fad_host_view_type;
2147 
2148  const size_type num_rows = global_num_rows;
2149  const size_type num_cols = global_num_cols;
2150  const size_type fad_size = global_fad_size;
2151 
2152  // Create and fill view
2153  ViewType v;
2154  host_view_type h_v;
2157  v = ViewType ("view", fad_size+1, num_rows, num_cols);
2158  h_v = Kokkos::create_mirror_view(v);
2159  for (size_type i=0; i<num_rows; ++i) {
2160  for (size_type j=0; j<num_cols; ++j) {
2161  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2162  for (size_type k=0; k<fad_size; k++)
2163  h_v(k,i,j) = f.dx(k);
2164  h_v(fad_size,i,j) = f.val();
2165  }
2166  }
2167  }
2168  else {
2169  v = ViewType ("view", num_rows, num_cols, fad_size+1);
2170  h_v = Kokkos::create_mirror_view(v);
2171  for (size_type i=0; i<num_rows; ++i) {
2172  for (size_type j=0; j<num_cols; ++j) {
2173  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2174  for (size_type k=0; k<fad_size; k++)
2175  h_v(i,j,k) = f.dx(k);
2176  h_v(i,j,fad_size) = f.val();
2177  }
2178  }
2179  }
2180  Kokkos::deep_copy(v, h_v);
2181 
2182  // Create unmanaged view
2183  FadViewType v_fad;//( v.data(), num_rows, num_cols, fad_size+1);
2184  fad_host_view_type h_v_fad;
2185 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2186  v_fad = FadViewType ( v.data(), num_rows, num_cols);
2187  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols);
2188 #else
2189  v_fad = FadViewType ( v.data(), num_rows, num_cols, fad_size+1);
2190  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols, fad_size+1);
2191 #endif
2192 
2193  // Copy back -- can't use create_mirror_view() because v_fad is unmanaged
2194  Kokkos::deep_copy(h_v_fad, v_fad);
2195 
2196  // Check
2197  success = true;
2198  for (size_type i=0; i<num_rows; ++i) {
2199  for (size_type j=0; j<num_cols; ++j) {
2200  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2201  success = success && checkFads(f, h_v_fad(i,j), out);
2202  }
2203  }
2204 }
2205 
2207  Kokkos_View_Fad, UnmanagedConst, FadType, Layout, Device )
2208 {
2209  // For LayoutContiguous or LayoutNatural, strip out the layout they are templated on
2210  typedef typename Kokkos::inner_layout<Layout>::type TestLayout;
2211 
2212  typedef typename FadType::value_type scalar_type;
2213  typedef Kokkos::View<scalar_type***,TestLayout,Device> ViewType;
2214  typedef Kokkos::View<const scalar_type***,TestLayout,Device> ConstViewType;
2215  typedef Kokkos::View<FadType**,TestLayout,Device,Kokkos::MemoryUnmanaged> FadViewType;
2216  typedef Kokkos::View<const FadType**,TestLayout,Device,Kokkos::MemoryUnmanaged> ConstFadViewType;
2217  typedef typename ViewType::size_type size_type;
2218  typedef typename ViewType::HostMirror host_view_type;
2219  typedef typename FadViewType::HostMirror fad_host_view_type;
2220 
2221  const size_type num_rows = global_num_rows;
2222  const size_type num_cols = global_num_cols;
2223  const size_type fad_size = global_fad_size;
2224 
2225  // Create and fill view
2226  ViewType v;
2227  host_view_type h_v;
2230  v = ViewType ("view", fad_size+1, num_rows, num_cols);
2231  h_v = Kokkos::create_mirror_view(v);
2232  for (size_type i=0; i<num_rows; ++i) {
2233  for (size_type j=0; j<num_cols; ++j) {
2234  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2235  for (size_type k=0; k<fad_size; k++)
2236  h_v(k,i,j) = f.dx(k);
2237  h_v(fad_size,i,j) = f.val();
2238  }
2239  }
2240  }
2241  else {
2242  v = ViewType ("view", num_rows, num_cols, fad_size+1);
2243  h_v = Kokkos::create_mirror_view(v);
2244  for (size_type i=0; i<num_rows; ++i) {
2245  for (size_type j=0; j<num_cols; ++j) {
2246  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2247  for (size_type k=0; k<fad_size; k++)
2248  h_v(i,j,k) = f.dx(k);
2249  h_v(i,j,fad_size) = f.val();
2250  }
2251  }
2252  }
2253  Kokkos::deep_copy(v, h_v);
2254  ConstViewType v_const = v;
2255 
2256  // Create unmanaged view
2257 
2258  ConstFadViewType v_fad;//( v.data(), num_rows, num_cols, fad_size+1);
2259  fad_host_view_type h_v_fad;
2260 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2261  v_fad = ConstFadViewType ( v_const.data(), num_rows, num_cols);
2262  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols);
2263 #else
2264  v_fad = ConstFadViewType ( v_const.data(), num_rows, num_cols, fad_size+1);
2265  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols, fad_size+1);
2266 #endif
2267 
2268  // Copy back -- can't use create_mirror_view() because v_fad is unmanaged
2269  Kokkos::deep_copy(h_v_fad, v_fad);
2270 
2271  // Check
2272  success = true;
2273  for (size_type i=0; i<num_rows; ++i) {
2274  for (size_type j=0; j<num_cols; ++j) {
2275  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2276  success = success && checkFads(f, h_v_fad(i,j), out);
2277  }
2278  }
2279 }
2280 
2282  Kokkos_View_Fad, UnmanagedConst2, FadType, Layout, Device )
2283 {
2284  // For LayoutContiguous or LayoutNatural, strip out the layout they are templated on
2285  typedef typename Kokkos::inner_layout<Layout>::type TestLayout;
2286  typedef typename FadType::value_type scalar_type;
2287  typedef Kokkos::View<scalar_type***,TestLayout,Device> ViewType;
2288  typedef Kokkos::View<const scalar_type***,TestLayout,Device> ConstViewType;
2289  typedef Kokkos::View<FadType**,TestLayout,Device> FadViewType;
2290  typedef Kokkos::View<const FadType**,TestLayout,Device> ConstFadViewType;
2291  typedef typename ViewType::size_type size_type;
2292  typedef typename ViewType::HostMirror host_view_type;
2293  typedef typename FadViewType::HostMirror fad_host_view_type;
2294 
2295  const size_type num_rows = global_num_rows;
2296  const size_type num_cols = global_num_cols;
2297  const size_type fad_size = global_fad_size;
2298 
2299  // Create and fill view
2300  ViewType v;
2301  host_view_type h_v;
2304  v = ViewType ("view", fad_size+1, num_rows, num_cols);
2305  h_v = Kokkos::create_mirror_view(v);
2306  for (size_type i=0; i<num_rows; ++i) {
2307  for (size_type j=0; j<num_cols; ++j) {
2308  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2309  for (size_type k=0; k<fad_size; k++)
2310  h_v(k,i,j) = f.dx(k);
2311  h_v(fad_size,i,j) = f.val();
2312  }
2313  }
2314  }
2315  else {
2316  v = ViewType ("view", num_rows, num_cols, fad_size+1);
2317  h_v = Kokkos::create_mirror_view(v);
2318  for (size_type i=0; i<num_rows; ++i) {
2319  for (size_type j=0; j<num_cols; ++j) {
2320  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2321  for (size_type k=0; k<fad_size; k++)
2322  h_v(i,j,k) = f.dx(k);
2323  h_v(i,j,fad_size) = f.val();
2324  }
2325  }
2326  }
2327  Kokkos::deep_copy(v, h_v);
2328  ConstViewType v_const = v;
2329 
2330  // Create unmanaged view
2331  ConstFadViewType v_fad;
2332  fad_host_view_type h_v_fad;
2333 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2334  v_fad = ConstFadViewType (v_const.data(), num_rows, num_cols);
2335  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols);
2336 #else
2337  v_fad = ConstFadViewType (v_const.data(), num_rows, num_cols, fad_size+1);
2338  h_v_fad = fad_host_view_type ("host_view_fad", num_rows, num_cols, fad_size+1);
2339 #endif
2340 
2341  // Copy back -- can't use create_mirror_view() because v_fad is unmanaged
2342  Kokkos::deep_copy(h_v_fad, v_fad);
2343 
2344  // Check
2345  success = true;
2346  for (size_type i=0; i<num_rows; ++i) {
2347  for (size_type j=0; j<num_cols; ++j) {
2348  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2349  success = success && checkFads(f, h_v_fad(i,j), out);
2350  }
2351  }
2352 }
2353 
2354 // This test checks we can allocate a view
2355 // with SFad without specifying the fad size in the constructor
2357  Kokkos_View_Fad, SFadNoSizeArg, FadType, Layout, Device )
2358 {
2359  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
2360  typedef typename ViewType::size_type size_type;
2361  typedef typename ViewType::HostMirror host_view_type;
2362 
2363  const size_type num_rows = global_num_rows;
2364  const size_type num_cols = global_num_cols;
2365  const size_type fad_size = global_fad_size;
2366 
2367  // Create and fill view
2368  ViewType v("view", num_rows, num_cols);
2369  host_view_type h_v = Kokkos::create_mirror_view(v);
2370  for (size_type i=0; i<num_rows; ++i) {
2371  for (size_type j=0; j<num_cols; ++j) {
2372  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2373  h_v(i,j) = f;
2374  }
2375  }
2376  Kokkos::deep_copy(v, h_v);
2377 
2378  // Copy back
2379  Kokkos::deep_copy(h_v, v);
2380 
2381  // Check
2382  success = true;
2383  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(v), fad_size+1, out, success);
2384  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_v), fad_size+1, out, success);
2385  for (size_type i=0; i<num_rows; ++i) {
2386  for (size_type j=0; j<num_cols; ++j) {
2387  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2388  success = success && checkFads(f, h_v(i,j), out);
2389  }
2390  }
2391 }
2392 
2394  Kokkos_View_Fad, Partition, FadType, Layout, Device )
2395 {
2396 #if !defined(SACADO_VIEW_CUDA_HIERARCHICAL) && !defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
2397  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
2398  typedef typename ViewType::size_type size_type;
2399  typedef typename ViewType::HostMirror host_view_type;
2400 
2401  const size_type num_rows = global_num_rows;
2402  const size_type num_cols = global_num_cols;
2403  const size_type fad_size = global_fad_size;
2404 
2405  // Create and fill view
2406  ViewType v;
2407 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2408  v = ViewType ("view", num_rows, num_cols);
2409 #else
2410  v = ViewType ("view", num_rows, num_cols, fad_size+1);
2411 #endif
2412  host_view_type h_v = Kokkos::create_mirror_view(v);
2413 
2414  for (size_type i=0; i<num_rows; ++i) {
2415  for (size_type j=0; j<num_cols; ++j) {
2416  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2417  h_v(i,j) = f;
2418  }
2419  }
2420  Kokkos::deep_copy(v, h_v);
2421 
2422  // Copy back
2423  Kokkos::deep_copy(h_v, v);
2424 
2425  // Partition derivative array of h_v into 2, first one starting at index 0,
2426  // the second at 1
2427  const size_type stride = 2;
2428  auto h_v1 = Kokkos::partition<2>(h_v, 0, stride);
2429  auto h_v2 = Kokkos::partition<2>(h_v, 1, stride);
2430 
2431  // Check
2432  const size_type fad_size_1 = (fad_size + stride - 0 - 1) / stride;
2433  const size_type fad_size_2 = (fad_size + stride - 1 - 1) / stride;
2434  success = true;
2435  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_v1), fad_size_1+1, out, success);
2436  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_v2), fad_size_2+1, out, success);
2437  for (size_type i=0; i<num_rows; ++i) {
2438  for (size_type j=0; j<num_cols; ++j) {
2439  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2440  Sacado::Fad::DFad<double> f1( fad_size_1, f.val() );
2441  Sacado::Fad::DFad<double> f2( fad_size_2, f.val() );
2442  for (unsigned int k=0; k<fad_size_1; ++k)
2443  if (2*k < fad_size) f1.fastAccessDx(k) = f.dx(2*k);
2444  for (unsigned int k=0; k<fad_size_2; ++k)
2445  if (2*k+1 < fad_size) f2.fastAccessDx(k) = f.dx(2*k+1);
2446  success = success && checkFads(f1, h_v1(i,j), out);
2447  success = success && checkFads(f2, h_v2(i,j), out);
2448  }
2449  }
2450 #endif
2451 }
2452 
2454  Kokkos_View_Fad, AssignLayoutContiguousToLayoutStride, FadType, Layout, Device )
2455 {
2456  typedef Kokkos::View<FadType**,Kokkos::LayoutContiguous<Layout>,Device> ContViewType;
2457  typedef Kokkos::View<FadType**,Kokkos::LayoutStride,Device> StrideViewType;
2458  typedef typename ContViewType::size_type size_type;
2459  typedef typename ContViewType::HostMirror cont_host_view_type;
2460  typedef typename StrideViewType::HostMirror stride_host_view_type;
2461 
2462  const size_type num_rows = global_num_rows;
2463  const size_type num_cols = global_num_cols;
2464  const size_type fad_size = global_fad_size;
2465 
2466  // Create and fill view
2467  ContViewType v;
2468 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2469  v = ContViewType ("view", num_rows, num_cols);
2470 #else
2471  v = ContViewType ("view", num_rows, num_cols, fad_size+1);
2472 #endif
2473  cont_host_view_type h_v = Kokkos::create_mirror_view(v);
2474 
2475  for (size_type i=0; i<num_rows; ++i) {
2476  for (size_type j=0; j<num_cols; ++j) {
2477  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2478  h_v(i,j) = f;
2479  }
2480  }
2481  Kokkos::deep_copy(v, h_v);
2482 
2483  // Assign to LayoutStride view
2484  StrideViewType vs = v;
2485 
2486  // Copy back
2487  // Note: don't create h_vs through create_mirror_view and deep_copy
2488  // since Kokkos doesn't support deep_copy of non-contiguous views
2489  //stride_host_view_type h_vs = Kokkos::create_mirror_view(vs);
2490  //Kokkos::deep_copy(h_vs, vs);
2491  stride_host_view_type h_vs = h_v;
2492 
2493  // Check
2494  success = true;
2495  TEUCHOS_TEST_EQUALITY(h_vs.extent(0), num_rows, out, success);
2496  TEUCHOS_TEST_EQUALITY(h_vs.extent(1), num_cols, out, success);
2497  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(h_vs), fad_size+1, out, success);
2498  for (size_type i=0; i<num_rows; ++i) {
2499  for (size_type j=0; j<num_cols; ++j) {
2500  FadType f = generate_fad<FadType>(num_rows, num_cols, fad_size, i, j);
2501  success = success && checkFads(f, h_vs(i,j), out);
2502  }
2503  }
2504 }
2505 
2507  Kokkos_View_Fad, CommonViewAllocMixedSpec, FadType, Layout, Device )
2508 {
2509  typedef Kokkos::View<FadType**,Kokkos::LayoutContiguous<Layout>,Device> ContViewType;
2510  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
2511  typedef typename ContViewType::size_type size_type;
2512 
2513  const size_type num_rows = global_num_rows;
2514  const size_type num_cols = global_num_cols;
2515  const size_type fad_size = global_fad_size;
2516 
2517  // Create contiguous view
2518  ContViewType v1;
2519 #if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
2520  v1 = ContViewType ("view", num_rows, num_cols);
2521 #else
2522  v1 = ContViewType ("view", num_rows, num_cols, fad_size+1);
2523 #endif
2524 
2525  // Create non-contiguous view using commen_view_alloc_prop
2526  auto cprop = Kokkos::common_view_alloc_prop(v1);
2527  ViewType v2(Kokkos::view_alloc("v2",cprop), num_rows, num_cols);
2528 
2529  // Check dimensions are correct for v2
2530  success = true;
2531  TEUCHOS_TEST_EQUALITY(v2.extent(0), num_rows, out, success);
2532  TEUCHOS_TEST_EQUALITY(v2.extent(1), num_cols, out, success);
2533  TEUCHOS_TEST_EQUALITY(Kokkos::dimension_scalar(v2), fad_size+1, out, success);
2534 }
2535 
2536 #else
2537 
2539  Kokkos_View_Fad, ShmemSize, FadType, Layout, Device )
2540 {
2541  typedef Kokkos::View<FadType**,Layout,Device> ViewType;
2542  typedef typename ViewType::size_type size_type;
2543 
2544  const size_type num_rows = global_num_rows;
2545  const size_type num_cols = global_num_cols;
2546 
2547  // Compute shared memory size for View
2548  const size_type shmem_size =
2549  ViewType::shmem_size(num_rows, num_cols);
2550 
2551  // Check
2552  static const size_type align = 8;
2553  static const size_type mask = align - 1;
2554  const size_type shmem_size_expected =
2555  (( sizeof(FadType) * global_num_rows * global_num_cols + mask ) & ~mask) + sizeof(typename ViewType::traits::value_type);
2556  TEUCHOS_TEST_EQUALITY(shmem_size, shmem_size_expected, out, success);
2557 }
2558 
2560  Kokkos_View_Fad, Unmanaged, FadType, Layout, Device ) {}
2561 
2563  Kokkos_View_Fad, Unmanaged2, FadType, Layout, Device ) {}
2564 
2566  Kokkos_View_Fad, UnmanagedConst, FadType, Layout, Device ) {}
2567 
2569  Kokkos_View_Fad, UnmanagedConst2, FadType, Layout, Device ) {}
2570 
2572  Kokkos_View_Fad, SFadNoSizeArg, FadType, Layout, Device ) {}
2573 
2575  Kokkos_View_Fad, Partition, FadType, Layout, Device ) {}
2576 
2578  Kokkos_View_Fad, AssignLayoutContiguousToLayoutStride, FadType, Layout, Device ) {}
2579 
2581  Kokkos_View_Fad, CommonViewAllocMixedSpec, FadType, Layout, Device ) {}
2582 
2583 #endif
2584 
2585 #define VIEW_FAD_TESTS_FLD( F, L, D ) \
2586  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Size, F, L, D ) \
2587  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DeepCopy, F, L, D ) \
2588  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DeepCopy_ConstantScalar, F, L, D ) \
2589  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DeepCopy_ConstantZero, F, L, D ) \
2590  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DeepCopy_ConstantFad, F, L, D ) \
2591  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DeepCopy_ConstantFadFull, F, L, D ) \
2592  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, LocalDeepCopy, F, L, D ) \
2593  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, LocalDeepCopyTeam, F, L, D ) \
2594  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, ScalarAssign, F, L, D ) \
2595  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, ValueAssign, F, L, D ) \
2596  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Resize, F, L, D ) \
2597  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Unmanaged, F, L, D ) \
2598  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Unmanaged2, F, L, D ) \
2599  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, UnmanagedConst, F, L, D ) \
2600  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, UnmanagedConst2, F, L, D ) \
2601  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Multiply, F, L, D ) \
2602  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, MultiplyUpdate, F, L, D ) \
2603  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, MultiplyConst, F, L, D ) \
2604  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, MultiplyMixed, F, L, D ) \
2605  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Rank8, F, L, D ) \
2606  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Roger, F, L, D ) \
2607  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, AtomicAdd, F, L, D ) \
2608  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, AssignDifferentStrides, F, L, D ) \
2609  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, ScalarValue, F, L, D ) \
2610  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DynRankDimensionScalar, F, L, D ) \
2611  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DynRankAssignStatic0, F, L, D ) \
2612  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DynRankAssignStatic1, F, L, D ) \
2613  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DynRankAssignStatic2, F, L, D ) \
2614  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DynRankMultiply, F, L, D ) \
2615  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, SubdynrankviewCol, F, L, D ) \
2616  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, SubdynrankviewRow, F, L, D ) \
2617  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, SubdynrankviewScalar, F, L, D ) \
2618  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Subview, F, L, D ) \
2619  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Subview2, F, L, D ) \
2620  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, ShmemSize, F, L, D ) \
2621  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, ConstViewAssign, F, L, D )
2622 
2623 #define VIEW_FAD_TESTS_SFLD( F, L, D ) \
2624  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, SFadNoSizeArg, F, L, D )
2625 
2626 #define VIEW_FAD_TESTS_FDI( F, D ) \
2627  using Kokkos::LayoutLeft; \
2628  using Kokkos::LayoutRight; \
2629  VIEW_FAD_TESTS_FLD( F, LayoutLeft, D ) \
2630  VIEW_FAD_TESTS_FLD( F, LayoutRight, D ) \
2631  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, AssignLayoutContiguousToLayoutStride, F, LayoutLeft, D ) \
2632  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, AssignLayoutContiguousToLayoutStride, F, LayoutRight, D ) \
2633  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, CommonViewAllocMixedSpec, F, LayoutLeft, D ) \
2634  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, CommonViewAllocMixedSpec, F, LayoutRight, D )
2635 
2636 #define VIEW_FAD_TESTS_SFDI( F, D ) \
2637  using Kokkos::LayoutLeft; \
2638  using Kokkos::LayoutRight; \
2639  VIEW_FAD_TESTS_SFLD( F, LayoutLeft, D ) \
2640  VIEW_FAD_TESTS_SFLD( F, LayoutRight, D )
2641 
2642 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
2645 #define VIEW_FAD_TESTS_FDC( F, D ) \
2646  VIEW_FAD_TESTS_FLD( F, LeftContiguous, D ) \
2647  VIEW_FAD_TESTS_FLD( F, RightContiguous, D ) \
2648  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Partition, F, LeftContiguous, D ) \
2649  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Partition, F, RightContiguous, D )
2650 
2651 #define VIEW_FAD_TESTS_SFDC( F, D ) \
2652  VIEW_FAD_TESTS_SFLD( F, LeftContiguous, D ) \
2653  VIEW_FAD_TESTS_SFLD( F, RightContiguous, D )
2654 #else
2655 #define VIEW_FAD_TESTS_FDC( F, D ) /* */
2656 #define VIEW_FAD_TESTS_SFDC( F, D ) /* */
2657 #endif
2658 
2659 #define VIEW_FAD_TESTS_FD( F, D ) \
2660  VIEW_FAD_TESTS_FDI( F, D ) \
2661  VIEW_FAD_TESTS_FDC( F, D )
2662 
2663 #define VIEW_FAD_TESTS_SFD( F, D ) \
2664  VIEW_FAD_TESTS_SFDI( F, D ) \
2665  VIEW_FAD_TESTS_SFDC( F, D )
2666 
2667 // We've unified the implementation for the different Fad variants, so
2668 // there is no reason to test ELRFad, CacheFad, and ELRCacheFad.
2672 
2673 /*
2674 typedef Sacado::ELRFad::DFad<double> ELRDFadType;
2675 typedef Sacado::ELRFad::SLFad<double,2*global_fad_size> ELRSLFadType;
2676 typedef Sacado::ELRFad::SFad<double,global_fad_size> ELRSFadType;
2677 
2678 typedef Sacado::CacheFad::DFad<double> CacheDFadType;
2679 typedef Sacado::CacheFad::SLFad<double,2*global_fad_size> CacheSLFadType;
2680 typedef Sacado::CacheFad::SFad<double,global_fad_size> CacheSFadType;
2681 
2682 typedef Sacado::ELRCacheFad::DFad<double> ELRCacheDFadType;
2683 typedef Sacado::ELRCacheFad::SLFad<double,2*global_fad_size> ELRCacheSLFadType;
2684 typedef Sacado::ELRCacheFad::SFad<double,global_fad_size> ELRCacheSFadType;
2685 */
2686 
2687 // We can't use DFad unless we use the View specialization
2688 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC) && SACADO_TEST_DFAD
2689 #define VIEW_FAD_TESTS_D( D ) \
2690  VIEW_FAD_TESTS_FD( SFadType, D ) \
2691  VIEW_FAD_TESTS_FD( SLFadType, D ) \
2692  VIEW_FAD_TESTS_FD( DFadType, D ) \
2693  VIEW_FAD_TESTS_SFD( SFadType, D )
2694 
2695 #if 0
2696  VIEW_FAD_TESTS_FD( ELRSFadType, D ) \
2697  VIEW_FAD_TESTS_FD( ELRSLFadType, D ) \
2698  VIEW_FAD_TESTS_FD( ELRDFadType, D ) \
2699  VIEW_FAD_TESTS_FD( CacheSFadType, D ) \
2700  VIEW_FAD_TESTS_FD( CacheSLFadType, D ) \
2701  VIEW_FAD_TESTS_FD( CacheDFadType, D ) \
2702  VIEW_FAD_TESTS_FD( ELRCacheSFadType, D ) \
2703  VIEW_FAD_TESTS_FD( ELRCacheSLFadType, D ) \
2704  VIEW_FAD_TESTS_FD( ELRCacheDFadType, D ) \
2705  VIEW_FAD_TESTS_SFD( SFadType, D ) \
2706  VIEW_FAD_TESTS_SFD( ELRSFadType, D ) \
2707  VIEW_FAD_TESTS_SFD( CacheSFadType, D ) \
2708  VIEW_FAD_TESTS_SFD( ELRCacheSFadType, D )
2709 #endif
2710 
2711 #else
2712 
2713 #define VIEW_FAD_TESTS_D( D ) \
2714  VIEW_FAD_TESTS_FD( SFadType, D ) \
2715  VIEW_FAD_TESTS_FD( SLFadType, D ) \
2716  VIEW_FAD_TESTS_SFD( SFadType, D )
2717 
2718 #if 0
2719  VIEW_FAD_TESTS_FD( ELRSFadType, D ) \
2720  VIEW_FAD_TESTS_FD( ELRSLFadType, D ) \
2721  VIEW_FAD_TESTS_FD( CacheSFadType, D ) \
2722  VIEW_FAD_TESTS_FD( CacheSLFadType, D ) \
2723  VIEW_FAD_TESTS_FD( ELRCacheSFadType, D ) \
2724  VIEW_FAD_TESTS_FD( ELRCacheSLFadType, D ) \
2725  VIEW_FAD_TESTS_SFD( SFadType, D ) \
2726  VIEW_FAD_TESTS_SFD( ELRSFadType, D ) \
2727  VIEW_FAD_TESTS_SFD( CacheSFadType, D ) \
2728  VIEW_FAD_TESTS_SFD( ELRCacheSFadType, D )
2729 #endif
2730 
2731 #endif
SACADO_INLINE_FUNCTION ScalarType< T >::type scalarValue(const T &x)
A simple template function for invoking ScalarValue&lt;&gt;
InputViewType::size_type size_type
team_policy_type::member_type team_handle
Kokkos::TeamPolicy< execution_space > team_policy_type
static const size_type stride
ViewType::value_type::value_type ScalarType
ViewType::execution_space execution_space
Kokkos::TeamPolicy< execution_space > team_policy_type
static const size_type stride
ViewType::execution_space execution_space
KOKKOS_INLINE_FUNCTION void operator()(const team_handle &team) const
void f()
ViewType::size_type size_type
const ViewType m_v
Kokkos::LayoutContiguous< Kokkos::LayoutRight > RightContiguous
Kokkos::RangePolicy< execution_space > range_policy_type
ViewType::value_type ValueType
Kokkos::LayoutContiguous< Kokkos::LayoutLeft > LeftContiguous
static const size_type stride
#define TEUCHOS_TEST_FLOATING_EQUALITY(v1, v2, tol, out, success)
KOKKOS_INLINE_FUNCTION void operator()(const team_handle &team) const
static const bool value
const int global_fad_size
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
const ScalarViewType m_s
#define VIEW_FAD_TESTS_SFD(F, D)
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
Sacado::Fad::DFad< double > FadType
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
#define VIEW_FAD_TESTS_FD(F, D)
expr true
const InputViewType1 m_v1
static const size_type stride
ViewType::size_type size_type
InputViewType::execution_space execution_space
static void apply(const InputViewType v1, const OutputViewType v2, const size_type col)
team_policy_type::member_type team_handle
bool checkFads(const FadType1 &x, const FadType2 &x2, Teuchos::FancyOStream &out, double tol=1.0e-15)
const ScalarViewType m_s
const int global_num_rows
Kokkos::RangePolicy< execution_space > range_policy_type
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
ViewType::execution_space execution_space
const ScalarType m_s
const ViewType m_v
scalar generate_fad(const size_t n0, const size_t n1, const size_t n2, const size_t n3, const int fad_size, const size_t i0, const size_t i1, const size_t i2, const size_t i3, const int i_fad)
Sacado::Fad::SFad< double, fad_dim > SFadType
const OutputViewType m_v2
expr val()
TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(Kokkos_View_FadFad, DeepCopy, FadFadType, Layout, Device)
AssignRank2Rank1Kernel(const InputViewType v1, const OutputViewType v2, const size_type col)
expr expr1 expr1 expr1 c expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr1 c expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr1 c *expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr1 c expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr2 expr1 expr2 expr1 expr1 expr1 c
team_policy_type::member_type team_handle
const int global_num_cols
KOKKOS_INLINE_FUNCTION void operator()(const team_handle &team) const
ValueAssignKernel(const ViewType &v, const ScalarViewType &s)
ViewType::size_type size_type
GeneralFad< DynamicStorage< T > > DFad
team_policy_type::member_type team_handle
static void apply(const InputViewType1 v1, const InputViewType2 v2, const OutputViewType v3, const bool update=false)
InputViewType1::execution_space execution_space
#define GLOBAL_FAD_SIZE
Kokkos::TeamPolicy< execution_space > team_policy_type
static void apply(const ViewType &v, const ScalarViewType &s)
Sacado::Fad::SLFad< double, fad_dim > SLFadType
void g()
Sacado::Fad::DFad< double > DFadType
static void apply(const ViewType &v, const ScalarViewType &s)
Kokkos::TeamPolicy< execution_space > team_policy_type
ScalarAssignKernel(const ViewType &v, const ScalarType &s)
KOKKOS_INLINE_FUNCTION void operator()(const team_handle &team) const
static const bool value
int value
AtomicAddKernel(const ViewType &v, const ScalarViewType &s)
expr expr expr fastAccessDx(i)) FAD_UNARYOP_MACRO(exp
#define TEUCHOS_TEST_EQUALITY(v1, v2, out, success)
KOKKOS_INLINE_FUNCTION void operator()(const team_handle &team) const
team_policy_type::member_type team_handle
Kokkos::ThreadLocalScalarType< ViewType >::type local_scalar_type
const double tol
static void apply(const ViewType &v, const ScalarType &s)
Kokkos::ThreadLocalScalarType< ViewType >::type local_scalar_type
MultiplyKernel(const InputViewType1 v1, const InputViewType2 v2, const OutputViewType v3, const bool update)
const InputViewType m_v1
InputViewType1::size_type size_type
Kokkos::RangePolicy< execution_space > range_policy_type
Kokkos::RangePolicy< execution_space > range_policy_type
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
GeneralFad< StaticFixedStorage< T, Num > > SFad
const InputViewType2 m_v2
Kokkos::RangePolicy< execution_space > range_policy_type
Kokkos::TeamPolicy< execution_space > team_policy_type
const OutputViewType m_v3