Finite Volume Solver  prototype
A framework to build finite volume solvers for the AG Klein at the Freie Universität Berlin.
AMReX/cutcell/boundary_condition/ReflectiveBoundary.hpp
Go to the documentation of this file.
1 // Copyright (c) 2019 Maikel Nadolski
2 //
3 // Permission is hereby granted, free of charge, to any person obtaining a copy
4 // of this software and associated documentation files (the "Software"), to deal
5 // in the Software without restriction, including without limitation the rights
6 // to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
7 // copies of the Software, and to permit persons to whom the Software is
8 // furnished to do so, subject to the following conditions:
9 //
10 // The above copyright notice and this permission notice shall be included in
11 // all copies or substantial portions of the Software.
12 //
13 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
14 // IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
15 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
16 // AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
17 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
18 // OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
19 // SOFTWARE.
20 
21 #ifndef FUB_AMREX_CUTCELL_BOUNDARY_CONDITION_REFLECTIVE_BOUNDARY_HPP
22 #define FUB_AMREX_CUTCELL_BOUNDARY_CONDITION_REFLECTIVE_BOUNDARY_HPP
23 
25 
26 #include "fub/AMReX/ForEachFab.hpp"
28 #include "fub/Execution.hpp"
30 
31 namespace fub::amrex::cutcell {
32 
33 /// \ingroup BoundaryCondition
34 ///
35 template <typename Tag, typename Equation> class ReflectiveBoundary {
36 public:
37  ReflectiveBoundary(Tag, const Equation& equation, Direction dir, int side);
38 
39  void FillBoundary(::amrex::MultiFab& mf, const GriddingAlgorithm& gridding,
40  int level);
41 
42  void FillBoundary(::amrex::MultiFab& mf, const GriddingAlgorithm& gridding,
43  int level, Direction dir);
44 
45 private:
48  int side_;
49 };
50 
51 template <typename Tag, typename Equation>
53  const Equation& equation,
54  Direction dir, int side)
55  : equation_{Local<Tag, Equation>{equation}}, dir_{dir}, side_{side} {}
56 
57 template <typename Tag, typename Equation>
59  ::amrex::MultiFab& mf, const GriddingAlgorithm& gridding, int level,
60  Direction dir) {
61  if (dir == dir_) {
62  FillBoundary(mf, gridding, level);
63  }
64 }
65 
66 template <typename Tag, typename Equation>
68  ::amrex::MultiFab& mf, const GriddingAlgorithm& grid, int level) {
69  const ::amrex::Geometry& geom = grid.GetPatchHierarchy().GetGeometry(level);
70  const int ngrow = mf.nGrow(int(dir_));
71  ::amrex::Box grown_box = geom.growNonPeriodicDomain(ngrow);
72  ::amrex::BoxList boundaries =
73  ::amrex::complementIn(grown_box, ::amrex::BoxList{geom.Domain()});
74  if (boundaries.isEmpty()) {
75  return;
76  }
77  static constexpr int Rank = Equation::Rank();
78  static constexpr std::size_t sRank = static_cast<std::size_t>(Rank);
79  const Eigen::Matrix<double, Rank, 1> unit = UnitVector<Rank>(dir_);
80  const ::amrex::MultiFab& alphas =
81  grid.GetPatchHierarchy().GetEmbeddedBoundary(level)->getVolFrac();
82  ForEachFab(Tag(), mf, [&](const ::amrex::MFIter& mfi) {
83  Complete<Equation> state(*equation_);
84  Complete<Equation> reflected(*equation_);
85  Complete<Equation> zeros(*equation_);
86  ::amrex::FArrayBox& fab = mf[mfi];
87  const ::amrex::FArrayBox& alpha = alphas[mfi];
88  for (const ::amrex::Box& boundary : boundaries) {
89  ::amrex::Box shifted =
90  ::amrex::shift(boundary, int(dir_), GetSign(side_) * ngrow);
91  if (!geom.Domain().intersects(shifted)) {
92  continue;
93  }
94  ::amrex::Box box_to_fill = mfi.growntilebox() & boundary;
95  if (!box_to_fill.isEmpty()) {
96  auto states =
97  MakeView<Complete<Equation>>(fab, *equation_, mfi.growntilebox());
98  auto box = AsIndexBox<Rank>(box_to_fill);
99  ForEachIndex(box, [&](auto... is) {
100  std::array<std::ptrdiff_t, sRank> dest{is...};
101  std::array<std::ptrdiff_t, sRank> src =
102  ReflectIndex(dest, box, dir_, side_);
103  ::amrex::IntVect iv{
104  AMREX_D_DECL(int(src[0]), int(src[1]), int(src[2]))};
105  ::amrex::IntVect dest_iv{
106  AMREX_D_DECL(int(dest[0]), int(dest[1]), int(dest[2]))};
107  if (alpha(dest_iv) > 0.0 && alpha(iv) > 0.0) {
108  Load(state, states, src);
109  FUB_ASSERT(state.density > 0.0);
110  Reflect(reflected, state, unit, *equation_);
111  Store(states, reflected, dest);
112  } else {
113  Store(states, zeros, dest);
114  }
115  });
116  }
117  }
118  });
119 }
120 
121 template <typename Tag, typename Equation>
123  ->ReflectiveBoundary<Tag, Equation>;
124 
125 } // namespace fub::amrex::cutcell
126 
127 #endif
#define FUB_ASSERT(x)
Definition: assert.hpp:39
const ::amrex::Geometry & GetGeometry(int level) const
Returns a Geometry object for a specified level.
void FillBoundary(::amrex::MultiFab &mf, const GriddingAlgorithm &gridding, int level)
Definition: AMReX/boundary_condition/ReflectiveBoundary.hpp:97
This class modifies and initializes a cutcell::PatchLevel in a cutcell::PatchHierarchy.
Definition: AMReX/cutcell/GriddingAlgorithm.hpp:56
const PatchHierarchy & GetPatchHierarchy() const noexcept
Definition: AMReX/cutcell/boundary_condition/ReflectiveBoundary.hpp:35
ReflectiveBoundary(Tag, const Equation &equation, Direction dir, int side)
Definition: AMReX/cutcell/boundary_condition/ReflectiveBoundary.hpp:52
int side_
Definition: AMReX/cutcell/boundary_condition/ReflectiveBoundary.hpp:48
void FillBoundary(::amrex::MultiFab &mf, const GriddingAlgorithm &gridding, int level)
Definition: AMReX/cutcell/boundary_condition/ReflectiveBoundary.hpp:67
Direction dir_
Definition: AMReX/cutcell/boundary_condition/ReflectiveBoundary.hpp:47
Local< Tag, Equation > equation_
Definition: AMReX/cutcell/boundary_condition/ReflectiveBoundary.hpp:46
ReflectiveBoundary(const Equation &, Direction, int) -> ReflectiveBoundary< execution::SequentialTag, Equation >
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
void ForEachIndex(const ::amrex::Box &box, F function)
Definition: ForEachIndex.hpp:29
std::decay_t< decltype(std::declval< T >().GetEquation())> Equation
A template typedef to detect the member function.
Definition: Meta.hpp:59
Index< 1 > ReflectIndex(Index< 1 > i, const IndexBox< 1 > &domain, Direction dir, int side)
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
void Reflect(Complete< IdealGasMix< 1 >> &reflected, const Complete< IdealGasMix< 1 >> &state, const Eigen::Matrix< double, 1, 1 > &normal, const IdealGasMix< 1 > &gas)
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
typename detail::LocalType< Tag, T >::type Local
Definition: Execution.hpp:56