Finite Volume Solver  prototype
A framework to build finite volume solvers for the AG Klein at the Freie Universität Berlin.
Reconstruction.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_RECONSTRUCTION_HPP
22 #define FUB_AMREX_RECONSTRUCTION_HPP
23 
24 #include "fub/CompleteFromCons.hpp"
25 #include "fub/Execution.hpp"
26 #include "fub/State.hpp"
27 
28 #include "fub/AMReX/ForEachFab.hpp"
30 
31 #include <AMReX_MultiFab.H>
32 
33 namespace fub::amrex {
34 
35 template <typename Tag, typename Equation> class Reconstruction {
36 public:
37  Reconstruction(const Equation& eq) : Reconstruction(Tag(), eq) {}
38 
39  Reconstruction(Tag, const Equation& eq) : rec_{} {
42  }
43 
44  void CompleteFromCons(::amrex::MultiFab& dest, const ::amrex::MultiFab& src) {
45  ForEachFab(Tag(), dest, [&](const ::amrex::MFIter& mfi) {
46  Equation& eq = rec_->equation_;
47  View<Complete<Equation>> complete =
48  MakeView<Complete<Equation>>(dest[mfi], eq, mfi.growntilebox());
49  View<const Conservative<Equation>> conservative =
50  MakeView<const Conservative<Equation>>(src[mfi], eq,
51  mfi.growntilebox());
52  rec_->CompleteFromCons(Tag(), complete, conservative);
53  });
54  }
55 
56  void CompleteFromCons(IntegratorContext& context, int level, Duration) {
57  ::amrex::MultiFab& dest = context.GetScratch(level);
58  const ::amrex::MultiFab& src = context.GetScratch(level);
59  CompleteFromCons(dest, src);
60  }
61 
62 private:
64 };
65 
66 template <typename Equation>
68  ->Reconstruction<execution::OpenMpSimdTag, Equation>;
69 
70 template <typename Tag, typename Equation>
71 Reconstruction(Tag, const Equation& eq)->Reconstruction<Tag, Equation>;
72 
74  void CompleteFromCons([[maybe_unused]] IntegratorContext& context,
75  [[maybe_unused]] int level,
76  [[maybe_unused]] Duration time_step_size) {}
77 };
78 
79 } // namespace fub::amrex
80 
81 #endif
This class is used by the HypebrolicSplitLevelIntegrator and delegates AMR related tasks to the AMReX...
Definition: AMReX/IntegratorContext.hpp:49
::amrex::MultiFab & GetScratch(int level)
Returns the MultiFab associated with level data with ghost cells on the specifed level number and dir...
Definition: Reconstruction.hpp:35
Reconstruction(const Equation &eq)
Definition: Reconstruction.hpp:37
void CompleteFromCons(::amrex::MultiFab &dest, const ::amrex::MultiFab &src)
Definition: Reconstruction.hpp:44
Local< Tag, CompleteFromConsFn< Equation > > rec_
Definition: Reconstruction.hpp:63
Reconstruction(Tag, const Equation &eq)
Definition: Reconstruction.hpp:39
void CompleteFromCons(IntegratorContext &context, int level, Duration)
Definition: Reconstruction.hpp:56
void ForEachFab(Tag, const ::amrex::FabArrayBase &fabarray, F function)
Iterate through all local FArrayBox objects in a MultiFab.
Definition: ForEachFab.hpp:34
The amrex namespace.
Definition: AverageState.hpp:33
Reconstruction(const Equation &eq) -> Reconstruction< execution::OpenMpSimdTag, Equation >
std::decay_t< decltype(std::declval< T >().GetEquation())> Equation
A template typedef to detect the member function.
Definition: Meta.hpp:59
std::chrono::duration< double > Duration
Definition: Duration.hpp:31
typename detail::LocalType< Tag, T >::type Local
Definition: Execution.hpp:56
Definition: State.hpp:403
Definition: CompleteFromCons.hpp:141
Definition: Reconstruction.hpp:73
void CompleteFromCons([[maybe_unused]] IntegratorContext &context, [[maybe_unused]] int level, [[maybe_unused]] Duration time_step_size)
Definition: Reconstruction.hpp:74
Definition: Execution.hpp:39