Skip to content

Commit c8bb901

Browse files
committed
Kernalize UpdateNewVariables.
1 parent 9434d33 commit c8bb901

File tree

2 files changed

+26
-3
lines changed

2 files changed

+26
-3
lines changed

swm_AMReX/swm_mini_app_kernels.h

+20
Original file line numberDiff line numberDiff line change
@@ -19,4 +19,24 @@ void UpdateIntermediateVariablesKernel( const int i, const int j, const int k, c
1919
h(i,j,k) = p(i,j,k) + 0.25*(u(i-1,j,k)*u(i-1,j,k) + u(i,j,k)*u(i,j,k) + v(i,j-1,k)*v(i,j-1,k) + v(i,j,k)*v(i,j,k));
2020
}
2121

22+
23+
AMREX_GPU_DEVICE AMREX_FORCE_INLINE
24+
void UpdateNewVariablesKernel( const int i, const int j, const int k,
25+
const double tdtsdx, const double tdtsdy, const double tdts8,
26+
const amrex::Array4<amrex::Real const>& p_old_array,
27+
const amrex::Array4<amrex::Real const>& u_old_array,
28+
const amrex::Array4<amrex::Real const>& v_old_array,
29+
const amrex::Array4<amrex::Real const>& cu_array,
30+
const amrex::Array4<amrex::Real const>& cv_array,
31+
const amrex::Array4<amrex::Real const>& h_array,
32+
const amrex::Array4<amrex::Real const>& z_array,
33+
const amrex::Array4<amrex::Real>& p_new_array,
34+
const amrex::Array4<amrex::Real>& u_new_array,
35+
const amrex::Array4<amrex::Real>& v_new_array)
36+
{
37+
u_new_array(i,j,k) = u_old_array(i,j,k) + tdts8 * (z_array(i,j-1,k)+z_array(i,j,k)) * (cv_array(i,j-1,k) + cv_array(i,j,k) + cv_array(i+1,j-1,k) + cv_array(i+1,j,k)) - tdtsdx * (h_array(i+1,j,k) - h_array(i,j,k));
38+
v_new_array(i,j,k) = v_old_array(i,j,k) - tdts8 * (z_array(i-1,j,k)+z_array(i,j,k)) * (cu_array(i-1,j,k) + cu_array(i-1,j+1,k) + cu_array(i,j,k) + cu_array(i,j+1,k)) - tdtsdy * (h_array(i,j+1,k) - h_array(i,j,k));
39+
p_new_array(i,j,k) = p_old_array(i,j,k) - tdtsdx * (cu_array(i,j,k) - cu_array(i-1,j,k)) - tdtsdy * (cv_array(i,j,k) - cv_array(i,j-1,k));
40+
}
41+
2242
#endif // SWM_MINI_APP_KERNELS_H_

swm_AMReX/swm_mini_app_utils.cpp

+6-3
Original file line numberDiff line numberDiff line change
@@ -389,9 +389,12 @@ void UpdateNewVariables(const double dx, const double dy, const double tdt, cons
389389

390390
amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
391391
{
392-
u_new_array(i,j,k) = u_old_array(i,j,k) + tdts8 * (z_array(i,j-1,k)+z_array(i,j,k)) * (cv_array(i,j-1,k) + cv_array(i,j,k) + cv_array(i+1,j-1,k) + cv_array(i+1,j,k)) - tdtsdx * (h_array(i+1,j,k) - h_array(i,j,k));
393-
v_new_array(i,j,k) = v_old_array(i,j,k) - tdts8 * (z_array(i-1,j,k)+z_array(i,j,k)) * (cu_array(i-1,j,k) + cu_array(i-1,j+1,k) + cu_array(i,j,k) + cu_array(i,j+1,k)) - tdtsdy * (h_array(i,j+1,k) - h_array(i,j,k));
394-
p_new_array(i,j,k) = p_old_array(i,j,k) - tdtsdx * (cu_array(i,j,k) - cu_array(i-1,j,k)) - tdtsdy * (cv_array(i,j,k) - cv_array(i,j-1,k));
392+
UpdateNewVariablesKernel(i, j, k,
393+
tdtsdx, tdtsdy, tdts8,
394+
p_old_array, u_old_array, v_old_array,
395+
cu_array, cv_array, h_array, z_array,
396+
p_new_array, u_new_array, v_new_array);
397+
395398
});
396399
}
397400

0 commit comments

Comments
 (0)