17 template <
typename ExecSpace>
19 static const bool value =
false;
22 #ifdef KOKKOS_ENABLE_CUDA
25 static const bool value =
true;
29 template <
typename scalar>
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,
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)
44 const scalar x_fad = 1.0 + scalar(fad_size) / scalar(i_fad+1);
48 template <
typename WgbView,
typename WbsView,
typename FluxView,
49 typename SrcView,
typename ResidualView>
50 void init_fad(
const WgbView& wgb,
const WbsView& wbs,
const FluxView& flux,
51 const SrcView& src,
const ResidualView& residual)
53 typedef typename ResidualView::non_const_value_type::value_type scalar;
55 const int ncells = wgb.extent(0);
56 const int num_basis = wgb.extent(1);
57 const int num_points = wgb.extent(2);
58 const int ndim = wgb.extent(3);
59 const int N = Kokkos::dimension_scalar(residual)-1;
61 auto wgb_h = Kokkos::create_mirror_view(wgb);
62 auto wbs_h = Kokkos::create_mirror_view(wbs);
63 auto flux_h = Kokkos::create_mirror_view(flux);
64 auto src_h = Kokkos::create_mirror_view(src);
65 for (
int cell=0; cell<ncells; ++cell) {
66 for (
int basis=0; basis<num_basis; ++basis) {
67 for (
int qp=0; qp<num_points; ++qp) {
68 for (
int dim=0; dim<ndim; ++dim) {
69 wgb_h(cell,basis,qp,dim) =
70 generate_fad<scalar>(ncells,num_basis,num_points,ndim,
N,cell,basis,qp,dim,
N);
72 wbs_h(cell,basis,qp) =
73 generate_fad<scalar>(ncells,num_basis,num_points,1,
N,cell,basis,qp,0,
N);
76 for (
int qp=0; qp<num_points; ++qp) {
77 for (
int dim=0; dim<ndim; ++dim) {
78 for (
int i=0;
i<
N; ++
i)
79 flux_h(cell,qp,dim).fastAccessDx(
i) =
80 generate_fad<scalar>(ncells,1,num_points,ndim,
N,cell,0,qp,dim,
i);
81 flux_h(cell,qp,dim).val() =
82 generate_fad<scalar>(ncells,1,num_points,ndim,
N,cell,0,qp,dim,
N);
84 for (
int i=0;
i<
N; ++
i)
85 src_h(cell,qp).fastAccessDx(
i) =
86 generate_fad<scalar>(ncells,1,num_points,1,
N,cell,0,qp,0,
i);
87 src_h(cell,qp).val() =
88 generate_fad<scalar>(ncells,1,num_points,1,
N,cell,0,qp,0,
N);
92 Kokkos::deep_copy( wgb, wgb_h );
93 Kokkos::deep_copy( wbs, wbs_h );
94 Kokkos::deep_copy( flux, flux_h );
95 Kokkos::deep_copy( src, src_h );
97 Kokkos::deep_copy(
typename ResidualView::array_type(residual), 0.0);
100 template <
typename WgbView,
typename WbsView,
typename FluxView,
101 typename SrcView,
typename ResidualView>
102 void init_array(
const WgbView& wgb,
const WbsView& wbs,
const FluxView& flux,
103 const SrcView& src,
const ResidualView& residual)
105 typedef typename ResidualView::non_const_value_type scalar;
107 const int ncells = wgb.extent(0);
108 const int num_basis = wgb.extent(1);
109 const int num_points = wgb.extent(2);
110 const int ndim = wgb.extent(3);
111 const int N = residual.extent(2)-1;
113 auto wgb_h = Kokkos::create_mirror_view(wgb);
114 auto wbs_h = Kokkos::create_mirror_view(wbs);
115 auto flux_h = Kokkos::create_mirror_view(flux);
116 auto src_h = Kokkos::create_mirror_view(src);
117 for (
int cell=0; cell<ncells; ++cell) {
118 for (
int basis=0; basis<num_basis; ++basis) {
119 for (
int qp=0; qp<num_points; ++qp) {
120 for (
int dim=0; dim<ndim; ++dim) {
121 wgb_h(cell,basis,qp,dim) =
122 generate_fad<scalar>(ncells,num_basis,num_points,ndim,
N,cell,basis,qp,dim,
N);
124 wbs_h(cell,basis,qp) =
125 generate_fad<scalar>(ncells,num_basis,num_points,1,
N,cell,basis,qp,0,
N);
128 for (
int qp=0; qp<num_points; ++qp) {
129 for (
int dim=0; dim<ndim; ++dim) {
130 for (
int i=0;
i<
N; ++
i)
131 flux_h(cell,qp,dim,
i) =
132 generate_fad<scalar>(ncells,1,num_points,ndim,
N,cell,0,qp,dim,
i);
133 flux_h(cell,qp,dim,N) =
134 generate_fad<scalar>(ncells,1,num_points,ndim,
N,cell,0,qp,dim,
N);
136 for (
int i=0;
i<
N; ++
i)
138 generate_fad<scalar>(ncells,1,num_points,1,
N,cell,0,qp,0,
i);
140 generate_fad<scalar>(ncells,1,num_points,1,
N,cell,0,qp,0,
N);
144 Kokkos::deep_copy( wgb, wgb_h );
145 Kokkos::deep_copy( wbs, wbs_h );
146 Kokkos::deep_copy( flux, flux_h );
147 Kokkos::deep_copy( src, src_h );
149 Kokkos::deep_copy(residual, 0.0);
152 template <
typename View1,
typename View2>
154 check(
const View1& v_gold,
const View2& v,
const double tol)
157 typename View1::HostMirror v_gold_h = Kokkos::create_mirror_view(v_gold);
158 typename View2::HostMirror v_h = Kokkos::create_mirror_view(v);
159 Kokkos::deep_copy(v_gold_h, v_gold);
160 Kokkos::deep_copy(v_h, v);
162 typedef typename View1::value_type value_type;
164 const size_t n0 = v_gold_h.extent(0);
165 const size_t n1 = v_gold_h.extent(1);
166 const size_t n2 = v_gold_h.extent(2);
169 for (
size_t i0 = 0 ; i0 < n0 ; ++i0 ) {
170 for (
size_t i1 = 0 ; i1 < n1 ; ++i1 ) {
171 for (
size_t i2 = 0 ; i2 < n2 ; ++i2 ) {
172 value_type x_gold = v_gold_h(i0,i1,i2);
173 value_type
x = v_h(i0,i1,i2);
175 std::cout <<
"Comparison failed! x_gold("
176 << i0 <<
"," << i1 <<
"," << i2 <<
") = "
177 << x_gold <<
" , x = " << x
188 template <
typename View1,
typename View2>
190 check(
const View1& v_gold,
const View2& v,
const double tol)
193 typename View1::HostMirror v_gold_h = Kokkos::create_mirror_view(v_gold);
194 typename View2::HostMirror v_h = Kokkos::create_mirror_view(v);
195 Kokkos::deep_copy(v_gold_h, v_gold);
196 Kokkos::deep_copy(v_h, v);
198 typedef typename View1::value_type value_type;
200 const size_t n0 = v_gold_h.extent(0);
201 const size_t n1 = v_gold_h.extent(1);
202 const size_t n2 = v_gold_h.extent(2);
205 for (
size_t i0 = 0 ; i0 < n0 ; ++i0 ) {
206 for (
size_t i1 = 0 ; i1 < n1 ; ++i1 ) {
207 for (
size_t i2 = 0 ; i2 < n2 ; ++i2 ) {
208 value_type x_gold = v_gold_h(i0,i1,i2);
209 value_type
x = (i2 == n2-1) ? v_h(i0,i1).val() : v_h(i0,i1).dx(i2);
211 std::cout <<
"Comparison failed! x_gold("
212 << i0 <<
"," << i1 <<
"," << i2 <<
") = "
213 << x_gold <<
" , x = " << x
224 template<
typename FluxView,
typename WgbView,
typename SrcView,
226 Kokkos::View<double***,typename FluxView::execution_space>
228 const FluxView& flux,
const WgbView& wgb,
const SrcView& src,
232 typedef typename FluxView::execution_space execution_space;
234 const size_t num_cells = wgb.extent(0);
235 const int num_basis = wgb.extent(1);
236 const int num_points = wgb.extent(2);
237 const int num_dim = wgb.extent(3);
238 const int N = Kokkos::dimension_scalar(flux)-1;
240 Kokkos::View<double***,typename FluxView::execution_space> residual(
241 "",num_cells,num_basis,N+1);
243 Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>( 0,num_cells ),
244 KOKKOS_LAMBDA (
const size_t cell)
246 double value, value2;
249 for (
int basis=0; basis<num_basis; ++basis) {
250 value = value2 = 0.0;
251 for (
int qp=0; qp<num_points; ++qp) {
252 for (
int dim=0; dim<num_dim; ++dim)
253 value += flux(cell,qp,dim).
val()*wgb(cell,basis,qp,dim);
254 value2 += src(cell,qp).val()*wbs(cell,basis,qp);
256 residual(cell,basis,N) = value+value2;
260 for (
int k=0; k<
N; ++k) {
261 for (
int basis=0; basis<num_basis; ++basis) {
262 value = value2 = 0.0;
263 for (
int qp=0; qp<num_points; ++qp) {
264 for (
int dim=0; dim<num_dim; ++dim)
266 flux(cell,qp,dim).dx(k)*wgb(cell,basis,qp,dim);
268 src(cell,qp).dx(k)*wbs(cell,basis,qp);
270 residual(cell,basis,k) = value+value2;
278 template<
typename FluxView,
typename WgbView,
typename SrcView,
280 Kokkos::View<double***,typename FluxView::execution_space>
282 const FluxView& flux,
const WgbView& wgb,
const SrcView& src,
286 typedef typename FluxView::execution_space execution_space;
288 const size_t num_cells = wgb.extent(0);
289 const int num_basis = wgb.extent(1);
290 const int num_points = wgb.extent(2);
291 const int num_dim = wgb.extent(3);
292 const int N = flux.extent(3)-1;
294 Kokkos::View<double***,typename FluxView::execution_space> residual(
295 "",num_cells,num_basis,N+1);
297 Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>( 0,num_cells ),
298 KOKKOS_LAMBDA (
const size_t cell)
300 double value, value2;
301 for (
int k=0; k<=
N; ++k) {
302 for (
int basis=0; basis<num_basis; ++basis) {
303 value = value2 = 0.0;
304 for (
int qp=0; qp<num_points; ++qp) {
305 for (
int dim=0; dim<num_dim; ++dim)
306 value += flux(cell,qp,dim,k)*wgb(cell,basis,qp,dim);
307 value2 += src(cell,qp,k)*wbs(cell,basis,qp);
309 residual(cell,basis,k) = value+value2;
317 template<
typename FluxView,
typename WgbView,
typename SrcView,
318 typename WbsView,
typename ResidualView>
320 const SrcView& src,
const WbsView& wbs,
321 const ResidualView& residual)
327 const double tol = 1.0e-14;
328 check(residual_gold, residual, tol);
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
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
Sacado::Fad::DFad< double > DFadType
void check_residual(const FluxView &flux, const WgbView &wgb, const SrcView &src, const WbsView &wbs, const ResidualView &residual)