ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERFPhysBCFunct_v Class Reference

#include <ERF_PhysBCFunct.H>

Collaboration diagram for ERFPhysBCFunct_v:

Public Member Functions

 ERFPhysBCFunct_v (const int lev, const amrex::Geometry &geom, const amrex::Vector< amrex::BCRec > &domain_bcs_type, const amrex::Gpu::DeviceVector< amrex::BCRec > &domain_bcs_type_d, amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_max > bc_extdir_vals, amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_max > bc_neumann_vals, std::unique_ptr< amrex::MultiFab > &z_phys_nd, const bool use_real_bcs)
 
 ~ERFPhysBCFunct_v ()
 
void operator() (amrex::MultiFab &mf, int icomp, int ncomp, amrex::IntVect const &nghost, const amrex::Real time, int bccomp)
 
void impose_lateral_yvel_bcs (const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, int bccomp)
 
void impose_vertical_yvel_bcs (const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, const amrex::Array4< amrex::Real const > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > dxInv, int bccomp)
 

Private Attributes

int m_lev
 
amrex::Geometry m_geom
 
amrex::Vector< amrex::BCRec > m_domain_bcs_type
 
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
 
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_maxm_bc_extdir_vals
 
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_maxm_bc_neumann_vals
 
amrex::MultiFab * m_z_phys_nd
 
bool m_use_real_bcs
 

Constructor & Destructor Documentation

◆ ERFPhysBCFunct_v()

ERFPhysBCFunct_v::ERFPhysBCFunct_v ( const int  lev,
const amrex::Geometry &  geom,
const amrex::Vector< amrex::BCRec > &  domain_bcs_type,
const amrex::Gpu::DeviceVector< amrex::BCRec > &  domain_bcs_type_d,
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_max bc_extdir_vals,
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_max bc_neumann_vals,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
const bool  use_real_bcs 
)
inline
140  : m_lev(lev),
141  m_geom(geom), m_domain_bcs_type(domain_bcs_type),
142  m_domain_bcs_type_d(domain_bcs_type_d),
143  m_bc_extdir_vals(bc_extdir_vals),
144  m_bc_neumann_vals(bc_neumann_vals),
145  m_z_phys_nd(z_phys_nd.get()),
146  m_use_real_bcs(use_real_bcs)
147  {}
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:175
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:177
int m_lev
Definition: ERF_PhysBCFunct.H:174
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:178
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:181
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:179
amrex::Vector< amrex::BCRec > m_domain_bcs_type
Definition: ERF_PhysBCFunct.H:176
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:180

◆ ~ERFPhysBCFunct_v()

ERFPhysBCFunct_v::~ERFPhysBCFunct_v ( )
inline
149 {}

Member Function Documentation

◆ impose_lateral_yvel_bcs()

void ERFPhysBCFunct_v::impose_lateral_yvel_bcs ( const amrex::Array4< amrex::Real > &  dest_arr,
const amrex::Box &  bx,
const amrex::Box &  domain,
int  bccomp 
)
17 {
18  BL_PROFILE_VAR("impose_lateral_yvel_bcs()",impose_lateral_yvel_bcs);
19  const auto& dom_lo = lbound(domain);
20  const auto& dom_hi = ubound(domain);
21 
22  // Based on BCRec for the domain, we need to make BCRec for this Box
23  // bccomp is used as starting index for m_domain_bcs_type
24  // 0 is used as starting index for bcrs
25  int ncomp = 1;
26  Vector<BCRec> bcrs(ncomp);
27  setBC(enclosedCells(bx), domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
28 
29  // xlo: ori = 0
30  // ylo: ori = 1
31  // zlo: ori = 2
32  // xhi: ori = 3
33  // yhi: ori = 4
34  // zhi: ori = 5
35 
36  Gpu::DeviceVector<BCRec> bcrs_d(ncomp);
37  Gpu::copyAsync(Gpu::hostToDevice, bcrs.begin(), bcrs.end(), bcrs_d.begin());
38  const BCRec* bc_ptr = bcrs_d.data();
39 
40  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>, AMREX_SPACEDIM+NVAR_max> l_bc_extdir_vals_d;
41 
42  for (int i = 0; i < ncomp; i++)
43  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++)
44  l_bc_extdir_vals_d[i][ori] = m_bc_extdir_vals[bccomp+i][ori];
45 
46  GeometryData const& geomdata = m_geom.data();
47  bool is_periodic_in_x = geomdata.isPeriodic(0);
48  bool is_periodic_in_y = geomdata.isPeriodic(1);
49 
50  // First do all ext_dir bcs
51  if (!is_periodic_in_x)
52  {
53  // Populate ghost cells on lo-x and hi-x domain boundaries
54  Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
55  Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
56  ParallelFor(
57  bx_xlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
58  int iflip = dom_lo.x - 1- i;
59  if (bc_ptr[n].lo(0) == ERFBCType::ext_dir) {
60  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][0];
61  } else if (bc_ptr[n].lo(0) == ERFBCType::foextrap) {
62  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
63  } else if (bc_ptr[n].lo(0) == ERFBCType::open) {
64  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
65  } else if (bc_ptr[n].lo(0) == ERFBCType::reflect_even) {
66  dest_arr(i,j,k) = dest_arr(iflip,j,k);
67  } else if (bc_ptr[n].lo(0) == ERFBCType::reflect_odd) {
68  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
69  }
70  },
71  bx_xhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
72  int iflip = 2*dom_hi.x + 1 - i;
73  if (bc_ptr[n].hi(0) == ERFBCType::ext_dir) {
74  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][3];
75  } else if (bc_ptr[n].hi(0) == ERFBCType::foextrap) {
76  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
77  } else if (bc_ptr[n].hi(0) == ERFBCType::open) {
78  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
79  } else if (bc_ptr[n].hi(0) == ERFBCType::reflect_even) {
80  dest_arr(i,j,k) = dest_arr(iflip,j,k);
81  } else if (bc_ptr[n].hi(0) == ERFBCType::reflect_odd) {
82  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
83  }
84  }
85  );
86  }
87 
88  if (!is_periodic_in_y)
89  {
90  // Populate ghost cells on lo-y and hi-y domain boundaries
91  Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
92  Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+2);
93  Box bx_ylo_face(bx); bx_ylo_face.setSmall(1,dom_lo.y ); bx_ylo_face.setBig(1,dom_lo.y );
94  Box bx_yhi_face(bx); bx_yhi_face.setSmall(1,dom_hi.y+1); bx_yhi_face.setBig(1,dom_hi.y+1);
95 
96  ParallelFor(
97  bx_ylo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
98  {
99  int jflip = dom_lo.y-j;
100  if (bc_ptr[n].lo(1) == ERFBCType::ext_dir) {
101  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][1];
102  } else if (bc_ptr[n].lo(1) == ERFBCType::foextrap) {
103  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
104  } else if (bc_ptr[n].lo(1) == ERFBCType::open) {
105  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
106  } else if (bc_ptr[n].lo(1) == ERFBCType::reflect_even) {
107  dest_arr(i,j,k) = dest_arr(i,jflip,k);
108  } else if (bc_ptr[n].lo(1) == ERFBCType::reflect_odd) {
109  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
110  } else if (bc_ptr[n].lo(1) == ERFBCType::neumann_int) {
111  dest_arr(i,j,k) = (4.0*dest_arr(i,dom_lo.y+1,k) - dest_arr(i,dom_lo.y+2,k))/3.0;
112  }
113  },
114  // We only set the values on the domain faces themselves if EXT_DIR
115  bx_ylo_face, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
116  {
117  if (bc_ptr[n].lo(1) == ERFBCType::ext_dir) {
118  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][1];
119  } else if (bc_ptr[n].lo(1) == ERFBCType::neumann_int) {
120  dest_arr(i,j,k) = (4.0*dest_arr(i,dom_lo.y+1,k) - dest_arr(i,dom_lo.y+2,k))/3.0;
121  }
122  }
123  );
124  ParallelFor(
125  bx_yhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
126  {
127  int jflip = 2*(dom_hi.y + 1) - j;
128  if (bc_ptr[n].hi(1) == ERFBCType::ext_dir) {
129  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][4];
130  } else if (bc_ptr[n].hi(1) == ERFBCType::foextrap) {
131  dest_arr(i,j,k) = dest_arr(i,dom_hi.y+1,k);
132  } else if (bc_ptr[n].hi(1) == ERFBCType::open) {
133  dest_arr(i,j,k) = dest_arr(i,dom_hi.y+1,k);
134  } else if (bc_ptr[n].hi(1) == ERFBCType::reflect_even) {
135  dest_arr(i,j,k) = dest_arr(i,jflip,k);
136  } else if (bc_ptr[n].hi(1) == ERFBCType::reflect_odd) {
137  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
138  } else if (bc_ptr[n].hi(1) == ERFBCType::neumann_int) {
139  dest_arr(i,j,k) = (4.0*dest_arr(i,dom_hi.y,k) - dest_arr(i,dom_hi.y-1,k))/3.0;
140  }
141  },
142  // We only set the values on the domain faces themselves if EXT_DIR
143  bx_yhi_face, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
144  {
145  if (bc_ptr[n].hi(1) == ERFBCType::ext_dir) {
146  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][4];
147  } else if (bc_ptr[n].hi(1) == ERFBCType::neumann_int) {
148  dest_arr(i,j,k) = (4.0*dest_arr(i,dom_hi.y,k) - dest_arr(i,dom_hi.y-1,k))/3.0;
149  }
150  }
151  );
152  }
153  Gpu::streamSynchronize();
154 }
#define NVAR_max
Definition: IndexDefines.H:27
void impose_lateral_yvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, int bccomp)
Definition: BoundaryConditions_yvel.cpp:14
@ open
Definition: IndexDefines.H:157
@ reflect_odd
Definition: IndexDefines.H:148
@ foextrap
Definition: IndexDefines.H:151
@ ext_dir
Definition: IndexDefines.H:152
@ neumann_int
Definition: IndexDefines.H:156
@ reflect_even
Definition: IndexDefines.H:150

◆ impose_vertical_yvel_bcs()

void ERFPhysBCFunct_v::impose_vertical_yvel_bcs ( const amrex::Array4< amrex::Real > &  dest_arr,
const amrex::Box &  bx,
const amrex::Box &  domain,
const amrex::Array4< amrex::Real const > &  z_nd,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM >  dxInv,
int  bccomp 
)
172 {
173  BL_PROFILE_VAR("impose_vertical_yvel_bcs()",impose_vertical_yvel_bcs);
174  const auto& dom_lo = lbound(domain);
175  const auto& dom_hi = ubound(domain);
176 
177  // Based on BCRec for the domain, we need to make BCRec for this Box
178  // bccomp is used as starting index for m_domain_bcs_type
179  // 0 is used as starting index for bcrs
180  int ncomp = 1;
181  Vector<BCRec> bcrs(ncomp);
182  setBC(enclosedCells(bx), domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
183 
184  // xlo: ori = 0
185  // ylo: ori = 1
186  // zlo: ori = 2
187  // xhi: ori = 3
188  // yhi: ori = 4
189  // zhi: ori = 5
190 
191  Gpu::DeviceVector<BCRec> bcrs_d(ncomp);
192  Gpu::copyAsync(Gpu::hostToDevice, bcrs.begin(), bcrs.end(), bcrs_d.begin());
193  const BCRec* bc_ptr = bcrs_d.data();
194 
195  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>, AMREX_SPACEDIM+NVAR_max> l_bc_extdir_vals_d;
196 
197  for (int i = 0; i < ncomp; i++)
198  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++)
199  l_bc_extdir_vals_d[i][ori] = m_bc_extdir_vals[bccomp+i][ori];
200 
201  GeometryData const& geomdata = m_geom.data();
202 
203  {
204  // Populate ghost cells on lo-z and hi-z domain boundaries
205  Box bx_zlo(bx); bx_zlo.setBig (2,dom_lo.z-1);
206  Box bx_zhi(bx); bx_zhi.setSmall(2,dom_hi.z+1);
207  ParallelFor(
208  bx_zlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
209  int kflip = dom_lo.z - 1 - k;
210  if (bc_ptr[n].lo(2) == ERFBCType::ext_dir) {
211  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][2];
212  } else if (bc_ptr[n].lo(2) == ERFBCType::foextrap) {
213  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
214  } else if (bc_ptr[n].lo(2) == ERFBCType::open) {
215  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
216  } else if (bc_ptr[n].lo(2) == ERFBCType::reflect_even) {
217  dest_arr(i,j,k) = dest_arr(i,j,kflip);
218  } else if (bc_ptr[n].lo(2) == ERFBCType::reflect_odd) {
219  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
220  }
221  },
222  bx_zhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
223  int kflip = 2*dom_hi.z + 1 - k;
224  if (bc_ptr[n].hi(2) == ERFBCType::ext_dir) {
225  dest_arr(i,j,k) = l_bc_extdir_vals_d[n][5];
226  } else if (bc_ptr[n].hi(2) == ERFBCType::foextrap) {
227  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
228  } else if (bc_ptr[n].hi(2) == ERFBCType::open) {
229  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
230  } else if (bc_ptr[n].hi(2) == ERFBCType::reflect_even) {
231  dest_arr(i,j,k) = dest_arr(i,j,kflip);
232  } else if (bc_ptr[n].hi(2) == ERFBCType::reflect_odd) {
233  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
234  }
235  }
236  );
237  }
238 
239  if (m_z_phys_nd) {
240  const auto& bx_lo = lbound(bx);
241  const auto& bx_hi = ubound(bx);
242  // Neumann conditions (d<var>/dn = 0) must be aware of the surface normal with terrain.
243  // An additional source term arises from d<var>/dx & d<var>/dy & met_h_xi/eta/zeta.
244  //=====================================================================================
245  // Only modify scalars, U, or V
246  // Loop over each component
247  for (int n = 0; n < ncomp; n++) {
248  // Hit for Neumann condition at kmin
249  if(bcrs[n].lo(2) == ERFBCType::foextrap) {
250  // Loop over ghost cells in bottom XY-plane (valid box)
251  Box xybx = bx;
252  xybx.setBig(2,-1);
253  xybx.setSmall(2,bx.smallEnd()[2]);
254  int k0 = 0;
255 
256  // Get the dz cell size
257  Real dz = geomdata.CellSize(2);
258 
259  // Fill all the Neumann srcs with terrain
260  ParallelFor(xybx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
261  {
262  // Clip indices for ghost-cells
263  int ii = amrex::min(amrex::max(i,dom_lo.x),dom_hi.x);
264  int jj = amrex::min(amrex::max(j,dom_lo.y),dom_hi.y);
265 
266  // Get metrics
267  Real met_h_xi = Compute_h_xi_AtJface (ii, jj, k0, dxInv, z_phys_nd);
268  Real met_h_eta = Compute_h_eta_AtJface (ii, jj, k0, dxInv, z_phys_nd);
269  Real met_h_zeta = Compute_h_zeta_AtJface(ii, jj, k0, dxInv, z_phys_nd);
270 
271  // GradX at IJK location inside domain -- this relies on the assumption that we have
272  // used foextrap for cell-centered quantities outside the domain to define the gradient as zero
273  Real GradVarx, GradVary;
274  if (i < dom_lo.x-1 || i > dom_hi.x+1)
275  GradVarx = 0.0;
276  else if (i+1 > bx_hi.x)
277  GradVarx = dxInv[0] * (dest_arr(i ,j,k0) - dest_arr(i-1,j,k0));
278  else if (i-1 < bx_lo.x)
279  GradVarx = dxInv[0] * (dest_arr(i+1,j,k0) - dest_arr(i ,j,k0));
280  else
281  GradVarx = 0.5 * dxInv[0] * (dest_arr(i+1,j,k0) - dest_arr(i-1,j,k0));
282 
283  // GradY at IJK location inside domain -- this relies on the assumption that we have
284  // used foextrap for cell-centered quantities outside the domain to define the gradient as zero
285  if (j < dom_lo.y-1 || j > dom_hi.y+1)
286  GradVary = 0.0;
287  else if (j+1 > bx_hi.y)
288  GradVary = dxInv[1] * (dest_arr(i,j ,k0) - dest_arr(i,j-1,k0));
289  else if (j-1 < bx_lo.y)
290  GradVary = dxInv[1] * (dest_arr(i,j+1,k0) - dest_arr(i,j ,k0));
291  else
292  GradVary = 0.5 * dxInv[1] * (dest_arr(i,j+1,k0) - dest_arr(i,j-1,k0));
293 
294  // Prefactor
295  Real met_fac = met_h_zeta / ( met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1. );
296 
297  // Accumulate in bottom ghost cell (EXTRAP already populated)
298  dest_arr(i,j,k) -= dz * met_fac * ( met_h_xi * GradVarx + met_h_eta * GradVary );
299  });
300  } // foextrap
301  } // ncomp
302  } //m_z_phys_nd
303  Gpu::streamSynchronize();
304 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: TerrainMetrics.H:136
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: TerrainMetrics.H:122
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: TerrainMetrics.H:150
void impose_vertical_yvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, const amrex::Array4< amrex::Real const > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > dxInv, int bccomp)
Definition: BoundaryConditions_yvel.cpp:167
Here is the call graph for this function:

◆ operator()()

void ERFPhysBCFunct_v::operator() ( amrex::MultiFab &  mf,
int  icomp,
int  ncomp,
amrex::IntVect const &  nghost,
const amrex::Real  time,
int  bccomp 
)
138 {
139  BL_PROFILE("ERFPhysBCFunct_v::()");
140 
141  if (m_geom.isAllPeriodic()) return;
142 
143  const auto& domain = m_geom.Domain();
144  const auto dxInv = m_geom.InvCellSizeArray();
145 
146  Box gdomainy = surroundingNodes(domain,1);
147  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
148  if (m_geom.isPeriodic(i)) {
149  gdomainy.grow(i, nghost[i]);
150  }
151  }
152 
153 #ifdef AMREX_USE_OMP
154 #pragma omp parallel if (Gpu::notInLaunchRegion())
155 #endif
156  {
157  for (MFIter mfi(mf,false); mfi.isValid(); ++mfi)
158  {
159  //
160  // This is the box we pass to the different routines
161  // NOTE -- this is the full grid NOT the tile box
162  //
163  Box bx = mfi.validbox();
164 
165  //
166  // These are the boxes we use to test on relative to the domain
167  //
168  Box ybx1 = surroundingNodes(bx,1); ybx1.grow(IntVect(nghost[0],nghost[1],0));
169  Box ybx2 = surroundingNodes(bx,1); ybx2.grow(nghost);
170 
171  Array4<const Real> z_nd_arr;
172 
173  if (m_z_phys_nd)
174  {
175  z_nd_arr = m_z_phys_nd->const_array(mfi);
176  }
177 
178  if (!gdomainy.contains(ybx2))
179  {
180  const Array4<Real> vely_arr = mf.array(mfi);;
181 
182  if (!m_use_real_bcs)
183  {
184  impose_lateral_yvel_bcs(vely_arr,ybx1,domain,bccomp);
185  }
186 
187  impose_vertical_yvel_bcs(vely_arr,ybx2,domain,z_nd_arr,dxInv,bccomp);
188  }
189 
190  } // MFIter
191  } // OpenMP
192 } // operator()

Member Data Documentation

◆ m_bc_extdir_vals

amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NVAR_max> ERFPhysBCFunct_v::m_bc_extdir_vals
private

◆ m_bc_neumann_vals

amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NVAR_max> ERFPhysBCFunct_v::m_bc_neumann_vals
private

◆ m_domain_bcs_type

amrex::Vector<amrex::BCRec> ERFPhysBCFunct_v::m_domain_bcs_type
private

◆ m_domain_bcs_type_d

amrex::Gpu::DeviceVector<amrex::BCRec> ERFPhysBCFunct_v::m_domain_bcs_type_d
private

◆ m_geom

amrex::Geometry ERFPhysBCFunct_v::m_geom
private

◆ m_lev

int ERFPhysBCFunct_v::m_lev
private

◆ m_use_real_bcs

bool ERFPhysBCFunct_v::m_use_real_bcs
private

◆ m_z_phys_nd

amrex::MultiFab* ERFPhysBCFunct_v::m_z_phys_nd
private

The documentation for this class was generated from the following files: