8 using namespace amrex::literals;
10 #ifndef HYDRO_EBMOL_EDGE_STATE_K_H_
11 #define HYDRO_EBMOL_EDGE_STATE_K_H_
32 int order,
const bool is_velocity) noexcept
35 #if (AMREX_SPACEDIM==2)
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;
45 if ( extlo && i <= domlo)
47 qs = q(domlo-1,j,k,n);
49 else if ( exthi && i >= domhi+1)
51 qs = q(domhi+1,j,k,n);
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);
64 amrex::Real yf = fcx(i,j,k,0);
65 #if (AMREX_SPACEDIM==3)
66 amrex::Real zf = fcx(i,j,k,1);
70 amrex::Real yc = ccc(i,j,k,1);,
71 amrex::Real zc = ccc(i,j,k,2););
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;);
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));
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););
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););
94 const auto& slopes_eb_hi =
95 amrex_lim_slopes_extdir_eb(i, j, k, n, q, ccc, vfrac,
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),
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];
109 amrex::Real qpls = q(i ,j,k,n) - delta_x * slopes_eb_hi[0]
110 + delta_y * slopes_eb_hi[1];
116 yc = ccc(i-1,j,k,1);,
117 zc = ccc(i-1,j,k,2););
124 const auto& slopes_eb_lo =
125 amrex_lim_slopes_extdir_eb(i-1, j, k, n, q, ccc, vfrac,
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),
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];
138 amrex::Real qmns = q(i-1,j,k,n) + delta_x * slopes_eb_lo[0]
139 + delta_y * slopes_eb_lo[1];
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);
146 if ( (i==domain_ilo) && (d_bcrec[n].lo(0) == amrex::BCType::foextrap || d_bcrec[n].lo(0) == amrex::BCType::hoextrap) )
148 if ( umac(i,j,k) >= 0. && n==
XVEL && is_velocity ) qpls =
amrex::min(qpls,amrex::Real(0.0));
151 if ( (i==domain_ihi+1) && (d_bcrec[n].hi(0) == amrex::BCType::foextrap || d_bcrec[n].hi(0) == amrex::BCType::hoextrap) )
153 if ( umac(i,j,k) <= 0. && n==
XVEL && is_velocity ) qmns =
amrex::max(qmns,amrex::Real(0.0));
167 qs = amrex::Real(0.5)*(qmns+qpls);
186 int order,
const bool is_velocity) noexcept
188 #if (AMREX_SPACEDIM==2)
192 const int domain_ilo = domain.smallEnd(0);
193 const int domain_ihi = domain.bigEnd(0);
196 amrex::Real yf = fcx(i,j,k,0);
197 #if (AMREX_SPACEDIM==3)
198 amrex::Real zf = fcx(i,j,k,1);
202 amrex::Real yc = ccc(i,j,k,1);,
203 amrex::Real zc = ccc(i,j,k,2););
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;);
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));
213 const auto& slopes_eb_hi = amrex_lim_slopes_eb(i, j, k, n, q, ccc, vfrac,
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];
222 amrex::Real qpls = q(i ,j,k,n) - delta_x * slopes_eb_hi[0]
223 + delta_y * slopes_eb_hi[1];
229 yc = ccc(i-1,j,k,1);,
230 zc = ccc(i-1,j,k,2););
237 const auto& slopes_eb_lo = amrex_lim_slopes_eb(i-1, j, k, n, q, ccc, vfrac,
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];
246 amrex::Real qmns = q(i-1,j,k,n) + delta_x * slopes_eb_lo[0]
247 + delta_y * slopes_eb_lo[1];
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);
254 if ( (i==domain_ilo) && (d_bcrec[n].lo(0) == amrex::BCType::foextrap || d_bcrec[n].lo(0) == amrex::BCType::hoextrap) )
256 if ( umac(i,j,k) >= 0. && n==
XVEL && is_velocity ) qpls =
amrex::min(qpls,amrex::Real(0.0));
259 if ( (i==domain_ihi+1) && (d_bcrec[n].hi(0) == amrex::BCType::foextrap || d_bcrec[n].hi(0) == amrex::BCType::hoextrap) )
261 if ( umac(i,j,k) <= 0. && n==
XVEL && is_velocity ) qmns =
amrex::max(qmns,amrex::Real(0.0));
277 qs = amrex::Real(0.5)*(qmns+qpls);
295 int order,
const bool is_velocity) noexcept
298 #if (AMREX_SPACEDIM==2)
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;
308 if ( extlo && j <= domlo)
310 qs = q(i,domlo-1,k,n);
312 else if (exthi && j >= domhi+1)
314 qs = q(i,domhi+1,k,n);
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);
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););
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););
341 amrex::Real xf = fcy(i,j,k,0);
342 #if (AMREX_SPACEDIM==3)
343 amrex::Real zf = fcy(i,j,k,1);
347 amrex::Real yc = ccc(i,j,k,1);,
348 amrex::Real zc = ccc(i,j,k,2););
351 amrex::Real delta_y = amrex::Real(0.5) + yc;,
352 amrex::Real delta_z = zf - zc;);
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));
358 const auto& slopes_eb_hi =
359 amrex_lim_slopes_extdir_eb(i, j, k, n, q, ccc, vfrac,
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),
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];
372 amrex::Real qpls = q(i ,j,k,n) + delta_x * slopes_eb_hi[0]
373 - delta_y * slopes_eb_hi[1];
379 yc = ccc(i,j-1,k,1);,
380 zc = ccc(i,j-1,k,2););
383 delta_y = amrex::Real(0.5) - yc;,
387 const auto& slopes_eb_lo =
388 amrex_lim_slopes_extdir_eb(i, j-1, k, n, q, ccc, vfrac,
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),
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];
402 amrex::Real qmns = q(i,j-1,k,n) + delta_x * slopes_eb_lo[0]
403 + delta_y * slopes_eb_lo[1];
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);
408 if ( (j==domain_jlo) && (d_bcrec[n].lo(1) == amrex::BCType::foextrap || d_bcrec[n].lo(1) == amrex::BCType::hoextrap) )
410 if ( vmac(i,j,k) >= 0. && n==
YVEL && is_velocity ) qpls =
amrex::min(qpls,amrex::Real(0.0));
413 if ( (j==domain_jhi+1) && (d_bcrec[n].hi(1) == amrex::BCType::foextrap || d_bcrec[n].hi(1) == amrex::BCType::hoextrap) )
415 if ( vmac(i,j,k) <= 0. && n==
YVEL && is_velocity ) qmns =
amrex::max(qmns,amrex::Real(0.0));
431 qs = amrex::Real(0.5)*(qmns+qpls);
452 int order,
const bool is_velocity) noexcept
454 #if (AMREX_SPACEDIM==2)
458 const int domain_jlo = domain.smallEnd(1);
459 const int domain_jhi = domain.bigEnd(1);
462 amrex::Real xf = fcy(i,j,k,0);
463 #if (AMREX_SPACEDIM==3)
464 amrex::Real zf = fcy(i,j,k,1);
468 amrex::Real yc = ccc(i,j,k,1);,
469 amrex::Real zc = ccc(i,j,k,2););
472 amrex::Real delta_y = amrex::Real(0.5) + yc;,
473 amrex::Real delta_z = zf - zc;);
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));
479 const auto& slopes_eb_hi = amrex_lim_slopes_eb(i, j, k, n, q, ccc, vfrac,
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];
488 amrex::Real qpls = q(i ,j,k,n) + delta_x * slopes_eb_hi[0]
489 - delta_y * slopes_eb_hi[1];
495 yc = ccc(i,j-1,k,1);,
496 zc = ccc(i,j-1,k,2););
499 delta_y = amrex::Real(0.5) - yc;,
503 const auto& slopes_eb_lo = amrex_lim_slopes_eb(i, j-1, k, n, q, ccc, vfrac,
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];
512 amrex::Real qmns = q(i,j-1,k,n) + delta_x * slopes_eb_lo[0]
513 + delta_y * slopes_eb_lo[1];
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);
520 if ( (j==domain_jlo) && (d_bcrec[n].lo(1) == amrex::BCType::foextrap || d_bcrec[n].lo(1) == amrex::BCType::hoextrap) )
522 if ( vmac(i,j,k) >= 0. && n==
YVEL && is_velocity ) qpls =
amrex::min(qpls,amrex::Real(0.0));
525 if ( (j==domain_jhi+1) && (d_bcrec[n].hi(1) == amrex::BCType::foextrap || d_bcrec[n].hi(1) == amrex::BCType::hoextrap) )
527 if ( vmac(i,j,k) <= 0. && n==
YVEL && is_velocity ) qmns =
amrex::max(qmns,amrex::Real(0.0));
543 qs = amrex::Real(0.5)*(qmns+qpls);
552 #if (AMREX_SPACEDIM==3)
555 amrex::Real hydro_ebmol_zedge_state_extdir (
int i,
int j,
int k,
int n,
566 int order,
const bool is_velocity) noexcept
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;
574 if ( extlo && k <= domlo)
576 qs = q(i,j,domlo-1,n);
578 else if ( exthi && k >= domhi+1)
580 qs = q(i,j,domhi+1,n);
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);
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););
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););
607 amrex::Real xf = fcz(i,j,k,0);
608 amrex::Real yf = fcz(i,j,k,1);
610 amrex::Real xc = ccc(i,j,k,0);
611 amrex::Real yc = ccc(i,j,k,1);
612 amrex::Real zc = ccc(i,j,k,2);
614 amrex::Real delta_x = xf - xc;
615 amrex::Real delta_y = yf - yc;
616 amrex::Real delta_z = amrex::Real(0.5) + zc;
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));
622 const auto& slopes_eb_hi =
623 amrex_lim_slopes_extdir_eb(i, j, k, n, q, ccc, vfrac,
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),
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];
643 delta_z = amrex::Real(0.5) - zc;
646 const auto& slopes_eb_lo =
647 amrex_lim_slopes_extdir_eb(i, j, k-1, n, q, ccc, vfrac,
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),
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];
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);
663 if ( (k==domain_klo) && (d_bcrec[n].lo(2) == amrex::BCType::foextrap || d_bcrec[n].lo(2) == amrex::BCType::hoextrap) )
665 if ( wmac(i,j,k) >= 0. && n==
ZVEL && is_velocity ) qpls =
amrex::min(qpls,amrex::Real(0.0));
668 if ( (k==domain_khi+1) && (d_bcrec[n].hi(2) == amrex::BCType::foextrap || d_bcrec[n].hi(2) == amrex::BCType::hoextrap) )
670 if ( wmac(i,j,k) <= 0. && n==
ZVEL && is_velocity ) qmns =
amrex::max(qmns,amrex::Real(0.0));
684 qs = amrex::Real(0.5)*(qmns+qpls);
695 amrex::Real hydro_ebmol_zedge_state (
int i,
int j,
int k,
int n,
706 int order,
const bool is_velocity) noexcept
708 const int domain_klo = domain.smallEnd(2);
709 const int domain_khi = domain.bigEnd(2);
711 amrex::Real xf = fcz(i,j,k,0);
712 amrex::Real yf = fcz(i,j,k,1);
714 amrex::Real xc = ccc(i,j,k,0);
715 amrex::Real yc = ccc(i,j,k,1);
716 amrex::Real zc = ccc(i,j,k,2);
718 amrex::Real delta_x = xf - xc;
719 amrex::Real delta_y = yf - yc;
720 amrex::Real delta_z = amrex::Real(0.5) + zc;
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));
726 const auto& slopes_eb_hi = amrex_lim_slopes_eb(i, j, k, n, q, ccc, vfrac,
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];
742 delta_z = amrex::Real(0.5) - zc;
745 const auto& slopes_eb_lo = amrex_lim_slopes_eb(i, j, k-1, n, q, ccc, vfrac,
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];
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);
757 if ( (k==domain_klo) && (d_bcrec[n].lo(2) == amrex::BCType::foextrap || d_bcrec[n].lo(2) == amrex::BCType::hoextrap) )
759 if ( wmac(i,j,k) >= 0. && n==
ZVEL && is_velocity ) qpls =
amrex::min(qpls,amrex::Real(0.0));
762 if ( (k==domain_khi+1) && (d_bcrec[n].hi(2) == amrex::BCType::foextrap || d_bcrec[n].hi(2) == amrex::BCType::hoextrap) )
764 if ( wmac(i,j,k) <= 0. && n==
ZVEL && is_velocity ) qmns =
amrex::max(qmns,amrex::Real(0.0));
780 qs = amrex::Real(0.5)*(qmns+qpls);
#define AMREX_FORCE_INLINE
#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