ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
MomentumToVelocity.cpp File Reference
#include <AMReX.H>
#include <AMReX_MultiFab.H>
#include <Utils.H>
Include dependency graph for MomentumToVelocity.cpp:

Functions

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)
 

Function Documentation

◆ MomentumToVelocity()

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 
)

Convert momentum to velocity by dividing by density averaged onto faces

Parameters
[out]xvelx-component of velocity
[out]yvely-component of velocity
[out]zvelz-component of velocity
[in]densitydensity at cell centers
[in]xmom_inx-component of momentum
[in]ymom_iny-component of momentum
[in]zmom_inz-component of momentum
[in]domainDomain at this level
[in]domain_bcs_type_hhost vector for domain boundary conditions
30 {
31  BL_PROFILE_VAR("MomentumToVelocity()",MomentumToVelocity);
32 
33  const BCRec* bc_ptr_h = domain_bcs_type_h.data();
34 
35 #ifdef _OPENMP
36 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
37 #endif
38  for ( MFIter mfi(density,TilingIfNotGPU()); mfi.isValid(); ++mfi)
39  {
40  // We need velocity in the interior ghost cells (init == real)
41  Box bx = mfi.tilebox();
42 
43  const Box& tbx = surroundingNodes(bx,0);
44  const Box& tby = surroundingNodes(bx,1);
45  const Box& tbz = surroundingNodes(bx,2);
46 
47  // Conserved variables on cell centers -- we use this for density
48  const Array4<const Real>& dens_arr = density.array(mfi);
49 
50  // Momentum on faces
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);
54 
55  // Velocity on faces
56  const Array4<Real>& velx = xvel.array(mfi);
57  const Array4<Real>& vely = yvel.array(mfi);
58  const Array4<Real>& velz = zvel.array(mfi);
59 
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));
63  },
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));
66  },
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));
69  });
70 
71  if ( (bx.smallEnd(0) == domain.smallEnd(0)) &&
72  (bc_ptr_h[BCVars::cons_bc].lo(0) == ERFBCType::ext_dir) ) {
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);
75  });
76  }
77  if ( (bx.bigEnd(0) == domain.bigEnd(0)) &&
78  (bc_ptr_h[BCVars::cons_bc].hi(0) == ERFBCType::ext_dir) ) {
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);
81  });
82  }
83  if ( (bx.smallEnd(1) == domain.smallEnd(1)) &&
84  (bc_ptr_h[BCVars::cons_bc].lo(1) == ERFBCType::ext_dir) ) {
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);
87  });
88  }
89  if ( (bx.bigEnd(1) == domain.bigEnd(1)) &&
90  (bc_ptr_h[BCVars::cons_bc].hi(1) == ERFBCType::ext_dir) ) {
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);
93  });
94  }
95  } // end MFIter
96 }
#define Rho_comp
Definition: IndexDefines.H:12
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: MomentumToVelocity.cpp:25
@ cons_bc
Definition: IndexDefines.H:47
@ ext_dir
Definition: IndexDefines.H:152
@ xvel
Definition: IndexDefines.H:100
@ zvel
Definition: IndexDefines.H:102
@ yvel
Definition: IndexDefines.H:101

Referenced by ERF::AverageDownTo(), ERF::FillCoarsePatch(), ERF::FillIntermediatePatch(), and ERF::FillPatch().

Here is the caller graph for this function: