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)
 
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 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, 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 
)
420 {
421  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
422  {
423  const amrex::Box& bx = mfi.growntilebox(nghost);
424  if (bx.ok())
425  {
426  auto dstFab = dst.array(mfi);
427  const auto maskFab = imask.const_array(mfi);
428  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
429  {
430  dstFab(i,j,k) *= (1-maskFab(i,j,k));
431  });
432  }
433  }
434 }

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

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 
)

◆ 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,
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
281 {
282  // RHS computation
283  amrex::IntVect iv = bx_xlo.type();
284  amrex::Real ioff = (iv[0]==1) ? 0.0 : 0.5;
285  amrex::Real joff = (iv[1]==1) ? 0.0 : 0.5;
286  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
287  {
288  // Corners with x boxes
289  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
290  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
291  amrex::Real x_end = ProbLo[0] + (width + ioff) * dx[0];
292  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
293  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
294  amrex::Real xi = (x_end - x) / (x_end - ProbLo[0]);
295  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
296  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
297  amrex::Real eta = std::min(eta_lo,eta_hi);
298  amrex::Real Factor = std::max(xi*xi,eta*eta);
299 
300  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
301  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
302  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
303  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
304  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
305  amrex::Real delta = arr_xlo(i ,j ,k,n) - d;
306  amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
307  amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
308  amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
309  amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
310  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
311  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
312  rhs_arr(i,j,k,n+icomp) += Temp;
313  },
314  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
315  {
316  // Corners with x boxes
317  amrex::Real x = ProbLo[0] + (i + ioff) * dx[0];
318  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
319  amrex::Real x_strt = ProbHi[0] - (width + ioff) * dx[0];
320  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
321  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
322  amrex::Real xi = (x - x_strt) / (ProbHi[0] - x_strt);
323  amrex::Real eta_lo = (y < y_end ) ? (y_end - y) / (y_end - ProbLo[1]) : 0.0;
324  amrex::Real eta_hi = (y > y_strt) ? (y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
325  amrex::Real eta = std::min(eta_lo,eta_hi);
326  amrex::Real Factor = std::max(xi*xi,eta*eta);
327 
328  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
329  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
330  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
331  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
332  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
333  amrex::Real delta = arr_xhi(i ,j ,k,n) - d;
334  amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
335  amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
336  amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
337  amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
338  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
339  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
340  rhs_arr(i,j,k,n+icomp) += Temp;
341  });
342 
343  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
344  {
345  // No corners for y boxes
346  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
347  amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
348  amrex::Real eta = (y_end - y) / (y_end - ProbLo[1]);
349  amrex::Real Factor = eta*eta;
350 
351  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
352  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
353  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
354  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
355  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
356  amrex::Real delta = arr_ylo(i ,j ,k,n) - d;
357  amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
358  amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
359  amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
360  amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
361  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
362  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
363  rhs_arr(i,j,k,n+icomp) += Temp;
364  },
365  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
366  {
367  // No corners for y boxes
368  amrex::Real y = ProbLo[1] + (j + joff) * dx[1];
369  amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
370  amrex::Real eta = (y - y_strt) / (ProbHi[1] - y_strt);
371  amrex::Real Factor = eta*eta;
372 
373  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
374  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
375  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
376  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
377  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
378  amrex::Real delta = arr_yhi(i ,j ,k,n) - d;
379  amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
380  amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
381  amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
382  amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
383  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
384  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
385  rhs_arr(i,j,k,n+icomp) += Temp;
386  });
387 }
amrex::Real Real
Definition: ERF_ShocInterface.H:16

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
187 {
188  int Spec_z_x = set_width_x;
189  int Spec_z_y = set_width_y;
190  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
191  {
192  // Corners with x boxes
193  int j_lo = std::min(j-dom_lo.y,width-1);
194  int j_hi = std::min(dom_hi.y-j,width-1);
195  int jj = std::min(j_lo,j_hi);
196  int n_ind_y = jj + 1;
197  int n_ind_x = i - dom_lo.x + 1;
198  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
199  amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
200  rhs_arr(i,j,k,n+icomp) = tend;
201  }
202  },
203  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
204  {
205  // Corners with x boxes
206  int j_lo = std::min(j-dom_lo.y,width-1);
207  int j_hi = std::min(dom_hi.y-j,width-1);
208  int jj = std::min(j_lo,j_hi);
209  int n_ind_y = jj + 1;
210  int n_ind_x = dom_hi.x - i + 1;
211  if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
212  amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
213  rhs_arr(i,j,k,n+icomp) = tend;
214  }
215  });
216 
217  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
218  {
219  // No corners for y boxes
220  int n_ind = j - dom_lo.y + 1;
221  if (n_ind <= Spec_z_y) {
222  amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
223  rhs_arr(i,j,k,n+icomp) = tend;
224  }
225  },
226  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
227  {
228  // No corners for y boxes
229  int n_ind = dom_hi.y - j + 1;
230  if (n_ind <= Spec_z_y) {
231  amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
232  rhs_arr(i,j,k,n+icomp) = tend;
233  }
234  });
235 }
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 
)

◆ WeatherDataInterpolation()

void WeatherDataInterpolation ( const amrex::Real  time)