REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_BoundaryConditions_cons.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 ] valid_bx valid box
10 * @param[in ] domain domain box
11 * @param[in ] dxInv pm or pn
12 * @param[in ] mskr land-sea mask on rho-points
13 * @param[in ] msku land-sea mask on u-points
14 * @param[in ] mskv land-sea mask on v-points
15 * @param[in ] calc_arr data to use in the RHS of calculations
16 * @param[in ] icomp component to update
17 * @param[in ] ncomp number of components to update, starting from icomp
18 * @param[in ] time current time
19 * @param[in ] bccomp index into both domain_bcs_type_bcr and bc_extdir_vals for icomp=0
20 * @param[in ] n_not_fill perimter of cells in x and y where BCs are not applied for non-ext_dir conditions
21 */
22void REMORAPhysBCFunct::impose_cons_bcs (const Array4<Real>& dest_arr, const Box& bx, const Box& valid_bx, const Box& domain,
23 const GpuArray<Real,AMREX_SPACEDIM> /*dxInv*/, const Array4<const Real>& mskr,
24 const Array4<const Real>& msku, const Array4<const Real>& mskv,
25 const Array4<const Real>& calc_arr,
26 int icomp, int ncomp, Real /*time*/, int bccomp, int n_not_fill)
27{
28 BL_PROFILE_VAR("impose_cons_bcs()",impose_cons_bcs);
29 const auto& dom_lo = amrex::lbound(domain);
30 const auto& dom_hi = amrex::ubound(domain);
31
32 // Based on BCRec for the domain, we need to make BCRec for this Box
33 // bccomp is used as starting index for m_domain_bcs_type
34 // 0 is used as starting index for bcrs
35 Vector<BCRec> bcrs(ncomp);
36 amrex::setBC(bx, domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
37
38 // xlo: ori = 0
39 // ylo: ori = 1
40 // zlo: ori = 2
41 // xhi: ori = 3
42 // yhi: ori = 4
43 // zhi: ori = 5
44
45 amrex::Gpu::DeviceVector<BCRec> bcrs_d(ncomp);
46#ifdef AMREX_USE_GPU
47 Gpu::htod_memcpy_async(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
48#else
49 std::memcpy(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
50#endif
51 const amrex::BCRec* bc_ptr = bcrs_d.data();
52
53 GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NCONS+8> l_bc_extdir_vals_d;
54
55 for (int i = 0; i < ncomp; i++) {
56 for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++) {
57 l_bc_extdir_vals_d[i][ori] = m_bc_extdir_vals[bccomp+i][ori];
58 }
59 }
60
61 GeometryData const& geomdata = m_geom.data();
62 bool is_periodic_in_x = geomdata.isPeriodic(0);
63 bool is_periodic_in_y = geomdata.isPeriodic(1);
64 const Real eps= 1.0e-20_rt;
65
66 // First do all ext_dir bcs
67 if (!is_periodic_in_x)
68 {
69 Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
70 bx_xlo.setSmall(1,std::max(valid_bx.smallEnd(1)-1,dom_lo.y)); bx_xlo.setBig(1,std::min(valid_bx.bigEnd(1)+1,dom_hi.y));
71 Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
72 bx_xhi.setSmall(1,std::max(valid_bx.smallEnd(1)-1,dom_lo.y)); bx_xhi.setBig(1,std::min(valid_bx.bigEnd(1)+1,dom_hi.y));
73 ParallelFor(
74 bx_xlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
75 if (bc_ptr[n].lo(0) == REMORABCType::ext_dir) {
76 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][0] * mskr(i,j,0);
77 } else if (bc_ptr[n].lo(0) == REMORABCType::orlanski_rad) {
78 Real grad_lo = (calc_arr(dom_lo.x ,j ,k,icomp+n) - calc_arr(dom_lo.x ,j-1,k,icomp+n)) * mskv(i,j,0);
79 Real grad_lo_jp1 = (calc_arr(dom_lo.x ,j+1,k,icomp+n) - calc_arr(dom_lo.x ,j ,k,icomp+n)) * mskv(i,j,0);
80 Real dTdt = calc_arr(dom_lo.x,j,k,icomp+n) - dest_arr(dom_lo.x ,j,k,icomp+n);
81 Real dTdx = dest_arr(dom_lo.x,j,k,icomp+n) - dest_arr(dom_lo.x+1,j,k,icomp+n);
82 if (dTdt*dTdx < 0.0_rt) dTdt = 0.0_rt;
83 Real dTde = (dTdt * (grad_lo+grad_lo_jp1) > 0.0_rt) ? grad_lo : grad_lo_jp1;
84 Real cff = std::max(dTdx*dTdx+dTde*dTde,eps);
85 Real Cx = dTdt * dTdx;
86 dest_arr(i,j,k,icomp+n) = (cff * calc_arr(dom_lo.x-1,j,k,icomp+n) + Cx * dest_arr(dom_lo.x,j,k,icomp+n)) * mskr(i,j,0) / (cff+Cx);
87 }
88 },
89 bx_xhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
90 if (bc_ptr[n].hi(0) == REMORABCType::ext_dir) {
91 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][3] * mskr(i,j,0);
92 } else if (bc_ptr[n].hi(0) == REMORABCType::orlanski_rad) {
93 Real grad_hi = (calc_arr(dom_hi.x ,j ,k,icomp+n) - calc_arr(dom_hi.x ,j-1,k,icomp+n)) * mskv(i,j,0);
94 Real grad_hi_jp1 = (calc_arr(dom_hi.x ,j+1,k,icomp+n) - calc_arr(dom_hi.x ,j ,k,icomp+n)) * mskv(i,j,0);
95 Real dTdt = calc_arr(dom_hi.x,j,k,icomp+n) - dest_arr(dom_hi.x ,j,k,icomp+n);
96 Real dTdx = dest_arr(dom_hi.x,j,k,icomp+n) - dest_arr(dom_hi.x-1,j,k,icomp+n);
97 if (dTdt * dTdx < 0.0_rt) dTdt = 0.0_rt;
98 Real dTde = (dTdt * (grad_hi + grad_hi_jp1) > 0.0_rt) ? grad_hi : grad_hi_jp1;
99 Real cff = std::max(dTdx*dTdx + dTde*dTde,eps);
100 Real Cx = dTdt * dTdx;
101 dest_arr(i,j,k,icomp+n) = (cff * calc_arr(dom_hi.x+1,j,k,icomp+n) + Cx * dest_arr(dom_hi.x,j,k,icomp+n)) * mskr(i,j,0) / (cff+Cx);
102 }
103 }
104 );
105 }
106
107 if (!is_periodic_in_y)
108 {
109 Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
110 bx_ylo.setSmall(0,std::max(valid_bx.smallEnd(0)-1,dom_lo.x)); bx_ylo.setBig(0,std::min(valid_bx.bigEnd(0)+1,dom_hi.x));
111 Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
112 bx_yhi.setSmall(0,std::max(valid_bx.smallEnd(0)-1,dom_lo.x)); bx_yhi.setBig(0,std::min(valid_bx.bigEnd(0)+1,dom_hi.x));
113 ParallelFor(
114 bx_ylo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
115 if (bc_ptr[n].lo(1) == REMORABCType::ext_dir) {
116 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][1] * mskr(i,j,0);
117 } else if (bc_ptr[n].lo(1) == REMORABCType::orlanski_rad) {
118 Real grad_lo = (calc_arr(i ,dom_lo.y, k,icomp+n) - calc_arr(i-1,dom_lo.y ,k,icomp+n)) * msku(i,j,0);
119 Real grad_lo_ip1 = (calc_arr(i+1,dom_lo.y ,k,icomp+n) - calc_arr(i ,dom_lo.y ,k,icomp+n)) * msku(i,j,0);
120 Real dTdt = calc_arr(i,dom_lo.y,k,icomp+n) - dest_arr(i,dom_lo.y ,k,icomp+n);
121 Real dTde = dest_arr(i,dom_lo.y,k,icomp+n) - dest_arr(i,dom_lo.y+1,k,icomp+n);
122 if (dTdt * dTde < 0.0_rt) dTdt = 0.0_rt;
123 Real dTdx = (dTdt * (grad_lo + grad_lo_ip1) > 0.0_rt) ? grad_lo : grad_lo_ip1;
124 Real cff = std::max(dTdx*dTdx + dTde*dTde, eps);
125 Real Ce = dTdt*dTde;
126 dest_arr(i,j,k,icomp+n) = (cff * calc_arr(i,dom_lo.y-1,k,icomp+n) + Ce * dest_arr(i,dom_lo.y,k,icomp+n)) * mskr(i,j,0) / (cff+Ce);
127 }
128 },
129 bx_yhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
130 if (bc_ptr[n].hi(1) == REMORABCType::ext_dir) {
131 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][4] * mskr(i,j,0);
132 } else if (bc_ptr[n].hi(1) == REMORABCType::orlanski_rad) {
133 Real grad_hi = (calc_arr(i ,dom_hi.y ,k,icomp+n) - calc_arr(i-1,dom_hi.y ,k,icomp+n)) * msku(i,j,0);
134 Real grad_hi_ip1 = (calc_arr(i+1,dom_hi.y ,k,icomp+n) - calc_arr(i ,dom_hi.y ,k,icomp+n)) * msku(i,j,0);
135 Real dTdt = calc_arr(i,dom_hi.y,k,icomp+n) - dest_arr(i,dom_hi.y ,k,icomp+n);
136 Real dTde = dest_arr(i,dom_hi.y,k,icomp+n) - dest_arr(i,dom_hi.y-1,k,icomp+n);
137 if (dTdt * dTde < 0.0_rt) dTdt = 0.0_rt;
138 Real dTdx = (dTdt * (grad_hi + grad_hi_ip1) > 0.0_rt) ? grad_hi : grad_hi_ip1;
139 Real cff = std::max(dTdx*dTdx + dTde*dTde, eps);
140 Real Ce = dTdt*dTde;
141 dest_arr(i,j,k,icomp+n) = (cff*calc_arr(i,dom_hi.y+1,k,icomp+n) + Ce*dest_arr(i,dom_hi.y,k,icomp+n)) * mskr(i,j,0) / (cff+Ce);
142 }
143 }
144 );
145 }
146
147 {
148 Box bx_zlo(bx); bx_zlo.setBig (2,dom_lo.z-1);
149 Box bx_zhi(bx); bx_zhi.setSmall(2,dom_hi.z+1);
150 ParallelFor(
151 bx_zlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
152 if (bc_ptr[n].lo(2) == REMORABCType::ext_dir) {
153 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][2] * mskr(i,j,0);
154 }
155 },
156 bx_zhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
157 if (bc_ptr[n].hi(2) == REMORABCType::ext_dir) {
158 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][5] * mskr(i,j,0);
159 }
160 }
161 );
162 }
163
164 Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1-n_not_fill);
165 bx_xlo.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
166 bx_xlo.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
167 Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1+n_not_fill);
168 bx_xhi.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
169 bx_xhi.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
170 Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1-n_not_fill);
171 bx_ylo.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
172 bx_ylo.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
173 Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1+n_not_fill);
174 bx_yhi.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
175 bx_yhi.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
176 // Calculate intersections for corners before adjusting to exclude them
177 Box xlo_ylo = bx_xlo & bx_ylo;
178 Box xhi_ylo = bx_xhi & bx_ylo;
179 Box xlo_yhi = bx_xlo & bx_yhi;
180 Box xhi_yhi = bx_xhi & bx_yhi;
181// bx_xlo.setSmall(1,valid_bx.smallEnd(1)); bx_xlo.setBig(1,valid_bx.bigEnd(1));
182// bx_xhi.setSmall(1,valid_bx.smallEnd(1)); bx_xhi.setBig(1,valid_bx.bigEnd(1));
183// bx_ylo.setSmall(0,valid_bx.smallEnd(0)); bx_ylo.setBig(0,valid_bx.bigEnd(0));
184// bx_yhi.setSmall(0,valid_bx.smallEnd(0)); bx_yhi.setBig(0,valid_bx.bigEnd(0));
185 // Next do ghost cells in x-direction but not reaching out in y
186 // The corners we miss here will be covered in the y-loop below or by periodicity
187 if (!is_periodic_in_x or bccomp==BCVars::foextrap_bc)
188 {
189 // Populate ghost cells on lo-x and hi-x domain boundaries
190 ParallelFor(bx_xlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
191 int iflip = dom_lo.x - 1 - i;
192 int inner = (bc_ptr[n].lo(0) == REMORABCType::orlanski_rad) ? 1 : 0;
193 if (bc_ptr[n].lo(0) == REMORABCType::foextrap || bc_ptr[n].lo(0) == REMORABCType::clamped || bc_ptr[n].lo(0) == REMORABCType::chapman || bc_ptr[n].lo(0) == REMORABCType::orlanski_rad ||
194 bc_ptr[n].lo(0) == REMORABCType::orlanski_rad_nudge) {
195 dest_arr(i,j,k,icomp+n) = dest_arr(dom_lo.x-n_not_fill-inner,j,k,icomp+n);
196 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_even) {
197 dest_arr(i,j,k,icomp+n) = dest_arr(iflip,j,k,icomp+n);
198 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_odd) {
199 dest_arr(i,j,k,icomp+n) = -dest_arr(iflip,j,k,icomp+n);
200 }
201 },
202 bx_xhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
203 int iflip = 2*dom_hi.x + 1 - i;
204 int inner = (bc_ptr[n].hi(0) == REMORABCType::orlanski_rad) ? 1 : 0;
205 if (bc_ptr[n].hi(0) == REMORABCType::foextrap || bc_ptr[n].hi(0) == REMORABCType::clamped || bc_ptr[n].hi(0) == REMORABCType::chapman || bc_ptr[n].hi(0) == REMORABCType::orlanski_rad ||
206 bc_ptr[n].hi(0) == REMORABCType::orlanski_rad_nudge) {
207 dest_arr(i,j,k,icomp+n) = dest_arr(dom_hi.x+n_not_fill+inner,j,k,icomp+n);
208 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_even) {
209 dest_arr(i,j,k,icomp+n) = dest_arr(iflip,j,k,icomp+n);
210 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_odd) {
211 dest_arr(i,j,k,icomp+n) = -dest_arr(iflip,j,k,icomp+n);
212 }
213 }
214 );
215 }
216
217 if (!is_periodic_in_y or bccomp==BCVars::foextrap_bc)
218 {
219 // Populate ghost cells on lo-y and hi-y domain boundaries
220 ParallelFor(
221 bx_ylo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
222 int jflip = dom_lo.y - 1 - j;
223 int inner = (bc_ptr[n].lo(1) == REMORABCType::orlanski_rad) ? 1 : 0;
224 if (bc_ptr[n].lo(1) == REMORABCType::foextrap || bc_ptr[n].lo(1) == REMORABCType::clamped || bc_ptr[n].lo(1) == REMORABCType::chapman || bc_ptr[n].lo(1) == REMORABCType::orlanski_rad ||
225 bc_ptr[n].lo(1) == REMORABCType::orlanski_rad_nudge) {
226 dest_arr(i,j,k,icomp+n) = dest_arr(i,dom_lo.y-n_not_fill-inner,k,icomp+n);
227 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_even) {
228 dest_arr(i,j,k,icomp+n) = dest_arr(i,jflip,k,icomp+n);
229 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_odd) {
230 dest_arr(i,j,k,icomp+n) = -dest_arr(i,jflip,k,icomp+n);
231 }
232 },
233 bx_yhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
234 int jflip = 2*dom_hi.y + 1 - j;
235 int inner = (bc_ptr[n].hi(1) == REMORABCType::orlanski_rad) ? 1 : 0;
236 if (bc_ptr[n].hi(1) == REMORABCType::foextrap || bc_ptr[n].hi(1) == REMORABCType::clamped || bc_ptr[n].hi(1) == REMORABCType::chapman || bc_ptr[n].lo(1) == REMORABCType::orlanski_rad ||
237 bc_ptr[n].hi(1) == REMORABCType::orlanski_rad_nudge) {
238 dest_arr(i,j,k,icomp+n) = dest_arr(i,dom_hi.y+n_not_fill+inner,k,icomp+n);
239 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_even) {
240 dest_arr(i,j,k,icomp+n) = dest_arr(i,jflip,k,icomp+n);
241 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_odd) {
242 dest_arr(i,j,k,icomp+n) = -dest_arr(i,jflip,k,icomp+n);
243 }
244 }
245 );
246 }
247
248 {
249 Box bx_zlo(bx); bx_zlo.setBig (2,std::max(dom_lo.z-1,bx.smallEnd(2)));
250 Box bx_zhi(bx); bx_zhi.setSmall(2,std::min(dom_hi.z+1,bx.bigEnd(2)));
251 // Populate ghost cells on lo-z and hi-z domain boundaries
252
253 if (bx_zlo.ok()) {
254 ParallelFor(bx_zlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
255 {
256 int kflip = dom_lo.z - 1 - i;
257 if (bc_ptr[n].lo(2) == REMORABCType::foextrap) {
258 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,dom_lo.z,icomp+n);
259 } else if (bc_ptr[n].lo(2) == REMORABCType::reflect_even) {
260 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,kflip,icomp+n);
261 } else if (bc_ptr[n].lo(2) == REMORABCType::reflect_odd) {
262 dest_arr(i,j,k,icomp+n) = -dest_arr(i,j,kflip,icomp+n);
263 }
264 });
265 }
266
267 if (bx_zlo.ok()) {
268 ParallelFor(bx_zhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
269 {
270 int kflip = 2*dom_hi.z + 1 - i;
271 if (bc_ptr[n].hi(2) == REMORABCType::foextrap) {
272 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,dom_hi.z,icomp+n);
273 } else if (bc_ptr[n].hi(2) == REMORABCType::reflect_even) {
274 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,kflip,icomp+n);
275 } else if (bc_ptr[n].hi(2) == REMORABCType::reflect_odd) {
276 dest_arr(i,j,k,icomp+n) = -dest_arr(i,j,kflip,icomp+n);
277 }
278 });
279 }
280 }
281 if ((!is_periodic_in_x && !is_periodic_in_y) or bccomp==BCVars::foextrap_bc) {
282 // If we've applied boundary conditions to either side, update the corner
283 if (!xlo_ylo.isEmpty()) {
284 ParallelFor(xlo_ylo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
285 {
286 if (!(bc_ptr[n].lo(0) == REMORABCType::clamped || bc_ptr[n].lo(0) == REMORABCType::flather ||
287 bc_ptr[n].lo(0) == REMORABCType::chapman || bc_ptr[n].lo(0) == REMORABCType::orlanski_rad_nudge)
288 && !(bc_ptr[n].lo(1) == REMORABCType::clamped || bc_ptr[n].lo(1) == REMORABCType::flather ||
289 bc_ptr[n].lo(1) == REMORABCType::chapman || bc_ptr[n].lo(1) == REMORABCType::orlanski_rad_nudge)) {
290 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_lo.y,k,icomp+n)
291 + dest_arr(dom_lo.x,j,k,icomp+n));
292 }
293 });
294 }
295 if (!xlo_yhi.isEmpty()) {
296 ParallelFor(xlo_yhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
297 {
298 if (!(bc_ptr[n].lo(0) == REMORABCType::clamped || bc_ptr[n].lo(0) == REMORABCType::flather ||
299 bc_ptr[n].lo(0) == REMORABCType::chapman || bc_ptr[n].lo(0) == REMORABCType::orlanski_rad_nudge)
300 && !(bc_ptr[n].hi(1) == REMORABCType::clamped || bc_ptr[n].hi(1) == REMORABCType::flather ||
301 bc_ptr[n].hi(1) == REMORABCType::chapman || bc_ptr[n].hi(1) == REMORABCType::orlanski_rad_nudge)) {
302 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_hi.y,k,icomp+n)
303 + dest_arr(dom_lo.x,j,k,icomp+n));
304 }
305 });
306 }
307 if (!xhi_ylo.isEmpty()) {
308 ParallelFor(xhi_ylo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
309 {
310 if (!(bc_ptr[n].hi(0) == REMORABCType::clamped || bc_ptr[n].hi(0) == REMORABCType::flather ||
311 bc_ptr[n].hi(0) == REMORABCType::chapman || bc_ptr[n].hi(0) == REMORABCType::orlanski_rad_nudge)
312 && !(bc_ptr[n].lo(1) == REMORABCType::clamped || bc_ptr[n].lo(1) == REMORABCType::flather ||
313 bc_ptr[n].lo(1) == REMORABCType::chapman || bc_ptr[n].lo(1) == REMORABCType::orlanski_rad_nudge)) {
314 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_lo.y,k,icomp+n)
315 + dest_arr(dom_hi.x,j,k,icomp+n));
316 }
317 });
318 }
319 if (!xhi_yhi.isEmpty()) {
320 ParallelFor(xhi_yhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
321 {
322 if (!(bc_ptr[n].hi(0) == REMORABCType::clamped || bc_ptr[n].hi(0) == REMORABCType::flather ||
323 bc_ptr[n].hi(0) == REMORABCType::chapman || bc_ptr[n].hi(0) == REMORABCType::orlanski_rad_nudge)
324 && !(bc_ptr[n].hi(1) == REMORABCType::clamped || bc_ptr[n].hi(1) == REMORABCType::flather ||
325 bc_ptr[n].hi(1) == REMORABCType::chapman || bc_ptr[n].hi(1) == REMORABCType::orlanski_rad_nudge)) {
326 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_hi.y,k,icomp+n)
327 + dest_arr(dom_hi.x,j,k,icomp+n));
328 }
329 });
330 }
331 }
332
333 Gpu::streamSynchronize();
334}
#define NCONS
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NCONS+8 > m_bc_extdir_vals
amrex::Vector< amrex::BCRec > m_domain_bcs_type
void impose_cons_bcs(const amrex::Array4< amrex::Real > &mf, const amrex::Box &bx, const amrex::Box &valid_bx, const amrex::Box &domain, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > dxInv, const amrex::Array4< const amrex::Real > &mskr, const amrex::Array4< const amrex::Real > &msku, const amrex::Array4< const amrex::Real > &mskv, const amrex::Array4< const amrex::Real > &calc_arr, int icomp, int ncomp, amrex::Real time, int bccomp, int n_not_fill)
apply scalar type boundary conditions
amrex::Geometry m_geom