33 const BCRec* bc_ptr_h = domain_bcs_type_h.data();
36 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
38 for ( MFIter mfi(density,TilingIfNotGPU()); mfi.isValid(); ++mfi)
41 Box bx = mfi.tilebox();
43 const Box& tbx = surroundingNodes(bx,0);
44 const Box& tby = surroundingNodes(bx,1);
45 const Box& tbz = surroundingNodes(bx,2);
48 const Array4<const Real>& dens_arr = density.array(mfi);
51 Array4<Real const>
const& momx = xmom_in.const_array(mfi);
52 Array4<Real const>
const& momy = ymom_in.const_array(mfi);
53 Array4<Real const>
const& momz = zmom_in.const_array(mfi);
56 const Array4<Real>& velx =
xvel.array(mfi);
57 const Array4<Real>& vely =
yvel.array(mfi);
58 const Array4<Real>& velz =
zvel.array(mfi);
60 ParallelFor(tbx, tby, tbz,
61 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
62 velx(i,j,k) = momx(i,j,k) * 2.0 / (dens_arr(i,j,k,
Rho_comp) + dens_arr(i-1,j,k,
Rho_comp));
64 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
65 vely(i,j,k) = momy(i,j,k) * 2.0 / (dens_arr(i,j,k,
Rho_comp) + dens_arr(i,j-1,k,
Rho_comp));
67 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
68 velz(i,j,k) = momz(i,j,k) * 2.0 / (dens_arr(i,j,k,
Rho_comp) + dens_arr(i,j,k-1,
Rho_comp));
71 if ( (bx.smallEnd(0) == domain.smallEnd(0)) &&
73 ParallelFor(makeSlab(tbx,0,domain.smallEnd(0)), [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
74 velx(i,j,k) = momx(i,j,k) / dens_arr(i-1,j,k,
Rho_comp);
77 if ( (bx.bigEnd(0) == domain.bigEnd(0)) &&
79 ParallelFor(makeSlab(tbx,0,domain.bigEnd(0)+1), [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
80 velx(i,j,k) = momx(i,j,k) / dens_arr(i,j,k,
Rho_comp);
83 if ( (bx.smallEnd(1) == domain.smallEnd(1)) &&
85 ParallelFor(makeSlab(tby,1,domain.smallEnd(1)), [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
86 vely(i,j,k) = momy(i,j,k) / dens_arr(i,j-1,k,
Rho_comp);
89 if ( (bx.bigEnd(1) == domain.bigEnd(1)) &&
91 ParallelFor(makeSlab(tby,1,domain.bigEnd(1)+1), [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
92 vely(i,j,k) = momy(i,j,k) / dens_arr(i,j,k,
Rho_comp);
#define Rho_comp
Definition: ERF_IndexDefines.H:36
void MomentumToVelocity(MultiFab &xvel, MultiFab &yvel, MultiFab &zvel, const MultiFab &density, const MultiFab &xmom_in, const MultiFab &ymom_in, const MultiFab &zmom_in, const Box &domain, const Vector< BCRec > &domain_bcs_type_h)
Definition: ERF_MomentumToVelocity.cpp:25
@ cons_bc
Definition: ERF_IndexDefines.H:76
@ ext_dir
Definition: ERF_IndexDefines.H:180
@ xvel
Definition: ERF_IndexDefines.H:130
@ zvel
Definition: ERF_IndexDefines.H:132
@ yvel
Definition: ERF_IndexDefines.H:131