AMReX-Hydro
AMReX-based hydro routines for low Mach number flows
hydro_ebmol_edge_state_K.H
Go to the documentation of this file.
1 /**
2  * \file hydro_ebmol_edge_state_K.H
3  * \addtogroup EBMOL
4  * @{
5  *
6  */
7 
8 using namespace amrex::literals;
9 
10 #ifndef HYDRO_EBMOL_EDGE_STATE_K_H_
11 #define HYDRO_EBMOL_EDGE_STATE_K_H_
12 
13 #include <AMReX_EB_Slopes_K.H>
14 #include <AMReX_BCRec.H>
15 #include <hydro_constants.H>
16 #include <hydro_bcs_K.H>
17 
18 namespace EBMOL {
19 
21 amrex::Real hydro_ebmol_xedge_state_extdir ( AMREX_D_DECL(int i, int j, int k), int n,
30  amrex::BCRec const* const d_bcrec,
31  amrex::Box const& domain,
32  int order, const bool is_velocity) noexcept
33 {
34 
35 #if (AMREX_SPACEDIM==2)
36  const int k = 0;
37 #endif
38 
39  amrex::Real qs;
40  int domlo = domain.smallEnd(0);
41  int domhi = domain.bigEnd(0);
42  bool extlo = d_bcrec[n].lo(0) == amrex::BCType::ext_dir;
43  bool exthi = d_bcrec[n].hi(0) == amrex::BCType::ext_dir;
44 
45  if ( extlo && i <= domlo)
46  {
47  qs = q(domlo-1,j,k,n);
48  }
49  else if ( exthi && i >= domhi+1)
50  {
51  qs = q(domhi+1,j,k,n);
52  }
53  else
54  {
55  const int domain_ilo = domain.smallEnd(0);
56  const int domain_ihi = domain.bigEnd(0);
57  const int domain_jlo = domain.smallEnd(1);
58  const int domain_jhi = domain.bigEnd(1);
59 #if (AMREX_SPACEDIM == 3)
60  const int domain_klo = domain.smallEnd(2);
61  const int domain_khi = domain.bigEnd(2);
62 #endif
63 
64  amrex::Real yf = fcx(i,j,k,0); // local (y,z) of centroid of z-face we are extrapolating to
65 #if (AMREX_SPACEDIM==3)
66  amrex::Real zf = fcx(i,j,k,1);
67 #endif
68 
69  AMREX_D_TERM(amrex::Real xc = ccc(i,j,k,0);, // centroid of cell (i,j,k)
70  amrex::Real yc = ccc(i,j,k,1);,
71  amrex::Real zc = ccc(i,j,k,2););
72 
73  AMREX_D_TERM(amrex::Real delta_x = amrex::Real(0.5) + xc;,
74  amrex::Real delta_y = yf - yc;,
75  amrex::Real delta_z = zf - zc;);
76 
77  amrex::Real cc_qmax = amrex::max(q(i,j,k,n),q(i-1,j,k,n));
78  amrex::Real cc_qmin = amrex::min(q(i,j,k,n),q(i-1,j,k,n));
79 
80  AMREX_D_TERM(bool extdir_or_ho_ilo = (d_bcrec[n].lo(0) == amrex::BCType::ext_dir) ||
81  (d_bcrec[n].lo(0) == amrex::BCType::hoextrap);,
82  bool extdir_or_ho_jlo = (d_bcrec[n].lo(1) == amrex::BCType::ext_dir) ||
83  (d_bcrec[n].lo(1) == amrex::BCType::hoextrap);,
84  bool extdir_or_ho_klo = (d_bcrec[n].lo(2) == amrex::BCType::ext_dir) ||
85  (d_bcrec[n].lo(2) == amrex::BCType::hoextrap););
86 
87  AMREX_D_TERM(bool extdir_or_ho_ihi = (d_bcrec[n].hi(0) == amrex::BCType::ext_dir) ||
88  (d_bcrec[n].hi(0) == amrex::BCType::hoextrap);,
89  bool extdir_or_ho_jhi = (d_bcrec[n].hi(1) == amrex::BCType::ext_dir) ||
90  (d_bcrec[n].hi(1) == amrex::BCType::hoextrap);,
91  bool extdir_or_ho_khi = (d_bcrec[n].hi(2) == amrex::BCType::ext_dir) ||
92  (d_bcrec[n].hi(2) == amrex::BCType::hoextrap););
93 
94  const auto& slopes_eb_hi =
95  amrex_lim_slopes_extdir_eb(i, j, k, n, q, ccc, vfrac,
96  AMREX_D_DECL(fcx,fcy,fcz), flag,
97  AMREX_D_DECL(extdir_or_ho_ilo, extdir_or_ho_jlo, extdir_or_ho_klo),
98  AMREX_D_DECL(extdir_or_ho_ihi, extdir_or_ho_jhi, extdir_or_ho_khi),
99  AMREX_D_DECL(domain_ilo, domain_jlo, domain_klo),
100  AMREX_D_DECL(domain_ihi, domain_jhi, domain_khi),
101  order);
102 
103 
104 #if (AMREX_SPACEDIM==3)
105  amrex::Real qpls = q(i ,j,k,n) - delta_x * slopes_eb_hi[0]
106  + delta_y * slopes_eb_hi[1]
107  + delta_z * slopes_eb_hi[2];
108 #else
109  amrex::Real qpls = q(i ,j,k,n) - delta_x * slopes_eb_hi[0]
110  + delta_y * slopes_eb_hi[1];
111 #endif
112 
113  qpls = amrex::max(amrex::min(qpls, cc_qmax), cc_qmin);
114 
115  AMREX_D_TERM(xc = ccc(i-1,j,k,0);, // centroid of cell (i-1,j,k)
116  yc = ccc(i-1,j,k,1);,
117  zc = ccc(i-1,j,k,2););
118 
119  AMREX_D_TERM( delta_x = amrex::Real(0.5) - xc;,
120  delta_y = yf - yc;,
121  delta_z = zf - zc;);
122 
123  // Compute slopes of component "n" of q
124  const auto& slopes_eb_lo =
125  amrex_lim_slopes_extdir_eb(i-1, j, k, n, q, ccc, vfrac,
126  AMREX_D_DECL(fcx,fcy,fcz), flag,
127  AMREX_D_DECL(extdir_or_ho_ilo, extdir_or_ho_jlo, extdir_or_ho_klo),
128  AMREX_D_DECL(extdir_or_ho_ihi, extdir_or_ho_jhi, extdir_or_ho_khi),
129  AMREX_D_DECL(domain_ilo, domain_jlo, domain_klo),
130  AMREX_D_DECL(domain_ihi, domain_jhi, domain_khi),
131  order);
132 
133 #if (AMREX_SPACEDIM==3)
134  amrex::Real qmns = q(i-1,j,k,n) + delta_x * slopes_eb_lo[0]
135  + delta_y * slopes_eb_lo[1]
136  + delta_z * slopes_eb_lo[2];
137 #else
138  amrex::Real qmns = q(i-1,j,k,n) + delta_x * slopes_eb_lo[0]
139  + delta_y * slopes_eb_lo[1];
140 #endif
141 
142  qmns = amrex::max(amrex::min(qmns, cc_qmax), cc_qmin);
143 
144  HydroBC::SetXEdgeBCs(i, j, k, 0, q, qmns, qpls, d_bcrec[n].lo(0), domain_ilo, d_bcrec[n].hi(0), domain_ihi, is_velocity);
145 
146  if ( (i==domain_ilo) && (d_bcrec[n].lo(0) == amrex::BCType::foextrap || d_bcrec[n].lo(0) == amrex::BCType::hoextrap) )
147  {
148  if ( umac(i,j,k) >= 0. && n==XVEL && is_velocity ) qpls = amrex::min(qpls,amrex::Real(0.0));
149  qmns = qpls;
150  }
151  if ( (i==domain_ihi+1) && (d_bcrec[n].hi(0) == amrex::BCType::foextrap || d_bcrec[n].hi(0) == amrex::BCType::hoextrap) )
152  {
153  if ( umac(i,j,k) <= 0. && n==XVEL && is_velocity ) qmns = amrex::max(qmns,amrex::Real(0.0));
154  qpls = qmns;
155  }
156 
157  if (umac(i,j,k) > small_vel)
158  {
159  qs = qmns;
160  }
161  else if (umac(i,j,k) < - small_vel)
162  {
163  qs = qpls;
164  }
165  else
166  {
167  qs = amrex::Real(0.5)*(qmns+qpls);
168  }
169  }
170 
171  return qs;
172 }
173 
175 amrex::Real hydro_ebmol_xedge_state ( AMREX_D_DECL(int i, int j, int k), int n,
184  amrex::BCRec const* const d_bcrec,
185  amrex::Box const& domain,
186  int order, const bool is_velocity) noexcept
187 {
188 #if (AMREX_SPACEDIM==2)
189  const int k = 0;
190 #endif
191 
192  const int domain_ilo = domain.smallEnd(0);
193  const int domain_ihi = domain.bigEnd(0);
194 
195  // local (y,z) of centroid of x-face we are extrapolating to
196  amrex::Real yf = fcx(i,j,k,0);
197 #if (AMREX_SPACEDIM==3)
198  amrex::Real zf = fcx(i,j,k,1);
199 #endif
200 
201  AMREX_D_TERM(amrex::Real xc = ccc(i,j,k,0);, // centroid of cell (i,j,k)
202  amrex::Real yc = ccc(i,j,k,1);,
203  amrex::Real zc = ccc(i,j,k,2););
204 
205  AMREX_D_TERM(amrex::Real delta_x = amrex::Real(0.5) + xc;,
206  amrex::Real delta_y = yf - yc;,
207  amrex::Real delta_z = zf - zc;);
208 
209  amrex::Real cc_qmax = amrex::max(q(i,j,k,n),q(i-1,j,k,n));
210  amrex::Real cc_qmin = amrex::min(q(i,j,k,n),q(i-1,j,k,n));
211 
212  // Compute slopes of component "n" of q
213  const auto& slopes_eb_hi = amrex_lim_slopes_eb(i, j, k, n, q, ccc, vfrac,
214  AMREX_D_DECL(fcx,fcy,fcz),
215  flag, order);
216 
217 #if (AMREX_SPACEDIM==3)
218  amrex::Real qpls = q(i ,j,k,n) - delta_x * slopes_eb_hi[0]
219  + delta_y * slopes_eb_hi[1]
220  + delta_z * slopes_eb_hi[2];
221 #else
222  amrex::Real qpls = q(i ,j,k,n) - delta_x * slopes_eb_hi[0]
223  + delta_y * slopes_eb_hi[1];
224 #endif
225 
226  qpls = amrex::max(amrex::min(qpls, cc_qmax), cc_qmin);
227 
228  AMREX_D_TERM(xc = ccc(i-1,j,k,0);, // centroid of cell (i-1,j,k)
229  yc = ccc(i-1,j,k,1);,
230  zc = ccc(i-1,j,k,2););
231 
232  AMREX_D_TERM(delta_x = amrex::Real(0.5) - xc;,
233  delta_y = yf - yc;,
234  delta_z = zf - zc;);
235 
236  // Compute slopes of component "n" of q
237  const auto& slopes_eb_lo = amrex_lim_slopes_eb(i-1, j, k, n, q, ccc, vfrac,
238  AMREX_D_DECL(fcx,fcy,fcz),
239  flag, order);
240 
241 #if (AMREX_SPACEDIM==3)
242  amrex::Real qmns = q(i-1,j,k,n) + delta_x * slopes_eb_lo[0]
243  + delta_y * slopes_eb_lo[1]
244  + delta_z * slopes_eb_lo[2];
245 #else
246  amrex::Real qmns = q(i-1,j,k,n) + delta_x * slopes_eb_lo[0]
247  + delta_y * slopes_eb_lo[1];
248 #endif
249 
250  qmns = amrex::max(amrex::min(qmns, cc_qmax), cc_qmin);
251 
252  HydroBC::SetXEdgeBCs(i, j, k, 0, q, qmns, qpls, d_bcrec[n].lo(0), domain_ilo, d_bcrec[n].hi(0), domain_ihi, is_velocity);
253 
254  if ( (i==domain_ilo) && (d_bcrec[n].lo(0) == amrex::BCType::foextrap || d_bcrec[n].lo(0) == amrex::BCType::hoextrap) )
255  {
256  if ( umac(i,j,k) >= 0. && n==XVEL && is_velocity ) qpls = amrex::min(qpls,amrex::Real(0.0));
257  qmns = qpls;
258  }
259  if ( (i==domain_ihi+1) && (d_bcrec[n].hi(0) == amrex::BCType::foextrap || d_bcrec[n].hi(0) == amrex::BCType::hoextrap) )
260  {
261  if ( umac(i,j,k) <= 0. && n==XVEL && is_velocity ) qmns = amrex::max(qmns,amrex::Real(0.0));
262  qpls = qmns;
263  }
264 
265  amrex::Real qs;
266 
267  if (umac(i,j,k) > small_vel)
268  {
269  qs = qmns;
270  }
271  else if (umac(i,j,k) < - small_vel)
272  {
273  qs = qpls;
274  }
275  else
276  {
277  qs = amrex::Real(0.5)*(qmns+qpls);
278  }
279 
280  return qs;
281 }
282 
284 amrex::Real hydro_ebmol_yedge_state_extdir ( AMREX_D_DECL(int i, int j, int k), int n,
293  amrex::BCRec const* const d_bcrec,
294  amrex::Box const& domain,
295  int order, const bool is_velocity) noexcept
296 {
297 
298 #if (AMREX_SPACEDIM==2)
299  const int k = 0;
300 #endif
301 
302  amrex::Real qs;
303  int domlo = domain.smallEnd(1);
304  int domhi = domain.bigEnd(1);
305  bool extlo = d_bcrec[n].lo(1) == amrex::BCType::ext_dir;
306  bool exthi = d_bcrec[n].hi(1) == amrex::BCType::ext_dir;
307 
308  if ( extlo && j <= domlo)
309  {
310  qs = q(i,domlo-1,k,n);
311  }
312  else if (exthi && j >= domhi+1)
313  {
314  qs = q(i,domhi+1,k,n);
315  }
316  else
317  {
318  const int domain_ilo = domain.smallEnd(0);
319  const int domain_ihi = domain.bigEnd(0);
320  const int domain_jlo = domain.smallEnd(1);
321  const int domain_jhi = domain.bigEnd(1);
322 #if (AMREX_SPACEDIM == 3)
323  const int domain_klo = domain.smallEnd(2);
324  const int domain_khi = domain.bigEnd(2);
325 #endif
326 
327  AMREX_D_TERM(bool extdir_or_ho_ilo = (d_bcrec[n].lo(0) == amrex::BCType::ext_dir) ||
328  (d_bcrec[n].lo(0) == amrex::BCType::hoextrap);,
329  bool extdir_or_ho_jlo = (d_bcrec[n].lo(1) == amrex::BCType::ext_dir) ||
330  (d_bcrec[n].lo(1) == amrex::BCType::hoextrap);,
331  bool extdir_or_ho_klo = (d_bcrec[n].lo(2) == amrex::BCType::ext_dir) ||
332  (d_bcrec[n].lo(2) == amrex::BCType::hoextrap););
333 
334  AMREX_D_TERM(bool extdir_or_ho_ihi = (d_bcrec[n].hi(0) == amrex::BCType::ext_dir) ||
335  (d_bcrec[n].hi(0) == amrex::BCType::hoextrap);,
336  bool extdir_or_ho_jhi = (d_bcrec[n].hi(1) == amrex::BCType::ext_dir) ||
337  (d_bcrec[n].hi(1) == amrex::BCType::hoextrap);,
338  bool extdir_or_ho_khi = (d_bcrec[n].hi(2) == amrex::BCType::ext_dir) ||
339  (d_bcrec[n].hi(2) == amrex::BCType::hoextrap););
340 
341  amrex::Real xf = fcy(i,j,k,0); // local (x,z) of centroid of z-face we are extrapolating to
342 #if (AMREX_SPACEDIM==3)
343  amrex::Real zf = fcy(i,j,k,1);
344 #endif
345 
346  AMREX_D_TERM(amrex::Real xc = ccc(i,j,k,0);, // centroid of cell (i,j,k)
347  amrex::Real yc = ccc(i,j,k,1);,
348  amrex::Real zc = ccc(i,j,k,2););
349 
350  AMREX_D_TERM(amrex::Real delta_x = xf - xc;,
351  amrex::Real delta_y = amrex::Real(0.5) + yc;,
352  amrex::Real delta_z = zf - zc;);
353 
354  amrex::Real cc_qmax = amrex::max(q(i,j,k,n),q(i,j-1,k,n));
355  amrex::Real cc_qmin = amrex::min(q(i,j,k,n),q(i,j-1,k,n));
356 
357  // Compute slopes of component "n" of q
358  const auto& slopes_eb_hi =
359  amrex_lim_slopes_extdir_eb(i, j, k, n, q, ccc, vfrac,
360  AMREX_D_DECL(fcx,fcy,fcz), flag,
361  AMREX_D_DECL(extdir_or_ho_ilo, extdir_or_ho_jlo, extdir_or_ho_klo),
362  AMREX_D_DECL(extdir_or_ho_ihi, extdir_or_ho_jhi, extdir_or_ho_khi),
363  AMREX_D_DECL(domain_ilo, domain_jlo, domain_klo),
364  AMREX_D_DECL(domain_ihi, domain_jhi, domain_khi),
365  order);
366 
367 #if (AMREX_SPACEDIM==3)
368  amrex::Real qpls = q(i ,j,k,n) + delta_x * slopes_eb_hi[0]
369  - delta_y * slopes_eb_hi[1]
370  + delta_z * slopes_eb_hi[2];
371 #else
372  amrex::Real qpls = q(i ,j,k,n) + delta_x * slopes_eb_hi[0]
373  - delta_y * slopes_eb_hi[1];
374 #endif
375 
376  qpls = amrex::max(amrex::min(qpls, cc_qmax), cc_qmin);
377 
378  AMREX_D_TERM(xc = ccc(i,j-1,k,0);, // centroid of cell (i,j-1,k)
379  yc = ccc(i,j-1,k,1);,
380  zc = ccc(i,j-1,k,2););
381 
382  AMREX_D_TERM(delta_x = xf - xc;,
383  delta_y = amrex::Real(0.5) - yc;,
384  delta_z = zf - zc;);
385 
386  // Compute slopes of component "n" of q
387  const auto& slopes_eb_lo =
388  amrex_lim_slopes_extdir_eb(i, j-1, k, n, q, ccc, vfrac,
389  AMREX_D_DECL(fcx,fcy,fcz), flag,
390  AMREX_D_DECL(extdir_or_ho_ilo, extdir_or_ho_jlo, extdir_or_ho_klo),
391  AMREX_D_DECL(extdir_or_ho_ihi, extdir_or_ho_jhi, extdir_or_ho_khi),
392  AMREX_D_DECL(domain_ilo, domain_jlo, domain_klo),
393  AMREX_D_DECL(domain_ihi, domain_jhi, domain_khi),
394  order);
395 
396 
397 #if (AMREX_SPACEDIM==3)
398  amrex::Real qmns = q(i,j-1,k,n) + delta_x * slopes_eb_lo[0]
399  + delta_y * slopes_eb_lo[1]
400  + delta_z * slopes_eb_lo[2];
401 #else
402  amrex::Real qmns = q(i,j-1,k,n) + delta_x * slopes_eb_lo[0]
403  + delta_y * slopes_eb_lo[1];
404 #endif
405 
406  HydroBC::SetYEdgeBCs(i, j, k, n, q, qmns, qpls, d_bcrec[n].lo(1), domain_jlo, d_bcrec[n].hi(1), domain_jhi, is_velocity);
407 
408  if ( (j==domain_jlo) && (d_bcrec[n].lo(1) == amrex::BCType::foextrap || d_bcrec[n].lo(1) == amrex::BCType::hoextrap) )
409  {
410  if ( vmac(i,j,k) >= 0. && n==YVEL && is_velocity ) qpls = amrex::min(qpls,amrex::Real(0.0));
411  qmns = qpls;
412  }
413  if ( (j==domain_jhi+1) && (d_bcrec[n].hi(1) == amrex::BCType::foextrap || d_bcrec[n].hi(1) == amrex::BCType::hoextrap) )
414  {
415  if ( vmac(i,j,k) <= 0. && n==YVEL && is_velocity ) qmns = amrex::max(qmns,amrex::Real(0.0));
416  qpls = qmns;
417  }
418 
419  qmns = amrex::max(amrex::min(qmns, cc_qmax), cc_qmin);
420 
421  if (vmac(i,j,k) > small_vel)
422  {
423  qs = qmns;
424  }
425  else if (vmac(i,j,k) < - small_vel)
426  {
427  qs = qpls;
428  }
429  else
430  {
431  qs = amrex::Real(0.5)*(qmns+qpls);
432  }
433  }
434 
435  return qs;
436 }
437 
438 
439 
441 amrex::Real hydro_ebmol_yedge_state ( AMREX_D_DECL(int i, int j, int k), int n,
450  amrex::BCRec const* const d_bcrec,
451  amrex::Box const& domain,
452  int order, const bool is_velocity) noexcept
453 {
454 #if (AMREX_SPACEDIM==2)
455  const int k = 0;
456 #endif
457 
458  const int domain_jlo = domain.smallEnd(1);
459  const int domain_jhi = domain.bigEnd(1);
460 
461  // local (x,z) of centroid of z-face we are extrapolating to
462  amrex::Real xf = fcy(i,j,k,0);
463 #if (AMREX_SPACEDIM==3)
464  amrex::Real zf = fcy(i,j,k,1);
465 #endif
466 
467  AMREX_D_TERM(amrex::Real xc = ccc(i,j,k,0);, // centroid of cell (i,j,k)
468  amrex::Real yc = ccc(i,j,k,1);,
469  amrex::Real zc = ccc(i,j,k,2););
470 
471  AMREX_D_TERM(amrex::Real delta_x = xf - xc;,
472  amrex::Real delta_y = amrex::Real(0.5) + yc;,
473  amrex::Real delta_z = zf - zc;);
474 
475  amrex::Real cc_qmax = amrex::max(q(i,j,k,n),q(i,j-1,k,n));
476  amrex::Real cc_qmin = amrex::min(q(i,j,k,n),q(i,j-1,k,n));
477 
478  // Compute slopes of component "n" of q
479  const auto& slopes_eb_hi = amrex_lim_slopes_eb(i, j, k, n, q, ccc, vfrac,
480  AMREX_D_DECL(fcx,fcy,fcz),
481  flag, order);
482 
483 #if (AMREX_SPACEDIM==3)
484  amrex::Real qpls = q(i ,j,k,n) + delta_x * slopes_eb_hi[0]
485  - delta_y * slopes_eb_hi[1]
486  + delta_z * slopes_eb_hi[2];
487 #else
488  amrex::Real qpls = q(i ,j,k,n) + delta_x * slopes_eb_hi[0]
489  - delta_y * slopes_eb_hi[1];
490 #endif
491 
492  qpls = amrex::max(amrex::min(qpls, cc_qmax), cc_qmin);
493 
494  AMREX_D_TERM(xc = ccc(i,j-1,k,0);, // centroid of cell (i,j-1,k)
495  yc = ccc(i,j-1,k,1);,
496  zc = ccc(i,j-1,k,2););
497 
498  AMREX_D_TERM(delta_x = xf - xc;,
499  delta_y = amrex::Real(0.5) - yc;,
500  delta_z = zf - zc;);
501 
502  // Compute slopes of component "n" of q
503  const auto& slopes_eb_lo = amrex_lim_slopes_eb(i, j-1, k, n, q, ccc, vfrac,
504  AMREX_D_DECL(fcx,fcy,fcz),
505  flag, order);
506 
507 #if (AMREX_SPACEDIM==3)
508  amrex::Real qmns = q(i,j-1,k,n) + delta_x * slopes_eb_lo[0]
509  + delta_y * slopes_eb_lo[1]
510  + delta_z * slopes_eb_lo[2];
511 #else
512  amrex::Real qmns = q(i,j-1,k,n) + delta_x * slopes_eb_lo[0]
513  + delta_y * slopes_eb_lo[1];
514 #endif
515 
516  qmns = amrex::max(amrex::min(qmns, cc_qmax), cc_qmin);
517 
518  HydroBC::SetYEdgeBCs(i, j, k, n, q, qmns, qpls, d_bcrec[n].lo(1), domain_jlo, d_bcrec[n].hi(1), domain_jhi, is_velocity);
519 
520  if ( (j==domain_jlo) && (d_bcrec[n].lo(1) == amrex::BCType::foextrap || d_bcrec[n].lo(1) == amrex::BCType::hoextrap) )
521  {
522  if ( vmac(i,j,k) >= 0. && n==YVEL && is_velocity ) qpls = amrex::min(qpls,amrex::Real(0.0));
523  qmns = qpls;
524  }
525  if ( (j==domain_jhi+1) && (d_bcrec[n].hi(1) == amrex::BCType::foextrap || d_bcrec[n].hi(1) == amrex::BCType::hoextrap) )
526  {
527  if ( vmac(i,j,k) <= 0. && n==YVEL && is_velocity ) qmns = amrex::max(qmns,amrex::Real(0.0));
528  qpls = qmns;
529  }
530 
531  amrex::Real qs;
532 
533  if (vmac(i,j,k) > small_vel)
534  {
535  qs = qmns;
536  }
537  else if (vmac(i,j,k) < - small_vel)
538  {
539  qs = qpls;
540  }
541  else
542  {
543  qs = amrex::Real(0.5)*(qmns+qpls);
544  }
545 
546  return qs;
547 }
548 
549 
550 
551 
552 #if (AMREX_SPACEDIM==3)
553 
555 amrex::Real hydro_ebmol_zedge_state_extdir ( int i, int j, int k, int n,
564  amrex::BCRec const* const d_bcrec,
565  amrex::Box const& domain,
566  int order, const bool is_velocity) noexcept
567 {
568  amrex::Real qs;
569  int domlo = domain.smallEnd(2);
570  int domhi = domain.bigEnd(2);
571  bool extlo = d_bcrec[n].lo(2) == amrex::BCType::ext_dir;
572  bool exthi = d_bcrec[n].hi(2) == amrex::BCType::ext_dir;
573 
574  if ( extlo && k <= domlo)
575  {
576  qs = q(i,j,domlo-1,n);
577  }
578  else if ( exthi && k >= domhi+1)
579  {
580  qs = q(i,j,domhi+1,n);
581  }
582  else
583  {
584  const int domain_ilo = domain.smallEnd(0);
585  const int domain_ihi = domain.bigEnd(0);
586  const int domain_jlo = domain.smallEnd(1);
587  const int domain_jhi = domain.bigEnd(1);
588 #if (AMREX_SPACEDIM == 3)
589  const int domain_klo = domain.smallEnd(2);
590  const int domain_khi = domain.bigEnd(2);
591 #endif
592 
593  AMREX_D_TERM(bool extdir_or_ho_ilo = (d_bcrec[n].lo(0) == amrex::BCType::ext_dir) ||
594  (d_bcrec[n].lo(0) == amrex::BCType::hoextrap);,
595  bool extdir_or_ho_jlo = (d_bcrec[n].lo(1) == amrex::BCType::ext_dir) ||
596  (d_bcrec[n].lo(1) == amrex::BCType::hoextrap);,
597  bool extdir_or_ho_klo = (d_bcrec[n].lo(2) == amrex::BCType::ext_dir) ||
598  (d_bcrec[n].lo(2) == amrex::BCType::hoextrap););
599 
600  AMREX_D_TERM(bool extdir_or_ho_ihi = (d_bcrec[n].hi(0) == amrex::BCType::ext_dir) ||
601  (d_bcrec[n].hi(0) == amrex::BCType::hoextrap);,
602  bool extdir_or_ho_jhi = (d_bcrec[n].hi(1) == amrex::BCType::ext_dir) ||
603  (d_bcrec[n].hi(1) == amrex::BCType::hoextrap);,
604  bool extdir_or_ho_khi = (d_bcrec[n].hi(2) == amrex::BCType::ext_dir) ||
605  (d_bcrec[n].hi(2) == amrex::BCType::hoextrap););
606 
607  amrex::Real xf = fcz(i,j,k,0); // local (x,y) of centroid of z-face we are extrapolating to
608  amrex::Real yf = fcz(i,j,k,1);
609 
610  amrex::Real xc = ccc(i,j,k,0); // centroid of cell (i,j,k)
611  amrex::Real yc = ccc(i,j,k,1);
612  amrex::Real zc = ccc(i,j,k,2);
613 
614  amrex::Real delta_x = xf - xc;
615  amrex::Real delta_y = yf - yc;
616  amrex::Real delta_z = amrex::Real(0.5) + zc;
617 
618  amrex::Real cc_qmax = amrex::max(q(i,j,k,n),q(i,j,k-1,n));
619  amrex::Real cc_qmin = amrex::min(q(i,j,k,n),q(i,j,k-1,n));
620 
621  // Compute slopes of component "n" of q
622  const auto& slopes_eb_hi =
623  amrex_lim_slopes_extdir_eb(i, j, k, n, q, ccc, vfrac,
624  AMREX_D_DECL(fcx,fcy,fcz), flag,
625  AMREX_D_DECL(extdir_or_ho_ilo, extdir_or_ho_jlo, extdir_or_ho_klo),
626  AMREX_D_DECL(extdir_or_ho_ihi, extdir_or_ho_jhi, extdir_or_ho_khi),
627  AMREX_D_DECL(domain_ilo, domain_jlo, domain_klo),
628  AMREX_D_DECL(domain_ihi, domain_jhi, domain_khi),
629  order);
630 
631  amrex::Real qpls = q(i,j,k ,n) + delta_x * slopes_eb_hi[0]
632  + delta_y * slopes_eb_hi[1]
633  - delta_z * slopes_eb_hi[2];
634 
635  qpls = amrex::max(amrex::min(qpls, cc_qmax), cc_qmin);
636 
637  xc = ccc(i,j,k-1,0); // centroid of cell (i,j,k-1)
638  yc = ccc(i,j,k-1,1);
639  zc = ccc(i,j,k-1,2);
640 
641  delta_x = xf - xc;
642  delta_y = yf - yc;
643  delta_z = amrex::Real(0.5) - zc;
644 
645  // Compute slopes of component "n" of q
646  const auto& slopes_eb_lo =
647  amrex_lim_slopes_extdir_eb(i, j, k-1, n, q, ccc, vfrac,
648  AMREX_D_DECL(fcx,fcy,fcz), flag,
649  AMREX_D_DECL(extdir_or_ho_ilo, extdir_or_ho_jlo, extdir_or_ho_klo),
650  AMREX_D_DECL(extdir_or_ho_ihi, extdir_or_ho_jhi, extdir_or_ho_khi),
651  AMREX_D_DECL(domain_ilo, domain_jlo, domain_klo),
652  AMREX_D_DECL(domain_ihi, domain_jhi, domain_khi),
653  order);
654 
655  amrex::Real qmns = q(i,j,k-1,n) + delta_x * slopes_eb_lo[0]
656  + delta_y * slopes_eb_lo[1]
657  + delta_z * slopes_eb_lo[2];
658 
659  qmns = amrex::max(amrex::min(qmns, cc_qmax), cc_qmin);
660 
661  HydroBC::SetZEdgeBCs(i, j, k, n, q, qmns, qpls, d_bcrec[n].lo(2), domain_klo, d_bcrec[n].hi(2), domain_khi, is_velocity);
662 
663  if ( (k==domain_klo) && (d_bcrec[n].lo(2) == amrex::BCType::foextrap || d_bcrec[n].lo(2) == amrex::BCType::hoextrap) )
664  {
665  if ( wmac(i,j,k) >= 0. && n==ZVEL && is_velocity ) qpls = amrex::min(qpls,amrex::Real(0.0));
666  qmns = qpls;
667  }
668  if ( (k==domain_khi+1) && (d_bcrec[n].hi(2) == amrex::BCType::foextrap || d_bcrec[n].hi(2) == amrex::BCType::hoextrap) )
669  {
670  if ( wmac(i,j,k) <= 0. && n==ZVEL && is_velocity ) qmns = amrex::max(qmns,amrex::Real(0.0));
671  qpls = qmns;
672  }
673 
674  if (wmac(i,j,k) > small_vel)
675  {
676  qs = qmns;
677  }
678  else if (wmac(i,j,k) < -small_vel)
679  {
680  qs = qpls;
681  }
682  else
683  {
684  qs = amrex::Real(0.5)*(qmns+qpls);
685  }
686  }
687 
688  return qs;
689 }
690 
691 
692 
693 
695 amrex::Real hydro_ebmol_zedge_state ( int i, int j, int k, int n,
704  amrex::BCRec const* const d_bcrec,
705  amrex::Box const& domain,
706  int order, const bool is_velocity) noexcept
707 {
708  const int domain_klo = domain.smallEnd(2);
709  const int domain_khi = domain.bigEnd(2);
710 
711  amrex::Real xf = fcz(i,j,k,0); // local (x,y) of centroid of z-face we are extrapolating to
712  amrex::Real yf = fcz(i,j,k,1);
713 
714  amrex::Real xc = ccc(i,j,k,0); // centroid of cell (i,j,k)
715  amrex::Real yc = ccc(i,j,k,1);
716  amrex::Real zc = ccc(i,j,k,2);
717 
718  amrex::Real delta_x = xf - xc;
719  amrex::Real delta_y = yf - yc;
720  amrex::Real delta_z = amrex::Real(0.5) + zc;
721 
722  amrex::Real cc_qmax = amrex::max(q(i,j,k,n),q(i,j,k-1,n));
723  amrex::Real cc_qmin = amrex::min(q(i,j,k,n),q(i,j,k-1,n));
724 
725  // Compute slopes of component "n" of q
726  const auto& slopes_eb_hi = amrex_lim_slopes_eb(i, j, k, n, q, ccc, vfrac,
727  AMREX_D_DECL(fcx,fcy,fcz),
728  flag, order);
729 
730  amrex::Real qpls = q(i,j,k ,n) + delta_x * slopes_eb_hi[0]
731  + delta_y * slopes_eb_hi[1]
732  - delta_z * slopes_eb_hi[2];
733 
734  qpls = amrex::max(amrex::min(qpls, cc_qmax), cc_qmin);
735 
736  xc = ccc(i,j,k-1,0); // centroid of cell (i,j,k-1)
737  yc = ccc(i,j,k-1,1);
738  zc = ccc(i,j,k-1,2);
739 
740  delta_x = xf - xc;
741  delta_y = yf - yc;
742  delta_z = amrex::Real(0.5) - zc;
743 
744  // Compute slopes of component "n" of q
745  const auto& slopes_eb_lo = amrex_lim_slopes_eb(i, j, k-1, n, q, ccc, vfrac,
746  AMREX_D_DECL(fcx,fcy,fcz),
747  flag, order);
748 
749  amrex::Real qmns = q(i,j,k-1,n) + delta_x * slopes_eb_lo[0]
750  + delta_y * slopes_eb_lo[1]
751  + delta_z * slopes_eb_lo[2];
752 
753  qmns = amrex::max(amrex::min(qmns, cc_qmax), cc_qmin);
754 
755  HydroBC::SetZEdgeBCs(i, j, k, n, q, qmns, qpls, d_bcrec[n].lo(2), domain_klo, d_bcrec[n].hi(2), domain_khi, is_velocity);
756 
757  if ( (k==domain_klo) && (d_bcrec[n].lo(2) == amrex::BCType::foextrap || d_bcrec[n].lo(2) == amrex::BCType::hoextrap) )
758  {
759  if ( wmac(i,j,k) >= 0. && n==ZVEL && is_velocity ) qpls = amrex::min(qpls,amrex::Real(0.0));
760  qmns = qpls;
761  }
762  if ( (k==domain_khi+1) && (d_bcrec[n].hi(2) == amrex::BCType::foextrap || d_bcrec[n].hi(2) == amrex::BCType::hoextrap) )
763  {
764  if ( wmac(i,j,k) <= 0. && n==ZVEL && is_velocity ) qmns = amrex::max(qmns,amrex::Real(0.0));
765  qpls = qmns;
766  }
767 
768  amrex::Real qs;
769 
770  if (wmac(i,j,k) > small_vel)
771  {
772  qs = qmns;
773  }
774  else if (wmac(i,j,k) < -small_vel)
775  {
776  qs = qpls;
777  }
778  else
779  {
780  qs = amrex::Real(0.5)*(qmns+qpls);
781  }
782 
783 
784  return qs;
785 }
786 
787 #endif
788 
789 }
790 
791 #endif
792 /** @} */
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
#define AMREX_D_TERM(a, b, c)
#define AMREX_D_DECL(a, b, c)
#define YVEL
Definition: hydro_constants.H:29
static constexpr amrex::Real small_vel
Definition: hydro_constants.H:37
#define XVEL
Definition: hydro_constants.H:28
#define ZVEL
Definition: hydro_constants.H:30
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real hydro_ebmol_yedge_state_extdir(AMREX_D_DECL(int i, int j, int k), int n, amrex::Array4< amrex::Real const > const &q, amrex::Array4< amrex::Real const > const &vmac, AMREX_D_DECL(amrex::Array4< amrex::Real const > const &fcx, amrex::Array4< amrex::Real const > const &fcy, amrex::Array4< amrex::Real const > const &fcz), amrex::Array4< amrex::Real const > const &ccc, amrex::Array4< amrex::Real const > const &vfrac, amrex::Array4< amrex::EBCellFlag const > const &flag, amrex::BCRec const *const d_bcrec, amrex::Box const &domain, int order, const bool is_velocity) noexcept
Definition: hydro_ebmol_edge_state_K.H:284
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real hydro_ebmol_yedge_state(AMREX_D_DECL(int i, int j, int k), int n, amrex::Array4< amrex::Real const > const &q, amrex::Array4< amrex::Real const > const &vmac, AMREX_D_DECL(amrex::Array4< amrex::Real const > const &fcx, amrex::Array4< amrex::Real const > const &fcy, amrex::Array4< amrex::Real const > const &fcz), amrex::Array4< amrex::Real const > const &ccc, amrex::Array4< amrex::Real const > const &vfrac, amrex::Array4< amrex::EBCellFlag const > const &flag, amrex::BCRec const *const d_bcrec, amrex::Box const &domain, int order, const bool is_velocity) noexcept
Definition: hydro_ebmol_edge_state_K.H:441
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real hydro_ebmol_xedge_state(AMREX_D_DECL(int i, int j, int k), int n, amrex::Array4< amrex::Real const > const &q, amrex::Array4< amrex::Real const > const &umac, AMREX_D_DECL(amrex::Array4< amrex::Real const > const &fcx, amrex::Array4< amrex::Real const > const &fcy, amrex::Array4< amrex::Real const > const &fcz), amrex::Array4< amrex::Real const > const &ccc, amrex::Array4< amrex::Real const > const &vfrac, amrex::Array4< amrex::EBCellFlag const > const &flag, amrex::BCRec const *const d_bcrec, amrex::Box const &domain, int order, const bool is_velocity) noexcept
Definition: hydro_ebmol_edge_state_K.H:175
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real hydro_ebmol_xedge_state_extdir(AMREX_D_DECL(int i, int j, int k), int n, amrex::Array4< amrex::Real const > const &q, amrex::Array4< amrex::Real const > const &umac, AMREX_D_DECL(amrex::Array4< amrex::Real const > const &fcx, amrex::Array4< amrex::Real const > const &fcy, amrex::Array4< amrex::Real const > const &fcz), amrex::Array4< amrex::Real const > const &ccc, amrex::Array4< amrex::Real const > const &vfrac, amrex::Array4< amrex::EBCellFlag const > const &flag, amrex::BCRec const *const d_bcrec, amrex::Box const &domain, int order, const bool is_velocity) noexcept
Definition: hydro_ebmol_edge_state_K.H:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void SetXEdgeBCs(int i, int j, int k, int n, const amrex::Array4< const amrex::Real > &s, amrex::Real &lo, amrex::Real &hi, int bclo, int domlo, int bchi, int domhi, bool is_velocity)
Boundary condition effects.
Definition: hydro_bcs_K.H:45
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void SetYEdgeBCs(int i, int j, int k, int n, const amrex::Array4< const amrex::Real > &s, amrex::Real &lo, amrex::Real &hi, int bclo, int domlo, int bchi, int domhi, bool is_velocity)
Boundary condition effects.
Definition: hydro_bcs_K.H:112
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & min(const T &a, const T &b) noexcept