REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_BoundaryConditions_xvel.cpp
Go to the documentation of this file.
1#include "AMReX_PhysBCFunct.H"
3
4using namespace amrex;
5
6/**
7 * @param[inout] dest_arr data on which to apply BCs
8 * @param[in ] bx box to update on
9 * @param[in ] domain domain box
10 * @param[in ] dxInv pm or pn
11 * @param[in ] msku land-sea mask on u-points
12 * @param[in ] calc_arr data to use in the RHS of calculations
13 * @param[in ] time current time
14 * @param[in ] bccomp index into both domain_bcs_type_bcr and bc_extdir_vals for icomp=0
15 */
16void REMORAPhysBCFunct::impose_xvel_bcs (const Array4<Real>& dest_arr, const Box& bx, const Box& domain,
17 const GpuArray<Real,AMREX_SPACEDIM> /*dxInv*/, const Array4<const Real>& msku,
18 const Array4<const Real>& calc_arr,
19 Real /*time*/, int bccomp)
20{
21 BL_PROFILE_VAR("impose_xvel_bcs()",impose_xvel_bcs);
22 const auto& dom_lo = amrex::lbound(domain);
23 const auto& dom_hi = amrex::ubound(domain);
24
25 // Based on BCRec for the domain, we need to make BCRec for this Box
26 // bccomp is used as starting index for m_domain_bcs_type
27 // 0 is used as starting index for bcrs
28 int ncomp = 1;
29 Vector<BCRec> bcrs(ncomp);
30 amrex::setBC(bx, domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
31
32 // xlo: ori = 0
33 // ylo: ori = 1
34 // zlo: ori = 2
35 // xhi: ori = 3
36 // yhi: ori = 4
37 // zhi: ori = 5
38
39 amrex::Gpu::DeviceVector<BCRec> bcrs_d(ncomp);
40#ifdef AMREX_USE_GPU
41 Gpu::htod_memcpy_async(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
42#else
43 std::memcpy(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
44#endif
45 const amrex::BCRec* bc_ptr = bcrs_d.data();
46
47 GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NCONS+8> l_bc_extdir_vals_d;
48
49 for (int i = 0; i < ncomp; i++)
50 for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++)
51 l_bc_extdir_vals_d[i][ori] = m_bc_extdir_vals[bccomp+i][ori];
52
53 GeometryData const& geomdata = m_geom.data();
54 bool is_periodic_in_x = geomdata.isPeriodic(0);
55 bool is_periodic_in_y = geomdata.isPeriodic(1);
56 const Real eps= 1.0e-20_rt;
57
58
59 // First do all ext_dir bcs
60 if (!is_periodic_in_x or bccomp==BCVars::foextrap_bc)
61 {
62 Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
63 Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+2);
64 Box bx_xlo_face(bx); bx_xlo_face.setSmall(0,dom_lo.x ); bx_xlo_face.setBig(0,dom_lo.x );
65 Box bx_xhi_face(bx); bx_xhi_face.setSmall(0,dom_hi.x+1); bx_xhi_face.setBig(0,dom_hi.x+1);
66 ParallelFor(
67 // We only set the values on the domain faces themselves if EXT_DIR or actual outflow
68 grow(bx_xlo_face,IntVect(0,-1,0)), ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
69 if (bc_ptr[n].lo(0) == REMORABCType::ext_dir) {
70 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][0]*msku(i,j,0);
71 } else if (bc_ptr[n].lo(0) == REMORABCType::foextrap) {
72 dest_arr(i,j,k) = dest_arr(dom_lo.x+1,j,k)*msku(i,j,0);
73 } else if (bc_ptr[n].lo(0) == REMORABCType::orlanski_rad) {
74 Real grad_lo_ip1 = calc_arr(dom_lo.x+1,j ,k) - calc_arr(dom_lo.x+1,j-1,k);
75 Real grad_lo_ijp1 = calc_arr(dom_lo.x+1,j+1,k) - calc_arr(dom_lo.x+1,j ,k);
76 Real dUdt = calc_arr(dom_lo.x+1,j,k) - dest_arr(dom_lo.x+1,j,k);
77 Real dUdx = dest_arr(dom_lo.x+1,j,k) - dest_arr(dom_lo.x+2,j,k);
78 if (dUdt * dUdx < 0.0_rt) dUdt = 0.0_rt;
79 Real dUde = (dUdt * (grad_lo_ip1 + grad_lo_ijp1) > 0.0_rt) ? grad_lo_ip1 : grad_lo_ijp1;
80 Real cff = std::max(dUdx*dUdx+dUde*dUde,eps);
81 Real Cx = dUdt*dUdx;
82 dest_arr(i,j,k) = (cff * calc_arr(i,j,k) + Cx * dest_arr(dom_lo.x+1,j,k)) * msku(i,j,0) / (cff+Cx);
83 }
84 });
85 ParallelFor(
86 grow(bx_xlo,IntVect(0,-1,0)), ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
87 int inner = (bc_ptr[n].lo(0) == REMORABCType::foextrap) ? 1 : 0;
88 int iflip = dom_lo.x - i;
89 if (bc_ptr[n].lo(0) == REMORABCType::ext_dir) {
90 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][0]*msku(i,j,0);
91 } else if (bc_ptr[n].lo(0) == REMORABCType::foextrap || bc_ptr[n].lo(0) == REMORABCType::clamped ||
92 bc_ptr[n].lo(0) == REMORABCType::orlanski_rad) {
93 dest_arr(i,j,k) = dest_arr(dom_lo.x+inner,j,k)*msku(i,j,0);
94 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_even) {
95 dest_arr(i,j,k) = dest_arr(iflip,j,k)*msku(i,j,0);
96 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_odd) {
97 dest_arr(i,j,k) = -dest_arr(iflip,j,k)*msku(i,j,0);
98 }
99 });
100 // We only set the values on the domain faces themselves if EXT_DIR or actual outflow
101 ParallelFor(
102 grow(bx_xhi_face,IntVect(0,-1,0)), ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
103 if (bc_ptr[n].hi(0) == REMORABCType::ext_dir) {
104 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][3]*msku(i,j,0);
105 } else if (bc_ptr[n].hi(0) == REMORABCType::foextrap) {
106 dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k)*msku(i,j,0);
107 } else if (bc_ptr[n].hi(0) == REMORABCType::orlanski_rad) {
108 Real grad_hi = calc_arr(dom_hi.x ,j ,k) - calc_arr(dom_hi.x ,j-1,k);
109 Real grad_hi_jp1 = calc_arr(dom_hi.x ,j+1,k) - calc_arr(dom_hi.x ,j ,k);
110 Real dUdt = calc_arr(dom_hi.x,j,k) - dest_arr(dom_hi.x ,j,k);
111 Real dUdx = dest_arr(dom_hi.x,j,k) - dest_arr(dom_hi.x-1,j,k);
112 if (dUdt * dUdx < 0.0_rt) dUdt = 0.0_rt;
113 Real dUde = (dUdt * (grad_hi + grad_hi_jp1) > 0.0_rt) ? grad_hi : grad_hi_jp1;
114 Real cff = std::max(dUdx*dUdx+dUde*dUde,eps);
115 Real Cx = dUdt * dUdx;
116 dest_arr(i,j,k) = (cff * calc_arr(dom_hi.x+1,j,k) + Cx * dest_arr(dom_hi.x,j,k)) * msku(i,j,0) / (cff + Cx);
117 }
118 });
119 ParallelFor(
120 grow(bx_xhi,IntVect(0,-1,0)), ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
121 int iflip = 2*(dom_hi.x + 1) - i;
122 int inner = (bc_ptr[n].hi(0) == REMORABCType::foextrap) ? 1 : 0;
123 if (bc_ptr[n].hi(0) == REMORABCType::ext_dir) {
124 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][3]*msku(i,j,0);
125 } else if (bc_ptr[n].hi(0) == REMORABCType::foextrap || bc_ptr[n].hi(0) == REMORABCType::clamped ||
126 bc_ptr[n].hi(0) == REMORABCType::orlanski_rad) {
127 dest_arr(i,j,k) = dest_arr(dom_hi.x+1-inner,j,k)*msku(i,j,0);
128 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_even) {
129 dest_arr(i,j,k) = dest_arr(iflip,j,k)*msku(i,j,0);
130 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_odd) {
131 dest_arr(i,j,k) = -dest_arr(iflip,j,k)*msku(i,j,0);
132 }
133 });
134 } // not periodic in x
135
136 if (!is_periodic_in_y or bccomp==BCVars::foextrap_bc)
137 {
138 // Populate ghost cells on lo-y and hi-y domain boundaries
139 Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
140 Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
141 ParallelFor(
142 grow(bx_ylo,IntVect(-1,0,0)), ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
143 int jflip = dom_lo.y - 1 - j;
144 if (bc_ptr[n].lo(1) == REMORABCType::ext_dir) {
145 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][1]*msku(i,j,0);
146 } else if (bc_ptr[n].lo(1) == REMORABCType::foextrap || bc_ptr[n].lo(1) == REMORABCType::clamped) {
147 dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k)*msku(i,j,0);
148 } else if (bc_ptr[n].lo(1) == REMORABCType::orlanski_rad) {
149 Real grad_lo = calc_arr(i+1,dom_lo.y ,k) - calc_arr(i ,dom_lo.y ,k);
150 Real grad_lo_im1 = calc_arr(i ,dom_lo.y ,k) - calc_arr(i-1,dom_lo.y ,k);
151 Real dUdt = calc_arr(i,dom_lo.y,k) - dest_arr(i,dom_lo.y ,k);
152 Real dUde = dest_arr(i,dom_lo.y,k) - dest_arr(i,dom_lo.y+1,k);
153 if (dUdt * dUde < 0.0_rt) dUdt = 0.0_rt;
154 Real dUdx = (dUdt * (grad_lo_im1 + grad_lo) > 0.0_rt) ? grad_lo_im1 : grad_lo;
155 Real cff = std::max(dUdx * dUdx + dUde * dUde, eps);
156 Real Ce = dUdt * dUde;
157 dest_arr(i,j,k) = (cff * calc_arr(i,dom_lo.y-1,k) + Ce * dest_arr(i,dom_lo.y,k)) * msku(i,j,0) / (cff + Ce);
158 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_even) {
159 dest_arr(i,j,k) = dest_arr(i,jflip,k)*msku(i,j,0);
160 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_odd) {
161 dest_arr(i,j,k) = -dest_arr(i,jflip,k)*msku(i,j,0);
162 }
163 },
164 grow(bx_yhi,IntVect(-1,0,0)), ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
165 int jflip = 2*dom_hi.y + 1 - j;
166 if (bc_ptr[n].hi(1) == REMORABCType::ext_dir) {
167 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][4]*msku(i,j,0);
168 } else if (bc_ptr[n].hi(1) == REMORABCType::foextrap || bc_ptr[n].hi(1) == REMORABCType::clamped) {
169 dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k)*msku(i,j,0);
170 } else if (bc_ptr[n].hi(1) == REMORABCType::orlanski_rad) {
171 Real grad_hi = calc_arr(i+1,dom_hi.y ,k) - calc_arr(i ,dom_hi.y ,k);
172 Real grad_hi_im1 = calc_arr(i ,dom_hi.y ,k) - calc_arr(i-1,dom_hi.y ,k);
173 Real dUdt = calc_arr(i,dom_hi.y,k) - dest_arr(i,dom_hi.y ,k);
174 Real dUde = dest_arr(i,dom_hi.y,k) - dest_arr(i,dom_hi.y-1,k);
175 if (dUdt * dUde < 0.0_rt) dUdt = 0.0_rt;
176 Real dUdx = (dUdt * (grad_hi_im1 + grad_hi) > 0.0_rt) ? grad_hi_im1 : grad_hi;
177 Real cff = std::max(dUdx*dUdx+dUde*dUde,eps);
178 Real Ce = dUdt * dUde;
179 dest_arr(i,j,k) = (cff * calc_arr(i,dom_hi.y+1,k) + Ce * dest_arr(i,dom_hi.y,k)) * msku(i,j,0) / (cff+Ce);
180 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_even) {
181 dest_arr(i,j,k) = dest_arr(i,jflip,k)*msku(i,j,0);
182 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_odd) {
183 dest_arr(i,j,k) = -dest_arr(i,jflip,k)*msku(i,j,0);
184 }
185 }
186 );
187 } // not periodic in y
188
189 {
190 // Populate ghost cells on lo-z and hi-z domain boundaries
191 Box bx_zlo(bx); bx_zlo.setBig (2,dom_lo.z-1);
192 Box bx_zhi(bx); bx_zhi.setSmall(2,dom_hi.z+1);
193 ParallelFor(
194 bx_zlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
195 int kflip = dom_lo.z - 1 - k;
196 if (bc_ptr[n].lo(2) == REMORABCType::ext_dir) {
197 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][2]*msku(i,j,0);
198 } else if (bc_ptr[n].lo(2) == REMORABCType::foextrap) {
199 dest_arr(i,j,k) = dest_arr(i,j,dom_lo.z)*msku(i,j,0);
200 } else if (bc_ptr[n].lo(2) == REMORABCType::reflect_even) {
201 dest_arr(i,j,k) = dest_arr(i,j,kflip)*msku(i,j,0);
202 } else if (bc_ptr[n].lo(2) == REMORABCType::reflect_odd) {
203 dest_arr(i,j,k) = -dest_arr(i,j,kflip)*msku(i,j,0);
204 }
205 },
206 bx_zhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
207 int kflip = 2*dom_hi.z + 1 - k;
208 if (bc_ptr[n].hi(2) == REMORABCType::ext_dir) {
209 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][5]*msku(i,j,0);
210 } else if (bc_ptr[n].hi(2) == REMORABCType::foextrap) {
211 dest_arr(i,j,k) = dest_arr(i,j,dom_hi.z)*msku(i,j,0);
212 } else if (bc_ptr[n].hi(2) == REMORABCType::reflect_even) {
213 dest_arr(i,j,k) = dest_arr(i,j,kflip)*msku(i,j,0);
214 } else if (bc_ptr[n].hi(2) == REMORABCType::reflect_odd) {
215 dest_arr(i,j,k) = -dest_arr(i,j,kflip)*msku(i,j,0);
216 }
217 }
218 );
219 } // z
220
221 if ((!is_periodic_in_x or bccomp==BCVars::foextrap_bc) and
222 (!is_periodic_in_y or bccomp==BCVars::foextrap_bc)) {
223 Box xlo(bx); xlo.setBig (0,dom_lo.x );
224 Box xhi(bx); xhi.setSmall(0,dom_hi.x+1);
225 Box ylo(bx); ylo.setBig (1,dom_lo.y-1);
226 Box yhi(bx); yhi.setSmall(1,dom_hi.y+1);
227 Box xlo_ylo = xlo & ylo;
228 Box xlo_yhi = xlo & yhi;
229 Box xhi_ylo = xhi & ylo;
230 Box xhi_yhi = xhi & yhi;
231
232 const bool clamp_west = m_domain_bcs_type[bccomp].lo(0) == REMORABCType::clamped;
233 const bool clamp_east = m_domain_bcs_type[bccomp].hi(0) == REMORABCType::clamped;
234 const bool clamp_south = m_domain_bcs_type[bccomp].lo(1) == REMORABCType::clamped;
235 const bool clamp_north = m_domain_bcs_type[bccomp].hi(1) == REMORABCType::clamped;
236
237 if (!clamp_west && !clamp_south) {
238 ParallelFor(xlo_ylo, [=] AMREX_GPU_DEVICE (int i, int j, int k)
239 {
240 dest_arr(i,j,k) = 0.5 * (dest_arr(i,dom_lo.y,k) + dest_arr(dom_lo.x+1,j,k));
241 });
242 }
243 if (!clamp_west && !clamp_north) {
244 ParallelFor(xlo_yhi, [=] AMREX_GPU_DEVICE (int i, int j, int k)
245 {
246 dest_arr(i,j,k) = 0.5 * (dest_arr(i,dom_hi.y,k) + dest_arr(dom_lo.x+1,j,k));
247 });
248 }
249 if (!clamp_east && !clamp_south) {
250 ParallelFor(xhi_ylo, [=] AMREX_GPU_DEVICE (int i, int j, int k)
251 {
252 dest_arr(i,j,k) = 0.5 * (dest_arr(i,dom_lo.y,k) + dest_arr(dom_hi.x,j,k));
253 });
254 }
255 if (!clamp_east && !clamp_north) {
256 ParallelFor(xhi_yhi, [=] AMREX_GPU_DEVICE (int i, int j, int k)
257 {
258 dest_arr(i,j,k) = 0.5 * (dest_arr(i,dom_hi.y,k) + dest_arr(dom_hi.x,j,k));
259 });
260 }
261 }
262
263
264 Gpu::streamSynchronize();
265}
#define NCONS
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NCONS+8 > m_bc_extdir_vals
void impose_xvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > dxInv, const amrex::Array4< const amrex::Real > &msku, const amrex::Array4< const amrex::Real > &calc_arr, amrex::Real time, int bccomp)
apply x-velocity type boundary conditions
amrex::Vector< amrex::BCRec > m_domain_bcs_type
amrex::Geometry m_geom