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