ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Interpolation_UPW.H
Go to the documentation of this file.
1 #ifndef ERF_INTERPOLATE_UPW_H_
2 #define ERF_INTERPOLATE_UPW_H_
3 
4 #include "ERF_DataStruct.H"
5 
6 /**
7  * Interpolation operators used for 2nd order centered scheme
8  */
9 struct CENTERED2
10 {
11  AMREX_GPU_HOST_DEVICE
12  AMREX_FORCE_INLINE
13  CENTERED2 (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  // Data to interpolate on
28  amrex::Real s = m_phi(i , j , k , qty_index);
29  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
30 
31  // Interpolate lo
32  val_lo = Evaluate(s,sm1);
33  }
34 
35  AMREX_GPU_DEVICE
36  AMREX_FORCE_INLINE
37  void
38  InterpolateInY (const int& i,
39  const int& j,
40  const int& k,
41  const int& qty_index,
42  amrex::Real& val_lo,
43  amrex::Real /*upw_lo*/) const
44  {
45  // Data to interpolate on
46  amrex::Real s = m_phi(i , j , k , qty_index);
47  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
48 
49  // Interpolate lo
50  val_lo = Evaluate(s,sm1);
51  }
52 
53  AMREX_GPU_DEVICE
54  AMREX_FORCE_INLINE
55  void
56  InterpolateInZ (const int& i,
57  const int& j,
58  const int& k,
59  const int& qty_index,
60  amrex::Real& val_lo,
61  amrex::Real /*upw_lo*/) const
62  {
63  // Data to interpolate on
64  amrex::Real s = m_phi(i , j , k , qty_index);
65  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
66 
67  // Interpolate lo
68  val_lo = Evaluate(s,sm1);
69  }
70 
71  AMREX_GPU_DEVICE
72  AMREX_FORCE_INLINE
75  const amrex::Real& sm1) const
76  {
77  // Averages and diffs
78  amrex::Real a1 = (s + sm1);
79 
80  // Interpolated value
81  return ( g1*a1 );
82  }
83 
84  int GetUpwindCellNumber() const {return 1;}
85 
86 private:
87  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
88  static constexpr amrex::Real g1=(myhalf);
89 };
90 
91 /**
92  * Interpolation operators used for 3rd order upwind scheme
93  */
94 struct UPWIND3
95 {
96  AMREX_GPU_HOST_DEVICE
97  AMREX_FORCE_INLINE
98  UPWIND3 (const amrex::Array4<const amrex::Real>& phi,
99  const amrex::Real upw_frac)
100  : m_phi(phi)
101  , m_upw_frac(upw_frac)
102  {}
103 
104  AMREX_GPU_DEVICE
105  AMREX_FORCE_INLINE
106  void
107  InterpolateInX (const int& i,
108  const int& j,
109  const int& k,
110  const int& qty_index,
111  amrex::Real& val_lo,
112  amrex::Real upw_lo) const
113  {
114  // Data to interpolate on
115  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
116  amrex::Real s = m_phi(i , j , k , qty_index);
117  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
118  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
119 
120  // Upwinding flags
121  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
122 
123  // Add blending
124  upw_lo *= m_upw_frac;
125 
126  // Interpolate lo
127  val_lo = Evaluate(sp1,s,sm1,sm2,upw_lo);
128  }
129 
130  AMREX_GPU_DEVICE
131  AMREX_FORCE_INLINE
132  void
133  InterpolateInY (const int& i,
134  const int& j,
135  const int& k,
136  const int& qty_index,
137  amrex::Real& val_lo,
138  amrex::Real upw_lo) const
139  {
140  // Data to interpolate on
141  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
142  amrex::Real s = m_phi(i , j , k , qty_index);
143  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
144  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
145 
146  // Upwinding flags
147  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
148 
149  // Add blending
150  upw_lo *= m_upw_frac;
151 
152  // Interpolate lo
153  val_lo = Evaluate(sp1,s,sm1,sm2,upw_lo);
154  }
155 
156  AMREX_GPU_DEVICE
157  AMREX_FORCE_INLINE
158  void
159  InterpolateInZ (const int& i,
160  const int& j,
161  const int& k,
162  const int& qty_index,
163  amrex::Real& val_lo,
164  amrex::Real upw_lo) const
165  {
166  // Data to interpolate on
167  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
168  amrex::Real s = m_phi(i , j , k , qty_index);
169  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
170  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
171 
172  // Upwinding flags
173  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
174 
175  // Add blending
176  upw_lo *= m_upw_frac;
177 
178  // Interpolate lo
179  val_lo = Evaluate(sp1,s,sm1,sm2,upw_lo);
180  }
181 
182  AMREX_GPU_DEVICE
183  AMREX_FORCE_INLINE
185  Evaluate (const amrex::Real& sp1,
186  const amrex::Real& s,
187  const amrex::Real& sm1,
188  const amrex::Real& sm2,
189  const amrex::Real& upw) const
190  {
191  // Averages and diffs
192  amrex::Real a1 = (s + sm1);
193  amrex::Real d1 = (s - sm1);
194  amrex::Real a2 = (sp1 + sm2);
195  amrex::Real d2 = (sp1 - sm2);
196 
197  // Interpolated value
198  return ( g1*a1 - g2*a2 + upw * g2 * (d2 - three*d1) );
199  }
200 
201  int GetUpwindCellNumber() const {return 2;}
202 
203  AMREX_GPU_HOST_DEVICE
204  AMREX_FORCE_INLINE
205  void SetUpwinding(amrex::Real upw_frac) {m_upw_frac = upw_frac;}
206 
207 private:
208  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
209  amrex::Real m_upw_frac; // Upwinding fraction
210  static constexpr amrex::Real g1=(amrex::Real(7.0)/amrex::Real(12.0));
211  static constexpr amrex::Real g2=(amrex::Real(1)/amrex::Real(12.0));
212 };
213 
214 
215 /**
216  * Interpolation operators used for 3rd order upwind scheme
217  */
218 struct UPWIND3SL
219 {
220  AMREX_GPU_HOST_DEVICE
221  AMREX_FORCE_INLINE
222  UPWIND3SL (const amrex::Array4<const amrex::Real>& phi,
223  const amrex::Real /*upw_frac*/)
224  : m_phi(phi)
225  {}
226 
227  AMREX_GPU_DEVICE
228  AMREX_FORCE_INLINE
229  void
230  InterpolateInX (const int& i,
231  const int& j,
232  const int& k,
233  const int& qty_index,
234  amrex::Real& val_lo,
235  amrex::Real upw_lo) const
236  {
237  // Upwinding flags
238  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
239 
240  // Data to interpolate on
241  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
242  amrex::Real s = m_phi(i , j , k , qty_index);
243  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
244  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
245 
246  // Left and right fluxes
247  amrex::Real Fhi = Evaluate(s ,sm1,sm2);
248  amrex::Real Flo = Evaluate(sm1,s ,sp1);
249 
250  // Numerical flux
251  val_lo = (amrex::Real(1) + upw_lo)/amrex::Real(2) * Fhi + (amrex::Real(1) - upw_lo)/amrex::Real(2) * Flo;
252  }
253 
254  AMREX_GPU_DEVICE
255  AMREX_FORCE_INLINE
256  void
257  InterpolateInY (const int& i,
258  const int& j,
259  const int& k,
260  const int& qty_index,
261  amrex::Real& val_lo,
262  amrex::Real upw_lo) const
263  {
264  // Upwinding flags
265  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
266 
267  // Data to interpolate on
268  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
269  amrex::Real s = m_phi(i , j , k , qty_index);
270  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
271  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
272 
273  // Left and right fluxes
274  amrex::Real Fhi = Evaluate(s ,sm1,sm2);
275  amrex::Real Flo = Evaluate(sm1,s ,sp1);
276 
277  // Numerical flux
278  val_lo = (amrex::Real(1) + upw_lo)/amrex::Real(2) * Fhi + (amrex::Real(1) - upw_lo)/amrex::Real(2) * Flo;
279  }
280 
281  AMREX_GPU_DEVICE
282  AMREX_FORCE_INLINE
283  void
284  InterpolateInZ (const int& i,
285  const int& j,
286  const int& k,
287  const int& qty_index,
288  amrex::Real& val_lo,
289  amrex::Real upw_lo) const
290  {
291  // Upwinding flags
292  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
293 
294  // Data to interpolate on
295  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
296  amrex::Real s = m_phi(i , j , k , qty_index);
297  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
298  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
299 
300  // Left and right fluxes
301  amrex::Real Fhi = Evaluate(s ,sm1,sm2);
302  amrex::Real Flo = Evaluate(sm1,s ,sp1);
303 
304  // Numerical flux
305  val_lo = (amrex::Real(1) + upw_lo)/amrex::Real(2) * Fhi + (amrex::Real(1) - upw_lo)/amrex::Real(2) * Flo;
306  }
307 
308  AMREX_GPU_DEVICE
309  AMREX_FORCE_INLINE
312  const amrex::Real& sm1,
313  const amrex::Real& sm2) const
314  {
315  amrex::Real rmh = (s - sm1) / ((sm1 - sm2) + eps);
316  amrex::Real K = l2 + amrex::Real(2)*l2*rmh;
317  amrex::Real phi = amrex::max(amrex::Real(0), amrex::min(amrex::Real(2)*rmh, amrex::min(amrex::Real(2), K)));
318  return (sm1 + myhalf * phi * (sm1 - sm2));
319  }
320 
321  int GetUpwindCellNumber() const {return 2;}
322 
323 private:
324  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
325  static constexpr amrex::Real eps=(amrex::Real(1.0e-16));
326  static constexpr amrex::Real l1 =(amrex::Real(1)/amrex::Real(6.0));
327  static constexpr amrex::Real l2 =(amrex::Real(1)/three);
328 };
329 
330 
331 /**
332  * Interpolation operators used for 4th order centered scheme
333  */
334 struct CENTERED4
335 {
336  AMREX_GPU_HOST_DEVICE
337  AMREX_FORCE_INLINE
338  CENTERED4 (const amrex::Array4<const amrex::Real>& phi,
339  const amrex::Real /*upw_frac*/)
340  : m_phi(phi) {}
341 
342  AMREX_GPU_DEVICE
343  AMREX_FORCE_INLINE
344  void
345  InterpolateInX (const int& i,
346  const int& j,
347  const int& k,
348  const int& qty_index,
349  amrex::Real& val_lo,
350  amrex::Real /*upw_lo*/) const
351  {
352  // Data to interpolate on
353  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
354  amrex::Real s = m_phi(i , j , k , qty_index);
355  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
356  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
357 
358  // Interpolate lo
359  val_lo = Evaluate(sp1,s,sm1,sm2);
360  }
361 
362  AMREX_GPU_DEVICE
363  AMREX_FORCE_INLINE
364  void
365  InterpolateInY (const int& i,
366  const int& j,
367  const int& k,
368  const int& qty_index,
369  amrex::Real& val_lo,
370  amrex::Real /*upw_lo*/) const
371  {
372  // Data to interpolate on
373  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
374  amrex::Real s = m_phi(i , j , k , qty_index);
375  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
376  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
377 
378  // Interpolate lo
379  val_lo = Evaluate(sp1,s,sm1,sm2);
380  }
381 
382  AMREX_GPU_DEVICE
383  AMREX_FORCE_INLINE
384  void
385  InterpolateInZ (const int& i,
386  const int& j,
387  const int& k,
388  const int& qty_index,
389  amrex::Real& val_lo,
390  amrex::Real /*upw_lo*/) const
391  {
392  // Data to interpolate on
393  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
394  amrex::Real s = m_phi(i , j , k , qty_index);
395  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
396  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
397 
398  // Interpolate lo
399  val_lo = Evaluate(sp1,s,sm1,sm2);
400  }
401 
402  AMREX_GPU_DEVICE
403  AMREX_FORCE_INLINE
405  Evaluate (const amrex::Real& sp1,
406  const amrex::Real& s,
407  const amrex::Real& sm1,
408  const amrex::Real& sm2) const
409  {
410  // Averages and diffs
411  amrex::Real a1 = (s + sm1);
412  amrex::Real a2 = (sp1 + sm2);
413 
414  // Interpolated value
415  return ( g1*a1 - g2*a2 );
416  }
417 
418  int GetUpwindCellNumber() const {return 2;}
419 
420 private:
421  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
422  static constexpr amrex::Real g1=(amrex::Real(7.0)/amrex::Real(12.0));
423  static constexpr amrex::Real g2=(amrex::Real(1)/amrex::Real(12.0));
424 };
425 
426 /**
427  * Interpolation operators used for 5th order upwind scheme
428  */
429 struct UPWIND5
430 {
431  AMREX_GPU_HOST_DEVICE
432  AMREX_FORCE_INLINE
433  UPWIND5 (const amrex::Array4<const amrex::Real>& phi,
434  const amrex::Real upw_frac)
435  : m_phi(phi)
436  , m_upw_frac(upw_frac)
437  {}
438 
439  AMREX_GPU_DEVICE
440  AMREX_FORCE_INLINE
441  void
442  InterpolateInX (const int& i,
443  const int& j,
444  const int& k,
445  const int& qty_index,
446  amrex::Real& val_lo,
447  amrex::Real upw_lo) const
448  {
449  // Data to interpolate on
450  amrex::Real sp2 = m_phi(i+2, j , k , qty_index);
451  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
452  amrex::Real s = m_phi(i , j , k , qty_index);
453  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
454  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
455  amrex::Real sm3 = m_phi(i-3, j , k , qty_index);
456 
457  // Upwinding flags
458  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
459 
460  // Add blending
461  upw_lo *= m_upw_frac;
462 
463  // Interpolate lo
464  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
465  }
466 
467  AMREX_GPU_DEVICE
468  AMREX_FORCE_INLINE
469  void
470  InterpolateInY (const int& i,
471  const int& j,
472  const int& k,
473  const int& qty_index,
474  amrex::Real& val_lo,
475  amrex::Real upw_lo) const
476  {
477  // Data to interpolate on
478  amrex::Real sp2 = m_phi(i , j+2, k , qty_index);
479  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
480  amrex::Real s = m_phi(i , j , k , qty_index);
481  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
482  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
483  amrex::Real sm3 = m_phi(i , j-3, k , qty_index);
484 
485  // Upwinding flags
486  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
487 
488  // Add blending
489  upw_lo *= m_upw_frac;
490 
491  // Interpolate lo
492  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
493  }
494 
495  AMREX_GPU_DEVICE
496  AMREX_FORCE_INLINE
497  void
498  InterpolateInZ (const int& i,
499  const int& j,
500  const int& k,
501  const int& qty_index,
502  amrex::Real& val_lo,
503  amrex::Real upw_lo) const
504  {
505  // Data to interpolate on
506  amrex::Real sp2 = m_phi(i , j , k+2, qty_index);
507  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
508  amrex::Real s = m_phi(i , j , k , qty_index);
509  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
510  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
511  amrex::Real sm3 = m_phi(i , j , k-3, qty_index);
512 
513  // Upwinding flags
514  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
515 
516  // Add blending
517  upw_lo *= m_upw_frac;
518 
519  // Interpolate lo
520  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
521  }
522 
523  AMREX_GPU_DEVICE
524  AMREX_FORCE_INLINE
526  Evaluate (const amrex::Real& sp2,
527  const amrex::Real& sp1,
528  const amrex::Real& s,
529  const amrex::Real& sm1,
530  const amrex::Real& sm2,
531  const amrex::Real& sm3,
532  const amrex::Real& upw) const
533  {
534  // Averages and diffs
535  amrex::Real a1 = (s + sm1);
536  amrex::Real a2 = (sp1 + sm2);
537  amrex::Real a3 = (sp2 + sm3);
538  amrex::Real d1 = (s - sm1);
539  amrex::Real d2 = (sp1 - sm2);
540  amrex::Real d3 = (sp2 - sm3);
541 
542  // Interpolated value
543  return ( g1*a1 - g2*a2 + g3*a3 - upw * g3 * (d3 - amrex::Real(5.0)*d2 + amrex::Real(10.0)*d1) );
544  }
545 
546  int GetUpwindCellNumber() const {return 3;}
547 
548  AMREX_GPU_HOST_DEVICE
549  AMREX_FORCE_INLINE
550  void SetUpwinding(amrex::Real upw_frac) {m_upw_frac = upw_frac;}
551 
552 private:
553  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
554  amrex::Real m_upw_frac; // Upwinding fraction
555  static constexpr amrex::Real g1=(amrex::Real(37.0)/amrex::Real(60.0));
556  static constexpr amrex::Real g2=(amrex::Real(2)/amrex::Real(15.0));
557  static constexpr amrex::Real g3=(amrex::Real(1)/amrex::Real(60.0));
558 };
559 
560 /**
561  * Interpolation operators used for 6th order centered scheme
562  */
563 struct CENTERED6
564 {
565  AMREX_GPU_HOST_DEVICE
566  AMREX_FORCE_INLINE
567  CENTERED6 (const amrex::Array4<const amrex::Real>& phi,
568  const amrex::Real /*upw_frac*/)
569  : m_phi(phi) {}
570 
571  AMREX_GPU_DEVICE
572  AMREX_FORCE_INLINE
573  void
574  InterpolateInX (const int& i,
575  const int& j,
576  const int& k,
577  const int& qty_index,
578  amrex::Real& val_lo,
579  amrex::Real /*upw_lo*/) const
580  {
581  // Data to interpolate on
582  amrex::Real sp2 = m_phi(i+2, j , k , qty_index);
583  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
584  amrex::Real s = m_phi(i , j , k , qty_index);
585  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
586  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
587  amrex::Real sm3 = m_phi(i-3, j , k , qty_index);
588 
589  // Interpolate lo
590  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
591  }
592 
593  AMREX_GPU_DEVICE
594  AMREX_FORCE_INLINE
595  void
596  InterpolateInY (const int& i,
597  const int& j,
598  const int& k,
599  const int& qty_index,
600  amrex::Real& val_lo,
601  amrex::Real /*upw_lo*/) const
602  {
603  // Data to interpolate on
604  amrex::Real sp2 = m_phi(i , j+2, k , qty_index);
605  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
606  amrex::Real s = m_phi(i , j , k , qty_index);
607  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
608  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
609  amrex::Real sm3 = m_phi(i , j-3, k , qty_index);
610 
611  // Interpolate lo
612  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
613  }
614 
615  AMREX_GPU_DEVICE
616  AMREX_FORCE_INLINE
617  void
618  InterpolateInZ (const int& i,
619  const int& j,
620  const int& k,
621  const int& qty_index,
622  amrex::Real& val_lo,
623  amrex::Real /*upw_lo*/) const
624  {
625  // Data to interpolate on
626  amrex::Real sp2 = m_phi(i , j , k+2, qty_index);
627  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
628  amrex::Real s = m_phi(i , j , k , qty_index);
629  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
630  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
631  amrex::Real sm3 = m_phi(i , j , k-3, qty_index);
632 
633  // Interpolate lo
634  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
635  }
636 
637  AMREX_GPU_DEVICE
638  AMREX_FORCE_INLINE
640  Evaluate (const amrex::Real& sp2,
641  const amrex::Real& sp1,
642  const amrex::Real& s,
643  const amrex::Real& sm1,
644  const amrex::Real& sm2,
645  const amrex::Real& sm3) const
646  {
647  // Averages and diffs
648  amrex::Real a1 = (s + sm1);
649  amrex::Real a2 = (sp1 + sm2);
650  amrex::Real a3 = (sp2 + sm3);
651 
652  // Interpolated value
653  return ( g1*a1 - g2*a2 + g3*a3 );
654  }
655 
656  int GetUpwindCellNumber() const {return 3;}
657 
658 private:
659  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
660  static constexpr amrex::Real g1=(amrex::Real(37.0)/amrex::Real(60.0));
661  static constexpr amrex::Real g2=(amrex::Real(2)/amrex::Real(15.0));
662  static constexpr amrex::Real g3=(amrex::Real(1)/amrex::Real(60.0));
663 };
664 
665 /**
666  * Interpolation operators used for all central/upwind schemes
667  */
668 struct UPWINDALL
669 {
670  AMREX_GPU_HOST_DEVICE
671  AMREX_FORCE_INLINE
672  UPWINDALL (const amrex::Array4<const amrex::Real>& phi,
673  const amrex::Real upw_frac)
674  : m_phi(phi)
675  , m_upw_frac(upw_frac)
676  {}
677 
678  AMREX_GPU_DEVICE
679  AMREX_FORCE_INLINE
680  void
681  InterpolateInZ (const int& i,
682  const int& j,
683  const int& k,
684  const int& qty_index,
685  amrex::Real& val_lo,
686  amrex::Real upw_lo,
687  const AdvType adv_type) const
688  {
689  if (adv_type == AdvType::Centered_2nd) {
690  val_lo = myhalf * ( m_phi(i,j,k,qty_index) + m_phi(i,j,k-1,qty_index) );
691  return;
692  } else {
693  // Data to interpolate on
694  amrex::Real sp2 = (adv_type == AdvType::Upwind_5th || adv_type == AdvType::Centered_6th ) ? m_phi(i , j , k+2, qty_index): amrex::Real(0);
695  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
696  amrex::Real s = m_phi(i , j , k , qty_index);
697  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
698  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
699  amrex::Real sm3 = (adv_type == AdvType::Upwind_5th || adv_type == AdvType::Centered_6th ) ? m_phi(i , j , k-3, qty_index) : amrex::Real(0);
700 
701  // Upwinding flags
702  if (upw_lo != amrex::Real(0)) upw_lo = (upw_lo > 0) ? amrex::Real(1) : -amrex::Real(1);
703 
704  // Add blending
705  upw_lo *= m_upw_frac;
706 
707  // Interpolate lo
708  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo,adv_type);
709  }
710  }
711 
712  AMREX_GPU_DEVICE
713  AMREX_FORCE_INLINE
715  Evaluate (const amrex::Real& sp2,
716  const amrex::Real& sp1,
717  const amrex::Real& s,
718  const amrex::Real& sm1,
719  const amrex::Real& sm2,
720  const amrex::Real& sm3,
721  const amrex::Real& upw,
722  const AdvType adv_type) const
723  {
724  // Averages and diffs
725  amrex::Real a1 = (s + sm1);
726  amrex::Real a2 = (sp1 + sm2);
727  amrex::Real d1 = (s - sm1);
728  amrex::Real d2 = (sp1 - sm2);
729  amrex::Real a3 = (sp2 + sm3);
730  amrex::Real d3 = (sp2 - sm3);
731 
732  // Interpolated value
733  amrex::Real iv(amrex::Real(0));
734  if (adv_type == AdvType::Centered_2nd) {
735  iv = myhalf * a1;
736  } else if (adv_type == AdvType::Upwind_3rd) {
737  iv = g1_3_4*a1 - g2_3_4*a2 + upw * g2_3_4 * (d2 - three*d1);
738  } else if (adv_type == AdvType::Centered_4th) {
739  iv = g1_3_4*a1 - g2_3_4*a2;
740  } else if (adv_type == AdvType::Upwind_5th) {
741  iv = g1_5_6*a1 - g2_5_6*a2 + g3_5_6*a3 - upw * g3_5_6 * (d3 - amrex::Real(5.0)*d2 + amrex::Real(10.0)*d1);
742  } else if (adv_type == AdvType::Centered_6th) {
743  iv = g1_5_6*a1 - g2_5_6*a2 + g3_5_6*a3;
744  }
745  return ( iv );
746  }
747 
748 private:
749  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
750  amrex::Real m_upw_frac; // Upwinding fraction
751  static constexpr amrex::Real g1_3_4=( amrex::Real(7.0)/amrex::Real(12.0));
752  static constexpr amrex::Real g2_3_4=( amrex::Real(1)/amrex::Real(12.0));
753  static constexpr amrex::Real g1_5_6=(amrex::Real(37.0)/amrex::Real(60.0));
754  static constexpr amrex::Real g2_5_6=( amrex::Real(2)/amrex::Real(15.0));
755  static constexpr amrex::Real g3_5_6=( amrex::Real(1)/amrex::Real(60.0));
756 };
757 #endif
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
AdvType
Definition: ERF_IndexDefines.H:255
@ Centered_4th
@ Centered_6th
@ Centered_2nd
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_UPW.H:10
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &s, const amrex::Real &sm1) const
Definition: ERF_Interpolation_UPW.H:74
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) const
Definition: ERF_Interpolation_UPW.H:38
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:84
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) const
Definition: ERF_Interpolation_UPW.H:56
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:87
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:88
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) const
Definition: ERF_Interpolation_UPW.H:20
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CENTERED2(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_UPW.H:13
Definition: ERF_Interpolation_UPW.H:335
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) const
Definition: ERF_Interpolation_UPW.H:385
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) const
Definition: ERF_Interpolation_UPW.H:365
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:421
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:422
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:418
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) const
Definition: ERF_Interpolation_UPW.H:345
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sp1, const amrex::Real &s, const amrex::Real &sm1, const amrex::Real &sm2) const
Definition: ERF_Interpolation_UPW.H:405
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CENTERED4(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_UPW.H:338
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:423
Definition: ERF_Interpolation_UPW.H:564
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) const
Definition: ERF_Interpolation_UPW.H:574
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:661
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CENTERED6(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_UPW.H:567
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) const
Definition: ERF_Interpolation_UPW.H:596
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:656
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:660
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:659
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) const
Definition: ERF_Interpolation_UPW.H:618
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sp2, const amrex::Real &sp1, const amrex::Real &s, const amrex::Real &sm1, const amrex::Real &sm2, const amrex::Real &sm3) const
Definition: ERF_Interpolation_UPW.H:640
static constexpr amrex::Real g3
Definition: ERF_Interpolation_UPW.H:662
Definition: ERF_Interpolation_UPW.H:219
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_UPW.H:230
static constexpr amrex::Real eps
Definition: ERF_Interpolation_UPW.H:325
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:321
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &s, const amrex::Real &sm1, const amrex::Real &sm2) const
Definition: ERF_Interpolation_UPW.H:311
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE UPWIND3SL(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_UPW.H:222
static constexpr amrex::Real l1
Definition: ERF_Interpolation_UPW.H:326
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_UPW.H:284
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_UPW.H:257
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:324
static constexpr amrex::Real l2
Definition: ERF_Interpolation_UPW.H:327
Definition: ERF_Interpolation_UPW.H:95
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sp1, const amrex::Real &s, const amrex::Real &sm1, const amrex::Real &sm2, const amrex::Real &upw) const
Definition: ERF_Interpolation_UPW.H:185
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:211
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void SetUpwinding(amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:205
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_UPW.H:159
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:208
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:210
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE UPWIND3(const amrex::Array4< const amrex::Real > &phi, const amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:98
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_UPW.H:133
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:201
amrex::Real m_upw_frac
Definition: ERF_Interpolation_UPW.H:209
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_UPW.H:107
Definition: ERF_Interpolation_UPW.H:430
amrex::Real m_upw_frac
Definition: ERF_Interpolation_UPW.H:554
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE UPWIND5(const amrex::Array4< const amrex::Real > &phi, const amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:433
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sp2, const amrex::Real &sp1, const amrex::Real &s, const amrex::Real &sm1, const amrex::Real &sm2, const amrex::Real &sm3, const amrex::Real &upw) const
Definition: ERF_Interpolation_UPW.H:526
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:553
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:556
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_UPW.H:498
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void SetUpwinding(amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:550
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_UPW.H:470
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_UPW.H:442
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:546
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:555
static constexpr amrex::Real g3
Definition: ERF_Interpolation_UPW.H:557
Definition: ERF_Interpolation_UPW.H:669
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:749
static constexpr amrex::Real g1_5_6
Definition: ERF_Interpolation_UPW.H:753
static constexpr amrex::Real g1_3_4
Definition: ERF_Interpolation_UPW.H:751
amrex::Real m_upw_frac
Definition: ERF_Interpolation_UPW.H:750
static constexpr amrex::Real g2_3_4
Definition: ERF_Interpolation_UPW.H:752
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Evaluate(const amrex::Real &sp2, const amrex::Real &sp1, const amrex::Real &s, const amrex::Real &sm1, const amrex::Real &sm2, const amrex::Real &sm3, const amrex::Real &upw, const AdvType adv_type) const
Definition: ERF_Interpolation_UPW.H:715
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE UPWINDALL(const amrex::Array4< const amrex::Real > &phi, const amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:672
static constexpr amrex::Real g3_5_6
Definition: ERF_Interpolation_UPW.H:755
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 AdvType adv_type) const
Definition: ERF_Interpolation_UPW.H:681
static constexpr amrex::Real g2_5_6
Definition: ERF_Interpolation_UPW.H:754