Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
advection/common.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 
10 #pragma once
11 
12 const int fad_dim = 50;
16 
17 template <typename ExecSpace>
18 struct is_cuda_space {
19  static const bool value = false;
20 };
21 
22 #ifdef KOKKOS_ENABLE_CUDA
23 template <>
24 struct is_cuda_space<Kokkos::Cuda> {
25  static const bool value = true;
26 };
27 #endif
28 
29 template <typename scalar>
30 scalar
31 generate_fad(const size_t n0, const size_t n1,
32  const size_t n2, const size_t n3, const int fad_size,
33  const size_t i0, const size_t i1,
34  const size_t i2, const size_t i3,
35  const int i_fad)
36 {
37  const scalar x0 = 10.0 + scalar(n0) / scalar(i0+1);
38  const scalar x1 = 100.0 + scalar(n1) / scalar(i1+1);
39  const scalar x2 = 1000.0 + scalar(n2) / scalar(i2+1);
40  const scalar x3 = 10000.0 + scalar(n3) / scalar(i3+1);
41  const scalar x = x0 + x1 + x2 + x3;
42  if (i_fad == fad_size)
43  return x;
44  const scalar x_fad = 1.0 + scalar(fad_size) / scalar(i_fad+1);
45  return x + x_fad;
46 }
47 
48 template <typename V1, typename V2, typename V3, typename V4, typename V5>
49 void init_fad(const V1& v1, const V2& v2, const V3& v3, const V4& v4,
50  const V5& v5)
51 {
52  typedef typename V1::non_const_value_type::value_type scalar;
53 
54  const int ncells = v1.extent(0);
55  const int num_basis = v1.extent(1);
56  const int num_points = v1.extent(2);
57  const int ndim = v1.extent(3);
58  const int N = Kokkos::dimension_scalar(v1)-1;
59 
60  // Kokkos::deep_copy(typename V1::array_type(v1), 1.0);
61  // Kokkos::deep_copy(typename V2::array_type(v2), 2.0);
62  // Kokkos::deep_copy(typename V3::array_type(v3), 3.0);
63  // Kokkos::deep_copy(typename V4::array_type(v4), 4.0);
64 
65  auto v1_h = Kokkos::create_mirror_view(v1);
66  auto v2_h = Kokkos::create_mirror_view(v2);
67  auto v3_h = Kokkos::create_mirror_view(v3);
68  auto v4_h = Kokkos::create_mirror_view(v4);
69  for (int cell=0; cell<ncells; ++cell) {
70  for (int basis=0; basis<num_basis; ++basis) {
71  for (int qp=0; qp<num_points; ++qp) {
72  for (int dim=0; dim<ndim; ++dim) {
73  for (int i=0; i<N; ++i)
74  v1_h(cell,basis,qp,dim).fastAccessDx(i) =
75  generate_fad<scalar>(ncells,num_basis,num_points,ndim,N,cell,basis,qp,dim,i);
76  v1_h(cell,basis,qp,dim).val() =
77  generate_fad<scalar>(ncells,num_basis,num_points,ndim,N,cell,basis,qp,dim,N);
78  }
79  for (int i=0; i<N; ++i)
80  v2_h(cell,basis,qp).fastAccessDx(i) =
81  generate_fad<scalar>(ncells,num_basis,num_points,1,N,cell,basis,qp,0,i);
82  v2_h(cell,basis,qp).val() =
83  generate_fad<scalar>(ncells,num_basis,num_points,1,N,cell,basis,qp,0,N);
84  }
85  }
86  for (int qp=0; qp<num_points; ++qp) {
87  for (int dim=0; dim<ndim; ++dim) {
88  for (int i=0; i<N; ++i)
89  v3_h(cell,qp,dim).fastAccessDx(i) =
90  generate_fad<scalar>(ncells,1,num_points,ndim,N,cell,0,qp,dim,i);
91  v3_h(cell,qp,dim).val() =
92  generate_fad<scalar>(ncells,1,num_points,ndim,N,cell,0,qp,dim,N);
93  }
94  for (int i=0; i<N; ++i)
95  v4_h(cell,qp).fastAccessDx(i) =
96  generate_fad<scalar>(ncells,1,num_points,1,N,cell,0,qp,0,i);
97  v4_h(cell,qp).val() =
98  generate_fad<scalar>(ncells,1,num_points,1,N,cell,0,qp,0,N);
99  }
100  }
101 
102  Kokkos::deep_copy( v1, v1_h );
103  Kokkos::deep_copy( v2, v2_h );
104  Kokkos::deep_copy( v3, v3_h );
105  Kokkos::deep_copy( v4, v4_h );
106 
107  Kokkos::deep_copy(typename V5::array_type(v5), 0.0);
108 }
109 
110 template <typename V1, typename V2, typename V3, typename V4, typename V5>
111 void init_array(const V1& v1, const V2& v2, const V3& v3, const V4& v4,
112  const V5& v5)
113 {
114  typedef typename V1::non_const_value_type scalar;
115 
116  const int ncells = v1.extent(0);
117  const int num_basis = v1.extent(1);
118  const int num_points = v1.extent(2);
119  const int ndim = v1.extent(3);
120  const int N = v1.extent(4)-1;
121 
122  // Kokkos::deep_copy(typename V1::array_type(v1), 1.0);
123  // Kokkos::deep_copy(typename V2::array_type(v2), 2.0);
124  // Kokkos::deep_copy(typename V3::array_type(v3), 3.0);
125  // Kokkos::deep_copy(typename V4::array_type(v4), 4.0);
126 
127  auto v1_h = Kokkos::create_mirror_view(v1);
128  auto v2_h = Kokkos::create_mirror_view(v2);
129  auto v3_h = Kokkos::create_mirror_view(v3);
130  auto v4_h = Kokkos::create_mirror_view(v4);
131  for (int cell=0; cell<ncells; ++cell) {
132  for (int basis=0; basis<num_basis; ++basis) {
133  for (int qp=0; qp<num_points; ++qp) {
134  for (int dim=0; dim<ndim; ++dim) {
135  for (int i=0; i<N; ++i)
136  v1_h(cell,basis,qp,dim,i) =
137  generate_fad<scalar>(ncells,num_basis,num_points,ndim,N,cell,basis,qp,dim,i);
138  v1_h(cell,basis,qp,dim,N) =
139  generate_fad<scalar>(ncells,num_basis,num_points,ndim,N,cell,basis,qp,dim,N);
140  }
141  for (int i=0; i<N; ++i)
142  v2_h(cell,basis,qp,i) =
143  generate_fad<scalar>(ncells,num_basis,num_points,1,N,cell,basis,qp,0,i);
144  v2_h(cell,basis,qp,N) =
145  generate_fad<scalar>(ncells,num_basis,num_points,1,N,cell,basis,qp,0,N);
146  }
147  }
148  for (int qp=0; qp<num_points; ++qp) {
149  for (int dim=0; dim<ndim; ++dim) {
150  for (int i=0; i<N; ++i)
151  v3_h(cell,qp,dim,i) =
152  generate_fad<scalar>(ncells,1,num_points,ndim,N,cell,0,qp,dim,i);
153  v3_h(cell,qp,dim,N) =
154  generate_fad<scalar>(ncells,1,num_points,ndim,N,cell,0,qp,dim,N);
155  }
156  for (int i=0; i<N; ++i)
157  v4_h(cell,qp,i) =
158  generate_fad<scalar>(ncells,1,num_points,1,N,cell,0,qp,0,i);
159  v4_h(cell,qp,N) =
160  generate_fad<scalar>(ncells,1,num_points,1,N,cell,0,qp,0,N);
161  }
162  }
163 
164  Kokkos::deep_copy( v1, v1_h );
165  Kokkos::deep_copy( v2, v2_h );
166  Kokkos::deep_copy( v3, v3_h );
167  Kokkos::deep_copy( v4, v4_h );
168 
169  Kokkos::deep_copy(typename V5::array_type(v5), 0.0);
170 }
171 
172 template <typename View1, typename View2>
174 check(const View1& v_gold, const View2& v, const double tol)
175 {
176  // Copy to host
177  typename View1::HostMirror v_gold_h = Kokkos::create_mirror_view(v_gold);
178  typename View2::HostMirror v_h = Kokkos::create_mirror_view(v);
179  Kokkos::deep_copy(v_gold_h, v_gold);
180  Kokkos::deep_copy(v_h, v);
181 
182  typedef typename View1::value_type value_type;
183 
184  const size_t n0 = v_gold_h.extent(0);
185  const size_t n1 = v_gold_h.extent(1);
186  const size_t n2 = v_gold_h.extent(2);
187 
188  bool success = true;
189  for ( size_t i0 = 0 ; i0 < n0 ; ++i0 ) {
190  for ( size_t i1 = 0 ; i1 < n1 ; ++i1 ) {
191  for ( size_t i2 = 0 ; i2 < n2 ; ++i2 ) {
192  value_type x_gold = v_gold_h(i0,i1,i2);
193  value_type x = v_h(i0,i1,i2);
194  if (std::abs(x_gold-x) > tol*std::abs(x_gold)) {
195  std::cout << "Comparison failed! x_gold("
196  << i0 << "," << i1 << "," << i2 << ") = "
197  << x_gold << " , x = " << x
198  << std::endl;
199  success = false;
200  }
201  }
202  }
203  }
204 
205  return success;
206 }
207 
208 template <typename View1, typename View2>
210 check(const View1& v_gold, const View2& v, const double tol)
211 {
212  // Copy to host
213  typename View1::HostMirror v_gold_h = Kokkos::create_mirror_view(v_gold);
214  typename View2::HostMirror v_h = Kokkos::create_mirror_view(v);
215  Kokkos::deep_copy(v_gold_h, v_gold);
216  Kokkos::deep_copy(v_h, v);
217 
218  typedef typename View1::value_type value_type;
219 
220  const size_t n0 = v_gold_h.extent(0);
221  const size_t n1 = v_gold_h.extent(1);
222  const size_t n2 = v_gold_h.extent(2);
223 
224  bool success = true;
225  for ( size_t i0 = 0 ; i0 < n0 ; ++i0 ) {
226  for ( size_t i1 = 0 ; i1 < n1 ; ++i1 ) {
227  for ( size_t i2 = 0 ; i2 < n2 ; ++i2 ) {
228  value_type x_gold = v_gold_h(i0,i1,i2);
229  value_type x = (i2 == n2-1) ? v_h(i0,i1).val() : v_h(i0,i1).dx(i2);
230  if (std::abs(x_gold-x) > tol*std::abs(x_gold)) {
231  std::cout << "Comparison failed! x_gold("
232  << i0 << "," << i1 << "," << i2 << ") = "
233  << x_gold << " , x = " << x
234  << std::endl;
235  success = false;
236  }
237  }
238  }
239  }
240 
241  return success;
242 }
243 
244 template<typename FluxView, typename WgbView, typename SrcView,
245  typename WbsView>
246 Kokkos::View<double***,typename FluxView::execution_space>
248  const FluxView& flux, const WgbView& wgb, const SrcView& src,
249  const WbsView& wbs,
250  typename std::enable_if< Kokkos::is_view_fad<FluxView>::value>::type* = 0)
251 {
252  typedef typename FluxView::execution_space execution_space;
253 
254  const size_t num_cells = wgb.extent(0);
255  const int num_basis = wgb.extent(1);
256  const int num_points = wgb.extent(2);
257  const int num_dim = wgb.extent(3);
258  const int N = Kokkos::dimension_scalar(wgb)-1;
259 
260  Kokkos::View<double***,typename FluxView::execution_space> residual(
261  "",num_cells,num_basis,N+1);
262 
263  Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>( 0,num_cells ),
264  KOKKOS_LAMBDA (const size_t cell)
265  {
266  double value, value2;
267 
268  // Value
269  for (int basis=0; basis<num_basis; ++basis) {
270  value = value2 = 0.0;
271  for (int qp=0; qp<num_points; ++qp) {
272  for (int dim=0; dim<num_dim; ++dim)
273  value += flux(cell,qp,dim).val()*wgb(cell,basis,qp,dim).val();
274  value2 += src(cell,qp).val()*wbs(cell,basis,qp).val();
275  }
276  residual(cell,basis,N) = value+value2;
277  }
278 
279  // Derivatives
280  for (int k=0; k<N; ++k) {
281  for (int basis=0; basis<num_basis; ++basis) {
282  value = value2 = 0.0;
283  for (int qp=0; qp<num_points; ++qp) {
284  for (int dim=0; dim<num_dim; ++dim)
285  value +=
286  flux(cell,qp,dim).val()*wgb(cell,basis,qp,dim).dx(k) +
287  flux(cell,qp,dim).dx(k)*wgb(cell,basis,qp,dim).val();
288  value2 +=
289  src(cell,qp).val()*wbs(cell,basis,qp).dx(k) +
290  src(cell,qp).dx(k)*wbs(cell,basis,qp).val();
291  }
292  residual(cell,basis,k) = value+value2;
293  }
294  }
295  });
296 
297  return residual;
298 }
299 
300 template<typename FluxView, typename WgbView, typename SrcView,
301  typename WbsView>
302 Kokkos::View<double***,typename FluxView::execution_space>
304  const FluxView& flux, const WgbView& wgb, const SrcView& src,
305  const WbsView& wbs,
306  typename std::enable_if< !Kokkos::is_view_fad<FluxView>::value>::type* = 0)
307 {
308  typedef typename FluxView::execution_space execution_space;
309 
310  const size_t num_cells = wgb.extent(0);
311  const int num_basis = wgb.extent(1);
312  const int num_points = wgb.extent(2);
313  const int num_dim = wgb.extent(3);
314  const int N = wgb.extent(4)-1;
315 
316  Kokkos::View<double***,typename FluxView::execution_space> residual(
317  "",num_cells,num_basis,N+1);
318 
319  Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>( 0,num_cells ),
320  KOKKOS_LAMBDA (const size_t cell)
321  {
322  double value, value2;
323 
324  // Value
325  for (int basis=0; basis<num_basis; ++basis) {
326  value = value2 = 0.0;
327  for (int qp=0; qp<num_points; ++qp) {
328  for (int dim=0; dim<num_dim; ++dim)
329  value += flux(cell,qp,dim,N)*wgb(cell,basis,qp,dim,N);
330  value2 += src(cell,qp,N)*wbs(cell,basis,qp,N);
331  }
332  residual(cell,basis,N) = value+value2;
333  }
334 
335  // Derivatives
336  for (int k=0; k<N; ++k) {
337  for (int basis=0; basis<num_basis; ++basis) {
338  value = value2 = 0.0;
339  for (int qp=0; qp<num_points; ++qp) {
340  for (int dim=0; dim<num_dim; ++dim)
341  value +=
342  flux(cell,qp,dim,N)*wgb(cell,basis,qp,dim,k) +
343  flux(cell,qp,dim,k)*wgb(cell,basis,qp,dim,N);
344  value2 +=
345  src(cell,qp,N)*wbs(cell,basis,qp,k) +
346  src(cell,qp,k)*wbs(cell,basis,qp,N);
347  }
348  residual(cell,basis,k) = value+value2;
349  }
350  }
351  });
352 
353  return residual;
354 }
355 
356 template<typename FluxView, typename WgbView, typename SrcView,
357  typename WbsView, typename ResidualView>
358 void check_residual(const FluxView& flux, const WgbView& wgb,
359  const SrcView& src, const WbsView& wbs,
360  const ResidualView& residual)
361 {
362  // Generate gold residual
363  auto residual_gold = compute_gold_residual(flux, wgb, src, wbs);
364 
365  // Compare residual and residual_gold
366  const double tol = 1.0e-14;
367  check(residual_gold, residual, tol);
368 }
abs(expr.val())
std::enable_if< !Kokkos::is_view_fad< View2 >::value, bool >::type check(const View1 &v_gold, const View2 &v, const double tol)
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
expr val()
Kokkos::View< double ***, typename FluxView::execution_space > compute_gold_residual(const FluxView &flux, const WgbView &wgb, const SrcView &src, const WbsView &wbs, typename std::enable_if< Kokkos::is_view_fad< FluxView >::value >::type *=0)
void init_array(const V1 &v1, const V2 &v2, const V3 &v3, const V4 &v4, const V5 &v5)
void init_fad(const V1 &v1, const V2 &v2, const V3 &v3, const V4 &v4, const V5 &v5)
Sacado::Fad::SLFad< double, fad_dim > SLFadType
const int N
Sacado::Fad::DFad< double > DFadType
void check_residual(const FluxView &flux, const WgbView &wgb, const SrcView &src, const WbsView &wbs, const ResidualView &residual)
static const bool value
const int fad_dim
int value
const double tol