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