ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_Interpolation_WENO_Z.H
Go to the documentation of this file.
1 #ifndef ERF_INTERPOLATE_WENO_Z_H_
2 #define ERF_INTERPOLATE_WENO_Z_H_
3 
4 #include "ERF_DataStruct.H"
5 
6 /**
7  * Interpolation operators used for WENO_Z-3 scheme
8  */
9 struct WENO_Z3
10 {
11  WENO_Z3 (const amrex::Array4<const amrex::Real>& phi,
12  const amrex::Real /*upw_frac*/)
13  : m_phi(phi) {}
14 
15  AMREX_GPU_DEVICE
16  AMREX_FORCE_INLINE
17  void
18  InterpolateInX (const int& i,
19  const int& j,
20  const int& k,
21  const int& qty_index,
22  amrex::Real& val_lo,
23  amrex::Real upw_lo) const
24  {
25  // Data to interpolate on
26  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
27  amrex::Real s = m_phi(i , j , k , qty_index);
28  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
29  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
30 
31  if (upw_lo > tol) {
32  val_lo = Evaluate(sm2,sm1,s );
33  } else if (upw_lo < -tol) {
34  val_lo = Evaluate(sp1,s ,sm1);
35  } else {
36  val_lo = 0.5 * (s + sm1);
37  }
38  }
39 
40  AMREX_GPU_DEVICE
41  AMREX_FORCE_INLINE
42  void
43  InterpolateInY (const int& i,
44  const int& j,
45  const int& k,
46  const int& qty_index,
47  amrex::Real& val_lo,
48  amrex::Real upw_lo) const
49  {
50  // Data to interpolate on
51  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
52  amrex::Real s = m_phi(i , j , k , qty_index);
53  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
54  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
55 
56  if (upw_lo > tol) {
57  val_lo = Evaluate(sm2,sm1,s );
58  } else if (upw_lo < -tol) {
59  val_lo = Evaluate(sp1,s ,sm1);
60  } else {
61  val_lo = 0.5 * (s + sm1);
62  }
63  }
64 
65  AMREX_GPU_DEVICE
66  AMREX_FORCE_INLINE
67  void
68  InterpolateInZ (const int& i,
69  const int& j,
70  const int& k,
71  const int& qty_index,
72  amrex::Real& val_lo,
73  amrex::Real upw_lo) const
74  {
75  // Data to interpolate on
76  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
77  amrex::Real s = m_phi(i , j , k , qty_index);
78  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
79  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
80 
81  if (upw_lo > tol) {
82  val_lo = Evaluate(sm2,sm1,s );
83  } else if (upw_lo < -tol) {
84  val_lo = Evaluate(sp1,s ,sm1);
85  } else {
86  val_lo = 0.5 * (s + sm1);
87  }
88  }
89 
90  AMREX_GPU_DEVICE
91  AMREX_FORCE_INLINE
92  amrex::Real
93  Evaluate (const amrex::Real& sm1,
94  const amrex::Real& s ,
95  const amrex::Real& sp1) const
96  {
97  // Smoothing factors
98  amrex::Real b1 = (s - sm1) * (s - sm1);
99  amrex::Real b2 = (sp1 - s) * (sp1 - s);
100 
101  // Weight factors
102  amrex::Real t5 = std::abs(b2 - b1);
103  amrex::Real w1 = g1 * ( 1.0 + (t5*t5) / ((eps + b1) * (eps + b1)) );
104  amrex::Real w2 = g2 * ( 1.0 + (t5*t5) / ((eps + b2) * (eps + b2)) );
105 
106  // Weight factor norm
107  amrex::Real wsum = w1 + w2;
108 
109  // Taylor expansions
110  amrex::Real v1 = -sm1 + 3.0 * s;
111  amrex::Real v2 = s + sp1;
112 
113  // Interpolated value
114  return ( (w1 * v1 + w2 * v2) / (2.0 * wsum) );
115  }
116 
117 private:
118  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
119  const amrex::Real eps=1.0e-6;
120  const amrex::Real tol=1.0e-12;
121  static constexpr amrex::Real g1=(1.0/3.0);
122  static constexpr amrex::Real g2=(2.0/3.0);
123 };
124 
125 /**
126  * Interpolation operators used for WENO_MZQ3 scheme
127  */
128 struct WENO_MZQ3
129 {
130  WENO_MZQ3 (const amrex::Array4<const amrex::Real>& phi,
131  const amrex::Real /*upw_frac*/)
132  : m_phi(phi) {}
133 
134  AMREX_GPU_DEVICE
135  AMREX_FORCE_INLINE
136  void
137  InterpolateInX (const int& i,
138  const int& j,
139  const int& k,
140  const int& qty_index,
141  amrex::Real& val_lo,
142  amrex::Real upw_lo) const
143  {
144  // Data to interpolate on
145  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
146  amrex::Real s = m_phi(i , j , k , qty_index);
147  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
148  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
149 
150  if (upw_lo > tol) {
151  val_lo = Evaluate(sm2,sm1,s );
152  } else if (upw_lo < -tol) {
153  val_lo = Evaluate(sp1,s ,sm1);
154  } else {
155  val_lo = 0.5 * (s + sm1);
156  }
157  }
158 
159  AMREX_GPU_DEVICE
160  AMREX_FORCE_INLINE
161  void
162  InterpolateInY (const int& i,
163  const int& j,
164  const int& k,
165  const int& qty_index,
166  amrex::Real& val_lo,
167  amrex::Real upw_lo) const
168  {
169  // Data to interpolate on
170  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
171  amrex::Real s = m_phi(i , j , k , qty_index);
172  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
173  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
174 
175  if (upw_lo > tol) {
176  val_lo = Evaluate(sm2,sm1,s );
177  } else if (upw_lo < -tol) {
178  val_lo = Evaluate(sp1,s ,sm1);
179  } else {
180  val_lo = 0.5 * (s + sm1);
181  }
182  }
183 
184  AMREX_GPU_DEVICE
185  AMREX_FORCE_INLINE
186  void
187  InterpolateInZ (const int& i,
188  const int& j,
189  const int& k,
190  const int& qty_index,
191  amrex::Real& val_lo,
192  amrex::Real upw_lo) const
193  {
194  // Data to interpolate on
195  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
196  amrex::Real s = m_phi(i , j , k , qty_index);
197  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
198  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
199 
200  if (upw_lo > tol) {
201  val_lo = Evaluate(sm2,sm1,s );
202  } else if (upw_lo < -tol) {
203  val_lo = Evaluate(sp1,s ,sm1);
204  } else {
205  val_lo = 0.5 * (s + sm1);
206  }
207  }
208 
209  AMREX_GPU_DEVICE
210  AMREX_FORCE_INLINE
211  amrex::Real
212  Evaluate (const amrex::Real& sm1,
213  const amrex::Real& s ,
214  const amrex::Real& sp1) const
215  {
216  // Smoothing factors
217  amrex::Real b1 = (s - sm1) * (s - sm1);
218  amrex::Real b2 = (sp1 - s) * (sp1 - s);
219  amrex::Real b3 = ( (13.0 / 12.0) * ((sm1 - 2.0*s + sp1)*(sm1 - 2.0*s + sp1)) ) + ( ((sp1 - sm1)*(sp1 - sm1)) / 4.0 );
220 
221  // Weight factors
222  amrex::Real t5 = ( std::abs(b3 - b1) + std::abs(b3 - b2) ) / 32.0;
223  amrex::Real a1 = g1 * ( 1.0 + (t5*t5) / ((eps + b1) * (eps + b1)) );
224  amrex::Real a2 = g2 * ( 1.0 + (t5*t5) / ((eps + b2) * (eps + b2)) );
225  amrex::Real a3 = g3 * ( 1.0 + (t5*t5) / ((eps + b3) * (eps + b3)) );
226  amrex::Real asum = a1 + a2 + a3;
227  amrex::Real w1 = a1 / asum;
228  amrex::Real w2 = a2 / asum;
229  amrex::Real w3 = a3 / asum;
230 
231  // Taylor expansions
232  amrex::Real v1 = (-sm1 + 3.0 * s) / 2.0;
233  amrex::Real v2 = ( s + sp1) / 2.0;
234  amrex::Real v3 = (-sm1 + 5.0 * s + 2.0 * sp1) / 6.0;
235 
236  // Interpolated value
237  return ( (w3 / g3) * (v3 - g1 * v1 - g2 * v2) + w1 * v1 + w2 * v2 );
238  }
239 
240 private:
241  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
242  const amrex::Real eps=1.0e-6;
243  const amrex::Real tol=1.0e-12;
244  static constexpr amrex::Real g1=(1.0/3.0);
245  static constexpr amrex::Real g2=(1.0/3.0);
246  static constexpr amrex::Real g3=(1.0/3.0);
247 };
248 
249 /**
250  * Interpolation operators used for WENO_Z-5 scheme
251  */
252 struct WENO_Z5
253 {
254  WENO_Z5 (const amrex::Array4<const amrex::Real>& phi,
255  const amrex::Real /*upw_frac*/)
256  : m_phi(phi) {}
257 
258  AMREX_GPU_DEVICE
259  AMREX_FORCE_INLINE
260  void
261  InterpolateInX (const int& i,
262  const int& j,
263  const int& k,
264  const int& qty_index,
265  amrex::Real& val_lo,
266  amrex::Real upw_lo) const
267  {
268  // Data to interpolate on
269  amrex::Real sp2 = m_phi(i+2, j , k , qty_index);
270  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
271  amrex::Real s = m_phi(i , j , k , qty_index);
272  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
273  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
274  amrex::Real sm3 = m_phi(i-3, j , k , qty_index);
275 
276  if (upw_lo > tol) {
277  val_lo = Evaluate(sm3,sm2,sm1,s ,sp1);
278  } else if (upw_lo < -tol) {
279  val_lo = Evaluate(sp2,sp1,s,sm1,sm2);
280  } else {
281  val_lo = 0.5 * (s + sm1);
282  }
283  }
284 
285  AMREX_GPU_DEVICE
286  AMREX_FORCE_INLINE
287  void
288  InterpolateInY (const int& i,
289  const int& j,
290  const int& k,
291  const int& qty_index,
292  amrex::Real& val_lo,
293  amrex::Real upw_lo) const
294  {
295  // Data to interpolate on
296  amrex::Real sp2 = m_phi(i , j+2, k , qty_index);
297  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
298  amrex::Real s = m_phi(i , j , k , qty_index);
299  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
300  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
301  amrex::Real sm3 = m_phi(i , j-3, k , qty_index);
302 
303  if (upw_lo > tol) {
304  val_lo = Evaluate(sm3,sm2,sm1,s ,sp1);
305  } else if (upw_lo < -tol) {
306  val_lo = Evaluate(sp2,sp1,s,sm1,sm2);
307  } else {
308  val_lo = 0.5 * (s + sm1);
309  }
310  }
311 
312  AMREX_GPU_DEVICE
313  AMREX_FORCE_INLINE
314  void
315  InterpolateInZ (const int& i,
316  const int& j,
317  const int& k,
318  const int& qty_index,
319  amrex::Real& val_lo,
320  amrex::Real upw_lo) const
321  {
322  // Data to interpolate on
323  amrex::Real sp2 = m_phi(i , j , k+2, qty_index);
324  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
325  amrex::Real s = m_phi(i , j , k , qty_index);
326  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
327  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
328  amrex::Real sm3 = m_phi(i , j , k-3, qty_index);
329 
330  if (upw_lo > tol) {
331  val_lo = Evaluate(sm3,sm2,sm1,s ,sp1);
332  } else if (upw_lo < -tol) {
333  val_lo = Evaluate(sp2,sp1,s,sm1,sm2);
334  } else {
335  val_lo = 0.5 * (s + sm1);
336  }
337  }
338 
339  AMREX_GPU_DEVICE
340  AMREX_FORCE_INLINE
341  amrex::Real
342  Evaluate (const amrex::Real& sm2,
343  const amrex::Real& sm1,
344  const amrex::Real& s ,
345  const amrex::Real& sp1,
346  const amrex::Real& sp2) const
347  {
348  // Smoothing factors
349  amrex::Real b1 = c1 * (sm2 - 2.0 * sm1 + s) * (sm2 - 2.0 * sm1 + s) +
350  0.25 * (sm2 - 4.0 * sm1 + 3.0 * s) * (sm2 - 4.0 * sm1 + 3.0 * s);
351  amrex::Real b2 = c1 * (sm1 - 2.0 * s + sp1) * (sm1 - 2.0 * s + sp1) +
352  0.25 * (sm1 - sp1) * (sm1 - sp1);
353  amrex::Real b3 = c1 * (s - 2.0 * sp1 + sp2) * (s - 2.0 * sp1 + sp2) +
354  0.25 * (3.0 * s - 4.0 * sp1 + sp2) * (3.0 * s - 4.0 * sp1 + sp2);
355 
356  // Weight factors
357  amrex::Real t5 = std::abs(b3 - b1);
358  amrex::Real w1 = g1 * ( 1.0 + (t5*t5) / ((eps + b1) * (eps + b1)) );
359  amrex::Real w2 = g2 * ( 1.0 + (t5*t5) / ((eps + b2) * (eps + b2)) );
360  amrex::Real w3 = g3 * ( 1.0 + (t5*t5) / ((eps + b3) * (eps + b3)) );
361 
362  // Weight factor norm
363  amrex::Real wsum = w1 + w2 + w3;
364 
365  // Taylor expansions
366  amrex::Real v1 = 2.0 * sm2 - 7.0 * sm1 + 11.0 * s;
367  amrex::Real v2 = -sm1 + 5.0 * s + 2.0 * sp1;
368  amrex::Real v3 = 2.0 * s + 5.0 * sp1 - sp2;
369 
370  // Interpolated value
371  return ( (w1 * v1 + w2 * v2 + w3 * v3) / (6.0 * wsum) );
372  }
373 
374 private:
375  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
376  const amrex::Real eps=1.0e-6;
377  const amrex::Real tol=1.0e-12;
378  static constexpr amrex::Real c1=(13.0/12.0);
379  static constexpr amrex::Real g1=(1.0/10.0);
380  static constexpr amrex::Real g2=(3.0/5.0);
381  static constexpr amrex::Real g3=(3.0/10.0);
382 };
383 
384 /**
385  * Interpolation operators used for WENO_Z-7 scheme
386  */
387 struct WENO_Z7
388 {
389  WENO_Z7 (const amrex::Array4<const amrex::Real>& phi,
390  const amrex::Real /*upw_frac*/)
391  : m_phi(phi) {}
392 
393  AMREX_GPU_DEVICE
394  AMREX_FORCE_INLINE
395  void
396  InterpolateInX (const int& i,
397  const int& j,
398  const int& k,
399  const int& qty_index,
400  amrex::Real& val_lo,
401  amrex::Real upw_lo) const
402  {
403  // Data to interpolate on
404  amrex::Real sp3 = m_phi(i+3, j , k , qty_index);
405  amrex::Real sp2 = m_phi(i+2, j , k , qty_index);
406  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
407  amrex::Real s = m_phi(i , j , k , qty_index);
408  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
409  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
410  amrex::Real sm3 = m_phi(i-3, j , k , qty_index);
411  amrex::Real sm4 = m_phi(i-4, j , k , qty_index);
412 
413  if (upw_lo > tol) {
414  val_lo = Evaluate(sm4,sm3,sm2,sm1,s ,sp1,sp2);
415  } else if (upw_lo < -tol) {
416  val_lo = Evaluate(sp3,sp2,sp1,s ,sm1,sm2,sm3);
417  } else {
418  val_lo = 0.5 * (s + sm1);
419  }
420  }
421 
422  AMREX_GPU_DEVICE
423  AMREX_FORCE_INLINE
424  void
425  InterpolateInY (const int& i,
426  const int& j,
427  const int& k,
428  const int& qty_index,
429  amrex::Real& val_lo,
430  amrex::Real upw_lo) const
431  {
432  // Data to interpolate on
433  amrex::Real sp3 = m_phi(i , j+3, k , qty_index);
434  amrex::Real sp2 = m_phi(i , j+2, k , qty_index);
435  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
436  amrex::Real s = m_phi(i , j , k , qty_index);
437  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
438  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
439  amrex::Real sm3 = m_phi(i , j-3, k , qty_index);
440  amrex::Real sm4 = m_phi(i , j-4, k , qty_index);
441 
442  if (upw_lo > tol) {
443  val_lo = Evaluate(sm4,sm3,sm2,sm1,s ,sp1,sp2);
444  } else if (upw_lo < -tol) {
445  val_lo = Evaluate(sp3,sp2,sp1,s ,sm1,sm2,sm3);
446  } else {
447  val_lo = 0.5 * (s + sm1);
448  }
449  }
450 
451  AMREX_GPU_DEVICE
452  AMREX_FORCE_INLINE
453  void
454  InterpolateInZ (const int& i,
455  const int& j,
456  const int& k,
457  const int& qty_index,
458  amrex::Real& val_lo,
459  amrex::Real upw_lo) const
460  {
461  // Data to interpolate on
462  amrex::Real sp3 = m_phi(i , j , k+2, qty_index);
463  amrex::Real sp2 = m_phi(i , j , k+2, qty_index);
464  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
465  amrex::Real s = m_phi(i , j , k , qty_index);
466  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
467  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
468  amrex::Real sm3 = m_phi(i , j , k-3, qty_index);
469  amrex::Real sm4 = m_phi(i , j , k-3, qty_index);
470 
471  if (upw_lo > tol) {
472  val_lo = Evaluate(sm4,sm3,sm2,sm1,s ,sp1,sp2);
473  } else if (upw_lo < -tol) {
474  val_lo = Evaluate(sp3,sp2,sp1,s ,sm1,sm2,sm3);
475  } else {
476  val_lo = 0.5 * (s + sm1);
477  }
478  }
479 
480  AMREX_GPU_DEVICE
481  AMREX_FORCE_INLINE
482  amrex::Real
483  Evaluate (const amrex::Real& sm3,
484  const amrex::Real& sm2,
485  const amrex::Real& sm1,
486  const amrex::Real& s ,
487  const amrex::Real& sp1,
488  const amrex::Real& sp2,
489  const amrex::Real& sp3) const
490  {
491  // Smoothing factors
492  amrex::Real b1 = ( sm3 * sm3 * 6649./2880.0
493  - sm3 * sm2 * 2623./160.0
494  + sm3 * sm1 * 9449./480.0
495  - sm3 * s * 11389./1440.0
496  + sm2 * sm2 * 28547./960.0
497  - sm2 * sm1 * 35047./480.0
498  + sm2 * s * 14369./480.0
499  + sm1 * sm1 * 44747./960.0
500  - sm1 * s * 6383./160.0
501  + s * s * 25729./2880.0 );
502  amrex::Real b2 = ( sm2 * sm2 * 3169/2880.0
503  - sm2 * sm1 * 3229/480.0
504  + sm2 * s * 3169/480.0
505  - sm2 * sp1 * 2989/1440.0
506  + sm1 * sm1 * 11147/960.0
507  - sm1 * s * 11767/480.0
508  + sm1 * sp1 * 1283/160.0
509  + s * s * 13667/960.0
510  - s * sp1 * 5069/480.0
511  + sp1 * sp1 * 6649/2880.0 );
512  amrex::Real b3 = ( sm1 * sm1 * 6649./2880.0
513  - sm1 * s * 5069./480.0
514  + sm1 * sp1 * 1283./160.0
515  - sm1 * sp2 * 2989./1440.0
516  + s * s * 13667./960.0
517  - s * sp1 * 11767./480.0
518  + s * sp2 * 3169./480.0
519  + sp1 * sp1 * 11147./960.0
520  - sp1 * sp2 * 3229./480.0
521  + sp2 * sp2 * 3169./2880.0 );
522  amrex::Real b4 = ( s * s * 25729./2880.
523  - s * sp1 * 6383./160.
524  + s * sp2 * 14369./480.
525  - s * sp3 * 11389./1440.
526  + sp1 * sp1 * 44747./960.
527  - sp1 * sp2 * 35047./480.
528  + sp1 * sp3 * 9449./480.
529  + sp2 * sp2 * 28547./960.
530  - sp2 * sp3 * 2623./160.
531  + sp3 * sp3 * 6649./2880. );
532 
533  // Weight factors
534  amrex::Real t5 = std::abs(b1 - b2 - b3 + b4);
535  amrex::Real w1 = g1 * ( 1.0 + (t5*t5) / ((eps + b1) * (eps + b1)) );
536  amrex::Real w2 = g2 * ( 1.0 + (t5*t5) / ((eps + b2) * (eps + b2)) );
537  amrex::Real w3 = g3 * ( 1.0 + (t5*t5) / ((eps + b3) * (eps + b3)) );
538  amrex::Real w4 = g4 * ( 1.0 + (t5*t5) / ((eps + b3) * (eps + b3)) );
539 
540  // Weight factor norm
541  amrex::Real wsum = w1 + w2 + w3 + w4;
542 
543  // Taylor expansions
544  amrex::Real v1 = (-0.3125)*sm3 + ( 1.3125)*sm2 + (-2.1875)*sm1 + ( 2.1875)*s;
545  amrex::Real v2 = ( 0.0625)*sm2 + (-0.3125)*sm1 + ( 0.9375)*s + ( 0.3125)*sp1;
546  amrex::Real v3 = (-0.0625)*sm1 + ( 0.5625)*s + ( 0.5625)*sp1 + (-0.0625)*sp2;
547  amrex::Real v4 = ( 0.3125)*s + ( 0.9375)*sp1 + (-0.3125)*sp2 + ( 0.0625)*sp3;
548 
549  // Interpolated value
550  return ( (w1 * v1 + w2 * v2 + w3 * v3 + w4 * v4) / (wsum) );
551  }
552 
553 private:
554  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
555  const amrex::Real eps=1.0e-40;
556  const amrex::Real tol=1.0e-12;
557  static constexpr amrex::Real g1=( 1.0/64.0);
558  static constexpr amrex::Real g2=(21.0/64.0);
559  static constexpr amrex::Real g3=(35.0/64.0);
560  static constexpr amrex::Real g4=( 7.0/64.0);
561 };
562 #endif
Definition: ERF_Interpolation_WENO_Z.H:129
static constexpr amrex::Real g2
Definition: ERF_Interpolation_WENO_Z.H:245
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_WENO_Z.H:241
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:187
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sm1, const amrex::Real &s, const amrex::Real &sp1) const
Definition: ERF_Interpolation_WENO_Z.H:212
static constexpr amrex::Real g3
Definition: ERF_Interpolation_WENO_Z.H:246
WENO_MZQ3(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_WENO_Z.H:130
const amrex::Real eps
Definition: ERF_Interpolation_WENO_Z.H:242
static constexpr amrex::Real g1
Definition: ERF_Interpolation_WENO_Z.H:244
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:137
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:162
const amrex::Real tol
Definition: ERF_Interpolation_WENO_Z.H:243
Definition: ERF_Interpolation_WENO_Z.H:10
const amrex::Real eps
Definition: ERF_Interpolation_WENO_Z.H:119
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:68
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_WENO_Z.H:118
WENO_Z3(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_WENO_Z.H:11
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sm1, const amrex::Real &s, const amrex::Real &sp1) const
Definition: ERF_Interpolation_WENO_Z.H:93
static constexpr amrex::Real g2
Definition: ERF_Interpolation_WENO_Z.H:122
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:18
static constexpr amrex::Real g1
Definition: ERF_Interpolation_WENO_Z.H:121
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:43
const amrex::Real tol
Definition: ERF_Interpolation_WENO_Z.H:120
Definition: ERF_Interpolation_WENO_Z.H:253
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sm2, const amrex::Real &sm1, const amrex::Real &s, const amrex::Real &sp1, const amrex::Real &sp2) const
Definition: ERF_Interpolation_WENO_Z.H:342
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:315
static constexpr amrex::Real g3
Definition: ERF_Interpolation_WENO_Z.H:381
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:288
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:261
static constexpr amrex::Real g1
Definition: ERF_Interpolation_WENO_Z.H:379
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_WENO_Z.H:375
static constexpr amrex::Real c1
Definition: ERF_Interpolation_WENO_Z.H:378
const amrex::Real eps
Definition: ERF_Interpolation_WENO_Z.H:376
static constexpr amrex::Real g2
Definition: ERF_Interpolation_WENO_Z.H:380
const amrex::Real tol
Definition: ERF_Interpolation_WENO_Z.H:377
WENO_Z5(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_WENO_Z.H:254
Definition: ERF_Interpolation_WENO_Z.H:388
static constexpr amrex::Real g1
Definition: ERF_Interpolation_WENO_Z.H:557
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sm3, const amrex::Real &sm2, const amrex::Real &sm1, const amrex::Real &s, const amrex::Real &sp1, const amrex::Real &sp2, const amrex::Real &sp3) const
Definition: ERF_Interpolation_WENO_Z.H:483
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:396
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_WENO_Z.H:554
const amrex::Real eps
Definition: ERF_Interpolation_WENO_Z.H:555
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:454
WENO_Z7(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_WENO_Z.H:389
static constexpr amrex::Real g4
Definition: ERF_Interpolation_WENO_Z.H:560
static constexpr amrex::Real g2
Definition: ERF_Interpolation_WENO_Z.H:558
const amrex::Real tol
Definition: ERF_Interpolation_WENO_Z.H:556
static constexpr amrex::Real g3
Definition: ERF_Interpolation_WENO_Z.H:559
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_WENO_Z.H:425