condition_action_sequence_2.hpp
1 #ifndef CONDITION_ACTION_SEQUENCE2_HPP
2 #define CONDITION_ACTION_SEQUENCE2_HPP 1
3 
4 #include "condition_action_sequence.hpp"
6 #include "../common/LagrangianHLLC.hpp"
7 
10 {
11 public:
12 
14  class Action2
15  {
16  public:
17 
31  virtual void operator()
32  (const Edge& edge,
33  const size_t index,
34  const Tessellation& tess,
35  const Vector2D& edge_velocity,
36  const vector<ComputationalCell>& cells,
37  const EquationOfState& eos,
38  const bool aux,
39  const pair<ComputationalCell, ComputationalCell> & edge_values,
40  Extensive &res,
41  double time,
42  TracerStickerNames const& tracerstickernames) const = 0;
43 
46  virtual void Reset(void) const {}
47 
48  virtual ~Action2(void);
49  private:
50  // Action2(const Action2&);
51  };
52 
59  (const vector<pair<const ConditionActionSequence::Condition*, const ConditionActionSequence::Action*> >& sequence,
60  const vector<pair<const ConditionActionSequence::Condition*, const ConditionActionSequence2::Action2*> >& sequence2,
61  SpatialReconstruction const& interp);
62 
64 
65  vector<Extensive> operator()
66  (const Tessellation& tess,
67  const vector<Vector2D>& edge_velocities,
68  const vector<ComputationalCell>& cells,
69  const vector<Extensive>& extensives,
70  const CacheData& cd,
71  const EquationOfState& eos,
72  const double time,
73  const double dt,
74  TracerStickerNames const& tracerstickernames) const;
75 
76 private:
77  const vector<pair<const ConditionActionSequence::Condition*, const ConditionActionSequence::Action*> > sequence_;
78  const vector<pair<const ConditionActionSequence::Condition*, const Action2*> > sequence2_;
79  const SpatialReconstruction & interp_;
80  mutable vector<pair<ComputationalCell, ComputationalCell> > edge_values_;
81 };
82 
85 {
86 public:
87 
91  explicit RegularFlux2(const RiemannSolver& rs);
92 
93  void operator()
94  (const Edge& edge,
95  const size_t index,
96  const Tessellation& tess,
97  const Vector2D& edge_velocity,
98  const vector<ComputationalCell>& cells,
99  const EquationOfState& eos,
100  const bool aux,
101  const pair<ComputationalCell, ComputationalCell> & edge_values,
102  Extensive &res, double time,
103  TracerStickerNames const& tracerstickernames) const;
104 
105 private:
106 
107  const RiemannSolver& rs_;
108 };
109 
110 
113 {
114 public:
115 
119  explicit RigidWallFlux2(const RiemannSolver& rs);
120 
121  void operator()
122  (const Edge& edge,
123  const size_t index,
124  const Tessellation& tess,
125  const Vector2D& edge_velocity,
126  const vector<ComputationalCell>& cells,
127  const EquationOfState& eos,
128  const bool aux,
129  const pair<ComputationalCell, ComputationalCell> & edge_values,
130  Extensive &res, double time,
131  TracerStickerNames const& tracerstickernames) const;
132 
133 private:
134  const RiemannSolver& rs_;
135 };
136 
139 {
140 public:
141 
146  Ratchet(const RiemannSolver& rs, const bool in);
147 
148  void operator()
149  (const Edge& edge,
150  const size_t index,
151  const Tessellation& tess,
152  const Vector2D& edge_velocity,
153  const vector<ComputationalCell>& cells,
154  const EquationOfState& eos,
155  const bool aux,
156  const pair<ComputationalCell, ComputationalCell> & edge_values,
157  Extensive &res, double time,
158  TracerStickerNames const& tracerstickernames) const;
159 
160 private:
161  const bool in_;
162  const RigidWallFlux2 wall_;
163  const FreeFlowFlux free_;
164 };
165 
168 {
169 public:
170 
173  {
174  public:
188  virtual bool operator()(const Edge& edge,
189  const size_t index,
190  const Tessellation& tess,
191  const Vector2D& edge_velocity,
192  const vector<ComputationalCell>& cells,
193  const EquationOfState& eos,
194  const bool aux,
195  const pair<ComputationalCell, ComputationalCell> & edge_values,
196  double time,
197  TracerStickerNames const& tracerstickernames) const = 0;
198 
199  virtual ~LagrangianCriteria();
200  };
201 
207  LagrangianFlux(const LagrangianHLLC& rs,const LagrangianHLLC& rs2,LagrangianCriteria const& criteria);
208 
209  void operator()
210  (const Edge& edge,
211  const size_t index,
212  const Tessellation& tess,
213  const Vector2D& edge_velocity,
214  const vector<ComputationalCell>& cells,
215  const EquationOfState& eos,
216  const bool aux,
217  const pair<ComputationalCell, ComputationalCell> & edge_values,
218  Extensive &res, double time,
219  TracerStickerNames const& tracerstickernames) const;
220 
223  void Reset(void) const;
224 
226  mutable vector<double> ws_;
228  mutable vector<double> edge_vel_;
230  mutable vector<bool> Lag_calc_;
231 private:
232  const LagrangianHLLC& rs_;
233  const LagrangianHLLC& rs2_;
234  const LagrangianCriteria& criteria_;
235 };
236 
239 {
240 public:
241  WallsMassFlux();
242 
243  bool operator()(const Edge& edge,
244  const size_t index,
245  const Tessellation& tess,
246  const Vector2D& edge_velocity,
247  const vector<ComputationalCell>& cells,
248  const EquationOfState& eos,
249  const bool aux,
250  const pair<ComputationalCell, ComputationalCell> & edge_values,
251  double time,
252  TracerStickerNames const& tracerstickernames) const;
253 };
254 
255 #endif // CONDITION_ACTION_SEQUENCE2_HPP
Extensive variables.
Definition: extensive.hpp:18
LagrangianHLLC Riemann solver for an Eulerian grid.
Calculates flux assuming rigid wall boundary conditions.
Spatial reconstruction of the primitive functions.
virtual void Reset(void) const
Return instance to initial state.
Abstract class for tessellation.
Condition on when to apply mass transfer fix.
vector< bool > Lag_calc_
Was this edge calculated Lagrangialy.
virtual void operator()(const Edge &edge, const size_t index, const Tessellation &tess, const Vector2D &edge_velocity, const vector< ComputationalCell > &cells, const EquationOfState &eos, const bool aux, const pair< ComputationalCell, ComputationalCell > &edge_values, Extensive &res, double time, TracerStickerNames const &tracerstickernames) const =0
Calculates flux.
Interface between two cells.
Definition: Edge.hpp:13
Base class for Riemann solver.
Allows matter to flow in only one direction.
Base class for flux calculator.
Base class for equation of state.
Second order flux calculator based on a series of conditions and actions.
Container for cache data.
Definition: cache_data.hpp:14
Estimate flux assuming free flow boundary conditions.
Calculates flux between two regular bulk cells.
Class for keeping the names of the tracers and stickers.
A flux scheme that minimises mass transfer between cells.
Abstract class for interpolation of the hydrodynamic variables.
ConditionActionSequence2(const vector< pair< const ConditionActionSequence::Condition *, const ConditionActionSequence::Action *> > &sequence, const vector< pair< const ConditionActionSequence::Condition *, const ConditionActionSequence2::Action2 *> > &sequence2, SpatialReconstruction const &interp)
Class constructor.
2D Mathematical vector
Definition: geometry.hpp:15
vector< double > edge_vel_
Velocity of the edges.
Criteria for having mass flux at outer edges of domain.
vector< double > ws_
Velocity of the interfaces.