Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_MFInterp_2D_C.H
Go to the documentation of this file.
1#ifndef AMREX_MF_INTERP_2D_C_H_
2#define AMREX_MF_INTERP_2D_C_H_
3
4#include <AMReX_Geometry.H>
5
6namespace amrex {
7
9void mf_cell_cons_lin_interp_limit_minmax_llslope (int i, int j, int, Array4<Real> const& slope,
10 Array4<Real const> const& u, int scomp, int ncomp,
11 Box const& domain, IntVect const& ratio, BCRec const* bc) noexcept
12{
13 Real sfx = Real(1.0);
14 Real sfy = Real(1.0);
15
16 Real sx = Real(0.0);
17 Real sy = Real(0.0);
18
19 Real dcx = Real(0.0);
20 Real dcy = Real(0.0);
21
22 for (int ns = 0; ns < ncomp; ++ns) {
23 int nu = ns + scomp;
24
25 // x-direction
26 if (ratio[0] > 1) {
27 dcx = mf_compute_slopes_x(i, j, 0, u, nu, domain, bc[ns]);
28 Real df = Real(2.0) * (u(i+1,j,0,nu) - u(i ,j,0,nu));
29 Real db = Real(2.0) * (u(i ,j,0,nu) - u(i-1,j,0,nu));
30 sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
31 sx = std::copysign(Real(1.),dcx)*amrex::min(sx,std::abs(dcx));
32 slope(i,j,0,ns ) = dcx;
33 } else {
34 slope(i,j,0,ns ) = Real(0.0);
35 }
36
37 // y-direction
38 if (ratio[1] > 1) {
39 dcy = mf_compute_slopes_y(i, j, 0, u, nu, domain, bc[ns]);
40 Real df = Real(2.0) * (u(i,j+1,0,nu) - u(i,j ,0,nu));
41 Real db = Real(2.0) * (u(i,j ,0,nu) - u(i,j-1,0,nu));
42 sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
43 sy = std::copysign(Real(1.),dcy)*amrex::min(sy,std::abs(dcy));
44 slope(i,j,0,ns+ ncomp) = dcy;
45 } else {
46 slope(i,j,0,ns+ ncomp) = Real(0.0);
47 }
48
49 // adjust limited slopes to prevent new min/max for this component
50 Real alpha = Real(1.0);
51 if (sx != Real(0.0) || sy != Real(0.0)) {
52 Real dumax = std::abs(sx) * Real(ratio[0]-1)/Real(2*ratio[0])
53 + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1]);
54 Real umax = u(i,j,0,nu);
55 Real umin = u(i,j,0,nu);
56 int ilim = ratio[0] > 1 ? 1 : 0;
57 int jlim = ratio[1] > 1 ? 1 : 0;
58 for (int joff = -jlim; joff <= jlim; ++joff) {
59 for (int ioff = -ilim; ioff <= ilim; ++ioff) {
60 umin = amrex::min(umin, u(i+ioff,j+joff,0,nu));
61 umax = amrex::max(umax, u(i+ioff,j+joff,0,nu));
62 }}
63 if (dumax * alpha > (umax - u(i,j,0,nu))) {
64 alpha = (umax - u(i,j,0,nu)) / dumax;
65 }
66 if (dumax * alpha > (u(i,j,0,nu) - umin)) {
67 alpha = (u(i,j,0,nu) - umin) / dumax;
68 }
69 }
70 sx *= alpha;
71 sy *= alpha;
72
73 // for each direction, compute minimum of the ratio of limited to unlimited slopes
74 if (dcx != Real(0.0)) {
75 sfx = amrex::min(sfx, sx / dcx);
76 }
77 if (dcy != Real(0.0)) {
78 sfy = amrex::min(sfy, sy / dcy);
79 }
80 }
81
82 // multiply unlimited slopes by the minimum of the ratio of limited to unlimited slopes
83 for (int ns = 0; ns < ncomp; ++ns) {
84 slope(i,j,0,ns ) *= sfx;
85 slope(i,j,0,ns+ ncomp) *= sfy;
86 }
87}
88
90void mf_cell_cons_lin_interp_llslope (int i, int j, int, Array4<Real> const& slope,
91 Array4<Real const> const& u, int scomp, int ncomp,
92 Box const& domain, IntVect const& ratio, BCRec const* bc) noexcept
93{
94 Real sfx = Real(1.0);
95 Real sfy = Real(1.0);
96
97 for (int ns = 0; ns < ncomp; ++ns) {
98 int nu = ns + scomp;
99
100 // x-direction
101 if (ratio[0] > 1) {
102 Real dc = mf_compute_slopes_x(i, j, 0, u, nu, domain, bc[ns]);
103 Real df = Real(2.0) * (u(i+1,j,0,nu) - u(i ,j,0,nu));
104 Real db = Real(2.0) * (u(i ,j,0,nu) - u(i-1,j,0,nu));
105 Real sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
106 sx = std::copysign(Real(1.),dc)*amrex::min(sx,std::abs(dc));
107 if (dc != Real(0.0)) {
108 sfx = amrex::min(sfx, sx / dc);
109 }
110 slope(i,j,0,ns ) = dc;
111 } else {
112 slope(i,j,0,ns ) = Real(0.0);
113 }
114
115 // y-direction
116 if (ratio[1] > 1) {
117 Real dc = mf_compute_slopes_y(i, j, 0, u, nu, domain, bc[ns]);
118 Real df = Real(2.0) * (u(i,j+1,0,nu) - u(i,j ,0,nu));
119 Real db = Real(2.0) * (u(i,j ,0,nu) - u(i,j-1,0,nu));
120 Real sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
121 sy = std::copysign(Real(1.),dc)*amrex::min(sy,std::abs(dc));
122 if (dc != Real(0.0)) {
123 sfy = amrex::min(sfy, sy / dc);
124 }
125 slope(i,j,0,ns+ ncomp) = dc;
126 } else {
127 slope(i,j,0,ns+ ncomp) = Real(0.0);
128 }
129 }
130
131 for (int ns = 0; ns < ncomp; ++ns) {
132 slope(i,j,0,ns ) *= sfx;
133 slope(i,j,0,ns+ ncomp) *= sfy;
134 }
135}
136
138void mf_cell_cons_lin_interp_mcslope (int i, int j, int /*k*/, int ns, Array4<Real> const& slope,
139 Array4<Real const> const& u, int scomp, int ncomp,
140 Box const& domain, IntVect const& ratio,
141 BCRec const* bc) noexcept
142{
143 int nu = ns + scomp;
144
145 Real sx = Real(0.);
146 Real sy = Real(0.);
147
148 // x-direction
149 if (ratio[0] > 1) {
150 Real dc = mf_compute_slopes_x(i, j, 0, u, nu, domain, bc[ns]);
151 Real df = Real(2.0) * (u(i+1,j,0,nu) - u(i ,j,0,nu));
152 Real db = Real(2.0) * (u(i ,j,0,nu) - u(i-1,j,0,nu));
153 sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
154 sx = std::copysign(Real(1.),dc)*amrex::min(sx,std::abs(dc));
155 }
156
157 // y-direction
158 if (ratio[1] > 1) {
159 Real dc = mf_compute_slopes_y(i, j, 0, u, nu, domain, bc[ns]);
160 Real df = Real(2.0) * (u(i,j+1,0,nu) - u(i,j ,0,nu));
161 Real db = Real(2.0) * (u(i,j ,0,nu) - u(i,j-1,0,nu));
162 sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
163 sy = std::copysign(Real(1.),dc)*amrex::min(sy,std::abs(dc));
164 }
165
166 Real alpha = Real(1.0);
167 if (sx != Real(0.0) || sy != Real(0.0)) {
168 Real dumax = std::abs(sx) * Real(ratio[0]-1)/Real(2*ratio[0])
169 + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1]);
170 Real umax = u(i,j,0,nu);
171 Real umin = u(i,j,0,nu);
172 int ilim = ratio[0] > 1 ? 1 : 0;
173 int jlim = ratio[1] > 1 ? 1 : 0;
174 for (int joff = -jlim; joff <= jlim; ++joff) {
175 for (int ioff = -ilim; ioff <= ilim; ++ioff) {
176 umin = amrex::min(umin, u(i+ioff,j+joff,0,nu));
177 umax = amrex::max(umax, u(i+ioff,j+joff,0,nu));
178 }}
179 if (dumax * alpha > (umax - u(i,j,0,nu))) {
180 alpha = (umax - u(i,j,0,nu)) / dumax;
181 }
182 if (dumax * alpha > (u(i,j,0,nu) - umin)) {
183 alpha = (u(i,j,0,nu) - umin) / dumax;
184 }
185 }
186
187 slope(i,j,0,ns ) = sx * alpha;
188 slope(i,j,0,ns+ ncomp) = sy * alpha;
189}
190
192void mf_cell_cons_lin_interp (int i, int j, int /*k*/, int ns, Array4<Real> const& fine, int fcomp,
193 Array4<Real const> const& slope, Array4<Real const> const& crse,
194 int ccomp, int ncomp, IntVect const& ratio) noexcept
195{
196 const int ic = amrex::coarsen(i, ratio[0]);
197 const int jc = amrex::coarsen(j, ratio[1]);
198 const Real xoff = (static_cast<Real>(i - ic*ratio[0]) + Real(0.5)) / Real(ratio[0]) - Real(0.5);
199 const Real yoff = (static_cast<Real>(j - jc*ratio[1]) + Real(0.5)) / Real(ratio[1]) - Real(0.5);
200 fine(i,j,0,fcomp+ns) = crse(ic,jc,0,ccomp+ns)
201 + xoff * slope(ic,jc,0,ns)
202 + yoff * slope(ic,jc,0,ns+ncomp);
203}
204
206void mf_cell_cons_lin_interp_mcslope_sph (int i, int j, int ns, Array4<Real> const& slope,
207 Array4<Real const> const& u, int scomp, int ncomp,
208 Box const& domain, IntVect const& ratio,
209 BCRec const* bc, Real drf, Real rlo,
210 Real dtf, Real tlo) noexcept
211{
212 int nu = ns + scomp;
213
214 Real sx = Real(0.0);
215 Real sy = Real(0.0);
216
217 // x-direction
218 if (ratio[0] > 1) {
219 Real dc = mf_compute_slopes_x(i, j, 0, u, nu, domain, bc[ns]);
220 Real df = Real(2.0) * (u(i+1,j,0,nu) - u(i ,j,0,nu));
221 Real db = Real(2.0) * (u(i ,j,0,nu) - u(i-1,j,0,nu));
222 sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
223 sx = std::copysign(Real(1.),dc)*amrex::min(sx,std::abs(dc));
224 }
225
226 // y-direction
227 if (ratio[1] > 1) {
228 Real dc = mf_compute_slopes_y(i, j, 0, u, nu, domain, bc[ns]);
229 Real df = Real(2.0) * (u(i,j+1,0,nu) - u(i,j ,0,nu));
230 Real db = Real(2.0) * (u(i,j ,0,nu) - u(i,j-1,0,nu));
231 sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
232 sy = std::copysign(Real(1.),dc)*amrex::min(sy,std::abs(dc));
233 }
234
235 Real alpha = Real(1.0);
236 if (sx != Real(0.0) || sy != Real(0.0)) {
237
238 // x-direction
239 const Real drc = drf * ratio[0];
240 const Real rcm = i * drc + rlo;
241 const Real rcp = (i+1) * drc + rlo;
242 Real vcm = rcm*rcm*rcm;
243 Real vcp = rcp*rcp*rcp;
244 Real rfm = i*ratio[0] * drf + rlo;
245 Real rfp = (i*ratio[0] + 1) * drf + rlo;
246 Real vfm = rfm*rfm*rfm;
247 Real vfp = rfp*rfp*rfp;
248 Real xlo = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
249 rfm = ((i+1)*ratio[0] - 1) * drf + rlo;
250 rfp = (i+1)*ratio[0] * drf + rlo;
251 vfm = rfm*rfm*rfm;
252 vfp = rfp*rfp*rfp;
253 Real xhi = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
254
255 // y-direction
256 const Real dtc = dtf * ratio[1];
257 const Real tcm = j * dtc + tlo;
258 const Real tcp = (j+1) * dtc + tlo;
259 vcm = std::cos(tcm);
260 vcp = std::cos(tcp);
261 Real tfm = j*ratio[1] * dtf + tlo;
262 Real tfp = (j*ratio[1] + 1) * dtf + tlo;
263 vfm = std::cos(tfm);
264 vfp = std::cos(tfp);
265 Real ylo = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
266 tfm = ((j+1)*ratio[1] - 1) * dtf + tlo;
267 tfp = (j+1)*ratio[1] * dtf + tlo;
268 vfm = std::cos(tfm);
269 vfp = std::cos(tfp);
270 Real yhi = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
271
272 Real dumax = amrex::max(sx*xlo, sx*xhi) + amrex::max(sy*ylo, sy*yhi);
273 Real dumin = -amrex::min(sx*xlo, sx*xhi) - amrex::min(sy*ylo, sy*yhi);;
274 Real umax = u(i,j,0,nu);
275 Real umin = u(i,j,0,nu);
276 int ilim = ratio[0] > 1 ? 1 : 0;
277 int jlim = ratio[1] > 1 ? 1 : 0;
278 for (int joff = -jlim; joff <= jlim; ++joff) {
279 for (int ioff = -ilim; ioff <= ilim; ++ioff) {
280 umin = amrex::min(umin, u(i+ioff,j+joff,0,nu));
281 umax = amrex::max(umax, u(i+ioff,j+joff,0,nu));
282 }}
283 if (dumax * alpha > (umax - u(i,j,0,nu))) {
284 alpha = (umax - u(i,j,0,nu)) / dumax;
285 }
286 if (dumin * alpha > (u(i,j,0,nu) - umin)) {
287 alpha = (u(i,j,0,nu) - umin) / dumin;
288 }
289 }
290
291 slope(i,j,0,ns ) = sx * alpha;
292 slope(i,j,0,ns+ ncomp) = sy * alpha;
293}
294
296void mf_cell_cons_lin_interp_sph (int i, int j, int ns, Array4<Real> const& fine, int fcomp,
298 int ccomp, int ncomp, IntVect const& ratio, Real drf, Real rlo,
299 Real dtf, Real tlo) noexcept
300{
301 // x-direction
302 const int ic = amrex::coarsen(i, ratio[0]);
303 const Real drc = drf * ratio[0];
304 const Real rcm = ic * drc + rlo;
305 const Real rcp = (ic+1) * drc + rlo;
306 const Real rfm = i * drf + rlo;
307 const Real rfp = (i +1) * drf + rlo;
308 Real vcm = rcm*rcm*rcm;
309 Real vcp = rcp*rcp*rcp;
310 Real vfm = rfm*rfm*rfm;
311 Real vfp = rfp*rfp*rfp;
312 const Real xoff = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
313
314 // y-direction
315 const int jc = amrex::coarsen(j, ratio[1]);
316 const Real dtc = dtf * ratio[1];
317 const Real tcm = jc * dtc + tlo;
318 const Real tcp = (jc+1) * dtc + tlo;
319 const Real tfm = j * dtf + tlo;
320 const Real tfp = (j +1) * dtf + tlo;
321 vcm = std::cos(tcm);
322 vcp = std::cos(tcp);
323 vfm = std::cos(tfm);
324 vfp = std::cos(tfp);
325 const Real yoff = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
326
327 fine(i,j,0,fcomp+ns) = crse(ic,jc,0,ccomp+ns)
328 + xoff * slope(ic,jc,0,ns)
329 + yoff * slope(ic,jc,0,ns+ncomp);
330}
331
333void mf_cell_cons_lin_interp_mcslope_rz (int i, int j, int ns, Array4<Real> const& slope,
334 Array4<Real const> const& u, int scomp, int ncomp,
335 Box const& domain, IntVect const& ratio,
336 BCRec const* bc, Real drf, Real rlo) noexcept
337{
338 int nu = ns + scomp;
339
340 Real sx = Real(0.0);
341 Real sy = Real(0.0);
342
343 // x-direction
344 if (ratio[0] > 1) {
345 Real dc = mf_compute_slopes_x(i, j, 0, u, nu, domain, bc[ns]);
346 Real df = Real(2.0) * (u(i+1,j,0,nu) - u(i ,j,0,nu));
347 Real db = Real(2.0) * (u(i ,j,0,nu) - u(i-1,j,0,nu));
348 sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
349 sx = std::copysign(Real(1.),dc)*amrex::min(sx,std::abs(dc));
350 }
351
352 // y-direction
353 if (ratio[1] > 1) {
354 Real dc = mf_compute_slopes_y(i, j, 0, u, nu, domain, bc[ns]);
355 Real df = Real(2.0) * (u(i,j+1,0,nu) - u(i,j ,0,nu));
356 Real db = Real(2.0) * (u(i,j ,0,nu) - u(i,j-1,0,nu));
357 sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
358 sy = std::copysign(Real(1.),dc)*amrex::min(sy,std::abs(dc));
359 }
360
361 Real alpha = Real(1.0);
362 if (sx != Real(0.0) || sy != Real(0.0)) {
363 const Real drc = drf * ratio[0];
364 const Real rcm = i * drc + rlo;
365 const Real rcp = (i+1) * drc + rlo;
366 const Real vcm = rcm*rcm;
367 const Real vcp = rcp*rcp;
368 Real rfm = i*ratio[0] * drf + rlo;
369 Real rfp = (i*ratio[0] + 1) * drf + rlo;
370 Real vfm = rfm*rfm;
371 Real vfp = rfp*rfp;
372 Real xlo = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
373 rfm = ((i+1)*ratio[0] - 1) * drf + rlo;
374 rfp = (i+1)*ratio[0] * drf + rlo;
375 vfm = rfm*rfm;
376 vfp = rfp*rfp;
377 Real xhi = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
378 Real dumax = amrex::max(sx*xlo, sx*xhi)
379 + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1]);
380 Real dumin = -amrex::min(sx*xlo, sx*xhi)
381 + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1]);
382 Real umax = u(i,j,0,nu);
383 Real umin = u(i,j,0,nu);
384 int ilim = ratio[0] > 1 ? 1 : 0;
385 int jlim = ratio[1] > 1 ? 1 : 0;
386 for (int joff = -jlim; joff <= jlim; ++joff) {
387 for (int ioff = -ilim; ioff <= ilim; ++ioff) {
388 umin = amrex::min(umin, u(i+ioff,j+joff,0,nu));
389 umax = amrex::max(umax, u(i+ioff,j+joff,0,nu));
390 }}
391 if (dumax * alpha > (umax - u(i,j,0,nu))) {
392 alpha = (umax - u(i,j,0,nu)) / dumax;
393 }
394 if (dumin * alpha > (u(i,j,0,nu) - umin)) {
395 alpha = (u(i,j,0,nu) - umin) / dumin;
396 }
397 }
398
399 slope(i,j,0,ns ) = sx * alpha;
400 slope(i,j,0,ns+ ncomp) = sy * alpha;
401}
402
404void mf_cell_cons_lin_interp_rz (int i, int j, int ns, Array4<Real> const& fine, int fcomp,
406 int ccomp, int ncomp, IntVect const& ratio, Real drf, Real rlo) noexcept
407{
408 const int ic = amrex::coarsen(i, ratio[0]);
409 const int jc = amrex::coarsen(j, ratio[1]);
410 const Real drc = drf * ratio[0];
411 const Real rcm = ic * drc + rlo;
412 const Real rcp = (ic+1) * drc + rlo;
413 const Real rfm = i * drf + rlo;
414 const Real rfp = (i +1) * drf + rlo;
415 const Real vcm = rcm*rcm;
416 const Real vcp = rcp*rcp;
417 const Real vfm = rfm*rfm;
418 const Real vfp = rfp*rfp;
419 const Real xoff = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
420 const Real yoff = (j - jc*ratio[1] + Real(0.5)) / Real(ratio[1]) - Real(0.5);
421 fine(i,j,0,fcomp+ns) = crse(ic,jc,0,ccomp+ns)
422 + xoff * slope(ic,jc,0,ns)
423 + yoff * slope(ic,jc,0,ns+ncomp);
424}
425
426template<typename T>
428void mf_cell_bilin_interp (int i, int j, int, int n, Array4<T> const& fine, int fcomp,
429 Array4<T const> const& crse, int ccomp, IntVect const& ratio) noexcept
430{
431 int ic = amrex::coarsen(i,ratio[0]);
432 int jc = amrex::coarsen(j,ratio[1]);
433 int ioff = i - ic*ratio[0];
434 int joff = j - jc*ratio[1];
435 int sx, sy;
436 Real wx, wy;
437 if (ioff*2 < ratio[0]) {
438 sx = -1;
439 wx = Real(ratio[0]+1+2*ioff) / Real(2*ratio[0]);
440 } else {
441 sx = 1;
442 wx = Real(3*ratio[0]-1-2*ioff) / Real(2*ratio[0]);
443 }
444 if (joff*2 < ratio[1]) {
445 sy = -1;
446 wy = Real(ratio[1]+1+2*joff) / Real(2*ratio[1]);
447 } else {
448 sy = 1;
449 wy = Real(3*ratio[1]-1-2*joff) / Real(2*ratio[1]);
450 }
451 fine(i,j,0,n+fcomp) =
452 crse(ic ,jc ,0,n+ccomp)* wx * wy +
453 crse(ic+sx,jc ,0,n+ccomp)*(Real(1.0)-wx)* wy +
454 crse(ic ,jc+sy,0,n+ccomp)* wx *(Real(1.0)-wy) +
455 crse(ic+sx,jc+sy,0,n+ccomp)*(Real(1.0)-wx)*(Real(1.0)-wy);
456}
457
459void mf_nodebilin_interp (int i, int j, int, int n, Array4<Real> const& fine, int fcomp,
460 Array4<Real const> const& crse, int ccomp, IntVect const& ratio) noexcept
461{
462 int ic = amrex::coarsen(i,ratio[0]);
463 int jc = amrex::coarsen(j,ratio[1]);
464 int ioff = i - ic*ratio[0];
465 int joff = j - jc*ratio[1];
466 Real rxinv = Real(1.0) / Real(ratio[0]);
467 Real ryinv = Real(1.0) / Real(ratio[1]);
468 if (ioff != 0 && joff != 0) {
469 // Node on a X-Y face
470 fine(i,j,0,n+fcomp) = rxinv * ryinv *
471 (crse(ic ,jc ,0,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * (ratio[1]-joff)) +
472 crse(ic+1,jc ,0,n+ccomp) * static_cast<Real>(( ioff) * (ratio[1]-joff)) +
473 crse(ic ,jc+1,0,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * ( joff)) +
474 crse(ic+1,jc+1,0,n+ccomp) * static_cast<Real>(( ioff) * ( joff)));
475 } else if (ioff != 0) {
476 // Node on X line
477 fine(i,j,0,n+fcomp) = rxinv*(static_cast<Real>(ratio[0]-ioff)*crse(ic ,jc,0,n+ccomp) +
478 static_cast<Real>( ioff)*crse(ic+1,jc,0,n+ccomp));
479 } else if (joff != 0) {
480 // Node on Y line
481 fine(i,j,0,n+fcomp) = ryinv*(static_cast<Real>(ratio[1]-joff)*crse(ic,jc ,0,n+ccomp) +
482 static_cast<Real>( joff)*crse(ic,jc+1,0,n+ccomp));
483 } else {
484 // Node coincident with coarse node
485 fine(i,j,0,n+fcomp) = crse(ic,jc,0,n+ccomp);
486 }
487}
488
490void mf_cell_quadratic_calcslope (int i, int j, int /*k*/, int n,
491 Array4<Real const> const& crse, int ccomp,
492 Array4<Real> const& slope,
493 Box const& domain,
494 BCRec const* bc) noexcept
495{
496 int nu = ccomp + n;
497 Real sx = mf_compute_slopes_x(i, j, 0, crse, nu, domain, bc[n]);
498 Real sy = mf_compute_slopes_y(i, j, 0, crse, nu, domain, bc[n]);
499 Real sxx = mf_cell_quadratic_compute_slopes_xx (i, j, 0, crse, nu, domain, bc[n]);
500 Real syy = mf_cell_quadratic_compute_slopes_yy (i, j, 0, crse, nu, domain, bc[n]);
501 Real sxy = mf_cell_quadratic_compute_slopes_xy (i, j, 0, crse, nu, domain, bc[n]);
502
503 slope(i,j,0,5*n ) = sx; // x
504 slope(i,j,0,5*n+1) = sy; // y
505 slope(i,j,0,5*n+2) = sxx; // x^2
506 slope(i,j,0,5*n+3) = syy; // y^2
507 slope(i,j,0,5*n+4) = sxy; // xy
508}
509
511void mf_cell_quadratic_interp (int i, int j, int /*k*/, int n,
512 Array4<Real> const& fine, int fcomp,
513 Array4<Real const> const& crse, int ccomp,
515 IntVect const& ratio) noexcept
516{
517 int ic = amrex::coarsen(i, ratio[0]);
518 int jc = amrex::coarsen(j, ratio[1]);
519 int irx = i - ic*ratio[0]; // = abs(i % ratio[0])
520 int jry = j - jc*ratio[1]; // = abs(j % ratio[1])
521
522 // Compute offsets.
523 Real xoff = ( Real(irx) + 0.5_rt ) / Real(ratio[0]) - 0.5_rt;
524 Real yoff = ( Real(jry) + 0.5_rt ) / Real(ratio[1]) - 0.5_rt;
525
526 fine(i,j,0,fcomp+n) = crse(ic,jc,0,ccomp+n)
527 + xoff * slope(ic,jc,0,5*n ) // x
528 + yoff * slope(ic,jc,0,5*n+1) // y
529 + 0.5_rt * xoff * xoff * slope(ic,jc,0,5*n+2) // x^2
530 + 0.5_rt * yoff * yoff * slope(ic,jc,0,5*n+3) // y^2
531 + xoff * yoff * slope(ic,jc,0,5*n+4); // xy
532}
533
535void mf_cell_quadratic_interp_rz (int i, int j, int /*k*/, int n,
536 Array4<Real> const& fine, int fcomp,
537 Array4<Real const> const& crse, int ccomp,
539 IntVect const& ratio,
540 GeometryData const& cs_geomdata,
541 GeometryData const& fn_geomdata) noexcept
542{
543 int ic = amrex::coarsen(i, ratio[0]);
544 int jc = amrex::coarsen(j, ratio[1]);
545 int jry = j - jc*ratio[1]; // = abs(j % ratio[1])
546
547 // Extract geometry data for radial direction.
548 Real fn_dr = fn_geomdata.CellSize(0);
549 Real fn_rlo = fn_geomdata.ProbLo(0);
550 Real cs_dr = cs_geomdata.CellSize(0);
551 Real cs_rlo = cs_geomdata.ProbLo(0);
552
553 // Compute radial mesh.
554 Real fn_rm = fn_rlo + Real(i) * fn_dr;
555 Real fn_rp = fn_rlo + Real(i+1) * fn_dr;
556 Real cs_rm = cs_rlo + Real(ic) * cs_dr;
557 Real cs_rp = cs_rlo + Real(ic+1) * cs_dr;
558
559 // Compute radial cell edge coords.
560 // Note: omit the 0.5 prefactor here -- we'll add it once at the last step
561 Real fn_vm = fn_rm*fn_rm;
562 Real fn_vp = fn_rp*fn_rp;
563 Real fcen = fn_vm + fn_vp;
564 Real cs_vm = cs_rm*cs_rm;
565 Real cs_vp = cs_rp*cs_rp;
566 Real ccen = cs_vm + cs_vp;
567
568 // Compute r offset
569 Real xoff = 0.5_rt * ( fcen - ccen ) / ( cs_vp - cs_vm ) ;
570
571 // Compute z offset (same formula as Cartesian)
572 Real yoff = ( Real(jry) + 0.5_rt ) / Real(ratio[1]) - 0.5_rt;
573
574 fine(i,j,0,fcomp+n) = crse(ic,jc,0,ccomp+n)
575 + xoff * slope(ic,jc,0,5*n ) // x
576 + yoff * slope(ic,jc,0,5*n+1) // y
577 + 0.5_rt * xoff * xoff * slope(ic,jc,0,5*n+2) // x^2
578 + 0.5_rt * yoff * yoff * slope(ic,jc,0,5*n+3) // y^2
579 + xoff * yoff * slope(ic,jc,0,5*n+4); // xy
580}
581
582} // namespace amrex
583
584#endif
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real > slope
Definition AMReX_InterpFaceRegister.cpp:91
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
Boundary Condition Records. Necessary information and functions for computing boundary conditions.
Definition AMReX_BCRec.H:17
Definition AMReX_Amr.cpp:49
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_cell_cons_lin_interp_mcslope_sph(int i, int ns, Array4< Real > const &slope, Array4< Real const > const &u, int scomp, int, Box const &domain, IntVect const &ratio, BCRec const *bc, Real drf, Real rlo) noexcept
Definition AMReX_MFInterp_1D_C.H:114
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_cell_bilin_interp(int i, int, int, int n, Array4< T > const &fine, int fcomp, Array4< T const > const &crse, int ccomp, IntVect const &ratio) noexcept
Definition AMReX_MFInterp_1D_C.H:186
__host__ __device__ void mf_cell_quadratic_calcslope(int i, int j, int, int n, Array4< Real const > const &crse, int ccomp, Array4< Real > const &slope, Box const &domain, BCRec const *bc) noexcept
Definition AMReX_MFInterp_2D_C.H:490
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_cell_cons_lin_interp_llslope(int i, int, int, Array4< Real > const &slope, Array4< Real const > const &u, int scomp, int ncomp, Box const &domain, IntVect const &, BCRec const *bc) noexcept
Definition AMReX_MFInterp_1D_C.H:39
__host__ __device__ Real mf_cell_quadratic_compute_slopes_xx(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:89
__host__ __device__ Real mf_compute_slopes_y(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:37
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_cell_cons_lin_interp_limit_minmax_llslope(int i, int, int, Array4< Real > const &slope, Array4< Real const > const &u, int scomp, int ncomp, Box const &domain, IntVect const &ratio, BCRec const *bc) noexcept
Definition AMReX_MFInterp_1D_C.H:7
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:1322
BoxND< 3 > Box
Definition AMReX_BaseFwd.H:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_cell_cons_lin_interp_mcslope(int i, int, int, int ns, Array4< Real > const &slope, Array4< Real const > const &u, int scomp, int, Box const &domain, IntVect const &ratio, BCRec const *bc) noexcept
Definition AMReX_MFInterp_1D_C.H:66
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:21
__host__ __device__ void mf_cell_quadratic_interp_rz(int i, int j, int, int n, Array4< Real > const &fine, int fcomp, Array4< Real const > const &crse, int ccomp, Array4< Real const > const &slope, IntVect const &ratio, GeometryData const &cs_geomdata, GeometryData const &fn_geomdata) noexcept
Definition AMReX_MFInterp_2D_C.H:535
__host__ __device__ void mf_cell_cons_lin_interp_rz(int i, int j, int ns, Array4< Real > const &fine, int fcomp, Array4< Real const > const &slope, Array4< Real const > const &crse, int ccomp, int ncomp, IntVect const &ratio, Real drf, Real rlo) noexcept
Definition AMReX_MFInterp_2D_C.H:404
__host__ __device__ Real mf_cell_quadratic_compute_slopes_yy(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:110
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_cell_cons_lin_interp(int i, int, int, int ns, Array4< Real > const &fine, int fcomp, Array4< Real const > const &slope, Array4< Real const > const &crse, int ccomp, int, IntVect const &ratio) noexcept
Definition AMReX_MFInterp_1D_C.H:102
IntVectND< 3 > IntVect
Definition AMReX_BaseFwd.H:30
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:35
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_nodebilin_interp(int i, int, int, int n, Array4< Real > const &fine, int fcomp, Array4< Real const > const &crse, int ccomp, IntVect const &ratio) noexcept
Definition AMReX_MFInterp_1D_C.H:206
__host__ __device__ Real mf_compute_slopes_x(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:11
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_cell_cons_lin_interp_sph(int i, int ns, Array4< Real > const &fine, int fcomp, Array4< Real const > const &slope, Array4< Real const > const &crse, int ccomp, int, IntVect const &ratio, Real drf, Real rlo) noexcept
Definition AMReX_MFInterp_1D_C.H:165
__host__ __device__ void mf_cell_cons_lin_interp_mcslope_rz(int i, int j, int ns, Array4< Real > const &slope, Array4< Real const > const &u, int scomp, int ncomp, Box const &domain, IntVect const &ratio, BCRec const *bc, Real drf, Real rlo) noexcept
Definition AMReX_MFInterp_2D_C.H:333
__host__ __device__ Real mf_cell_quadratic_compute_slopes_xy(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:152
__host__ __device__ void mf_cell_quadratic_interp(int i, int j, int, int n, Array4< Real > const &fine, int fcomp, Array4< Real const > const &crse, int ccomp, Array4< Real const > const &slope, IntVect const &ratio) noexcept
Definition AMReX_MFInterp_2D_C.H:511
Definition AMReX_Array4.H:61
Definition AMReX_Geometry.H:25