Alsvinn  0.5.3
The fast FVM simulator with UQ support
Public Member Functions | List of all members
alsfvm::grid::Grid Class Reference

#include <alsfvm/grid/Grid.hpp>

Public Member Functions

 Grid (rvec3 origin, rvec3 top, ivec3 dimensions, const std::array< boundary::Type, 6 > &boundaryConditions=boundary::allPeriodic())
 
 Grid (rvec3 origin, rvec3 top, ivec3 dimensions, const std::array< boundary::Type, 6 > &boundaryConditions, const ivec3 &globalPosition, const ivec3 &globalSize)
 
 Grid (rvec3 origin, rvec3 top, ivec3 dimensions, const std::array< boundary::Type, 6 > &boundaryConditions, const ivec3 &globalPosition, const ivec3 &globalSize, const rvec3 &cellLengths)
 
 Grid (rvec3 origin, rvec3 top, ivec3 dimensions, const std::array< boundary::Type, 6 > &boundaryConditions, const ivec3 &globalPosition, const ivec3 &globalSize, const rvec3 &cellLengths, const std::vector< rvec3 > &cellMidpoints)
 
rvec3 getOrigin () const
 
rvec3 getTop () const
 
ivec3 getDimensions () const
 
size_t getActiveDimension () const
 
rvec3 getCellLengths () const
 
const std::vector< rvec3 > & getCellMidpoints () const
 
boundary::Type getBoundaryCondition (int side) const
 
std::array< boundary::Type, 6 > getBoundaryConditions () const
 
std::vector< boundary::TypegetActiveBoundaryConditions () const
 
ivec3 getGlobalPosition () const
 
ivec3 getGlobalSize () const
 Get the total size (in number of cells) of the larger grid. More...
 

Detailed Description

Holds the information about the grid

Note
We only support regular, uniform cartesian grid

Constructor & Destructor Documentation

◆ Grid() [1/4]

alsfvm::grid::Grid::Grid ( rvec3  origin,
rvec3  top,
ivec3  dimensions,
const std::array< boundary::Type, 6 > &  boundaryConditions = boundary::allPeriodic() 
)

Constructs the Grid

Parameters
originthe origin point of the grid (the smallest point in lexicographical order)
topthe top right corner of the grid (maximum point in lexicographical order)
dimensionsthe dimensions of the grid (in number of cells in each direction)
boundaryConditionsfor each side, list the boundary conditions.
Index Spatial side 1D Spatial side 2D Spatial side 3D
0 left left left
1 right right right
2 < not used > bottom bottom
3 < not used > top top
4 < not used > < not used > front
5 < not used > < not used > back

◆ Grid() [2/4]

alsfvm::grid::Grid::Grid ( rvec3  origin,
rvec3  top,
ivec3  dimensions,
const std::array< boundary::Type, 6 > &  boundaryConditions,
const ivec3 globalPosition,
const ivec3 globalSize 
)

Constructs the Grid

Parameters
originthe origin point of the grid (the smallest point in lexicographical order)
topthe top right corner of the grid (maximum point in lexicographical order)
dimensionsthe dimensions of the grid (in number of cells in each direction)
boundaryConditionsfor each side, list the boundary conditions.
Index Spatial side 1D Spatial side 2D Spatial side 3D
0 left left left
1 right right right
2 < not used > bottom bottom
3 < not used > top top
4 < not used > < not used > front
5 < not used > < not used > back
globalPositionthe global position of the current grid in the large grid (used for MPI)
globalSizethe total size of the grid

Constructs the Grid

Parameters
originthe origin point of the grid (the smallest point in lexicographical order)
topthe top right corner of the grid (maximum point in lexicographical order)
dimensionsthe dimensions of the grid (in number of cells in each direction)

◆ Grid() [3/4]

alsfvm::grid::Grid::Grid ( rvec3  origin,
rvec3  top,
ivec3  dimensions,
const std::array< boundary::Type, 6 > &  boundaryConditions,
const ivec3 globalPosition,
const ivec3 globalSize,
const rvec3 cellLengths 
)

Constructs the Grid

This is the "least dummy proof version", since it lets the user specify the cellLengths. This should only be used for domain decomposition in MPI for instance. Unless you know what you are doing, don't use this version.

Parameters
originthe origin point of the grid (the smallest point in lexicographical order)
topthe top right corner of the grid (maximum point in lexicographical order)
dimensionsthe dimensions of the grid (in number of cells in each direction)
boundaryConditionsfor each side, list the boundary conditions.
Index Spatial side 1D Spatial side 2D Spatial side 3D
0 left left left
1 right right right
2 < not used > bottom bottom
3 < not used > top top
4 < not used > < not used > front
5 < not used > < not used > back
globalPositionthe global position of the current grid in the large grid (used for MPI)
globalSizethe total size of the grid
cellLengthsthe cell lengths in each direction
Note
The user is responsible for cellLengths being compatible with the rest of the parameters.

◆ Grid() [4/4]

alsfvm::grid::Grid::Grid ( rvec3  origin,
rvec3  top,
ivec3  dimensions,
const std::array< boundary::Type, 6 > &  boundaryConditions,
const ivec3 globalPosition,
const ivec3 globalSize,
const rvec3 cellLengths,
const std::vector< rvec3 > &  cellMidpoints 
)

Constructs the Grid

This is the "least dummy proof version", since it lets the user specify the cellLengths. This should only be used for domain decomposition in MPI for instance. Unless you know what you are doing, don't use this version.

Parameters
originthe origin point of the grid (the smallest point in lexicographical order)
topthe top right corner of the grid (maximum point in lexicographical order)
dimensionsthe dimensions of the grid (in number of cells in each direction)
boundaryConditionsfor each side, list the boundary conditions.
Index Spatial side 1D Spatial side 2D Spatial side 3D
0 left left left
1 right right right
2 < not used > bottom bottom
3 < not used > top top
4 < not used > < not used > front
5 < not used > < not used > back
globalPositionthe global position of the current grid in the large grid (used for MPI)
globalSizethe total size of the grid
cellLengthsthe cell lengths in each direction
cellMidpointsare the cell midpoints with respect to a larger grid, and indexed according to globalPosition
Note
The user is responsible for cellLengths being compatible with the rest of the parameters.

Member Function Documentation

◆ getActiveBoundaryConditions()

std::vector< boundary::Type > alsfvm::grid::Grid::getActiveBoundaryConditions ( ) const

Gets the boundary conditions for each side

Index Spatial side 1D Spatial side 2D Spatial side 3D
0 left left left
1 right right right
2 < not used > bottom bottom
3 < not used > top top
4 < not used > < not used > front
5 < not used > < not used > back

Will have length 2 * getActiveDimensions

◆ getActiveDimension()

size_t alsfvm::grid::Grid::getActiveDimension ( ) const

Gets the number of active dimensions (eg. 1D, 2D, 3D)

Returns
the number of active space dimensions

◆ getBoundaryCondition()

boundary::Type alsfvm::grid::Grid::getBoundaryCondition ( int  side) const

Gets the boundary conditions for the given side

Index Spatial side 1D Spatial side 2D Spatial side 3D
0 left left left
1 right right right
2 < not used > bottom bottom
3 < not used > top top
4 < not used > < not used > front
5 < not used > < not used > back

◆ getBoundaryConditions()

std::array< boundary::Type, 6 > alsfvm::grid::Grid::getBoundaryConditions ( ) const

Gets the boundary conditions for each side

Index Spatial side 1D Spatial side 2D Spatial side 3D
0 left left left
1 right right right
2 < not used > bottom bottom
3 < not used > top top
4 < not used > < not used > front
5 < not used > < not used > back

◆ getCellLengths()

rvec3 alsfvm::grid::Grid::getCellLengths ( ) const

Gets the cell lengths in each direction

◆ getCellMidpoints()

const std::vector< rvec3 > & alsfvm::grid::Grid::getCellMidpoints ( ) const

Returns the midpoints of the grid.

Note
this vector is indexed by
size_t index = z*nx*ny + y*nx + x;

◆ getDimensions()

ivec3 alsfvm::grid::Grid::getDimensions ( ) const

Gets the dimensions

Returns
the dimensions (number of cells in each direction)

◆ getGlobalPosition()

ivec3 alsfvm::grid::Grid::getGlobalPosition ( ) const

Gets the global position index of the grid, that is, the grid is contained in a virtual larger grid. This is used for MPI parallelization to know which part of the grid we are working on

◆ getGlobalSize()

ivec3 alsfvm::grid::Grid::getGlobalSize ( ) const

Get the total size (in number of cells) of the larger grid.

◆ getOrigin()

rvec3 alsfvm::grid::Grid::getOrigin ( ) const

Gets the origin point

Returns
the origin point

◆ getTop()

rvec3 alsfvm::grid::Grid::getTop ( ) const

Gets the top point

Returns
the top point

The documentation for this class was generated from the following files: