hydrodynamics_2d.hpp
Go to the documentation of this file.
1 
6 #ifndef HYDRODYNAMICS_2D_HPP
7 #define HYDRODYNAMICS_2D_HPP 1
8 
10 #include "../common/equation_of_state.hpp"
11 #include "point_motion.hpp"
13 #include "SourceTerm.hpp"
14 #include "../../misc/utils.hpp"
15 #include "../../misc/lazy_list.hpp"
16 #include "../common/hydrodynamics.hpp"
17 #include "../../misc/universal_error.hpp"
18 #include "../../tessellation/ConvexHull.hpp"
19 #include "../../tessellation/polygon_overlap_area.hpp"
20 #include <boost/scoped_ptr.hpp>
21 #include "../../tessellation/EdgeLengthCorrect.hpp"
22 #include "physical_geometry.hpp"
23 #include "time_step_function.hpp"
24 #include "cache_data.hpp"
25 #include "../common/riemann_solver.hpp"
26 
33 Primitive RotatePrimitive(const Vector2D& normaldir,
34  const Vector2D& paraldir,
35  const Primitive& p);
36 
44  const Vector2D& n,
45  const Vector2D& p);
46 
52 int get_other_index(const Edge& edge,
53  const int index);
54 
65 vector<Primitive> InitialiseCells
66 (SpatialDistribution const& density,
67  SpatialDistribution const& pressure,
68  SpatialDistribution const& xvelocity,
69  SpatialDistribution const& yvelocity,
70  EquationOfState const& eos,
71  Tessellation const& tess, bool CMvalue = true);
72 
77 vector<Conserved> CalcConservedIntensive
78 (vector<Primitive> const& cells);
79 
86 vector<Conserved> CalcConservedExtensive
87 (const vector<Conserved>& cons_int,
88  const Tessellation& tess,
89  const PhysicalGeometry& pg);
90 
97 double TimeStepForCell(Primitive const& cell,
98  double width, vector<Vector2D> const& face_velocites);
99 
108 vector<int> MoveMeshPoints(vector<Vector2D> const& pointvelocity,
109  double dt, Tessellation& tessellation, bool reorder,
110  vector<Vector2D> oldpoints = vector<Vector2D>());
120 vector<int> MoveMeshPoints(vector<Vector2D> const& pointvelocity,
121  double dt, Tessellation& tessellation, Tessellation const& vproc, bool reorder,
122  vector<Vector2D> oldpoints = vector<Vector2D>());
123 
130 vector<Conserved> calc_conserved_intensive
131 (const Tessellation& tess,
132  const vector<Conserved>& extensive,
133  const PhysicalGeometry& pg);
134 
140 void UpdateConservedIntensive(Tessellation const& tessellation,
141  vector<Conserved> const& conservedextensive,
142  vector<Conserved>& conservedintensive);
143 
153 Conserved FluxInBulk(Vector2D const& normaldir,
154  Vector2D const& paraldir,
155  Primitive const& left,
156  Primitive const& right,
157  Vector2D const& edge_velocity,
158  RiemannSolver const& rs);
159 
174 (const Tessellation& tess,
175  const PhysicalGeometry& pg,
176  const CacheData& cd,
177  const vector<ComputationalCell>& cells,
178  const vector<Extensive>& fluxes,
179  const vector<Vector2D>& point_velocities,
180  const SourceTerm& source,
181  double t,
182  double dt,
183  vector<Extensive>& extensives,
184  TracerStickerNames const& tracerstickernames);
185 
190 vector<Vector2D> get_all_mesh_points
191 (Tessellation const& tess);
192 
198 vector<Primitive> make_eos_consistent
199 (vector<Primitive> const& vp,
200  EquationOfState const& eos);
201 
208 vector<double> GetForceEnergy(Tessellation const& tess,
209  vector<double> const& g);
210 
218 vector<vector<double> > calc_extensive_tracer
219 (const vector<vector<double> > & intensive_tracer,
220  const Tessellation& tess,
221  const vector<Primitive>& cells,
222  const PhysicalGeometry& pg);
223 
231 (vector<vector<double> > const&tracer,
232  Tessellation const& tess, vector<Primitive> const& cells,
233  vector<vector<double> > &result);
234 
243 void GetPointToRemove(Tessellation const& tess, Vector2D const& point,
244  double R, vector<int> & PointToRemove, int Inner);
245 
256 void FixAdvection(vector<Conserved>& extensive,
257  vector<Conserved> const& intensive,
258  Tessellation const& tessold,
259  Tessellation const& tessnew,
260  vector<Vector2D> const& facevelocity,
261  double dt, vector<Vector2D> const& pointvelocity);
262 
270 double determine_time_step(double hydro_time_step,
271  double external_dt,
272  double current_time,
273  double end_time);
274 
275 #endif // HYDRODYNAMICS_2D_HPP
Set of conserved variables (extensive)
double TimeStepForCell(Primitive const &cell, double width, vector< Vector2D > const &face_velocites)
Calculates the time step for a cell.
double determine_time_step(double hydro_time_step, double external_dt, double current_time, double end_time)
Determines the time step.
Abstract class for tessellation.
void GetPointToRemove(Tessellation const &tess, Vector2D const &point, double R, vector< int > &PointToRemove, int Inner)
Makes a list of points to remove.
vector< Conserved > CalcConservedIntensive(vector< Primitive > const &cells)
Calculates the intensive conserved variables.
vector< vector< double > > calc_extensive_tracer(const vector< vector< double > > &intensive_tracer, const Tessellation &tess, const vector< Primitive > &cells, const PhysicalGeometry &pg)
Calculates extensive tracers.
Interface between two cells.
Definition: Edge.hpp:13
void FixAdvection(vector< Conserved > &extensive, vector< Conserved > const &intensive, Tessellation const &tessold, Tessellation const &tessnew, vector< Vector2D > const &facevelocity, double dt, vector< Vector2D > const &pointvelocity)
Applies a correction to the extensive variables due to the change in volume during time step...
vector< Primitive > InitialiseCells(SpatialDistribution const &density, SpatialDistribution const &pressure, SpatialDistribution const &xvelocity, SpatialDistribution const &yvelocity, EquationOfState const &eos, Tessellation const &tess, bool CMvalue=true)
Initialize computational cells.
vector< Conserved > calc_conserved_intensive(const Tessellation &tess, const vector< Conserved > &extensive, const PhysicalGeometry &pg)
Calculates the intensive conserved variables.
Spatial distribution for initial conditions.
Abstract class for external forces.
Definition: SourceTerm.hpp:17
vector< Primitive > make_eos_consistent(vector< Primitive > const &vp, EquationOfState const &eos)
Changes the energy and sound speed so they would satisfy the equation of state.
Base class for Riemann solver.
int get_other_index(const Edge &edge, const int index)
Given an edge and an index of one neighbor, returns the index of another neighbor.
vector< Vector2D > get_all_mesh_points(Tessellation const &tess)
Returns the position of all mesh generating points.
vector< double > GetForceEnergy(Tessellation const &tess, vector< double > const &g)
Returns the energies due to external potentials.
Base class for equation of state.
void ExternalForceContribution(const Tessellation &tess, const PhysicalGeometry &pg, const CacheData &cd, const vector< ComputationalCell > &cells, const vector< Extensive > &fluxes, const vector< Vector2D > &point_velocities, const SourceTerm &source, double t, double dt, vector< Extensive > &extensives, TracerStickerNames const &tracerstickernames)
Adds force contribution to the extensive conserved variables.
Abstract class for source terms.
Container for cache data.
Definition: cache_data.hpp:14
Abstract class for motion of the mesh generating points.
void reorder(vector< T > &v, vector< size_t > const &order)
Reorder a vector according to an index vector (obtained from the &#39;ordered&#39; function) ...
Class for keeping the names of the tracers and stickers.
vector< Conserved > CalcConservedExtensive(const vector< Conserved > &cons_int, const Tessellation &tess, const PhysicalGeometry &pg)
Calculates the extensive conserved variables.
Conserved RotateFluxBack(const Conserved &c, const Vector2D &n, const Vector2D &p)
Rotates flux from the edge frame back to the lab frame.
Conserved FluxInBulk(Vector2D const &normaldir, Vector2D const &paraldir, Primitive const &left, Primitive const &right, Vector2D const &edge_velocity, RiemannSolver const &rs)
Calculates the flux in the bulk of the fluid.
Physical geometry of the grid.
Primitive hydrodynamic variables.
Abstract class for interpolation of the hydrodynamic variables.
void UpdateConservedIntensive(Tessellation const &tessellation, vector< Conserved > const &conservedextensive, vector< Conserved > &conservedintensive)
Updates the intensive conserved variables.
Geometrical cache data for optimisation.
Spatial distribution for initial conditions.
2D Mathematical vector
Definition: geometry.hpp:15
vector< int > MoveMeshPoints(vector< Vector2D > const &pointvelocity, double dt, Tessellation &tessellation, bool reorder, vector< Vector2D > oldpoints=vector< Vector2D >())
Move mesh points.
Base class for physical geometry.
Abstract class for time step function.
void MakeTracerExtensive(vector< vector< double > > const &tracer, Tessellation const &tess, vector< Primitive > const &cells, vector< vector< double > > &result)
Calculates the extensive tracer.
Primitive RotatePrimitive(const Vector2D &normaldir, const Vector2D &paraldir, const Primitive &p)
Rotates primitive variables to align with edge.