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