OGS
MeshLib::Mesh Class Reference

Detailed Description

A basic mesh.

Definition at line 33 of file Mesh.h.

#include <Mesh.h>

Inheritance diagram for MeshLib::Mesh:
[legend]
Collaboration diagram for MeshLib::Mesh:
[legend]

Public Member Functions

 Mesh (std::string name, std::vector< Node * > nodes, std::vector< Element * > elements, bool const compute_element_neighbors=false, Properties const &properties=Properties())
 Mesh (const Mesh &mesh)
 Copy constructor.
 Mesh (Mesh &&mesh)
Meshoperator= (const Mesh &mesh)=delete
Meshoperator= (Mesh &&mesh)=delete
virtual ~Mesh ()
 Destructor.
void shallowClean ()
void addElement (Element *elem)
 Add an element to the mesh.
unsigned getDimension () const
 Returns the dimension of the mesh (determined by the maximum dimension over all elements).
const NodegetNode (std::size_t idx) const
 Get the node with the given index.
const ElementgetElement (std::size_t idx) const
 Get the element with the given index.
std::size_t getNumberOfElements () const
 Get the number of elements.
std::size_t getNumberOfNodes () const
 Get the number of nodes.
const std::string getName () const
 Get name of the mesh.
std::vector< Node * > const & getNodes () const
 Get the nodes-vector for the mesh.
std::vector< Element * > const & getElements () const
 Get the element-vector for the mesh.
void resetElementIDs ()
 Resets the IDs of all mesh-elements to their position in the element vector.
void resetNodeIDs ()
 Resets the IDs of all mesh-nodes to their position in the node vector.
void setName (const std::string &name)
 Changes the name of the mesh.
std::size_t getID () const
 Get id of the mesh.
std::size_t computeNumberOfBaseNodes () const
 Get the number of base nodes.
bool hasNonlinearElement () const
 Check if the mesh contains any nonlinear element.
std::vector< Element const * > const & getElementsConnectedToNode (std::size_t node_id) const
std::vector< Element const * > const & getElementsConnectedToNode (Node const &node) const
PropertiesgetProperties ()
Properties const & getProperties () const
bool isAxiallySymmetric () const
void setAxiallySymmetric (bool is_axial_symmetric)

Protected Member Functions

void calcEdgeLengthRange ()
 Set the minimum and maximum length over the edges of the mesh.
void setDimension ()
 Sets the dimension of the mesh.
void setElementNeighbors ()

Protected Attributes

std::size_t const _id
unsigned _mesh_dimension
std::pair< double, double > _node_distance
 The minimal and maximal distance of nodes within an element over all elements in the mesh.
std::string _name
std::vector< Node * > _nodes
std::vector< Element * > _elements
Properties _properties
std::vector< std::vector< Element const * > > _elements_connected_to_nodes
bool _is_axially_symmetric = false
bool const _compute_element_neighbors

Friends

class ApplicationUtils::NodeWiseMeshPartitioner
void removeMeshNodes (Mesh &mesh, const std::vector< std::size_t > &nodes)

Constructor & Destructor Documentation

◆ Mesh() [1/3]

MeshLib::Mesh::Mesh ( std::string name,
std::vector< Node * > nodes,
std::vector< Element * > elements,
bool const compute_element_neighbors = false,
Properties const & properties = Properties() )

Constructor using a mesh name and an array of nodes and elements

Parameters
nameMesh name.
nodesA vector of mesh nodes.
elementsAn array of mesh elements.
compute_element_neighborsswitch to compute element neighbors or not
propertiesMesh properties.

Definition at line 52 of file Mesh.cpp.

61 _node_distance(std::numeric_limits<double>::max(), 0),
62 _name(std::move(name)),
63 _nodes(std::move(nodes)),
64 _elements(std::move(elements)),
65 _properties(properties),
66 _compute_element_neighbors(compute_element_neighbors)
67{
68 this->resetNodeIDs();
69 this->resetElementIDs();
70 this->setDimension();
71
73
75 {
76 this->setElementNeighbors();
77 }
78}
static std::size_t global_mesh_counter
Mesh counter used to uniquely identify meshes by id.
Definition Mesh.cpp:29
Properties _properties
Definition Mesh.h:151
std::size_t const _id
Definition Mesh.h:144
bool const _compute_element_neighbors
Definition Mesh.h:156
std::vector< std::vector< Element const * > > _elements_connected_to_nodes
Definition Mesh.h:153
unsigned _mesh_dimension
Definition Mesh.h:145
std::string _name
Definition Mesh.h:148
void resetNodeIDs()
Resets the IDs of all mesh-nodes to their position in the node vector.
Definition Mesh.cpp:149
std::vector< Element * > _elements
Definition Mesh.h:150
void setDimension()
Sets the dimension of the mesh.
Definition Mesh.cpp:167
void resetElementIDs()
Resets the IDs of all mesh-elements to their position in the element vector.
Definition Mesh.cpp:158
void setElementNeighbors()
Definition Mesh.cpp:194
std::pair< double, double > _node_distance
The minimal and maximal distance of nodes within an element over all elements in the mesh.
Definition Mesh.h:147
std::vector< Node * > _nodes
Definition Mesh.h:149
std::vector< std::vector< Element const * > > findElementsConnectedToNodes(Mesh const &mesh)
Definition Mesh.cpp:35

References _compute_element_neighbors, _elements, _elements_connected_to_nodes, _id, _mesh_dimension, _name, _node_distance, _nodes, _properties, MeshLib::findElementsConnectedToNodes(), global_mesh_counter, resetElementIDs(), resetNodeIDs(), setDimension(), and setElementNeighbors().

Referenced by Mesh(), Mesh(), MeshLib::NodePartitionedMesh::NodePartitionedMesh(), MeshLib::NodePartitionedMesh::NodePartitionedMesh(), ApplicationUtils::NodeWiseMeshPartitioner, operator=(), operator=(), and removeMeshNodes.

◆ Mesh() [2/3]

MeshLib::Mesh::Mesh ( const Mesh & mesh)

Copy constructor.

Definition at line 80 of file Mesh.cpp.

82 _mesh_dimension(mesh.getDimension()),
83 _node_distance(mesh._node_distance.first, mesh._node_distance.second),
84 _name(mesh.getName()),
85 _nodes(mesh.getNumberOfNodes()),
86 _elements(mesh.getNumberOfElements()),
87 _properties(mesh._properties),
88 _compute_element_neighbors(mesh._compute_element_neighbors)
89{
90 const std::vector<Node*>& nodes(mesh.getNodes());
91 const std::size_t nNodes(nodes.size());
92 for (unsigned i = 0; i < nNodes; ++i)
93 {
94 _nodes[i] = new Node(*nodes[i]);
95 }
96
97 const std::vector<Element*>& elements(mesh.getElements());
98 const std::size_t nElements(elements.size());
99 for (unsigned i = 0; i < nElements; ++i)
100 {
101 _elements[i] = elements[i]->clone();
102 for (auto const& [j, node_id] :
103 elements[i]->nodes() | views::ids | ranges::views::enumerate)
104 {
105 _elements[i]->setNode(static_cast<unsigned>(j), _nodes[node_id]);
106 }
107 }
108
109 if (_mesh_dimension == 0)
110 {
111 this->setDimension();
112 }
113
116 {
117 this->setElementNeighbors();
118 }
119}
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:216

References Mesh(), _compute_element_neighbors, _elements, _elements_connected_to_nodes, _id, _mesh_dimension, _name, _node_distance, _nodes, _properties, MeshLib::findElementsConnectedToNodes(), getDimension(), getElements(), getName(), getNodes(), getNumberOfElements(), getNumberOfNodes(), global_mesh_counter, MeshLib::views::ids, MeshLib::Node, setDimension(), and setElementNeighbors().

◆ Mesh() [3/3]

MeshLib::Mesh::Mesh ( Mesh && mesh)
default

References Mesh().

◆ ~Mesh()

MeshLib::Mesh::~Mesh ( )
virtual

Destructor.

Definition at line 129 of file Mesh.cpp.

130{
131 const std::size_t nElements(_elements.size());
132 for (std::size_t i = 0; i < nElements; ++i)
133 {
134 delete _elements[i];
135 }
136
137 const std::size_t nNodes(_nodes.size());
138 for (std::size_t i = 0; i < nNodes; ++i)
139 {
140 delete _nodes[i];
141 }
142}

References _elements, and _nodes.

Member Function Documentation

◆ addElement()

void MeshLib::Mesh::addElement ( Element * elem)

Add an element to the mesh.

Definition at line 144 of file Mesh.cpp.

145{
146 _elements.push_back(elem);
147}

References _elements.

◆ calcEdgeLengthRange()

void MeshLib::Mesh::calcEdgeLengthRange ( )
protected

Set the minimum and maximum length over the edges of the mesh.

◆ computeNumberOfBaseNodes()

std::size_t MeshLib::Mesh::computeNumberOfBaseNodes ( ) const

Get the number of base nodes.

Definition at line 228 of file Mesh.cpp.

229{
230 return std::count_if(begin(_nodes), end(_nodes),
231 [this](auto const* const node) {
232 return isBaseNode(
233 *node,
234 _elements_connected_to_nodes[node->getID()]);
235 });
236}
bool isBaseNode(Node const &node, std::vector< Element const * > const &elements_connected_to_node)
Definition Mesh.cpp:336

References _elements_connected_to_nodes, _nodes, and MeshLib::isBaseNode().

Referenced by MeshLib::NodePartitionedMesh::NodePartitionedMesh(), MeshLib::computeNumberOfRegularBaseNodesAtRank(), MeshLib::computeNumberOfRegularHigherOrderNodesAtRank(), and ApplicationUtils::setNumberOfNodesInPartitions().

◆ getDimension()

unsigned MeshLib::Mesh::getDimension ( ) const
inline

Returns the dimension of the mesh (determined by the maximum dimension over all elements).

Definition at line 79 of file Mesh.h.

79{ return _mesh_dimension; }

References _mesh_dimension.

Referenced by Mesh(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::PhaseFieldProcess(), ProcessLib::LIE::PostProcessTool::PostProcessTool(), ProcessLib::SurfaceFlux::SurfaceFlux(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::ThermoMechanicsProcess(), MeshToolsLib::MeshGenerator::AddFaultToVoxelGrid::addFaultToVoxelGrid(), MeshToolsLib::addLayerToMesh(), MeshGeoToolsLib::GeoMapper::advancedMapOnMesh(), MeshToolsLib::RasterDataToMesh::checkMesh(), MeshView::contextMenuEvent(), MeshToolsLib::convertMeshToGeo(), ProcessLib::createBoundaryCondition(), MeshToolsLib::createBoundaryElements(), ProcessLib::ComponentTransport::createComponentTransportProcess(), ProcessLib::createConstraintDirichletBoundaryCondition(), ProcessLib::createDirichletBoundaryCondition(), ProcessLib::createDirichletBoundaryConditionWithinTimeInterval(), MeshToolsLib::createFlippedMesh(), ProcessLib::createHCNonAdvectiveFreeComponentFlowBoundaryCondition(), ProcessLib::HT::createHTProcess(), LayeredMeshGenerator::createLayers(), ProcessLib::LiquidFlow::createLiquidFlowProcess(), ProcessLib::createNeumannBoundaryCondition(), ProcessLib::createPrimaryVariableConstraintDirichletBoundaryCondition(), ProcessLib::createPythonBoundaryCondition(), ProcessLib::createPythonSourceTerm(), LayeredVolume::createRasterLayers(), MeshToolsLib::MeshLayerMapper::createRasterLayers(), ProcessLib::createReleaseNodalForce(), ProcessLib::RichardsComponentTransport::createRichardsComponentTransportProcess(), ProcessLib::createRobinBoundaryCondition(), ProcessLib::createSolutionDependentDirichletBoundaryCondition(), ProcessLib::createSourceTerm(), MeshToolsLib::MeshLayerMapper::createStaticLayers(), ProcessLib::ThermoRichardsFlow::createThermoRichardsFlowProcess(), ProcessLib::createTimeDecayDirichletBoundaryCondition(), ProcessLib::createVariableDependentNeumannBoundaryCondition(), ProcessLib::createWellboreCompensateNeumannBoundaryCondition(), ProcessLib::WellboreSimulator::createWellboreSimulatorProcess(), MeshToolsLib::MeshValidation::detectHoles(), ProcessLib::LIE::anonymous_namespace{MeshUtils.cpp}::findFracutreIntersections(), MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh(), ProcessLib::LIE::getFractureMatrixDataInMesh(), MeshToolsLib::MeshSurfaceExtraction::getMeshSurface(), MeshToolsLib::MeshSurfaceExtraction::getSurfaceAreaForNodes(), getSurfaceIntegratedValuesForNodes(), MeshToolsLib::MeshSurfaceExtraction::getSurfaceNodes(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::initializeConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::initializeConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::initializeConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::initializeConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::initializeConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::initializeConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::initializeConcreteProcess(), MeshLib::is2DMeshOnRotatedVerticalPlane(), MeshToolsLib::MeshLayerMapper::layerMapping(), MeshGeoToolsLib::GeoMapper::mapOnMesh(), MeshToolsLib::MeshLayerMapper::mapToStaticValue(), setDimension(), MeshToolsLib::Mesh2MeshPropertyInterpolation::setPropertiesForMesh(), testIfMeshesAreEqual(), FileIO::SHPInterface::write2dMeshToSHP(), and FileIO::TetGenInterface::writeTetGenSmesh().

◆ getElement()

◆ getElements()

std::vector< Element * > const & MeshLib::Mesh::getElements ( ) const
inline

Get the element-vector for the mesh.

Definition at line 100 of file Mesh.h.

100{ return _elements; }

References _elements.

Referenced by Mesh(), ProcessLib::LIE::PostProcessTool::PostProcessTool(), ProcessLib::SurfaceFlux::SurfaceFlux(), ProcessLib::VolumetricSourceTerm::VolumetricSourceTerm(), LayeredVolume::addLayerBoundaries(), LayeredVolume::addLayerToMesh(), MeshToolsLib::addLayerToMesh(), MeshToolsLib::MeshLayerMapper::addLayerToMesh(), MeshModel::addMeshObject(), MeshGeoToolsLib::appendLinesAlongPolylines(), MaterialPropertyLib::checkMaterialSpatialDistributionMap(), MaterialPropertyLib::checkMPLPhasesForSinglePhaseFlow(), ProcessLib::ComponentTransport::checkMPLProperties(), ProcessLib::RichardsComponentTransport::anonymous_namespace{CreateRichardsComponentTransportProcess.cpp}::checkMPLProperties(), computePointCloudNodes(), MeshToolsLib::convertMeshToGeo(), MeshToolsLib::convertToLinearMesh(), ProcessLib::createAnchorTerm(), MeshToolsLib::createBoundaryElements(), ProcessLib::ComponentTransport::createComponentTransportProcess(), ProcessLib::createDeactivatedSubdomainMesh(), ProcessLib::createEmbeddedAnchor(), MeshToolsLib::createFlippedMesh(), ProcessLib::HT::createHTProcess(), ProcessLib::LIE::HydroMechanics::createHydroMechanicsProcess(), MeshToolsLib::createIntegrationPointProperties(), ProcessLib::LiquidFlow::createLiquidFlowProcess(), NumLib::createNumericalStabilization(), createPropertyVector(), MeshToolsLib::createQuadraticOrderMesh(), ParameterLib::createRandomFieldMeshElementParameter(), MeshToolsLib::MeshLayerMapper::createRasterLayers(), MeshToolsLib::MeshLayerMapper::createStaticLayers(), ApplicationUtils::distributeElementsIntoPartitions(), extractBoundaries(), extractBoundaryMeshes(), MeshLib::findElementsConnectedToNodes(), ProcessLib::HeatTransportBHE::getBHEDataInMesh(), ProcessLib::LIE::getFractureMatrixDataInMesh(), MeshLib::IO::getLocalTopologyType(), MeshLib::getMeshElementsForMaterialIDs(), MeshToolsLib::MeshSurfaceExtraction::getMeshSurface(), MeshToolsLib::MeshInformation::getNumberOfElementTypes(), MeshToolsLib::MeshSurfaceExtraction::getSurfaceNodes(), anonymous_namespace{IdentifySubdomainMesh.cpp}::identifySubdomainMeshElements(), ProcessLib::ComponentTransport::ComponentTransportProcess::initializeConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::initializeConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::initializeConcreteProcess(), ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::HT::HTProcess::initializeConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::initializeConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::initializeConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::initializeConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::initializeConcreteProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::initializeConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::initializeConcreteProcess(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >::initializeConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::initializeConcreteProcess(), ProcessLib::WellboreSimulator::WellboreSimulatorProcess::initializeConcreteProcess(), MeshToolsLib::Mesh2MeshPropertyInterpolation::interpolatePropertiesForMesh(), MeshLib::is2DMeshOnRotatedVerticalPlane(), MeshToolsLib::MeshGenerator::AddFaultToVoxelGrid::isVoxelGrid(), main(), ApplicationUtils::markDuplicateGhostCells(), anonymous_namespace{AddFaultToVoxelGrid.cpp}::markFaults(), MeshToolsLib::mergeMeshToBulkMesh(), ApplicationUtils::partitionMesh(), MeshToolsLib::RasterDataToMesh::projectToElements(), MeshToolsLib::removeElements(), removeLineElements(), MeshToolsLib::removeNodes(), reorderNonlinearNodes(), MeshGeoToolsLib::resetMeshElementProperty(), MeshToolsLib::ElementValueModification::setByElementType(), setMaterialIDs(), ElementTreeModel::setMesh(), MeshLib::MeshElementGrid::sortElementsInGridCells(), MeshToolsLib::MeshValidation::testElementGeometry(), MeshLib::transformMeshToNodePartitionedMesh(), MeshLib::IO::transformToXDMFTopology(), FileIO::TetGenInterface::write2dElements(), FileIO::SHPInterface::write2dMeshToSHP(), FileIO::TetGenInterface::write3dElements(), and MeshToolsLib::zeroMeshFieldDataByMaterialIDs().

◆ getElementsConnectedToNode() [1/2]

std::vector< MeshLib::Element const * > const & MeshLib::Mesh::getElementsConnectedToNode ( Node const & node) const

Definition at line 252 of file Mesh.cpp.

254{
255 return _elements_connected_to_nodes[node.getID()];
256}

References _elements_connected_to_nodes, and MathLib::Point3dWithID::getID().

◆ getElementsConnectedToNode() [2/2]

◆ getID()

std::size_t MeshLib::Mesh::getID ( ) const
inline

◆ getName()

const std::string MeshLib::Mesh::getName ( ) const
inline

Get name of the mesh.

Definition at line 94 of file Mesh.h.

94{ return _name; }

References _name.

Referenced by ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::HydroMechanicsProcess(), Mesh(), ProcessLib::LIE::PostProcessTool::PostProcessTool(), ProcessLib::SolutionDependentDirichletBoundaryCondition::SolutionDependentDirichletBoundaryCondition(), ProcessLib::addBulkMeshPropertyToSubMesh(), addGlobalIDsToMesh(), MeshToolsLib::addLayerToMesh(), MeshModel::addMeshObject(), addPrimaryVariablesToMesh(), ProcessLib::Output::addProcess(), MeshGeoToolsLib::appendLinesAlongPolylines(), anonymous_namespace{SubmeshResiduumOutputConfig.cpp}::checkBulkIDMappingsPresent(), anonymous_namespace{SubmeshResiduumOutputConfig.cpp}::checkNonOverlappingCover(), ProcessLib::checkParametersOfDirichletBoundaryCondition(), MeshLib::computeGhostBaseNodeGlobalNodeIDsOfSubDomainPartition(), anonymous_namespace{SubmeshResiduumOutputConfig.cpp}::computeNonOverlappingBulkMeshCoverBySubmeshes(), MeshLib::computeNumberOfRegularNodes(), MeshLib::computeRegularBaseNodeGlobalNodeIDsOfSubDomainPartition(), anonymous_namespace{convertMeshToGeo.cpp}::convertMeshNodesToGeoPoints(), MeshToolsLib::convertMeshToGeo(), MeshToolsLib::convertSurfaceToMesh(), ProcessLib::createHCNonAdvectiveFreeComponentFlowBoundaryCondition(), ProcessLib::HydroMechanics::createHydroMechanicsProcess(), ProcessLib::createNeumannBoundaryCondition(), MeshToolsLib::createQuadraticOrderMesh(), ProcessLib::createReleaseNodalForce(), ProcessLib::createRobinBoundaryCondition(), ProcessLib::createSourceTerm(), ProcessLib::createVariableDependentNeumannBoundaryCondition(), ProcessLib::createWellboreCompensateNeumannBoundaryCondition(), ProcessLib::Output::doOutputAlways(), ProcessLib::Output::doOutputNonlinearIteration(), MeshView::exportToShapefile(), MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh(), MeshToolsLib::MeshSurfaceExtraction::getMeshSurface(), MeshLib::getNumberOfGlobalNodes(), NumLib::MeshComponentMap::getSubset(), MeshLib::is2DMeshOnRotatedVerticalPlane(), ParameterLib::isDefinedOnSameMesh(), main(), MeshView::openRasterDataToMeshDialog(), outputAABB(), ProcessLib::outputMeshVtk(), parseOutputMeshConfig(), ElementTreeModel::setMesh(), transferPropertiesFromPartitionedMeshToUnpartitionedMesh(), MeshLib::transformMeshToNodePartitionedMesh(), anonymous_namespace{IdentifySubdomainMesh.cpp}::updateOrCheckExistingSubdomainProperty(), MeshLib::vtkStandardNewMacro(), and FileIO::SHPInterface::write2dMeshToSHP().

◆ getNode()

◆ getNodes()

std::vector< Node * > const & MeshLib::Mesh::getNodes ( ) const
inline

Get the nodes-vector for the mesh.

Definition at line 97 of file Mesh.h.

97{ return _nodes; }

References _nodes.

Referenced by MeshGeoToolsLib::BoundaryElementsAtPoint::BoundaryElementsAtPoint(), ProcessLib::ConstraintDirichletBoundaryCondition::ConstraintDirichletBoundaryCondition(), Mesh(), ProcessLib::LIE::PostProcessTool::PostProcessTool(), ProcessLib::PrimaryVariableConstraintDirichletBoundaryCondition::PrimaryVariableConstraintDirichletBoundaryCondition(), ProcessLib::SolutionDependentDirichletBoundaryCondition::SolutionDependentDirichletBoundaryCondition(), ProcessLib::SurfaceFlux::SurfaceFlux(), LayeredVolume::addLayerToMesh(), MeshToolsLib::addLayerToMesh(), MeshToolsLib::MeshLayerMapper::addLayerToMesh(), adjustExtent(), MeshGeoToolsLib::appendLinesAlongPolylines(), MeshLib::calculateNodesConnectedByElements(), MeshLib::computeGhostBaseNodeGlobalNodeIDsOfSubDomainPartition(), MeshLib::computeNumberOfRegularNodes(), MeshLib::computeRegularBaseNodeGlobalNodeIDsOfSubDomainPartition(), anonymous_namespace{convertMeshToGeo.cpp}::convertMeshNodesToGeoPoints(), MeshToolsLib::convertToLinearMesh(), ProcessLib::createBoundaryDOFTable(), ProcessLib::ComponentTransport::createComponentTransportProcess(), MeshToolsLib::createFlippedMesh(), ProcessLib::HeatConduction::createHeatConductionProcess(), ProcessLib::HT::createHTProcess(), ProcessLib::LiquidFlow::createLiquidFlowProcess(), createMeshFromElements(), MeshToolsLib::createQuadraticOrderMesh(), ProcessLib::createReleaseNodalForce(), ProcessLib::createSourceTerm(), MeshToolsLib::MeshLayerMapper::createStaticLayers(), MeshLib::NodeAdjacencyTable::createTable(), ProcessLib::createVariableDependentNeumannBoundaryCondition(), ProcessLib::createWellboreCompensateNeumannBoundaryCondition(), ApplicationUtils::determineAndAppendGhostNodesToPartitions(), ProcessLib::extractInnerAndOuterNodes(), MeshLib::findElementsConnectedToNodes(), MeshToolsLib::MeshInformation::getBoundingBox(), ProcessLib::getEssentialBCValuesLocal(), getMeshExtent(), MeshLib::getNumberOfGlobalNodes(), MeshToolsLib::MeshSurfaceExtraction::getSurfaceAreaForNodes(), getSurfaceIntegratedValuesForNodes(), anonymous_namespace{IdentifySubdomainMesh.cpp}::identifySubdomainMeshElements(), anonymous_namespace{IdentifySubdomainMesh.cpp}::identifySubdomainMeshNodes(), MeshToolsLib::MeshLayerMapper::layerMapping(), MeshToolsLib::MeshLayerMapper::mapToStaticValue(), anonymous_namespace{AddFaultToVoxelGrid.cpp}::markFaults(), MeshToolsLib::mergeMeshToBulkMesh(), TranslateDataDialog::moveMesh(), ApplicationUtils::partitionMesh(), MeshToolsLib::RasterDataToMesh::projectToNodes(), MeshToolsLib::removeElements(), MeshToolsLib::removeNodes(), reorderNonlinearNodes(), MeshGeoToolsLib::resetMeshElementProperty(), rotateMesh(), swapNodeCoordinateAxes(), MeshLib::IO::transformGeometry(), MeshLib::IO::transformToXDMFGeometry(), and FileIO::TetGenInterface::writeTetGenSmesh().

◆ getNumberOfElements()

std::size_t MeshLib::Mesh::getNumberOfElements ( ) const
inline

Get the number of elements.

Definition at line 88 of file Mesh.h.

88{ return _elements.size(); }

References _elements.

Referenced by ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::HydroMechanicsProcess(), Mesh(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::SmallDeformationProcess(), CreateStructuredGridDialog::accept(), addGlobalIDsToMesh(), addIntegrationPointData(), MeshToolsLib::addLayerToMesh(), MeshToolsLib::MeshLayerMapper::addLayerToMesh(), MeshLib::addPropertyToMesh(), FileIO::SwmmInterface::addResultsToMesh(), addSecondaryVariableResiduals(), MeshGeoToolsLib::appendLinesAlongPolylines(), anonymous_namespace{SubmeshResiduumOutputConfig.cpp}::checkMatchingElementCounts(), anonymous_namespace{SubmeshResiduumOutputConfig.cpp}::checkNonOverlappingCover(), anonymous_namespace{SubmeshResiduumOutputConfig.cpp}::computeNonOverlappingBulkMeshCoverBySubmeshes(), MeshToolsLib::convertMeshToGeo(), MeshToolsLib::createCellProperties(), ProcessLib::createConstraintDirichletBoundaryCondition(), ProcessLib::createDirichletBoundaryCondition(), ProcessLib::createDirichletBoundaryConditionWithinTimeInterval(), ChemistryLib::PhreeqcIOData::createEquilibriumReactants(), MeshToolsLib::createFlippedMesh(), ProcessLib::createHCNonAdvectiveFreeComponentFlowBoundaryCondition(), ChemistryLib::PhreeqcIOData::createKineticReactants(), ProcessLib::createNeumannBoundaryCondition(), ProcessLib::createPrimaryVariableConstraintDirichletBoundaryCondition(), ProcessLib::createPythonBoundaryCondition(), ProcessLib::createPythonSourceTerm(), ProcessLib::createReleaseNodalForce(), ProcessLib::createRobinBoundaryCondition(), MeshToolsLib::createSfcMeshProperties(), ProcessLib::createSolutionDependentDirichletBoundaryCondition(), MeshToolsLib::MeshLayerMapper::createStaticLayers(), ProcessLib::createTimeDecayDirichletBoundaryCondition(), ProcessLib::createVariableDependentNeumannBoundaryCondition(), ProcessLib::HeatTransportBHE::getBHEDataInMesh(), ProcessLib::LIE::getFractureMatrixDataInMesh(), MeshLib::getOrCreateMeshProperty(), anonymous_namespace{IdentifySubdomainMesh.cpp}::identifySubdomainMeshElements(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::initializeConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), main(), anonymous_namespace{AddFaultToVoxelGrid.cpp}::markFaults(), MeshAnalysisDialog::on_startButton_pressed(), setMaterialIDs(), ElementTreeModel::setMesh(), MeshToolsLib::Mesh2MeshPropertyInterpolation::setPropertiesForMesh(), MeshToolsLib::MeshValidation::testElementGeometry(), testIfMeshesAreEqual(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::updateElementLevelSets(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::updateElementLevelSets(), FileIO::SHPInterface::write2dMeshToSHP(), and FileIO::TetGenInterface::writeTetGenSmesh().

◆ getNumberOfNodes()

std::size_t MeshLib::Mesh::getNumberOfNodes ( ) const
inline

Get the number of nodes.

Definition at line 91 of file Mesh.h.

91{ return _nodes.size(); }

References _nodes.

Referenced by Mesh(), MeshLib::NodePartitionedMesh::NodePartitionedMesh(), ProcessLib::SolutionDependentDirichletBoundaryCondition::SolutionDependentDirichletBoundaryCondition(), addGlobalIDsToMesh(), LayeredVolume::addLayerBoundaries(), LayeredVolume::addLayerToMesh(), MeshToolsLib::MeshLayerMapper::addLayerToMesh(), MeshLib::addPropertyToMesh(), FileIO::SwmmInterface::addResultsToMesh(), addSecondaryVariableNodes(), ProcessLib::checkParametersOfDirichletBoundaryCondition(), MeshLib::computeNumberOfRegularHigherOrderNodesAtRank(), anonymous_namespace{convertMeshToGeo.cpp}::convertMeshNodesToGeoPoints(), ProcessLib::createConstraintDirichletBoundaryCondition(), ProcessLib::createDirichletBoundaryCondition(), ProcessLib::createDirichletBoundaryConditionWithinTimeInterval(), ProcessLib::createHCNonAdvectiveFreeComponentFlowBoundaryCondition(), ProcessLib::createNeumannBoundaryCondition(), MeshToolsLib::createNodeProperties(), NumLib::MeshComponentMap::createParallelMeshComponentMap(), ProcessLib::createPrimaryVariableConstraintDirichletBoundaryCondition(), ProcessLib::createPythonBoundaryCondition(), ProcessLib::createPythonSourceTerm(), MeshToolsLib::MeshLayerMapper::createRasterLayers(), ProcessLib::createReleaseNodalForce(), ProcessLib::createRobinBoundaryCondition(), MeshToolsLib::createSfcMeshProperties(), ProcessLib::createSolutionDependentDirichletBoundaryCondition(), MeshToolsLib::MeshLayerMapper::createStaticLayers(), MeshLib::NodeAdjacencyTable::createTable(), ProcessLib::createTimeDecayDirichletBoundaryCondition(), ProcessLib::createVariableDependentNeumannBoundaryCondition(), ProcessLib::extractInnerAndOuterNodes(), MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh(), MeshToolsLib::MeshSurfaceExtraction::getMeshSurface(), MeshLib::getOrCreateMeshProperty(), MeshToolsLib::MeshSurfaceExtraction::getSurfaceAreaForNodes(), MeshToolsLib::MeshSurfaceExtraction::getSurfaceNodes(), anonymous_namespace{IdentifySubdomainMesh.cpp}::identifySubdomainMeshElements(), anonymous_namespace{IdentifySubdomainMesh.cpp}::identifySubdomainMeshNodes(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), MeshToolsLib::MeshLayerMapper::layerMapping(), main(), MeshLib::views::meshLocations(), MeshAnalysisDialog::on_startButton_pressed(), ElementTreeModel::setMesh(), ApplicationUtils::setNumberOfNodesInPartitions(), and testIfMeshesAreEqual().

◆ getProperties() [1/2]

Properties & MeshLib::Mesh::getProperties ( )
inline

Definition at line 125 of file Mesh.h.

125{ return _properties; }

References _properties.

Referenced by MeshLib::ElementStatus::ElementStatus(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::HydroMechanicsProcess(), ProcessLib::LIE::PostProcessTool::PostProcessTool(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::SmallDeformationProcess(), ProcessLib::SolutionDependentDirichletBoundaryCondition::SolutionDependentDirichletBoundaryCondition(), ProcessLib::SurfaceFlux::SurfaceFlux(), CreateStructuredGridDialog::accept(), MeshElementRemovalDialog::accept(), ProcessLib::addBulkMeshPropertyToSubMesh(), addGlobalIDsToMesh(), MeshToolsLib::addLayerToMesh(), MeshLib::addPropertyToMesh(), MeshElementRemovalDialog::addScalarArrays(), MeshLib::bulkElementIDs(), MeshLib::bulkNodeIDs(), anonymous_namespace{SubmeshResiduumOutputConfig.cpp}::checkBulkIDMappingsPresent(), ProcessLib::checkParametersOfDirichletBoundaryCondition(), MeshToolsLib::ElementValueModification::condense(), MeshLib::VtkMeshConverter::convertScalarArrays(), MeshToolsLib::convertToLinearMesh(), MeshToolsLib::createFlippedMesh(), ParameterLib::createGroupBasedParameter(), ProcessLib::createHCNonAdvectiveFreeComponentFlowBoundaryCondition(), ParameterLib::createMeshElementParameter(), ParameterLib::createMeshNodeParameter(), MeshToolsLib::createQuadraticOrderMesh(), ParameterLib::createRandomFieldMeshElementParameter(), ProcessLib::createSourceTerm(), determineIntegrationOrder(), ProcessLib::extractElementsAlongOuterNodes(), extractMatGroup(), MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh(), MeshToolsLib::MeshSurfaceExtraction::getMeshSurface(), MeshLib::getOrCreateMeshProperty(), getSurfaceIntegratedValuesForNodes(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::initializeConcreteProcess(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >::initializeConcreteProcess(), MeshLib::IO::Legacy::MeshIO::loadMeshFromFile(), main(), MeshToolsLib::MeshGenerator::VoxelGridFromMesh::mapArray(), MeshLib::materialIDs(), MeshToolsLib::mergeMeshToBulkMesh(), MeshElementRemovalDialog::on_scalarArrayComboBox_currentIndexChanged(), ProcessLib::Output::prepareSubmesh(), MeshToolsLib::RasterDataToMesh::projectToElements(), MeshToolsLib::RasterDataToMesh::projectToNodes(), FileIO::GMSH::readGMSHMesh(), MeshToolsLib::removeElements(), MeshToolsLib::removeNodes(), MeshToolsLib::ElementValueModification::replace(), MeshGeoToolsLib::resetMeshElementProperty(), MeshLib::scaleMeshPropertyVector(), MeshToolsLib::ElementValueModification::setByElementType(), setMaterialIDs(), ElementTreeModel::setMesh(), MeshToolsLib::Mesh2MeshPropertyInterpolation::setPropertiesForMesh(), transferPropertiesFromPartitionedMeshToUnpartitionedMesh(), MeshLib::IO::transformAttributes(), MeshLib::transformMeshToNodePartitionedMesh(), anonymous_namespace{IdentifySubdomainMesh.cpp}::updateOrCheckExistingSubdomainProperty(), FileIO::TetGenInterface::write2dElements(), FileIO::SHPInterface::write2dMeshToSHP(), FileIO::TetGenInterface::write3dElements(), MeshToolsLib::MeshInformation::writePropertyVectorInformation(), and MeshToolsLib::zeroMeshFieldDataByMaterialIDs().

◆ getProperties() [2/2]

Properties const & MeshLib::Mesh::getProperties ( ) const
inline

Definition at line 126 of file Mesh.h.

126{ return _properties; }

References _properties.

◆ hasNonlinearElement()

bool MeshLib::Mesh::hasNonlinearElement ( ) const

Check if the mesh contains any nonlinear element.

Definition at line 238 of file Mesh.cpp.

239{
240 return std::any_of(
241 std::begin(_elements), std::end(_elements),
242 [](Element const* const e)
243 { return e->getNumberOfNodes() != e->getNumberOfBaseNodes(); });
244}

References _elements, MeshLib::Element::getNumberOfBaseNodes(), and MeshLib::Element::getNumberOfNodes().

◆ isAxiallySymmetric()

bool MeshLib::Mesh::isAxiallySymmetric ( ) const
inline

Definition at line 128 of file Mesh.h.

128{ return _is_axially_symmetric; }
bool _is_axially_symmetric
Definition Mesh.h:155

References _is_axially_symmetric.

Referenced by ProcessLib::SurfaceFlux::SurfaceFlux(), ProcessLib::VolumetricSourceTerm::VolumetricSourceTerm(), MeshGeoToolsLib::constructAdditionalMeshesFromGeoObjects(), ProcessLib::createBoundaryCondition(), ProcessLib::HydroMechanics::createHydroMechanicsProcess(), ProcessLib::LargeDeformation::createLargeDeformationProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::initializeConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::initializeConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::initializeConcreteProcess(), ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::HT::HTProcess::initializeConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::initializeConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::initializeConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::initializeConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::initializeConcreteProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::initializeConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::initializeConcreteProcess(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >::initializeConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::initializeConcreteProcess(), and ProcessLib::WellboreSimulator::WellboreSimulatorProcess::initializeConcreteProcess().

◆ operator=() [1/2]

Mesh & MeshLib::Mesh::operator= ( const Mesh & mesh)
delete

References Mesh().

◆ operator=() [2/2]

Mesh & MeshLib::Mesh::operator= ( Mesh && mesh)
delete

References Mesh().

◆ resetElementIDs()

void MeshLib::Mesh::resetElementIDs ( )

Resets the IDs of all mesh-elements to their position in the element vector.

Definition at line 158 of file Mesh.cpp.

159{
160 const std::size_t nElements(this->_elements.size());
161 for (unsigned i = 0; i < nElements; ++i)
162 {
163 _elements[i]->setID(i);
164 }
165}

References _elements.

Referenced by Mesh().

◆ resetNodeIDs()

void MeshLib::Mesh::resetNodeIDs ( )

Resets the IDs of all mesh-nodes to their position in the node vector.

Definition at line 149 of file Mesh.cpp.

150{
151 const std::size_t nNodes(_nodes.size());
152 for (std::size_t i = 0; i < nNodes; ++i)
153 {
154 _nodes[i]->setID(i);
155 }
156}

References _nodes.

Referenced by Mesh(), and reorderNonlinearNodes().

◆ setAxiallySymmetric()

void MeshLib::Mesh::setAxiallySymmetric ( bool is_axial_symmetric)
inline

Definition at line 129 of file Mesh.h.

129 {
130 _is_axially_symmetric = is_axial_symmetric;
131 }

References _is_axially_symmetric.

◆ setDimension()

void MeshLib::Mesh::setDimension ( )
protected

Sets the dimension of the mesh.

Definition at line 167 of file Mesh.cpp.

168{
169 const std::size_t nElements(_elements.size());
170 for (unsigned i = 0; i < nElements; ++i)
171 {
173 {
174 _mesh_dimension = _elements[i]->getDimension();
175 }
176 }
177}
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:79

References _elements, _mesh_dimension, and getDimension().

Referenced by Mesh(), and Mesh().

◆ setElementNeighbors()

void MeshLib::Mesh::setElementNeighbors ( )
protected

Fills in the neighbor-information for elements. Note: Using this implementation, an element e can only have neighbors that have the same dimensionality as e.

Definition at line 194 of file Mesh.cpp.

195{
196 std::vector<Element const*> neighbors;
197 for (auto element : _elements)
198 {
199 // create vector with all elements connected to current element
200 // (includes lots of doubles!)
201 const std::size_t nNodes(element->getNumberOfBaseNodes());
202 for (unsigned n(0); n < nNodes; ++n)
203 {
204 auto const& conn_elems(
205 _elements_connected_to_nodes[element->getNode(n)->getID()]);
206 neighbors.insert(neighbors.end(), conn_elems.begin(),
207 conn_elems.end());
208 }
209 std::sort(neighbors.begin(), neighbors.end());
210 auto const neighbors_new_end =
211 std::unique(neighbors.begin(), neighbors.end());
212
213 for (auto neighbor = neighbors.begin(); neighbor != neighbors_new_end;
214 ++neighbor)
215 {
216 std::optional<unsigned> const opposite_face_id =
217 element->addNeighbor(const_cast<Element*>(*neighbor));
218 if (opposite_face_id)
219 {
220 const_cast<Element*>(*neighbor)->setNeighbor(element,
221 *opposite_face_id);
222 }
223 }
224 neighbors.clear();
225 }
226}

References _elements, _elements_connected_to_nodes, and MeshLib::Element::setNeighbor().

Referenced by Mesh(), and Mesh().

◆ setName()

void MeshLib::Mesh::setName ( const std::string & name)
inline

Changes the name of the mesh.

Definition at line 109 of file Mesh.h.

109{ this->_name = name; }

References _name.

Referenced by VtkVisPipelineView::convertVTKToOGSMesh().

◆ shallowClean()

void MeshLib::Mesh::shallowClean ( )

Only cleans vector members _nodes and _elements, and does not touch the pointer entries of the vectors. This function can only be called in the case that the pointers of _nodes and _elements are shared with other instances of this class and are deleted by them as well.

Definition at line 123 of file Mesh.cpp.

124{
125 _elements.clear();
126 _nodes.clear();
127}

References _elements, and _nodes.

◆ ApplicationUtils::NodeWiseMeshPartitioner

Definition at line 39 of file Mesh.h.

References Mesh().

◆ removeMeshNodes

void removeMeshNodes ( Mesh & mesh,
const std::vector< std::size_t > & nodes )
friend

References Mesh().

Member Data Documentation

◆ _compute_element_neighbors

bool const MeshLib::Mesh::_compute_element_neighbors
protected

Definition at line 156 of file Mesh.h.

Referenced by Mesh(), and Mesh().

◆ _elements

std::vector<Element*> MeshLib::Mesh::_elements
protected

◆ _elements_connected_to_nodes

std::vector<std::vector<Element const*> > MeshLib::Mesh::_elements_connected_to_nodes
protected

◆ _id

std::size_t const MeshLib::Mesh::_id
protected

Definition at line 144 of file Mesh.h.

Referenced by Mesh(), Mesh(), and getID().

◆ _is_axially_symmetric

bool MeshLib::Mesh::_is_axially_symmetric = false
protected

Definition at line 155 of file Mesh.h.

Referenced by isAxiallySymmetric(), and setAxiallySymmetric().

◆ _mesh_dimension

unsigned MeshLib::Mesh::_mesh_dimension
protected

Definition at line 145 of file Mesh.h.

Referenced by Mesh(), Mesh(), getDimension(), and setDimension().

◆ _name

std::string MeshLib::Mesh::_name
protected

Definition at line 148 of file Mesh.h.

Referenced by Mesh(), Mesh(), getName(), and setName().

◆ _node_distance

std::pair<double, double> MeshLib::Mesh::_node_distance
protected

The minimal and maximal distance of nodes within an element over all elements in the mesh.

Definition at line 147 of file Mesh.h.

Referenced by Mesh(), and Mesh().

◆ _nodes

std::vector<Node*> MeshLib::Mesh::_nodes
protected

◆ _properties

Properties MeshLib::Mesh::_properties
protected

Definition at line 151 of file Mesh.h.

Referenced by Mesh(), Mesh(), getProperties(), and getProperties().


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