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

#include <ERF_PhysBCFunct.H>

Collaboration diagram for ERFPhysBCFunct_u:

Public Member Functions

 ERFPhysBCFunct_u (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, std::unique_ptr< amrex::MultiFab > &z_phys_nd, const bool use_real_bcs, amrex::Real *u_bc_data)
 
 ~ERFPhysBCFunct_u ()
 
void operator() (amrex::MultiFab &mf, int icomp, int ncomp, amrex::IntVect const &nghost, const amrex::Real time, int bccomp, bool do_fb)
 
void impose_lateral_xvel_bcs (const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, int bccomp)
 
void impose_vertical_xvel_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, 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
 
amrex::MultiFab * m_z_phys_nd
 
bool m_use_real_bcs
 
amrex::Real * m_u_bc_data
 

Constructor & Destructor Documentation

◆ ERFPhysBCFunct_u()

ERFPhysBCFunct_u::ERFPhysBCFunct_u ( 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,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
const bool  use_real_bcs,
amrex::Real *  u_bc_data 
)
inline
89  : m_lev(lev), m_geom(geom),
90  m_domain_bcs_type(domain_bcs_type),
91  m_domain_bcs_type_d(domain_bcs_type_d),
92  m_bc_extdir_vals(bc_extdir_vals),
93  m_bc_neumann_vals(bc_neumann_vals),
94  m_z_phys_nd(z_phys_nd.get()),
95  m_use_real_bcs(use_real_bcs),
96  m_u_bc_data(u_bc_data)
97  { }
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:125
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:128
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:127
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:131
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:129
amrex::Real * m_u_bc_data
Definition: ERF_PhysBCFunct.H:132
int m_lev
Definition: ERF_PhysBCFunct.H:124
amrex::Vector< amrex::BCRec > m_domain_bcs_type
Definition: ERF_PhysBCFunct.H:126
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:130

◆ ~ERFPhysBCFunct_u()

ERFPhysBCFunct_u::~ERFPhysBCFunct_u ( )
inline
99 {}

Member Function Documentation

◆ impose_lateral_xvel_bcs()

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

◆ impose_vertical_xvel_bcs()

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

◆ operator()()

void ERFPhysBCFunct_u::operator() ( amrex::MultiFab &  mf,
int  icomp,
int  ncomp,
amrex::IntVect const &  nghost,
const amrex::Real  time,
int  bccomp,
bool  do_fb 
)
89 {
90  BL_PROFILE("ERFPhysBCFunct_u::()");
91 
92  if (m_geom.isAllPeriodic()) return;
93 
94  const auto& domain = m_geom.Domain();
95  const auto dxInv = m_geom.InvCellSizeArray();
96 
97  Box gdomainx = surroundingNodes(domain,0);
98  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
99  if (m_geom.isPeriodic(i)) {
100  gdomainx.grow(i, nghost[i]);
101  }
102  }
103 
104  //
105  // We fill all of the interior and periodic ghost cells first, so we can fill
106  // those directly inside the lateral and vertical calls.
107  //
108  if (do_fb) {
109  mf.FillBoundary(m_geom.periodicity());
110  }
111 
112 #ifdef AMREX_USE_OMP
113 #pragma omp parallel if (Gpu::notInLaunchRegion())
114 #endif
115  {
116  for (MFIter mfi(mf,false); mfi.isValid(); ++mfi)
117  {
118  //
119  // This is the box we pass to the different routines
120  // NOTE -- this is the full grid NOT the tile box
121  //
122  Box bx = mfi.validbox();
123 
124  //
125  // These are the boxes we use to test on relative to the domain
126  //
127  Box xbx1 = surroundingNodes(bx,0); xbx1.grow(nghost);
128  if(xbx1.smallEnd(2) < domain.smallEnd(2)) xbx1.setSmall(2,domain.smallEnd(2));
129  if(xbx1.bigEnd(2) > domain.bigEnd(2)) xbx1.setBig(2,domain.bigEnd(2));
130 
131  Box xbx2 = surroundingNodes(bx,0); xbx2.grow(nghost);
132 
133  Array4<const Real> z_nd_arr;
134 
135  if (m_z_phys_nd)
136  {
137  z_nd_arr = m_z_phys_nd->const_array(mfi);
138  }
139 
140  if (!gdomainx.contains(xbx2))
141  {
142  const Array4<Real> velx_arr = mf.array(mfi);
143 
144  if (!m_use_real_bcs)
145  {
146  if (!gdomainx.contains(xbx1))
147  {
148  impose_lateral_xvel_bcs(velx_arr,xbx1,domain,bccomp);
149  }
150  }
151 
152  impose_vertical_xvel_bcs(velx_arr,xbx2,domain,z_nd_arr,dxInv,bccomp,time);
153  }
154  } // MFIter
155  } // OpenMP
156 } // operator()

Member Data Documentation

◆ m_bc_extdir_vals

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

◆ m_bc_neumann_vals

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

◆ m_domain_bcs_type

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

◆ m_domain_bcs_type_d

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

◆ m_geom

amrex::Geometry ERFPhysBCFunct_u::m_geom
private

◆ m_lev

int ERFPhysBCFunct_u::m_lev
private

◆ m_u_bc_data

amrex::Real* ERFPhysBCFunct_u::m_u_bc_data
private

◆ m_use_real_bcs

bool ERFPhysBCFunct_u::m_use_real_bcs
private

◆ m_z_phys_nd

amrex::MultiFab* ERFPhysBCFunct_u::m_z_phys_nd
private

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