AMReX-Hydro
AMReX-based hydro routines for low Mach number flows
hydro_godunov_ppm.H
Go to the documentation of this file.
1 /*
2  * \file hydro_godunov_ppm.H
3  *
4  * \addtogroup Godunov
5  * @{
6  */
7 
8 #ifndef HYDRO_PPM_GODUNOV_H
9 #define HYDRO_PPM_GODUNOV_H
10 
11 #include <AMReX_MultiFab.H>
12 #include <AMReX_BCRec.H>
13 #include <hydro_constants.H>
14 
15 namespace PPM {
16 
18 
19 struct nolimiter {
20 
21  explicit nolimiter() = default;
22 
23 
25  amrex::Real
26  sedge1(const amrex::Real sm2,
27  const amrex::Real sm1,
28  const amrex::Real s0,
29  const amrex::Real sp1,
30  const amrex::Real /*sp2*/) const
31  {
32 
33  amrex::Real d1 = sp1 - sm1;
34  amrex::Real d2 = s0 - sm2;
35 
36  amrex::Real sedge = m_half*(s0 + sm1) - m_sixth*(d1 - d2);
37  return amrex::min(amrex::max(sedge, amrex::min(s0, sm1)),amrex::max(s0,sm1));
38 
39  }
40 
42  amrex::Real
43  sedge2(const amrex::Real /*sm2*/,
44  const amrex::Real sm1,
45  const amrex::Real s0,
46  const amrex::Real sp1,
47  const amrex::Real sp2) const
48  {
49  amrex::Real d1 = sp2 - s0;
50  amrex::Real d2 = sp1 - sm1;
51 
52  amrex::Real sedge = m_half*(sp1 + s0) - m_sixth*(d1 - d2);
53  return amrex::min(amrex::max(sedge, amrex::min(s0, sp1)),amrex::max(s0,sp1));
54 
55 
56  }
57 
58  [[nodiscard]] AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE
60  sm_sp(const amrex::Real /*s0*/,
61  const amrex::Real sedge1,
62  const amrex::Real sedge2)
63  {
65  }
66 
67  const amrex::Real m_half{amrex::Real(0.5)};
68  const amrex::Real m_sixth{amrex::Real(1.0)/amrex::Real(6.0)};
69 };
70 
71 struct vanleer {
72 
73  explicit vanleer() = default;
74 
76  amrex::Real vanLeer(const amrex::Real a, const amrex::Real b, const amrex::Real c) const
77  {
78  constexpr amrex::Real small_qty_sq = amrex::Real(1.e-10) * amrex::Real(1.e-10);
79 
80  amrex::Real dsc = m_half*(b - c);
81  amrex::Real dsl = amrex::Real(amrex::Real(2.0))*(a - c);
82  amrex::Real dsr = amrex::Real(amrex::Real(2.0))*(b - a);
83  return (dsl*dsr > small_qty_sq) ?
84  std::copysign(amrex::Real(amrex::Real(1.0)), dsc)*amrex::min(amrex::Math::abs(dsc),amrex::min(amrex::Math::abs(dsl), amrex::Math::abs(dsr))) : amrex::Real(0.0);
85  }
86 
88  amrex::Real
89  sedge1(const amrex::Real sm2,
90  const amrex::Real sm1,
91  const amrex::Real s0,
92  const amrex::Real sp1,
93  const amrex::Real /*sp2*/) const
94  {
95 
96  amrex::Real d1 = vanLeer(s0,sp1,sm1);
97  amrex::Real d2 = vanLeer(sm1,s0,sm2);
98 
99  amrex::Real sedge = m_half*(s0 + sm1) - m_sixth*(d1 - d2);
100  return amrex::min(amrex::max(sedge, amrex::min(s0, sm1)),amrex::max(s0,sm1));
101 
102  }
103 
105  amrex::Real
106  sedge2(const amrex::Real /*sm2*/,
107  const amrex::Real sm1,
108  const amrex::Real s0,
109  const amrex::Real sp1,
110  const amrex::Real sp2) const
111  {
112 
113  amrex::Real d1 = vanLeer(sp1,sp2,s0);
114  amrex::Real d2 = vanLeer(s0,sp1,sm1);
115 
116  amrex::Real sedge = m_half*(sp1 + s0) - m_sixth*(d1 - d2);
117  return amrex::min(amrex::max(sedge, amrex::min(s0, sp1)),amrex::max(s0,sp1));
118 
119 
120  }
121 
122  [[nodiscard]] AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE
124  sm_sp(const amrex::Real s0,
125  const amrex::Real sedge1,
126  const amrex::Real sedge2)
127  {
128  amrex::Real sm_ = sedge1;
129  amrex::Real sp_ = sedge2;
130 
131  if ((sedge2-s0)*(s0-sedge1) < 0.0) {
132  sm_ = s0;
133  sp_ = s0;
134  } else if (amrex::Math::abs(sedge2-s0) >= amrex::Real(2.0)*amrex::Math::abs(sedge1-s0)) {
135  sp_ = amrex::Real(3.0)*s0 - amrex::Real(amrex::Real(2.0))*sedge1;
136  } else if (amrex::Math::abs(sedge1-s0) >= amrex::Real(2.0)*amrex::Math::abs(sedge2-s0)) {
137  sm_ = amrex::Real(3.0)*s0 - amrex::Real(amrex::Real(2.0))*sedge2;
138  }
139  return amrex::makeTuple(sm_,sp_);
140 
141  // this might be the correct way to do this in case the else if are both == 0.0
142  // could also change the >= to just > for each
143 // if ((sedge2-s0)*(s0-sedge1) < amrex::Real(0.0)) {
144 // return amrex::makeTuple(s0,s0);
145 // }
146 // amrex::Real sm_ = sedge1;
147 // amrex::Real sp_ = sedge2;
148 // if (amrex::Math::abs(sedge2-s0) >= amrex::Real(amrex::Real(2.0))*amrex::Math::abs(sedge1-s0)) {
149 // sp_ = amrex::Real(3.0)*s0 - amrex::Real(amrex::Real(2.0))*sedge1;
150 // }
151 // if (amrex::Math::abs(sedge1-s0) >= amrex::Real(amrex::Real(2.0))*amrex::Math::abs(sedge2-s0)) {
152 // sm_ = amrex::Real(3.0)*s0 - amrex::Real(amrex::Real(2.0))*sedge2;
153 // }
154 // return amrex::makeTuple(sm_,sp_);
155 
156  }
157 
158  const amrex::Real m_half{amrex::Real(0.5)};
159  const amrex::Real m_sixth{amrex::Real(1.0)/amrex::Real(6.0)};
160 };
161 
162 struct wenoz {
163 
164  explicit wenoz()= default;
165 
167  amrex::Real
168  sedge2(const amrex::Real sm2,
169  const amrex::Real sm1,
170  const amrex::Real s0,
171  const amrex::Real sp1,
172  const amrex::Real sp2) const
173  {
174  const amrex::Real beta1 =
175  amrex::Real(13.0) / amrex::Real(12.0) * (sm2 - amrex::Real(amrex::Real(2.0)) * sm1 + s0) * (sm2 - amrex::Real(amrex::Real(2.0)) * sm1 + s0) +
176  amrex::Real(0.25) * (sm2 - amrex::Real(4.0) * sm1 + amrex::Real(3.0) * s0) * (sm2 - amrex::Real(4.0) * sm1 + amrex::Real(3.0) * s0);
177  const amrex::Real beta2 =
178  amrex::Real(13.0) / amrex::Real(12.0) * (sm1 - amrex::Real(amrex::Real(2.0)) * s0 + sp1) * (sm1 - amrex::Real(amrex::Real(2.0)) * s0 + sp1) +
179  amrex::Real(0.25) * (sm1 - sp1) * (sm1 - sp1);
180  const amrex::Real beta3 =
181  amrex::Real(13.0) / amrex::Real(12.0) * (s0 - amrex::Real(amrex::Real(2.0)) * sp1 + sp2) * (s0 - amrex::Real(amrex::Real(2.0)) * sp1 + sp2) +
182  amrex::Real(0.25) * (amrex::Real(3.0) * s0 - amrex::Real(4.0) * sp1 + sp2) * (amrex::Real(3.0) * s0 - amrex::Real(4.0) * sp1 + sp2);
183 
184  const amrex::Real t5 = amrex::Math::abs(beta3 - beta1);
185  const amrex::Real omega1 = amrex::Real(0.1) * (amrex::Real(amrex::Real(1.0)) + t5 / (m_eps + beta1));
186  const amrex::Real omega2 = amrex::Real(0.6) * (amrex::Real(amrex::Real(1.0)) + t5 / (m_eps + beta2));
187  const amrex::Real omega3 = amrex::Real(0.3) * (amrex::Real(amrex::Real(1.0)) + t5 / (m_eps + beta3));
188 
189  const amrex::Real omega = omega1 + omega2 + omega3;
190 
191  const amrex::Real v_1 = amrex::Real(amrex::Real(2.0)) * sm2 - amrex::Real(7.0) * sm1 + amrex::Real(11.0) * s0;
192  const amrex::Real v_2 = -sm1 + amrex::Real(5.0) * s0 + amrex::Real(2.0) * sp1;
193  const amrex::Real v_3 = amrex::Real(amrex::Real(2.0)) * s0 + amrex::Real(5.0) * sp1 - sp2;
194 
195  return (omega1 * v_1 + omega2 * v_2 + omega3 * v_3) / (amrex::Real(6.0) * omega);
196  }
197 
199  amrex::Real
200  sedge1(const amrex::Real sm2,
201  const amrex::Real sm1,
202  const amrex::Real s0,
203  const amrex::Real sp1,
204  const amrex::Real sp2) const
205  {
206  return sedge2(sp2,sp1,s0,sm1,sm2); // NOLINT(readability-suspicious-call-argument)
207  }
208 
209  [[nodiscard]] AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE
211  sm_sp(const amrex::Real /*s0*/,
212  const amrex::Real sedge1,
213  const amrex::Real sedge2)
214  {
215  return amrex::makeTuple(sedge1, sedge2);
216  }
217 
218  amrex::Real m_eps{amrex::Real(1.0e-6)};
219 };
220 
221 
223 void SetXBCs ( const int i, const int j, const int k, const int n,
224  amrex::Real &sm, amrex::Real &sp,
225  amrex::Real &sedge1, amrex::Real &sedge2,
227  const int bclo, const int bchi,
228  const int domlo, const int domhi)
229 {
230  using namespace amrex;
231 
232  if (bclo == BCType::ext_dir || bclo == BCType::hoextrap)
233  {
234  if (i == domlo)
235  {
236  sedge2 = -amrex::Real(0.2)*s(domlo-1,j,k,n) + amrex::Real(0.75)*s(domlo,j, k, n)
237  +amrex::Real(amrex::Real(0.5))*s(domlo+1,j,k,n) - amrex::Real(0.05)*s(domlo+2,j,k,n);
238  sedge2 = amrex::max(sedge2, amrex::min(s(domlo+1,j,k,n), s(domlo,j,k,n)));
239  sedge2 = amrex::min(sedge2, amrex::max(s(domlo+1,j,k,n), s(domlo,j,k,n)));
240 
241  sm = s(domlo-1,j,k,n);
242  sp = sedge2;
243 
244  } else if (i == domlo+1) {
245 
246  sedge1 = -amrex::Real(0.2)*s(domlo-1,j,k,n) + amrex::Real(0.75)*s(domlo ,j,k,n)
247  + amrex::Real(amrex::Real(0.5))*s(domlo+1,j,k,n) - amrex::Real(0.05)*s(domlo+2,j,k,n);
248  sedge1 = amrex::max(sedge1, amrex::min(s(domlo+1,j,k,n), s(domlo,j,k,n)));
249  sedge1 = amrex::min(sedge1, amrex::max(s(domlo+1,j,k,n), s(domlo,j,k,n)));
250 
251  sp = sedge2;
252  sm = sedge1;
253 
254  if ( (sp - s(domlo+1,j,k,n))*(s(domlo+1,j,k,n) - sm) <= amrex::Real(0.0))
255  {
256  sp = s(domlo+1,j,k,n);
257  sm = s(domlo+1,j,k,n);
258  }
259  else if(amrex::Math::abs(sp - s(domlo+1,j,k,n)) >= amrex::Real(amrex::Real(2.0))*amrex::Math::abs(sm - s(domlo+1,j,k,n)))
260  sp = amrex::Real(3.0)*s(domlo+1,j,k,n) - amrex::Real(amrex::Real(2.0))*sm;
261  else if(amrex::Math::abs(sm - s(domlo+1,j,k,n)) >= amrex::Real(amrex::Real(2.0))*amrex::Math::abs(sp - s(domlo+1,j,k,n)))
262  sm = amrex::Real(3.0)*s(domlo+1,j,k,n) - amrex::Real(amrex::Real(2.0))*sp;
263  }
264  }
265 
266  if (bchi == BCType::ext_dir || bchi == BCType::hoextrap)
267  {
268  if (i == domhi)
269  {
270  sedge1 = -amrex::Real(0.2)*s(domhi+1,j,k,n) + amrex::Real(0.75)*s(domhi,j,k, n)
271  +amrex::Real(amrex::Real(0.5))*s(domhi-1,j,k,n) - amrex::Real(0.05)*s(domhi-2,j,k,n);
272  sedge1 = amrex::max(sedge1, amrex::min(s(domhi-1,j,k,n), s(domhi,j,k,n)));
273  sedge1 = amrex::min(sedge1, amrex::max(s(domhi-1,j,k,n), s(domhi,j,k,n)));
274 
275  sp = s(domhi+1,j,k,n);
276  sm = sedge1;
277 
278  } else if (i == domhi-1) {
279 
280  sedge2 = -amrex::Real(0.2)*s(domhi+1,j,k,n) + amrex::Real(0.75)*s(domhi ,j,k,n)
281  +amrex::Real(amrex::Real(0.5))*s(domhi-1,j,k,n) - amrex::Real(0.05)*s(domhi-2,j,k,n);
282  sedge2 = amrex::max(sedge2, amrex::min(s(domhi-1,j,k,n), s(domhi,j,k,n)));
283  sedge2 = amrex::min(sedge2, amrex::max(s(domhi-1,j,k,n), s(domhi,j,k,n)));
284 
285  sp = sedge2;
286  sm = sedge1;
287 
288  if( (sp - s(domhi-1,j,k,n))*(s(domhi-1,j,k,n) - sm) <= amrex::Real(0.0))
289  {
290  sp = s(domhi-1,j,k,n);
291  sm = s(domhi-1,j,k,n);
292  }
293  else if(amrex::Math::abs(sp - s(domhi-1,j,k,n)) >= 2.*amrex::Math::abs(sm - s(domhi-1,j,k,n)))
294  sp = amrex::Real(3.0)*s(domhi-1,j,k,n) - amrex::Real(amrex::Real(2.0))*sm;
295  else if(amrex::Math::abs(sm - s(domhi-1,j,k,n)) >= 2.*amrex::Math::abs(sp - s(domhi-1,j,k,n)))
296  sm = amrex::Real(3.0)*s(domhi-1,j,k,n) - amrex::Real(amrex::Real(2.0))*sp;
297  }
298  }
299 }
300 
302 void SetYBCs ( const int i, const int j, const int k, const int n,
303  amrex::Real &sm, amrex::Real &sp,
304  amrex::Real &sedge1, amrex::Real &sedge2,
306  const int bclo, const int bchi,
307  const int domlo, const int domhi)
308 {
309  using namespace amrex;
310 
311  if (bclo == BCType::ext_dir || bclo == BCType::hoextrap)
312  {
313  if (j == domlo)
314  {
315  sedge2 = -amrex::Real(0.2)*s(i,domlo-1,k,n) + amrex::Real(0.75)*s(i,domlo ,k,n)
316  +amrex::Real(amrex::Real(0.5))*s(i,domlo+1,k,n) - amrex::Real(0.05)*s(i,domlo+2,k,n);
317  sedge2 = amrex::max(sedge2, amrex::min(s(i,domlo+1,k,n), s(i,domlo,k,n)));
318  sedge2 = amrex::min(sedge2, amrex::max(s(i,domlo+1,k,n), s(i,domlo,k,n)));
319 
320  sm = s(i,domlo-1,k,n);
321  sp = sedge2;
322 
323  } else if (j == domlo+1) {
324 
325  sedge1 = -amrex::Real(0.2)*s(i,domlo-1,k,n) + amrex::Real(0.75)*s(i,domlo ,k,n)
326  +amrex::Real(amrex::Real(0.5))*s(i,domlo+1,k,n) - amrex::Real(0.05)*s(i,domlo+2,k,n);
327  sedge1 = amrex::max(sedge1, amrex::min(s(i,domlo+1,k,n), s(i,domlo,k,n)));
328  sedge1 = amrex::min(sedge1, amrex::max(s(i,domlo+1,k,n), s(i,domlo,k,n)));
329 
330  sp = sedge2;
331  sm = sedge1;
332 
333  if ( (sp - s(i,domlo+1,k,n))*(s(i,domlo+1,k,n) - sm) <= amrex::Real(0.0))
334  {
335  sp = s(i,domlo+1,k,n);
336  sm = s(i,domlo+1,k,n);
337  }
338  else if(amrex::Math::abs(sp - s(i,domlo+1,k,n)) >= amrex::Real(amrex::Real(2.0))*amrex::Math::abs(sm - s(i,domlo+1,k,n)))
339  sp = amrex::Real(3.0)*s(i,domlo+1,k,n) - amrex::Real(amrex::Real(2.0))*sm;
340  else if(amrex::Math::abs(sm - s(i,domlo+1,k,n)) >= amrex::Real(amrex::Real(2.0))*amrex::Math::abs(sp - s(i,domlo+1,k,n)))
341  sm = amrex::Real(3.0)*s(i,domlo+1,k,n) - amrex::Real(amrex::Real(2.0))*sp;
342  }
343  }
344 
345  if (bchi == BCType::ext_dir || bchi == BCType::hoextrap)
346  {
347  if (j == domhi)
348  {
349  sedge1 = -amrex::Real(0.2)*s(i,domhi+1,k,n) + amrex::Real(0.75)*s(i,domhi ,k,n)
350  +amrex::Real(amrex::Real(0.5))*s(i,domhi-1,k,n) - amrex::Real(0.05)*s(i,domhi-2,k,n);
351  sedge1 = amrex::max(sedge1, amrex::min(s(i,domhi-1,k,n), s(i,domhi,k,n)));
352  sedge1 = amrex::min(sedge1, amrex::max(s(i,domhi-1,k,n), s(i,domhi,k,n)));
353 
354  sp = s(i,domhi+1,k, n);
355  sm = sedge1;
356 
357  } else if (j == domhi-1) {
358 
359  sedge2 = -amrex::Real(0.2)*s(i,domhi+1,k,n) + amrex::Real(0.75)*s(i,domhi ,k,n)
360  +amrex::Real(amrex::Real(0.5))*s(i,domhi-1,k,n) - amrex::Real(0.05)*s(i,domhi-2,k,n);
361  sedge2 = amrex::max(sedge2, amrex::min(s(i,domhi-1,k,n), s(i,domhi,k,n)));
362  sedge2 = amrex::min(sedge2, amrex::max(s(i,domhi-1,k,n), s(i,domhi,k,n)));
363 
364  sp = sedge2;
365  sm = sedge1;
366 
367  if( (sp - s(i,domhi-1,k,n))*(s(i,domhi-1,k,n) - sm) <= amrex::Real(0.0)){
368  sp = s(i,domhi-1,k,n);
369  sm = s(i,domhi-1,k,n);
370  }
371  else if(amrex::Math::abs(sp - s(i,domhi-1,k,n)) >= 2.*amrex::Math::abs(sm - s(i,domhi-1,k,n)))
372  sp = amrex::Real(3.0)*s(i,domhi-1,k,n) - amrex::Real(amrex::Real(2.0))*sm;
373 
374  else if(amrex::Math::abs(sm - s(i,domhi-1,k,n)) >= 2.*amrex::Math::abs(sp - s(i,domhi-1,k,n)))
375  sm = amrex::Real(3.0)*s(i,domhi-1,k,n) - amrex::Real(amrex::Real(2.0))*sp;
376  }
377  }
378 }
379 
380 #if (AMREX_SPACEDIM==3)
382 void SetZBCs ( const int i, const int j, const int k, const int n,
383  amrex::Real &sm, amrex::Real &sp,
384  amrex::Real &sedge1, amrex::Real &sedge2,
386  const int bclo, const int bchi,
387  const int domlo, const int domhi)
388 {
389  using namespace amrex;
390 
391  if (bclo == BCType::ext_dir || bclo == BCType::hoextrap)
392  {
393 
394  if (k == domlo)
395  {
396  sedge2 = -amrex::Real(0.2)*s(i,j,domlo-1,n) + amrex::Real(0.75)*s(i,j,domlo ,n)
397  +amrex::Real(amrex::Real(0.5))*s(i,j,domlo+1,n) - amrex::Real(0.05)*s(i,j,domlo+2,n);
398  sedge2 = amrex::max(sedge2, amrex::min(s(i,j,domlo+1,n), s(i,j,domlo,n)));
399  sedge2 = amrex::min(sedge2, amrex::max(s(i,j,domlo+1,n), s(i,j,domlo,n)));
400 
401  sm = s(i,j,domlo-1,n);
402  sp = sedge2;
403 
404  } else if (k == domlo+1) {
405 
406  sedge1 = -amrex::Real(0.2)*s(i,j,domlo-1,n) + amrex::Real(0.75)*s(i,j,domlo ,n)
407  +amrex::Real(amrex::Real(0.5))*s(i,j,domlo+1,n) - amrex::Real(0.05)*s(i,j,domlo+2,n);
408  sedge1 = amrex::max(sedge1, amrex::min(s(i,j,domlo+1,n), s(i,j,domlo,n)));
409  sedge1 = amrex::min(sedge1, amrex::max(s(i,j,domlo+1,n), s(i,j,domlo,n)));
410 
411  sp = sedge2;
412  sm = sedge1;
413 
414  if ( (sp - s(i,j,domlo+1,n))*(s(i,j,domlo+1,n) - sm) <= 0. )
415  {
416  sp = s(i,j,domlo+1,n);
417  sm = s(i,j,domlo+1,n);
418  }
419  else if(amrex::Math::abs(sp - s(i,j,domlo+1,n)) >= 2.*amrex::Math::abs(sm - s(i,j,domlo+1,n)))
420  sp = amrex::Real(3.0)*s(i,j,domlo+1,n) - amrex::Real(amrex::Real(2.0))*sm;
421 
422  else if(amrex::Math::abs(sm - s(i,j,domlo+1,n)) >= 2.*amrex::Math::abs(sp - s(i,j,domlo+1,n)))
423  sm = amrex::Real(3.0)*s(i,j,domlo+1,n) - amrex::Real(amrex::Real(2.0))*sp;
424  }
425  }
426 
427  if (bchi == BCType::ext_dir || bchi == BCType::hoextrap)
428  {
429  if (k == domhi)
430  {
431  sedge1 = -amrex::Real(0.2)*s(i,j,domhi+1,n) + amrex::Real(0.75)*s(i,j,domhi ,n)
432  +amrex::Real(amrex::Real(0.5))*s(i,j,domhi-1,n) - amrex::Real(0.05)*s(i,j,domhi-2,n);
433  sedge1 = amrex::max(sedge1, amrex::min(s(i,j,domhi-1,n), s(i,j,domhi,n)));
434  sedge1 = amrex::min(sedge1, amrex::max(s(i,j,domhi-1,n), s(i,j,domhi,n)));
435 
436  sp = s(i,j,domhi+1,n);
437  sm = sedge1;
438 
439  } else if (k == domhi-1) {
440 
441  sedge2 = -amrex::Real(0.2)*s(i,j,domhi+1,n) + amrex::Real(0.75)*s(i,j,domhi ,n)
442  +amrex::Real(amrex::Real(0.5))*s(i,j,domhi-1,n) - amrex::Real(0.05)*s(i,j,domhi-2,n);
443  sedge2 = amrex::max(sedge2, amrex::min(s(i,j,domhi-1,n), s(i,j,domhi,n)));
444  sedge2 = amrex::min(sedge2, amrex::max(s(i,j,domhi-1,n), s(i,j,domhi,n)));
445 
446  sp = sedge2;
447  sm = sedge1;
448 
449  if ( (sp - s(i,j,domhi-1,n))*(s(i,j,domhi-1,n) - sm) <= 0. )
450  {
451  sp = s(i,j,domhi-1,n);
452  sm = s(i,j,domhi-1,n);
453  }
454  else if(amrex::Math::abs(sp - s(i,j,domhi-1,n)) >= 2.*amrex::Math::abs(sm - s(i,j,domhi-1,n)))
455  sp = amrex::Real(3.0)*s(i,j,domhi-1,n) - amrex::Real(amrex::Real(2.0))*sm;
456 
457  else if(amrex::Math::abs(sm - s(i,j,domhi-1,n)) >= 2.*amrex::Math::abs(sp - s(i,j,domhi-1,n)))
458  sm = amrex::Real(3.0)*s(i,j,domhi-1,n) - amrex::Real(amrex::Real(2.0))*sp;
459 
460  }
461  }
462 }
463 #endif
464 
465 
466 // Right now only ppm type 1 is supported on GPU
467 // This version is called before the MAC projection, when we use the cell-centered velocity
468 // for upwinding
469 template <typename Limiter>
471 void PredictVelOnXFace ( const int i, const int j, const int k, const int n,
472  const amrex::Real dtdx, const amrex::Real v_ad,
474  const amrex::Array4<amrex::Real> &Im,
475  const amrex::Array4<amrex::Real> &Ip,
476  const amrex::BCRec bc, const int domlo, const int domhi,
477  const Limiter& limiter)
478 {
479  constexpr amrex::Real m_half{amrex::Real(0.5)};
480  constexpr amrex::Real m_two3rds{amrex::Real(2.0/3.0)};
481 
482  amrex::Real sm2 = S(i-2,j,k,n);
483  amrex::Real sm1 = S(i-1,j,k,n);
484  amrex::Real s0 = S(i ,j,k,n);
485  amrex::Real sp1 = S(i+1,j,k,n);
486  amrex::Real sp2 = S(i+2,j,k,n);
487 
488  amrex::Real sedge1 = limiter.sedge1(sm2,sm1,s0,sp1,sp2);
489  amrex::Real sedge2 = limiter.sedge2(sm2,sm1,s0,sp1,sp2);
490 
491  amrex::Real sm, sp;
492  amrex::Tie(sm,sp) = limiter.sm_sp(s0,sedge1,sedge2);
493 
494  SetXBCs(i, j, k, n, sm, sp, sedge1, sedge2, S, bc.lo(0), bc.hi(0), domlo, domhi);
495 
496  amrex::Real s6 = amrex::Real(6.0)*s0 - amrex::Real(3.0)*(sm + sp);
497 
498  amrex::Real sigma = amrex::Math::abs(v_ad)*dtdx;
499 
500  if (v_ad > small_vel)
501  {
502  Ip(i,j,k,n) = sp - (m_half*sigma)*((sp-sm) - (amrex::Real(amrex::Real(1.0)) - m_two3rds*sigma)*s6);
503  Im(i,j,k,n) = S(i,j,k,n);
504  }
505  else if (v_ad < -small_vel)
506  {
507  Ip(i,j,k,n) = S(i,j,k,n);
508  Im(i,j,k,n) = sm + (m_half*sigma)*((sp-sm) + (amrex::Real(amrex::Real(1.0)) - m_two3rds*sigma)*s6);
509  } else
510  {
511  Ip(i,j,k,n) = S(i,j,k,n);
512  Im(i,j,k,n) = S(i,j,k,n);
513  }
514 }
515 
516 template <typename Limiter>
518 void PredictVelOnYFace ( const int i, const int j, const int k, const int n,
519  const amrex::Real dtdy, const amrex::Real v_ad,
521  const amrex::Array4<amrex::Real> &Im,
522  const amrex::Array4<amrex::Real> &Ip,
523  const amrex::BCRec bc, const int domlo, const int domhi,
524  const Limiter& limiter)
525 {
526  constexpr amrex::Real m_half{amrex::Real(0.5)};
527  constexpr amrex::Real m_two3rds{amrex::Real(2.0/3.0)};
528 
529  amrex::Real sm2 = S(i,j-2,k,n);
530  amrex::Real sm1 = S(i,j-1,k,n);
531  amrex::Real s0 = S(i,j ,k,n);
532  amrex::Real sp1 = S(i,j+1,k,n);
533  amrex::Real sp2 = S(i,j+2,k,n);
534 
535  amrex::Real sedge1 = limiter.sedge1(sm2,sm1,s0,sp1,sp2);
536  amrex::Real sedge2 = limiter.sedge2(sm2,sm1,s0,sp1,sp2);
537 
538  amrex::Real sm, sp;
539  amrex::Tie(sm,sp) = limiter.sm_sp(s0,sedge1,sedge2);
540 
541  SetYBCs(i, j, k, n, sm, sp, sedge1, sedge2, S, bc.lo(1), bc.hi(1), domlo, domhi);
542 
543  amrex::Real s6 = amrex::Real(6.0)*s0 - amrex::Real(3.0)*(sm + sp);
544 
545  amrex::Real sigma = amrex::Math::abs(v_ad)*dtdy;
546 
547  if (v_ad > small_vel)
548  {
549  Ip(i,j,k,n) = sp - (m_half*sigma)*((sp-sm) - (amrex::Real(1.0) - m_two3rds*sigma)*s6);
550  Im(i,j,k,n) = S(i,j,k,n);
551  }
552  else if (v_ad < -small_vel)
553  {
554  Ip(i,j,k,n) = S(i,j,k,n);
555  Im(i,j,k,n) = sm + (m_half*sigma)*((sp-sm) + (amrex::Real(1.0) - m_two3rds*sigma)*s6);
556  } else
557  {
558  Ip(i,j,k,n) = S(i,j,k,n);
559  Im(i,j,k,n) = S(i,j,k,n);
560  }
561 }
562 
563 #if (AMREX_SPACEDIM==3)
564 template <typename Limiter>
566 void PredictVelOnZFace ( const int i, const int j, const int k, const int n,
567  const amrex::Real dtdz, const amrex::Real v_ad,
569  const amrex::Array4<amrex::Real> &Im,
570  const amrex::Array4<amrex::Real> &Ip,
571  const amrex::BCRec bc, const int domlo, const int domhi,
572  const Limiter& limiter)
573 {
574  constexpr amrex::Real m_half{amrex::Real(0.5)};
575  constexpr amrex::Real m_two3rds{amrex::Real(2.0/3.0)};
576 
577  amrex::Real sm2 = S(i,j,k-2,n);
578  amrex::Real sm1 = S(i,j,k-1,n);
579  amrex::Real s0 = S(i,j,k ,n);
580  amrex::Real sp1 = S(i,j,k+1,n);
581  amrex::Real sp2 = S(i,j,k+2,n);
582 
583  amrex::Real sedge1 = limiter.sedge1(sm2,sm1,s0,sp1,sp2);
584  amrex::Real sedge2 = limiter.sedge2(sm2,sm1,s0,sp1,sp2);
585 
586  amrex::Real sm, sp;
587  amrex::Tie(sm,sp) = limiter.sm_sp(s0,sedge1,sedge2);
588 
589  SetZBCs(i, j, k, n, sm, sp, sedge1, sedge2, S, bc.lo(2), bc.hi(2), domlo, domhi);
590 
591  amrex::Real s6 = amrex::Real(6.0)*s0 - amrex::Real(3.0)*(sm + sp);
592 
593  amrex::Real sigma = amrex::Math::abs(v_ad)*dtdz;
594 
595  if (v_ad > small_vel)
596  {
597  Ip(i,j,k,n) = sp - (m_half*sigma)*((sp-sm) - (amrex::Real(1.0) - m_two3rds*sigma)*s6);
598  Im(i,j,k,n) = S(i,j,k,n);
599  }
600  else if (v_ad < -small_vel)
601  {
602  Ip(i,j,k,n) = S(i,j,k,n);
603  Im(i,j,k,n) = sm + (m_half*sigma)*((sp-sm) + (amrex::Real(1.0) - m_two3rds*sigma)*s6);
604  } else
605  {
606  Ip(i,j,k,n) = S(i,j,k,n);
607  Im(i,j,k,n) = S(i,j,k,n);
608  }
609 }
610 
611 #endif
612 
613 // Right now only ppm type 1 is supported on GPU
614 // This version is called after the MAC projection, when we use the MAC-projected velocity
615 // for upwinding
616 template <typename Limiter>
618 void PredictStateOnXFace ( const int i, const int j, const int k, const int n,
619  const amrex::Real dt, const amrex::Real dx,
620  amrex::Real& Im, amrex::Real& Ip,
622  const amrex::Array4<const amrex::Real> &vel_edge,
623  const amrex::BCRec bc,
624  const int domlo, const int domhi,
625  const Limiter& limiter)
626 {
627  {
628  using namespace amrex;
629 
630  constexpr amrex::Real m_half{amrex::Real(0.5)};
631  constexpr amrex::Real m_two3rds{amrex::Real(2.0/3.0)};
632 
633  amrex::Real sm2 = S(i-2,j,k,n);
634  amrex::Real sm1 = S(i-1,j,k,n);
635  amrex::Real s0 = S(i ,j,k,n);
636  amrex::Real sp1 = S(i+1,j,k,n);
637  amrex::Real sp2 = S(i+2,j,k,n);
638 
639  amrex::Real sedge1 = limiter.sedge1(sm2,sm1,s0,sp1,sp2);
640  amrex::Real sedge2 = limiter.sedge2(sm2,sm1,s0,sp1,sp2);
641 
642  amrex::Real sm, sp;
643  amrex::Tie(sm,sp) = limiter.sm_sp(s0,sedge1,sedge2);
644 
645  SetXBCs(i, j, k, n, sm, sp, sedge1, sedge2, S, bc.lo(0), bc.hi(0), domlo, domhi);
646 
647  amrex::Real s6 = 6.0*s0 - 3.0*(sm + sp);
648 
649  amrex::Real sigmap = amrex::Math::abs(vel_edge(i+1,j,k))*dt/dx;
650  amrex::Real sigmam = amrex::Math::abs(vel_edge(i ,j,k))*dt/dx;
651 
652  if (vel_edge(i+1,j,k) > small_vel)
653  Ip = sp - (m_half*sigmap)*((sp - sm) - (amrex::Real(amrex::Real(1.0)) -m_two3rds*sigmap)*s6);
654  else
655  Ip = S(i,j,k,n);
656 
657  if(vel_edge(i,j,k) < -small_vel)
658  Im = sm + (m_half*sigmam)*((sp-sm) + (amrex::Real(amrex::Real(1.0)) - m_two3rds*sigmam)*s6);
659  else
660  Im = S(i,j,k,n);
661  }
662 }
663 template <typename Limiter>
664 
666 void PredictStateOnYFace ( const int i, const int j, const int k, const int n,
667  const amrex::Real dt, const amrex::Real dx,
668  amrex::Real& Im, amrex::Real& Ip,
670  const amrex::Array4<const amrex::Real> &vel_edge,
671  const amrex::BCRec bc,
672  const int domlo, const int domhi,
673  const Limiter& limiter)
674 {
675  {
676  using namespace amrex;
677 
678  constexpr amrex::Real m_half{amrex::Real(0.5)};
679  constexpr amrex::Real m_two3rds{amrex::Real(2.0/3.0)};
680 
681  amrex::Real sm2 = S(i,j-2,k,n);
682  amrex::Real sm1 = S(i,j-1,k,n);
683  amrex::Real s0 = S(i,j ,k,n);
684  amrex::Real sp1 = S(i,j+1,k,n);
685  amrex::Real sp2 = S(i,j+2,k,n);
686 
687  amrex::Real sedge1 = limiter.sedge1(sm2,sm1,s0,sp1,sp2);
688  amrex::Real sedge2 = limiter.sedge2(sm2,sm1,s0,sp1,sp2);
689 
690  amrex::Real sm, sp;
691  amrex::Tie(sm,sp) = limiter.sm_sp(s0,sedge1,sedge2);
692 
693  SetYBCs(i, j, k, n, sm, sp, sedge1, sedge2, S, bc.lo(1), bc.hi(1), domlo, domhi);
694 
695  amrex::Real s6 = 6.0*s0- 3.0*(sm + sp);
696 
697  amrex::Real sigmap = amrex::Math::abs(vel_edge(i,j+1,k))*dt/dx;
698  amrex::Real sigmam = amrex::Math::abs(vel_edge(i,j ,k))*dt/dx;
699 
700  if (vel_edge(i,j+1,k) > small_vel)
701  Ip = sp - (m_half*sigmap)*((sp - sm) - (amrex::Real(amrex::Real(1.0)) -m_two3rds*sigmap)*s6);
702  else
703  Ip = S(i,j,k,n);
704 
705  if (vel_edge(i,j,k) < -small_vel)
706  Im = sm + (m_half*sigmam)*((sp-sm) + (amrex::Real(amrex::Real(1.0)) - m_two3rds*sigmam)*s6);
707  else
708  Im = S(i,j,k,n);
709  }
710 }
711 
712 
713 #if (AMREX_SPACEDIM==3)
714 template <typename Limiter>
716 void PredictStateOnZFace ( const int i, const int j, const int k, const int n,
717  const amrex::Real dt, const amrex::Real dx,
718  amrex::Real& Im, amrex::Real& Ip,
720  const amrex::Array4<const amrex::Real> &vel_edge,
721  const amrex::BCRec bc,
722  const int domlo, const int domhi,
723  const Limiter& limiter)
724 {
725  {
726  using namespace amrex;
727 
728  constexpr amrex::Real m_half{amrex::Real(0.5)};
729  constexpr amrex::Real m_two3rds{amrex::Real(2.0/3.0)};
730 
731  amrex::Real sm2 = S(i,j,k-2,n);
732  amrex::Real sm1 = S(i,j,k-1,n);
733  amrex::Real s0 = S(i,j,k ,n);
734  amrex::Real sp1 = S(i,j,k+1,n);
735  amrex::Real sp2 = S(i,j,k+2,n);
736 
737  amrex::Real sedge1 = limiter.sedge1(sm2,sm1,s0,sp1,sp2);
738  amrex::Real sedge2 = limiter.sedge2(sm2,sm1,s0,sp1,sp2);
739 
740  amrex::Real sm, sp;
741  amrex::Tie(sm,sp) = limiter.sm_sp(s0,sedge1,sedge2);
742 
743  SetZBCs(i, j, k, n, sm, sp, sedge1, sedge2, S, bc.lo(2), bc.hi(2), domlo, domhi);
744 
745  amrex::Real s6 = 6.0*s0- 3.0*(sm + sp);
746  amrex::Real sigmap = amrex::Math::abs(vel_edge(i,j,k+1))*dt/dx;
747  amrex::Real sigmam = amrex::Math::abs(vel_edge(i,j,k ))*dt/dx;
748 
749  if(vel_edge(i,j,k+1) > small_vel)
750  Ip = sp - (m_half*sigmap)*((sp-sm) - (amrex::Real(amrex::Real(1.0)) -m_two3rds*sigmap)*s6);
751  else
752  Ip = S(i,j,k,n);
753 
754  if(vel_edge(i,j,k) < -small_vel)
755  Im = sm + (m_half*sigmam)*((sp-sm) + (amrex::Real(amrex::Real(1.0)) - m_two3rds*sigmam)*s6);
756  else
757  Im = S(i,j,k,n);
758  }
759 }
760 #endif
761 
762 template <typename Limiter>
765  amrex::Array4<amrex::Real> const& Imy,
766  amrex::Array4<amrex::Real> const& Imz),
768  amrex::Array4<amrex::Real> const& Ipy,
769  amrex::Array4<amrex::Real> const& Ipz),
772  amrex::Geometry geom,
773  amrex::Real dt,
774  amrex::BCRec const* pbc,
775  const Limiter& limiter)
776 {
777  const amrex::Box& domain = geom.Domain();
778  const amrex::Dim3 dlo = amrex::lbound(domain);
779  const amrex::Dim3 dhi = amrex::ubound(domain);
780 
781  const auto dx = geom.CellSizeArray();
782  AMREX_D_TERM( amrex::Real l_dtdx = dt / dx[0];,
783  amrex::Real l_dtdy = dt / dx[1];,
784  amrex::Real l_dtdz = dt / dx[2];);
785 
786  amrex::ParallelFor(bx, AMREX_SPACEDIM,
787  [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
788  {
789  PredictVelOnXFace(i,j,k,n,l_dtdx,vel(i,j,k,0),q,Imx,Ipx,pbc[n],dlo.x,dhi.x,limiter);
790  PredictVelOnYFace(i,j,k,n,l_dtdy,vel(i,j,k,1),q,Imy,Ipy,pbc[n],dlo.y,dhi.y,limiter);
791 #if (AMREX_SPACEDIM==3)
792  PredictVelOnZFace(i,j,k,n,l_dtdz,vel(i,j,k,2),q,Imz,Ipz,pbc[n],dlo.z,dhi.z,limiter);
793 #endif
794  });
795 }
796 
797 template <typename Limiter>
800  amrex::Array4<amrex::Real> const& Imy,
801  amrex::Array4<amrex::Real> const& Imz),
803  amrex::Array4<amrex::Real> const& Ipy,
804  amrex::Array4<amrex::Real> const& Ipz),
809  amrex::Geometry geom,
810  amrex::Real l_dt,
811  amrex::BCRec const* pbc,
812  const int ncomp,
813  const Limiter& limiter)
814 {
815  const amrex::Box& domain = geom.Domain();
816  const amrex::Dim3 dlo = amrex::lbound(domain);
817  const amrex::Dim3 dhi = amrex::ubound(domain);
818 
819  AMREX_D_TERM (const auto dx = geom.CellSize(0);,
820  const auto dy = geom.CellSize(1);,
821  const auto dz = geom.CellSize(2););
822 
823  amrex::ParallelFor(bx, ncomp,
824  [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
825  {
826  PPM::PredictStateOnXFace(i, j, k, n, l_dt, dx, Imx(i,j,k,n), Ipx(i,j,k,n),
827  q, umac, pbc[n], dlo.x, dhi.x,limiter);
828  PPM::PredictStateOnYFace(i, j, k, n, l_dt, dy, Imy(i,j,k,n), Ipy(i,j,k,n),
829  q, vmac, pbc[n], dlo.y, dhi.y,limiter);
830 #if (AMREX_SPACEDIM==3)
831  PPM::PredictStateOnZFace(i, j, k, n, l_dt, dz, Imz(i,j,k,n), Ipz(i,j,k,n),
832  q, wmac, pbc[n], dlo.z, dhi.z,limiter);
833 #endif
834  });
835 }
836 
837 
838 }
839 
840 #endif
841 /** @} */
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
#define AMREX_GPU_HOST_DEVICE
#define AMREX_D_TERM(a, b, c)
#define AMREX_D_DECL(a, b, c)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE const int * hi() const &noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE const int * lo() const &noexcept
GpuArray< Real, AMREX_SPACEDIM > CellSizeArray() const noexcept
const Real * CellSize() const noexcept
const Box & Domain() const noexcept
static constexpr amrex::Real small_vel
Definition: hydro_constants.H:37
Definition: hydro_godunov_ppm.H:15
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void PredictVelOnXFace(const int i, const int j, const int k, const int n, const amrex::Real dtdx, const amrex::Real v_ad, const amrex::Array4< const amrex::Real > &S, const amrex::Array4< amrex::Real > &Im, const amrex::Array4< amrex::Real > &Ip, const amrex::BCRec bc, const int domlo, const int domhi, const Limiter &limiter)
Definition: hydro_godunov_ppm.H:471
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void PredictStateOnXFace(const int i, const int j, const int k, const int n, const amrex::Real dt, const amrex::Real dx, amrex::Real &Im, amrex::Real &Ip, const amrex::Array4< const amrex::Real > &S, const amrex::Array4< const amrex::Real > &vel_edge, const amrex::BCRec bc, const int domlo, const int domhi, const Limiter &limiter)
Definition: hydro_godunov_ppm.H:618
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void SetYBCs(const int i, const int j, const int k, const int n, amrex::Real &sm, amrex::Real &sp, amrex::Real &sedge1, amrex::Real &sedge2, const amrex::Array4< const amrex::Real > &s, const int bclo, const int bchi, const int domlo, const int domhi)
Definition: hydro_godunov_ppm.H:302
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void PredictVelOnYFace(const int i, const int j, const int k, const int n, const amrex::Real dtdy, const amrex::Real v_ad, const amrex::Array4< const amrex::Real > &S, const amrex::Array4< amrex::Real > &Im, const amrex::Array4< amrex::Real > &Ip, const amrex::BCRec bc, const int domlo, const int domhi, const Limiter &limiter)
Definition: hydro_godunov_ppm.H:518
void PredictStateOnFaces(amrex::Box const &bx, AMREX_D_DECL(amrex::Array4< amrex::Real > const &Imx, amrex::Array4< amrex::Real > const &Imy, amrex::Array4< amrex::Real > const &Imz), AMREX_D_DECL(amrex::Array4< amrex::Real > const &Ipx, amrex::Array4< amrex::Real > const &Ipy, amrex::Array4< amrex::Real > const &Ipz), AMREX_D_DECL(amrex::Array4< amrex::Real const > const &umac, amrex::Array4< amrex::Real const > const &vmac, amrex::Array4< amrex::Real const > const &wmac), amrex::Array4< amrex::Real const > const &q, amrex::Geometry geom, amrex::Real l_dt, amrex::BCRec const *pbc, const int ncomp, const Limiter &limiter)
Definition: hydro_godunov_ppm.H:798
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void SetXBCs(const int i, const int j, const int k, const int n, amrex::Real &sm, amrex::Real &sp, amrex::Real &sedge1, amrex::Real &sedge2, const amrex::Array4< const amrex::Real > &s, const int bclo, const int bchi, const int domlo, const int domhi)
Definition: hydro_godunov_ppm.H:223
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void PredictStateOnYFace(const int i, const int j, const int k, const int n, const amrex::Real dt, const amrex::Real dx, amrex::Real &Im, amrex::Real &Ip, const amrex::Array4< const amrex::Real > &S, const amrex::Array4< const amrex::Real > &vel_edge, const amrex::BCRec bc, const int domlo, const int domhi, const Limiter &limiter)
Definition: hydro_godunov_ppm.H:666
void PredictVelOnFaces(amrex::Box const &bx, AMREX_D_DECL(amrex::Array4< amrex::Real > const &Imx, amrex::Array4< amrex::Real > const &Imy, amrex::Array4< amrex::Real > const &Imz), AMREX_D_DECL(amrex::Array4< amrex::Real > const &Ipx, amrex::Array4< amrex::Real > const &Ipy, amrex::Array4< amrex::Real > const &Ipz), amrex::Array4< amrex::Real const > const &q, amrex::Array4< amrex::Real const > const &vel, amrex::Geometry geom, amrex::Real dt, amrex::BCRec const *pbc, const Limiter &limiter)
Definition: hydro_godunov_ppm.H:763
limiters
Definition: hydro_godunov_ppm.H:17
@ VanLeer
Definition: hydro_godunov_ppm.H:17
@ WENOZ
Definition: hydro_godunov_ppm.H:17
@ NoLimiter
Definition: hydro_godunov_ppm.H:17
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
constexpr AMREX_GPU_HOST_DEVICE GpuTuple< Args &... > Tie(Args &... args) noexcept
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & min(const T &a, const T &b) noexcept
constexpr AMREX_GPU_HOST_DEVICE GpuTuple< detail::tuple_decay_t< Ts >... > makeTuple(Ts &&... args)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 ubound(Array4< T > const &a) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 lbound(Array4< T > const &a) noexcept
std::enable_if_t< std::is_integral< T >::value > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition: hydro_godunov_ppm.H:19
const amrex::Real m_sixth
Definition: hydro_godunov_ppm.H:68
const amrex::Real m_half
Definition: hydro_godunov_ppm.H:67
static AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::GpuTuple< amrex::Real, amrex::Real > sm_sp(const amrex::Real, const amrex::Real sedge1, const amrex::Real sedge2)
Definition: hydro_godunov_ppm.H:60
nolimiter()=default
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sedge2(const amrex::Real, const amrex::Real sm1, const amrex::Real s0, const amrex::Real sp1, const amrex::Real sp2) const
Definition: hydro_godunov_ppm.H:43
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sedge1(const amrex::Real sm2, const amrex::Real sm1, const amrex::Real s0, const amrex::Real sp1, const amrex::Real) const
Definition: hydro_godunov_ppm.H:26
Definition: hydro_godunov_ppm.H:71
const amrex::Real m_sixth
Definition: hydro_godunov_ppm.H:159
vanleer()=default
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sedge1(const amrex::Real sm2, const amrex::Real sm1, const amrex::Real s0, const amrex::Real sp1, const amrex::Real) const
Definition: hydro_godunov_ppm.H:89
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real vanLeer(const amrex::Real a, const amrex::Real b, const amrex::Real c) const
Definition: hydro_godunov_ppm.H:76
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sedge2(const amrex::Real, const amrex::Real sm1, const amrex::Real s0, const amrex::Real sp1, const amrex::Real sp2) const
Definition: hydro_godunov_ppm.H:106
static AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::GpuTuple< amrex::Real, amrex::Real > sm_sp(const amrex::Real s0, const amrex::Real sedge1, const amrex::Real sedge2)
Definition: hydro_godunov_ppm.H:124
const amrex::Real m_half
Definition: hydro_godunov_ppm.H:158
Definition: hydro_godunov_ppm.H:162
amrex::Real m_eps
Definition: hydro_godunov_ppm.H:218
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sedge2(const amrex::Real sm2, const amrex::Real sm1, const amrex::Real s0, const amrex::Real sp1, const amrex::Real sp2) const
Definition: hydro_godunov_ppm.H:168
static AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::GpuTuple< amrex::Real, amrex::Real > sm_sp(const amrex::Real, const amrex::Real sedge1, const amrex::Real sedge2)
Definition: hydro_godunov_ppm.H:211
wenoz()=default
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sedge1(const amrex::Real sm2, const amrex::Real sm1, const amrex::Real s0, const amrex::Real sp1, const amrex::Real sp2) const
Definition: hydro_godunov_ppm.H:200