Finite Volume Solver  prototype
A framework to build finite volume solvers for the AG Klein at the Freie Universität Berlin.
cutcell/boundary_condition/IsentropicPressureExpansion.hpp
Go to the documentation of this file.
1 // Copyright (c) 2020 Maikel Nadolski
2 // Copyright (c) 2020 Christian Zenker
3 //
4 // Permission is hereby granted, free of charge, to any person obtaining a copy
5 // of this software and associated documentation files (the "Software"), to deal
6 // in the Software without restriction, including without limitation the rights
7 // to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
8 // copies of the Software, and to permit persons to whom the Software is
9 // furnished to do so, subject to the following conditions:
10 //
11 // The above copyright notice and this permission notice shall be included in
12 // all copies or substantial portions of the Software.
13 //
14 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
15 // IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
16 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
17 // AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
18 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
19 // OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
20 // SOFTWARE.
21 
22 #ifndef FUB_AMREX_CUTCELL_BOUNDARY_CONDITION_ISENTROPIC_PRESSURE_EXPANSION_HPP
23 #define FUB_AMREX_CUTCELL_BOUNDARY_CONDITION_ISENTROPIC_PRESSURE_EXPANSION_HPP
24 
28 
29 namespace fub::amrex::cutcell {
30 
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);
40  // src[std::size_t(dir)] = boundary;
41  return src;
42 }
43 
44 /// \ingroup BoundaryCondition
45 ///
46 /// \brief This boundary models an isentropic pressure expansion for the
47 /// one-dimensional ideal gas equations for mixtures.
48 template <typename EulerEquation> class IsentropicPressureExpansion {
49 public:
51  const EulerEquation& eq,
53 
54  void FillBoundary(::amrex::MultiFab& mf, const GriddingAlgorithm& gridding,
55  int level);
56 
57  void FillBoundary(::amrex::MultiFab& mf, const GriddingAlgorithm& gridding,
58  int level, Direction dir);
59 
60  void FillBoundary(::amrex::MultiFab& mf, const ::amrex::MultiFab& alphas,
61  const ::amrex::Geometry& geom);
62 
63 private:
64  EulerEquation equation_;
66 };
67 
68 template <typename EulerEquation>
69 void ExpandState(EulerEquation& eq, Complete<EulerEquation>& dest,
70  const Complete<EulerEquation>& src, double pressure_dest,
71  double efficiency) {
72  static constexpr int N = EulerEquation::Rank();
73  const Array<double, N, 1> old_velocity = euler::Velocity(eq, src);
74  dest = src;
76  const double h_before = euler::TotalEnthalpy(eq, dest);
77  euler::SetIsentropicPressure(eq, dest, dest, pressure_dest);
78  const double h_after = euler::TotalEnthalpy(eq, dest);
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));
83  Array<double, N, 1> u_border = old_velocity;
84  u_border[0] = u_border0;
85  euler::SetVelocity(eq, dest, u_border);
86 }
87 
88 template <typename EulerEquation>
90  const EulerEquation& eq,
92  : equation_{eq}, options_{options} {}
93 
94 template <typename EulerEquation>
96  ::amrex::MultiFab& mf, const GriddingAlgorithm& grid, int level) {
97  auto factory = grid.GetPatchHierarchy().GetEmbeddedBoundary(level);
98  const ::amrex::MultiFab& alphas = factory->getVolFrac();
99  FillBoundary(mf, alphas, grid.GetPatchHierarchy().GetGeometry(level));
100 }
101 
102 template <typename EulerEquation>
104  ::amrex::MultiFab& mf, const GriddingAlgorithm& gridding, int level,
105  Direction dir) {
106  if (dir == options_.direction) {
107  FillBoundary(mf, gridding, level);
108  }
109 }
110 
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()) {
120  return;
121  }
122  Complete<EulerEquation> state{equation_};
123  Complete<EulerEquation> expanded{equation_};
124  ForEachFab(execution::seq, mf, [&](const ::amrex::MFIter& mfi) {
125  ::amrex::FArrayBox& fab = mf[mfi];
126  const ::amrex::FArrayBox& alpha = alphas[mfi];
127  for (const ::amrex::Box& boundary : boundaries) {
128  ::amrex::Box shifted = ::amrex::shift(boundary, int(options_.direction),
129  GetSign(options_.side) * ngrow);
130  if (!geom.Domain().intersects(shifted)) {
131  continue;
132  }
133  ::amrex::Box box_to_fill = mfi.growntilebox() & boundary;
134  if (!box_to_fill.isEmpty()) {
135  auto states = MakeView<Complete<EulerEquation>>(fab, equation_,
136  mfi.growntilebox());
137  ForEachIndex(
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]))};
144  ::amrex::IntVect iv{
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);
151  }
152  });
153  }
154  }
155  });
156 }
157 
158 } // namespace fub::amrex::cutcell
159 
160 #endif
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