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