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 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, 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 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_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 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 &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)
 
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 
)
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 }
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);} } })

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

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 
)

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

◆ 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,
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 
)
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 }
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_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,
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 
)

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

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
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 }
const auto dx
Definition: ERF_InitCustomPertVels_ParticleTests.H:15
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: