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 4th order centered scheme
211  */
212 struct CENTERED4
213 {
214  CENTERED4 (const amrex::Array4<const amrex::Real>& phi,
215  const amrex::Real /*upw_frac*/)
216  : m_phi(phi) {}
217 
218  AMREX_GPU_DEVICE
219  AMREX_FORCE_INLINE
220  void
221  InterpolateInX (const int& i,
222  const int& j,
223  const int& k,
224  const int& qty_index,
225  amrex::Real& val_lo,
226  amrex::Real /*upw_lo*/) const
227  {
228  // Data to interpolate on
229  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
230  amrex::Real s = m_phi(i , j , k , qty_index);
231  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
232  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
233 
234  // Interpolate lo
235  val_lo = Evaluate(sp1,s,sm1,sm2);
236  }
237 
238  AMREX_GPU_DEVICE
239  AMREX_FORCE_INLINE
240  void
241  InterpolateInY (const int& i,
242  const int& j,
243  const int& k,
244  const int& qty_index,
245  amrex::Real& val_lo,
246  amrex::Real /*upw_lo*/) const
247  {
248  // Data to interpolate on
249  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
250  amrex::Real s = m_phi(i , j , k , qty_index);
251  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
252  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
253 
254  // Interpolate lo
255  val_lo = Evaluate(sp1,s,sm1,sm2);
256  }
257 
258  AMREX_GPU_DEVICE
259  AMREX_FORCE_INLINE
260  void
261  InterpolateInZ (const int& i,
262  const int& j,
263  const int& k,
264  const int& qty_index,
265  amrex::Real& val_lo,
266  amrex::Real /*upw_lo*/) const
267  {
268  // Data to interpolate on
269  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
270  amrex::Real s = m_phi(i , j , k , qty_index);
271  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
272  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
273 
274  // Interpolate lo
275  val_lo = Evaluate(sp1,s,sm1,sm2);
276  }
277 
278  AMREX_GPU_DEVICE
279  AMREX_FORCE_INLINE
280  amrex::Real
281  Evaluate (const amrex::Real& sp1,
282  const amrex::Real& s,
283  const amrex::Real& sm1,
284  const amrex::Real& sm2) const
285  {
286  // Averages and diffs
287  amrex::Real a1 = (s + sm1);
288  amrex::Real a2 = (sp1 + sm2);
289 
290  // Interpolated value
291  return ( g1*a1 - g2*a2 );
292  }
293 
294  int GetUpwindCellNumber() const {return 2;}
295 
296 private:
297  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
298  static constexpr amrex::Real g1=(7.0/12.0);
299  static constexpr amrex::Real g2=(1.0/12.0);
300 };
301 
302 /**
303  * Interpolation operators used for 5th order upwind scheme
304  */
305 struct UPWIND5
306 {
307  UPWIND5 (const amrex::Array4<const amrex::Real>& phi,
308  const amrex::Real upw_frac)
309  : m_phi(phi)
310  , m_upw_frac(upw_frac)
311  {}
312 
313  AMREX_GPU_DEVICE
314  AMREX_FORCE_INLINE
315  void
316  InterpolateInX (const int& i,
317  const int& j,
318  const int& k,
319  const int& qty_index,
320  amrex::Real& val_lo,
321  amrex::Real upw_lo) const
322  {
323  // Data to interpolate on
324  amrex::Real sp2 = m_phi(i+2, j , k , qty_index);
325  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
326  amrex::Real s = m_phi(i , j , k , qty_index);
327  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
328  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
329  amrex::Real sm3 = m_phi(i-3, j , k , qty_index);
330 
331  // Upwinding flags
332  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
333 
334  // Add blending
335  upw_lo *= m_upw_frac;
336 
337  // Interpolate lo
338  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
339  }
340 
341  AMREX_GPU_DEVICE
342  AMREX_FORCE_INLINE
343  void
344  InterpolateInY (const int& i,
345  const int& j,
346  const int& k,
347  const int& qty_index,
348  amrex::Real& val_lo,
349  amrex::Real upw_lo) const
350  {
351  // Data to interpolate on
352  amrex::Real sp2 = m_phi(i , j+2, k , qty_index);
353  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
354  amrex::Real s = m_phi(i , j , k , qty_index);
355  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
356  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
357  amrex::Real sm3 = m_phi(i , j-3, k , qty_index);
358 
359  // Upwinding flags
360  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
361 
362  // Add blending
363  upw_lo *= m_upw_frac;
364 
365  // Interpolate lo
366  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
367  }
368 
369  AMREX_GPU_DEVICE
370  AMREX_FORCE_INLINE
371  void
372  InterpolateInZ (const int& i,
373  const int& j,
374  const int& k,
375  const int& qty_index,
376  amrex::Real& val_lo,
377  amrex::Real upw_lo) const
378  {
379  // Data to interpolate on
380  amrex::Real sp2 = m_phi(i , j , k+2, qty_index);
381  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
382  amrex::Real s = m_phi(i , j , k , qty_index);
383  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
384  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
385  amrex::Real sm3 = m_phi(i , j , k-3, qty_index);
386 
387  // Upwinding flags
388  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
389 
390  // Add blending
391  upw_lo *= m_upw_frac;
392 
393  // Interpolate lo
394  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
395  }
396 
397  AMREX_GPU_DEVICE
398  AMREX_FORCE_INLINE
399  amrex::Real
400  Evaluate (const amrex::Real& sp2,
401  const amrex::Real& sp1,
402  const amrex::Real& s,
403  const amrex::Real& sm1,
404  const amrex::Real& sm2,
405  const amrex::Real& sm3,
406  const amrex::Real& upw) const
407  {
408  // Averages and diffs
409  amrex::Real a1 = (s + sm1);
410  amrex::Real a2 = (sp1 + sm2);
411  amrex::Real a3 = (sp2 + sm3);
412  amrex::Real d1 = (s - sm1);
413  amrex::Real d2 = (sp1 - sm2);
414  amrex::Real d3 = (sp2 - sm3);
415 
416  // Interpolated value
417  return ( g1*a1 - g2*a2 + g3*a3 - upw * g3 * (d3 - 5.0*d2 + 10.0*d1) );
418  }
419 
420  int GetUpwindCellNumber() const {return 3;}
421 
422  void SetUpwinding(amrex::Real upw_frac) {m_upw_frac = upw_frac;}
423 
424 private:
425  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
426  amrex::Real m_upw_frac; // Upwinding fraction
427  static constexpr amrex::Real g1=(37.0/60.0);
428  static constexpr amrex::Real g2=(2.0/15.0);
429  static constexpr amrex::Real g3=(1.0/60.0);
430 };
431 
432 /**
433  * Interpolation operators used for 6th order centered scheme
434  */
435 struct CENTERED6
436 {
437  CENTERED6 (const amrex::Array4<const amrex::Real>& phi,
438  const amrex::Real /*upw_frac*/)
439  : m_phi(phi) {}
440 
441  AMREX_GPU_DEVICE
442  AMREX_FORCE_INLINE
443  void
444  InterpolateInX (const int& i,
445  const int& j,
446  const int& k,
447  const int& qty_index,
448  amrex::Real& val_lo,
449  amrex::Real /*upw_lo*/) const
450  {
451  // Data to interpolate on
452  amrex::Real sp2 = m_phi(i+2, j , k , qty_index);
453  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
454  amrex::Real s = m_phi(i , j , k , qty_index);
455  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
456  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
457  amrex::Real sm3 = m_phi(i-3, j , k , qty_index);
458 
459  // Interpolate lo
460  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
461  }
462 
463  AMREX_GPU_DEVICE
464  AMREX_FORCE_INLINE
465  void
466  InterpolateInY (const int& i,
467  const int& j,
468  const int& k,
469  const int& qty_index,
470  amrex::Real& val_lo,
471  amrex::Real /*upw_lo*/) const
472  {
473  // Data to interpolate on
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 
481  // Interpolate lo
482  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
483  }
484 
485  AMREX_GPU_DEVICE
486  AMREX_FORCE_INLINE
487  void
488  InterpolateInZ (const int& i,
489  const int& j,
490  const int& k,
491  const int& qty_index,
492  amrex::Real& val_lo,
493  amrex::Real /*upw_lo*/) const
494  {
495  // Data to interpolate on
496  amrex::Real sp2 = m_phi(i , j , k+2, qty_index);
497  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
498  amrex::Real s = m_phi(i , j , k , qty_index);
499  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
500  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
501  amrex::Real sm3 = m_phi(i , j , k-3, qty_index);
502 
503  // Interpolate lo
504  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
505  }
506 
507  AMREX_GPU_DEVICE
508  AMREX_FORCE_INLINE
509  amrex::Real
510  Evaluate (const amrex::Real& sp2,
511  const amrex::Real& sp1,
512  const amrex::Real& s,
513  const amrex::Real& sm1,
514  const amrex::Real& sm2,
515  const amrex::Real& sm3) const
516  {
517  // Averages and diffs
518  amrex::Real a1 = (s + sm1);
519  amrex::Real a2 = (sp1 + sm2);
520  amrex::Real a3 = (sp2 + sm3);
521 
522  // Interpolated value
523  return ( g1*a1 - g2*a2 + g3*a3 );
524  }
525 
526  int GetUpwindCellNumber() const {return 3;}
527 
528 private:
529  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
530  static constexpr amrex::Real g1=(37.0/60.0);
531  static constexpr amrex::Real g2=(2.0/15.0);
532  static constexpr amrex::Real g3=(1.0/60.0);
533 };
534 
535 /**
536  * Interpolation operators used for all central/upwind schemes
537  */
538 struct UPWINDALL
539 {
540  UPWINDALL (const amrex::Array4<const amrex::Real>& phi,
541  const amrex::Real upw_frac)
542  : m_phi(phi)
543  , m_upw_frac(upw_frac)
544  {}
545 
546  AMREX_GPU_DEVICE
547  AMREX_FORCE_INLINE
548  void
549  InterpolateInZ (const int& i,
550  const int& j,
551  const int& k,
552  const int& qty_index,
553  amrex::Real& val_lo,
554  amrex::Real upw_lo,
555  const AdvType adv_type) const
556  {
557  if (adv_type == AdvType::Centered_2nd) {
558  val_lo = 0.5 * ( m_phi(i,j,k,qty_index) + m_phi(i,j,k-1,qty_index) );
559  return;
560  } else {
561  // Data to interpolate on
562  amrex::Real sp2 = (adv_type == AdvType::Upwind_5th || adv_type == AdvType::Centered_6th ) ? m_phi(i , j , k+2, qty_index): 0.;
563  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
564  amrex::Real s = m_phi(i , j , k , qty_index);
565  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
566  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
567  amrex::Real sm3 = (adv_type == AdvType::Upwind_5th || adv_type == AdvType::Centered_6th ) ? m_phi(i , j , k-3, qty_index) : 0.;
568 
569  // Upwinding flags
570  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
571 
572  // Add blending
573  upw_lo *= m_upw_frac;
574 
575  // Interpolate lo
576  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo,adv_type);
577  }
578  }
579 
580  AMREX_GPU_DEVICE
581  AMREX_FORCE_INLINE
582  amrex::Real
583  Evaluate (const amrex::Real& sp2,
584  const amrex::Real& sp1,
585  const amrex::Real& s,
586  const amrex::Real& sm1,
587  const amrex::Real& sm2,
588  const amrex::Real& sm3,
589  const amrex::Real& upw,
590  const AdvType adv_type) const
591  {
592  // Averages and diffs
593  amrex::Real a1 = (s + sm1);
594  amrex::Real a2 = (sp1 + sm2);
595  amrex::Real d1 = (s - sm1);
596  amrex::Real d2 = (sp1 - sm2);
597  amrex::Real a3 = (sp2 + sm3);
598  amrex::Real d3 = (sp2 - sm3);
599 
600  // Interpolated value
601  amrex::Real iv(0.0);
602  if (adv_type == AdvType::Centered_2nd) {
603  iv = 0.5 * a1;
604  } else if (adv_type == AdvType::Upwind_3rd) {
605  iv = g1_3_4*a1 - g2_3_4*a2 + upw * g2_3_4 * (d2 - 3.0*d1);
606  } else if (adv_type == AdvType::Centered_4th) {
607  iv = g1_3_4*a1 - g2_3_4*a2;
608  } else if (adv_type == AdvType::Upwind_5th) {
609  iv = g1_5_6*a1 - g2_5_6*a2 + g3_5_6*a3 - upw * g3_5_6 * (d3 - 5.0*d2 + 10.0*d1);
610  } else if (adv_type == AdvType::Centered_6th) {
611  iv = g1_5_6*a1 - g2_5_6*a2 + g3_5_6*a3;
612  }
613  return ( iv );
614  }
615 
616 private:
617  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
618  amrex::Real m_upw_frac; // Upwinding fraction
619  static constexpr amrex::Real g1_3_4=( 7.0/12.0);
620  static constexpr amrex::Real g2_3_4=( 1.0/12.0);
621  static constexpr amrex::Real g1_5_6=(37.0/60.0);
622  static constexpr amrex::Real g2_5_6=( 2.0/15.0);
623  static constexpr amrex::Real g3_5_6=( 1.0/60.0);
624 };
625 #endif
AdvType
Definition: ERF_IndexDefines.H:203
@ Centered_4th
@ Centered_6th
@ Centered_2nd
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:213
CENTERED4(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_UPW.H:214
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:261
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:241
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:297
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:298
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:294
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:221
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:281
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:299
Definition: ERF_Interpolation_UPW.H:436
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:444
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:531
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:466
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:526
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:530
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:529
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:488
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:510
CENTERED6(const amrex::Array4< const amrex::Real > &phi, const amrex::Real)
Definition: ERF_Interpolation_UPW.H:437
static constexpr amrex::Real g3
Definition: ERF_Interpolation_UPW.H:532
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:306
amrex::Real m_upw_frac
Definition: ERF_Interpolation_UPW.H:426
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:400
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:425
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:428
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:372
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:344
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:316
UPWIND5(const amrex::Array4< const amrex::Real > &phi, const amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:307
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:420
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:427
static constexpr amrex::Real g3
Definition: ERF_Interpolation_UPW.H:429
void SetUpwinding(amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:422
Definition: ERF_Interpolation_UPW.H:539
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:617
static constexpr amrex::Real g1_5_6
Definition: ERF_Interpolation_UPW.H:621
static constexpr amrex::Real g1_3_4
Definition: ERF_Interpolation_UPW.H:619
amrex::Real m_upw_frac
Definition: ERF_Interpolation_UPW.H:618
static constexpr amrex::Real g2_3_4
Definition: ERF_Interpolation_UPW.H:620
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:583
static constexpr amrex::Real g3_5_6
Definition: ERF_Interpolation_UPW.H:623
UPWINDALL(const amrex::Array4< const amrex::Real > &phi, const amrex::Real upw_frac)
Definition: ERF_Interpolation_UPW.H:540
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:549
static constexpr amrex::Real g2_5_6
Definition: ERF_Interpolation_UPW.H:622