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 //
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 FadType1, typename FadType2>
35 bool checkFads(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  TEUCHOS_TEST_EQUALITY(x.val(), x2.val(), out, success);
45 
46  // Check derivatives match
47  for (int i=0; i<x.size(); ++i)
48  TEUCHOS_TEST_FLOATING_EQUALITY(x.dx(i), x2.dx(i), tol, out, success);
49 
50  return success;
51 }
52 
53 template <typename FadType1, typename FadType2>
54 bool checkNestedFads(const FadType1& x, const FadType2& x2,
55  Teuchos::FancyOStream& out, double tol = 1.0e-15)
56 {
57  bool success = true;
58 
59  // Check sizes match
60  TEUCHOS_TEST_EQUALITY(x.size(), x2.size(), out, success);
61 
62  // Check values match
63  success = success && checkFads(x.val(), x2.val(), out, tol);
64 
65  // Check derivatives match
66  for (int i=0; i<x.size(); ++i)
67  success = success && checkFads(x.dx(i), x2.dx(i), out, tol);
68 
69  return success;
70 }
71 
72 template <typename fadfadtype, typename ordinal>
73 inline
74 fadfadtype generate_nested_fad( const ordinal num_rows,
75  const ordinal num_cols,
76  const ordinal outer_fad_size,
77  const ordinal inner_fad_size,
78  const ordinal row,
79  const ordinal col )
80 {
81  typedef typename fadfadtype::value_type fadtype;
82  typedef typename fadtype::value_type scalar;
83  fadfadtype x(outer_fad_size, scalar(0.0));
84  fadtype y(inner_fad_size, scalar(0.0));
85 
86  const scalar x_row = 1000.0 + scalar(num_rows) / scalar(row+1);
87  const scalar x_col = 100.0 + scalar(num_cols) / scalar(col+1);
88  y.val() = x_row + x_col;
89  for (ordinal j=0; j<inner_fad_size; ++j) {
90  const scalar y_fad = 1.0 + scalar(inner_fad_size) / scalar(j+1);
91  y.fastAccessDx(j) = x_row + x_col + y_fad;
92  }
93  x.val() = y;
94  for (ordinal i=0; i<outer_fad_size; ++i) {
95  const scalar x_fad = 10.0 + scalar(outer_fad_size) / scalar(i+1);
96  y.val() = x_fad;
97  for (ordinal j=0; j<inner_fad_size; ++j) {
98  const scalar y_fad = 1.0 + scalar(inner_fad_size) / scalar(j+1);
99  y.fastAccessDx(j) = x_row + x_col + x_fad + y_fad;
100  }
101  x.fastAccessDx(i) = y;
102  }
103  return x;
104 }
105 
106 const int global_num_rows = 11;
107 const int global_num_cols = 7;
108 const int global_outer_fad_size = 5;
109 const int global_inner_fad_size = 3;
110 
112  Kokkos_View_FadFad, DeepCopy, FadFadType, Layout, Device )
113 {
114  typedef Kokkos::View<FadFadType**,Layout,Device> ViewType;
115  typedef typename ViewType::size_type size_type;
116  typedef typename ViewType::HostMirror host_view_type;
117 
118  const size_type num_rows = global_num_rows;
119  const size_type num_cols = global_num_cols;
120  const size_type outer_fad_size = global_outer_fad_size;
121  const size_type inner_fad_size = global_inner_fad_size;
122 
123  // Create and fill view
124  ViewType v1("view1", num_rows, num_cols, outer_fad_size+1);
125  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
126  for (size_type i=0; i<num_rows; ++i)
127  for (size_type j=0; j<num_cols; ++j)
128  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
129  num_cols,
130  outer_fad_size,
131  inner_fad_size,
132  i, j);
133  Kokkos::deep_copy(v1, h_v1);
134 
135  // Deep copy
136  ViewType v2("view2", num_rows, num_cols, outer_fad_size+1);
137  Kokkos::deep_copy(v2, v1);
138 
139  // Copy back
140  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
141  Kokkos::deep_copy(h_v2, v2);
142 
143  // Check
144  success = true;
145  for (size_type i=0; i<num_rows; ++i) {
146  for (size_type j=0; j<num_cols; ++j) {
147  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
148  num_cols,
149  outer_fad_size,
150  inner_fad_size,
151  i, j);
152  success = success && checkNestedFads(f, h_v2(i,j), out);
153  }
154  }
155 }
156 
157 #ifdef HAVE_SACADO_KOKKOSCONTAINERS
158 
160  Kokkos_DynRankView_FadFad, DeepCopy, FadFadType, Layout, Device )
161 {
162  typedef Kokkos::DynRankView<FadFadType,Layout,Device> ViewType;
163  typedef typename ViewType::size_type size_type;
164  typedef typename ViewType::HostMirror host_view_type;
165 
166  const size_type num_rows = global_num_rows;
167  const size_type num_cols = global_num_cols;
168  const size_type outer_fad_size = global_outer_fad_size;
169  const size_type inner_fad_size = global_inner_fad_size;
170 
171  // Create and fill view
172  ViewType v1("view1", num_rows, num_cols, outer_fad_size+1);
173  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
174  for (size_type i=0; i<num_rows; ++i)
175  for (size_type j=0; j<num_cols; ++j)
176  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
177  num_cols,
178  outer_fad_size,
179  inner_fad_size,
180  i, j);
181  Kokkos::deep_copy(v1, h_v1);
182 
183  // Deep copy
184  ViewType v2("view2", num_rows, num_cols, outer_fad_size+1);
185  Kokkos::deep_copy(v2, v1);
186 
187  // Copy back
188  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
189  Kokkos::deep_copy(h_v2, v2);
190 
191  // Check
192  success = true;
193  for (size_type i=0; i<num_rows; ++i) {
194  for (size_type j=0; j<num_cols; ++j) {
195  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
196  num_cols,
197  outer_fad_size,
198  inner_fad_size,
199  i, j);
200  success = success && checkNestedFads(f, h_v2(i,j), out);
201  }
202  }
203 }
204 
205 // To test DynRankView - View interoperabitlity
206 // Deep copy of DynRankView to View
207 // Assignment operator
209  Kokkos_DynRankView_FadFad, Interop, FadFadType, Layout, Device )
210 {
211  typedef Kokkos::DynRankView<FadFadType,Layout,Device> DRViewType;
212  typedef typename DRViewType::size_type size_type;
213  typedef typename DRViewType::HostMirror host_view_type;
214 
215  typedef Kokkos::View<FadFadType**,Layout,Device> NoDynViewType;
216  typedef typename NoDynViewType::HostMirror host_nondynrankview_type;
217 
218  const size_type num_rows = global_num_rows;
219  const size_type num_cols = global_num_cols;
220  const size_type outer_fad_size = global_outer_fad_size;
221  const size_type inner_fad_size = global_inner_fad_size;
222 
223  // Create and fill view
224  DRViewType v1("drview1", num_rows, num_cols, outer_fad_size+1);
225  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
226 
227  NoDynViewType ndv2("nodview2", num_rows, num_cols, outer_fad_size+1);
228  host_nondynrankview_type h_ndv2 = Kokkos::create_mirror_view(ndv2);
229 
230  for (size_type i=0; i<num_rows; ++i)
231  for (size_type j=0; j<num_cols; ++j)
232  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
233  num_cols,
234  outer_fad_size,
235  inner_fad_size,
236  i, j);
237  Kokkos::deep_copy(v1, h_v1); //v1 unused here
238 
239  // Deep copy DynRankView to View on device
240  Kokkos::deep_copy(ndv2, h_v1);
241  // Assign View to DynRankView
242  DRViewType v2("drview2", num_rows, num_cols, outer_fad_size+1);
243  v2 = ndv2 ;
244 
245  // Copy back
246  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
247  Kokkos::deep_copy(h_v2, v2);
248 
249  // Check
250  success = true;
251  for (size_type i=0; i<num_rows; ++i) {
252  for (size_type j=0; j<num_cols; ++j) {
253  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
254  num_cols,
255  outer_fad_size,
256  inner_fad_size,
257  i, j);
258  success = success && checkNestedFads(f, h_v2(i,j), out);
259  }
260  }
261 }
262 
263 
264 // Deep copy of DynRankView to View
265 // Copy ctor
267  Kokkos_DynRankView_FadFad, Interop2, FadFadType, Layout, Device )
268 {
269  typedef Kokkos::DynRankView<FadFadType,Layout,Device> DRViewType;
270  typedef typename DRViewType::size_type size_type;
271  typedef typename DRViewType::HostMirror host_view_type;
272 
273  typedef Kokkos::View<FadFadType**,Layout,Device> NoDynViewType;
274  typedef typename NoDynViewType::HostMirror host_nondynrankview_type;
275 
276  const size_type num_rows = global_num_rows;
277  const size_type num_cols = global_num_cols;
278  const size_type outer_fad_size = global_outer_fad_size;
279  const size_type inner_fad_size = global_inner_fad_size;
280 
281  // Create and fill view
282  DRViewType v1("drview1", num_rows, num_cols, outer_fad_size+1);
283  host_view_type h_v1 = Kokkos::create_mirror_view(v1);
284 
285  NoDynViewType ndv2("nodview2", num_rows, num_cols, outer_fad_size+1);
286  host_nondynrankview_type h_ndv2 = Kokkos::create_mirror_view(ndv2);
287 
288  for (size_type i=0; i<num_rows; ++i)
289  for (size_type j=0; j<num_cols; ++j)
290  h_v1(i,j) = generate_nested_fad<FadFadType>(num_rows,
291  num_cols,
292  outer_fad_size,
293  inner_fad_size,
294  i, j);
295  Kokkos::deep_copy(v1, h_v1); //v1 unused here
296 
297  // Deep copy DynRankView to View on device
298  Kokkos::deep_copy(ndv2, h_v1);
299  // Copy construct DynRankView from View
300  DRViewType v2(ndv2) ;
301 
302  // Copy back
303  host_view_type h_v2 = Kokkos::create_mirror_view(v2);
304  Kokkos::deep_copy(h_v2, v2);
305 
306  // Check
307  success = true;
308  for (size_type i=0; i<num_rows; ++i) {
309  for (size_type j=0; j<num_cols; ++j) {
310  FadFadType f = generate_nested_fad<FadFadType>(num_rows,
311  num_cols,
312  outer_fad_size,
313  inner_fad_size,
314  i, j);
315  success = success && checkNestedFads(f, h_v2(i,j), out);
316  }
317  }
318 }
319 
320 
321 
322 #else
323 
325  Kokkos_DynRankView_FadFad, DeepCopy, FadFadType, Layout, Device ) {}
327  Kokkos_DynRankView_FadFad, Interop, FadFadType, Layout, Device ) {}
329  Kokkos_DynRankView_FadFad, Interop2, FadFadType, Layout, Device ) {}
330 
331 #endif
332 
333 #define VIEW_FAD_TESTS_FLD( F, L, D ) \
334  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_FadFad, DeepCopy, F, L, D ) \
335  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_DynRankView_FadFad, DeepCopy, F, L, D ) \
336  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_DynRankView_FadFad, Interop, F, L, D ) \
337  TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_DynRankView_FadFad, Interop2, F, L, D )
338 
339 #define VIEW_FAD_TESTS_FD( F, D ) \
340  using Kokkos::LayoutLeft; \
341  using Kokkos::LayoutRight; \
342  VIEW_FAD_TESTS_FLD( F, LayoutLeft, D) \
343  VIEW_FAD_TESTS_FLD( F, LayoutRight, D)
344 
345 // We've unified the implementation for the different Fad variants, so
346 // there is no reason to test ELRFad, CacheFad, and ELRCacheFad.
351 
352 // These tests are only relevant when we have the experimental view
353 // specialization
354 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
355 #if SACADO_TEST_DFAD
356 #define VIEW_FAD_TESTS_D( D ) \
357  VIEW_FAD_TESTS_FD( SFadType, D ) \
358  VIEW_FAD_TESTS_FD( SLFadType, D ) \
359  VIEW_FAD_TESTS_FD( DFadType, D )
360 #else
361 #define VIEW_FAD_TESTS_D( D ) \
362  VIEW_FAD_TESTS_FD( SFadType, D ) \
363  VIEW_FAD_TESTS_FD( SLFadType, D )
364 #endif
365 #else
366 #define VIEW_FAD_TESTS_D( D ) /* */
367 #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
Sacado::Fad::DFad< double > DFadType
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
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)