ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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, amrex::MultiFab &xvel, amrex::MultiFab &yvel, 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::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_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
94  : m_lev(lev), m_geom(geom),
95  m_domain_bcs_type(domain_bcs_type),
96  m_domain_bcs_type_d(domain_bcs_type_d),
97  m_bc_extdir_vals(bc_extdir_vals),
98  m_bc_neumann_vals(bc_neumann_vals),
99  m_z_phys_nd(z_phys_nd.get()),
100  m_use_real_bcs(use_real_bcs),
101  m_u_bc_data(u_bc_data)
102  { }
amrex::Geometry m_geom
Definition: ERF_PhysBCFunct.H:131
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF_PhysBCFunct.H:134
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
Definition: ERF_PhysBCFunct.H:133
bool m_use_real_bcs
Definition: ERF_PhysBCFunct.H:137
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF_PhysBCFunct.H:135
amrex::Real * m_u_bc_data
Definition: ERF_PhysBCFunct.H:138
int m_lev
Definition: ERF_PhysBCFunct.H:130
amrex::Vector< amrex::BCRec > m_domain_bcs_type
Definition: ERF_PhysBCFunct.H:132
amrex::MultiFab * m_z_phys_nd
Definition: ERF_PhysBCFunct.H:136

◆ ~ERFPhysBCFunct_u()

ERFPhysBCFunct_u::~ERFPhysBCFunct_u ( )
inline
104 {}

Member Function Documentation

◆ impose_lateral_xvel_bcs()

void ERFPhysBCFunct_u::impose_lateral_xvel_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_xvel_bcs()",impose_lateral_xvel_bcs);
21  const auto& dom_lo = lbound(domain);
22  const auto& dom_hi = ubound(domain);
23 
24  // xlo: ori = 0
25  // ylo: ori = 1
26  // zlo: ori = 2
27  // xhi: ori = 3
28  // yhi: ori = 4
29  // zhi: ori = 5
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(1);
35  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs);
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  Real* xvel_bc_ptr = m_u_bc_data;
55  Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
56  Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+2);
57  Box bx_xlo_face(bx); bx_xlo_face.setSmall(0,dom_lo.x ); bx_xlo_face.setBig(0,dom_lo.x );
58  Box bx_xhi_face(bx); bx_xhi_face.setSmall(0,dom_hi.x+1); bx_xhi_face.setBig(0,dom_hi.x+1);
59  ParallelFor(bx_xlo, bx_xlo_face,
60  [=] AMREX_GPU_DEVICE (int i, int j, int k)
61  {
62  int iflip = dom_lo.x - i;
63  if (bc_ptr[0].lo(0) == ERFBCType::ext_dir) {
64  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
65  } else if (bc_ptr[0].lo(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_lo.x,j,k) >= 0.) {
66  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
67  } else if (bc_ptr[0].lo(0) == ERFBCType::foextrap) {
68  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
69  } else if (bc_ptr[0].lo(0) == ERFBCType::open) {
70  dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
71  } else if (bc_ptr[0].lo(0) == ERFBCType::reflect_even) {
72  dest_arr(i,j,k) = dest_arr(iflip,j,k);
73  } else if (bc_ptr[0].lo(0) == ERFBCType::reflect_odd) {
74  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
75  } else if (bc_ptr[0].lo(0) == ERFBCType::neumann_int) {
76  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;
77  } else if (bc_ptr[0].lo(0) == ERFBCType::hoextrap) {
78  Real delta_i = (dom_lo.x - i);
79  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);
80  }
81  },
82  // We only set the values on the domain faces themselves if EXT_DIR or neumann_int
83  [=] AMREX_GPU_DEVICE (int i, int j, int k)
84  {
85  if (bc_ptr[0].lo(0) == ERFBCType::ext_dir) {
86  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
87  } else if (bc_ptr[0].lo(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_lo.x,j,k) >= 0.) {
88  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][0];
89  } else if (bc_ptr[0].lo(0) == ERFBCType::neumann_int) {
90  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;
91  }
92  }
93  );
94  ParallelFor(bx_xhi, bx_xhi_face,
95  [=] AMREX_GPU_DEVICE (int i, int j, int k)
96  {
97  int iflip = 2*(dom_hi.x + 1) - i;
98  if (bc_ptr[0].hi(0) == ERFBCType::ext_dir) {
99  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
100  } else if (bc_ptr[0].hi(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_hi.x+1,j,k) <= 0.) {
101  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
102  } else if (bc_ptr[0].hi(0) == ERFBCType::foextrap) {
103  dest_arr(i,j,k) = dest_arr(dom_hi.x+1,j,k);
104  } else if (bc_ptr[0].hi(0) == ERFBCType::open) {
105  dest_arr(i,j,k) = dest_arr(dom_hi.x+1,j,k);
106  } else if (bc_ptr[0].hi(0) == ERFBCType::reflect_even) {
107  dest_arr(i,j,k) = dest_arr(iflip,j,k);
108  } else if (bc_ptr[0].hi(0) == ERFBCType::reflect_odd) {
109  dest_arr(i,j,k) = -dest_arr(iflip,j,k);
110  } else if (bc_ptr[0].hi(0) == ERFBCType::neumann_int) {
111  dest_arr(i,j,k) = (4.0*dest_arr(dom_hi.x,j,k) - dest_arr(dom_hi.x-1,j,k))/3.0;
112  } else if (bc_ptr[0].hi(0) == ERFBCType::hoextrap) {
113  Real delta_i = (i - dom_hi.x - 1);
114  dest_arr(i,j,k) = (1.0 + delta_i)*dest_arr(dom_hi.x+1,j,k) - delta_i*dest_arr(dom_hi.x,j,k);
115  }
116  },
117  // We only set the values on the domain faces themselves if EXT_DIR or neumann_int
118  [=] AMREX_GPU_DEVICE (int i, int j, int k)
119  {
120  if (bc_ptr[0].hi(0) == ERFBCType::ext_dir) {
121  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
122  } else if (bc_ptr[0].hi(0) == ERFBCType::ext_dir_upwind && xvel_arr(dom_hi.x+1,j,k) <= 0.) {
123  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][3];
124  } else if (bc_ptr[0].hi(0) == ERFBCType::neumann_int) {
125  dest_arr(i,j,k) = (4.0*dest_arr(dom_hi.x,j,k) - dest_arr(dom_hi.x-1,j,k))/3.0;
126  }
127  }
128  );
129  } // not periodic in x
130 
131  if (!is_periodic_in_y)
132  {
133  // Populate ghost cells on lo-y and hi-y domain boundaries
134  Real* xvel_bc_ptr = m_u_bc_data;
135  Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
136  Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
137  ParallelFor(bx_ylo, bx_yhi,
138  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
139  int jflip = dom_lo.y - 1 - j;
140  if (bc_ptr[0].lo(1) == ERFBCType::ext_dir) {
141  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
142  } else if (bc_ptr[0].lo(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_lo.y,k) >= 0.) {
143  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][1];
144  } else if (bc_ptr[0].lo(1) == ERFBCType::foextrap) {
145  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
146  } else if (bc_ptr[0].lo(1) == ERFBCType::open) {
147  dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
148  } else if (bc_ptr[0].lo(1) == ERFBCType::reflect_even) {
149  dest_arr(i,j,k) = dest_arr(i,jflip,k);
150  } else if (bc_ptr[0].lo(1) == ERFBCType::reflect_odd) {
151  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
152  } else if (bc_ptr[0].lo(1) == ERFBCType::hoextrap) {
153  Real delta_j = (dom_lo.y - j);
154  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);
155  }
156  },
157  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
158  int jflip = 2*dom_hi.y + 1 - j;
159  if (bc_ptr[0].hi(1) == ERFBCType::ext_dir) {
160  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
161  } else if (bc_ptr[0].hi(1) == ERFBCType::ext_dir_upwind && yvel_arr(i,dom_hi.y+1,k) <= 0.) {
162  dest_arr(i,j,k) = (xvel_bc_ptr) ? xvel_bc_ptr[k] : l_bc_extdir_vals_d[0][4];
163  } else if (bc_ptr[0].hi(1) == ERFBCType::foextrap) {
164  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
165  } else if (bc_ptr[0].hi(1) == ERFBCType::open) {
166  dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
167  } else if (bc_ptr[0].hi(1) == ERFBCType::reflect_even) {
168  dest_arr(i,j,k) = dest_arr(i,jflip,k);
169  } else if (bc_ptr[0].hi(1) == ERFBCType::reflect_odd) {
170  dest_arr(i,j,k) = -dest_arr(i,jflip,k);
171  } else if (bc_ptr[0].hi(1) == ERFBCType::hoextrap) {
172  Real delta_j = (j - dom_hi.y);
173  dest_arr(i,j,k) = (1.0 + delta_j)*dest_arr(i,dom_hi.y,k) - delta_j*dest_arr(i,dom_hi.y-1,k);
174  }
175  }
176  );
177  } // not periodic in y
178 
179  Gpu::streamSynchronize();
180 }
void impose_lateral_xvel_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_BoundaryConditionsXvel.cpp:15
@ 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_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 
)
202 {
203  BL_PROFILE_VAR("impose_vertical_xvel_bcs()",impose_vertical_xvel_bcs);
204  const auto& dom_lo = lbound(domain);
205  const auto& dom_hi = ubound(domain);
206 
207  Box per_grown_domain(domain);
208  int growx = (m_geom.isPeriodic(0)) ? 1 : 0;
209  int growy = (m_geom.isPeriodic(1)) ? 1 : 0;
210  per_grown_domain.grow(IntVect(growx,growy,0));
211  const auto& perdom_lo = lbound(per_grown_domain);
212  const auto& perdom_hi = ubound(per_grown_domain);
213 
214  GeometryData const& geomdata = m_geom.data();
215 
216  // Based on BCRec for the domain, we need to make BCRec for this Box
217  // bccomp is used as starting index for m_domain_bcs_type
218  // 0 is used as starting index for bcrs
219  Vector<BCRec> bcrs(1);
220  setBC(enclosedCells(bx), domain, bccomp, 0, 1, m_domain_bcs_type, bcrs);
221 
222  Gpu::DeviceVector<BCRec> bcrs_d(1);
223  Gpu::copyAsync(Gpu::hostToDevice, bcrs.begin(), bcrs.end(), bcrs_d.begin());
224  const BCRec* bc_ptr = bcrs_d.data();
225 
226  GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,1> l_bc_extdir_vals_d;
227 
228  for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
229  l_bc_extdir_vals_d[0][ori] = m_bc_extdir_vals[bccomp][ori];
230  }
231 
232  {
233  // Populate ghost cells on lo-z and hi-z domain boundaries
234  Box bx_zlo(bx); bx_zlo.setBig (2,dom_lo.z-1);
235  Box bx_zhi(bx); bx_zhi.setSmall(2,dom_hi.z+1);
236  ParallelFor(bx_zlo, bx_zhi,
237  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
238  int kflip = dom_lo.z - 1 - k;
239  if (bc_ptr[0].lo(2) == ERFBCType::ext_dir) {
240  #ifdef ERF_USE_TERRAIN_VELOCITY
241  dest_arr(i,j,k) = prob->compute_terrain_velocity(time);
242  #else
243  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][2];
244  #endif
245  } else if (bc_ptr[0].lo(2) == ERFBCType::foextrap) {
246  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
247  } else if (bc_ptr[0].lo(2) == ERFBCType::open) {
248  dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z);
249  } else if (bc_ptr[0].lo(2) == ERFBCType::reflect_even) {
250  dest_arr(i,j,k) = dest_arr(i,j,kflip);
251  } else if (bc_ptr[0].lo(2) == ERFBCType::reflect_odd) {
252  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
253  } else if (bc_ptr[0].lo(2) == ERFBCType::hoextrap) {
254  Real delta_k = (dom_lo.z - k);
255  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);
256  }
257  },
258  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
259  int kflip = 2*dom_hi.z + 1 - k;
260  if (bc_ptr[0].hi(2) == ERFBCType::ext_dir) {
261  dest_arr(i,j,k) = l_bc_extdir_vals_d[0][5];
262  } else if (bc_ptr[0].hi(2) == ERFBCType::foextrap) {
263  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
264  } else if (bc_ptr[0].hi(2) == ERFBCType::open) {
265  dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z);
266  } else if (bc_ptr[0].hi(2) == ERFBCType::reflect_even) {
267  dest_arr(i,j,k) = dest_arr(i,j,kflip);
268  } else if (bc_ptr[0].hi(2) == ERFBCType::reflect_odd) {
269  dest_arr(i,j,k) = -dest_arr(i,j,kflip);
270  } else if (bc_ptr[0].hi(2) == ERFBCType::hoextrap){
271  Real delta_k = (k - dom_hi.z);
272  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);
273  }
274  }
275  );
276  } // z
277 
278  if (m_z_phys_nd) {
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_AtIface (ii, jj, k0, dxInv, z_phys_nd);
312  Real met_h_eta = Compute_h_eta_AtIface (ii, jj, k0, dxInv, z_phys_nd);
313  Real met_h_zeta = Compute_h_zeta_AtIface(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 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 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_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:110
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:96
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:125
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:192
Here is the call graph for this function:

◆ operator()()

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