ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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_ABLMost.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 
57 /*
58  * Convert velocity to momentum by multiplying by density averaged onto faces
59  */
60 void VelocityToMomentum (const amrex::MultiFab& xvel_in,
61  const amrex::IntVect & xvel_ngrow,
62  const amrex::MultiFab& yvel_in,
63  const amrex::IntVect & yvel_ngrow,
64  const amrex::MultiFab& zvel_in,
65  const amrex::IntVect & zvel_ngrow,
66  const amrex::MultiFab& cons_in,
67  amrex::MultiFab& xmom_out,
68  amrex::MultiFab& ymom_out,
69  amrex::MultiFab& zmom_out,
70  const amrex::Box& domain,
71  const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
72 
73 /*
74  * Convert (den_div u) to (den_mlt u) by multiplying by (den_mlt/den_div)
75  */
76 void
77 ConvertForProjection (const amrex::MultiFab& den_div, const amrex::MultiFab& den_mlt,
78  amrex::MultiFab& xmom, amrex::MultiFab& ymom, amrex::MultiFab& zmom,
79  const amrex::Box& domain, const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
80 
81 /*
82  * Compute boxes for looping over interior/exterior ghost cells
83  * for use by fillpatch, erf_slow_rhs_pre, and erf_slow_rhs_post
84  */
85 void realbdy_interior_bxs_xy (const amrex::Box& bx,
86  const amrex::Box& domain,
87  const int& width,
88  amrex::Box& bx_xlo,
89  amrex::Box& bx_xhi,
90  amrex::Box& bx_ylo,
91  amrex::Box& bx_yhi,
92  const int& set_width=0,
93  const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0),
94  const bool get_int_ng=false);
95 
96 /*
97  * Compute boxes for looping over set region cells
98  * for use by fillpatch, erf_slow_rhs_pre, and erf_slow_rhs_post
99  */
100 void realbdy_bc_bxs_xy (const amrex::Box& bx,
101  const amrex::Box& domain,
102  const int& set_width,
103  amrex::Box& bx_xlo,
104  amrex::Box& bx_xhi,
105  amrex::Box& bx_ylo,
106  amrex::Box& bx_yhi,
107  const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0));
108 
109 /*
110  * Compute relaxation region RHS with wrfbdy
111  */
112 void realbdy_compute_interior_ghost_rhs (const amrex::Real& bdy_time_interval,
113  const amrex::Real& start_bdy_time,
114  const amrex::Real& time,
115  const amrex::Real& delta_t,
116  int width,
117  int set_width,
118  const amrex::Geometry& geom,
119  amrex::Vector<amrex::MultiFab>& S_rhs,
120  amrex::Vector<amrex::MultiFab>& S_old_data,
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
131 fine_compute_interior_ghost_rhs (const amrex::Real& time,
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
148 Time_Avg_Vel_atCC (const amrex::Real& dt,
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  * Zero RHS in the set region
157  *
158  * @param[in] icomp component offset
159  * @param[in] num_var number of variables to loop
160  * @param[in] bx_xlo box for low x relaxation
161  * @param[in] bx_xhi box for high x relaxation
162  * @param[in] bx_ylo box for low y relaxation
163  * @param[in] bx_yhi box for high y relaxation
164  * @param[out] rhs_arr RHS array
165  */
166 AMREX_GPU_HOST
167 AMREX_FORCE_INLINE
168 void
169 realbdy_set_rhs_in_spec_region (const amrex::Real& dt,
170  const int& icomp,
171  const int& num_var,
172  const int& width,
173  const int& set_width_x,
174  const int& set_width_y,
175  const amrex::Dim3& dom_lo,
176  const amrex::Dim3& dom_hi,
177  const amrex::Box& bx_xlo,
178  const amrex::Box& bx_xhi,
179  const amrex::Box& bx_ylo,
180  const amrex::Box& bx_yhi,
181  const amrex::Array4<const amrex::Real>& arr_xlo,
182  const amrex::Array4<const amrex::Real>& arr_xhi,
183  const amrex::Array4<const amrex::Real>& arr_ylo,
184  const amrex::Array4<const amrex::Real>& arr_yhi,
185  const amrex::Array4<const amrex::Real>& data_arr,
186  const amrex::Array4<amrex::Real>& rhs_arr)
187 {
188  int Spec_z_x = set_width_x;
189  int Spec_z_y = set_width_y;
190  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
191  {
192  // Corners with x boxes
193  int j_lo = std::min(j-dom_lo.y,width-1);
194  int j_hi = std::min(dom_hi.y-j,width-1);
195  int jj = std::min(j_lo,j_hi);
196  int n_ind_y = jj + 1;
197  int n_ind_x = i - dom_lo.x + 1;
198  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
199  amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
200  rhs_arr(i,j,k,n+icomp) = tend;
201  }
202  },
203  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
204  {
205  // Corners with x boxes
206  int j_lo = std::min(j-dom_lo.y,width-1);
207  int j_hi = std::min(dom_hi.y-j,width-1);
208  int jj = std::min(j_lo,j_hi);
209  int n_ind_y = jj + 1;
210  int n_ind_x = dom_hi.x - i + 1;
211  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
212  amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
213  rhs_arr(i,j,k,n+icomp) = tend;
214  }
215  });
216 
217  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
218  {
219  // No corners for y boxes
220  int n_ind = j - dom_lo.y + 1;
221  if (n_ind <= Spec_z_y) {
222  amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
223  rhs_arr(i,j,k,n+icomp) = tend;
224  }
225  },
226  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
227  {
228  // No corners for y boxes
229  int n_ind = dom_hi.y - j + 1;
230  if (n_ind <= Spec_z_y) {
231  amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
232  rhs_arr(i,j,k,n+icomp) = tend;
233  }
234  });
235 }
236 
237 /**
238  * Compute the Laplacian RHS in the relaxation zone
239  *
240  * @param[in] delta_t time step
241  * @param[in] icomp component offset
242  * @param[in] num_var number of variables to loop
243  * @param[in] width width of wrf bdy file
244  * @param[in] set_width width the set region
245  * @param[in] dom_lo low bound of domain
246  * @param[in] dom_hi high bound of domain
247  * @param[in] F1 drift relaxation parameter
248  * @param[in] F2 Laplacian relaxation parameter
249  * @param[in] bx_xlo box for low x relaxation
250  * @param[in] bx_xhi box for high x relaxation
251  * @param[in] bx_ylo box for low y relaxation
252  * @param[in] bx_yhi box for high y relaxation
253  * @param[in] arr_xlo array for low x relaxation
254  * @param[in] arr_xhi array for high x relaxation
255  * @param[in] arr_ylo array for low y relaxation
256  * @param[in] arr_yhi array for high y relaxation
257  * @param[in] data_arr data array
258  * @param[out] rhs_arr RHS array
259  */
260 AMREX_GPU_HOST
261 AMREX_FORCE_INLINE
262 void
264  const int& num_var,
265  const int& width,
266  const int& set_width,
267  const amrex::Dim3& dom_lo,
268  const amrex::Dim3& dom_hi,
269  const amrex::Real& F1,
270  const amrex::Real& F2,
271  const amrex::Box& bx_xlo,
272  const amrex::Box& bx_xhi,
273  const amrex::Box& bx_ylo,
274  const amrex::Box& bx_yhi,
275  const amrex::Array4<const amrex::Real>& arr_xlo,
276  const amrex::Array4<const amrex::Real>& arr_xhi,
277  const amrex::Array4<const amrex::Real>& arr_ylo,
278  const amrex::Array4<const amrex::Real>& arr_yhi,
279  const amrex::Array4<const amrex::Real>& data_arr,
280  const amrex::Array4<amrex::Real>& rhs_arr)
281 {
282  // RHS computation
283  int Spec_z = set_width;
284  int Relax_z = width - Spec_z + 1;
285  amrex::Real num = amrex::Real(Spec_z + Relax_z);
286  amrex::Real denom = amrex::Real(Relax_z - 1);
287  amrex::Real SpecExp = -std::log(0.1) / amrex::Real(width - Spec_z);
288  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
289  {
290  // Corners with x boxes
291  int j_lo = std::min(j-dom_lo.y,width-1);
292  int j_hi = std::min(dom_hi.y-j,width-1);
293  int jj = std::min(j_lo,j_hi);
294  int n_ind = std::min(i-dom_lo.x,jj) + 1;
295  //AMREX_ASSERT(n_ind > Spec_z);
296  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
297  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
298  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
299  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
300  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
301  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
302  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
303  amrex::Real delta = arr_xlo(i ,j ,k,n) - d;
304  amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
305  amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
306  amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
307  amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
308  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
309  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
310  rhs_arr(i,j,k,n+icomp) += Temp;
311  },
312  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
313  {
314  // Corners with x boxes
315  int j_lo = std::min(j-dom_lo.y,width-1);
316  int j_hi = std::min(dom_hi.y-j,width-1);
317  int jj = std::min(j_lo,j_hi);
318  int n_ind = std::min(dom_hi.x-i,jj) + 1;
319  //AMREX_ASSERT(n_ind > Spec_z);
320  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
321  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
322  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
323  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
324  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
325  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
326  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
327  amrex::Real delta = arr_xhi(i ,j ,k,n) - d;
328  amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
329  amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
330  amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
331  amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
332  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
333  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
334  rhs_arr(i,j,k,n+icomp) += Temp;
335  });
336 
337  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
338  {
339  // No corners for y boxes
340  int n_ind = j - dom_lo.y + 1;
341  //AMREX_ASSERT(n_ind > Spec_z);
342  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
343  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
344  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
345  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
346  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
347  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
348  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
349  amrex::Real delta = arr_ylo(i ,j ,k,n) - d;
350  amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
351  amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
352  amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
353  amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
354  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
355  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
356  rhs_arr(i,j,k,n+icomp) += Temp;
357  },
358  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
359  {
360  // No corners for y boxes
361  int n_ind = dom_hi.y - j + 1;
362  //AMREX_ASSERT(n_ind > Spec_z);
363  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
364  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
365  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
366  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
367  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
368  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
369  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
370  amrex::Real delta = arr_yhi(i ,j ,k,n) - d;
371  amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
372  amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
373  amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
374  amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
375  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
376  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
377  rhs_arr(i,j,k,n+icomp) += Temp;
378  });
379 }
380 
381 /*
382  * Effectively a Multiply for a MultiFab and an iMultiFab mask
383  */
384 AMREX_GPU_HOST
385 AMREX_FORCE_INLINE
386 void
387 ApplyMask (amrex::MultiFab& dst,
388  const amrex::iMultiFab& imask,
389  const int nghost = 0)
390 {
391  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
392  {
393  const amrex::Box& bx = mfi.growntilebox(nghost);
394  if (bx.ok())
395  {
396  auto dstFab = dst.array(mfi);
397  const auto maskFab = imask.const_array(mfi);
398  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
399  {
400  dstFab(i,j,k) *= maskFab(i,j,k);
401  });
402  }
403  }
404 }
405 
406 AMREX_GPU_HOST
407 AMREX_FORCE_INLINE
408 void
409 ApplyInvertedMask (amrex::MultiFab& dst,
410  const amrex::iMultiFab& imask,
411  const int nghost = 0)
412 {
413  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
414  {
415  const amrex::Box& bx = mfi.growntilebox(nghost);
416  if (bx.ok())
417  {
418  auto dstFab = dst.array(mfi);
419  const auto maskFab = imask.const_array(mfi);
420  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
421  {
422  dstFab(i,j,k) *= (1-maskFab(i,j,k));
423  });
424  }
425  }
426 }
427 
428 void
429 thinbody_wall_dist (std::unique_ptr<amrex::MultiFab>& wdist,
430  amrex::Vector<amrex::IntVect>& xfaces,
431  amrex::Vector<amrex::IntVect>& yfaces,
432  amrex::Vector<amrex::IntVect>& zfaces,
433  const amrex::Geometry& geomdata,
434  std::unique_ptr<amrex::MultiFab>& z_phys_cc);
435 #endif
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_set_rhs_in_spec_region(const amrex::Real &dt, const int &icomp, const int &num_var, const int &width, const int &set_width_x, const int &set_width_y, const amrex::Dim3 &dom_lo, const amrex::Dim3 &dom_hi, 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:169
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)
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)
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 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)
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 int &set_width=0, const amrex::IntVect &ng_vect=amrex::IntVect(0, 0, 0), const bool get_int_ng=false)
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_compute_interior_ghost_rhs(const amrex::Real &bdy_time_interval, const amrex::Real &start_bdy_time, const amrex::Real &time, const amrex::Real &delta_t, int width, int set_width, const amrex::Geometry &geom, amrex::Vector< amrex::MultiFab > &S_rhs, amrex::Vector< amrex::MultiFab > &S_old_data, 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 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 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:409
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:387
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)
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_compute_laplacian_relaxation(const int &icomp, const int &num_var, const int &width, const int &set_width, const amrex::Dim3 &dom_lo, const amrex::Dim3 &dom_hi, const amrex::Real &F1, const amrex::Real &F2, 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:263
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:152
@ zmom
Definition: ERF_IndexDefines.H:153
@ xmom
Definition: ERF_IndexDefines.H:151
@ xvel
Definition: ERF_IndexDefines.H:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142