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