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, const amrex::Real time)
 
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, const amrex::Real time)
 

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
218  : m_lev(lev),
219  m_geom(geom), m_domain_bcs_type(domain_bcs_type),
220  m_domain_bcs_type_d(domain_bcs_type_d),
221  m_bc_extdir_vals(bc_extdir_vals),
222  m_bc_neumann_vals(bc_neumann_vals),
223  m_terrain_type(terrain_type),
224  m_z_phys_nd(z_phys_nd.get()),
225  m_use_real_bcs(use_real_bcs),
226  m_w_bc_data(w_bc_data)
227  {
228  m_mapfac_u = mapfac_lev[MapFacType::u_x].get();
229  m_mapfac_v = mapfac_lev[MapFacType::v_y].get();
230  }
@ 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:273
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:272
int m_lev
Definition: ERF_PhysBCFunct.H:271
amrex::Real * m_w_bc_data
Definition: ERF_PhysBCFunct.H:282
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:281
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:274
amrex::MultiFab * m_mapfac_v
Definition: ERF_PhysBCFunct.H:279
TerrainType m_terrain_type
Definition: ERF_PhysBCFunct.H:277
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:280
amrex::MultiFab * m_mapfac_u
Definition: ERF_PhysBCFunct.H:278
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:276
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:275

◆ ~ERFPhysBCFunct_w()

ERFPhysBCFunct_w::~ERFPhysBCFunct_w ( )
inline
232 {}

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,
const amrex::Real  time 
)
26 {
27  BL_PROFILE_VAR("impose_lateral_zvel_bcs()",impose_lateral_zvel_bcs);
28  const auto& dom_lo = lbound(domain);
29  const auto& dom_hi = ubound(domain);
30 
31  // Based on BCRec for the domain, we need to make BCRec for this Box
32  // bccomp is used as starting index for m_domain_bcs_type
33  // 0 is used as starting index for bcrs
34  Vector<BCRec> bcrs_w(1);
35  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs_w);
36 
37  bool l_use_terrain_fitted_coords = ( (terrain_type == TerrainType::StaticFittedMesh) ||
38  (terrain_type == TerrainType::MovingFittedMesh) );
39 
40  // xlo: ori = 0
41  // ylo: ori = 1
42  // zlo: ori = 2
43  // xhi: ori = 3
44  // yhi: ori = 4
45  // zhi: ori = 5
46 
47  Gpu::DeviceVector<BCRec> bcrs_w_d(1);
48  Gpu::copyAsync(Gpu::hostToDevice, bcrs_w.begin(), bcrs_w.end(), bcrs_w_d.begin());
49  const BCRec* bc_ptr_w = bcrs_w_d.data();
50 
51  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,1> l_bc_extdir_vals_d;
52 
53  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
54  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp][ori];
55  }
56 
57  GeometryData const& geomdata = m_geom.data();
58  bool is_periodic_in_x = geomdata.isPeriodic(0);
59  bool is_periodic_in_y = geomdata.isPeriodic(1);
60 
61  FArrayBox dhdtfab;
62 
63  // First do all ext_dir bcs
64  if (!is_periodic_in_x)
65  {
66  Real* zvel_bc_ptr = m_w_bc_data;
67  Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
68  Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
69  ParallelFor(bx_xlo, bx_xhi,
70  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
71  int iflip = dom_lo.x - 1 - i;
72  if ( (bc_ptr_w[0].lo(0) == ERFBCType::ext_dir) ||
73  (bc_ptr_w[0].lo(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_lo.x,j,k) >= 0.) )
74  {
75  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
76  if (l_use_terrain_fitted_coords) {
77  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
78  xvel_arr,yvel_arr,
79  mf_u,mf_v,z_phys_nd,dxInv);
80  }
81  } else if (bc_ptr_w[0].lo(0) == ERFBCType::foextrap) {
82  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
83  } else if (bc_ptr_w[0].lo(0) == ERFBCType::open) {
84  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
85  } else if (bc_ptr_w[0].lo(0) == ERFBCType::reflect_even) {
86  dest_arr(i,j,k) = dest_arr(iflip,j,k);
87  } else if (bc_ptr_w[0].lo(0) == ERFBCType::reflect_odd) {
88  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
89  }
90  },
91  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
92  int iflip = 2*dom_hi.x + 1 - i;
93  if ( (bc_ptr_w[0].hi(0) == ERFBCType::ext_dir) ||
94  (bc_ptr_w[0].hi(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_hi.x+1,j,k) <= 0.) )
95  {
96  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
97  if (l_use_terrain_fitted_coords) {
98  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
99  xvel_arr,yvel_arr,
100  mf_u,mf_v,z_phys_nd,dxInv);
101  }
102  } else if (bc_ptr_w[0].hi(0) == ERFBCType::foextrap) {
103  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
104  } else if (bc_ptr_w[0].hi(0) == ERFBCType::open) {
105  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
106  } else if (bc_ptr_w[0].hi(0) == ERFBCType::reflect_even) {
107  dest_arr(i,j,k) = dest_arr(iflip,j,k);
108  } else if (bc_ptr_w[0].hi(0) == ERFBCType::reflect_odd) {
109  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
110  }
111  }
112  );
113  }
114 
115  // First do all ext_dir bcs
116  if (!is_periodic_in_y)
117  {
118  Real* zvel_bc_ptr = m_w_bc_data;
119  Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
120  Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
121  ParallelFor(bx_ylo, bx_yhi,
122  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
123  int jflip = dom_lo.y - 1 - j;
124  if ( (bc_ptr_w[0].lo(1) == ERFBCType::ext_dir) ||
125  (bc_ptr_w[0].lo(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_lo.y,k) >= 0.) )
126  {
127  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
128  if (l_use_terrain_fitted_coords) {
129  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
130  xvel_arr,yvel_arr,
131  mf_u,mf_v,z_phys_nd,dxInv);
132  }
133  } else if (bc_ptr_w[0].lo(1) == ERFBCType::foextrap) {
134  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
135  } else if (bc_ptr_w[0].lo(1) == ERFBCType::open) {
136  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
137  } else if (bc_ptr_w[0].lo(1) == ERFBCType::reflect_even) {
138  dest_arr(i,j,k) = dest_arr(i,jflip,k);
139  } else if (bc_ptr_w[0].lo(1) == ERFBCType::reflect_odd) {
140  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
141  }
142  },
143  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
144  int jflip = 2*dom_hi.y + 1 - j;
145  if ( (bc_ptr_w[0].hi(1) == ERFBCType::ext_dir) ||
146  (bc_ptr_w[0].hi(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_hi.y+1,k) <= 0.) )
147  {
148  dest_arr(i,j,k) = (zvel_bc_ptr) ? zvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
149  if (l_use_terrain_fitted_coords) {
150  dest_arr(i,j,k) = WFromOmega(i,j,k,dest_arr(i,j,k),
151  xvel_arr,yvel_arr,
152  mf_u,mf_v,z_phys_nd,dxInv);
153  }
154  } else if (bc_ptr_w[0].hi(1) == ERFBCType::foextrap) {
155  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
156  } else if (bc_ptr_w[0].hi(1) == ERFBCType::open) {
157  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
158  } else if (bc_ptr_w[0].hi(1) == ERFBCType::reflect_even) {
159  dest_arr(i,j,k) = dest_arr(i,jflip,k);
160  } else if (bc_ptr_w[0].hi(1) == ERFBCType::reflect_odd) {
161  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
162  }
163  });
164  }
165  Gpu::streamSynchronize();
166 }
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, const amrex::Real time)
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,
const amrex::Real  time 
)
193 {
194  BL_PROFILE_VAR("impose_vertical_zvel_bcs()",impose_vertical_zvel_bcs);
195  const auto& dom_lo = lbound(domain);
196  const auto& dom_hi = ubound(domain);
197 
198  // xlo: ori = 0
199  // ylo: ori = 1
200  // zlo: ori = 2
201  // xhi: ori = 3
202  // yhi: ori = 4
203  // zhi: ori = 5
204 
205  // Based on BCRec for the domain, we need to make BCRec for this Box
206  // bccomp is used as starting index for m_domain_bcs_type
207  // 0 is used as starting index for bcrs
208  Vector<BCRec> bcrs_u(1), bcrs_v(1), bcrs_w(1);
209  setBC(enclosedCells(bx), domain, bccomp_u, 0, 1, m_domain_bcs_type, bcrs_u);
210  setBC(enclosedCells(bx), domain, bccomp_v, 0, 1, m_domain_bcs_type, bcrs_v);
211  setBC(enclosedCells(bx), domain, bccomp_w, 0, 1, m_domain_bcs_type, bcrs_w);
212 
213  // We use these for the asserts below
214  const BCRec* bc_ptr_u_h = bcrs_u.data();
215  const BCRec* bc_ptr_v_h = bcrs_v.data();
216  const BCRec* bc_ptr_w_h = bcrs_w.data();
217 
218  bool l_use_terrain_fitted_coords = ( (terrain_type == TerrainType::StaticFittedMesh) ||
219  (terrain_type == TerrainType::MovingFittedMesh) );
220  bool l_moving_terrain = (terrain_type == TerrainType::MovingFittedMesh);
221 
222  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,1> l_bc_extdir_vals_d;
223 
224  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
225  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp_w][ori];
226  }
227 
228  // *******************************************************
229  // Bottom boundary
230  // *******************************************************
231 
232  // *******************************************************
233  // Moving terrain
234  // *******************************************************
235  if (l_moving_terrain)
236  {
237  //************************************************************
238  // NOTE: z_t depends on the time interval in which it is
239  // evaluated so we can't arbitrarily define it at a
240  // given time, we must specify an interval
241  //************************************************************
242 
243  // Static terrain
244  } else if (l_use_terrain_fitted_coords) {
245 
246  if (m_lev == 0) {
247  AMREX_ALWAYS_ASSERT( (bc_ptr_u_h[0].lo(2) == ERFBCType::ext_dir && bc_ptr_v_h[0].lo(2) == ERFBCType::ext_dir) ||
248  (bc_ptr_u_h[0].lo(2) != ERFBCType::ext_dir && bc_ptr_v_h[0].lo(2) != ERFBCType::ext_dir) );
249  } else {
250  // If we do not reach to the top or bottom boundary then the z-vel should be
251  // filled by interpolation from the coarser grid using ERF_FillPatcher.
252  }
253  if (bx.smallEnd(2) == dom_lo.z) {
254  ParallelFor(makeSlab(bx,2,dom_lo.z), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
255  dest_arr(i,j,k) = WFromOmega(i,j,k,l_bc_extdir_vals_d[0][2],
256  xvel_arr,yvel_arr,
257  mf_u,mf_v,z_phys_nd,dxInv);
258  });
259  }
260 
261  // No terrain
262  } else {
263  if (bx.smallEnd(2) == dom_lo.z) {
264  ParallelFor(makeSlab(bx,2,dom_lo.z), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
265  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][2];
266  });
267  }
268  }
269 
270  // *******************************************************
271  // Top boundary
272  // *******************************************************
273 
274  // NOTE: if we set SlipWall at top, that generates ERFBCType::ext_dir which sets w=0 here
275  // NOTE: if we set Outflow at top, that generates ERFBCType::foextrap which doesn't touch w here
276  if (bx.bigEnd(2) == dom_hi.z+1) {
277  if (bc_ptr_w_h[0].hi(2) == ERFBCType::ext_dir) {
278  ParallelFor(makeSlab(bx,2,dom_hi.z+1), [=] AMREX_GPU_DEVICE (int i, int j, int k)
279  {
280  if (l_use_terrain_fitted_coords) {
281  dest_arr(i,j,k) = WFromOmega(i,j,k,l_bc_extdir_vals_d[0][5],
282  xvel_arr,yvel_arr,
283  mf_u,mf_v,z_phys_nd,dxInv);
284  } else {
285  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][5];
286  }
287  });
288  } else if (bc_ptr_w_h[0].hi(2) == ERFBCType::neumann_int) {
289  ParallelFor(makeSlab(bx,2,dom_hi.z+1), [=] AMREX_GPU_DEVICE (int i, int j, int k)
290  {
291  dest_arr(i,j,k) = (4.0*dest_arr(i,j,dom_hi.z) - dest_arr(i,j,dom_hi.z-1))/3.0;
292  });
293  }
294  }
295  Gpu::streamSynchronize();
296 }
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, const amrex::Real time)
Definition: ERF_BoundaryConditionsZvel.cpp:182
@ 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,time);
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, time);
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: