Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Fad_Fad_KokkosTests.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Sacado Package
4 //
5 // Copyright 2006 NTESS and the Sacado contributors.
6 // SPDX-License-Identifier: LGPL-2.1-or-later
7 // *****************************************************************************
8 // @HEADER
9 
11 
12 #include "Sacado.hpp"
14 
15 template <typename FadType1, typename FadType2>
16 bool checkFads(const FadType1& x, const FadType2& x2,
17  Teuchos::FancyOStream& out, double tol = 1.0e-15)
18 {
19  bool success = true;
20 
21  // Check sizes match
22  TEUCHOS_TEST_EQUALITY(x.size(), x2.size(), out, success);
23 
24  // Check values match
25  TEUCHOS_TEST_EQUALITY(x.val(), x2.val(), out, success);
26 
27  // Check derivatives match
28  for (int i=0; i<x.size(); ++i)
29  TEUCHOS_TEST_FLOATING_EQUALITY(x.dx(i), x2.dx(i), tol, out, success);
30 
31  return success;
32 }
33 
34 template <typename FadType1, typename FadType2>
35 bool checkNestedFads(const FadType1& x, const FadType2& x2,
36  Teuchos::FancyOStream& out, double tol = 1.0e-15)
37 {
38  bool success = true;
39 
40  // Check sizes match
41  TEUCHOS_TEST_EQUALITY(x.size(), x2.size(), out, success);
42 
43  // Check values match
44  success = success && checkFads(x.val(), x2.val(), out, tol);
45 
46  // Check derivatives match
47  for (int i=0; i<x.size(); ++i)
48  success = success && checkFads(x.dx(i), x2.dx(i), out, tol);
49 
50  return success;
51 }
52 
53 template <typename fadfadtype, typename ordinal>
54 inline
55 fadfadtype generate_nested_fad( const ordinal num_rows,
56  const ordinal num_cols,
57  const ordinal outer_fad_size,
58  const ordinal inner_fad_size,
59  const ordinal row,
60  const ordinal col )
61 {
62  typedef typename fadfadtype::value_type fadtype;
63  typedef typename fadtype::value_type scalar;
64  fadfadtype x(outer_fad_size, scalar(0.0));
65  fadtype y(inner_fad_size, scalar(0.0));
66 
67  const scalar x_row = 1000.0 + scalar(num_rows) / scalar(row+1);
68  const scalar x_col = 100.0 + scalar(num_cols) / scalar(col+1);
69  y.val() = x_row + x_col;
70  for (ordinal j=0; j<inner_fad_size; ++j) {
71  const scalar y_fad = 1.0 + scalar(inner_fad_size) / scalar(j+1);
72  y.fastAccessDx(j) = x_row + x_col + y_fad;
73  }
74  x.val() = y;
75  for (ordinal i=0; i<outer_fad_size; ++i) {
76  const scalar x_fad = 10.0 + scalar(outer_fad_size) / scalar(i+1);
77  y.val() = x_fad;
78  for (ordinal j=0; j<inner_fad_size; ++j) {
79  const scalar y_fad = 1.0 + scalar(inner_fad_size) / scalar(j+1);
80  y.fastAccessDx(j) = x_row + x_col + x_fad + y_fad;
81  }
82  x.fastAccessDx(i) = y;
83  }
84  return x;
85 }
86 
87 const int global_num_rows = 11;
88 const int global_num_cols = 7;
89 const int global_outer_fad_size = 5;
90 const int global_inner_fad_size = 3;
91 
93  Kokkos_View_FadFad, DeepCopy, FadFadType, Layout, Device )
94 {
95  typedef Kokkos::View<FadFadType**,Layout,Device> ViewType;
96  typedef typename ViewType::size_type size_type;
97  typedef typename ViewType::HostMirror host_view_type;
98 
99  const size_type num_rows = global_num_rows;
100  const size_type num_cols = global_num_cols;
101  const size_type outer_fad_size = global_outer_fad_size;
102  const size_type inner_fad_size = global_inner_fad_size;
103 
104  // Create and fill view
105  ViewType v1("view1", num_rows, num_cols, outer_fad_size+1);
106  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
107  for (size_type i=0; i<num_rows; ++i)
108  for (size_type j=0; j<num_cols; ++j)
109  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
110  num_cols,
111  outer_fad_size,
112  inner_fad_size,
113  i, j);
114  Kokkos::deep_copy(v1, h_v1);
115 
116  // Deep copy
117  ViewType v2("view2", num_rows, num_cols, outer_fad_size+1);
118  Kokkos::deep_copy(v2, v1);
119 
120  // Copy back
121  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
122  Kokkos::deep_copy(h_v2, v2);
123 
124  // Check
125  success = true;
126  for (size_type i=0; i<num_rows; ++i) {
127  for (size_type j=0; j<num_cols; ++j) {
128  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
129  num_cols,
130  outer_fad_size,
131  inner_fad_size,
132  i, j);
133  success = success && checkNestedFads(f, h_v2(i,j), out);
134  }
135  }
136 }
137 
138 #ifdef HAVE_SACADO_KOKKOS
139 
141  Kokkos_DynRankView_FadFad, DeepCopy, FadFadType, Layout, Device )
142 {
143  typedef Kokkos::DynRankView<FadFadType,Layout,Device> ViewType;
144  typedef typename ViewType::size_type size_type;
145  typedef typename ViewType::HostMirror host_view_type;
146 
147  const size_type num_rows = global_num_rows;
148  const size_type num_cols = global_num_cols;
149  const size_type outer_fad_size = global_outer_fad_size;
150  const size_type inner_fad_size = global_inner_fad_size;
151 
152  // Create and fill view
153  ViewType v1("view1", num_rows, num_cols, outer_fad_size+1);
154  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
155  for (size_type i=0; i<num_rows; ++i)
156  for (size_type j=0; j<num_cols; ++j)
157  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
158  num_cols,
159  outer_fad_size,
160  inner_fad_size,
161  i, j);
162  Kokkos::deep_copy(v1, h_v1);
163 
164  // Deep copy
165  ViewType v2("view2", num_rows, num_cols, outer_fad_size+1);
166  Kokkos::deep_copy(v2, v1);
167 
168  // Copy back
169  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
170  Kokkos::deep_copy(h_v2, v2);
171 
172  // Check
173  success = true;
174  for (size_type i=0; i<num_rows; ++i) {
175  for (size_type j=0; j<num_cols; ++j) {
176  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
177  num_cols,
178  outer_fad_size,
179  inner_fad_size,
180  i, j);
181  success = success && checkNestedFads(f, h_v2(i,j), out);
182  }
183  }
184 }
185 
186 // To test DynRankView - View interoperabitlity
187 // Deep copy of DynRankView to View
188 // Assignment operator
190  Kokkos_DynRankView_FadFad, Interop, FadFadType, Layout, Device )
191 {
192  typedef Kokkos::DynRankView<FadFadType,Layout,Device> DRViewType;
193  typedef typename DRViewType::size_type size_type;
194  typedef typename DRViewType::HostMirror host_view_type;
195 
196  typedef Kokkos::View<FadFadType**,Layout,Device> NoDynViewType;
197  typedef typename NoDynViewType::HostMirror host_nondynrankview_type;
198 
199  const size_type num_rows = global_num_rows;
200  const size_type num_cols = global_num_cols;
201  const size_type outer_fad_size = global_outer_fad_size;
202  const size_type inner_fad_size = global_inner_fad_size;
203 
204  // Create and fill view
205  DRViewType v1("drview1", num_rows, num_cols, outer_fad_size+1);
206  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
207 
208  NoDynViewType ndv2("nodview2", num_rows, num_cols, outer_fad_size+1);
209  host_nondynrankview_type h_ndv2 = Kokkos::create_mirror_view(ndv2);
210 
211  for (size_type i=0; i<num_rows; ++i)
212  for (size_type j=0; j<num_cols; ++j)
213  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
214  num_cols,
215  outer_fad_size,
216  inner_fad_size,
217  i, j);
218  Kokkos::deep_copy(v1, h_v1); //v1 unused here
219 
220  // Deep copy DynRankView to View on device
221  Kokkos::deep_copy(ndv2, h_v1);
222  // Assign View to DynRankView
223  DRViewType v2("drview2", num_rows, num_cols, outer_fad_size+1);
224  v2 = ndv2 ;
225 
226  // Copy back
227  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
228  Kokkos::deep_copy(h_v2, v2);
229 
230  // Check
231  success = true;
232  for (size_type i=0; i<num_rows; ++i) {
233  for (size_type j=0; j<num_cols; ++j) {
234  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
235  num_cols,
236  outer_fad_size,
237  inner_fad_size,
238  i, j);
239  success = success && checkNestedFads(f, h_v2(i,j), out);
240  }
241  }
242 }
243 
244 
245 // Deep copy of DynRankView to View
246 // Copy ctor
248  Kokkos_DynRankView_FadFad, Interop2, FadFadType, Layout, Device )
249 {
250  typedef Kokkos::DynRankView<FadFadType,Layout,Device> DRViewType;
251  typedef typename DRViewType::size_type size_type;
252  typedef typename DRViewType::HostMirror host_view_type;
253 
254  typedef Kokkos::View<FadFadType**,Layout,Device> NoDynViewType;
255  typedef typename NoDynViewType::HostMirror host_nondynrankview_type;
256 
257  const size_type num_rows = global_num_rows;
258  const size_type num_cols = global_num_cols;
259  const size_type outer_fad_size = global_outer_fad_size;
260  const size_type inner_fad_size = global_inner_fad_size;
261 
262  // Create and fill view
263  DRViewType v1("drview1", num_rows, num_cols, outer_fad_size+1);
264  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
265 
266  NoDynViewType ndv2("nodview2", num_rows, num_cols, outer_fad_size+1);
267  host_nondynrankview_type h_ndv2 = Kokkos::create_mirror_view(ndv2);
268 
269  for (size_type i=0; i<num_rows; ++i)
270  for (size_type j=0; j<num_cols; ++j)
271  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
272  num_cols,
273  outer_fad_size,
274  inner_fad_size,
275  i, j);
276  Kokkos::deep_copy(v1, h_v1); //v1 unused here
277 
278  // Deep copy DynRankView to View on device
279  Kokkos::deep_copy(ndv2, h_v1);
280  // Copy construct DynRankView from View
281  DRViewType v2(ndv2) ;
282 
283  // Copy back
284  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
285  Kokkos::deep_copy(h_v2, v2);
286 
287  // Check
288  success = true;
289  for (size_type i=0; i<num_rows; ++i) {
290  for (size_type j=0; j<num_cols; ++j) {
291  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
292  num_cols,
293  outer_fad_size,
294  inner_fad_size,
295  i, j);
296  success = success && checkNestedFads(f, h_v2(i,j), out);
297  }
298  }
299 }
300 
301 
302 
303 #else
304 
306  Kokkos_DynRankView_FadFad, DeepCopy, FadFadType, Layout, Device ) {}
308  Kokkos_DynRankView_FadFad, Interop, FadFadType, Layout, Device ) {}
310  Kokkos_DynRankView_FadFad, Interop2, FadFadType, Layout, Device ) {}
311 
312 #endif
313 
314 #define VIEW_FAD_TESTS_FLD( F, L, D ) \
315  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_FadFad, DeepCopy, F, L, D ) \
316  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_DynRankView_FadFad, DeepCopy, F, L, D ) \
317  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_DynRankView_FadFad, Interop, F, L, D ) \
318  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_DynRankView_FadFad, Interop2, F, L, D )
319 
320 #define VIEW_FAD_TESTS_FD( F, D ) \
321  using Kokkos::LayoutLeft; \
322  using Kokkos::LayoutRight; \
323  VIEW_FAD_TESTS_FLD( F, LayoutLeft, D) \
324  VIEW_FAD_TESTS_FLD( F, LayoutRight, D)
325 
326 // We've unified the implementation for the different Fad variants, so
327 // there is no reason to test ELRFad, CacheFad, and ELRCacheFad.
332 
333 // These tests are only relevant when we have the experimental view
334 // specialization
335 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
336 #if SACADO_TEST_DFAD
337 #define VIEW_FAD_TESTS_D( D ) \
338  VIEW_FAD_TESTS_FD( SFadType, D ) \
339  VIEW_FAD_TESTS_FD( SLFadType, D ) \
340  VIEW_FAD_TESTS_FD( DFadType, D )
341 #else
342 #define VIEW_FAD_TESTS_D( D ) \
343  VIEW_FAD_TESTS_FD( SFadType, D ) \
344  VIEW_FAD_TESTS_FD( SLFadType, D )
345 #endif
346 #else
347 #define VIEW_FAD_TESTS_D( D ) /* */
348 #endif
void f()
#define TEUCHOS_TEST_FLOATING_EQUALITY(v1, v2, tol, out, success)
bool checkFads(const FadType1 &x, const FadType2 &x2, Teuchos::FancyOStream &out, double tol=1.0e-15)
const int global_num_rows
Sacado::Fad::SFad< double, fad_dim > SFadType
TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(Kokkos_View_FadFad, DeepCopy, FadFadType, Layout, Device)
const int global_num_cols
Sacado::Fad::SLFad< double, fad_dim > SLFadType
Sacado::Fad::DFad< double > DFadType
const int global_inner_fad_size
#define TEUCHOS_TEST_EQUALITY(v1, v2, out, success)
Sacado::Fad::SFad< double, global_inner_fad_size > InnerFadType
const double tol
bool checkNestedFads(const FadType1 &x, const FadType2 &x2, Teuchos::FancyOStream &out, double tol=1.0e-15)
const int global_outer_fad_size
fadfadtype generate_nested_fad(const ordinal num_rows, const ordinal num_cols, const ordinal outer_fad_size, const ordinal inner_fad_size, const ordinal row, const ordinal col)
const double y