22 #ifndef FUB_AMREX_CUTCELL_BOUNDARY_CONDITION_ISENTROPIC_PRESSURE_EXPANSION_HPP
23 #define FUB_AMREX_CUTCELL_BOUNDARY_CONDITION_ISENTROPIC_PRESSURE_EXPANSION_HPP
31 inline std::array<std::ptrdiff_t, 2>
32 MapToSrc(
const std::array<std::ptrdiff_t, 2>& dest,
33 const ::amrex::Geometry& geom,
int side,
Direction dir) {
34 const int boundary = (side == 0) * geom.Domain().smallEnd(
int(dir)) +
35 (side == 1) * geom.Domain().bigEnd(
int(dir));
36 const auto distance = dest[std::size_t(dir)] - boundary;
37 const int sign = int((distance > 0) - (distance < 0));
38 std::array<std::ptrdiff_t, 2> src{dest};
39 src[std::size_t(dir)] -=
static_cast<int>(2 * distance - sign);
51 const EulerEquation& eq,
60 void FillBoundary(::amrex::MultiFab& mf, const ::amrex::MultiFab& alphas,
61 const ::amrex::Geometry& geom);
68 template <
typename EulerEquation>
72 static constexpr
int N = EulerEquation::Rank();
79 const double h_diff = h_before - h_after;
80 auto Sign = [](
double x) {
return (x >= 0) - (x < 0); };
81 const double e_kin = efficiency * h_diff * 2.0 + old_velocity[0] * old_velocity[0];
82 const double u_border0 = Sign(e_kin) * std::sqrt(std::abs(e_kin));
84 u_border[0] = u_border0;
88 template <
typename EulerEquation>
90 const EulerEquation& eq,
92 : equation_{eq}, options_{options} {}
94 template <
typename EulerEquation>
98 const ::amrex::MultiFab& alphas = factory->getVolFrac();
102 template <
typename EulerEquation>
106 if (dir == options_.direction) {
107 FillBoundary(mf, gridding, level);
111 template <
typename EulerEquation>
113 ::amrex::MultiFab& mf, const ::amrex::MultiFab& alphas,
114 const ::amrex::Geometry& geom) {
115 const int ngrow = mf.nGrow(
int(options_.direction));
116 ::amrex::Box grown_box = geom.growNonPeriodicDomain(ngrow);
117 ::amrex::BoxList boundaries =
118 ::amrex::complementIn(grown_box, ::amrex::BoxList{geom.Domain()});
119 if (boundaries.isEmpty()) {
125 ::amrex::FArrayBox& fab = mf[mfi];
126 const ::amrex::FArrayBox& alpha = alphas[mfi];
128 ::amrex::Box shifted = ::amrex::shift(boundary,
int(options_.direction),
129 GetSign(options_.side) * ngrow);
130 if (!geom.Domain().intersects(shifted)) {
133 ::amrex::Box box_to_fill = mfi.growntilebox() & boundary;
134 if (!box_to_fill.isEmpty()) {
135 auto states = MakeView<Complete<EulerEquation>>(fab, equation_,
138 AsIndexBox<EulerEquation::Rank()>(box_to_fill), [&](
auto... is) {
139 Index<EulerEquation::Rank()> dest{is...};
140 Index<EulerEquation::Rank()> src =
141 MapToSrc(dest, geom, options_.side, options_.direction);
142 ::amrex::IntVect src_iv{
143 AMREX_D_DECL(
int(src[0]),
int(src[1]),
int(src[2]))};
145 AMREX_D_DECL(
int(dest[0]),
int(dest[1]),
int(dest[2]))};
146 if (alpha(iv) > 0.0 && alpha(src_iv) > 0.0) {
147 Load(state, states, src);
148 ExpandState(equation_, expanded, state, options_.outer_pressure,
149 options_.efficiency);
150 Store(states, expanded, dest);
This class modifies and initializes a cutcell::PatchLevel in a cutcell::PatchHierarchy.
Definition: AMReX/cutcell/GriddingAlgorithm.hpp:56
const PatchHierarchy & GetPatchHierarchy() const noexcept
This boundary models an isentropic pressure expansion for the one-dimensional ideal gas equations for...
Definition: cutcell/boundary_condition/IsentropicPressureExpansion.hpp:48
IsentropicPressureExpansion(const EulerEquation &eq, const fub::amrex::IsentropicPressureBoundaryOptions &options)
Definition: cutcell/boundary_condition/IsentropicPressureExpansion.hpp:89
void FillBoundary(::amrex::MultiFab &mf, const GriddingAlgorithm &gridding, int level)
Definition: cutcell/boundary_condition/IsentropicPressureExpansion.hpp:95
EulerEquation equation_
Definition: cutcell/boundary_condition/IsentropicPressureExpansion.hpp:64
fub::amrex::IsentropicPressureBoundaryOptions options_
Definition: cutcell/boundary_condition/IsentropicPressureExpansion.hpp:65
const std::shared_ptr<::amrex::EBFArrayBoxFactory > & GetEmbeddedBoundary(int level) const
Returns the number of cycles done for a specific level.
const ::amrex::Geometry & GetGeometry(int level) const
Returns the number of cycles done for a specific level.
void ForEachFab(Tag, const ::amrex::FabArrayBase &fabarray, F function)
Iterate through all local FArrayBox objects in a MultiFab.
Definition: ForEachFab.hpp:34
Definition: FillCutCellData.hpp:30
std::array< std::ptrdiff_t, 2 > MapToSrc(const std::array< std::ptrdiff_t, 2 > &dest, const ::amrex::Geometry &geom, int side, Direction dir)
Definition: cutcell/boundary_condition/IsentropicPressureExpansion.hpp:32
void ExpandState(EulerEquation &eq, Complete< EulerEquation > &dest, const Complete< EulerEquation > &src, double pressure_dest, double efficiency)
Definition: cutcell/boundary_condition/IsentropicPressureExpansion.hpp:69
void ForEachIndex(const ::amrex::Box &box, F function)
Definition: ForEachIndex.hpp:29
IndexBox< Rank > AsIndexBox(const ::amrex::Box &box)
Definition: ViewFArrayBox.hpp:51
constexpr struct fub::euler::TotalEnthalpyFn TotalEnthalpy
constexpr struct fub::euler::SetIsentropicPressureFn SetIsentropicPressure
constexpr struct fub::euler::VelocityFn Velocity
constexpr struct fub::euler::SetVelocityFn SetVelocity
constexpr SequentialTag seq
Definition: Execution.hpp:31
std::conditional_t< N==1||M==1, Eigen::Array< T, N, M >, Eigen::Array< T, N, M, Eigen::RowMajor > > Array
Definition: Eigen.hpp:50
void Load(State &state, const BasicView< const State, Layout, Rank > &view, const std::array< std::ptrdiff_t, State::Equation::Rank()> &index)
Definition: State.hpp:640
int GetSign(int side)
Definition: AnyBoundaryCondition.hpp:104
Direction
This is a type safe type to denote a dimensional split direction.
Definition: Direction.hpp:30
std::array< std::ptrdiff_t, static_cast< std::size_t >(Rank)> Index
Definition: PatchDataView.hpp:34
void Store(const BasicView< Conservative< Eq >, Layout, Eq::Rank()> &view, const Conservative< Eq > &state, const std::array< std::ptrdiff_t, Eq::Rank()> &index)
Definition: State.hpp:663
IndexBox< Rank > Box(const BasicView< State, Layout, Rank > &view)
Definition: State.hpp:486
Definition: boundary_condition/IsentropicPressureBoundary.hpp:38