ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Utils.H
Go to the documentation of this file.
1 #ifndef ERF_UTILS_H_
2 #define ERF_UTILS_H_
3 
4 #include "AMReX.H"
5 #include "AMReX_MultiFab.H"
6 #include "AMReX_BCRec.H"
7 #include "ERF_DataStruct.H"
8 #include "ERF_IndexDefines.H"
9 #include "ERF_SurfaceLayer.H"
10 #include "ERF_FillPatcher.H"
11 
12 /*
13  * Create a new BoxArray in which every grid touches the bottom boundary
14  */
15 void ChopGrids2D (amrex::BoxArray& ba, const amrex::Box& domain, int target_size);
16 
17 /*
18  * Create a new BoxArray with exactly the number of boxes as MPI ranks.
19  * This is designed to work only on the base level which covers the entire domain.
20  */
21 amrex::BoxArray
22 ERFPostProcessBaseGrids (const amrex::Box& domain, bool decompose_in_z);
23 
24 /*
25  * Create the Jacobian for the metric transformation when use_terrain is true
26  */
27 void make_J (const amrex::Geometry& geom,
28  amrex::MultiFab& z_phys_nd,
29  amrex::MultiFab& detJ_cc);
30 
31 void make_areas (const amrex::Geometry& geom,
32  amrex::MultiFab& z_phys_nd,
33  amrex::MultiFab& ax,
34  amrex::MultiFab& ay,
35  amrex::MultiFab& az);
36 
37 /*
38  * Average z_phys_nd on nodes to cell centers
39  */
40 void make_zcc (const amrex::Geometry& geom,
41  amrex::MultiFab& z_phys_nd,
42  amrex::MultiFab& z_phys_cc);
43 
44 /*
45  * Convert momentum to velocity by dividing by density averaged onto faces
46  */
47 void MomentumToVelocity (amrex::MultiFab& xvel_out,
48  amrex::MultiFab& yvel_out,
49  amrex::MultiFab& zvel_out,
50  const amrex::MultiFab& cons_in,
51  const amrex::MultiFab& xmom_in,
52  const amrex::MultiFab& ymom_in,
53  const amrex::MultiFab& zmom_in,
54  const amrex::Box& domain,
55  const amrex::Vector<amrex::BCRec>& domain_bcs_type_h,
56  const amrex::MultiFab* c_vfrac = nullptr // optional
57  );
58 
59 /*
60  * Convert velocity to momentum by multiplying by density averaged onto faces
61  */
62 void VelocityToMomentum (const amrex::MultiFab& xvel_in,
63  const amrex::IntVect & xvel_ngrow,
64  const amrex::MultiFab& yvel_in,
65  const amrex::IntVect & yvel_ngrow,
66  const amrex::MultiFab& zvel_in,
67  const amrex::IntVect & zvel_ngrow,
68  const amrex::MultiFab& cons_in,
69  amrex::MultiFab& xmom_out,
70  amrex::MultiFab& ymom_out,
71  amrex::MultiFab& zmom_out,
72  const amrex::Box& domain,
73  const amrex::Vector<amrex::BCRec>& domain_bcs_type_h,
74  const amrex::MultiFab* c_vfrac = nullptr // optional
75  );
76 
77 /*
78  * Convert (den_div u) to (den_mlt u) by multiplying by (den_mlt/den_div)
79  */
80 void
81 ConvertForProjection (const amrex::MultiFab& den_div, const amrex::MultiFab& den_mlt,
82  amrex::MultiFab& xmom, amrex::MultiFab& ymom, amrex::MultiFab& zmom,
83  const amrex::Box& domain, const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
84 
85 /*
86  * \brief Enforces solvability by scaling outflow to match with inflow.
87  *
88  */
90  amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& vels_vec,
91  amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& area_vec,
92  const amrex::Geometry& geom);
93 
94 /*
95  * Compute boxes for looping over interior/exterior ghost cells
96  * for use by fillpatch, erf_slow_rhs_pre, and erf_slow_rhs_post
97  */
98 void realbdy_interior_bxs_xy (const amrex::Box& bx,
99  const amrex::Box& domain,
100  const int& width,
101  amrex::Box& bx_xlo,
102  amrex::Box& bx_xhi,
103  amrex::Box& bx_ylo,
104  amrex::Box& bx_yhi,
105  const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0),
106  const bool get_int_ng=false);
107 
108 /*
109  * Compute relaxation region RHS with wrfbdy
110  */
112  const amrex::Real& delta_t,
113  const amrex::Real& start_bdy_time,
114  const amrex::Real& final_bdy_time,
115  const amrex::Real& bdy_time_interval,
116  const amrex::Real& nudge_factor,
117  int width,
118  bool do_upwind,
119  const amrex::Geometry& geom,
120  amrex::Vector<amrex::MultiFab>& S_rhs,
121  amrex::Vector<amrex::MultiFab>& S_cur_data,
122  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xlo,
123  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xhi,
124  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_ylo,
125  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_yhi);
126 
127 /*
128  * Compute relaxation region RHS at fine-crse interface
129  */
130 void
132  const amrex::Real& delta_t,
133  const int& width,
134  const int& set_width,
135  const amrex::Geometry& geom,
136  ERFFillPatcher* FPr_c,
137  ERFFillPatcher* FPr_u,
138  ERFFillPatcher* FPr_v,
139  ERFFillPatcher* FPr_w,
140  amrex::Vector<amrex::BCRec>& domain_bcs_type,
141  amrex::Vector<amrex::MultiFab>& S_rhs_f,
142  amrex::Vector<amrex::MultiFab>& S_data_f);
143 
144 /*
145  * Accumulate time averaged velocity fields
146  */
147 void
149  amrex::Real& t_avg_cnt,
150  amrex::MultiFab* vel_t_avg,
151  amrex::MultiFab& xvel,
152  amrex::MultiFab& yvel,
153  amrex::MultiFab& zvel);
154 
155 /**
156  * Compute the nudging relaxation
157  *
158  * @param[in] delta_t time step
159  * @param[in] icomp component offset
160  * @param[in] num_var number of variables to loop
161  * @param[in] width width of wrf bdy file
162  * @param[in] dom_lo low bound of domain
163  * @param[in] dom_hi high bound of domain
164  * @param[in] F1 drift relaxation parameter
165  * @param[in] F2 Laplacian relaxation parameter
166  * @param[in] bx_xlo box for low x relaxation
167  * @param[in] bx_xhi box for high x relaxation
168  * @param[in] bx_ylo box for low y relaxation
169  * @param[in] bx_yhi box for high y relaxation
170  * @param[in] arr_xlo array for low x relaxation
171  * @param[in] arr_xhi array for high x relaxation
172  * @param[in] arr_ylo array for low y relaxation
173  * @param[in] arr_yhi array for high y relaxation
174  * @param[in] data_arr data array
175  * @param[out] rhs_arr RHS array
176  */
177 AMREX_GPU_HOST
178 AMREX_FORCE_INLINE
179 void
180 realbdy_compute_relaxation (const int& icomp,
181  const int& num_var,
182  const int& width,
183  const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& dx,
184  const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbLo,
185  const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbHi,
186  const amrex::Real& F1,
187  const amrex::Box& domain_cc,
188  const amrex::Box& bx_xlo,
189  const amrex::Box& bx_xhi,
190  const amrex::Box& bx_ylo,
191  const amrex::Box& bx_yhi,
192  const amrex::Array4<const amrex::Real>& arr_xlo,
193  const amrex::Array4<const amrex::Real>& arr_xhi,
194  const amrex::Array4<const amrex::Real>& arr_ylo,
195  const amrex::Array4<const amrex::Real>& arr_yhi,
196  const amrex::Array4<const amrex::Real>& u_xlo,
197  const amrex::Array4<const amrex::Real>& u_xhi,
198  const amrex::Array4<const amrex::Real>& v_xlo,
199  const amrex::Array4<const amrex::Real>& v_xhi,
200  const amrex::Array4<const amrex::Real>& v_ylo,
201  const amrex::Array4<const amrex::Real>& v_yhi,
202  const amrex::Array4<const amrex::Real>& data_arr,
203  const amrex::Array4<amrex::Real>& rhs_arr,
204  const bool& do_upwind)
205 {
206  amrex::IntVect iv = bx_xlo.type();
207  amrex::Real ioff = (iv[0]==1) ? 0.0 : 0.5;
208  amrex::Real joff = (iv[1]==1) ? 0.0 : 0.5;
209  const auto& dom_cc_lo = lbound(domain_cc);
210  const auto& dom_cc_hi = ubound(domain_cc);
211  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
212  {
213  // Corners with x boxes
214  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
215  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
216  amrex::Real x_end = ProbLo[0] + width * dx[0];
217  amrex::Real y_end = ProbLo[1] + width * dx[1];
218  amrex::Real y_strt = ProbHi[1] - width * dx[1];
219  amrex::Real xi = (x_end - x) / (x_end - ProbLo[0]);
220  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
221  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
222  amrex::Real eta = std::max(eta_lo,eta_hi);
223  amrex::Real Factor = std::max(xi*xi,eta*eta);
224 
225  amrex::Real delta = arr_xlo(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
226  amrex::Real Temp = Factor*F1*delta;
227  if (!do_upwind) {
228  rhs_arr(i,j,k,n+icomp) += Temp;
229  } else {
230  int jju = std::min(std::max(j,dom_cc_lo.y),dom_cc_hi.y);
231  int iiv = std::min(std::max(i,dom_cc_lo.x),dom_cc_hi.x);
232  if ( (u_xlo(dom_cc_lo.x,jju,k) >= 0.0) ||
233  ((j == dom_cc_lo.y ) && (v_xlo(iiv,dom_cc_lo.y ,k) >= 0.0)) ||
234  ((j == dom_cc_hi.y+iv[1]) && (v_xlo(iiv,dom_cc_hi.y+1,k) <= 0.0)) ) {
235  rhs_arr(i,j,k,n+icomp) += Temp;
236  }
237  }
238  },
239  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
240  {
241  // Corners with x boxes
242  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
243  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
244  amrex::Real x_strt = ProbHi[0] - width * dx[0];
245  amrex::Real y_strt = ProbHi[1] - width * dx[1];
246  amrex::Real y_end = ProbLo[1] + width * dx[1];
247  amrex::Real xi = (x - x_strt) / (ProbHi[0] - x_strt);
248  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
249  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
250  amrex::Real eta = std::max(eta_lo,eta_hi);
251  amrex::Real Factor = std::max(xi*xi,eta*eta);
252 
253  amrex::Real delta = arr_xhi(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
254  amrex::Real Temp = Factor*F1*delta;
255  if (!do_upwind) {
256  rhs_arr(i,j,k,n+icomp) += Temp;
257  } else {
258  int jju = std::min(std::max(j,dom_cc_lo.y),dom_cc_hi.y);
259  int iiv = std::min(std::max(i,dom_cc_lo.x),dom_cc_hi.x);
260  if ( (u_xhi(dom_cc_hi.x+1,jju,k) <= 0.0) ||
261  ((j == dom_cc_lo.y ) && (v_xhi(iiv,dom_cc_lo.y ,k) >= 0.0)) ||
262  ((j == dom_cc_hi.y+iv[1]) && (v_xhi(iiv,dom_cc_hi.y+1,k) <= 0.0)) ) {
263  rhs_arr(i,j,k,n+icomp) += Temp;
264  }
265  }
266  });
267 
268  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
269  {
270  // No corners for y boxes
271  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
272  amrex::Real y_end = ProbLo[1] + width * dx[1];
273  amrex::Real eta = (y_end - y) / (y_end - ProbLo[1]);
274  amrex::Real Factor = eta*eta;
275 
276  amrex::Real delta = arr_ylo(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
277  amrex::Real Temp = Factor*F1*delta;
278  if (!do_upwind) {
279  rhs_arr(i,j,k,n+icomp) += Temp;
280  } else {
281  int iiv = std::min(std::max(i,dom_cc_lo.x+width),dom_cc_hi.x-width);
282  if (v_ylo(iiv,dom_cc_lo.y,k) >= 0.0) {
283  rhs_arr(i,j,k,n+icomp) += Temp;
284  }
285  }
286  },
287  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
288  {
289  // No corners for y boxes
290  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
291  amrex::Real y_strt = ProbHi[1] - width * dx[1];
292  amrex::Real eta = (y - y_strt) / (ProbHi[1] - y_strt);
293  amrex::Real Factor = eta*eta;
294 
295  amrex::Real delta = arr_yhi(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
296  amrex::Real Temp = Factor*F1*delta;
297  if (!do_upwind) {
298  rhs_arr(i,j,k,n+icomp) += Temp;
299  } else {
300  int iiv = std::min(std::max(i,dom_cc_lo.x+width),dom_cc_hi.x-width);
301  if (v_yhi(iiv,dom_cc_hi.y+1,k) >= 0.0) {
302  rhs_arr(i,j,k,n+icomp) += Temp;
303  }
304  }
305  });
306 }
307 
308 /*
309  * Effectively a Multiply for a MultiFab and an iMultiFab mask
310  */
311 AMREX_GPU_HOST
312 AMREX_FORCE_INLINE
313 void
314 ApplyMask (amrex::MultiFab& dst,
315  const amrex::iMultiFab& imask,
316  const int nghost = 0)
317 {
318  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
319  {
320  const amrex::Box& bx = mfi.growntilebox(nghost);
321  if (bx.ok())
322  {
323  auto dstFab = dst.array(mfi);
324  const auto maskFab = imask.const_array(mfi);
325  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
326  {
327  dstFab(i,j,k) *= maskFab(i,j,k);
328  });
329  }
330  }
331 }
332 
333 AMREX_GPU_HOST
334 AMREX_FORCE_INLINE
335 void
336 ApplyInvertedMask (amrex::MultiFab& dst,
337  const amrex::iMultiFab& imask,
338  const int nghost = 0)
339 {
340  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
341  {
342  const amrex::Box& bx = mfi.growntilebox(nghost);
343  if (bx.ok())
344  {
345  auto dstFab = dst.array(mfi);
346  const auto maskFab = imask.const_array(mfi);
347  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
348  {
349  dstFab(i,j,k) *= (1-maskFab(i,j,k));
350  });
351  }
352  }
353 }
354 
355 void
356 thinbody_wall_dist (std::unique_ptr<amrex::MultiFab>& wdist,
357  amrex::Vector<amrex::IntVect>& xfaces,
358  amrex::Vector<amrex::IntVect>& yfaces,
359  amrex::Vector<amrex::IntVect>& zfaces,
360  const amrex::Geometry& geomdata,
361  std::unique_ptr<amrex::MultiFab>& z_phys_cc);
362 
363 
365 
366 AMREX_GPU_HOST_DEVICE
367 AMREX_FORCE_INLINE
371  amrex::Real z1, amrex::Real p1,
372  amrex::Real z2, amrex::Real p2)
373 {
374  return p0 * ( (z - z1) * (z - z2) ) / ( (z0 - z1) * (z0 - z2) )
375  + p1 * ( (z - z0) * (z - z2) ) / ( (z1 - z0) * (z1 - z2) )
376  + p2 * ( (z - z0) * (z - z1) ) / ( (z2 - z0) * (z2 - z1) );
377 }
378 
379 #endif
const auto dx
Definition: ERF_InitCustomPertVels_ParticleTests.H:15
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+0.5) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
amrex::Real Real
Definition: ERF_ShocInterface.H:19
void MomentumToVelocity(amrex::MultiFab &xvel_out, amrex::MultiFab &yvel_out, amrex::MultiFab &zvel_out, const amrex::MultiFab &cons_in, const amrex::MultiFab &xmom_in, const amrex::MultiFab &ymom_in, const amrex::MultiFab &zmom_in, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h, const amrex::MultiFab *c_vfrac=nullptr)
void VelocityToMomentum(const amrex::MultiFab &xvel_in, const amrex::IntVect &xvel_ngrow, const amrex::MultiFab &yvel_in, const amrex::IntVect &yvel_ngrow, const amrex::MultiFab &zvel_in, const amrex::IntVect &zvel_ngrow, const amrex::MultiFab &cons_in, amrex::MultiFab &xmom_out, amrex::MultiFab &ymom_out, amrex::MultiFab &zmom_out, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h, const amrex::MultiFab *c_vfrac=nullptr)
void realbdy_compute_interior_ghost_rhs(const amrex::Real &total_time, const amrex::Real &delta_t, const amrex::Real &start_bdy_time, const amrex::Real &final_bdy_time, const amrex::Real &bdy_time_interval, const amrex::Real &nudge_factor, int width, bool do_upwind, const amrex::Geometry &geom, amrex::Vector< amrex::MultiFab > &S_rhs, amrex::Vector< amrex::MultiFab > &S_cur_data, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_xlo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_xhi, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_ylo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_yhi)
void WeatherDataInterpolation(const amrex::Real time)
void fine_compute_interior_ghost_rhs(const amrex::Real &time, const amrex::Real &delta_t, const int &width, const int &set_width, const amrex::Geometry &geom, ERFFillPatcher *FPr_c, ERFFillPatcher *FPr_u, ERFFillPatcher *FPr_v, ERFFillPatcher *FPr_w, amrex::Vector< amrex::BCRec > &domain_bcs_type, amrex::Vector< amrex::MultiFab > &S_rhs_f, amrex::Vector< amrex::MultiFab > &S_data_f)
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_compute_relaxation(const int &icomp, const int &num_var, const int &width, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dx, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbLo, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbHi, const amrex::Real &F1, const amrex::Box &domain_cc, const amrex::Box &bx_xlo, const amrex::Box &bx_xhi, const amrex::Box &bx_ylo, const amrex::Box &bx_yhi, const amrex::Array4< const amrex::Real > &arr_xlo, const amrex::Array4< const amrex::Real > &arr_xhi, const amrex::Array4< const amrex::Real > &arr_ylo, const amrex::Array4< const amrex::Real > &arr_yhi, const amrex::Array4< const amrex::Real > &u_xlo, const amrex::Array4< const amrex::Real > &u_xhi, const amrex::Array4< const amrex::Real > &v_xlo, const amrex::Array4< const amrex::Real > &v_xhi, const amrex::Array4< const amrex::Real > &v_ylo, const amrex::Array4< const amrex::Real > &v_yhi, const amrex::Array4< const amrex::Real > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr, const bool &do_upwind)
Definition: ERF_Utils.H:180
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real quad_interp_1d(amrex::Real z, amrex::Real z0, amrex::Real p0, amrex::Real z1, amrex::Real p1, amrex::Real z2, amrex::Real p2)
Definition: ERF_Utils.H:369
void thinbody_wall_dist(std::unique_ptr< amrex::MultiFab > &wdist, amrex::Vector< amrex::IntVect > &xfaces, amrex::Vector< amrex::IntVect > &yfaces, amrex::Vector< amrex::IntVect > &zfaces, const amrex::Geometry &geomdata, std::unique_ptr< amrex::MultiFab > &z_phys_cc)
amrex::BoxArray ERFPostProcessBaseGrids(const amrex::Box &domain, bool decompose_in_z)
void make_zcc(const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc)
void realbdy_interior_bxs_xy(const amrex::Box &bx, const amrex::Box &domain, const int &width, amrex::Box &bx_xlo, amrex::Box &bx_xhi, amrex::Box &bx_ylo, amrex::Box &bx_yhi, const amrex::IntVect &ng_vect=amrex::IntVect(0, 0, 0), const bool get_int_ng=false)
void Time_Avg_Vel_atCC(const amrex::Real &dt, amrex::Real &t_avg_cnt, amrex::MultiFab *vel_t_avg, amrex::MultiFab &xvel, amrex::MultiFab &yvel, amrex::MultiFab &zvel)
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask(amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
Definition: ERF_Utils.H:336
void make_J(const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &detJ_cc)
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyMask(amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
Definition: ERF_Utils.H:314
void ConvertForProjection(const amrex::MultiFab &den_div, const amrex::MultiFab &den_mlt, amrex::MultiFab &xmom, amrex::MultiFab &ymom, amrex::MultiFab &zmom, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h)
void ChopGrids2D(amrex::BoxArray &ba, const amrex::Box &domain, int target_size)
void enforceInOutSolvability(int lev, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &vels_vec, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &area_vec, const amrex::Geometry &geom)
void make_areas(const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &ax, amrex::MultiFab &ay, amrex::MultiFab &az)
Definition: ERF_FillPatcher.H:9
@ ymom
Definition: ERF_IndexDefines.H:160
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
@ xvel
Definition: ERF_IndexDefines.H:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40