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_SurfaceLayer.H"
#include "ERF_FillPatcher.H"
#include "ERF_ReadBndryPlanes.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 cons_to_prim (const amrex::MultiFab &cons_state, amrex::MultiFab &S_prim, int ng)
 
void make_qt (const amrex::MultiFab &cons_state, amrex::MultiFab &qt, int n_qstate_into_total)
 
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, 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 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 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 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)
 
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 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 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_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 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_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)
 
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)
 
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 WeatherDataInterpolation (const amrex::Real time)
 
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)
 

Function Documentation

◆ ApplyInvertedMask()

AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask ( amrex::MultiFab &  dst,
const amrex::iMultiFab &  imask,
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 }
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);})

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

Here is the call graph for this function:
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 
)
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 }

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ChopGrids2D()

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

◆ compute_influx_outflux_bdy()

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 
)

◆ cons_to_prim()

void cons_to_prim ( const amrex::MultiFab &  cons_state,
amrex::MultiFab &  S_prim,
int  ng 
)

◆ ConvertForProjection()

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 
)

◆ enforceInOutSolvability()

void enforceInOutSolvability ( int  lev,
amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &  vels_vec,
amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &  area_vec,
const amrex::Geometry &  geom 
)

◆ enforceInOutSolvability_bdy()

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 
)

◆ 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_qt()

void make_qt ( const amrex::MultiFab &  cons_state,
amrex::MultiFab &  qt,
int  n_qstate_into_total 
)

◆ 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,
const amrex::MultiFab *  c_vfrac = nullptr 
)

◆ quad_interp_1d()

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 
)
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 }
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40

Referenced by compute_gradp_interpz().

Here is the caller graph for this function:

◆ realbdy_bc_bxs_xy()

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

◆ realbdy_compute_interior_ghost_rhs()

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 
)

◆ realbdy_compute_relaxation()

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 
)

Compute the nudging relaxation

Parameters
[in]delta_ttime step
[in]icompcomponent offset
[in]num_varnumber of variables to loop
[in]widthwidth of wrf bdy file
[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
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 }
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by realbdy_compute_interior_ghost_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ realbdy_interior_bxs_xy()

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 
)

◆ thinbody_wall_dist()

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 
)

◆ 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,
const amrex::MultiFab *  c_vfrac = nullptr 
)

◆ WeatherDataInterpolation()

void WeatherDataInterpolation ( const amrex::Real  time)

Referenced by ERF::Evolve(), and ERF::init_stuff().

Here is the caller graph for this function: