Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ViewFactoryTests.cpp
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 
13 
14 #include "Sacado.hpp"
15 
16 #include "Kokkos_ViewFactory.hpp"
17 
18 TEUCHOS_UNIT_TEST(view_factory, dyn_rank_views)
19 {
21  using Kokkos::View;
22  using Kokkos::DynRankView;
26  using Kokkos::dimension_scalar;
27  using Kokkos::view_alloc;
28  using Kokkos::WithoutInitializing;
29  const unsigned derivative_dim_plus_one = 7;
30 
31 
32  // Test constructing a View pod using Kokkos::view_alloc with deduce_value_type
33  {
34  // Typedef View
35  typedef View<double**, Kokkos::DefaultExecutionSpace> view_type;
36 
37  // Create two rank 2 Views that will be used for deducing types and Fad dims
38  view_type v1("v1", 10, 4);
39  view_type v2("v2", 10, 4);
40 
41  // Get common type of the Views
42  using CommonValueType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::value_type;
43  using ScalarArrayType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::scalar_array_type;
44  // Create an instance of this returned type to pass to ViewCtorProp via view_alloc function
45  auto cvt_for_ctorprop = Kokkos::common_view_alloc_prop(v1, v2);
46 
47  // Create a view with the common type and the max fad_dim of the views passed to deduce_value_type
48  typedef View< CommonValueType** > ViewCommonType;
49  ViewCommonType vct1( Kokkos::view_alloc("vct1", cvt_for_ctorprop), 10, 4 ); // fad_dim deduced and comes from the cvt_for_ctorprop
50 
51  TEST_EQUALITY(vct1.extent(0), v1.extent(0));
52  TEST_EQUALITY(vct1.extent(1), v1.extent(1));
53  TEST_EQUALITY(vct1.extent(2), v1.extent(2));
54  TEST_EQUALITY( Kokkos::dimension_scalar(vct1), 0);
56  bool check_eq_scalar_double = std::is_same < double, ScalarArrayType >::value;
57  TEST_EQUALITY(check_eq_kokkos_type, true);
58  TEST_EQUALITY(check_eq_scalar_double, true);
59  }
60  // Test constructing a View of Fad using Kokkos::view_alloc with deduce_value_type
61  {
62  // Typedef View
63  typedef View<FadType**, Kokkos::DefaultExecutionSpace> view_type;
64 
65  // Create two rank 2 Views that will be used for deducing types and Fad dims
66  view_type v1("v1", 10, 4, derivative_dim_plus_one );
67  view_type v2("v2", 10, 4, derivative_dim_plus_one );
68 
69  // Get common type of the Views
70  using CommonValueType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::value_type;
71  using ScalarArrayType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::scalar_array_type;
72  // Create an instance of this returned type to pass to ViewCtorProp via view_alloc function
73  auto cvt_for_ctorprop = Kokkos::common_view_alloc_prop(v1, v2);
74 
75  // Create a view with the common type and the max fad_dim of the views passed to deduce_value_type
76  typedef View< CommonValueType** > ViewCommonType;
77  ViewCommonType vct1( Kokkos::view_alloc("vct1", cvt_for_ctorprop), 10, 4 ); // fad_dim deduced and comes from the cvt_for_ctorprop
78 
79  TEST_EQUALITY(dimension_scalar(vct1), derivative_dim_plus_one);
80  TEST_EQUALITY(vct1.extent(0), v1.extent(0));
81  TEST_EQUALITY(vct1.extent(1), v1.extent(1));
82  TEST_EQUALITY(vct1.extent(2), v1.extent(2));
84  bool check_eq_fad_type = std::is_same < CommonValueType, FadType >::value;
85  bool check_eq_scalar_double = std::is_same < double, ScalarArrayType >::value;
86  TEST_EQUALITY(check_neq_kokkos_type, false);
87  TEST_EQUALITY(check_eq_fad_type, true);
88  TEST_EQUALITY(check_eq_scalar_double, true);
89  }
90  // Test constructing a View from mix of View and Viewof Fads using Kokkos::view_alloc with deduce_value_type
91  {
92  // Typedef View
93  typedef View<FadType**, Kokkos::DefaultExecutionSpace> view_of_fad_type;
94  typedef View<double**, Kokkos::DefaultExecutionSpace> view_of_pod_type;
95 
96  // Create two rank 2 Views that will be used for deducing types and Fad dims
97  view_of_fad_type v1("v1", 10, 4, derivative_dim_plus_one );
98  view_of_pod_type v2("v2", 10, 4);
99 
100  // Get common type of the Views
101  using CommonValueType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::value_type;
102  using ScalarArrayType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::scalar_array_type;
103  // Create an instance of this returned type to pass to ViewCtorProp via view_alloc function
104  auto cvt_for_ctorprop = Kokkos::common_view_alloc_prop(v1, v2);
105 
106  // Create a view with the common type and the max fad_dim of the views passed to deduce_value_type
107  typedef View< CommonValueType** > ViewCommonType;
108  ViewCommonType vct1( Kokkos::view_alloc("vct1", cvt_for_ctorprop), 10, 4 ); // fad_dim deduced and comes from the cvt_for_ctorprop
109 
110  TEST_EQUALITY(dimension_scalar(vct1), derivative_dim_plus_one);
111  TEST_EQUALITY(vct1.extent(0), v1.extent(0));
112  TEST_EQUALITY(vct1.extent(1), v1.extent(1));
113  TEST_EQUALITY(vct1.extent(2), v1.extent(2));
115  bool check_eq_fad_type = std::is_same < CommonValueType, FadType >::value;
116  bool check_eq_scalar_double = std::is_same < double, ScalarArrayType >::value;
117  TEST_EQUALITY(check_neq_kokkos_type, false);
118  TEST_EQUALITY(check_eq_fad_type, true);
119  TEST_EQUALITY(check_eq_scalar_double, true);
120  }
121  // Test constructing a DynRankView using Kokkos::view_alloc with deduce_value_type
122  {
123  // Typedef View
124  typedef DynRankView<FadType, Kokkos::DefaultExecutionSpace> view_type;
125 
126  // Create two rank 2 Views that will be used for deducing types and Fad dims
127  view_type v1("v1", 10, 4, derivative_dim_plus_one );
128  view_type v2("v2", 10, 4, derivative_dim_plus_one );
129 
130  // Get common type of the Views
131  using CommonValueType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::value_type;
132  using ScalarArrayType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::scalar_array_type;
133  // Create an instance of this returned type to pass to ViewCtorProp via view_alloc function
134  auto cvt_for_ctorprop = Kokkos::common_view_alloc_prop(v1, v2);
135 
136  // Create a view with the common type and the max fad_dim of the views passed to deduce_value_type
137  typedef DynRankView< CommonValueType > ViewCommonType;
138  ViewCommonType vct1( Kokkos::view_alloc("vct1", cvt_for_ctorprop), 10, 4 ); // fad_dim deduced and comes from the cvt_for_ctorprop
139 
140  TEST_EQUALITY(dimension_scalar(vct1), derivative_dim_plus_one);
141  TEST_EQUALITY(vct1.extent(0), v1.extent(0));
142  TEST_EQUALITY(vct1.extent(1), v1.extent(1));
143  TEST_EQUALITY(vct1.extent(2), v1.extent(2));
144  TEST_EQUALITY(Kokkos::rank(vct1), 2);
146  bool check_eq_fad_type = std::is_same < CommonValueType, FadType >::value;
147  bool check_eq_scalar_double = std::is_same < double, ScalarArrayType >::value;
148  TEST_EQUALITY(check_neq_kokkos_type, false);
149  TEST_EQUALITY(check_eq_fad_type, true);
150  TEST_EQUALITY(check_eq_scalar_double, true);
151  }
152  // Test constructing a DynRankView from mix of DynRankView and DynRankView of Fads using Kokkos::view_alloc with deduce_value_type
153  {
154  // Typedef View
155  typedef DynRankView<FadType, Kokkos::DefaultExecutionSpace> view_of_fad_type;
156  typedef DynRankView<double, Kokkos::DefaultExecutionSpace> view_of_pod_type;
157 
158  // Create two rank 2 Views that will be used for deducing types and Fad dims
159  view_of_fad_type v1("v1", 10, 4, derivative_dim_plus_one );
160  view_of_pod_type v2("v2", 10, 4);
161 
162  // Get common type of the Views
163  using CommonValueType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::value_type;
164  using ScalarArrayType = typename decltype( Kokkos::common_view_alloc_prop( v1, v2 ) )::scalar_array_type;
165  // Create an instance of this returned type to pass to ViewCtorProp via view_alloc function
166  auto cvt_for_ctorprop = Kokkos::common_view_alloc_prop(v1, v2);
167 
168  // Create a view with the common type and the max fad_dim of the views passed to deduce_value_type
169  typedef DynRankView< CommonValueType > ViewCommonType;
170  ViewCommonType vct1( Kokkos::view_alloc("vct1", cvt_for_ctorprop), 10, 4 ); // fad_dim deduced and comes from the cvt_for_ctorprop
171 
172  TEST_EQUALITY(dimension_scalar(vct1), derivative_dim_plus_one);
173  TEST_EQUALITY(vct1.extent(0), v1.extent(0));
174  TEST_EQUALITY(vct1.extent(1), v1.extent(1));
175  TEST_EQUALITY(vct1.extent(2), v1.extent(2));
176  TEST_EQUALITY(Kokkos::rank(vct1), 2);
178  bool check_eq_fad_type = std::is_same < CommonValueType, FadType >::value;
179  bool check_eq_scalar_double = std::is_same < double, ScalarArrayType >::value;
180  TEST_EQUALITY(check_neq_kokkos_type, false);
181  TEST_EQUALITY(check_eq_fad_type, true);
182  TEST_EQUALITY(check_eq_scalar_double, true);
183  }
184 
185 
186 
187 
188  // Test a DynRankView from a DynRankView
189  {
190  DynRankView<FadType> a("a",10,4,13,derivative_dim_plus_one);
191  TEST_EQUALITY(dimension_scalar(a),derivative_dim_plus_one);
192  TEST_EQUALITY(a.rank(),3);
193 
194  auto b = createDynRankView(a,"b",5,3,8);
195  TEST_EQUALITY(dimension_scalar(b),derivative_dim_plus_one);
196  TEST_EQUALITY(b.rank(),3);
197 
198  auto c = createDynRankView(a,view_alloc("c",WithoutInitializing),5,3,8);
199  TEST_EQUALITY(dimension_scalar(c),derivative_dim_plus_one);
200  TEST_EQUALITY(c.rank(),3);
201 
202  using d_type = Kokkos::DynRankView<FadType,Kokkos::LayoutRight>;
203  d_type d = createDynRankViewWithType<d_type>(a,"d",5,3,8);
204  TEST_EQUALITY(dimension_scalar(d),derivative_dim_plus_one);
205  TEST_EQUALITY(d.rank(),3);
206  }
207 
208  // Test a DynRankView from a View
209  {
210  View<FadType*> a("a",8,derivative_dim_plus_one);
211  TEST_EQUALITY(dimension_scalar(a),derivative_dim_plus_one);
212 
213  auto b = createDynRankView(a,"b",5,3,8);
214  TEST_EQUALITY(dimension_scalar(b),derivative_dim_plus_one);
215  TEST_EQUALITY(b.rank(),3);
216 
217  auto c = createDynRankView(a,view_alloc("c",WithoutInitializing),5,3,8);
218  TEST_EQUALITY(dimension_scalar(c),derivative_dim_plus_one);
219  TEST_EQUALITY(c.rank(),3);
220 
221  using d_type = Kokkos::DynRankView<FadType,Kokkos::LayoutRight>;
222  d_type d = createDynRankViewWithType<d_type>(a,"d",5,3,8);
223  TEST_EQUALITY(dimension_scalar(d),derivative_dim_plus_one);
224  TEST_EQUALITY(d.rank(),3);
225  }
226 
227  // Test a View from a View
228  {
229  View<FadType*> a("a",8,derivative_dim_plus_one);
230  TEST_EQUALITY(dimension_scalar(a),derivative_dim_plus_one);
231 
232  using b_type = Kokkos::View<FadType***>;
233  b_type b = createViewWithType<b_type>(a,"b",5,3,8);
234  TEST_EQUALITY(dimension_scalar(b),derivative_dim_plus_one);
235 
236  b_type c = createViewWithType<b_type>(a,view_alloc("c",WithoutInitializing),5,3,8);
237  TEST_EQUALITY(dimension_scalar(c),derivative_dim_plus_one);
238 
239  using d_type = Kokkos::View<FadType***,Kokkos::LayoutRight>;
240  d_type d = createViewWithType<d_type>(a,"d",5,3,8);
241  TEST_EQUALITY(dimension_scalar(d),derivative_dim_plus_one);
242  }
243 
244  // Test a View from a DynRankView
245  {
246  DynRankView<FadType> a("a",10,4,13,derivative_dim_plus_one);
247  TEST_EQUALITY(dimension_scalar(a),derivative_dim_plus_one);
248  TEST_EQUALITY(a.rank(),3);
249 
250  using b_type = Kokkos::View<FadType***>;
251  b_type b = createViewWithType<b_type>(a,"b",5,3,8);
252  TEST_EQUALITY(dimension_scalar(b),derivative_dim_plus_one);
253 
254  b_type c = createViewWithType<b_type>(a,view_alloc("c",WithoutInitializing),5,3,8);
255  TEST_EQUALITY(dimension_scalar(c),derivative_dim_plus_one);
256 
257  using d_type = Kokkos::View<FadType***,Kokkos::LayoutRight>;
258  d_type d = createViewWithType<d_type>(a,"d",5,3,8);
259  TEST_EQUALITY(dimension_scalar(d),derivative_dim_plus_one);
260  }
261 
262  // Test creation of a Fad DynRankView from a double DynRankView
263  {
264  DynRankView<double> a("a",10,4,13);
265  TEST_EQUALITY(dimension_scalar(a),0);
266  TEST_EQUALITY(a.rank(),3);
267 
268  using b_type = Kokkos::DynRankView<FadType,Kokkos::LayoutRight>;
269  b_type b = createDynRankViewWithType<b_type>(a,"b",5,3,8);
270  TEST_EQUALITY(dimension_scalar(b),1);
271  TEST_EQUALITY(b.rank(),3);
272  }
273 
274  // Test a double DynRankView from a double DynRankView
275  {
276  DynRankView<double> a("a",10,4,13);
277  TEST_EQUALITY(dimension_scalar(a),0);
278  TEST_EQUALITY(a.rank(),3);
279 
280  auto b = createDynRankView(a,"b",5,3,8);
281  TEST_EQUALITY(dimension_scalar(b),0);
282  TEST_EQUALITY(b.rank(),3);
283  }
284 
285  // Test double rank 0
286  {
287  DynRankView<double> a("a",10,4,13);
288  TEST_EQUALITY(dimension_scalar(a),0);
289  TEST_EQUALITY(a.rank(),3);
290 
291  auto b = createDynRankView(a,"b");
292  TEST_EQUALITY(dimension_scalar(b),0);
293  TEST_EQUALITY(b.rank(),0);
294  }
295 
296  // Test Fad rank 0
297  {
298  DynRankView<FadType> a("a",10,4,13,derivative_dim_plus_one);
299  TEST_EQUALITY(dimension_scalar(a),derivative_dim_plus_one);
300  TEST_EQUALITY(a.rank(),3);
301 
302  auto b = createDynRankView(a,"b");
303  TEST_EQUALITY(dimension_scalar(b),derivative_dim_plus_one);
304  TEST_EQUALITY(b.rank(),0);
305  }
306 
307  // Test unmanaged view of double
308  {
309  Kokkos::View<double*> a("a",5*3);
310  using b_type = Kokkos::View<double**,Kokkos::MemoryUnmanaged>;
311  b_type b = createViewWithType<b_type>(a,a.data(),5,3);
312  TEST_EQUALITY(b.extent(0),5);
313  TEST_EQUALITY(b.extent(1),3);
314  TEST_EQUALITY(dimension_scalar(b),0);
315  }
316 
317  // Test unmanaged view of Fad
318  {
319  Kokkos::View<FadType*> a("a",5*3,derivative_dim_plus_one);
320  using b_type = Kokkos::View<FadType**,Kokkos::MemoryUnmanaged>;
321  b_type b = createViewWithType<b_type>(a,a.data(),5,3);
322  TEST_EQUALITY(b.extent(0),5);
323  TEST_EQUALITY(b.extent(1),3);
324  TEST_EQUALITY(dimension_scalar(b),derivative_dim_plus_one);
325  }
326 
327  // Test LayoutStride view of double
328  {
329  Kokkos::DynRankView<double> a("a",10,13);
330  auto b = Kokkos::subview(a, std::make_pair(4,8), std::make_pair(5,11));
331  auto c = createDynRankView(b,"c",5,3);
332  using b_type = decltype(b);
333  using c_type = decltype(c);
334  using b_layout = typename b_type::array_layout;
335  using c_layout = typename c_type::array_layout;
336  using default_layout = typename b_type::device_type::execution_space::array_layout;
337  const bool is_b_layout_stride =
339  const bool is_c_default_layout =
341  TEST_EQUALITY(is_b_layout_stride,true);
342  TEST_EQUALITY(is_c_default_layout,true);
343  TEST_EQUALITY(c.rank(),2);
344  TEST_EQUALITY(c.extent(0),5);
345  TEST_EQUALITY(c.extent(1),3);
346  TEST_EQUALITY(dimension_scalar(b),0);
347  }
348 
349  // Test LayoutStride view of Fad
350  {
351  Kokkos::DynRankView<FadType> a("a",10,13,derivative_dim_plus_one);
352  auto b = Kokkos::subview(a, std::make_pair(4,8), std::make_pair(5,11));
353  auto c = createDynRankView(b,"c",5,3);
354  using b_type = decltype(b);
355  using c_type = decltype(c);
356  using b_layout = typename b_type::array_layout;
357  using c_layout = typename c_type::array_layout;
358  using default_layout = typename b_type::device_type::execution_space::array_layout;
359  const bool is_b_layout_stride =
361  const bool is_c_default_layout =
363  TEST_EQUALITY(is_b_layout_stride,true);
364  TEST_EQUALITY(is_c_default_layout,true);
365  TEST_EQUALITY(c.rank(),2);
366  TEST_EQUALITY(c.extent(0),5);
367  TEST_EQUALITY(c.extent(1),3);
368  TEST_EQUALITY(dimension_scalar(b),derivative_dim_plus_one);
369  }
370 
371 }
372 
373 int main( int argc, char* argv[] ) {
374  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
375 
376  Kokkos::initialize();
377 
379 
380  Kokkos::finalize();
381 
382  return res;
383 }
std::enable_if< is_view< InputViewType >::value||is_dyn_rank_view< InputViewType >::value, typename Impl::ResultDynRankView< InputViewType >::type >::type createDynRankView(const InputViewType &a, const CtorProp &prop, const Dims...dims)
Wrapper to simplify use of Sacado ViewFactory.
TEUCHOS_UNIT_TEST(Conversion, IsConvertible)
static int runUnitTestsFromMain(int argc, char *argv[])
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
int main()
Definition: ad_example.cpp:171
std::enable_if< is_view< InputViewType >::value||is_dyn_rank_view< InputViewType >::value, ResultViewType >::type createDynRankViewWithType(const InputViewType &a, const CtorProp &prop, const Dims...dims)
Wrapper to simplify use of Sacado ViewFactory.
std::enable_if< is_view< InputViewType >::value||is_dyn_rank_view< InputViewType >::value, ResultViewType >::type createViewWithType(const InputViewType &a, const CtorProp &prop, const Dims...dims)
Wrapper to simplify use of Sacado ViewFactory.
int value
#define TEST_EQUALITY(v1, v2)