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+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 *v_bc_data)
 
 ~ERFPhysBCFunct_v ()
 
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_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+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_v_bc_data
 

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+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 *  v_bc_data 
)
inline
146  : m_lev(lev),
147  m_geom(geom), m_domain_bcs_type(domain_bcs_type),
148  m_domain_bcs_type_d(domain_bcs_type_d),
149  m_bc_extdir_vals(bc_extdir_vals),
150  m_bc_neumann_vals(bc_neumann_vals),
151  m_z_phys_nd(z_phys_nd.get()),
152  m_use_real_bcs(use_real_bcs),
153  m_v_bc_data(v_bc_data)
154  { }
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:183
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:185
amrex::Real * m_v_bc_data
Definition: ERF_PhysBCFunct.H:190
int m_lev
Definition: ERF_PhysBCFunct.H:182
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:187
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:189
amrex::Vector< amrex::BCRec > m_domain_bcs_type
Definition: ERF_PhysBCFunct.H:184
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:186
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:188

◆ ~ERFPhysBCFunct_v()

ERFPhysBCFunct_v::~ERFPhysBCFunct_v ( )
inline
156 {}

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  Vector<BCRec> bcrs(1);
26  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs);
27 
28  // xlo: ori = 0
29  // ylo: ori = 1
30  // zlo: ori = 2
31  // xhi: ori = 3
32  // yhi: ori = 4
33  // zhi: ori = 5
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  // Populate ghost cells on lo-x and hi-x domain boundaries
53  Real* yvel_bc_ptr = m_v_bc_data;
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(bx_xlo, bx_xhi,
57  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
58  int iflip = dom_lo.x - 1- i;
59  if (bc_ptr[0].lo(0) == ERFBCType::ext_dir) {
60  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
61  } else if (bc_ptr[0].lo(0) == ERFBCType::foextrap) {
62  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
63  } else if (bc_ptr[0].lo(0) == ERFBCType::open) {
64  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
65  } else if (bc_ptr[0].lo(0) == ERFBCType::reflect_even) {
66  dest_arr(i,j,k) = dest_arr(iflip,j,k);
67  } else if (bc_ptr[0].lo(0) == ERFBCType::reflect_odd) {
68  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
69  }
70  },
71  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
72  int iflip = 2*dom_hi.x + 1 - i;
73  if (bc_ptr[0].hi(0) == ERFBCType::ext_dir) {
74  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
75  } else if (bc_ptr[0].hi(0) == ERFBCType::foextrap) {
76  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
77  } else if (bc_ptr[0].hi(0) == ERFBCType::open) {
78  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
79  } else if (bc_ptr[0].hi(0) == ERFBCType::reflect_even) {
80  dest_arr(i,j,k) = dest_arr(iflip,j,k);
81  } else if (bc_ptr[0].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  Real* yvel_bc_ptr = m_v_bc_data;
92  Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
93  Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+2);
94  Box bx_ylo_face(bx); bx_ylo_face.setSmall(1,dom_lo.y ); bx_ylo_face.setBig(1,dom_lo.y );
95  Box bx_yhi_face(bx); bx_yhi_face.setSmall(1,dom_hi.y+1); bx_yhi_face.setBig(1,dom_hi.y+1);
96  ParallelFor(bx_ylo, bx_ylo_face,
97  [=] AMREX_GPU_DEVICE (int i, int j, int k)
98  {
99  int jflip = dom_lo.y-j;
100  if (bc_ptr[0].lo(1) == ERFBCType::ext_dir) {
101  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
102  } else if (bc_ptr[0].lo(1) == ERFBCType::foextrap) {
103  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
104  } else if (bc_ptr[0].lo(1) == ERFBCType::open) {
105  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
106  } else if (bc_ptr[0].lo(1) == ERFBCType::reflect_even) {
107  dest_arr(i,j,k) = dest_arr(i,jflip,k);
108  } else if (bc_ptr[0].lo(1) == ERFBCType::reflect_odd) {
109  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
110  } else if (bc_ptr[0].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  [=] AMREX_GPU_DEVICE (int i, int j, int k)
116  {
117  if (bc_ptr[0].lo(1) == ERFBCType::ext_dir) {
118  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
119  } else if (bc_ptr[0].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(bx_yhi, bx_yhi_face,
125  [=] AMREX_GPU_DEVICE (int i, int j, int k)
126  {
127  int jflip = 2*(dom_hi.y + 1) - j;
128  if (bc_ptr[0].hi(1) == ERFBCType::ext_dir) {
129  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
130  } else if (bc_ptr[0].hi(1) == ERFBCType::foextrap) {
131  dest_arr(i,j,k) = dest_arr(i,dom_hi.y+1,k);
132  } else if (bc_ptr[0].hi(1) == ERFBCType::open) {
133  dest_arr(i,j,k) = dest_arr(i,dom_hi.y+1,k);
134  } else if (bc_ptr[0].hi(1) == ERFBCType::reflect_even) {
135  dest_arr(i,j,k) = dest_arr(i,jflip,k);
136  } else if (bc_ptr[0].hi(1) == ERFBCType::reflect_odd) {
137  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
138  } else if (bc_ptr[0].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  [=] AMREX_GPU_DEVICE (int i, int j, int k)
144  {
145  if (bc_ptr[0].hi(1) == ERFBCType::ext_dir) {
146  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
147  } else if (bc_ptr[0].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 }
void impose_lateral_yvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, int bccomp)
Definition: ERF_BoundaryConditionsYvel.cpp:14
@ 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_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  Box per_grown_domain(domain);
178  int growx = (m_geom.isPeriodic(0)) ? 1 : 0;
179  int growy = (m_geom.isPeriodic(1)) ? 1 : 0;
180  per_grown_domain.grow(IntVect(growx,growy,0));
181  const auto& perdom_lo = lbound(per_grown_domain);
182  const auto& perdom_hi = ubound(per_grown_domain);
183 
184  // Based on BCRec for the domain, we need to make BCRec for this Box
185  // bccomp is used as starting index for m_domain_bcs_type
186  // 0 is used as starting index for bcrs
187  Vector<BCRec> bcrs(1);
188  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs);
189 
190  // xlo: ori = 0
191  // ylo: ori = 1
192  // zlo: ori = 2
193  // xhi: ori = 3
194  // yhi: ori = 4
195  // zhi: ori = 5
196 
197  Gpu::DeviceVector<BCRec> bcrs_d(1);
198  Gpu::copyAsync(Gpu::hostToDevice, bcrs.begin(), bcrs.end(), bcrs_d.begin());
199  const BCRec* bc_ptr = bcrs_d.data();
200 
201  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>, 1> l_bc_extdir_vals_d;
202 
203  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
204  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp][ori];
205  }
206 
207  GeometryData const& geomdata = m_geom.data();
208 
209  {
210  // Populate ghost cells on lo-z and hi-z domain boundaries
211  Box bx_zlo(bx); bx_zlo.setBig (2,dom_lo.z-1);
212  Box bx_zhi(bx); bx_zhi.setSmall(2,dom_hi.z+1);
213  ParallelFor(bx_zlo, bx_zhi,
214  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
215  int kflip = dom_lo.z - 1 - k;
216  if (bc_ptr[0].lo(2) == ERFBCType::ext_dir) {
217  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][2];
218  } else if (bc_ptr[0].lo(2) == ERFBCType::foextrap) {
219  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
220  } else if (bc_ptr[0].lo(2) == ERFBCType::open) {
221  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
222  } else if (bc_ptr[0].lo(2) == ERFBCType::reflect_even) {
223  dest_arr(i,j,k) = dest_arr(i,j,kflip);
224  } else if (bc_ptr[0].lo(2) == ERFBCType::reflect_odd) {
225  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
226  }
227  },
228  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
229  int kflip = 2*dom_hi.z + 1 - k;
230  if (bc_ptr[0].hi(2) == ERFBCType::ext_dir) {
231  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][5];
232  } else if (bc_ptr[0].hi(2) == ERFBCType::foextrap) {
233  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
234  } else if (bc_ptr[0].hi(2) == ERFBCType::open) {
235  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
236  } else if (bc_ptr[0].hi(2) == ERFBCType::reflect_even) {
237  dest_arr(i,j,k) = dest_arr(i,j,kflip);
238  } else if (bc_ptr[0].hi(2) == ERFBCType::reflect_odd) {
239  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
240  }
241  }
242  );
243  }
244 
245  if (m_z_phys_nd) {
246 
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_AtJface (ii, jj, k0, dxInv, z_phys_nd);
280  Real met_h_eta = Compute_h_eta_AtJface (ii, jj, k0, dxInv, z_phys_nd);
281  Real met_h_zeta = Compute_h_zeta_AtJface(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 cell-centered 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 cell-centered 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_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: ERF_TerrainMetrics.H:144
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: ERF_TerrainMetrics.H:130
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: ERF_TerrainMetrics.H:158
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: ERF_BoundaryConditionsYvel.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,
bool  do_fb 
)
161 {
162  BL_PROFILE("ERFPhysBCFunct_v::()");
163 
164  if (m_geom.isAllPeriodic()) return;
165 
166  const auto& domain = m_geom.Domain();
167  const auto dxInv = m_geom.InvCellSizeArray();
168 
169  Box gdomainy = surroundingNodes(domain,1);
170  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
171  if (m_geom.isPeriodic(i)) {
172  gdomainy.grow(i, nghost[i]);
173  }
174  }
175 
176  //
177  // We fill all of the interior and periodic ghost cells first, so we can fill
178  // those directly inside the lateral and vertical calls.
179  //
180  if (do_fb) {
181  mf.FillBoundary(m_geom.periodicity());
182  }
183 
184 #ifdef AMREX_USE_OMP
185 #pragma omp parallel if (Gpu::notInLaunchRegion())
186 #endif
187  {
188  for (MFIter mfi(mf,false); mfi.isValid(); ++mfi)
189  {
190  //
191  // This is the box we pass to the different routines
192  // NOTE -- this is the full grid NOT the tile box
193  //
194  Box bx = mfi.validbox();
195 
196  //
197  // These are the boxes we use to test on relative to the domain
198  //
199  Box ybx1 = surroundingNodes(bx,1); ybx1.grow(nghost);
200  if (ybx1.smallEnd(2) < domain.smallEnd(2)) ybx1.setSmall(2,domain.smallEnd(2));
201  if (ybx1.bigEnd(2) > domain.bigEnd(2)) ybx1.setBig(2,domain.bigEnd(2));
202 
203  Box ybx2 = surroundingNodes(bx,1); ybx2.grow(nghost);
204 
205  Array4<const Real> z_nd_arr;
206 
207  if (m_z_phys_nd)
208  {
209  z_nd_arr = m_z_phys_nd->const_array(mfi);
210  }
211 
212  if (!gdomainy.contains(ybx2))
213  {
214  const Array4<Real> vely_arr = mf.array(mfi);
215 
216  if (!m_use_real_bcs)
217  {
218  impose_lateral_yvel_bcs(vely_arr,ybx1,domain,bccomp);
219  }
220 
221  impose_vertical_yvel_bcs(vely_arr,ybx2,domain,z_nd_arr,dxInv,bccomp);
222  }
223 
224  } // MFIter
225  } // OpenMP
226 } // operator()

Member Data Documentation

◆ m_bc_extdir_vals

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

◆ m_bc_neumann_vals

amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NBCVAR_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_v_bc_data

amrex::Real* ERFPhysBCFunct_v::m_v_bc_data
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: