ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_ConvertForProjection.cpp File Reference
#include <AMReX.H>
#include <AMReX_MultiFab.H>
#include <ERF_Utils.H>
Include dependency graph for ERF_ConvertForProjection.cpp:

Functions

void ConvertForProjection (const MultiFab &den_div, const MultiFab &den_mlt, MultiFab &xmom, MultiFab &ymom, MultiFab &zmom, const Box &domain, const Vector< BCRec > &domain_bcs_type_h)
 

Function Documentation

◆ ConvertForProjection()

void ConvertForProjection ( const MultiFab &  den_div,
const MultiFab &  den_mlt,
MultiFab &  xmom,
MultiFab &  ymom,
MultiFab &  zmom,
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
28 {
29  BL_PROFILE_VAR("ConvertForProjection()",onvertForProjection);
30 
31  const BCRec* bc_ptr_h = domain_bcs_type_h.data();
32 
33 #ifdef _OPENMP
34 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
35 #endif
36  for ( MFIter mfi(den_div,TilingIfNotGPU()); mfi.isValid(); ++mfi)
37  {
38  // We need velocity in the interior ghost cells (init == real)
39  Box bx = mfi.tilebox();
40 
41  Box tbx = surroundingNodes(bx,0);
42  Box tby = surroundingNodes(bx,1);
43  Box tbz = surroundingNodes(bx,2);
44 
45  if ( (bx.smallEnd(0) == domain.smallEnd(0)) &&
46  ( (bc_ptr_h[BCVars::cons_bc].lo(0) == ERFBCType::ext_dir) ||
47  (bc_ptr_h[BCVars::cons_bc].lo(0) == ERFBCType::ext_dir_upwind) ) )
48  {
49  tbx.growLo(0,-1);
50  }
51  if ( (bx.bigEnd(0) == domain.bigEnd(0)) &&
52  ( (bc_ptr_h[BCVars::cons_bc].hi(0) == ERFBCType::ext_dir) ||
53  (bc_ptr_h[BCVars::cons_bc].hi(0) == ERFBCType::ext_dir_upwind) ) )
54  {
55  tbx.growHi(0,-1);
56  }
57  if ( (bx.smallEnd(1) == domain.smallEnd(1)) &&
58  ( (bc_ptr_h[BCVars::cons_bc].lo(1) == ERFBCType::ext_dir) ||
59  (bc_ptr_h[BCVars::cons_bc].lo(1) == ERFBCType::ext_dir_upwind) ) )
60  {
61  tby.growLo(1,-1);
62  }
63  if ( (bx.bigEnd(1) == domain.bigEnd(1)) &&
64  ( (bc_ptr_h[BCVars::cons_bc].hi(1) == ERFBCType::ext_dir) ||
65  (bc_ptr_h[BCVars::cons_bc].hi(1) == ERFBCType::ext_dir_upwind) ) )
66  {
67  tby.growHi(1,-1);
68  }
69 
70  // Conserved variables on cell centers -- we use this for density
71  const Array4<const Real>& den_div_arr = den_div.const_array(mfi);
72  const Array4<const Real>& den_mlt_arr = den_mlt.const_array(mfi);
73 
74  // Momentum on faces
75  Array4<Real> const& momx = xmom.array(mfi);
76  Array4<Real> const& momy = ymom.array(mfi);
77  Array4<Real> const& momz = zmom.array(mfi);
78 
79  ParallelFor(tbx, tby, tbz,
80  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
81  momx(i,j,k) *= ( den_mlt_arr(i,j,k,Rho_comp) + den_mlt_arr(i-1,j,k,Rho_comp) )
82  / ( den_div_arr(i,j,k,Rho_comp) + den_div_arr(i-1,j,k,Rho_comp) );
83  },
84  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
85  momy(i,j,k) *= ( den_mlt_arr(i,j,k,Rho_comp) + den_mlt_arr(i,j-1,k,Rho_comp) )
86  / ( den_div_arr(i,j,k,Rho_comp) + den_div_arr(i,j-1,k,Rho_comp) );
87  },
88  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
89  momz(i,j,k) *= ( den_mlt_arr(i,j,k,Rho_comp) + den_mlt_arr(i,j,k-1,Rho_comp) )
90  / ( den_div_arr(i,j,k,Rho_comp) + den_div_arr(i,j,k-1,Rho_comp) );
91  });
92 
93  if (bx.smallEnd(0) == domain.smallEnd(0)) {
94  if (bc_ptr_h[BCVars::cons_bc].lo(0) == ERFBCType::ext_dir)
95  {
96  ParallelFor(makeSlab(tbx,0,domain.smallEnd(0)), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
97  momx(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
98  });
99  }
100  else if (bc_ptr_h[BCVars::cons_bc].lo(0) == ERFBCType::ext_dir_upwind)
101  {
102  ParallelFor(makeSlab(tbx,0,domain.smallEnd(0)), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
103  if (momx(i,j,k) >= 0.) {
104  momx(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
105  } else {
106  momx(i,j,k) *= ( den_mlt_arr(i,j,k,Rho_comp) + den_mlt_arr(i-1,j,k,Rho_comp) )
107  / ( den_div_arr(i,j,k,Rho_comp) + den_div_arr(i-1,j,k,Rho_comp) );
108  }
109  });
110  }
111  }
112 
113  if (bx.bigEnd(0) == domain.bigEnd(0)) {
114  if (bc_ptr_h[BCVars::cons_bc].hi(0) == ERFBCType::ext_dir)
115  {
116  ParallelFor(makeSlab(tbx,0,domain.bigEnd(0)+1), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
117  momx(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
118  });
119  }
120  else if (bc_ptr_h[BCVars::cons_bc].hi(0) == ERFBCType::ext_dir_upwind)
121  {
122  ParallelFor(makeSlab(tbx,0,domain.smallEnd(0)), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
123  if (momx(i,j,k) <= 0.) {
124  momx(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
125  } else {
126  momx(i,j,k) *= ( den_mlt_arr(i,j,k,Rho_comp) + den_mlt_arr(i-1,j,k,Rho_comp) )
127  / ( den_div_arr(i,j,k,Rho_comp) + den_div_arr(i-1,j,k,Rho_comp) );
128  }
129  });
130  }
131  }
132 
133  if (bx.smallEnd(1) == domain.smallEnd(1)) {
134  if (bc_ptr_h[BCVars::cons_bc].lo(1) == ERFBCType::ext_dir)
135  {
136  ParallelFor(makeSlab(tby,1,domain.smallEnd(1)), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
137  momy(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
138  });
139  }
140  else if (bc_ptr_h[BCVars::cons_bc].lo(1) == ERFBCType::ext_dir_upwind)
141  {
142  ParallelFor(makeSlab(tby,1,domain.smallEnd(1)), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
143  if (momy(i,j,k) >= 0.) {
144  momy(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
145  } else {
146  momy(i,j,k) *= ( den_mlt_arr(i,j,k,Rho_comp) + den_mlt_arr(i-1,j,k,Rho_comp) )
147  / ( den_div_arr(i,j,k,Rho_comp) + den_div_arr(i-1,j,k,Rho_comp) );
148  }
149  });
150  }
151  }
152 
153  if (bx.bigEnd(1) == domain.bigEnd(1)) {
154  if (bc_ptr_h[BCVars::cons_bc].hi(1) == ERFBCType::ext_dir)
155  {
156  ParallelFor(makeSlab(tby,1,domain.bigEnd(1)+1), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
157  momy(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
158  });
159  }
160  else if (bc_ptr_h[BCVars::cons_bc].hi(1) == ERFBCType::ext_dir_upwind)
161  {
162  ParallelFor(makeSlab(tby,1,domain.bigEnd(1)+1), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
163  if (momy(i,j,k) <= 0.) {
164  momy(i,j,k) *= den_mlt_arr(i-1,j,k,Rho_comp) / den_div_arr(i-1,j,k,Rho_comp) ;
165  } else {
166  momy(i,j,k) *= ( den_mlt_arr(i,j,k,Rho_comp) + den_mlt_arr(i-1,j,k,Rho_comp) )
167  / ( den_div_arr(i,j,k,Rho_comp) + den_div_arr(i-1,j,k,Rho_comp) );
168  }
169  });
170  }
171  }
172 
173 
174  } // end MFIter
175 }
#define Rho_comp
Definition: ERF_IndexDefines.H:36
@ cons_bc
Definition: ERF_IndexDefines.H:76
@ ext_dir
Definition: ERF_IndexDefines.H:191
@ ext_dir_upwind
Definition: ERF_IndexDefines.H:199
@ ymom
Definition: ERF_IndexDefines.H:152
@ zmom
Definition: ERF_IndexDefines.H:153
@ xmom
Definition: ERF_IndexDefines.H:151

Referenced by ERF::project_velocities().

Here is the caller graph for this function: