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

Function Documentation

◆ ApplyInvertedMask()

AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask ( amrex::MultiFab &  dst,
const amrex::iMultiFab &  imask,
const int  nghost = 0 
)
424 {
425  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
426  {
427  const amrex::Box& bx = mfi.growntilebox(nghost);
428  if (bx.ok())
429  {
430  auto dstFab = dst.array(mfi);
431  const auto maskFab = imask.const_array(mfi);
432  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
433  {
434  dstFab(i,j,k) *= (1-maskFab(i,j,k));
435  });
436  }
437  }
438 }

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 
)
402 {
403  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
404  {
405  const amrex::Box& bx = mfi.growntilebox(nghost);
406  if (bx.ok())
407  {
408  auto dstFab = dst.array(mfi);
409  const auto maskFab = imask.const_array(mfi);
410  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
411  {
412  dstFab(i,j,k) *= maskFab(i,j,k);
413  });
414  }
415  }
416 }

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 
)

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

◆ 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,
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 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::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
285 {
286  // RHS computation
287  amrex::IntVect iv = bx_xlo.type();
288  amrex::Real ioff = (iv[0]==1) ? 0.0 : 0.5;
289  amrex::Real joff = (iv[1]==1) ? 0.0 : 0.5;
290  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
291  {
292  // Corners with x boxes
293  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
294  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
295  amrex::Real x_end = ProbLo[0] + (width + ioff) * dx[0];
296  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
297  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
298  amrex::Real xi = (x_end - x) / (x_end - ProbLo[0]);
299  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
300  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
301  amrex::Real eta = std::min(eta_lo,eta_hi);
302  amrex::Real Factor = std::max(xi*xi,eta*eta);
303 
304  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
305  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
306  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
307  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
308  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
309  amrex::Real delta = arr_xlo(i ,j ,k,n) - d;
310  amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
311  amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
312  amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
313  amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
314  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
315  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
316  rhs_arr(i,j,k,n+icomp) += Temp;
317  },
318  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
319  {
320  // Corners with x boxes
321  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
322  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
323  amrex::Real x_strt = ProbHi[0] - (width + ioff) * dx[0];
324  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
325  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
326  amrex::Real xi = (x - x_strt) / (ProbHi[0] - x_strt);
327  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
328  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
329  amrex::Real eta = std::min(eta_lo,eta_hi);
330  amrex::Real Factor = std::max(xi*xi,eta*eta);
331 
332  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
333  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
334  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
335  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
336  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
337  amrex::Real delta = arr_xhi(i ,j ,k,n) - d;
338  amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
339  amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
340  amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
341  amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
342  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
343  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
344  rhs_arr(i,j,k,n+icomp) += Temp;
345  });
346 
347  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
348  {
349  // No corners for y boxes
350  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
351  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
352  amrex::Real eta = (y_end - y) / (y_end - ProbLo[1]);
353  amrex::Real Factor = eta*eta;
354 
355  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
356  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
357  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
358  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
359  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
360  amrex::Real delta = arr_ylo(i ,j ,k,n) - d;
361  amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
362  amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
363  amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
364  amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
365  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
366  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
367  rhs_arr(i,j,k,n+icomp) += Temp;
368  },
369  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
370  {
371  // No corners for y boxes
372  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
373  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
374  amrex::Real eta = (y - y_strt) / (ProbHi[1] - y_strt);
375  amrex::Real Factor = eta*eta;
376 
377  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
378  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
379  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
380  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
381  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
382  amrex::Real delta = arr_yhi(i ,j ,k,n) - d;
383  amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
384  amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
385  amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
386  amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
387  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
388  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
389  rhs_arr(i,j,k,n+icomp) += Temp;
390  });
391 }
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::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
191 {
192  int Spec_z_x = set_width_x;
193  int Spec_z_y = set_width_y;
194  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
195  {
196  // Corners with x boxes
197  int j_lo = std::min(j-dom_lo.y,width-1);
198  int j_hi = std::min(dom_hi.y-j,width-1);
199  int jj = std::min(j_lo,j_hi);
200  int n_ind_y = jj + 1;
201  int n_ind_x = i - dom_lo.x + 1;
202  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
203  amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
204  rhs_arr(i,j,k,n+icomp) = tend;
205  }
206  },
207  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
208  {
209  // Corners with x boxes
210  int j_lo = std::min(j-dom_lo.y,width-1);
211  int j_hi = std::min(dom_hi.y-j,width-1);
212  int jj = std::min(j_lo,j_hi);
213  int n_ind_y = jj + 1;
214  int n_ind_x = dom_hi.x - i + 1;
215  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
216  amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
217  rhs_arr(i,j,k,n+icomp) = tend;
218  }
219  });
220 
221  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
222  {
223  // No corners for y boxes
224  int n_ind = j - dom_lo.y + 1;
225  if (n_ind <= Spec_z_y) {
226  amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
227  rhs_arr(i,j,k,n+icomp) = tend;
228  }
229  },
230  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
231  {
232  // No corners for y boxes
233  int n_ind = dom_hi.y - j + 1;
234  if (n_ind <= Spec_z_y) {
235  amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
236  rhs_arr(i,j,k,n+icomp) = tend;
237  }
238  });
239 }
const auto & dom_hi
Definition: ERF_DiffSetup.H:10
const auto & dom_lo
Definition: ERF_DiffSetup.H:9

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::Advance().

Here is the caller graph for this function: