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