Finite Volume Solver  prototype
A framework to build finite volume solvers for the AG Klein at the Freie Universität Berlin.
MultiBlockIntegratorContext.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_MultiBlock_INTEGRATOR_CONTEXT_HPP
22 #define FUB_AMREX_MultiBlock_INTEGRATOR_CONTEXT_HPP
23 
26 
28 
29 #include <memory>
30 #include <vector>
31 
32 namespace fub::amrex {
33 
34 /// \ingroup IntegratorContext
36 public:
38  std::vector<IntegratorContext> tubes,
39  std::vector<cutcell::IntegratorContext> plena,
40  std::vector<BlockConnection> connectivity);
41 
43  FlameMasterReactor reactor, std::vector<IntegratorContext> tubes,
44  std::vector<cutcell::IntegratorContext> plena,
45  std::vector<BlockConnection> connectivity,
46  std::vector<std::shared_ptr<PressureValve>> valves);
47 
48  /// @{
49  /// \name Member Accessors
50  [[nodiscard]] span<IntegratorContext> Tubes() noexcept;
51  [[nodiscard]] span<const IntegratorContext> Tubes() const noexcept;
52 
53  [[nodiscard]] span<cutcell::IntegratorContext> Plena() noexcept;
54  [[nodiscard]] span<const cutcell::IntegratorContext> Plena() const noexcept;
55 
56  [[nodiscard]] const std::shared_ptr<MultiBlockGriddingAlgorithm>&
57  GetGriddingAlgorithm() const noexcept;
58 
59  /// \brief Returns the current time level for data at the specified refinement
60  /// level and direction.
61  [[nodiscard]] Duration GetTimePoint(int level = 0) const;
62 
63  /// \brief Returns the current number of cycles for data at the specified
64  /// refinement level and direction.
65  [[nodiscard]] std::ptrdiff_t GetCycles(int level = 0) const;
66  /// @}
67 
68  /// \brief Returns the current boundary condition for the specified level.
69  [[nodiscard]] MultiBlockBoundary& GetBoundaryCondition(int level);
70 
71  [[nodiscard]] MPI_Comm GetMpiCommunicator() const noexcept;
72 
73  /// \brief Returns a shared pointer to the counter registry.
74  [[nodiscard]] const std::shared_ptr<CounterRegistry>&
75  GetCounterRegistry() const noexcept;
76 
77  /// @{
78  /// \name Observers
79 
80  /// \brief Returns true if the data exists for the specified level number.
81  bool LevelExists(int level) const noexcept;
82 
83  /// \brief Returns the refinement ratio in the specified direction.
84  int GetRatioToCoarserLevel(int level, Direction dir) const noexcept;
85 
86  /// \brief Returns the refinement ratio for all directions.
87  ::amrex::IntVect GetRatioToCoarserLevel(int level) const noexcept;
88  /// @}
89 
90  /// @{
91  /// \name Modifiers
92 
93  void CopyDataToScratch(int level);
94  void CopyScratchToData(int level);
95 
96  /// \brief Replaces the underlying gridding algorithm with the specified one.
98  std::shared_ptr<MultiBlockGriddingAlgorithm> gridding);
99 
100  /// \brief Whenever the gridding algorithm changes the data hierarchy this
101  /// function will regrid all distributed helper variables managed by the
102  /// context.
103  ///
104  /// \param[in] level The level number of the coarsest level which changed its
105  /// shape. Regrid all levels finer than level.
106  void ResetHierarchyConfiguration(int level = 0);
107 
108  /// \brief Sets the cycle count for a specific level number and direction.
109  void SetCycles(std::ptrdiff_t cycle, int level);
110 
111  /// \brief Sets the time point for a specific level number and direction.
112  void SetTimePoint(Duration t, int level);
113  /// @}
114 
115  /// @{
116  /// \name Member functions relevant for the level integrator algorithm.
117 
119 
121 
122  /// \brief On each first subcycle this will regrid the data if neccessary.
123  int PreAdvanceLevel(int level_num, Duration dt, std::pair<int, int> subcycle);
124 
125  /// \brief Increases the internal time stamps and cycle counters for the
126  /// specified level number and direction.
128  std::pair<int, int> subcycle);
129 
130  /// \brief Applies the boundary condition for the scratch space on level
131  /// `level` in direcition `dir`.
132  ///
133  /// \param level The refinement level on which the boundary condition shall
134  /// be used.
135  void ApplyBoundaryCondition(int level, Direction dir);
136 
137  /// \brief Fills the ghost layer of the scratch data and interpolates in the
138  /// coarse fine layer.
139  void FillGhostLayerTwoLevels(int level, int coarse);
140 
141  /// \brief Fills the ghost layer of the scratch data and does nothing in the
142  /// coarse fine layer.
143  void FillGhostLayerSingleLevel(int level);
144 
145  /// \brief Returns a estimate for a stable time step size which can be taken
146  /// for specified level number in direction dir.
148 
149  /// \brief Fill the flux MultiFab with numeric fluxes based on current states
150  /// in scratch.
151  void ComputeNumericFluxes(int level, Duration dt, Direction dir);
152 
153  /// \brief Apply a conservative time update for each conservative variable on
154  /// the specified level number and direction.
155  void UpdateConservatively(int level, Duration dt, Direction dir);
156 
157  /// \brief Reconstruct complete state variables from conservative ones.
158  void CompleteFromCons(int level, Duration dt);
159 
160  /// \brief Accumulate fluxes on the coarse fine interfaces for a specified
161  /// fine level number.
162  void AccumulateCoarseFineFluxes(int level, double time_scale, Direction dir);
163 
164  /// \brief Replace the coarse fluxes by accumulated fine fluxes on the coarse
165  /// fine interfaces.
166  void ApplyFluxCorrection(int fine, int coarse, Duration dt);
167 
168  /// \brief Resets all accumulates fluxes to zero.
169  void ResetCoarseFineFluxes(int fine, int coarse);
170 
171  /// \brief Coarsen scratch data from a fine level number to a coarse level
172  /// number.
173  void CoarsenConservatively(int fine, int coarse);
174  ///@}
175 
176 private:
177  std::vector<IntegratorContext> tubes_;
178  std::vector<cutcell::IntegratorContext> plena_;
180 };
181 
182 } // namespace fub::amrex
183 
184 #endif
A class mimicking the IdealGasMix / Reactor / ReactorNet interface of Cantera, but with FlameMaster c...
Definition: FlameMasterReactor.hpp:159
This class is used by the HypebrolicSplitLevelIntegrator and delegates AMR related tasks to the AMReX...
Definition: AMReX/IntegratorContext.hpp:49
Definition: MultiBlockBoundary.hpp:64
Definition: MultiBlockGriddingAlgorithm.hpp:39
Definition: MultiBlockIntegratorContext.hpp:35
void ComputeNumericFluxes(int level, Duration dt, Direction dir)
Fill the flux MultiFab with numeric fluxes based on current states in scratch.
void ApplyFluxCorrection(int fine, int coarse, Duration dt)
Replace the coarse fluxes by accumulated fine fluxes on the coarse fine interfaces.
void ApplyBoundaryCondition(int level, Direction dir)
Applies the boundary condition for the scratch space on level level in direcition dir.
void FillGhostLayerSingleLevel(int level)
Fills the ghost layer of the scratch data and does nothing in the coarse fine layer.
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...
MPI_Comm GetMpiCommunicator() const noexcept
void PreAdvanceHierarchy()
On each first subcycle this will regrid the data if neccessary.
span< IntegratorContext > Tubes() noexcept
Returns the current time level for data at the specified refinement level and direction.
std::vector< cutcell::IntegratorContext > plena_
Definition: MultiBlockIntegratorContext.hpp:178
Duration GetTimePoint(int level=0) const
Returns the current time level for data at the specified refinement level and direction.
std::ptrdiff_t GetCycles(int level=0) const
Returns the current number of cycles for data at the specified refinement level and direction.
void CopyDataToScratch(int level)
Replaces the underlying gridding algorithm with the specified one.
void PostAdvanceHierarchy()
On each first subcycle this will regrid the data if neccessary.
bool LevelExists(int level) const noexcept
Returns true if the data exists for the specified level number.
void SetTimePoint(Duration t, int level)
Sets the time point for a specific level number and direction.
const std::shared_ptr< MultiBlockGriddingAlgorithm > & GetGriddingAlgorithm() const noexcept
Returns the current time level for data at the specified refinement level and direction.
void CopyScratchToData(int level)
Replaces the underlying gridding algorithm with the specified one.
MultiBlockBoundary & GetBoundaryCondition(int level)
Returns the current boundary condition for the specified level.
void SetCycles(std::ptrdiff_t cycle, int level)
Sets the cycle count for a specific level number and direction.
std::vector< IntegratorContext > tubes_
Definition: MultiBlockIntegratorContext.hpp:177
int GetRatioToCoarserLevel(int level, Direction dir) const noexcept
Returns the refinement ratio in the specified direction.
void AccumulateCoarseFineFluxes(int level, double time_scale, Direction dir)
Accumulate fluxes on the coarse fine interfaces for a specified fine level number.
void ResetHierarchyConfiguration(std::shared_ptr< MultiBlockGriddingAlgorithm > gridding)
Replaces the underlying gridding algorithm with the specified one.
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.
std::shared_ptr< MultiBlockGriddingAlgorithm > gridding_
Definition: MultiBlockIntegratorContext.hpp:179
const std::shared_ptr< CounterRegistry > & GetCounterRegistry() const noexcept
Returns a shared pointer to the counter registry.
void CoarsenConservatively(int fine, int coarse)
Coarsen scratch data from a fine level number to a coarse level number.
MultiBlockIntegratorContext(FlameMasterReactor reactor, std::vector< IntegratorContext > tubes, std::vector< cutcell::IntegratorContext > plena, std::vector< BlockConnection > connectivity, std::vector< std::shared_ptr< PressureValve >> valves)
void UpdateConservatively(int level, Duration dt, Direction dir)
Apply a conservative time update for each conservative variable on the specified level number and dir...
int PreAdvanceLevel(int level_num, Duration dt, std::pair< int, int > subcycle)
On each first subcycle this will regrid the data if neccessary.
void ResetCoarseFineFluxes(int fine, int coarse)
Resets all accumulates fluxes to zero.
void CompleteFromCons(int level, Duration dt)
Reconstruct complete state variables from conservative ones.
span< cutcell::IntegratorContext > Plena() noexcept
Returns the current time level for data at the specified refinement level and direction.
MultiBlockIntegratorContext(FlameMasterReactor reactor, std::vector< IntegratorContext > tubes, std::vector< cutcell::IntegratorContext > plena, std::vector< BlockConnection > connectivity)
void FillGhostLayerTwoLevels(int level, int coarse)
Fills the ghost layer of the scratch data and interpolates in the coarse fine layer.
A span is a view over a contiguous sequence of objects, the storage of which is owned by some other o...
Definition: span.hpp:81
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: TimeStepError.hpp:71