My Project
|
Public Member Functions | |
SimMeshDataAccessor (const GridView &gridView, Dune::PartitionSet< partitions > dunePartition) | |
Construct a SimMeshDataAccessor working on a specific GridView and specialize to a Dune::PartitionSet<>. | |
bool | polyhedralCellPresent () const |
Checks for cells that have polyhedral type within the current partition of cells. | |
void | countEntities () |
Count the vertices, cells and corners. | |
template<typename T > | |
long | writeGridPoints (T *x_inout, T *y_inout, T *z_inout, long max_size=0) const |
Write the positions of vertices - directly to the pointers given in parameters. | |
template<typename VectType > | |
long | writeGridPoints (VectType &x_inout, VectType &y_inout, VectType &z_inout) const |
Write the positions of vertices - directly to the pointers given in parameters. | |
template<typename T > | |
long | writeGridPoints_AOS (T *xyz_inout, long max_size=0) const |
Write the positions of vertices - directly to the pointers given in parameters as Array of Structures x,y,z,x,y,z,x,y,z,... | |
template<typename VectType > | |
long | writeGridPoints_AOS (VectType &xyz_inout) const |
Write the positions of vertices - directly to the pointers given in parameters as Array of Structures x,y,z,x,y,z,x,y,z,... | |
template<typename T > | |
long | writeGridPoints_SOA (T *xyz_inout, long max_size=0) const |
Write the positions of vertices - directly to the pointers given in parameters as Structure of Arrays: x,x,x,...,y,y,y,...,z,z,z,... | |
template<typename VectType > | |
long | writeGridPoints_SOA (VectType &xyz_inout) const |
Write the positions of vertices - directly to the pointers given in parameters as Structure of Arrays: x,x,x,...,y,y,y,...,z,z,z,... | |
template<typename Integer > | |
long | writeConnectivity (Integer *connectivity_inout, ConnectivityVertexOrder whichOrder, long max_size=0) const |
Write the connectivity array - directly to the pointer given in parameter 1 Reorders the indices as selected either in DUNE order or VTK order. | |
template<typename VectType > | |
long | writeConnectivity (VectType &connectivity_inout, ConnectivityVertexOrder whichOrder) const |
Write the connectivity array - directly to a VectType object given in parameter 1 Reorders the indices as selected either in DUNE order or VTK order. | |
template<typename Integer > | |
long | writeOffsetsCells (Integer *offsets_inout, long max_size=0) const |
Write the offsets values - directly to the pointer given in parameter 1. | |
template<typename VectType > | |
long | writeOffsetsCells (VectType &offsets_inout) const |
Write the offsets values - directly to a VectType object given in parameter 1. | |
template<typename Integer > | |
long | writeCellTypes (Integer *types_inout, long max_size=0) const |
Write the cell types values - directly to the pointer given in parameter 1. | |
template<typename VectType > | |
long | writeCellTypes (VectType &types_inout) const |
Write the cell types values - directly to the VectType object given in parameter 1. | |
std::string | getPartitionTypeString () const |
Dune::PartitionSet< partitions > | getPartition (void) |
void | printGridDetails (std::ostream &outstr) const |
int | getNCells () const |
int | getNVertices () const |
int | getNCorners () const |
std::string | getError () const |
void | clearError () |
bool | hasError () const |
|
explicit |
Construct a SimMeshDataAccessor working on a specific GridView and specialize to a Dune::PartitionSet<>.
gridView | The gridView |
PartitionSet<> | the set of cells from which to extract geometric data |
The PartitionSet of the data can be specified from one of: Dune::Partitions::all Dune::Partitions::interior Dune::Partitions::border Dune::Partitions::overlap Dune::Partitions::front Dune::Partitions::ghost Dune::Partitions::interiorBorder Dune::Partitions::interiorBorderOverlap Dune::Partitions::interiorBorderOverlapFront Dune::Partitions::all
N.B. To visualise 'field' data on the extracted grid mesh then the field variable should contain at least as many vlaues as the mesh has cells (ncells_) or vertices (nvertices_) depending on if data is cell centred or vertex centred, respectively.
As we are templated on the Dune::PartitionSet<partitions>, values for ncorners_, nvertices_ and ncells_ cannot change
This class does not work with grids containing polyhedral cells (well, it has not been tested with this kind of grid data). The user should call polyhedralCellPresent() to test if polyhedral cells are present and decide what they want to do before copying data using the data accessor methods.
void Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::countEntities | ( | ) |
Count the vertices, cells and corners.
Count all the vertices ( the Dune::Partitions::all partition ) as then we do not need to renumber the vertices as all the subsets use references to the full set.
bool Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::polyhedralCellPresent | ( | ) | const |
Checks for cells that have polyhedral type within the current partition of cells.
Returns true if a polyhedral sell is found. If this is the case then this partition is not going to be available for visualisation as this class does not yet handle polyhedral cells.
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeCellTypes | ( | Integer * | types_inout, |
long | max_size = 0 ) const |
Write the cell types values - directly to the pointer given in parameter 1.
types_inout | is the array to be filled with the cell types (VTK defined values) |
max_size | is used to check that the space available in the input pointer parameter will fit the number of cell offset values written. |
Returns number of cells type values written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeCellTypes | ( | VectType & | types_inout | ) | const |
Write the cell types values - directly to the VectType object given in parameter 1.
types_inout | is the array to be filled with the cell types (VTK defined values) The object VectType must have a size() and data() method (e.g. a std::vector<T>) |
Returns number of cells type values written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeConnectivity | ( | Integer * | connectivity_inout, |
ConnectivityVertexOrder | whichOrder, | ||
long | max_size = 0 ) const |
Write the connectivity array - directly to the pointer given in parameter 1 Reorders the indices as selected either in DUNE order or VTK order.
connectivity_inout | is the array to be filled with connectivity indexes (i.e. the index into the vertex array) |
whichOrder,is | the order that verticies are traversed to create a cell (VTK or DUNE) |
max_size | is used to check that the space available in the input pointer parameter will fit the number of corner values written. Returns the number of corner indices written. |
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeConnectivity | ( | VectType & | connectivity_inout, |
ConnectivityVertexOrder | whichOrder ) const |
Write the connectivity array - directly to a VectType object given in parameter 1 Reorders the indices as selected either in DUNE order or VTK order.
connectivity_inout | is the array to be filled with connectivity indexes (i.e. the index into the vertex array) The object VectType must have a size() and data() method (e.g. a std::vector<T>) |
whichOrder,is | the order that verticies are traversed to create a cell (VTK or DUNE) |
max_size | is used to check that the space available in the input pointer parameter will fit the number of corner values written. Returns the number of corner indices written. |
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeGridPoints | ( | T * | x_inout, |
T * | y_inout, | ||
T * | z_inout, | ||
long | max_size = 0 ) const |
Write the positions of vertices - directly to the pointers given in parameters.
x_inout | to be filled with x coordinate verticies |
y_inout | to be filled with y coordinate verticies |
y_inout | to be filled with z coordinate verticies |
max_size | the maximum number of elements of type T that can be written to the input pointer memory regions. |
Returns the number of vertices written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeGridPoints | ( | VectType & | x_inout, |
VectType & | y_inout, | ||
VectType & | z_inout ) const |
Write the positions of vertices - directly to the pointers given in parameters.
x_inout | to be filled with x coordinate verticies |
y_inout | to be filled with y coordinate verticies |
y_inout | to be filled with z coordinate verticies |
All parameters must have a size() and data() method (e.g. a std::vector<T>) and the current size() must be big enough
Returns the number of vertices written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeGridPoints_AOS | ( | T * | xyz_inout, |
long | max_size = 0 ) const |
Write the positions of vertices - directly to the pointers given in parameters as Array of Structures x,y,z,x,y,z,x,y,z,...
xyz_inout | is the array to be filled with x,y,z coordinate verticies. |
max_size | is the maximum number x,y,z structures with elements of type T that can be written to the input pointer memory regions. |
Returns the number of vertices written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeGridPoints_AOS | ( | VectType & | xyz_inout | ) | const |
Write the positions of vertices - directly to the pointers given in parameters as Array of Structures x,y,z,x,y,z,x,y,z,...
xyz_inout | is the array to be filled with x,y,z coordinate verticies. The object VectType must have a size() and data() method (e.g. a std::vector<T>) |
Returns the number of vertices written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeGridPoints_SOA | ( | T * | xyz_inout, |
long | max_size = 0 ) const |
Write the positions of vertices - directly to the pointers given in parameters as Structure of Arrays: x,x,x,...,y,y,y,...,z,z,z,...
xyz_inout | is the array to be filled with x,y,z coordinate verticies. |
max_size | number of verticies (x,...y,...z,... structures) with elements of type T that can be written to the input pointer memory regions. |
Returns the number of vertices written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeGridPoints_SOA | ( | VectType & | xyz_inout | ) | const |
Write the positions of vertices - directly to the pointers given in parameters as Structure of Arrays: x,x,x,...,y,y,y,...,z,z,z,...
xyz_inout | is the array to be filled with x,y,z coordinate verticies. The object VectType must have a size() and data() method (e.g. a std::vector<T>) |
Returns the number of vertices written
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeOffsetsCells | ( | Integer * | offsets_inout, |
long | max_size = 0 ) const |
Write the offsets values - directly to the pointer given in parameter 1.
offsets_inout | is the array to be filled with offsets into the connectivity array (i.e. the index into the connectivity array to determine the vertices used for the particular cell) |
max_size | is used to check that the space available in the input pointer parameter will fit the number of cell offset values written. |
Returns number of offset values written + 1
long Opm::GridDataOutput::SimMeshDataAccessor< GridView, partitions >::writeOffsetsCells | ( | VectType & | offsets_inout | ) | const |
Write the offsets values - directly to a VectType object given in parameter 1.
offsets_inout | is the array to be filled with offsets into the connectivity array (i.e. the index into the connectivity array to determine the vertices used for the particular cell). The object VectType must have a size() and data() method (e.g. a std::vector<T>) |
Returns number of offset values written + 1