Finite Volume Solver  prototype
A framework to build finite volume solvers for the AG Klein at the Freie Universität Berlin.
AMReX/IntegratorContext.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_INTEGRATOR_CONTEXT_HPP
22 #define FUB_AMREX_INTEGRATOR_CONTEXT_HPP
23 
25 #include "fub/HyperbolicMethod.hpp"
26 #include "fub/TimeStepError.hpp"
27 #include "fub/ext/outcome.hpp"
28 
29 #include <AMReX_FluxRegister.H>
30 #include <AMReX_MultiFab.H>
31 
32 #include <array>
33 #include <memory>
34 #include <vector>
35 
36 namespace fub::amrex {
37 
38 /// \defgroup IntegratorContext Integrator Contexts
39 /// \brief This group summarizes all classes and functions that are realted to integrator contexts.
40 
41 class IntegratorContext;
42 
44 
45 /// \ingroup IntegratorContext
46 ///
47 /// \brief This class is used by the HypebrolicSplitLevelIntegrator and
48 /// delegates AMR related tasks to the AMReX library.
50 public:
51  /// @{
52  /// \name Constructors, Assignment Operators and Desctructor
53 
54  /// \brief Constructs a context object from given a gridding algorithm and a
55  /// numerical method.
56  IntegratorContext(std::shared_ptr<GriddingAlgorithm> gridding,
57  HyperbolicMethod method);
58 
59  IntegratorContext(std::shared_ptr<GriddingAlgorithm> gridding,
60  HyperbolicMethod method, int cell_gcw, int face_gcw);
61 
62  /// \brief Deeply copies a context and all its distributed data for all MPI
63  /// ranks.
65 
66  /// \brief Deeply copies a context and all its distributed data for all MPI
67  /// ranks.
69 
70  IntegratorContext(IntegratorContext&&) noexcept = default;
71 
72  IntegratorContext& operator=(IntegratorContext&&) noexcept = default;
73 
74  virtual ~IntegratorContext() = default;
75  /// @}
76 
77  /// @{
78  /// \name Member Accessors
79 
80  [[nodiscard]] int Rank() const noexcept;
81 
82  /// \brief Returns a shared pointer to the underlying GriddingAlgorithm which
83  /// owns the simulation.
84  [[nodiscard]] const std::shared_ptr<GriddingAlgorithm>&
85  GetGriddingAlgorithm() const noexcept;
86 
87  /// \brief Returns a reference to const PatchHierarchy which is a member of
88  /// the GriddingAlgorithm.
89  [[nodiscard]] const PatchHierarchy& GetPatchHierarchy() const noexcept;
90 
91  /// \brief Returns a reference to PatchHierarchy which is a member of the
92  /// GriddingAlgorithm.
93  [[nodiscard]] PatchHierarchy& GetPatchHierarchy() noexcept;
94 
95  /// \brief Returns the MPI communicator which is associated with this context.
96  [[nodiscard]] MPI_Comm GetMpiCommunicator() const noexcept;
97  /// @}
98 
99  /// @{
100  /// \name Access Level-specific data
101 
102  /// \brief Returns the MultiFab associated with level data on the specifed
103  /// level number.
104  [[nodiscard]] ::amrex::MultiFab& GetData(int level);
105  [[nodiscard]] const ::amrex::MultiFab& GetData(int level) const;
106 
107  /// \brief Returns the MultiFab associated with level data with ghost cells on
108  /// the specifed level number and direction.
109  [[nodiscard]] ::amrex::MultiFab& GetScratch(int level);
110  [[nodiscard]] const ::amrex::MultiFab& GetScratch(int level) const;
111 
112  /// \brief Returns the MultiFab associated with flux data on the specifed
113  /// level number and direction.
114  [[nodiscard]] ::amrex::MultiFab& GetFluxes(int level, Direction dir);
115  [[nodiscard]] const ::amrex::MultiFab& GetFluxes(int level,
116  Direction dir) const;
117 
118  /// \brief Returns the current time level for data at the specified refinement
119  /// level and direction.
120  [[nodiscard]] Duration GetTimePoint(int level = 0) const;
121 
122  /// \brief Returns the current number of cycles for data at the specified
123  /// refinement level and direction.
124  [[nodiscard]] std::ptrdiff_t GetCycles(int level = 0) const;
125 
126  /// \brief Returns the geometry object for the specified refinement level.
127  [[nodiscard]] const ::amrex::Geometry& GetGeometry(int level) const;
128 
129  /// \brief Returns a shared pointer to the counter registry.
130  [[nodiscard]] const std::shared_ptr<CounterRegistry>&
131  GetCounterRegistry() const noexcept;
132 
133  /// @}
134 
135  /// @{
136  /// \name Observers
137 
138  /// \brief Returns true if the data exists for the specified level number.
139  [[nodiscard]] bool LevelExists(int level) const noexcept;
140 
141  /// \brief Returns the refinement ratio in the specified direction.
142  [[nodiscard]] int GetRatioToCoarserLevel(int level, Direction dir) const
143  noexcept;
144 
145  /// \brief Returns the refinement ratio for all directions.
146  [[nodiscard]] ::amrex::IntVect GetRatioToCoarserLevel(int level) const
147  noexcept;
148  /// @}
149 
150  /// @{
151  /// \name Modifiers
152 
153  /// \brief Replaces the underlying gridding algorithm with the specified one.
154  virtual void
155  ResetHierarchyConfiguration(std::shared_ptr<GriddingAlgorithm> gridding);
156 
157  /// \brief Whenever the gridding algorithm changes the data hierarchy this
158  /// function will regrid all distributed helper variables managed by the
159  /// context.
160  ///
161  /// \param[in] level The level number of the coarsest level which changed its
162  /// shape. Regrid all levels finer than level.
163  virtual void ResetHierarchyConfiguration(int level = 0);
164 
165  /// \brief Sets the cycle count for a specific level number and direction.
166  void SetCycles(std::ptrdiff_t cycle, int level);
167 
168  /// \brief Sets the time point for a specific level number and direction.
169  void SetTimePoint(Duration t, int level);
170  /// @}
171 
172  /// @{
173  /// \name Member functions relevant for the level integrator algorithm.
174 
175  /// \brief Updates time point and cycle counter for the patch hierarchy.
177 
178  /// \brief On each first subcycle this will regrid the data if neccessary.
179  ///
180  /// \return Returns the coarsest level that changed. This function returns
181  /// the maximum number of levels if no change happened.
182  int PreAdvanceLevel(int level_num, Duration dt, std::pair<int, int> subcycle);
183 
184  /// \brief Increases the internal time stamps and cycle counters for the
185  /// specified level number and direction.
187  std::pair<int, int> subcycle);
188 
189  void CopyDataToScratch(int level_num);
190  void CopyScratchToData(int level_num);
191 
192  /// \brief Applies the boundary condition for the scratch space on level
193  /// `level` in direcition `dir`.
194  ///
195  /// \param level The refinement level on which the boundary condition shall
196  /// be used.
197  void ApplyBoundaryCondition(int level, Direction dir);
198  void ApplyBoundaryCondition(int level, Direction dir,
200 
201  /// \brief Fills the ghost layer of the scratch data and interpolates in the
202  /// coarse fine layer.
203  void FillGhostLayerTwoLevels(int level, AnyBoundaryCondition& fbc, int coarse,
204  AnyBoundaryCondition& cbc);
205  void FillGhostLayerTwoLevels(int level, int coarse);
206 
207  /// \brief Fills the ghost layer of the scratch data and does nothing in the
208  /// coarse fine layer.
210  void FillGhostLayerSingleLevel(int level);
211 
212  /// \brief Returns a estimate for a stable time step size which can be taken
213  /// for specified level number in direction dir.
214  [[nodiscard]] Duration ComputeStableDt(int level, Direction dir);
215 
216  /// \brief Fill the flux MultiFab with numeric fluxes based on current states
217  /// in scratch.
218  void ComputeNumericFluxes(int level, Duration dt, Direction dir);
219 
220  /// \brief Apply a conservative time update for each conservative variable on
221  /// the specified level number and direction.
222  void UpdateConservatively(int level, Duration dt, Direction dir);
223 
224  /// \brief Reconstruct complete state variables from conservative ones.
225  void CompleteFromCons(int level, Duration dt);
226 
227  /// \brief Accumulate fluxes on the coarse fine interfaces for a specified
228  /// fine level number.
229  void AccumulateCoarseFineFluxes(int level, double time_scale, Direction dir);
230 
231  /// \brief Replace the coarse fluxes by accumulated fine fluxes on the coarse
232  /// fine interfaces.
233  void ApplyFluxCorrection(int fine, int coarse, Duration dt);
234 
235  /// \brief Resets all accumulates fluxes to zero.
236  void ResetCoarseFineFluxes(int fine, int coarse);
237 
238  /// \brief Coarsen scratch data from a fine level number to a coarse level
239  /// number.
240  void CoarsenConservatively(int fine, int coarse);
241  ///@}
242 
243  bool count_per_level = false;
244 
245 private:
246  /// \brief This class holds auxiliary data on each refinement level.
247  struct LevelData {
248  LevelData() = default;
249  LevelData(const LevelData& other) = delete;
250  LevelData& operator=(const LevelData& other) = delete;
251  LevelData(LevelData&&) noexcept = default;
252  LevelData& operator=(LevelData&&) noexcept;
253  ~LevelData() noexcept = default;
254 
255  /// Scratch space with ghost cell widths
256  ::amrex::MultiFab scratch;
257 
258  /// These arrays will store the fluxes for each patch level which is present
259  /// in the patch hierarchy. These will need to be rebuilt if the
260  /// PatchHierarchy changes.
261  std::array<::amrex::MultiFab, AMREX_SPACEDIM> fluxes;
262 
263  /// FluxRegister accumulate fluxes on coarse fine interfaces between
264  /// refinement level. These will need to be rebuilt whenever the hierarchy
265  /// changes.
266  ::amrex::FluxRegister coarse_fine;
267 
268  Duration time_point{};
269  Duration regrid_time_point{};
270  std::ptrdiff_t cycles{};
271  };
272 
275  std::shared_ptr<GriddingAlgorithm> gridding_;
276  std::vector<LevelData> data_;
278 };
279 
280 } // namespace fub::amrex
281 
282 #endif
This is a polymorphic value type that wraps any BoundaryCondition object.
Definition: AnyBoundaryCondition.hpp:55
Definition: AMReX/Geometry.hpp:30
This class modifies and initializes a PatchLevel in a PatchHierarchy.
Definition: AMReX/GriddingAlgorithm.hpp:60
This class is used by the HypebrolicSplitLevelIntegrator and delegates AMR related tasks to the AMReX...
Definition: AMReX/IntegratorContext.hpp:49
IntegratorContext(std::shared_ptr< GriddingAlgorithm > gridding, HyperbolicMethod method, int cell_gcw, int face_gcw)
Constructs a context object from given a gridding algorithm and a numerical method.
std::ptrdiff_t GetCycles(int level=0) const
Returns the current number of cycles for data at the specified refinement level and direction.
HyperbolicMethod method_
Definition: AMReX/IntegratorContext.hpp:277
MPI_Comm GetMpiCommunicator() const noexcept
Returns the MPI communicator which is associated with this context.
void CopyScratchToData(int level_num)
Updates time point and cycle counter for the patch hierarchy.
void SetTimePoint(Duration t, int level)
Sets the time point for a specific level number and direction.
void ApplyBoundaryCondition(int level, Direction dir)
Applies the boundary condition for the scratch space on level level in direcition dir.
const std::shared_ptr< GriddingAlgorithm > & GetGriddingAlgorithm() const noexcept
Returns a shared pointer to the underlying GriddingAlgorithm which owns the simulation.
::amrex::MultiFab & GetScratch(int level)
Returns the MultiFab associated with level data with ghost cells on the specifed level number and dir...
void PostAdvanceHierarchy()
Updates time point and cycle counter for the patch hierarchy.
void CompleteFromCons(int level, Duration dt)
Reconstruct complete state variables from conservative ones.
int cell_ghost_cell_width_
Definition: AMReX/IntegratorContext.hpp:273
const std::shared_ptr< CounterRegistry > & GetCounterRegistry() const noexcept
Returns a shared pointer to the counter registry.
int face_ghost_cell_width_
Definition: AMReX/IntegratorContext.hpp:274
::amrex::MultiFab & GetFluxes(int level, Direction dir)
Returns the MultiFab associated with flux data on the specifed level number and direction.
Result< void, TimeStepTooLarge > PostAdvanceLevel(int level_num, Duration dt, std::pair< int, int > subcycle)
Increases the internal time stamps and cycle counters for the specified level number and direction.
void SetCycles(std::ptrdiff_t cycle, int level)
Sets the cycle count for a specific level number and direction.
IntegratorContext(IntegratorContext &&) noexcept=default
Constructs a context object from given a gridding algorithm and a numerical method.
IntegratorContext & operator=(const IntegratorContext &)
Deeply copies a context and all its distributed data for all MPI ranks.
void CoarsenConservatively(int fine, int coarse)
Coarsen scratch data from a fine level number to a coarse level number.
void AccumulateCoarseFineFluxes(int level, double time_scale, Direction dir)
Accumulate fluxes on the coarse fine interfaces for a specified fine level number.
Duration GetTimePoint(int level=0) const
Returns the current time level for data at the specified refinement level and direction.
void CopyDataToScratch(int level_num)
Updates time point and cycle counter for the patch hierarchy.
int PreAdvanceLevel(int level_num, Duration dt, std::pair< int, int > subcycle)
On each first subcycle this will regrid the data if neccessary.
Duration ComputeStableDt(int level, Direction dir)
Returns a estimate for a stable time step size which can be taken for specified level number in direc...
std::vector< LevelData > data_
Definition: AMReX/IntegratorContext.hpp:276
void ApplyFluxCorrection(int fine, int coarse, Duration dt)
Replace the coarse fluxes by accumulated fine fluxes on the coarse fine interfaces.
const PatchHierarchy & GetPatchHierarchy() const noexcept
Returns a reference to const PatchHierarchy which is a member of the GriddingAlgorithm.
void ResetCoarseFineFluxes(int fine, int coarse)
Resets all accumulates fluxes to zero.
void UpdateConservatively(int level, Duration dt, Direction dir)
Apply a conservative time update for each conservative variable on the specified level number and dir...
::amrex::MultiFab & GetData(int level)
Returns the MultiFab associated with level data on the specifed level number.
IntegratorContext(std::shared_ptr< GriddingAlgorithm > gridding, HyperbolicMethod method)
Constructs a context object from given a gridding algorithm and a numerical method.
IntegratorContext(const IntegratorContext &)
Deeply copies a context and all its distributed data for all MPI ranks.
int Rank() const noexcept
Returns a shared pointer to the underlying GriddingAlgorithm which owns the simulation.
void FillGhostLayerSingleLevel(int level, AnyBoundaryCondition &bc)
Fills the ghost layer of the scratch data and does nothing in the coarse fine layer.
virtual void ResetHierarchyConfiguration(std::shared_ptr< GriddingAlgorithm > gridding)
Replaces the underlying gridding algorithm with the specified one.
bool LevelExists(int level) const noexcept
Returns true if the data exists for the specified level number.
const ::amrex::Geometry & GetGeometry(int level) const
Returns the geometry object for the specified refinement level.
bool count_per_level
Definition: AMReX/IntegratorContext.hpp:243
int GetRatioToCoarserLevel(int level, Direction dir) const noexcept
Returns the refinement ratio in the specified direction.
void FillGhostLayerTwoLevels(int level, AnyBoundaryCondition &fbc, int coarse, AnyBoundaryCondition &cbc)
Fills the ghost layer of the scratch data and interpolates in the coarse fine layer.
void ComputeNumericFluxes(int level, Duration dt, Direction dir)
Fill the flux MultiFab with numeric fluxes based on current states in scratch.
std::shared_ptr< GriddingAlgorithm > gridding_
Definition: AMReX/IntegratorContext.hpp:275
The PatchHierarchy holds simulation data on multiple refinement levels.
Definition: AMReX/PatchHierarchy.hpp:156
The amrex namespace.
Definition: AverageState.hpp:33
std::chrono::duration< double > Duration
Definition: Duration.hpp:31
Direction
This is a type safe type to denote a dimensional split direction.
Definition: Direction.hpp:30
boost::outcome_v2::result< T, E > Result
Definition: outcome.hpp:32
Definition: HyperbolicMethod.hpp:150
Definition: TimeStepError.hpp:71
This class holds auxiliary data on each refinement level.
Definition: AMReX/IntegratorContext.hpp:247
LevelData(const LevelData &other)=delete
LevelData & operator=(const LevelData &other)=delete
LevelData(LevelData &&) noexcept=default