Various manipulations of hydrodynamic variables. More...
#include "spatial_distribution2d.hpp"
#include "../common/equation_of_state.hpp"
#include "point_motion.hpp"
#include "spatial_reconstruction.hpp"
#include "SourceTerm.hpp"
#include "../../misc/utils.hpp"
#include "../../misc/lazy_list.hpp"
#include "../common/hydrodynamics.hpp"
#include "../../misc/universal_error.hpp"
#include "../../tessellation/ConvexHull.hpp"
#include "../../tessellation/polygon_overlap_area.hpp"
#include <boost/scoped_ptr.hpp>
#include "../../tessellation/EdgeLengthCorrect.hpp"
#include "physical_geometry.hpp"
#include "time_step_function.hpp"
#include "cache_data.hpp"
#include "../common/riemann_solver.hpp"
Go to the source code of this file.
Functions | |
Primitive | RotatePrimitive (const Vector2D &normaldir, const Vector2D ¶ldir, const Primitive &p) |
Rotates primitive variables to align with edge. More... | |
Conserved | RotateFluxBack (const Conserved &c, const Vector2D &n, const Vector2D &p) |
Rotates flux from the edge frame back to the lab frame. More... | |
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. More... | |
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. More... | |
vector< Conserved > | CalcConservedIntensive (vector< Primitive > const &cells) |
Calculates the intensive conserved variables. More... | |
vector< Conserved > | CalcConservedExtensive (const vector< Conserved > &cons_int, const Tessellation &tess, const PhysicalGeometry &pg) |
Calculates the extensive conserved variables. More... | |
double | TimeStepForCell (Primitive const &cell, double width, vector< Vector2D > const &face_velocites) |
Calculates the time step for a cell. More... | |
vector< int > | MoveMeshPoints (vector< Vector2D > const &pointvelocity, double dt, Tessellation &tessellation, bool reorder, vector< Vector2D > oldpoints=vector< Vector2D >()) |
Move mesh points. More... | |
vector< int > | MoveMeshPoints (vector< Vector2D > const &pointvelocity, double dt, Tessellation &tessellation, Tessellation const &vproc, bool reorder, vector< Vector2D > oldpoints=vector< Vector2D >()) |
Move mesh points. More... | |
vector< Conserved > | calc_conserved_intensive (const Tessellation &tess, const vector< Conserved > &extensive, const PhysicalGeometry &pg) |
Calculates the intensive conserved variables. More... | |
void | UpdateConservedIntensive (Tessellation const &tessellation, vector< Conserved > const &conservedextensive, vector< Conserved > &conservedintensive) |
Updates the intensive conserved variables. More... | |
Conserved | FluxInBulk (Vector2D const &normaldir, Vector2D const ¶ldir, Primitive const &left, Primitive const &right, Vector2D const &edge_velocity, RiemannSolver const &rs) |
Calculates the flux in the bulk of the fluid. More... | |
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. More... | |
vector< Vector2D > | get_all_mesh_points (Tessellation const &tess) |
Returns the position of all mesh generating points. More... | |
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. More... | |
vector< double > | GetForceEnergy (Tessellation const &tess, vector< double > const &g) |
Returns the energies due to external potentials. More... | |
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. More... | |
void | MakeTracerExtensive (vector< vector< double > > const &tracer, Tessellation const &tess, vector< Primitive > const &cells, vector< vector< double > > &result) |
Calculates the extensive tracer. More... | |
void | GetPointToRemove (Tessellation const &tess, Vector2D const &point, double R, vector< int > &PointToRemove, int Inner) |
Makes a list of points to remove. More... | |
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. More... | |
double | determine_time_step (double hydro_time_step, double external_dt, double current_time, double end_time) |
Determines the time step. More... | |
Various manipulations of hydrodynamic variables.
Definition in file hydrodynamics_2d.hpp.
vector<Conserved> calc_conserved_intensive | ( | const Tessellation & | tess, |
const vector< Conserved > & | extensive, | ||
const PhysicalGeometry & | pg | ||
) |
Calculates the intensive conserved variables.
tess | Tessellation |
extensive | Extensive conserved variables |
pg | Physical geometry |
Definition at line 296 of file hydrodynamics_2d.cpp.
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.
intensive_tracer | Intensive tracers |
tess | Tessellation |
cells | List of primitive variables |
pg | Physical geometry |
Definition at line 453 of file hydrodynamics_2d.cpp.
vector<Conserved> CalcConservedExtensive | ( | const vector< Conserved > & | cons_int, |
const Tessellation & | tess, | ||
const PhysicalGeometry & | pg | ||
) |
Calculates the extensive conserved variables.
cons_int | Conserved intensive variables |
tess | Tessellation |
pg | Physical geometry |
Definition at line 190 of file hydrodynamics_2d.cpp.
Calculates the intensive conserved variables.
cells | Hydrodynamical cells |
Definition at line 133 of file hydrodynamics_2d.cpp.
double determine_time_step | ( | double | hydro_time_step, |
double | external_dt, | ||
double | current_time, | ||
double | end_time | ||
) |
Determines the time step.
hydro_time_step | Time step derived from hydrodynamics |
external_dt | Time step suggested by user |
current_time | Current simulation time |
end_time | Termination time |
Definition at line 566 of file hydrodynamics_2d.cpp.
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.
tess | Tessellation |
pg | Physical geometry |
cd | Cache data |
cells | Computational cells |
fluxes | Fluxes |
point_velocities | Velocities of the mesh generating points |
source | Source term |
t | Time |
dt | Time step |
extensives | Extensive variables |
tracerstickernames | The names of the tracers and stickers |
Definition at line 349 of file hydrodynamics_2d.cpp.
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.
This method calculates the change in extensive by calculating the volume swept by an edge and multiplying it by the intensive variables of the respective cell.
extensive | Extensive variables |
intensive | Intensive variables |
tessold | Old tessellation |
tessnew | New tessellation |
facevelocity | Face velocity |
dt | Time step |
pointvelocity | Velocities of mesh generating points |
Definition at line 496 of file hydrodynamics_2d.cpp.
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.
normaldir | A unit vector normal to the interface |
paraldir | A unit vector parallel to the interface |
left | Primitive variables on the left side of the interface |
right | Primitive variables on the right side of the interface |
edge_velocity | Velocity of the interface |
rs | Riemann solver |
Definition at line 334 of file hydrodynamics_2d.cpp.
vector<Vector2D> get_all_mesh_points | ( | Tessellation const & | tess | ) |
Returns the position of all mesh generating points.
tess | Tessellation |
Definition at line 375 of file hydrodynamics_2d.cpp.
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.
edge | Voronoi edge |
index | Index of one neighbor |
Definition at line 7 of file hydrodynamics_2d.cpp.
vector<double> GetForceEnergy | ( | Tessellation const & | tess, |
vector< double > const & | g | ||
) |
Returns the energies due to external potentials.
tess | Tessellation |
g | TBA |
Definition at line 393 of file hydrodynamics_2d.cpp.
void GetPointToRemove | ( | Tessellation const & | tess, |
Vector2D const & | point, | ||
double | R, | ||
vector< int > & | PointToRemove, | ||
int | Inner | ||
) |
Makes a list of points to remove.
tess | Tessellation |
point | TBA |
R | Radius? |
PointToRemove | output |
Inner | TBA |
Definition at line 476 of file hydrodynamics_2d.cpp.
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.
density | Density distribution |
pressure | Pressure distribution |
xvelocity | Distribution of the x component of the velocity |
yvelocity | Distribution of the y component of the velocity |
eos | Equation of state |
tess | Tessellation |
CMvalue | Determines whether to evaluate the spatial distributions in the mesh generating points or the center of mass |
Definition at line 95 of file hydrodynamics_2d.cpp.
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.
vp | Primitive variables |
eos | Equation of state |
Definition at line 384 of file hydrodynamics_2d.cpp.
void MakeTracerExtensive | ( | vector< vector< double > > const & | tracer, |
Tessellation const & | tess, | ||
vector< Primitive > const & | cells, | ||
vector< vector< double > > & | result | ||
) |
Calculates the extensive tracer.
tracer | Intensive tracer |
tess | Tessellation |
cells | Fluid elements |
result | List of extensive tracers |
Definition at line 464 of file hydrodynamics_2d.cpp.
vector<int> MoveMeshPoints | ( | vector< Vector2D > const & | pointvelocity, |
double | dt, | ||
Tessellation & | tessellation, | ||
bool | reorder, | ||
vector< Vector2D > | oldpoints = vector< Vector2D >() |
||
) |
Move mesh points.
pointvelocity | Velocities of all mesh points |
dt | Time step |
tessellation | Tessellation |
oldpoints | Possible input for the location of the old points |
reorder | Should we reorder the points |
Definition at line 235 of file hydrodynamics_2d.cpp.
vector<int> MoveMeshPoints | ( | vector< Vector2D > const & | pointvelocity, |
double | dt, | ||
Tessellation & | tessellation, | ||
Tessellation const & | vproc, | ||
bool | reorder, | ||
vector< Vector2D > | oldpoints = vector< Vector2D >() |
||
) |
Move mesh points.
pointvelocity | Velocities of all mesh points |
dt | Time step |
tessellation | The local tessellation |
vproc | The processor tessellation |
oldpoints | Possible input for the location of the old points |
reorder | Should we reorder the points |
Definition at line 250 of file hydrodynamics_2d.cpp.
Rotates flux from the edge frame back to the lab frame.
c | Flux |
n | Normal direction |
p | Parallel direction |
Definition at line 324 of file hydrodynamics_2d.cpp.
Primitive RotatePrimitive | ( | const Vector2D & | normaldir, |
const Vector2D & | paraldir, | ||
const Primitive & | p | ||
) |
Rotates primitive variables to align with edge.
normaldir | Normal directions |
paraldir | Parallel direction |
p | Primitive variables |
Definition at line 314 of file hydrodynamics_2d.cpp.
double TimeStepForCell | ( | Primitive const & | cell, |
double | width, | ||
vector< Vector2D > const & | face_velocites | ||
) |
Calculates the time step for a cell.
cell | Computational cell |
width | Cell width |
face_velocites | Velocities of the edges of the cell |
Definition at line 197 of file hydrodynamics_2d.cpp.
void UpdateConservedIntensive | ( | Tessellation const & | tessellation, |
vector< Conserved > const & | conservedextensive, | ||
vector< Conserved > & | conservedintensive | ||
) |
Updates the intensive conserved variables.
tessellation | Tessellation |
conservedextensive | Extensive conserved variables |
conservedintensive | Intensive conserved variables |
Definition at line 303 of file hydrodynamics_2d.cpp.