ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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_ABLMost.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 &start_bdy_time, const amrex::Real &time, const amrex::Real &delta_t, 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 int &set_width, const amrex::Dim3 &dom_lo, const amrex::Dim3 &dom_hi, 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)
 

Function Documentation

◆ ApplyInvertedMask()

AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask ( amrex::MultiFab &  dst,
const amrex::iMultiFab &  imask,
const int  nghost = 0 
)
412 {
413  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
414  {
415  const amrex::Box& bx = mfi.growntilebox(nghost);
416  if (bx.ok())
417  {
418  auto dstFab = dst.array(mfi);
419  const auto maskFab = imask.const_array(mfi);
420  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
421  {
422  dstFab(i,j,k) *= (1-maskFab(i,j,k));
423  });
424  }
425  }
426 }

Referenced by add_thin_body_sources(), and ERF::project_velocities_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 
)
390 {
391  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
392  {
393  const amrex::Box& bx = mfi.growntilebox(nghost);
394  if (bx.ok())
395  {
396  auto dstFab = dst.array(mfi);
397  const auto maskFab = imask.const_array(mfi);
398  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
399  {
400  dstFab(i,j,k) *= maskFab(i,j,k);
401  });
402  }
403  }
404 }

Referenced by ERF::FillIntermediatePatch(), and ERF::project_velocities_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 &  start_bdy_time,
const amrex::Real &  time,
const amrex::Real &  delta_t,
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 int &  set_width,
const amrex::Dim3 &  dom_lo,
const amrex::Dim3 &  dom_hi,
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  int Spec_z = set_width;
284  int Relax_z = width - Spec_z + 1;
285  amrex::Real num = amrex::Real(Spec_z + Relax_z);
286  amrex::Real denom = amrex::Real(Relax_z - 1);
287  amrex::Real SpecExp = -std::log(0.1) / amrex::Real(width - Spec_z);
288  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
289  {
290  // Corners with x boxes
291  int j_lo = std::min(j-dom_lo.y,width-1);
292  int j_hi = std::min(dom_hi.y-j,width-1);
293  int jj = std::min(j_lo,j_hi);
294  int n_ind = std::min(i-dom_lo.x,jj) + 1;
295  //AMREX_ASSERT(n_ind > Spec_z);
296  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
297  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
298  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
299  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
300  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
301  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
302  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
303  amrex::Real delta = arr_xlo(i ,j ,k,n) - d;
304  amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
305  amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
306  amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
307  amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
308  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
309  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
310  rhs_arr(i,j,k,n+icomp) += Temp;
311  },
312  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
313  {
314  // Corners with x boxes
315  int j_lo = std::min(j-dom_lo.y,width-1);
316  int j_hi = std::min(dom_hi.y-j,width-1);
317  int jj = std::min(j_lo,j_hi);
318  int n_ind = std::min(dom_hi.x-i,jj) + 1;
319  //AMREX_ASSERT(n_ind > Spec_z);
320  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
321  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
322  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
323  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
324  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
325  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
326  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
327  amrex::Real delta = arr_xhi(i ,j ,k,n) - d;
328  amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
329  amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
330  amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
331  amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
332  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
333  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
334  rhs_arr(i,j,k,n+icomp) += Temp;
335  });
336 
337  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
338  {
339  // No corners for y boxes
340  int n_ind = j - dom_lo.y + 1;
341  //AMREX_ASSERT(n_ind > Spec_z);
342  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
343  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
344  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
345  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
346  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
347  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
348  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
349  amrex::Real delta = arr_ylo(i ,j ,k,n) - d;
350  amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
351  amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
352  amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
353  amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
354  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
355  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
356  rhs_arr(i,j,k,n+icomp) += Temp;
357  },
358  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
359  {
360  // No corners for y boxes
361  int n_ind = dom_hi.y - j + 1;
362  //AMREX_ASSERT(n_ind > Spec_z);
363  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
364  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
365  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
366  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
367  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
368  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
369  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
370  amrex::Real delta = arr_yhi(i ,j ,k,n) - d;
371  amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
372  amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
373  amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
374  amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
375  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
376  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
377  rhs_arr(i,j,k,n+icomp) += Temp;
378  });
379 }

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 }

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 
)

Referenced by ERF::Advance(), ERF::advance_dycore(), ERF::AverageDownTo(), ERF::FillCoarsePatch(), and ERF::FillIntermediatePatch().

Here is the caller graph for this function: