ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERFPhysBCFunct_w Class Reference

#include <ERF_PhysBCFunct.H>

Collaboration diagram for ERFPhysBCFunct_w:

Public Member Functions

 ERFPhysBCFunct_w (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+NBCVAR_max > bc_extdir_vals, amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > bc_neumann_vals, const TerrainType &terrain_type, std::unique_ptr< amrex::MultiFab > &z_phys_nd, const bool use_real_bcs, amrex::Real *w_bc_data)
 
 ~ERFPhysBCFunct_w ()
 
void operator() (amrex::MultiFab &mf, amrex::MultiFab &xvel, amrex::MultiFab &yvel, amrex::IntVect const &nghost, const amrex::Real time, int bccomp, bool do_fb)
 
void impose_lateral_zvel_bcs (const amrex::Array4< amrex::Real > &dest_arr, const amrex::Array4< amrex::Real const > &xvel_arr, const amrex::Array4< amrex::Real const > &yvel_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, TerrainType terrain_type, int bccomp)
 
void impose_vertical_zvel_bcs (const amrex::Array4< amrex::Real > &dest_arr, const amrex::Array4< amrex::Real const > &xvel_arr, const amrex::Array4< amrex::Real const > &yvel_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_u, int bccomp_v, int bccomp_w, TerrainType terrain_type)
 

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+NBCVAR_maxm_bc_extdir_vals
 
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_maxm_bc_neumann_vals
 
TerrainType m_terrain_type
 
amrex::MultiFab * m_z_phys_nd
 
bool m_use_real_bcs
 
amrex::Real * m_w_bc_data
 

Constructor & Destructor Documentation

◆ ERFPhysBCFunct_w()

ERFPhysBCFunct_w::ERFPhysBCFunct_w ( 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+NBCVAR_max bc_extdir_vals,
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max bc_neumann_vals,
const TerrainType &  terrain_type,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
const bool  use_real_bcs,
amrex::Real *  w_bc_data 
)
inline
211  : m_lev(lev),
212  m_geom(geom), m_domain_bcs_type(domain_bcs_type),
213  m_domain_bcs_type_d(domain_bcs_type_d),
214  m_bc_extdir_vals(bc_extdir_vals),
215  m_bc_neumann_vals(bc_neumann_vals),
216  m_terrain_type(terrain_type),
217  m_z_phys_nd(z_phys_nd.get()),
218  m_use_real_bcs(use_real_bcs),
219  m_w_bc_data(w_bc_data)
220  { }
amrex::Vector< amrex::BCRec > m_domain_bcs_type
Definition: ERF_PhysBCFunct.H:257
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:256
int m_lev
Definition: ERF_PhysBCFunct.H:255
amrex::Real * m_w_bc_data
Definition: ERF_PhysBCFunct.H:264
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:263
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:258
TerrainType m_terrain_type
Definition: ERF_PhysBCFunct.H:261
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:262
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:260
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:259

◆ ~ERFPhysBCFunct_w()

ERFPhysBCFunct_w::~ERFPhysBCFunct_w ( )
inline
222 {}

Member Function Documentation

◆ impose_lateral_zvel_bcs()

void ERFPhysBCFunct_w::impose_lateral_zvel_bcs ( const amrex::Array4< amrex::Real > &  dest_arr,
const amrex::Array4< amrex::Real const > &  xvel_arr,
const amrex::Array4< amrex::Real const > &  yvel_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,
TerrainType  terrain_type,
int  bccomp 
)
23 {
24  BL_PROFILE_VAR("impose_lateral_zvel_bcs()",impose_lateral_zvel_bcs);
25  const auto& dom_lo = lbound(domain);
26  const auto& dom_hi = ubound(domain);
27 
28  // Based on BCRec for the domain, we need to make BCRec for this Box
29  // bccomp is used as starting index for m_domain_bcs_type
30  // 0 is used as starting index for bcrs
31  Vector<BCRec> bcrs_w(1);
32  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs_w);
33 
34  bool l_use_terrain_fitted_coords = ( (terrain_type == TerrainType::StaticFittedMesh) ||
35  (terrain_type == TerrainType::MovingFittedMesh) );
36 
37  // xlo: ori = 0
38  // ylo: ori = 1
39  // zlo: ori = 2
40  // xhi: ori = 3
41  // yhi: ori = 4
42  // zhi: ori = 5
43 
44  Gpu::DeviceVector<BCRec> bcrs_w_d(1);
45  Gpu::copyAsync(Gpu::hostToDevice, bcrs_w.begin(), bcrs_w.end(), bcrs_w_d.begin());
46  const BCRec* bc_ptr_w = bcrs_w_d.data();
47 
48  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,1> l_bc_extdir_vals_d;
49 
50  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
51  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp][ori];
52  }
53 
54  GeometryData const& geomdata = m_geom.data();
55  bool is_periodic_in_x = geomdata.isPeriodic(0);
56  bool is_periodic_in_y = geomdata.isPeriodic(1);
57 
58  FArrayBox dhdtfab;
59 
60  // First do all ext_dir bcs
61  if (!is_periodic_in_x)
62  {
63  Real* zvel_bc_ptr = m_w_bc_data;
64  Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
65  Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
66  ParallelFor(bx_xlo, bx_xhi,
67  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
68  int iflip = dom_lo.x - 1 - i;
69  if ( (bc_ptr_w[0].lo(0) == ERFBCType::ext_dir) ||
70  (bc_ptr_w[0].lo(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_lo.x,j,k) >= 0.) )
71  {
72  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
73  if (l_use_terrain_fitted_coords) {
74  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),xvel_arr,yvel_arr,z_phys_nd,dxInv);
75  }
76  } else if (bc_ptr_w[0].lo(0) == ERFBCType::foextrap) {
77  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
78  } else if (bc_ptr_w[0].lo(0) == ERFBCType::open) {
79  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
80  } else if (bc_ptr_w[0].lo(0) == ERFBCType::reflect_even) {
81  dest_arr(i,j,k) = dest_arr(iflip,j,k);
82  } else if (bc_ptr_w[0].lo(0) == ERFBCType::reflect_odd) {
83  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
84  }
85  },
86  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
87  int iflip = 2*dom_hi.x + 1 - i;
88  if ( (bc_ptr_w[0].hi(0) == ERFBCType::ext_dir) ||
89  (bc_ptr_w[0].hi(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_hi.x+1,j,k) <= 0.) )
90  {
91  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
92  if (l_use_terrain_fitted_coords) {
93  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),xvel_arr,yvel_arr,z_phys_nd,dxInv);
94  }
95  } else if (bc_ptr_w[0].hi(0) == ERFBCType::foextrap) {
96  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
97  } else if (bc_ptr_w[0].hi(0) == ERFBCType::open) {
98  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
99  } else if (bc_ptr_w[0].hi(0) == ERFBCType::reflect_even) {
100  dest_arr(i,j,k) = dest_arr(iflip,j,k);
101  } else if (bc_ptr_w[0].hi(0) == ERFBCType::reflect_odd) {
102  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
103  }
104  }
105  );
106  }
107 
108  // First do all ext_dir bcs
109  if (!is_periodic_in_y)
110  {
111  Real* zvel_bc_ptr = m_w_bc_data;
112  Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
113  Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
114  ParallelFor(bx_ylo, bx_yhi,
115  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
116  int jflip = dom_lo.y - 1 - j;
117  if ( (bc_ptr_w[0].lo(1) == ERFBCType::ext_dir) ||
118  (bc_ptr_w[0].lo(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_lo.y,k) >= 0.) )
119  {
120  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
121  if (l_use_terrain_fitted_coords) {
122  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),xvel_arr,yvel_arr,z_phys_nd,dxInv);
123  }
124  } else if (bc_ptr_w[0].lo(1) == ERFBCType::foextrap) {
125  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
126  } else if (bc_ptr_w[0].lo(1) == ERFBCType::open) {
127  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
128  } else if (bc_ptr_w[0].lo(1) == ERFBCType::reflect_even) {
129  dest_arr(i,j,k) = dest_arr(i,jflip,k);
130  } else if (bc_ptr_w[0].lo(1) == ERFBCType::reflect_odd) {
131  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
132  }
133  },
134  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
135  int jflip = 2*dom_hi.y + 1 - j;
136  if ( (bc_ptr_w[0].hi(1) == ERFBCType::ext_dir) ||
137  (bc_ptr_w[0].hi(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_hi.y+1,k) <= 0.) )
138  {
139  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
140  if (l_use_terrain_fitted_coords) {
141  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),xvel_arr,yvel_arr,z_phys_nd,dxInv);
142  }
143  } else if (bc_ptr_w[0].hi(1) == ERFBCType::foextrap) {
144  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
145  } else if (bc_ptr_w[0].hi(1) == ERFBCType::open) {
146  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
147  } else if (bc_ptr_w[0].hi(1) == ERFBCType::reflect_even) {
148  dest_arr(i,j,k) = dest_arr(i,jflip,k);
149  } else if (bc_ptr_w[0].hi(1) == ERFBCType::reflect_odd) {
150  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
151  }
152  });
153  }
154  Gpu::streamSynchronize();
155 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega(int i, int j, int k, amrex::Real omega, amrex::Real u, amrex::Real v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:414
void impose_lateral_zvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Array4< amrex::Real const > &xvel_arr, const amrex::Array4< amrex::Real const > &yvel_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, TerrainType terrain_type, int bccomp)
Definition: ERF_BoundaryConditionsZvel.cpp:16
@ open
Definition: ERF_IndexDefines.H:197
@ reflect_odd
Definition: ERF_IndexDefines.H:187
@ foextrap
Definition: ERF_IndexDefines.H:190
@ ext_dir
Definition: ERF_IndexDefines.H:191
@ ext_dir_upwind
Definition: ERF_IndexDefines.H:199
@ reflect_even
Definition: ERF_IndexDefines.H:189
Here is the call graph for this function:

◆ impose_vertical_zvel_bcs()

void ERFPhysBCFunct_w::impose_vertical_zvel_bcs ( const amrex::Array4< amrex::Real > &  dest_arr,
const amrex::Array4< amrex::Real const > &  xvel_arr,
const amrex::Array4< amrex::Real const > &  yvel_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_u,
int  bccomp_v,
int  bccomp_w,
TerrainType  terrain_type 
)
179 {
180  BL_PROFILE_VAR("impose_vertical_zvel_bcs()",impose_vertical_zvel_bcs);
181  const auto& dom_lo = lbound(domain);
182  const auto& dom_hi = ubound(domain);
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  // Based on BCRec for the domain, we need to make BCRec for this Box
192  // bccomp is used as starting index for m_domain_bcs_type
193  // 0 is used as starting index for bcrs
194  Vector<BCRec> bcrs_u(1), bcrs_v(1), bcrs_w(1);
195  setBC(enclosedCells(bx), domain, bccomp_u, 0, 1, m_domain_bcs_type, bcrs_u);
196  setBC(enclosedCells(bx), domain, bccomp_v, 0, 1, m_domain_bcs_type, bcrs_v);
197  setBC(enclosedCells(bx), domain, bccomp_w, 0, 1, m_domain_bcs_type, bcrs_w);
198 
199  // We use these for the asserts below
200  const BCRec* bc_ptr_u_h = bcrs_u.data();
201  const BCRec* bc_ptr_v_h = bcrs_v.data();
202  const BCRec* bc_ptr_w_h = bcrs_w.data();
203 
204  bool l_use_terrain_fitted_coords = ( (terrain_type == TerrainType::StaticFittedMesh) ||
205  (terrain_type == TerrainType::MovingFittedMesh) );
206  bool l_moving_terrain = (terrain_type == TerrainType::MovingFittedMesh);
207 
208  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,1> l_bc_extdir_vals_d;
209 
210  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
211  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp_w][ori];
212  }
213 
214  // *******************************************************
215  // Bottom boundary
216  // *******************************************************
217 
218  // *******************************************************
219  // Moving terrain
220  // *******************************************************
221  if (l_moving_terrain)
222  {
223  //************************************************************
224  // NOTE: z_t depends on the time interval in which it is
225  // evaluated so we can't arbitrarily define it at a
226  // given time, we must specify an interval
227  //************************************************************
228 
229  // Static terrain
230  } else if (l_use_terrain_fitted_coords) {
231 
232  if (m_lev == 0) {
233  AMREX_ALWAYS_ASSERT( (bc_ptr_u_h[0].lo(2) == ERFBCType::ext_dir && bc_ptr_v_h[0].lo(2) == ERFBCType::ext_dir) ||
234  (bc_ptr_u_h[0].lo(2) != ERFBCType::ext_dir && bc_ptr_v_h[0].lo(2) != ERFBCType::ext_dir) );
235  } else {
236  // If we do not reach to the top or bottom boundary then the z-vel should be
237  // filled by interpolation from the coarser grid using ERF_FillPatcher.
238  }
239  if (bx.smallEnd(2) == dom_lo.z) {
240  ParallelFor(makeSlab(bx,2,dom_lo.z), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
241  dest_arr(i,j,k) = WFromOmega(i,j,k,l_bc_extdir_vals_d[0][2],xvel_arr,yvel_arr,z_phys_nd,dxInv);
242  });
243  }
244 
245  // No terrain
246  } else {
247  if (bx.smallEnd(2) == dom_lo.z) {
248  ParallelFor(makeSlab(bx,2,dom_lo.z), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
249  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][2];
250  });
251  }
252  }
253 
254  // *******************************************************
255  // Top boundary
256  // *******************************************************
257 
258  // NOTE: if we set SlipWall at top, that generates ERFBCType::ext_dir which sets w=0 here
259  // NOTE: if we set Outflow at top, that generates ERFBCType::foextrap which doesn't touch w here
260  if (bx.bigEnd(2) == dom_hi.z+1) {
261  if (bc_ptr_w_h[0].hi(2) == ERFBCType::ext_dir) {
262  ParallelFor(makeSlab(bx,2,dom_hi.z+1), [=] AMREX_GPU_DEVICE (int i, int j, int k)
263  {
264  if (l_use_terrain_fitted_coords) {
265  dest_arr(i,j,k) = WFromOmega(i,j,k,l_bc_extdir_vals_d[0][5],xvel_arr,yvel_arr,z_phys_nd,dxInv);
266  } else {
267  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][5];
268  }
269  });
270  } else if (bc_ptr_w_h[0].hi(2) == ERFBCType::neumann_int) {
271  ParallelFor(makeSlab(bx,2,dom_hi.z+1), [=] AMREX_GPU_DEVICE (int i, int j, int k)
272  {
273  dest_arr(i,j,k) = (4.0*dest_arr(i,j,dom_hi.z) - dest_arr(i,j,dom_hi.z-1))/3.0;
274  });
275  }
276  }
277  Gpu::streamSynchronize();
278 }
void impose_vertical_zvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Array4< amrex::Real const > &xvel_arr, const amrex::Array4< amrex::Real const > &yvel_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_u, int bccomp_v, int bccomp_w, TerrainType terrain_type)
Definition: ERF_BoundaryConditionsZvel.cpp:171
@ neumann_int
Definition: ERF_IndexDefines.H:196
Here is the call graph for this function:

◆ operator()()

void ERFPhysBCFunct_w::operator() ( amrex::MultiFab &  mf,
amrex::MultiFab &  xvel,
amrex::MultiFab &  yvel,
amrex::IntVect const &  nghost,
const amrex::Real  time,
int  bccomp,
bool  do_fb 
)
238 {
239  BL_PROFILE("ERFPhysBCFunct_w::()");
240 
241  int bccomp_u = BCVars::xvel_bc;
242  int bccomp_v = BCVars::yvel_bc;
243 
244  if (m_geom.isAllPeriodic()) return;
245 
246  const auto& domain = m_geom.Domain();
247  const auto dxInv = m_geom.InvCellSizeArray();
248 
249  Box gdomainz = surroundingNodes(domain,2);
250  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
251  if (m_geom.isPeriodic(i)) {
252  gdomainz.grow(i, nghost[i]);
253  }
254  }
255  //
256  // We want to make sure we impose the z-vels at k=0 if the box includes k=0
257  //
258  if (gdomainz.smallEnd(2) == 0) gdomainz.setSmall(2,1);
259 
260  //
261  // We fill all of the interior and periodic ghost cells first, so we can fill
262  // those directly inside the lateral and vertical calls.
263  //
264  if (do_fb) {
265  mf.FillBoundary(m_geom.periodicity());
266  }
267 
268 #ifdef AMREX_USE_OMP
269 #pragma omp parallel if (Gpu::notInLaunchRegion())
270 #endif
271  {
272  for (MFIter mfi(mf,false); mfi.isValid(); ++mfi)
273  {
274  //
275  // This is the box we pass to the different routines
276  // NOTE -- this is the full grid NOT the tile box
277  //
278  Box bx = mfi.validbox();
279 
280  //
281  // These are the boxes we use to test on relative to the domain
282  //
283  Box zbx = surroundingNodes(bx,2); zbx.grow(nghost);
284  if (zbx.smallEnd(2) < domain.smallEnd(2)) zbx.setSmall(2,domain.smallEnd(2));
285  if (zbx.bigEnd(2) > domain.bigEnd(2)) zbx.setBig(2,domain.bigEnd(2)+1);
286 
287  Array4<const Real> z_nd_arr;
288 
289  if (m_z_phys_nd)
290  {
291  z_nd_arr = m_z_phys_nd->const_array(mfi);
292  }
293 
294  //
295  // Recall that gdomainz.smallEnd(2) = 1 not 0!
296  //
297  if (!gdomainz.contains(zbx))
298  {
299  Array4<const Real> const& velx_arr = xvel.const_array(mfi);
300  Array4<const Real> const& vely_arr = yvel.const_array(mfi);
301  Array4< Real> const& velz_arr = mf.array(mfi);
302 
303  if (!m_use_real_bcs)
304  {
305  if (!gdomainz.contains(zbx))
306  {
307  impose_lateral_zvel_bcs(velz_arr,velx_arr,vely_arr,zbx,domain,z_nd_arr,
308  dxInv,m_terrain_type,bccomp_w);
309  }
310  }
311 
312  impose_vertical_zvel_bcs(velz_arr,velx_arr,vely_arr,zbx,domain,z_nd_arr,dxInv,
313  bccomp_u, bccomp_v, bccomp_w, m_terrain_type);
314  }
315  } // MFIter
316  } // OpenMP
317 } // operator()
@ yvel_bc
Definition: ERF_IndexDefines.H:88
@ xvel_bc
Definition: ERF_IndexDefines.H:87
@ xvel
Definition: ERF_IndexDefines.H:141
@ yvel
Definition: ERF_IndexDefines.H:142

Member Data Documentation

◆ m_bc_extdir_vals

amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NBCVAR_max> ERFPhysBCFunct_w::m_bc_extdir_vals
private

◆ m_bc_neumann_vals

amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NBCVAR_max> ERFPhysBCFunct_w::m_bc_neumann_vals
private

◆ m_domain_bcs_type

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

◆ m_domain_bcs_type_d

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

◆ m_geom

amrex::Geometry ERFPhysBCFunct_w::m_geom
private

◆ m_lev

int ERFPhysBCFunct_w::m_lev
private

◆ m_terrain_type

TerrainType ERFPhysBCFunct_w::m_terrain_type
private

◆ m_use_real_bcs

bool ERFPhysBCFunct_w::m_use_real_bcs
private

◆ m_w_bc_data

amrex::Real* ERFPhysBCFunct_w::m_w_bc_data
private

◆ m_z_phys_nd

amrex::MultiFab* ERFPhysBCFunct_w::m_z_phys_nd
private

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