ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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, amrex::MultiFab &xvel, amrex::MultiFab &yvel, 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::Array4< amrex::Real const > &xvel_arr, const amrex::Array4< amrex::Real const > &yvel_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
152  : m_lev(lev),
153  m_geom(geom), m_domain_bcs_type(domain_bcs_type),
154  m_domain_bcs_type_d(domain_bcs_type_d),
155  m_bc_extdir_vals(bc_extdir_vals),
156  m_bc_neumann_vals(bc_neumann_vals),
157  m_z_phys_nd(z_phys_nd.get()),
158  m_use_real_bcs(use_real_bcs),
159  m_v_bc_data(v_bc_data)
160  { }
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:190
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:192
amrex::Real * m_v_bc_data
Definition: ERF_PhysBCFunct.H:197
int m_lev
Definition: ERF_PhysBCFunct.H:189
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:194
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:196
amrex::Vector< amrex::BCRec > m_domain_bcs_type
Definition: ERF_PhysBCFunct.H:191
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:193
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:195

◆ ~ERFPhysBCFunct_v()

ERFPhysBCFunct_v::~ERFPhysBCFunct_v ( )
inline
162 {}

Member Function Documentation

◆ impose_lateral_yvel_bcs()

void ERFPhysBCFunct_v::impose_lateral_yvel_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,
int  bccomp 
)
19 {
20  BL_PROFILE_VAR("impose_lateral_yvel_bcs()",impose_lateral_yvel_bcs);
21  const auto& dom_lo = lbound(domain);
22  const auto& dom_hi = ubound(domain);
23 
24  // Based on BCRec for the domain, we need to make BCRec for this Box
25  // bccomp is used as starting index for m_domain_bcs_type
26  // 0 is used as starting index for bcrs
27  Vector<BCRec> bcrs(1);
28  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs);
29 
30  // xlo: ori = 0
31  // ylo: ori = 1
32  // zlo: ori = 2
33  // xhi: ori = 3
34  // yhi: ori = 4
35  // zhi: ori = 5
36 
37  Gpu::DeviceVector<BCRec> bcrs_d(1);
38  Gpu::copyAsync(Gpu::hostToDevice, bcrs.begin(), bcrs.end(), bcrs_d.begin());
39  const BCRec* bc_ptr = bcrs_d.data();
40 
41  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>, 1> l_bc_extdir_vals_d;
42 
43  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
44  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp][ori];
45  }
46 
47  GeometryData const& geomdata = m_geom.data();
48  bool is_periodic_in_x = geomdata.isPeriodic(0);
49  bool is_periodic_in_y = geomdata.isPeriodic(1);
50 
51  // First do all ext_dir bcs
52  if (!is_periodic_in_x)
53  {
54  // Populate ghost cells on lo-x and hi-x domain boundaries
55  Real* yvel_bc_ptr = m_v_bc_data;
56  Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
57  Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
58  ParallelFor(bx_xlo, bx_xhi,
59  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
60  int iflip = dom_lo.x - 1- i;
61  if (bc_ptr[0].lo(0) == ERFBCType::ext_dir) {
62  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
63  } else if (bc_ptr[0].lo(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_lo.x,i,j) >= 0.) {
64  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
65  } else if (bc_ptr[0].lo(0) == ERFBCType::foextrap) {
66  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
67  } else if (bc_ptr[0].lo(0) == ERFBCType::open) {
68  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
69  } else if (bc_ptr[0].lo(0) == ERFBCType::reflect_even) {
70  dest_arr(i,j,k) = dest_arr(iflip,j,k);
71  } else if (bc_ptr[0].lo(0) == ERFBCType::reflect_odd) {
72  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
73  } else if (bc_ptr[0].lo(0) == ERFBCType::hoextrap) {
74  Real delta_i = (dom_lo.x - i);
75  dest_arr(i,j,k) = (1.0 + delta_i)*dest_arr(dom_lo.x,j,k) - delta_i*dest_arr(dom_lo.x+1,j,k);
76  }
77  },
78  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
79  int iflip = 2*dom_hi.x + 1 - i;
80  if (bc_ptr[0].hi(0) == ERFBCType::ext_dir) {
81  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
82  } else if (bc_ptr[0].hi(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_hi.x+1,i,j) <= 0.) {
83  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
84  } else if (bc_ptr[0].hi(0) == ERFBCType::foextrap) {
85  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
86  } else if (bc_ptr[0].hi(0) == ERFBCType::open) {
87  dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
88  } else if (bc_ptr[0].hi(0) == ERFBCType::reflect_even) {
89  dest_arr(i,j,k) = dest_arr(iflip,j,k);
90  } else if (bc_ptr[0].hi(0) == ERFBCType::reflect_odd) {
91  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
92  } else if (bc_ptr[0].hi(0) == ERFBCType::hoextrap) {
93  Real delta_i = (i - dom_hi.x);
94  dest_arr(i,j,k) = (1.0 + delta_i)*dest_arr(dom_hi.x,j,k) - delta_i*dest_arr(dom_hi.x-1,j,k);
95  }
96  }
97  );
98  }
99 
100  if (!is_periodic_in_y)
101  {
102  // Populate ghost cells on lo-y and hi-y domain boundaries
103  Real* yvel_bc_ptr = m_v_bc_data;
104  Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
105  Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+2);
106  Box bx_ylo_face(bx); bx_ylo_face.setSmall(1,dom_lo.y ); bx_ylo_face.setBig(1,dom_lo.y );
107  Box bx_yhi_face(bx); bx_yhi_face.setSmall(1,dom_hi.y+1); bx_yhi_face.setBig(1,dom_hi.y+1);
108  ParallelFor(bx_ylo, bx_ylo_face,
109  [=] AMREX_GPU_DEVICE (int i, int j, int k)
110  {
111  int jflip = dom_lo.y-j;
112  if (bc_ptr[0].lo(1) == ERFBCType::ext_dir) {
113  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
114  } else if (bc_ptr[0].lo(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_lo.y,k) >= 0.) {
115  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
116  } else if (bc_ptr[0].lo(1) == ERFBCType::foextrap) {
117  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
118  } else if (bc_ptr[0].lo(1) == ERFBCType::open) {
119  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
120  } else if (bc_ptr[0].lo(1) == ERFBCType::reflect_even) {
121  dest_arr(i,j,k) = dest_arr(i,jflip,k);
122  } else if (bc_ptr[0].lo(1) == ERFBCType::reflect_odd) {
123  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
124  } else if (bc_ptr[0].lo(1) == ERFBCType::neumann_int) {
125  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;
126  } else if (bc_ptr[0].lo(1) == ERFBCType::hoextrap) {
127  Real delta_j = (dom_lo.y - j);
128  dest_arr(i,j,k) = (1.0 + delta_j)*dest_arr(i,dom_lo.y,k) - delta_j*dest_arr(i,dom_lo.y+1,k);
129  }
130  },
131  // We only set the values on the domain faces themselves if EXT_DIR or neumann_int
132  [=] AMREX_GPU_DEVICE (int i, int j, int k)
133  {
134  if (bc_ptr[0].lo(1) == ERFBCType::ext_dir) {
135  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
136  } else if (bc_ptr[0].lo(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_lo.y,k) >= 0.) {
137  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
138  } else if (bc_ptr[0].lo(1) == ERFBCType::neumann_int) {
139  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;
140  }
141  }
142  );
143  ParallelFor(bx_yhi, bx_yhi_face,
144  [=] AMREX_GPU_DEVICE (int i, int j, int k)
145  {
146  int jflip = 2*(dom_hi.y + 1) - j;
147  if (bc_ptr[0].hi(1) == ERFBCType::ext_dir) {
148  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
149  } else if (bc_ptr[0].hi(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_hi.y+1,k) <= 0.) {
150  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
151  } else if (bc_ptr[0].hi(1) == ERFBCType::foextrap) {
152  dest_arr(i,j,k) = dest_arr(i,dom_hi.y+1,k);
153  } else if (bc_ptr[0].hi(1) == ERFBCType::open) {
154  dest_arr(i,j,k) = dest_arr(i,dom_hi.y+1,k);
155  } else if (bc_ptr[0].hi(1) == ERFBCType::reflect_even) {
156  dest_arr(i,j,k) = dest_arr(i,jflip,k);
157  } else if (bc_ptr[0].hi(1) == ERFBCType::reflect_odd) {
158  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
159  } else if (bc_ptr[0].hi(1) == ERFBCType::neumann_int) {
160  dest_arr(i,j,k) = (4.0*dest_arr(i,dom_hi.y,k) - dest_arr(i,dom_hi.y-1,k))/3.0;
161  } else if (bc_ptr[0].hi(1) == ERFBCType::hoextrap) {
162  Real delta_j = (j - dom_hi.y - 1);
163  dest_arr(i,j,k) = (1.0 + delta_j)*dest_arr(i,dom_hi.y+1,k) - delta_j*dest_arr(i,dom_hi.y,k);
164  }
165  },
166  // We only set the values on the domain faces themselves if EXT_DIR or neumann_int
167  [=] AMREX_GPU_DEVICE (int i, int j, int k)
168  {
169  if (bc_ptr[0].hi(1) == ERFBCType::ext_dir) {
170  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
171  } else if (bc_ptr[0].hi(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_hi.y+1,k) <= 0.) {
172  dest_arr(i,j,k) = (yvel_bc_ptr) ? yvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
173  } else if (bc_ptr[0].hi(1) == ERFBCType::neumann_int) {
174  dest_arr(i,j,k) = (4.0*dest_arr(i,dom_hi.y,k) - dest_arr(i,dom_hi.y-1,k))/3.0;
175  }
176  }
177  );
178  }
179  Gpu::streamSynchronize();
180 }
void impose_lateral_yvel_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, int bccomp)
Definition: ERF_BoundaryConditionsYvel.cpp:14
@ open
Definition: ERF_IndexDefines.H:215
@ reflect_odd
Definition: ERF_IndexDefines.H:205
@ hoextrap
Definition: ERF_IndexDefines.H:216
@ foextrap
Definition: ERF_IndexDefines.H:208
@ ext_dir
Definition: ERF_IndexDefines.H:209
@ ext_dir_upwind
Definition: ERF_IndexDefines.H:217
@ neumann_int
Definition: ERF_IndexDefines.H:214
@ reflect_even
Definition: ERF_IndexDefines.H:207

◆ 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 
)
198 {
199  BL_PROFILE_VAR("impose_vertical_yvel_bcs()",impose_vertical_yvel_bcs);
200  const auto& dom_lo = lbound(domain);
201  const auto& dom_hi = ubound(domain);
202 
203  Box per_grown_domain(domain);
204  int growx = (m_geom.isPeriodic(0)) ? 1 : 0;
205  int growy = (m_geom.isPeriodic(1)) ? 1 : 0;
206  per_grown_domain.grow(IntVect(growx,growy,0));
207  const auto& perdom_lo = lbound(per_grown_domain);
208  const auto& perdom_hi = ubound(per_grown_domain);
209 
210  // Based on BCRec for the domain, we need to make BCRec for this Box
211  // bccomp is used as starting index for m_domain_bcs_type
212  // 0 is used as starting index for bcrs
213  Vector<BCRec> bcrs(1);
214  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs);
215 
216  // xlo: ori = 0
217  // ylo: ori = 1
218  // zlo: ori = 2
219  // xhi: ori = 3
220  // yhi: ori = 4
221  // zhi: ori = 5
222 
223  Gpu::DeviceVector<BCRec> bcrs_d(1);
224  Gpu::copyAsync(Gpu::hostToDevice, bcrs.begin(), bcrs.end(), bcrs_d.begin());
225  const BCRec* bc_ptr = bcrs_d.data();
226 
227  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>, 1> l_bc_extdir_vals_d;
228 
229  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
230  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp][ori];
231  }
232 
233  GeometryData const& geomdata = m_geom.data();
234 
235  {
236  // Populate ghost cells on lo-z and hi-z domain boundaries
237  Box bx_zlo(bx); bx_zlo.setBig (2,dom_lo.z-1);
238  Box bx_zhi(bx); bx_zhi.setSmall(2,dom_hi.z+1);
239  ParallelFor(bx_zlo, bx_zhi,
240  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
241  int kflip = dom_lo.z - 1 - k;
242  if (bc_ptr[0].lo(2) == ERFBCType::ext_dir) {
243  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][2];
244  } else if (bc_ptr[0].lo(2) == ERFBCType::foextrap) {
245  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
246  } else if (bc_ptr[0].lo(2) == ERFBCType::open) {
247  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
248  } else if (bc_ptr[0].lo(2) == ERFBCType::reflect_even) {
249  dest_arr(i,j,k) = dest_arr(i,j,kflip);
250  } else if (bc_ptr[0].lo(2) == ERFBCType::reflect_odd) {
251  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
252  } else if (bc_ptr[0].lo(2) == ERFBCType::hoextrap) {
253  Real delta_k = (dom_lo.z - k);
254  dest_arr(i,j,k) = (1.0 + delta_k)*dest_arr(i,j,dom_lo.z) - delta_k*dest_arr(i,j,dom_lo.z+1);
255  }
256  },
257  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
258  int kflip = 2*dom_hi.z + 1 - k;
259  if (bc_ptr[0].hi(2) == ERFBCType::ext_dir) {
260  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][5];
261  } else if (bc_ptr[0].hi(2) == ERFBCType::foextrap) {
262  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
263  } else if (bc_ptr[0].hi(2) == ERFBCType::open) {
264  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
265  } else if (bc_ptr[0].hi(2) == ERFBCType::reflect_even) {
266  dest_arr(i,j,k) = dest_arr(i,j,kflip);
267  } else if (bc_ptr[0].hi(2) == ERFBCType::reflect_odd) {
268  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
269  } else if (bc_ptr[0].hi(2) == ERFBCType::hoextrap){
270  Real delta_k = (k - dom_hi.z);
271  dest_arr(i,j,k) = (1.0 + delta_k)*dest_arr(i,j,dom_hi.z) - delta_k*dest_arr(i,j,dom_hi.z-1);
272  }
273  }
274  );
275  }
276 
277  if (m_z_phys_nd) {
278 
279  const auto& bx_lo = lbound(bx);
280  const auto& bx_hi = ubound(bx);
281 
282  const auto& zphys_lo = lbound(Box(z_phys_nd));
283  const auto& zphys_hi = ubound(Box(z_phys_nd));
284 
285  // Neumann conditions (d<var>/dn = 0) must be aware of the surface normal with terrain.
286  // An additional source term arises from d<var>/dx & d<var>/dy & met_h_xi/eta/zeta.
287  //=====================================================================================
288  // Only modify scalars, U, or V
289  // Loop over each component
290  // Hit for Neumann condition at kmin
291  if(bcrs[0].lo(2) == ERFBCType::foextrap) {
292  // Loop over ghost cells in bottom XY-plane (valid box)
293  Box xybx = bx;
294  xybx.setBig(2,-1);
295  xybx.setSmall(2,bx.smallEnd()[2]);
296  int k0 = 0;
297 
298  // Get the dz cell size
299  Real dz = geomdata.CellSize(2);
300 
301  // Fill all the Neumann srcs with terrain
302  ParallelFor(xybx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
303  {
304  // Clip indices for ghost-cells
305  int ii = amrex::min(amrex::max(i,perdom_lo.x),perdom_hi.x);
306  ii = amrex::min(amrex::max(ii,zphys_lo.x),zphys_hi.x);
307  int jj = amrex::min(amrex::max(j,perdom_lo.y),perdom_hi.y);
308  jj = amrex::min(amrex::max(jj,zphys_lo.y),zphys_hi.y);
309 
310  // Get metrics
311  Real met_h_xi = Compute_h_xi_AtJface (ii, jj, k0, dxInv, z_phys_nd);
312  Real met_h_eta = Compute_h_eta_AtJface (ii, jj, k0, dxInv, z_phys_nd);
313  Real met_h_zeta = Compute_h_zeta_AtJface(ii, jj, k0, dxInv, z_phys_nd);
314 
315  // GradX at IJK location inside domain -- this relies on the assumption that we have
316  // used foextrap for cell-centered quantities outside the domain to define the gradient as zero
317  Real GradVarx, GradVary;
318  if ( i < dom_lo.x-1 || i > dom_hi.x+1 || (i+1 > bx_hi.x && i-1 < bx_lo.x) ) {
319  GradVarx = 0.0;
320  } else if (i+1 > bx_hi.x) {
321  GradVarx = dxInv[0] * (dest_arr(i ,j,k0) - dest_arr(i-1,j,k0));
322  } else if (i-1 < bx_lo.x) {
323  GradVarx = dxInv[0] * (dest_arr(i+1,j,k0) - dest_arr(i ,j,k0));
324  } else {
325  GradVarx = 0.5 * dxInv[0] * (dest_arr(i+1,j,k0) - dest_arr(i-1,j,k0));
326  }
327 
328  // GradY at IJK location inside domain -- this relies on the assumption that we have
329  // used foextrap for cell-centered quantities outside the domain to define the gradient as zero
330  if ( j < dom_lo.y-1 || j > dom_hi.y+1 || (j+1 > bx_hi.y && j-1 < bx_lo.y) ) {
331  GradVary = 0.0;
332  } else if (j+1 > bx_hi.y) {
333  GradVary = dxInv[1] * (dest_arr(i,j ,k0) - dest_arr(i,j-1,k0));
334  } else if (j-1 < bx_lo.y) {
335  GradVary = dxInv[1] * (dest_arr(i,j+1,k0) - dest_arr(i,j ,k0));
336  } else {
337  GradVary = 0.5 * dxInv[1] * (dest_arr(i,j+1,k0) - dest_arr(i,j-1,k0));
338  }
339 
340  // Prefactor
341  Real met_fac = met_h_zeta / ( met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1. );
342 
343  // Accumulate in bottom ghost cell (EXTRAP already populated)
344  dest_arr(i,j,k) -= dz * met_fac * ( met_h_xi * GradVarx + met_h_eta * GradVary );
345  });
346  } // foextrap
347  } //m_z_phys_nd
348  Gpu::streamSynchronize();
349 }
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:153
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:139
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:167
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:193
Here is the call graph for this function:

◆ operator()()

void ERFPhysBCFunct_v::operator() ( amrex::MultiFab &  mf,
amrex::MultiFab &  xvel,
amrex::MultiFab &  yvel,
amrex::IntVect const &  nghost,
const amrex::Real  time,
int  bccomp,
bool  do_fb 
)
166 {
167  BL_PROFILE("ERFPhysBCFunct_v::()");
168 
169  if (m_geom.isAllPeriodic()) return;
170 
171  const auto& domain = m_geom.Domain();
172  const auto dxInv = m_geom.InvCellSizeArray();
173 
174  Box gdomainy = surroundingNodes(domain,1);
175  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
176  if (m_geom.isPeriodic(i)) {
177  gdomainy.grow(i, nghost[i]);
178  }
179  }
180 
181  //
182  // We fill all of the interior and periodic ghost cells first, so we can fill
183  // those directly inside the lateral and vertical calls.
184  //
185  if (do_fb) {
186  mf.FillBoundary(m_geom.periodicity());
187  }
188 
189 #ifdef AMREX_USE_OMP
190 #pragma omp parallel if (Gpu::notInLaunchRegion())
191 #endif
192  {
193  for (MFIter mfi(mf,false); mfi.isValid(); ++mfi)
194  {
195  //
196  // This is the box we pass to the different routines
197  // NOTE -- this is the full grid NOT the tile box
198  //
199  Box bx = mfi.validbox();
200 
201  //
202  // These are the boxes we use to test on relative to the domain
203  //
204  Box ybx1 = surroundingNodes(bx,1); ybx1.grow(nghost);
205  if (ybx1.smallEnd(2) < domain.smallEnd(2)) ybx1.setSmall(2,domain.smallEnd(2));
206  if (ybx1.bigEnd(2) > domain.bigEnd(2)) ybx1.setBig(2,domain.bigEnd(2));
207 
208  Box ybx2 = surroundingNodes(bx,1); ybx2.grow(nghost);
209 
210  Array4<const Real> z_nd_arr;
211 
212  if (m_z_phys_nd)
213  {
214  z_nd_arr = m_z_phys_nd->const_array(mfi);
215  }
216 
217  if (!gdomainy.contains(ybx2))
218  {
219  Array4< Real> const& dest_arr = mf.array(mfi);
220  Array4<const Real> const& velx_arr = xvel.const_array(mfi);
221  Array4<const Real> const& vely_arr = yvel.const_array(mfi);
222 
223  if (!m_use_real_bcs)
224  {
225  impose_lateral_yvel_bcs(dest_arr,velx_arr,vely_arr,ybx1,domain,bccomp);
226  }
227 
228  impose_vertical_yvel_bcs(dest_arr,ybx2,domain,z_nd_arr,dxInv,bccomp);
229  }
230 
231  } // MFIter
232  } // OpenMP
233 } // operator()
@ 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_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: