ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Utils.H File Reference
#include <AMReX.H>
#include <AMReX_MultiFab.H>
#include <AMReX_BCRec.H>
#include <ERF_DataStruct.H>
#include <ERF_IndexDefines.H>
#include <ERF_ABLMost.H>
#include <ERF_FillPatcher.H>
Include dependency graph for ERF_Utils.H:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

void ChopGrids2D (amrex::BoxArray &ba, const amrex::Box &domain, int target_size)
 
amrex::BoxArray ERFPostProcessBaseGrids (const amrex::Box &domain, bool decompose_in_z)
 
void make_J (const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &detJ_cc)
 
void make_areas (const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &ax, amrex::MultiFab &ay, amrex::MultiFab &az)
 
void make_zcc (const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc)
 
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 compute_interior_ghost_bxs_xy (const amrex::Box &bx, const amrex::Box &domain, const int &width, 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), const bool get_int_ng=false)
 
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 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 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 realbdy_set_rhs_in_spec_region (const amrex::Real &dt, 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::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)
 
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)
 
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyMask (amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
 
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask (amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
 

Function Documentation

◆ ApplyInvertedMask()

AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask ( amrex::MultiFab &  dst,
const amrex::iMultiFab &  imask,
const int  nghost = 0 
)
386 {
387  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
388  {
389  const amrex::Box& bx = mfi.growntilebox(nghost);
390  if (bx.ok())
391  {
392  auto dstFab = dst.array(mfi);
393  const auto maskFab = imask.const_array(mfi);
394  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
395  {
396  dstFab(i,j,k) *= (1-maskFab(i,j,k));
397  });
398  }
399  }
400 }

Referenced by add_thin_body_sources(), and ERF::project_velocities_tb().

Here is the caller graph for this function:

◆ ApplyMask()

AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyMask ( amrex::MultiFab &  dst,
const amrex::iMultiFab &  imask,
const int  nghost = 0 
)
365 {
366  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
367  {
368  const amrex::Box& bx = mfi.growntilebox(nghost);
369  if (bx.ok())
370  {
371  auto dstFab = dst.array(mfi);
372  const auto maskFab = imask.const_array(mfi);
373  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
374  {
375  dstFab(i,j,k) *= maskFab(i,j,k);
376  });
377  }
378  }
379 }

Referenced by ERF::FillIntermediatePatch(), and ERF::project_velocities_tb().

Here is the caller graph for this function:

◆ ChopGrids2D()

void ChopGrids2D ( amrex::BoxArray &  ba,
const amrex::Box &  domain,
int  target_size 
)

◆ compute_interior_ghost_bxs_xy()

void compute_interior_ghost_bxs_xy ( const amrex::Box &  bx,
const amrex::Box &  domain,
const int &  width,
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),
const bool  get_int_ng = false 
)

◆ ERFPostProcessBaseGrids()

amrex::BoxArray ERFPostProcessBaseGrids ( const amrex::Box &  domain,
bool  decompose_in_z 
)

◆ fine_compute_interior_ghost_rhs()

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 
)

◆ make_areas()

void make_areas ( const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::MultiFab &  ax,
amrex::MultiFab &  ay,
amrex::MultiFab &  az 
)

◆ make_J()

void make_J ( const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::MultiFab &  detJ_cc 
)

◆ make_zcc()

void make_zcc ( const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::MultiFab &  z_phys_cc 
)

◆ MomentumToVelocity()

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 
)

◆ realbdy_compute_interior_ghost_rhs()

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 
)

◆ realbdy_compute_laplacian_relaxation()

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 
)

Compute the Laplacian RHS in the relaxation zone

Parameters
[in]delta_ttime step
[in]icompcomponent offset
[in]num_varnumber of variables to loop
[in]widthwidth of wrf bdy file
[in]set_widthwidth the set region
[in]dom_lolow bound of domain
[in]dom_hihigh bound of domain
[in]F1drift relaxation parameter
[in]F2Laplacian relaxation parameter
[in]bx_xlobox for low x relaxation
[in]bx_xhibox for high x relaxation
[in]bx_ylobox for low y relaxation
[in]bx_yhibox for high y relaxation
[in]arr_xloarray for low x relaxation
[in]arr_xhiarray for high x relaxation
[in]arr_yloarray for low y relaxation
[in]arr_yhiarray for high y relaxation
[in]data_arrdata array
[out]rhs_arrRHS array
256 {
257  // RHS computation
258  int Spec_z = set_width;
259  int Relax_z = width - Spec_z + 1;
260  amrex::Real num = amrex::Real(Spec_z + Relax_z);
261  amrex::Real denom = amrex::Real(Relax_z - 1);
262  amrex::Real SpecExp = -std::log(0.1) / amrex::Real(width - Spec_z);
263  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
264  {
265  // Corners with x boxes
266  int j_lo = std::min(j-dom_lo.y,width-1);
267  int j_hi = std::min(dom_hi.y-j,width-1);
268  int jj = std::min(j_lo,j_hi);
269  int n_ind = std::min(i-dom_lo.x,jj) + 1;
270  AMREX_ASSERT(n_ind > Spec_z);
271  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
272  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
273  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
274  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
275  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
276  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
277  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
278  amrex::Real delta = arr_xlo(i ,j ,k,n) - d;
279  amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
280  amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
281  amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
282  amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
283  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
284  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
285  rhs_arr(i,j,k,n+icomp) += Temp;
286  },
287  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
288  {
289  // Corners with x boxes
290  int j_lo = std::min(j-dom_lo.y,width-1);
291  int j_hi = std::min(dom_hi.y-j,width-1);
292  int jj = std::min(j_lo,j_hi);
293  int n_ind = std::min(dom_hi.x-i,jj) + 1;
294  AMREX_ASSERT(n_ind > Spec_z);
295  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
296  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
297  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
298  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
299  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
300  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
301  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
302  amrex::Real delta = arr_xhi(i ,j ,k,n) - d;
303  amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
304  amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
305  amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
306  amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
307  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
308  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
309  rhs_arr(i,j,k,n+icomp) += Temp;
310  });
311 
312  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
313  {
314  // No corners for y boxes
315  int n_ind = j - dom_lo.y + 1;
316  AMREX_ASSERT(n_ind > Spec_z);
317  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
318  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
319  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
320  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
321  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
322  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
323  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
324  amrex::Real delta = arr_ylo(i ,j ,k,n) - d;
325  amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
326  amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
327  amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
328  amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
329  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
330  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
331  rhs_arr(i,j,k,n+icomp) += Temp;
332  },
333  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
334  {
335  // No corners for y boxes
336  int n_ind = dom_hi.y - j + 1;
337  AMREX_ASSERT(n_ind > Spec_z);
338  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
339  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
340  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
341  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
342  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
343  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
344  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
345  amrex::Real delta = arr_yhi(i ,j ,k,n) - d;
346  amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
347  amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
348  amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
349  amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
350  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
351  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
352  rhs_arr(i,j,k,n+icomp) += Temp;
353  });
354 }

Referenced by realbdy_compute_interior_ghost_rhs().

Here is the caller graph for this function:

◆ realbdy_set_rhs_in_spec_region()

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,
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 
)

Zero RHS in the set region

Parameters
[in]icompcomponent offset
[in]num_varnumber of variables to loop
[in]bx_xlobox for low x relaxation
[in]bx_xhibox for high x relaxation
[in]bx_ylobox for low y relaxation
[in]bx_yhibox for high y relaxation
[out]rhs_arrRHS array
165 {
166  int Spec_z = set_width;
167  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
168  {
169  // Corners with x boxes
170  int j_lo = std::min(j-dom_lo.y,width-1);
171  int j_hi = std::min(dom_hi.y-j,width-1);
172  int jj = std::min(j_lo,j_hi);
173  int n_ind = std::min(i-dom_lo.x,jj) + 1;
174  if (n_ind <= Spec_z) {
175  amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
176  rhs_arr(i,j,k,n+icomp) = tend;
177  }
178  },
179  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
180  {
181  // Corners with x boxes
182  int j_lo = std::min(j-dom_lo.y,width-1);
183  int j_hi = std::min(dom_hi.y-j,width-1);
184  int jj = std::min(j_lo,j_hi);
185  int n_ind = std::min(dom_hi.x-i,jj) + 1;
186  if (n_ind <= Spec_z) {
187  amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
188  rhs_arr(i,j,k,n+icomp) = tend;
189  }
190  });
191 
192  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
193  {
194  // No corners for y boxes
195  int n_ind = j - dom_lo.y + 1;
196  if (n_ind <= Spec_z) {
197  amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
198  rhs_arr(i,j,k,n+icomp) = tend;
199  }
200  },
201  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
202  {
203  // No corners for y boxes
204  int n_ind = dom_hi.y - j + 1;
205  if (n_ind <= Spec_z) {
206  amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
207  rhs_arr(i,j,k,n+icomp) = tend;
208  }
209  });
210 }

Referenced by realbdy_compute_interior_ghost_rhs().

Here is the caller graph for this function:

◆ Time_Avg_Vel_atCC()

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 
)

◆ VelocityToMomentum()

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 
)

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

Here is the caller graph for this function: