Mesh
-
class Mesh
Public Functions
-
Mesh(const Config *settings, std::vector<Vector> nodes, std::vector<std::vector<unsigned>> cells, std::vector<std::pair<BOUNDARY_TYPE, std::vector<unsigned>>> boundaries)
Constructor of mesh. Needs nodes, cells, and boundary descriptions as specified above. See LoadSU2MeshFromFile in io.cpp for setup information.
-
const std::vector<Vector> &GetNodes() const
Returns all node coordinates.
- Returns:
dimension: numNodes x dim
-
const std::vector<Vector> &GetCellMidPoints() const
Returns the mid point coordinates of each cell.
- Returns:
dimension: numCells x dim
-
const std::vector<std::vector<unsigned>> &GetCells() const
Returns all node IDs that construct up each cell.
- Returns:
dimension: numCells x numNodes
-
const std::vector<double> &GetCellAreas() const
Returns the cell area of each cell.
- Returns:
dimension: numCells
-
const std::vector<std::vector<unsigned>> &GetNeighbours() const
Returns the neighbor cell IDs for every cell.
- Returns:
dimension: numCells x numNodes
-
const std::vector<std::vector<Vector>> &GetNormals() const
Returns the edge length scaled normal vectors of each cell.
- Returns:
dimension: numCells x numNodes x dim
-
const std::vector<BOUNDARY_TYPE> &GetBoundaryTypes() const
Returns the boundary enum for each cell. BOUNDARY_TYPE::NONE is the default.
- Returns:
dimension: numCells
-
void SetBoundaryType(int idx_cell, BOUNDARY_TYPE boundary_type)
Set a boundary type for a cell. Use with caution!
- Parameters:
idx_cell – cell index
boundary_type – boundary_type to change the cell to
- Returns:
void
-
const std::vector<std::pair<double, double>> GetBounds() const
Returns the minimal and maximal coordinates of all nodes for each dimension.
- Returns:
dimension: dim
-
const std::vector<std::vector<Vector>> GetInterfaceMidPoints() const
Returns the the coordinates of midpoints of all interfaces of all cells.
- Returns:
dimension: numCells<numNeighborsPerCell<nDim>>
-
double GetDistanceToOrigin(unsigned idx_cell) const
Returns distance of a specified cells center to the coordinate systems origin.
- Returns:
dimension: scalar
-
unsigned GetCellOfKoordinate(const double x, const double y) const
Returns index of cell containing the coordinate (x,y)
- Returns:
cell_idx: unsigned
-
std::vector<unsigned> GetCellsofBall(const double x, const double y, const double r) const
Returns index of cells contained in the ball around the coordinate (x,y) with radius r.
- Returns:
cell_idxs: std::vector<unsigned>
-
std::vector<unsigned> GetCellsofRectangle(const std::vector<std::vector<double>> &cornercoordinates) const
Returns index of cells contained in the rectangle with specified corner coordinates.
- Returns:
cell_idxs: std::vector<unsigned>
-
void ComputeSlopes(unsigned nq, VectorVector &psiDerX, VectorVector &psiDerY, const VectorVector &psi) const
ComputeSlopes calculates the slope in every cell into x and y direction using the divergence theorem.
- Parameters:
nq – is number of quadrature points
psiDerX – is slope in x direction (gets computed. Slope is stored here)
psiDerY – is slope in y direction (gets computed. Slope is stored here)
psi – is solution for which slope is computed
-
void ComputeLimiter(unsigned nSys, const VectorVector &solDx, const VectorVector &solDy, const VectorVector &sol, VectorVector &limiter) const
Use gauss theorem and limiters. For unstructured mesh *.
- Parameters:
nq – is number of quadrature points
psiDerX – is slope in x direction (gets computed. Slope is stored here)
psiDerY – is slope in y direction (gets computed. Slope is stored here)
psi – is solution for which slope is computed
-
void ComputeLimiter1D(unsigned nSys, const VectorVector &sol, VectorVector &limiter) const
Use gauss theorem and limiters. For unstructured mesh *.
- Parameters:
nq – is number of quadrature points
psiDerX – is slope in x direction (gets computed. Slope is stored here)
psiDerY – is slope in y direction (gets computed. Slope is stored here)
psi – is solution for which slope is computed
-
void ComputeSlopes1D(unsigned nq, VectorVector &psiDerX, const VectorVector &psi) const
ComputeSlopes calculates for 1D meshes using finite difference formula in x direction.
- Parameters:
nq – is number of quadrature points
psiDerX – is slope in x direction (gets computed. Slope is stored here)
psi – is solution for which slope is computed
Protected Functions
-
void ComputeCellAreas()
Computes only the areas of the mesh cells. Write to _cellAreas.
-
void ComputeCellMidpoints()
Compute only the midpoints of the cells. Write to _cellMidPoints.
-
void ComputeConnectivity()
Computes _cellNeighbors and _nodeNeighbors, i.e. neighborship relation in mesh.
-
Vector ComputeOutwardFacingNormal(const Vector &nodeA, const Vector &nodeB, const Vector &cellCenter)
Computes outward facing normal of two neighboring nodes nodeA and nodeB with common cellCellcenter. Normals are scaled with their respective edge length.
- Parameters:
nodeA – first node
nodeB – neighboring node to nodeA
cellCenter – Center of the cell that has nodeA and nodeB as nodes.
- Returns:
outward facing normal
-
void ComputeBounds()
Computes the spatial bounds of a 2D domain.
-
Vector ComputeCellInterfaceMidpoints(const Vector &nodeA, const Vector &nodeB)
compute the midpoint of the edge between nodeA and nodeB
Protected Attributes
-
const unsigned _dim
spatial dimension of the mesh, i.e. 1D,2D,3D
-
const unsigned _numCells
number of cells in the mesh
-
const unsigned _numNodes
number of nodes in the mesh (for node centered view)
-
const unsigned _numNodesPerCell
number of nodes per cell
-
const unsigned _numBoundaries
number of boundary cells in the mesh
-
const unsigned _ghostCellID
Id of the ghost cell. (we use only one ghost cell). equal to _numCells and therefore has the ID of the last cell + 1.
-
std::vector<Vector> _nodes
nodes coordinates in the mesh. dimension:_numNodes<_dim>
-
std::vector<std::vector<unsigned>> _cells
node indices for each cell. dimension:_numCells<_numNodesPerCell>
-
std::vector<std::pair<BOUNDARY_TYPE, std::vector<unsigned>>> _boundaries
boundary cells in the mesh. Pair defines boundary type of the boundary nodes of the cell. numBoundaries<(1,numBoundaryNodes)>
-
std::vector<double> _cellAreas
cell areas of the mesh. dimension: numCells
-
std::vector<Vector> _cellMidPoints
cell midpoints of the mesh. dimension: numCells<dim>
-
std::vector<std::vector<unsigned>> _cellNeighbors
neighbors of each cell. dimension: numCells<numNodesPerCell>
-
std::vector<std::vector<Vector>> _cellInterfaceMidPoints
coordinates of the interface midpoints of all cells of the mesh dimension: numCells<numNeighborsPerCell<nDim>>. interfaces of each cell are in same order as _cellNeighbors
-
std::vector<std::vector<Vector>> _cellNormals
outward facing normals of each side of each cell. dimension: numCells<numNodesPerCell<dim>>, all normals are facing away from the cell center, and scaled with the edge length
-
std::vector<BOUNDARY_TYPE> _cellBoundaryTypes
Tags each cell with its boundary type. None means no boundary. dimension: numCells.
-
blaze::CompressedMatrix<bool> _nodeNeighbors
neighborshood relationship of nodes for (par-)metis
Private Functions
-
bool IsPointInsideCell(unsigned idx_cell, double x, double y) const
Function to check if a point is inside a polygon (cell)
-
Mesh(const Config *settings, std::vector<Vector> nodes, std::vector<std::vector<unsigned>> cells, std::vector<std::pair<BOUNDARY_TYPE, std::vector<unsigned>>> boundaries)