Block-Structured AMR Software Framework
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 
6 namespace amrex {
7 
9 void 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 
90 void 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 
138 void 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 
192 void 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 
206 void mf_cell_cons_lin_interp_mcslope_rz (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) noexcept
210 {
211  int nu = ns + scomp;
212 
213  Real sx = Real(0.0);
214  Real sy = Real(0.0);
215 
216  // x-direction
217  if (ratio[0] > 1) {
218  Real dc = mf_compute_slopes_x(i, j, 0, u, nu, domain, bc[ns]);
219  Real df = Real(2.0) * (u(i+1,j,0,nu) - u(i ,j,0,nu));
220  Real db = Real(2.0) * (u(i ,j,0,nu) - u(i-1,j,0,nu));
221  sx = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
222  sx = std::copysign(Real(1.),dc)*amrex::min(sx,std::abs(dc));
223  }
224 
225  // y-direction
226  if (ratio[1] > 1) {
227  Real dc = mf_compute_slopes_y(i, j, 0, u, nu, domain, bc[ns]);
228  Real df = Real(2.0) * (u(i,j+1,0,nu) - u(i,j ,0,nu));
229  Real db = Real(2.0) * (u(i,j ,0,nu) - u(i,j-1,0,nu));
230  sy = (df*db >= Real(0.0)) ? amrex::min(std::abs(df),std::abs(db)) : Real(0.);
231  sy = std::copysign(Real(1.),dc)*amrex::min(sy,std::abs(dc));
232  }
233 
234  Real alpha = Real(1.0);
235  if (sx != Real(0.0) || sy != Real(0.0)) {
236  const Real drc = drf * ratio[0];
237  const Real rcm = i * drc + rlo;
238  const Real rcp = (i+1) * drc + rlo;
239  const Real vcm = rcm*rcm;
240  const Real vcp = rcp*rcp;
241  Real rfm = i*ratio[0] * drf + rlo;
242  Real rfp = (i*ratio[0] + 1) * drf + rlo;
243  Real vfm = rfm*rfm;
244  Real vfp = rfp*rfp;
245  Real xlo = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
246  rfm = ((i+1)*ratio[0] - 1) * drf + rlo;
247  rfp = (i+1)*ratio[0] * drf + rlo;
248  vfm = rfm*rfm;
249  vfp = rfp*rfp;
250  Real xhi = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
251  Real dumax = amrex::max(sx*xlo, sx*xhi)
252  + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1]);
253  Real dumin = -amrex::min(sx*xlo, sx*xhi)
254  + std::abs(sy) * Real(ratio[1]-1)/Real(2*ratio[1]);
255  Real umax = u(i,j,0,nu);
256  Real umin = u(i,j,0,nu);
257  int ilim = ratio[0] > 1 ? 1 : 0;
258  int jlim = ratio[1] > 1 ? 1 : 0;
259  for (int joff = -jlim; joff <= jlim; ++joff) {
260  for (int ioff = -ilim; ioff <= ilim; ++ioff) {
261  umin = amrex::min(umin, u(i+ioff,j+joff,0,nu));
262  umax = amrex::max(umax, u(i+ioff,j+joff,0,nu));
263  }}
264  if (dumax * alpha > (umax - u(i,j,0,nu))) {
265  alpha = (umax - u(i,j,0,nu)) / dumax;
266  }
267  if (dumin * alpha > (u(i,j,0,nu) - umin)) {
268  alpha = (u(i,j,0,nu) - umin) / dumin;
269  }
270  }
271 
272  slope(i,j,0,ns ) = sx * alpha;
273  slope(i,j,0,ns+ ncomp) = sy * alpha;
274 }
275 
277 void mf_cell_cons_lin_interp_rz (int i, int j, int ns, Array4<Real> const& fine, int fcomp,
279  int ccomp, int ncomp, IntVect const& ratio, Real drf, Real rlo) noexcept
280 {
281  const int ic = amrex::coarsen(i, ratio[0]);
282  const int jc = amrex::coarsen(j, ratio[1]);
283  const Real drc = drf * ratio[0];
284  const Real rcm = ic * drc + rlo;
285  const Real rcp = (ic+1) * drc + rlo;
286  const Real rfm = i * drf + rlo;
287  const Real rfp = (i +1) * drf + rlo;
288  const Real vcm = rcm*rcm;
289  const Real vcp = rcp*rcp;
290  const Real vfm = rfm*rfm;
291  const Real vfp = rfp*rfp;
292  const Real xoff = Real(0.5) * ((vfm+vfp) - (vcm+vcp)) / (vcp - vcm);
293  const Real yoff = (j - jc*ratio[1] + Real(0.5)) / Real(ratio[1]) - Real(0.5);
294  fine(i,j,0,fcomp+ns) = crse(ic,jc,0,ccomp+ns)
295  + xoff * slope(ic,jc,0,ns)
296  + yoff * slope(ic,jc,0,ns+ncomp);
297 }
298 
299 template<typename T>
301 void mf_cell_bilin_interp (int i, int j, int, int n, Array4<T> const& fine, int fcomp,
302  Array4<T const> const& crse, int ccomp, IntVect const& ratio) noexcept
303 {
304  int ic = amrex::coarsen(i,ratio[0]);
305  int jc = amrex::coarsen(j,ratio[1]);
306  int ioff = i - ic*ratio[0];
307  int joff = j - jc*ratio[1];
308  int sx, sy;
309  Real wx, wy;
310  if (ioff*2 < ratio[0]) {
311  sx = -1;
312  wx = Real(ratio[0]+1+2*ioff) / Real(2*ratio[0]);
313  } else {
314  sx = 1;
315  wx = Real(3*ratio[0]-1-2*ioff) / Real(2*ratio[0]);
316  }
317  if (joff*2 < ratio[1]) {
318  sy = -1;
319  wy = Real(ratio[1]+1+2*joff) / Real(2*ratio[1]);
320  } else {
321  sy = 1;
322  wy = Real(3*ratio[1]-1-2*joff) / Real(2*ratio[1]);
323  }
324  fine(i,j,0,n+fcomp) =
325  crse(ic ,jc ,0,n+ccomp)* wx * wy +
326  crse(ic+sx,jc ,0,n+ccomp)*(Real(1.0)-wx)* wy +
327  crse(ic ,jc+sy,0,n+ccomp)* wx *(Real(1.0)-wy) +
328  crse(ic+sx,jc+sy,0,n+ccomp)*(Real(1.0)-wx)*(Real(1.0)-wy);
329 }
330 
332 void mf_nodebilin_interp (int i, int j, int, int n, Array4<Real> const& fine, int fcomp,
333  Array4<Real const> const& crse, int ccomp, IntVect const& ratio) noexcept
334 {
335  int ic = amrex::coarsen(i,ratio[0]);
336  int jc = amrex::coarsen(j,ratio[1]);
337  int ioff = i - ic*ratio[0];
338  int joff = j - jc*ratio[1];
339  Real rxinv = Real(1.0) / Real(ratio[0]);
340  Real ryinv = Real(1.0) / Real(ratio[1]);
341  if (ioff != 0 && joff != 0) {
342  // Node on a X-Y face
343  fine(i,j,0,n+fcomp) = rxinv * ryinv *
344  (crse(ic ,jc ,0,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * (ratio[1]-joff)) +
345  crse(ic+1,jc ,0,n+ccomp) * static_cast<Real>(( ioff) * (ratio[1]-joff)) +
346  crse(ic ,jc+1,0,n+ccomp) * static_cast<Real>((ratio[0]-ioff) * ( joff)) +
347  crse(ic+1,jc+1,0,n+ccomp) * static_cast<Real>(( ioff) * ( joff)));
348  } else if (ioff != 0) {
349  // Node on X line
350  fine(i,j,0,n+fcomp) = rxinv*(static_cast<Real>(ratio[0]-ioff)*crse(ic ,jc,0,n+ccomp) +
351  static_cast<Real>( ioff)*crse(ic+1,jc,0,n+ccomp));
352  } else if (joff != 0) {
353  // Node on Y line
354  fine(i,j,0,n+fcomp) = ryinv*(static_cast<Real>(ratio[1]-joff)*crse(ic,jc ,0,n+ccomp) +
355  static_cast<Real>( joff)*crse(ic,jc+1,0,n+ccomp));
356  } else {
357  // Node coincident with coarse node
358  fine(i,j,0,n+fcomp) = crse(ic,jc,0,n+ccomp);
359  }
360 }
361 
363 void mf_cell_quadratic_calcslope (int i, int j, int /*k*/, int n,
364  Array4<Real const> const& crse, int ccomp,
365  Array4<Real> const& slope,
366  Box const& domain,
367  BCRec const* bc) noexcept
368 {
369  int nu = ccomp + n;
370  Real sx = mf_compute_slopes_x(i, j, 0, crse, nu, domain, bc[n]);
371  Real sy = mf_compute_slopes_y(i, j, 0, crse, nu, domain, bc[n]);
372  Real sxx = mf_cell_quadratic_compute_slopes_xx (i, j, 0, crse, nu, domain, bc[n]);
373  Real syy = mf_cell_quadratic_compute_slopes_yy (i, j, 0, crse, nu, domain, bc[n]);
374  Real sxy = mf_cell_quadratic_compute_slopes_xy (i, j, 0, crse, nu, domain, bc[n]);
375 
376  slope(i,j,0,5*n ) = sx; // x
377  slope(i,j,0,5*n+1) = sy; // y
378  slope(i,j,0,5*n+2) = sxx; // x^2
379  slope(i,j,0,5*n+3) = syy; // y^2
380  slope(i,j,0,5*n+4) = sxy; // xy
381 }
382 
384 void mf_cell_quadratic_interp (int i, int j, int /*k*/, int n,
385  Array4<Real> const& fine, int fcomp,
386  Array4<Real const> const& crse, int ccomp,
387  Array4<Real const> const& slope,
388  IntVect const& ratio) noexcept
389 {
390  int ic = amrex::coarsen(i, ratio[0]);
391  int jc = amrex::coarsen(j, ratio[1]);
392  int irx = i - ic*ratio[0]; // = abs(i % ratio[0])
393  int jry = j - jc*ratio[1]; // = abs(j % ratio[1])
394 
395  // Compute offsets.
396  Real xoff = ( Real(irx) + 0.5_rt ) / Real(ratio[0]) - 0.5_rt;
397  Real yoff = ( Real(jry) + 0.5_rt ) / Real(ratio[1]) - 0.5_rt;
398 
399  fine(i,j,0,fcomp+n) = crse(ic,jc,0,ccomp+n)
400  + xoff * slope(ic,jc,0,5*n ) // x
401  + yoff * slope(ic,jc,0,5*n+1) // y
402  + 0.5_rt * xoff * xoff * slope(ic,jc,0,5*n+2) // x^2
403  + 0.5_rt * yoff * yoff * slope(ic,jc,0,5*n+3) // y^2
404  + xoff * yoff * slope(ic,jc,0,5*n+4); // xy
405 }
406 
408 void mf_cell_quadratic_interp_rz (int i, int j, int /*k*/, int n,
409  Array4<Real> const& fine, int fcomp,
410  Array4<Real const> const& crse, int ccomp,
411  Array4<Real const> const& slope,
412  IntVect const& ratio,
413  GeometryData const& cs_geomdata,
414  GeometryData const& fn_geomdata) noexcept
415 {
416  int ic = amrex::coarsen(i, ratio[0]);
417  int jc = amrex::coarsen(j, ratio[1]);
418  int jry = j - jc*ratio[1]; // = abs(j % ratio[1])
419 
420  // Extract geometry data for radial direction.
421  Real fn_dr = fn_geomdata.CellSize(0);
422  Real fn_rlo = fn_geomdata.ProbLo(0);
423  Real cs_dr = cs_geomdata.CellSize(0);
424  Real cs_rlo = cs_geomdata.ProbLo(0);
425 
426  // Compute radial mesh.
427  Real fn_rm = fn_rlo + Real(i) * fn_dr;
428  Real fn_rp = fn_rlo + Real(i+1) * fn_dr;
429  Real cs_rm = cs_rlo + Real(ic) * cs_dr;
430  Real cs_rp = cs_rlo + Real(ic+1) * cs_dr;
431 
432  // Compute radial cell edge coords.
433  // Note: omit the 0.5 prefactor here -- we'll add it once at the last step
434  Real fn_vm = fn_rm*fn_rm;
435  Real fn_vp = fn_rp*fn_rp;
436  Real fcen = fn_vm + fn_vp;
437  Real cs_vm = cs_rm*cs_rm;
438  Real cs_vp = cs_rp*cs_rp;
439  Real ccen = cs_vm + cs_vp;
440 
441  // Compute r offset
442  Real xoff = 0.5_rt * ( fcen - ccen ) / ( cs_vp - cs_vm ) ;
443 
444  // Compute z offset (same formula as Cartesian)
445  Real yoff = ( Real(jry) + 0.5_rt ) / Real(ratio[1]) - 0.5_rt;
446 
447  fine(i,j,0,fcomp+n) = crse(ic,jc,0,ccomp+n)
448  + xoff * slope(ic,jc,0,5*n ) // x
449  + yoff * slope(ic,jc,0,5*n+1) // y
450  + 0.5_rt * xoff * xoff * slope(ic,jc,0,5*n+2) // x^2
451  + 0.5_rt * yoff * yoff * slope(ic,jc,0,5*n+3) // y^2
452  + xoff * yoff * slope(ic,jc,0,5*n+4); // xy
453 }
454 
455 } // namespace amrex
456 
457 #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 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_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:206
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 constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
Definition: AMReX_Algorithm.H:35
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE 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 T abs(const GpuComplex< T > &a_z) noexcept
Return the absolute value of a complex number.
Definition: AMReX_GpuComplex.H:356
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_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
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 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:408
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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:277
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 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
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
Definition: AMReX_Array4.H:61
Definition: AMReX_Geometry.H:25