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 Box dest_arr_box = Box(dest_arr);
67 // First do all ext_dir bcs
68 if (!is_periodic_in_x)
69 {
70 Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
71 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));
72 Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
73 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));
74 ParallelFor(
75 bx_xlo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
76 if (bc_ptr[n].lo(0) == REMORABCType::ext_dir) {
77 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][0] * mskr(i,j,0);
78 } else if (bc_ptr[n].lo(0) == REMORABCType::orlanski_rad) {
79 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);
80 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);
81 Real dTdt = calc_arr(dom_lo.x,j,k,icomp+n) - dest_arr(dom_lo.x ,j,k,icomp+n);
82 Real dTdx = dest_arr(dom_lo.x,j,k,icomp+n) - dest_arr(dom_lo.x+1,j,k,icomp+n);
83 if (dTdt*dTdx < 0.0_rt) dTdt = 0.0_rt;
84 Real dTde = (dTdt * (grad_lo+grad_lo_jp1) > 0.0_rt) ? grad_lo : grad_lo_jp1;
85 Real cff = std::max(dTdx*dTdx+dTde*dTde,eps);
86 Real Cx = dTdt * dTdx;
87 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);
88 }
89 },
90 bx_xhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
91 if (bc_ptr[n].hi(0) == REMORABCType::ext_dir) {
92 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][3] * mskr(i,j,0);
93 } else if (bc_ptr[n].hi(0) == REMORABCType::orlanski_rad) {
94 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);
95 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);
96 Real dTdt = calc_arr(dom_hi.x,j,k,icomp+n) - dest_arr(dom_hi.x ,j,k,icomp+n);
97 Real dTdx = dest_arr(dom_hi.x,j,k,icomp+n) - dest_arr(dom_hi.x-1,j,k,icomp+n);
98 if (dTdt * dTdx < 0.0_rt) dTdt = 0.0_rt;
99 Real dTde = (dTdt * (grad_hi + grad_hi_jp1) > 0.0_rt) ? grad_hi : grad_hi_jp1;
100 Real cff = std::max(dTdx*dTdx + dTde*dTde,eps);
101 Real Cx = dTdt * dTdx;
102 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);
103 }
104 }
105 );
106 }
107
108 if (!is_periodic_in_y)
109 {
110 Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
111 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));
112 Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
113 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));
114 ParallelFor(
115 bx_ylo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
116 if (bc_ptr[n].lo(1) == REMORABCType::ext_dir) {
117 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][1] * mskr(i,j,0);
118 } else if (bc_ptr[n].lo(1) == REMORABCType::orlanski_rad) {
119 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);
120 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);
121 Real dTdt = calc_arr(i,dom_lo.y,k,icomp+n) - dest_arr(i,dom_lo.y ,k,icomp+n);
122 Real dTde = dest_arr(i,dom_lo.y,k,icomp+n) - dest_arr(i,dom_lo.y+1,k,icomp+n);
123 if (dTdt * dTde < 0.0_rt) dTdt = 0.0_rt;
124 Real dTdx = (dTdt * (grad_lo + grad_lo_ip1) > 0.0_rt) ? grad_lo : grad_lo_ip1;
125 Real cff = std::max(dTdx*dTdx + dTde*dTde, eps);
126 Real Ce = dTdt*dTde;
127 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);
128 }
129 },
130 bx_yhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
131 if (bc_ptr[n].hi(1) == REMORABCType::ext_dir) {
132 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][4] * mskr(i,j,0);
133 } else if (bc_ptr[n].hi(1) == REMORABCType::orlanski_rad) {
134 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);
135 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);
136 Real dTdt = calc_arr(i,dom_hi.y,k,icomp+n) - dest_arr(i,dom_hi.y ,k,icomp+n);
137 Real dTde = dest_arr(i,dom_hi.y,k,icomp+n) - dest_arr(i,dom_hi.y-1,k,icomp+n);
138 if (dTdt * dTde < 0.0_rt) dTdt = 0.0_rt;
139 Real dTdx = (dTdt * (grad_hi + grad_hi_ip1) > 0.0_rt) ? grad_hi : grad_hi_ip1;
140 Real cff = std::max(dTdx*dTdx + dTde*dTde, eps);
141 Real Ce = dTdt*dTde;
142 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);
143 }
144 }
145 );
146 }
147
148 {
149 Box bx_zlo(bx); bx_zlo.setBig (2,dom_lo.z-1);
150 Box bx_zhi(bx); bx_zhi.setSmall(2,dom_hi.z+1);
151 ParallelFor(
152 bx_zlo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
153 if (bc_ptr[n].lo(2) == REMORABCType::ext_dir) {
154 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][2] * mskr(i,j,0);
155 }
156 },
157 bx_zhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
158 if (bc_ptr[n].hi(2) == REMORABCType::ext_dir) {
159 dest_arr(i,j,k,icomp+n) = l_bc_extdir_vals_d[n][5] * mskr(i,j,0);
160 }
161 }
162 );
163 }
164
165 Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1-n_not_fill);
166 bx_xlo.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
167 bx_xlo.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
168 Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1+n_not_fill);
169 bx_xhi.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
170 bx_xhi.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
171 Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1-n_not_fill);
172 bx_ylo.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
173 bx_ylo.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
174 Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1+n_not_fill);
175 bx_yhi.setSmall(2,std::max(dom_lo.z,bx.smallEnd(2)));
176 bx_yhi.setBig (2,std::min(dom_hi.z,bx.bigEnd(2)));
177 // Calculate intersections for corners before adjusting to exclude them
178 Box xlo_ylo = bx_xlo & bx_ylo;
179 Box xhi_ylo = bx_xhi & bx_ylo;
180 Box xlo_yhi = bx_xlo & bx_yhi;
181 Box xhi_yhi = bx_xhi & bx_yhi;
182// bx_xlo.setSmall(1,valid_bx.smallEnd(1)); bx_xlo.setBig(1,valid_bx.bigEnd(1));
183// bx_xhi.setSmall(1,valid_bx.smallEnd(1)); bx_xhi.setBig(1,valid_bx.bigEnd(1));
184// bx_ylo.setSmall(0,valid_bx.smallEnd(0)); bx_ylo.setBig(0,valid_bx.bigEnd(0));
185// bx_yhi.setSmall(0,valid_bx.smallEnd(0)); bx_yhi.setBig(0,valid_bx.bigEnd(0));
186 // Next do ghost cells in x-direction but not reaching out in y
187 // The corners we miss here will be covered in the y-loop below or by periodicity
188 if (!is_periodic_in_x or bccomp==BCVars::foextrap_bc)
189 {
190 // Populate ghost cells on lo-x and hi-x domain boundaries
191 ParallelFor(bx_xlo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
192 int iflip = dom_lo.x - 1 - i;
193 int inner = (bc_ptr[n].lo(0) == REMORABCType::orlanski_rad) ? 1 : 0;
194 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 ||
195 bc_ptr[n].lo(0) == REMORABCType::orlanski_rad_nudge) {
196 dest_arr(i,j,k,icomp+n) = dest_arr(dom_lo.x-n_not_fill-inner,j,k,icomp+n);
197 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_even) {
198 dest_arr(i,j,k,icomp+n) = dest_arr(iflip,j,k,icomp+n);
199 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_odd) {
200 dest_arr(i,j,k,icomp+n) = -dest_arr(iflip,j,k,icomp+n);
201 }
202 },
203 bx_xhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
204 int iflip = 2*dom_hi.x + 1 - i;
205 int inner = (bc_ptr[n].hi(0) == REMORABCType::orlanski_rad) ? 1 : 0;
206 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 ||
207 bc_ptr[n].hi(0) == REMORABCType::orlanski_rad_nudge) {
208 dest_arr(i,j,k,icomp+n) = dest_arr(dom_hi.x+n_not_fill+inner,j,k,icomp+n);
209 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_even) {
210 dest_arr(i,j,k,icomp+n) = dest_arr(iflip,j,k,icomp+n);
211 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_odd) {
212 dest_arr(i,j,k,icomp+n) = -dest_arr(iflip,j,k,icomp+n);
213 }
214 }
215 );
216 }
217
218 if (!is_periodic_in_y or bccomp==BCVars::foextrap_bc)
219 {
220 // Populate ghost cells on lo-y and hi-y domain boundaries
221 ParallelFor(
222 bx_ylo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
223 int jflip = dom_lo.y - 1 - j;
224 int inner = (bc_ptr[n].lo(1) == REMORABCType::orlanski_rad) ? 1 : 0;
225 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 ||
226 bc_ptr[n].lo(1) == REMORABCType::orlanski_rad_nudge) {
227 dest_arr(i,j,k,icomp+n) = dest_arr(i,dom_lo.y-n_not_fill-inner,k,icomp+n);
228 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_even) {
229 dest_arr(i,j,k,icomp+n) = dest_arr(i,jflip,k,icomp+n);
230 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_odd) {
231 dest_arr(i,j,k,icomp+n) = -dest_arr(i,jflip,k,icomp+n);
232 }
233 },
234 bx_yhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
235 int jflip = 2*dom_hi.y + 1 - j;
236 int inner = (bc_ptr[n].hi(1) == REMORABCType::orlanski_rad) ? 1 : 0;
237 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 ||
238 bc_ptr[n].hi(1) == REMORABCType::orlanski_rad_nudge) {
239 dest_arr(i,j,k,icomp+n) = dest_arr(i,dom_hi.y+n_not_fill+inner,k,icomp+n);
240 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_even) {
241 dest_arr(i,j,k,icomp+n) = dest_arr(i,jflip,k,icomp+n);
242 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_odd) {
243 dest_arr(i,j,k,icomp+n) = -dest_arr(i,jflip,k,icomp+n);
244 }
245 }
246 );
247 }
248
249 {
250 Box bx_zlo(bx); bx_zlo.setBig (2,std::max(dom_lo.z-1,bx.smallEnd(2)));
251 Box bx_zhi(bx); bx_zhi.setSmall(2,std::min(dom_hi.z+1,bx.bigEnd(2)));
252 // Populate ghost cells on lo-z and hi-z domain boundaries
253
254 if (bx_zlo.ok()) {
255 ParallelFor(bx_zlo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
256 {
257 int kflip = dom_lo.z - 1 - i;
258 if (bc_ptr[n].lo(2) == REMORABCType::foextrap) {
259 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,dom_lo.z,icomp+n);
260 } else if (bc_ptr[n].lo(2) == REMORABCType::reflect_even) {
261 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,kflip,icomp+n);
262 } else if (bc_ptr[n].lo(2) == REMORABCType::reflect_odd) {
263 dest_arr(i,j,k,icomp+n) = -dest_arr(i,j,kflip,icomp+n);
264 }
265 });
266 }
267
268 if (bx_zlo.ok()) {
269 ParallelFor(bx_zhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
270 {
271 int kflip = 2*dom_hi.z + 1 - i;
272 if (bc_ptr[n].hi(2) == REMORABCType::foextrap) {
273 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,dom_hi.z,icomp+n);
274 } else if (bc_ptr[n].hi(2) == REMORABCType::reflect_even) {
275 dest_arr(i,j,k,icomp+n) = dest_arr(i,j,kflip,icomp+n);
276 } else if (bc_ptr[n].hi(2) == REMORABCType::reflect_odd) {
277 dest_arr(i,j,k,icomp+n) = -dest_arr(i,j,kflip,icomp+n);
278 }
279 });
280 }
281 }
282 if ((!is_periodic_in_x && !is_periodic_in_y) or bccomp==BCVars::foextrap_bc) {
283 // If we've applied boundary conditions to either side, update the corner
284 if (!xlo_ylo.isEmpty()) {
285 ParallelFor(xlo_ylo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
286 {
287 if (!(bc_ptr[n].lo(0) == REMORABCType::clamped || bc_ptr[n].lo(0) == REMORABCType::flather ||
288 bc_ptr[n].lo(0) == REMORABCType::chapman || bc_ptr[n].lo(0) == REMORABCType::orlanski_rad_nudge)
289 && !(bc_ptr[n].lo(1) == REMORABCType::clamped || bc_ptr[n].lo(1) == REMORABCType::flather ||
290 bc_ptr[n].lo(1) == REMORABCType::chapman || bc_ptr[n].lo(1) == REMORABCType::orlanski_rad_nudge)) {
291 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_lo.y,k,icomp+n)
292 + dest_arr(dom_lo.x,j,k,icomp+n));
293 }
294 });
295 }
296 if (!xlo_yhi.isEmpty()) {
297 ParallelFor(xlo_yhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
298 {
299 if (!(bc_ptr[n].lo(0) == REMORABCType::clamped || bc_ptr[n].lo(0) == REMORABCType::flather ||
300 bc_ptr[n].lo(0) == REMORABCType::chapman || bc_ptr[n].lo(0) == REMORABCType::orlanski_rad_nudge)
301 && !(bc_ptr[n].hi(1) == REMORABCType::clamped || bc_ptr[n].hi(1) == REMORABCType::flather ||
302 bc_ptr[n].hi(1) == REMORABCType::chapman || bc_ptr[n].hi(1) == REMORABCType::orlanski_rad_nudge)) {
303 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_hi.y,k,icomp+n)
304 + dest_arr(dom_lo.x,j,k,icomp+n));
305 }
306 });
307 }
308 if (!xhi_ylo.isEmpty()) {
309 ParallelFor(xhi_ylo & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
310 {
311 if (!(bc_ptr[n].hi(0) == REMORABCType::clamped || bc_ptr[n].hi(0) == REMORABCType::flather ||
312 bc_ptr[n].hi(0) == REMORABCType::chapman || bc_ptr[n].hi(0) == REMORABCType::orlanski_rad_nudge)
313 && !(bc_ptr[n].lo(1) == REMORABCType::clamped || bc_ptr[n].lo(1) == REMORABCType::flather ||
314 bc_ptr[n].lo(1) == REMORABCType::chapman || bc_ptr[n].lo(1) == REMORABCType::orlanski_rad_nudge)) {
315 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_lo.y,k,icomp+n)
316 + dest_arr(dom_hi.x,j,k,icomp+n));
317 }
318 });
319 }
320 if (!xhi_yhi.isEmpty()) {
321 ParallelFor(xhi_yhi & dest_arr_box, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
322 {
323 if (!(bc_ptr[n].hi(0) == REMORABCType::clamped || bc_ptr[n].hi(0) == REMORABCType::flather ||
324 bc_ptr[n].hi(0) == REMORABCType::chapman || bc_ptr[n].hi(0) == REMORABCType::orlanski_rad_nudge)
325 && !(bc_ptr[n].hi(1) == REMORABCType::clamped || bc_ptr[n].hi(1) == REMORABCType::flather ||
326 bc_ptr[n].hi(1) == REMORABCType::chapman || bc_ptr[n].hi(1) == REMORABCType::orlanski_rad_nudge)) {
327 dest_arr(i,j,k,icomp+n) = 0.5 * (dest_arr(i,dom_hi.y,k,icomp+n)
328 + dest_arr(dom_hi.x,j,k,icomp+n));
329 }
330 });
331 }
332 }
333
334 Gpu::streamSynchronize();
335}
#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