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 int &set_width=0, 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 &bdy_time_interval, const amrex::Real &time, const amrex::Real &delta_t, const amrex::Real &stop_time_elapsed, int width, int set_width, bool do_upwind, 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_x, const int &set_width_y, const amrex::Box &domain, 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 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 
)
482 {
483  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
484  {
485  const amrex::Box& bx = mfi.growntilebox(nghost);
486  if (bx.ok())
487  {
488  auto dstFab = dst.array(mfi);
489  const auto maskFab = imask.const_array(mfi);
490  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
491  {
492  dstFab(i,j,k) *= (1-maskFab(i,j,k));
493  });
494  }
495  }
496 }

Referenced by add_thin_body_sources(), and ERF::project_velocity_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 
)
460 {
461  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
462  {
463  const amrex::Box& bx = mfi.growntilebox(nghost);
464  if (bx.ok())
465  {
466  auto dstFab = dst.array(mfi);
467  const auto maskFab = imask.const_array(mfi);
468  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
469  {
470  dstFab(i,j,k) *= maskFab(i,j,k);
471  });
472  }
473  }
474 }

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

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 
)
516 {
517  return p0 * ( (z - z1) * (z - z2) ) / ( (z0 - z1) * (z0 - z2) )
518  + p1 * ( (z - z0) * (z - z2) ) / ( (z1 - z0) * (z1 - z2) )
519  + p2 * ( (z - z0) * (z - z1) ) / ( (z2 - z0) * (z2 - z1) );
520 }
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 bdy_time_interval,
const amrex::Real time,
const amrex::Real delta_t,
const amrex::Real stop_time_elapsed,
int  width,
int  set_width,
bool  do_upwind,
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_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]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
348 {
349  amrex::IntVect iv = bx_xlo.type();
350  amrex::Real ioff = (iv[0]==1) ? 0.0 : 0.5;
351  amrex::Real joff = (iv[1]==1) ? 0.0 : 0.5;
352  const auto& dom_cc_lo = lbound(domain_cc);
353  const auto& dom_cc_hi = ubound(domain_cc);
354  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
355  {
356  // Corners with x boxes
357  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
358  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
359  amrex::Real x_end = ProbLo[0] + (width + ioff) * dx[0];
360  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
361  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
362  amrex::Real xi = (x_end - x) / (x_end - ProbLo[0]);
363  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
364  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
365  amrex::Real eta = std::max(eta_lo,eta_hi);
366  amrex::Real Factor = std::max(xi*xi,eta*eta);
367 
368  amrex::Real delta = arr_xlo(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
369  amrex::Real Temp = Factor*F1*delta;
370  if (!do_upwind) {
371  rhs_arr(i,j,k,n+icomp) += Temp;
372  } else {
373  int jju = std::min(std::max(j,dom_cc_lo.y),dom_cc_hi.y);
374  int iiv = std::min(std::max(i,dom_cc_lo.x),dom_cc_hi.x);
375  if ( (u_xlo(dom_cc_lo.x,jju,k) >= 0.0) ||
376  ((j == dom_cc_lo.y ) && (v_xlo(iiv,dom_cc_lo.y ,k) >= 0.0)) ||
377  ((j == dom_cc_hi.y+iv[1]) && (v_xlo(iiv,dom_cc_hi.y+1,k) <= 0.0)) ) {
378  rhs_arr(i,j,k,n+icomp) += Temp;
379  }
380  }
381  },
382  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
383  {
384  // Corners with x boxes
385  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
386  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
387  amrex::Real x_strt = ProbHi[0] - (width + ioff) * dx[0];
388  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
389  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
390  amrex::Real xi = (x - x_strt) / (ProbHi[0] - x_strt);
391  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
392  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
393  amrex::Real eta = std::max(eta_lo,eta_hi);
394  amrex::Real Factor = std::max(xi*xi,eta*eta);
395 
396  amrex::Real delta = arr_xhi(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
397  amrex::Real Temp = Factor*F1*delta;
398  if (!do_upwind) {
399  rhs_arr(i,j,k,n+icomp) += Temp;
400  } else {
401  int jju = std::min(std::max(j,dom_cc_lo.y),dom_cc_hi.y);
402  int iiv = std::min(std::max(i,dom_cc_lo.x),dom_cc_hi.x);
403  if ( (u_xhi(dom_cc_hi.x+1,jju,k) <= 0.0) ||
404  ((j == dom_cc_lo.y ) && (v_xhi(iiv,dom_cc_lo.y ,k) >= 0.0)) ||
405  ((j == dom_cc_hi.y+iv[1]) && (v_xhi(iiv,dom_cc_hi.y+1,k) <= 0.0)) ) {
406  rhs_arr(i,j,k,n+icomp) += Temp;
407  }
408  }
409  });
410 
411  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
412  {
413  // No corners for y boxes
414  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
415  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
416  amrex::Real eta = (y_end - y) / (y_end - ProbLo[1]);
417  amrex::Real Factor = eta*eta;
418 
419  amrex::Real delta = arr_ylo(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
420  amrex::Real Temp = Factor*F1*delta;
421  if (!do_upwind) {
422  rhs_arr(i,j,k,n+icomp) += Temp;
423  } else {
424  int iiv = std::min(std::max(i,dom_cc_lo.x+width),dom_cc_hi.x-width);
425  if (v_ylo(iiv,dom_cc_lo.y,k) >= 0.0) {
426  rhs_arr(i,j,k,n+icomp) += Temp;
427  }
428  }
429  },
430  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
431  {
432  // No corners for y boxes
433  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
434  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
435  amrex::Real eta = (y - y_strt) / (ProbHi[1] - y_strt);
436  amrex::Real Factor = eta*eta;
437 
438  amrex::Real delta = arr_yhi(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
439  amrex::Real Temp = Factor*F1*delta;
440  if (!do_upwind) {
441  rhs_arr(i,j,k,n+icomp) += Temp;
442  } else {
443  int iiv = std::min(std::max(i,dom_cc_lo.x+width),dom_cc_hi.x-width);
444  if (v_yhi(iiv,dom_cc_hi.y+1,k) >= 0.0) {
445  rhs_arr(i,j,k,n+icomp) += Temp;
446  }
447  }
448  });
449 }
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by realbdy_compute_interior_ghost_rhs().

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 int &  set_width = 0,
const amrex::IntVect &  ng_vect = amrex::IntVect(0, 0, 0),
const bool  get_int_ng = false 
)

◆ 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_x,
const int &  set_width_y,
const amrex::Box &  domain,
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 
)

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
208 {
209  int Spec_z_x = set_width_x;
210  int Spec_z_y = set_width_y;
211  amrex::IntVect iv = bx_xlo.type();
212  const auto& dom_lo = lbound(domain);
213  const auto& dom_hi = ubound(domain);
214  const auto& dom_cc_lo = lbound(domain_cc);
215  const auto& dom_cc_hi = ubound(domain_cc);
216  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
217  {
218  // Corners with x boxes
219  int j_lo = std::min(j-dom_lo.y,width-1);
220  int j_hi = std::min(dom_hi.y-j,width-1);
221  int jj = std::min(j_lo,j_hi);
222  int n_ind_y = jj + 1;
223  int n_ind_x = i - dom_lo.x + 1;
224  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
225  amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
226  if (!do_upwind) {
227  rhs_arr(i,j,k,n+icomp) = tend;
228  } else {
229  int jju = std::min(std::max(j,dom_cc_lo.y),dom_cc_hi.y);
230  int iiv = std::min(std::max(i,dom_cc_lo.x),dom_cc_hi.x);
231  if ( (u_xlo(dom_cc_lo.x,jju,k) >= 0.0) ||
232  ((j == dom_cc_lo.y ) && (v_xlo(iiv,dom_cc_lo.y ,k) >= 0.0)) ||
233  ((j == dom_cc_hi.y+iv[1]) && (v_xlo(iiv,dom_cc_hi.y+1,k) <= 0.0)) ) {
234  rhs_arr(i,j,k,n+icomp) = tend;
235  }
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  int j_lo = std::min(j-dom_lo.y,width-1);
243  int j_hi = std::min(dom_hi.y-j,width-1);
244  int jj = std::min(j_lo,j_hi);
245  int n_ind_y = jj + 1;
246  int n_ind_x = dom_hi.x - i + 1;
247  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
248  amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
249  if (!do_upwind) {
250  rhs_arr(i,j,k,n+icomp) = tend;
251  } else {
252  int jju = std::min(std::max(j,dom_cc_lo.y),dom_cc_hi.y);
253  int iiv = std::min(std::max(i,dom_cc_lo.x),dom_cc_hi.x);
254  if ( (u_xhi(dom_cc_hi.x+1,jju,k) <= 0.0) ||
255  ((j == dom_cc_lo.y ) && (v_xhi(iiv,dom_cc_lo.y ,k) >= 0.0)) ||
256  ((j == dom_cc_hi.y+iv[1]) && (v_xhi(iiv,dom_cc_hi.y+1,k) <= 0.0)) ) {
257  rhs_arr(i,j,k,n+icomp) = tend;
258  }
259  }
260  }
261  });
262 
263  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
264  {
265  // No corners for y boxes
266  int n_ind = j - dom_lo.y + 1;
267  if (n_ind <= Spec_z_y) {
268  amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
269  if (!do_upwind) {
270  rhs_arr(i,j,k,n+icomp) = tend;
271  } else {
272  int iiv = std::min(std::max(i,dom_cc_lo.x+width),dom_cc_hi.x-width);
273  if (v_ylo(iiv,dom_cc_lo.y,k) >= 0.0) {
274  rhs_arr(i,j,k,n+icomp) = tend;
275  }
276  }
277  }
278  },
279  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
280  {
281  // No corners for y boxes
282  int n_ind = dom_hi.y - j + 1;
283  if (n_ind <= Spec_z_y) {
284  amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
285  if (!do_upwind) {
286  rhs_arr(i,j,k,n+icomp) = tend;
287  } else {
288  int iiv = std::min(std::max(i,dom_cc_lo.x+width),dom_cc_hi.x-width);
289  if (v_yhi(iiv,dom_cc_hi.y+1,k) >= 0.0) {
290  rhs_arr(i,j,k,n+icomp) = tend;
291  }
292  }
293  }
294  });
295 }
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1

Referenced by realbdy_compute_interior_ghost_rhs().

Here is the caller graph for this function:

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