Finite Volume Solver  prototype
A framework to build finite volume solvers for the AG Klein at the Freie Universität Berlin.
cutcell/boundary_condition/ConstantBoundary.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_BOUNDARY_CUTCELL_CONDITION_CONSTANT_BOUNDARY_HPP
22 #define FUB_AMREX_BOUNDARY_CUTCELL_CONDITION_CONSTANT_BOUNDARY_HPP
23 
25 
26 #include "fub/AMReX/ForEachFab.hpp"
28 #include "fub/ForEach.hpp"
29 
30 #include <AMReX_MultiFab.H>
31 
32 namespace fub::amrex::cutcell {
33 
34 /// \ingroup BoundaryCondition
35 ///
36 template <typename Equation> struct ConstantBoundary {
38  int side;
41 
42  void FillBoundary(::amrex::MultiFab& mf, const GriddingAlgorithm& grid,
43  int level) {
44  const ::amrex::Geometry& geom = grid.GetPatchHierarchy().GetGeometry(level);
45  const std::shared_ptr<::amrex::EBFArrayBoxFactory>& factory =
47  FillBoundary(mf, geom, factory);
48  }
49 
50  void FillBoundary(::amrex::MultiFab& mf, const GriddingAlgorithm& gridding,
51  int level, Direction dir) {
52  if (dir == this->dir) {
53  FillBoundary(mf, gridding, level);
54  }
55  }
56 
57  void
58  FillBoundary(::amrex::MultiFab& mf, const ::amrex::Geometry& geom,
59  const std::shared_ptr<::amrex::EBFArrayBoxFactory>& factory) {
60  int idir = static_cast<int>(dir);
61  if (!geom.isPeriodic(idir)) {
62  const int ngrow = mf.nGrow(idir);
63  ::amrex::Box fully_grown_box = geom.growNonPeriodicDomain(ngrow);
64  ::amrex::Box grown_dir_box = ::amrex::grow(fully_grown_box, idir, -ngrow);
65  ::amrex::BoxList boundaries = ::amrex::complementIn(
66  fully_grown_box, ::amrex::BoxList{grown_dir_box});
67  if (boundaries.isEmpty()) {
68  return;
69  }
70  const ::amrex::MultiFab& volfrac = factory->getVolFrac();
71  ForEachFab(fub::execution::openmp, mf, [&](const ::amrex::MFIter& mfi) {
72  ::amrex::FArrayBox& fab = mf[mfi];
73  const ::amrex::FArrayBox& alpha = volfrac[mfi];
74  for (const ::amrex::Box& boundary : boundaries) {
75  ::amrex::Box shifted =
76  ::amrex::shift(boundary, idir, GetSign(side) * ngrow);
77  if (!geom.Domain().intersects(shifted)) {
78  continue;
79  }
80  ::amrex::Box box_to_fill = mfi.growntilebox() & boundary;
81  if (!box_to_fill.isEmpty()) {
82  auto states =
83  MakeView<Complete<Equation>>(fab, equation, box_to_fill);
84  ForEachIndex(box_to_fill, [&](auto... is) {
85  ::amrex::IntVect iv{static_cast<int>(is)...};
86  if (alpha(iv) > 0.0) {
87  Store(states, state, {is...});
88  } else {
89  for (int comp = 0; comp < fab.nComp(); ++comp) {
90  fab(iv, comp) = 0.0;
91  }
92  }
93  });
94  }
95  }
96  });
97  }
98  }
99 };
100 
101 } // namespace fub::amrex::cutcell
102 
103 #endif
This class modifies and initializes a cutcell::PatchLevel in a cutcell::PatchHierarchy.
Definition: AMReX/cutcell/GriddingAlgorithm.hpp:56
const PatchHierarchy & GetPatchHierarchy() const noexcept
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
void ForEachIndex(const ::amrex::Box &box, F function)
Definition: ForEachIndex.hpp:29
constexpr OpenMpTag openmp
Definition: Execution.hpp:37
std::decay_t< decltype(std::declval< T >().GetEquation())> Equation
A template typedef to detect the member function.
Definition: Meta.hpp:59
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 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: cutcell/boundary_condition/ConstantBoundary.hpp:36
void FillBoundary(::amrex::MultiFab &mf, const ::amrex::Geometry &geom, const std::shared_ptr<::amrex::EBFArrayBoxFactory > &factory)
Definition: cutcell/boundary_condition/ConstantBoundary.hpp:58
Complete< Equation > state
Definition: cutcell/boundary_condition/ConstantBoundary.hpp:40
void FillBoundary(::amrex::MultiFab &mf, const GriddingAlgorithm &gridding, int level, Direction dir)
Definition: cutcell/boundary_condition/ConstantBoundary.hpp:50
Equation equation
Definition: cutcell/boundary_condition/ConstantBoundary.hpp:39
int side
Definition: cutcell/boundary_condition/ConstantBoundary.hpp:38
void FillBoundary(::amrex::MultiFab &mf, const GriddingAlgorithm &grid, int level)
Definition: cutcell/boundary_condition/ConstantBoundary.hpp:42
Direction dir
Definition: cutcell/boundary_condition/ConstantBoundary.hpp:37