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, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &mapfac_lev, 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 > &mf_u, const amrex::Array4< amrex::Real const > &mf_v, 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 > &mf_u, const amrex::Array4< amrex::Real const > &mf_v, 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_mapfac_u
 
amrex::MultiFab * m_mapfac_v
 
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,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  mapfac_lev,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
const bool  use_real_bcs,
amrex::Real *  w_bc_data 
)
inline
213  : m_lev(lev),
214  m_geom(geom), m_domain_bcs_type(domain_bcs_type),
215  m_domain_bcs_type_d(domain_bcs_type_d),
216  m_bc_extdir_vals(bc_extdir_vals),
217  m_bc_neumann_vals(bc_neumann_vals),
218  m_terrain_type(terrain_type),
219  m_z_phys_nd(z_phys_nd.get()),
220  m_use_real_bcs(use_real_bcs),
221  m_w_bc_data(w_bc_data)
222  {
223  m_mapfac_u = mapfac_lev[MapFacType::u_x].get();
224  m_mapfac_v = mapfac_lev[MapFacType::v_y].get();
225  }
@ v_y
Definition: ERF_DataStruct.H:23
@ u_x
Definition: ERF_DataStruct.H:22
amrex::Vector< amrex::BCRec > m_domain_bcs_type
Definition: ERF_PhysBCFunct.H:266
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:265
int m_lev
Definition: ERF_PhysBCFunct.H:264
amrex::Real * m_w_bc_data
Definition: ERF_PhysBCFunct.H:275
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:274
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:267
amrex::MultiFab * m_mapfac_v
Definition: ERF_PhysBCFunct.H:272
TerrainType m_terrain_type
Definition: ERF_PhysBCFunct.H:270
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:273
amrex::MultiFab * m_mapfac_u
Definition: ERF_PhysBCFunct.H:271
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:269
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:268

◆ ~ERFPhysBCFunct_w()

ERFPhysBCFunct_w::~ERFPhysBCFunct_w ( )
inline
227 {}

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 > &  mf_u,
const amrex::Array4< amrex::Real const > &  mf_v,
const amrex::Array4< amrex::Real const > &  z_nd,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM >  dxInv,
TerrainType  terrain_type,
int  bccomp 
)
25 {
26  BL_PROFILE_VAR("impose_lateral_zvel_bcs()",impose_lateral_zvel_bcs);
27  const auto& dom_lo = lbound(domain);
28  const auto& dom_hi = ubound(domain);
29 
30  // Based on BCRec for the domain, we need to make BCRec for this Box
31  // bccomp is used as starting index for m_domain_bcs_type
32  // 0 is used as starting index for bcrs
33  Vector<BCRec> bcrs_w(1);
34  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs_w);
35 
36  bool l_use_terrain_fitted_coords = ( (terrain_type == TerrainType::StaticFittedMesh) ||
37  (terrain_type == TerrainType::MovingFittedMesh) );
38 
39  // xlo: ori = 0
40  // ylo: ori = 1
41  // zlo: ori = 2
42  // xhi: ori = 3
43  // yhi: ori = 4
44  // zhi: ori = 5
45 
46  Gpu::DeviceVector<BCRec> bcrs_w_d(1);
47  Gpu::copyAsync(Gpu::hostToDevice, bcrs_w.begin(), bcrs_w.end(), bcrs_w_d.begin());
48  const BCRec* bc_ptr_w = bcrs_w_d.data();
49 
50  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,1> l_bc_extdir_vals_d;
51 
52  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
53  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp][ori];
54  }
55 
56  GeometryData const& geomdata = m_geom.data();
57  bool is_periodic_in_x = geomdata.isPeriodic(0);
58  bool is_periodic_in_y = geomdata.isPeriodic(1);
59 
60  FArrayBox dhdtfab;
61 
62  // First do all ext_dir bcs
63  if (!is_periodic_in_x)
64  {
65  Real* zvel_bc_ptr = m_w_bc_data;
66  Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
67  Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
68  ParallelFor(bx_xlo, bx_xhi,
69  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
70  int iflip = dom_lo.x - 1 - i;
71  if ( (bc_ptr_w[0].lo(0) == ERFBCType::ext_dir) ||
72  (bc_ptr_w[0].lo(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_lo.x,j,k) >= 0.) )
73  {
74  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
75  if (l_use_terrain_fitted_coords) {
76  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
77  xvel_arr,yvel_arr,
78  mf_u,mf_v,z_phys_nd,dxInv);
79  }
80  } else if (bc_ptr_w[0].lo(0) == ERFBCType::foextrap) {
81  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
82  } else if (bc_ptr_w[0].lo(0) == ERFBCType::open) {
83  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
84  } else if (bc_ptr_w[0].lo(0) == ERFBCType::reflect_even) {
85  dest_arr(i,j,k) = dest_arr(iflip,j,k);
86  } else if (bc_ptr_w[0].lo(0) == ERFBCType::reflect_odd) {
87  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
88  }
89  },
90  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
91  int iflip = 2*dom_hi.x + 1 - i;
92  if ( (bc_ptr_w[0].hi(0) == ERFBCType::ext_dir) ||
93  (bc_ptr_w[0].hi(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_hi.x+1,j,k) <= 0.) )
94  {
95  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
96  if (l_use_terrain_fitted_coords) {
97  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
98  xvel_arr,yvel_arr,
99  mf_u,mf_v,z_phys_nd,dxInv);
100  }
101  } else if (bc_ptr_w[0].hi(0) == ERFBCType::foextrap) {
102  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
103  } else if (bc_ptr_w[0].hi(0) == ERFBCType::open) {
104  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
105  } else if (bc_ptr_w[0].hi(0) == ERFBCType::reflect_even) {
106  dest_arr(i,j,k) = dest_arr(iflip,j,k);
107  } else if (bc_ptr_w[0].hi(0) == ERFBCType::reflect_odd) {
108  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
109  }
110  }
111  );
112  }
113 
114  // First do all ext_dir bcs
115  if (!is_periodic_in_y)
116  {
117  Real* zvel_bc_ptr = m_w_bc_data;
118  Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
119  Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
120  ParallelFor(bx_ylo, bx_yhi,
121  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
122  int jflip = dom_lo.y - 1 - j;
123  if ( (bc_ptr_w[0].lo(1) == ERFBCType::ext_dir) ||
124  (bc_ptr_w[0].lo(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_lo.y,k) >= 0.) )
125  {
126  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
127  if (l_use_terrain_fitted_coords) {
128  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
129  xvel_arr,yvel_arr,
130  mf_u,mf_v,z_phys_nd,dxInv);
131  }
132  } else if (bc_ptr_w[0].lo(1) == ERFBCType::foextrap) {
133  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
134  } else if (bc_ptr_w[0].lo(1) == ERFBCType::open) {
135  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
136  } else if (bc_ptr_w[0].lo(1) == ERFBCType::reflect_even) {
137  dest_arr(i,j,k) = dest_arr(i,jflip,k);
138  } else if (bc_ptr_w[0].lo(1) == ERFBCType::reflect_odd) {
139  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
140  }
141  },
142  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
143  int jflip = 2*dom_hi.y + 1 - j;
144  if ( (bc_ptr_w[0].hi(1) == ERFBCType::ext_dir) ||
145  (bc_ptr_w[0].hi(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_hi.y+1,k) <= 0.) )
146  {
147  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
148  if (l_use_terrain_fitted_coords) {
149  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
150  xvel_arr,yvel_arr,
151  mf_u,mf_v,z_phys_nd,dxInv);
152  }
153  } else if (bc_ptr_w[0].hi(1) == ERFBCType::foextrap) {
154  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
155  } else if (bc_ptr_w[0].hi(1) == ERFBCType::open) {
156  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
157  } else if (bc_ptr_w[0].hi(1) == ERFBCType::reflect_even) {
158  dest_arr(i,j,k) = dest_arr(i,jflip,k);
159  } else if (bc_ptr_w[0].hi(1) == ERFBCType::reflect_odd) {
160  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
161  }
162  });
163  }
164  Gpu::streamSynchronize();
165 }
const auto & dom_hi
Definition: ERF_DiffSetup.H:10
const auto & dom_lo
Definition: ERF_DiffSetup.H:9
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega(int &i, int &j, int &k, amrex::Real omega, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:465
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 > &mf_u, const amrex::Array4< amrex::Real const > &mf_v, 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:215
@ reflect_odd
Definition: ERF_IndexDefines.H:205
@ foextrap
Definition: ERF_IndexDefines.H:208
@ ext_dir
Definition: ERF_IndexDefines.H:209
@ ext_dir_upwind
Definition: ERF_IndexDefines.H:217
@ reflect_even
Definition: ERF_IndexDefines.H:207
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 > &  mf_u,
const amrex::Array4< amrex::Real const > &  mf_v,
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 
)
191 {
192  BL_PROFILE_VAR("impose_vertical_zvel_bcs()",impose_vertical_zvel_bcs);
193  const auto& dom_lo = lbound(domain);
194  const auto& dom_hi = ubound(domain);
195 
196  // xlo: ori = 0
197  // ylo: ori = 1
198  // zlo: ori = 2
199  // xhi: ori = 3
200  // yhi: ori = 4
201  // zhi: ori = 5
202 
203  // Based on BCRec for the domain, we need to make BCRec for this Box
204  // bccomp is used as starting index for m_domain_bcs_type
205  // 0 is used as starting index for bcrs
206  Vector<BCRec> bcrs_u(1), bcrs_v(1), bcrs_w(1);
207  setBC(enclosedCells(bx), domain, bccomp_u, 0, 1, m_domain_bcs_type, bcrs_u);
208  setBC(enclosedCells(bx), domain, bccomp_v, 0, 1, m_domain_bcs_type, bcrs_v);
209  setBC(enclosedCells(bx), domain, bccomp_w, 0, 1, m_domain_bcs_type, bcrs_w);
210 
211  // We use these for the asserts below
212  const BCRec* bc_ptr_u_h = bcrs_u.data();
213  const BCRec* bc_ptr_v_h = bcrs_v.data();
214  const BCRec* bc_ptr_w_h = bcrs_w.data();
215 
216  bool l_use_terrain_fitted_coords = ( (terrain_type == TerrainType::StaticFittedMesh) ||
217  (terrain_type == TerrainType::MovingFittedMesh) );
218  bool l_moving_terrain = (terrain_type == TerrainType::MovingFittedMesh);
219 
220  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,1> l_bc_extdir_vals_d;
221 
222  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
223  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp_w][ori];
224  }
225 
226  // *******************************************************
227  // Bottom boundary
228  // *******************************************************
229 
230  // *******************************************************
231  // Moving terrain
232  // *******************************************************
233  if (l_moving_terrain)
234  {
235  //************************************************************
236  // NOTE: z_t depends on the time interval in which it is
237  // evaluated so we can't arbitrarily define it at a
238  // given time, we must specify an interval
239  //************************************************************
240 
241  // Static terrain
242  } else if (l_use_terrain_fitted_coords) {
243 
244  if (m_lev == 0) {
245  AMREX_ALWAYS_ASSERT( (bc_ptr_u_h[0].lo(2) == ERFBCType::ext_dir && bc_ptr_v_h[0].lo(2) == ERFBCType::ext_dir) ||
246  (bc_ptr_u_h[0].lo(2) != ERFBCType::ext_dir && bc_ptr_v_h[0].lo(2) != ERFBCType::ext_dir) );
247  } else {
248  // If we do not reach to the top or bottom boundary then the z-vel should be
249  // filled by interpolation from the coarser grid using ERF_FillPatcher.
250  }
251  if (bx.smallEnd(2) == dom_lo.z) {
252  ParallelFor(makeSlab(bx,2,dom_lo.z), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
253  dest_arr(i,j,k) = WFromOmega(i,j,k,l_bc_extdir_vals_d[0][2],
254  xvel_arr,yvel_arr,
255  mf_u,mf_v,z_phys_nd,dxInv);
256  });
257  }
258 
259  // No terrain
260  } else {
261  if (bx.smallEnd(2) == dom_lo.z) {
262  ParallelFor(makeSlab(bx,2,dom_lo.z), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
263  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][2];
264  });
265  }
266  }
267 
268  // *******************************************************
269  // Top boundary
270  // *******************************************************
271 
272  // NOTE: if we set SlipWall at top, that generates ERFBCType::ext_dir which sets w=0 here
273  // NOTE: if we set Outflow at top, that generates ERFBCType::foextrap which doesn't touch w here
274  if (bx.bigEnd(2) == dom_hi.z+1) {
275  if (bc_ptr_w_h[0].hi(2) == ERFBCType::ext_dir) {
276  ParallelFor(makeSlab(bx,2,dom_hi.z+1), [=] AMREX_GPU_DEVICE (int i, int j, int k)
277  {
278  if (l_use_terrain_fitted_coords) {
279  dest_arr(i,j,k) = WFromOmega(i,j,k,l_bc_extdir_vals_d[0][5],
280  xvel_arr,yvel_arr,
281  mf_u,mf_v,z_phys_nd,dxInv);
282  } else {
283  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][5];
284  }
285  });
286  } else if (bc_ptr_w_h[0].hi(2) == ERFBCType::neumann_int) {
287  ParallelFor(makeSlab(bx,2,dom_hi.z+1), [=] AMREX_GPU_DEVICE (int i, int j, int k)
288  {
289  dest_arr(i,j,k) = (4.0*dest_arr(i,j,dom_hi.z) - dest_arr(i,j,dom_hi.z-1))/3.0;
290  });
291  }
292  }
293  Gpu::streamSynchronize();
294 }
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 > &mf_u, const amrex::Array4< amrex::Real const > &mf_v, 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:181
@ neumann_int
Definition: ERF_IndexDefines.H:214
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 
)
241 {
242  BL_PROFILE("ERFPhysBCFunct_w::()");
243 
244  //
245  // We fill all of the interior and periodic ghost cells first, so we can fill
246  // those directly inside the lateral and vertical calls.
247  // If triply periodic this is all we do
248  //
249  if (do_fb) {
250  mf.FillBoundary(m_geom.periodicity());
251  }
252 
253  if (m_geom.isAllPeriodic()) return;
254 
255  int bccomp_u = BCVars::xvel_bc;
256  int bccomp_v = BCVars::yvel_bc;
257 
258  const auto& domain = m_geom.Domain();
259  const auto dxInv = m_geom.InvCellSizeArray();
260 
261  Box gdomainz = surroundingNodes(domain,2);
262  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
263  if (m_geom.isPeriodic(i)) {
264  gdomainz.grow(i, nghost[i]);
265  }
266  }
267  //
268  // We want to make sure we impose the z-vels at k=0 if the box includes k=0
269  //
270  if (gdomainz.smallEnd(2) == 0) gdomainz.setSmall(2,1);
271 
272 #ifdef AMREX_USE_OMP
273 #pragma omp parallel if (Gpu::notInLaunchRegion())
274 #endif
275  {
276  for (MFIter mfi(mf,false); mfi.isValid(); ++mfi)
277  {
278  //
279  // This is the box we pass to the different routines
280  // NOTE -- this is the full grid NOT the tile box
281  //
282  Box bx = mfi.validbox();
283 
284  //
285  // These are the boxes we use to test on relative to the domain
286  //
287  Box zbx = surroundingNodes(bx,2); zbx.grow(nghost);
288  if (zbx.smallEnd(2) < domain.smallEnd(2)) zbx.setSmall(2,domain.smallEnd(2));
289  if (zbx.bigEnd(2) > domain.bigEnd(2)) zbx.setBig(2,domain.bigEnd(2)+1);
290 
291  Array4<const Real> z_nd_arr;
292  const Array4<const Real>& mf_u = m_mapfac_u->const_array(mfi);
293  const Array4<const Real>& mf_v = m_mapfac_v->const_array(mfi);
294 
295  if (m_z_phys_nd)
296  {
297  z_nd_arr = m_z_phys_nd->const_array(mfi);
298  }
299 
300  //
301  // Recall that gdomainz.smallEnd(2) = 1 not 0!
302  //
303  if (!gdomainz.contains(zbx))
304  {
305  Array4<const Real> const& velx_arr = xvel.const_array(mfi);
306  Array4<const Real> const& vely_arr = yvel.const_array(mfi);
307  Array4< Real> const& velz_arr = mf.array(mfi);
308 
309  if (!m_use_real_bcs)
310  {
311  if (!gdomainz.contains(zbx))
312  {
313  impose_lateral_zvel_bcs(velz_arr,velx_arr,vely_arr,zbx,domain,
314  mf_u,mf_v,z_nd_arr,dxInv,m_terrain_type,bccomp_w);
315  }
316  }
317 
318  impose_vertical_zvel_bcs(velz_arr,velx_arr,vely_arr,zbx,domain,mf_u,mf_v,
319  z_nd_arr,dxInv,bccomp_u, bccomp_v, bccomp_w, m_terrain_type);
320  }
321  } // MFIter
322  } // OpenMP
323 } // operator()
const Box zbx
Definition: ERF_DiffSetup.H:23
@ 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_mapfac_u

amrex::MultiFab* ERFPhysBCFunct_w::m_mapfac_u
private

Referenced by ERFPhysBCFunct_w().

◆ m_mapfac_v

amrex::MultiFab* ERFPhysBCFunct_w::m_mapfac_v
private

Referenced by ERFPhysBCFunct_w().

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