Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_MFInterp_3D_C.H
Go to the documentation of this file.
1#ifndef AMREX_MF_INTERP_3D_C_H_
2#define AMREX_MF_INTERP_3D_C_H_
3
4namespace amrex {
5
7void mf_cell_cons_lin_interp_limit_minmax_llslope (int i, int j, int k, Array4<Real> const& slope,
8 Array4<Real const> const& u, int scomp, int ncomp,
9 Box const& domain, IntVect const& ratio, BCRec const* bc) noexcept
10{
11 Real sfx = Real(1.0);
12 Real sfy = Real(1.0);
13 Real sfz = Real(1.0);
14
15 Real sx = Real(0.0);
16 Real sy = Real(0.0);
17 Real sz = Real(0.0);
18
19 Real dcx = Real(0.0);
20 Real dcy = Real(0.0);
21 Real dcz = Real(0.0);
22
23 for (int ns = 0; ns < ncomp; ++ns) {
24 int nu = ns + scomp;
25
26 // x-direction
27 if (ratio[0] > 1) {
28 dcx = mf_compute_slopes_x(i, j, k, u, nu, domain, bc[ns]);
29 Real df = Real(2.0) * (u(i+1,j,k,nu) - u(i ,j,k,nu));
30 Real db = Real(2.0) * (u(i ,j,k,nu) - u(i-1,j,k,nu));
31 sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
32 sx = std::copysign(Real(1.),dcx)*amrex::min(sx,std::abs(dcx));
33 slope(i,j,k,ns ) = dcx; // unlimited slope
34 } else {
35 slope(i,j,k,ns ) = Real(0.0);
36 }
37
38 // y-direction
39 if (ratio[1] > 1) {
40 dcy = mf_compute_slopes_y(i, j, k, u, nu, domain, bc[ns]);
41 Real df = Real(2.0) * (u(i,j+1,k,nu) - u(i,j ,k,nu));
42 Real db = Real(2.0) * (u(i,j ,k,nu) - u(i,j-1,k,nu));
43 sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
44 sy = std::copysign(Real(1.),dcy)*amrex::min(sy,std::abs(dcy));
45 slope(i,j,k,ns+ ncomp) = dcy; // unlimited slope
46 } else {
47 slope(i,j,k,ns+ ncomp) = Real(0.0);
48 }
49
50 // z-direction
51 if (ratio[2] > 1) {
52 dcz = mf_compute_slopes_z(i, j, k, u, nu, domain, bc[ns]);
53 Real df = Real(2.0) * (u(i,j,k+1,nu) - u(i,j,k ,nu));
54 Real db = Real(2.0) * (u(i,j,k ,nu) - u(i,j,k-1,nu));
55 sz = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
56 sz = std::copysign(Real(1.),dcz)*amrex::min(sz,std::abs(dcz));
57 slope(i,j,k,ns+2*ncomp) = dcz; // unlimited slope
58 } else {
59 slope(i,j,k,ns+2*ncomp) = Real(0.0);
60 }
61
62 // adjust limited slopes to prevent new min/max for this component
63 Real alpha = 1.0;
64 if (sx != Real(0.0) || sy != Real(0.0) || sz != Real(0.0)) {
65 Real dumax = std::abs(sx) * Real(ratio[0]-1)/Real(2*ratio[0])
66 + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1])
67 + std::abs(sz) * Real(ratio[2]-1)/Real(2*ratio[2]);
68 Real umax = u(i,j,k,nu);
69 Real umin = u(i,j,k,nu);
70 int ilim = ratio[0] > 1 ? 1 : 0;
71 int jlim = ratio[1] > 1 ? 1 : 0;
72 int klim = ratio[2] > 1 ? 1 : 0;
73 for (int koff = -klim; koff <= klim; ++koff) {
74 for (int joff = -jlim; joff <= jlim; ++joff) {
75 for (int ioff = -ilim; ioff <= ilim; ++ioff) {
76 umin = amrex::min(umin, u(i+ioff,j+joff,k+koff,nu));
77 umax = amrex::max(umax, u(i+ioff,j+joff,k+koff,nu));
78 }}}
79 if (dumax * alpha > (umax - u(i,j,k,nu))) {
80 alpha = (umax - u(i,j,k,nu)) / dumax;
81 }
82 if (dumax * alpha > (u(i,j,k,nu) - umin)) {
83 alpha = (u(i,j,k,nu) - umin) / dumax;
84 }
85 }
86 sx *= alpha;
87 sy *= alpha;
88 sz *= alpha;
89
90 // for each direction, compute minimum of the ratio of limited to unlimited slopes
91 if (dcx != Real(0.0)) {
92 sfx = amrex::min(sfx, sx / dcx);
93 }
94 if (dcy != Real(0.0)) {
95 sfy = amrex::min(sfy, sy / dcy);
96 }
97 if (dcz != Real(0.0)) {
98 sfz = amrex::min(sfz, sz / dcz);
99 }
100 }
101
102 // multiply unlimited slopes by the minimum of the ratio of limited to unlimited slopes
103 for (int ns = 0; ns < ncomp; ++ns) {
104 slope(i,j,k,ns ) *= sfx;
105 slope(i,j,k,ns+ ncomp) *= sfy;
106 slope(i,j,k,ns+2*ncomp) *= sfz;
107 }
108}
109
111void mf_cell_cons_lin_interp_llslope (int i, int j, int k, Array4<Real> const& slope,
112 Array4<Real const> const& u, int scomp, int ncomp,
113 Box const& domain, IntVect const& ratio, BCRec const* bc) noexcept
114{
115 Real sfx = Real(1.0);
116 Real sfy = Real(1.0);
117 Real sfz = Real(1.0);
118
119 for (int ns = 0; ns < ncomp; ++ns) {
120 int nu = ns + scomp;
121
122 // x-direction
123 if (ratio[0] > 1) {
124 Real dc = mf_compute_slopes_x(i, j, k, u, nu, domain, bc[ns]);
125 Real df = Real(2.0) * (u(i+1,j,k,nu) - u(i ,j,k,nu));
126 Real db = Real(2.0) * (u(i ,j,k,nu) - u(i-1,j,k,nu));
127 Real sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
128 sx = std::copysign(Real(1.),dc)*amrex::min(sx,std::abs(dc));
129 if (dc != Real(0.0)) {
130 sfx = amrex::min(sfx, sx / dc);
131 }
132 slope(i,j,k,ns ) = dc;
133 } else {
134 slope(i,j,k,ns ) = Real(0.0);
135 }
136
137 // y-direction
138 if (ratio[1] > 1) {
139 Real dc = mf_compute_slopes_y(i, j, k, u, nu, domain, bc[ns]);
140 Real df = Real(2.0) * (u(i,j+1,k,nu) - u(i,j ,k,nu));
141 Real db = Real(2.0) * (u(i,j ,k,nu) - u(i,j-1,k,nu));
142 Real sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
143 sy = std::copysign(Real(1.),dc)*amrex::min(sy,std::abs(dc));
144 if (dc != Real(0.0)) {
145 sfy = amrex::min(sfy, sy / dc);
146 }
147 slope(i,j,k,ns+ ncomp) = dc;
148 } else {
149 slope(i,j,k,ns+ ncomp) = Real(0.0);
150 }
151
152
153 // z-direction
154 if (ratio[2] > 1) {
155 Real dc = mf_compute_slopes_z(i, j, k, u, nu, domain, bc[ns]);
156 Real df = Real(2.0) * (u(i,j,k+1,nu) - u(i,j,k ,nu));
157 Real db = Real(2.0) * (u(i,j,k ,nu) - u(i,j,k-1,nu));
158 Real sz = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
159 sz = std::copysign(Real(1.),dc)*amrex::min(sz,std::abs(dc));
160 if (dc != Real(0.0)) {
161 sfz = amrex::min(sfz, sz / dc);
162 }
163 slope(i,j,k,ns+2*ncomp) = dc;
164 } else {
165 slope(i,j,k,ns+2*ncomp) = Real(0.0);
166 }
167
168 }
169
170 for (int ns = 0; ns < ncomp; ++ns) {
171 slope(i,j,k,ns ) *= sfx;
172 slope(i,j,k,ns+ ncomp) *= sfy;
173 slope(i,j,k,ns+2*ncomp) *= sfz;
174 }
175}
176
178void mf_cell_cons_lin_interp_mcslope (int i, int j, int k, int ns, Array4<Real> const& slope,
179 Array4<Real const> const& u, int scomp, int ncomp,
180 Box const& domain, IntVect const& ratio,
181 BCRec const* bc) noexcept
182{
183 int nu = ns + scomp;
184
185 Real sx = Real(0.);
186 Real sy = Real(0.);
187 Real sz = Real(0.);
188
189 // x-direction
190 if (ratio[0] > 1) {
191 Real dc = mf_compute_slopes_x(i, j, k, u, nu, domain, bc[ns]);
192 Real df = Real(2.0) * (u(i+1,j,k,nu) - u(i ,j,k,nu));
193 Real db = Real(2.0) * (u(i ,j,k,nu) - u(i-1,j,k,nu));
194 sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
195 sx = std::copysign(Real(1.),dc)*amrex::min(sx,std::abs(dc));
196 }
197
198 // y-direction
199 if (ratio[1] > 1) {
200 Real dc = mf_compute_slopes_y(i, j, k, u, nu, domain, bc[ns]);
201 Real df = Real(2.0) * (u(i,j+1,k,nu) - u(i,j ,k,nu));
202 Real db = Real(2.0) * (u(i,j ,k,nu) - u(i,j-1,k,nu));
203 sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
204 sy = std::copysign(Real(1.),dc)*amrex::min(sy,std::abs(dc));
205 }
206
207 // z-direction
208 if (ratio[2] > 1) {
209 Real dc = mf_compute_slopes_z(i, j, k, u, nu, domain, bc[ns]);
210 Real df = Real(2.0) * (u(i,j,k+1,nu) - u(i,j,k ,nu));
211 Real db = Real(2.0) * (u(i,j,k ,nu) - u(i,j,k-1,nu));
212 sz = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
213 sz = std::copysign(Real(1.),dc)*amrex::min(sz,std::abs(dc));
214 }
215
216 Real alpha = 1.0;
217 if (sx != Real(0.0) || sy != Real(0.0) || sz != Real(0.0)) {
218 Real dumax = std::abs(sx) * Real(ratio[0]-1)/Real(2*ratio[0])
219 + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1])
220 + std::abs(sz) * Real(ratio[2]-1)/Real(2*ratio[2]);
221 Real umax = u(i,j,k,nu);
222 Real umin = u(i,j,k,nu);
223 int ilim = ratio[0] > 1 ? 1 : 0;
224 int jlim = ratio[1] > 1 ? 1 : 0;
225 int klim = ratio[2] > 1 ? 1 : 0;
226 for (int koff = -klim; koff <= klim; ++koff) {
227 for (int joff = -jlim; joff <= jlim; ++joff) {
228 for (int ioff = -ilim; ioff <= ilim; ++ioff) {
229 umin = amrex::min(umin, u(i+ioff,j+joff,k+koff,nu));
230 umax = amrex::max(umax, u(i+ioff,j+joff,k+koff,nu));
231 }}}
232 if (dumax * alpha > (umax - u(i,j,k,nu))) {
233 alpha = (umax - u(i,j,k,nu)) / dumax;
234 }
235 if (dumax * alpha > (u(i,j,k,nu) - umin)) {
236 alpha = (u(i,j,k,nu) - umin) / dumax;
237 }
238 }
239
240 slope(i,j,k,ns ) = sx * alpha;
241 slope(i,j,k,ns+ ncomp) = sy * alpha;
242 slope(i,j,k,ns+2*ncomp) = sz * alpha;
243}
244
246void mf_cell_cons_lin_interp (int i, int j, int k, int ns, Array4<Real> const& fine, int fcomp,
247 Array4<Real const> const& slope, Array4<Real const> const& crse,
248 int ccomp, int ncomp, IntVect const& ratio) noexcept
249{
250 const int ic = amrex::coarsen(i, ratio[0]);
251 const int jc = amrex::coarsen(j, ratio[1]);
252 const int kc = amrex::coarsen(k, ratio[2]);
253 const Real xoff = (static_cast<Real>(i - ic*ratio[0]) + Real(0.5)) / Real(ratio[0]) - Real(0.5);
254 const Real yoff = (static_cast<Real>(j - jc*ratio[1]) + Real(0.5)) / Real(ratio[1]) - Real(0.5);
255 const Real zoff = (static_cast<Real>(k - kc*ratio[2]) + Real(0.5)) / Real(ratio[2]) - Real(0.5);
256 fine(i,j,k,fcomp+ns) = crse(ic,jc,kc,ccomp+ns)
257 + xoff * slope(ic,jc,kc,ns)
258 + yoff * slope(ic,jc,kc,ns+ncomp)
259 + zoff * slope(ic,jc,kc,ns+ncomp*2);
260}
261
262template<typename T>
264void mf_cell_bilin_interp (int i, int j, int k, int n, Array4<T> const& fine, int fcomp,
265 Array4<T const> const& crse, int ccomp, IntVect const& ratio) noexcept
266{
267 int ic = amrex::coarsen(i,ratio[0]);
268 int jc = amrex::coarsen(j,ratio[1]);
269 int kc = amrex::coarsen(k,ratio[2]);
270 int ioff = i - ic*ratio[0];
271 int joff = j - jc*ratio[1];
272 int koff = k - kc*ratio[2];
273 int sx, sy, sz;
274 Real wx, wy, wz;
275 if (ioff*2 < ratio[0]) {
276 sx = -1;
277 wx = Real(ratio[0]+1+2*ioff) / Real(2*ratio[0]);
278 } else {
279 sx = 1;
280 wx = Real(3*ratio[0]-1-2*ioff) / Real(2*ratio[0]);
281 }
282 if (joff*2 < ratio[1]) {
283 sy = -1;
284 wy = Real(ratio[1]+1+2*joff) / Real(2*ratio[1]);
285 } else {
286 sy = 1;
287 wy = Real(3*ratio[1]-1-2*joff) / Real(2*ratio[1]);
288 }
289 if (koff*2 < ratio[2]) {
290 sz = -1;
291 wz = Real(ratio[2]+1+2*koff) / Real(2*ratio[2]);
292 } else {
293 sz = 1;
294 wz = Real(3*ratio[2]-1-2*koff) / Real(2*ratio[2]);
295 }
296 fine(i,j,k,n+fcomp) =
297 crse(ic ,jc ,kc ,n+ccomp)* wx * wy * wz +
298 crse(ic+sx,jc ,kc ,n+ccomp)*(Real(1.0)-wx)* wy * wz +
299 crse(ic ,jc+sy,kc ,n+ccomp)* wx *(Real(1.0)-wy)* wz +
300 crse(ic+sx,jc+sy,kc ,n+ccomp)*(Real(1.0)-wx)*(Real(1.0)-wy)* wz +
301 crse(ic ,jc ,kc+sz,n+ccomp)* wx * wy *(Real(1.0)-wz) +
302 crse(ic+sx,jc ,kc+sz,n+ccomp)*(Real(1.0)-wx)* wy *(Real(1.0)-wz) +
303 crse(ic ,jc+sy,kc+sz,n+ccomp)* wx *(Real(1.0)-wy)*(Real(1.0)-wz) +
304 crse(ic+sx,jc+sy,kc+sz,n+ccomp)*(Real(1.0)-wx)*(Real(1.0)-wy)*(Real(1.0)-wz);
305}
306
308void mf_nodebilin_interp (int i, int j, int k, int n, Array4<Real> const& fine, int fcomp,
309 Array4<Real const> const& crse, int ccomp, IntVect const& ratio) noexcept
310{
311 int ic = amrex::coarsen(i,ratio[0]);
312 int jc = amrex::coarsen(j,ratio[1]);
313 int kc = amrex::coarsen(k,ratio[2]);
314 int ioff = i - ic*ratio[0];
315 int joff = j - jc*ratio[1];
316 int koff = k - kc*ratio[2];
317 Real rxinv = Real(1.0) / Real(ratio[0]);
318 Real ryinv = Real(1.0) / Real(ratio[1]);
319 Real rzinv = Real(1.0) / Real(ratio[2]);
320 if (ioff != 0 && joff != 0 && koff != 0) {
321 // Fine node at center of cell
322 fine(i,j,k,n+fcomp) = rxinv * ryinv * rzinv *
323 (crse(ic ,jc ,kc ,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * (ratio[1]-joff) * (ratio[2]-koff)) +
324 crse(ic+1,jc ,kc ,n+ccomp) * static_cast<Real>(( ioff) * (ratio[1]-joff) * (ratio[2]-koff)) +
325 crse(ic ,jc+1,kc ,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * ( joff) * (ratio[2]-koff)) +
326 crse(ic+1,jc+1,kc ,n+ccomp) * static_cast<Real>(( ioff) * ( joff) * (ratio[2]-koff)) +
327 crse(ic ,jc ,kc+1,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * (ratio[1]-joff) * ( koff)) +
328 crse(ic+1,jc ,kc+1,n+ccomp) * static_cast<Real>(( ioff) * (ratio[1]-joff) * ( koff)) +
329 crse(ic ,jc+1,kc+1,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * ( joff) * ( koff)) +
330 crse(ic+1,jc+1,kc+1,n+ccomp) * static_cast<Real>(( ioff) * ( joff) * ( koff)));
331 } else if (joff != 0 && koff != 0) {
332 // Node on a Y-Z face
333 fine(i,j,k,n+fcomp) = ryinv * rzinv *
334 (crse(ic,jc ,kc ,n+ccomp) * static_cast<Real>((ratio[1]-joff) * (ratio[2]-koff)) +
335 crse(ic,jc+1,kc ,n+ccomp) * static_cast<Real>(( joff) * (ratio[2]-koff)) +
336 crse(ic,jc ,kc+1,n+ccomp) * static_cast<Real>((ratio[1]-joff) * ( koff)) +
337 crse(ic,jc+1,kc+1,n+ccomp) * static_cast<Real>(( joff) * ( koff)));
338 } else if (ioff != 0 && koff != 0) {
339 // Node on a Z-X face
340 fine(i,j,k,n+fcomp) = rxinv * rzinv *
341 (crse(ic ,jc,kc ,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * (ratio[2]-koff)) +
342 crse(ic+1,jc,kc ,n+ccomp) * static_cast<Real>(( ioff) * (ratio[2]-koff)) +
343 crse(ic ,jc,kc+1,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * ( koff)) +
344 crse(ic+1,jc,kc+1,n+ccomp) * static_cast<Real>(( ioff) * ( koff)));
345 } else if (ioff != 0 && joff != 0) {
346 // Node on a X-Y face
347 fine(i,j,k,n+fcomp) = rxinv * ryinv *
348 (crse(ic ,jc ,kc,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * (ratio[1]-joff)) +
349 crse(ic+1,jc ,kc,n+ccomp) * static_cast<Real>(( ioff) * (ratio[1]-joff)) +
350 crse(ic ,jc+1,kc,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * ( joff)) +
351 crse(ic+1,jc+1,kc,n+ccomp) * static_cast<Real>(( ioff) * ( joff)));
352 } else if (ioff != 0) {
353 // Node on X line
354 fine(i,j,k,n+fcomp) = rxinv*(static_cast<Real>(ratio[0]-ioff)*crse(ic ,jc,kc,n+ccomp) +
355 static_cast<Real>( ioff)*crse(ic+1,jc,kc,n+ccomp));
356 } else if (joff != 0) {
357 // Node on Y line
358 fine(i,j,k,n+fcomp) = ryinv*(static_cast<Real>(ratio[1]-joff)*crse(ic,jc ,kc,n+ccomp) +
359 static_cast<Real>( joff)*crse(ic,jc+1,kc,n+ccomp));
360 } else if (koff != 0) {
361 // Node on Z line
362 fine(i,j,k,n+fcomp) = rzinv*(static_cast<Real>(ratio[2]-koff)*crse(ic,jc,kc ,n+ccomp) +
363 static_cast<Real>( koff)*crse(ic,jc,kc+1,n+ccomp));
364 } else {
365 // Node coincident with coarse node
366 fine(i,j,k,n+fcomp) = crse(ic,jc,kc,n+ccomp);
367 }
368}
369
371void mf_cell_quadratic_calcslope (int i, int j, int k, int n,
372 Array4<Real const> const& crse, int ccomp,
373 Array4<Real> const& slope,
374 Box const& domain,
375 BCRec const* bc) noexcept
376{
377 int nu = ccomp + n;
378 Real sx = mf_compute_slopes_x(i, j, k, crse, nu, domain, bc[n]);
379 Real sy = mf_compute_slopes_y(i, j, k, crse, nu, domain, bc[n]);
380 Real sz = mf_compute_slopes_z(i, j, k, crse, nu, domain, bc[n]);
381 Real sxx = mf_cell_quadratic_compute_slopes_xx (i, j, k, crse, nu, domain, bc[n]);
382 Real syy = mf_cell_quadratic_compute_slopes_yy (i, j, k, crse, nu, domain, bc[n]);
383 Real szz = mf_cell_quadratic_compute_slopes_zz (i, j, k, crse, nu, domain, bc[n]);
384 Real sxy = mf_cell_quadratic_compute_slopes_xy (i, j, k, crse, nu, domain, bc[n]);
385 Real sxz = mf_cell_quadratic_compute_slopes_xz (i, j, k, crse, nu, domain, bc[n]);
386 Real syz = mf_cell_quadratic_compute_slopes_yz (i, j, k, crse, nu, domain, bc[n]);
387
388 slope(i,j,k,9*n ) = sx; // x
389 slope(i,j,k,9*n+1) = sy; // y
390 slope(i,j,k,9*n+2) = sz; // z
391 slope(i,j,k,9*n+3) = sxx; // x^2
392 slope(i,j,k,9*n+4) = syy; // y^2
393 slope(i,j,k,9*n+5) = szz; // z^2
394 slope(i,j,k,9*n+6) = sxy; // xy
395 slope(i,j,k,9*n+7) = sxz; // xz
396 slope(i,j,k,9*n+8) = syz; // yz
397}
398
400void mf_cell_quadratic_interp (int i, int j, int k, int n,
401 Array4<Real> const& fine, int fcomp,
402 Array4<Real const> const& crse, int ccomp,
403 Array4<Real const> const& slope,
404 IntVect const& ratio) noexcept
405{
406 int ic = amrex::coarsen(i, ratio[0]);
407 int jc = amrex::coarsen(j, ratio[1]);
408 int kc = amrex::coarsen(k, ratio[2]);
409 int irx = i - ic*ratio[0]; // = abs(i % ratio[0])
410 int jry = j - jc*ratio[1]; // = abs(j % ratio[1])
411 int krz = k - kc*ratio[2]; // = abs(k % ratio[2])
412
413 // Compute offsets.
414 Real xoff = ( Real(irx) + 0.5_rt ) / Real(ratio[0]) - 0.5_rt;
415 Real yoff = ( Real(jry) + 0.5_rt ) / Real(ratio[1]) - 0.5_rt;
416 Real zoff = ( Real(krz) + 0.5_rt ) / Real(ratio[2]) - 0.5_rt;
417
418 fine(i,j,k,fcomp+n) = crse(ic,jc,kc,ccomp+n)
419 + xoff * slope(ic,jc,kc,9*n ) // x
420 + yoff * slope(ic,jc,kc,9*n+1) // y
421 + zoff * slope(ic,jc,kc,9*n+2) // z
422 + 0.5_rt * xoff * xoff * slope(ic,jc,kc,9*n+3) // x^2
423 + 0.5_rt * yoff * yoff * slope(ic,jc,kc,9*n+4) // y^2
424 + 0.5_rt * zoff * zoff * slope(ic,jc,kc,9*n+5) // z^2
425 + xoff * yoff * slope(ic,jc,kc,9*n+6) // xy
426 + xoff * zoff * slope(ic,jc,kc,9*n+7) // xz
427 + yoff * zoff * slope(ic,jc,kc,9*n+8); // yz
428}
429
430} // namespace amrex
431
432#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
Definition AMReX_Amr.cpp:49
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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_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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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:363
BoxND< AMREX_SPACEDIM > Box
Definition AMReX_BaseFwd.H:27
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:21
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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 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
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real mf_cell_quadratic_compute_slopes_yz(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:220
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real mf_cell_quadratic_compute_slopes_zz(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:131
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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:1304
IntVectND< AMREX_SPACEDIM > IntVect
Definition AMReX_BaseFwd.H:30
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real mf_cell_quadratic_compute_slopes_xz(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:186
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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 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_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:384
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real mf_compute_slopes_z(int i, int j, int k, Array4< Real const > const &u, int nu, Box const &domain, BCRec const &bc)
Definition AMReX_MFInterp_C.H:63