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  CENTERED2 (const amrex::Array4<const amrex::Real>& phi)
12  : m_phi(phi) {}
13 
14  AMREX_GPU_DEVICE
15  AMREX_FORCE_INLINE
16  void
17  InterpolateInX (const int& i,
18  const int& j,
19  const int& k,
20  const int& qty_index,
21  amrex::Real& val_lo,
22  amrex::Real /*upw_lo*/,
23  const amrex::Real /*upw_frac*/) 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*/,
42  const amrex::Real /*upw_frac*/) const
43  {
44  // Data to interpolate on
45  amrex::Real s = m_phi(i , j , k , qty_index);
46  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
47 
48  // Interpolate lo
49  val_lo = Evaluate(s,sm1);
50  }
51 
52  AMREX_GPU_DEVICE
53  AMREX_FORCE_INLINE
54  void
55  InterpolateInZ (const int& i,
56  const int& j,
57  const int& k,
58  const int& qty_index,
59  amrex::Real& val_lo,
60  amrex::Real /*upw_lo*/,
61  const amrex::Real /*upw_frac*/) 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
73  amrex::Real
74  Evaluate(const amrex::Real& s,
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=(0.5);
89 };
90 
91 /**
92  * Interpolation operators used for 3rd order upwind scheme
93  */
94 struct UPWIND3
95 {
96  UPWIND3 (const amrex::Array4<const amrex::Real>& phi)
97  : m_phi(phi) {}
98 
99  AMREX_GPU_DEVICE
100  AMREX_FORCE_INLINE
101  void
102  InterpolateInX (const int& i,
103  const int& j,
104  const int& k,
105  const int& qty_index,
106  amrex::Real& val_lo,
107  amrex::Real upw_lo,
108  const amrex::Real upw_frac) 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 *= 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,
135  const amrex::Real upw_frac) const
136  {
137  // Data to interpolate on
138  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
139  amrex::Real s = m_phi(i , j , k , qty_index);
140  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
141  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
142 
143  // Upwinding flags
144  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
145 
146  // Add blending
147  upw_lo *= upw_frac;
148 
149  // Interpolate lo
150  val_lo = Evaluate(sp1,s,sm1,sm2,upw_lo);
151  }
152 
153  AMREX_GPU_DEVICE
154  AMREX_FORCE_INLINE
155  void
156  InterpolateInZ (const int& i,
157  const int& j,
158  const int& k,
159  const int& qty_index,
160  amrex::Real& val_lo,
161  amrex::Real upw_lo,
162  const amrex::Real upw_frac) const
163  {
164  // Data to interpolate on
165  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
166  amrex::Real s = m_phi(i , j , k , qty_index);
167  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
168  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
169 
170  // Upwinding flags
171  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
172 
173  // Add blending
174  upw_lo *= upw_frac;
175 
176  // Interpolate lo
177  val_lo = Evaluate(sp1,s,sm1,sm2,upw_lo);
178  }
179 
180  AMREX_GPU_DEVICE
181  AMREX_FORCE_INLINE
182  amrex::Real
183  Evaluate (const amrex::Real& sp1,
184  const amrex::Real& s,
185  const amrex::Real& sm1,
186  const amrex::Real& sm2,
187  const amrex::Real& upw) const
188  {
189  // Averages and diffs
190  amrex::Real a1 = (s + sm1);
191  amrex::Real d1 = (s - sm1);
192  amrex::Real a2 = (sp1 + sm2);
193  amrex::Real d2 = (sp1 - sm2);
194 
195  // Interpolated value
196  return ( g1*a1 - g2*a2 + upw * g2 * (d2 - 3.0*d1) );
197  }
198 
199  int GetUpwindCellNumber() const {return 2;}
200 
201 private:
202  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
203  static constexpr amrex::Real g1=(7.0/12.0);
204  static constexpr amrex::Real g2=(1.0/12.0);
205 };
206 
207 
208 /**
209  * Interpolation operators used for 4th order centered scheme
210  */
211 struct CENTERED4
212 {
213  CENTERED4 (const amrex::Array4<const amrex::Real>& phi)
214  : m_phi(phi) {}
215 
216  AMREX_GPU_DEVICE
217  AMREX_FORCE_INLINE
218  void
219  InterpolateInX (const int& i,
220  const int& j,
221  const int& k,
222  const int& qty_index,
223  amrex::Real& val_lo,
224  amrex::Real /*upw_lo*/,
225  const amrex::Real /*upw_frac*/) const
226  {
227  // Data to interpolate on
228  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
229  amrex::Real s = m_phi(i , j , k , qty_index);
230  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
231  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
232 
233  // Interpolate lo
234  val_lo = Evaluate(sp1,s,sm1,sm2);
235  }
236 
237  AMREX_GPU_DEVICE
238  AMREX_FORCE_INLINE
239  void
240  InterpolateInY (const int& i,
241  const int& j,
242  const int& k,
243  const int& qty_index,
244  amrex::Real& val_lo,
245  amrex::Real /*upw_lo*/,
246  const amrex::Real /*upw_frac*/) 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*/,
267  const amrex::Real /*upw_frac*/) const
268  {
269  // Data to interpolate on
270  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
271  amrex::Real s = m_phi(i , j , k , qty_index);
272  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
273  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
274 
275  // Interpolate lo
276  val_lo = Evaluate(sp1,s,sm1,sm2);
277  }
278 
279  AMREX_GPU_DEVICE
280  AMREX_FORCE_INLINE
281  amrex::Real
282  Evaluate (const amrex::Real& sp1,
283  const amrex::Real& s,
284  const amrex::Real& sm1,
285  const amrex::Real& sm2) const
286  {
287  // Averages and diffs
288  amrex::Real a1 = (s + sm1);
289  amrex::Real a2 = (sp1 + sm2);
290 
291  // Interpolated value
292  return ( g1*a1 - g2*a2 );
293  }
294 
295  int GetUpwindCellNumber() const {return 2;}
296 
297 private:
298  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
299  static constexpr amrex::Real g1=(7.0/12.0);
300  static constexpr amrex::Real g2=(1.0/12.0);
301 };
302 
303 /**
304  * Interpolation operators used for 5th order upwind scheme
305  */
306 struct UPWIND5
307 {
308  UPWIND5 (const amrex::Array4<const amrex::Real>& phi)
309  : m_phi(phi) {}
310 
311  AMREX_GPU_DEVICE
312  AMREX_FORCE_INLINE
313  void
314  InterpolateInX (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,
320  const amrex::Real upw_frac) const
321  {
322  // Data to interpolate on
323  amrex::Real sp2 = m_phi(i+2, j , k , qty_index);
324  amrex::Real sp1 = m_phi(i+1, j , k , qty_index);
325  amrex::Real s = m_phi(i , j , k , qty_index);
326  amrex::Real sm1 = m_phi(i-1, j , k , qty_index);
327  amrex::Real sm2 = m_phi(i-2, j , k , qty_index);
328  amrex::Real sm3 = m_phi(i-3, j , k , qty_index);
329 
330  // Upwinding flags
331  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
332 
333  // Add blending
334  upw_lo *= upw_frac;
335 
336  // Interpolate lo
337  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
338  }
339 
340  AMREX_GPU_DEVICE
341  AMREX_FORCE_INLINE
342  void
343  InterpolateInY (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,
349  const amrex::Real upw_frac) 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 *= 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,
378  const amrex::Real upw_frac) const
379  {
380  // Data to interpolate on
381  amrex::Real sp2 = m_phi(i , j , k+2, qty_index);
382  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
383  amrex::Real s = m_phi(i , j , k , qty_index);
384  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
385  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
386  amrex::Real sm3 = m_phi(i , j , k-3, qty_index);
387 
388  // Upwinding flags
389  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
390 
391  // Add blending
392  upw_lo *= upw_frac;
393 
394  // Interpolate lo
395  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo);
396  }
397 
398  AMREX_GPU_DEVICE
399  AMREX_FORCE_INLINE
400  amrex::Real
401  Evaluate (const amrex::Real& sp2,
402  const amrex::Real& sp1,
403  const amrex::Real& s,
404  const amrex::Real& sm1,
405  const amrex::Real& sm2,
406  const amrex::Real& sm3,
407  const amrex::Real& upw) const
408  {
409  // Averages and diffs
410  amrex::Real a1 = (s + sm1);
411  amrex::Real a2 = (sp1 + sm2);
412  amrex::Real a3 = (sp2 + sm3);
413  amrex::Real d1 = (s - sm1);
414  amrex::Real d2 = (sp1 - sm2);
415  amrex::Real d3 = (sp2 - sm3);
416 
417  // Interpolated value
418  return ( g1*a1 - g2*a2 + g3*a3 - upw * g3 * (d3 - 5.0*d2 + 10.0*d1) );
419  }
420 
421  int GetUpwindCellNumber() const {return 3;}
422 
423 private:
424  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
425  static constexpr amrex::Real g1=(37.0/60.0);
426  static constexpr amrex::Real g2=(2.0/15.0);
427  static constexpr amrex::Real g3=(1.0/60.0);
428 };
429 
430 /**
431  * Interpolation operators used for 6th order centered scheme
432  */
433 struct CENTERED6
434 {
435  CENTERED6 (const amrex::Array4<const amrex::Real>& phi)
436  : m_phi(phi) {}
437 
438  AMREX_GPU_DEVICE
439  AMREX_FORCE_INLINE
440  void
441  InterpolateInX (const int& i,
442  const int& j,
443  const int& k,
444  const int& qty_index,
445  amrex::Real& val_lo,
446  amrex::Real /*upw_lo*/,
447  const amrex::Real /*upw_frac*/) 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  // Interpolate lo
458  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
459  }
460 
461  AMREX_GPU_DEVICE
462  AMREX_FORCE_INLINE
463  void
464  InterpolateInY (const int& i,
465  const int& j,
466  const int& k,
467  const int& qty_index,
468  amrex::Real& val_lo,
469  amrex::Real /*upw_lo*/,
470  const amrex::Real /*upw_frac*/) const
471  {
472  // Data to interpolate on
473  amrex::Real sp2 = m_phi(i , j+2, k , qty_index);
474  amrex::Real sp1 = m_phi(i , j+1, k , qty_index);
475  amrex::Real s = m_phi(i , j , k , qty_index);
476  amrex::Real sm1 = m_phi(i , j-1, k , qty_index);
477  amrex::Real sm2 = m_phi(i , j-2, k , qty_index);
478  amrex::Real sm3 = m_phi(i , j-3, k , qty_index);
479 
480  // Interpolate lo
481  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3);
482  }
483 
484  AMREX_GPU_DEVICE
485  AMREX_FORCE_INLINE
486  void
487  InterpolateInZ (const int& i,
488  const int& j,
489  const int& k,
490  const int& qty_index,
491  amrex::Real& val_lo,
492  amrex::Real /*upw_lo*/,
493  const amrex::Real /*upw_frac*/) 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  : m_phi(phi) {}
542 
543  AMREX_GPU_DEVICE
544  AMREX_FORCE_INLINE
545  void
546  InterpolateInZ (const int& i,
547  const int& j,
548  const int& k,
549  const int& qty_index,
550  amrex::Real& val_lo,
551  amrex::Real upw_lo,
552  const amrex::Real upw_frac,
553  const AdvType adv_type) const
554  {
555  if (adv_type == AdvType::Centered_2nd) {
556  val_lo = 0.5 * ( m_phi(i,j,k,qty_index) + m_phi(i,j,k-1,qty_index) );
557  return;
558  } else {
559  // Data to interpolate on
560  amrex::Real sp2 = (adv_type == AdvType::Upwind_5th || adv_type == AdvType::Centered_6th ) ? m_phi(i , j , k+2, qty_index): 0.;
561  amrex::Real sp1 = m_phi(i , j , k+1, qty_index);
562  amrex::Real s = m_phi(i , j , k , qty_index);
563  amrex::Real sm1 = m_phi(i , j , k-1, qty_index);
564  amrex::Real sm2 = m_phi(i , j , k-2, qty_index);
565  amrex::Real sm3 = (adv_type == AdvType::Upwind_5th || adv_type == AdvType::Centered_6th ) ? m_phi(i , j , k-3, qty_index) : 0.;
566 
567  // Upwinding flags
568  if (upw_lo != 0.) upw_lo = (upw_lo > 0) ? 1. : -1.;
569 
570  // Add blending
571  upw_lo *= upw_frac;
572 
573  // Interpolate lo
574  val_lo = Evaluate(sp2,sp1,s,sm1,sm2,sm3,upw_lo,adv_type);
575  }
576  }
577 
578  AMREX_GPU_DEVICE
579  AMREX_FORCE_INLINE
580  amrex::Real
581  Evaluate (const amrex::Real& sp2,
582  const amrex::Real& sp1,
583  const amrex::Real& s,
584  const amrex::Real& sm1,
585  const amrex::Real& sm2,
586  const amrex::Real& sm3,
587  const amrex::Real& upw,
588  const AdvType adv_type) const
589  {
590  // Averages and diffs
591  amrex::Real a1 = (s + sm1);
592  amrex::Real a2 = (sp1 + sm2);
593  amrex::Real d1 = (s - sm1);
594  amrex::Real d2 = (sp1 - sm2);
595  amrex::Real a3 = (sp2 + sm3);
596  amrex::Real d3 = (sp2 - sm3);
597 
598  // Interpolated value
599  amrex::Real iv(0.0);
600  if (adv_type == AdvType::Centered_2nd) {
601  iv = 0.5 * a1;
602  } else if (adv_type == AdvType::Upwind_3rd) {
603  iv = g1_3_4*a1 - g2_3_4*a2 + upw * g2_3_4 * (d2 - 3.0*d1);
604  } else if (adv_type == AdvType::Centered_4th) {
605  iv = g1_3_4*a1 - g2_3_4*a2;
606  } else if (adv_type == AdvType::Upwind_5th) {
607  iv = g1_5_6*a1 - g2_5_6*a2 + g3_5_6*a3 - upw * g3_5_6 * (d3 - 5.0*d2 + 10.0*d1);
608  } else if (adv_type == AdvType::Centered_6th) {
609  iv = g1_5_6*a1 - g2_5_6*a2 + g3_5_6*a3;
610  }
611  return ( iv );
612  }
613 
614 private:
615  amrex::Array4<const amrex::Real> m_phi; // Quantity to interpolate
616  static constexpr amrex::Real g1_3_4=( 7.0/12.0);
617  static constexpr amrex::Real g2_3_4=( 1.0/12.0);
618  static constexpr amrex::Real g1_5_6=(37.0/60.0);
619  static constexpr amrex::Real g2_5_6=( 2.0/15.0);
620  static constexpr amrex::Real g3_5_6=( 1.0/60.0);
621 };
622 #endif
AdvType
Definition: ERF_IndexDefines.H:191
@ Centered_4th
@ Centered_6th
@ Centered_2nd
Definition: ERF_Interpolation_UPW.H:10
CENTERED2(const amrex::Array4< const amrex::Real > &phi)
Definition: ERF_Interpolation_UPW.H:11
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 InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real, const amrex::Real) const
Definition: ERF_Interpolation_UPW.H:55
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:84
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 amrex::Real) const
Definition: ERF_Interpolation_UPW.H:36
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 amrex::Real) const
Definition: ERF_Interpolation_UPW.H:17
Definition: ERF_Interpolation_UPW.H:212
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 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 amrex::Real) const
Definition: ERF_Interpolation_UPW.H:240
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:298
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:299
CENTERED4(const amrex::Array4< const amrex::Real > &phi)
Definition: ERF_Interpolation_UPW.H:213
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:295
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 amrex::Real) const
Definition: ERF_Interpolation_UPW.H:219
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:282
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:300
Definition: ERF_Interpolation_UPW.H:434
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 amrex::Real) const
Definition: ERF_Interpolation_UPW.H:464
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:531
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 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)
Definition: ERF_Interpolation_UPW.H:435
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 amrex::Real) const
Definition: ERF_Interpolation_UPW.H:487
static constexpr amrex::Real g3
Definition: ERF_Interpolation_UPW.H:532
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 amrex::Real) const
Definition: ERF_Interpolation_UPW.H:441
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:183
static constexpr amrex::Real g2
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 amrex::Real upw_frac) const
Definition: ERF_Interpolation_UPW.H:129
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:202
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:203
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:199
UPWIND3(const amrex::Array4< const amrex::Real > &phi)
Definition: ERF_Interpolation_UPW.H:96
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 amrex::Real upw_frac) const
Definition: ERF_Interpolation_UPW.H:102
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 amrex::Real upw_frac) const
Definition: ERF_Interpolation_UPW.H:156
Definition: ERF_Interpolation_UPW.H:307
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:401
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:424
static constexpr amrex::Real g2
Definition: ERF_Interpolation_UPW.H:426
int GetUpwindCellNumber() const
Definition: ERF_Interpolation_UPW.H:421
static constexpr amrex::Real g1
Definition: ERF_Interpolation_UPW.H:425
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 amrex::Real upw_frac) const
Definition: ERF_Interpolation_UPW.H:372
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 amrex::Real upw_frac) const
Definition: ERF_Interpolation_UPW.H:314
static constexpr amrex::Real g3
Definition: ERF_Interpolation_UPW.H:427
UPWIND5(const amrex::Array4< const amrex::Real > &phi)
Definition: ERF_Interpolation_UPW.H:308
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 amrex::Real upw_frac) const
Definition: ERF_Interpolation_UPW.H:343
Definition: ERF_Interpolation_UPW.H:539
UPWINDALL(const amrex::Array4< const amrex::Real > &phi)
Definition: ERF_Interpolation_UPW.H:540
amrex::Array4< const amrex::Real > m_phi
Definition: ERF_Interpolation_UPW.H:615
static constexpr amrex::Real g1_5_6
Definition: ERF_Interpolation_UPW.H:618
static constexpr amrex::Real g1_3_4
Definition: ERF_Interpolation_UPW.H:616
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 amrex::Real upw_frac, const AdvType adv_type) const
Definition: ERF_Interpolation_UPW.H:546
static constexpr amrex::Real g2_3_4
Definition: ERF_Interpolation_UPW.H:617
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:581
static constexpr amrex::Real g3_5_6
Definition: ERF_Interpolation_UPW.H:620
static constexpr amrex::Real g2_5_6
Definition: ERF_Interpolation_UPW.H:619