OGS
MeshToolsLib Namespace Reference

Namespaces

namespace  anonymous_namespace{AngleSkewMetric.cpp}
 
namespace  anonymous_namespace{ConvertToLinearMesh.cpp}
 
namespace  BoundaryExtraction
 
namespace  details
 
namespace  MeshGenerator
 
namespace  MeshGenerators
 
namespace  ProjectPointOnMesh
 
namespace  RasterDataToMesh
 Adding pixel values from a raster onto nodes or cells of a mesh.
 

Classes

struct  AngleSkewMetric
 
struct  EdgeRatioMetric
 
class  ElementQualityInterface
 
class  ElementQualityMetric
 
class  ElementSizeMetric
 
class  ElementValueModification
 A set of methods for manipulating mesh element values. More...
 
class  Mesh2MeshPropertyInterpolation
 
class  MeshInformation
 A set of tools for extracting information from a mesh. More...
 
class  MeshLayerMapper
 Manipulating and adding prism element layers to an existing 2D mesh. More...
 
class  MeshRevision
 
class  MeshSurfaceExtraction
 A set of tools concerned with extracting nodes and elements from a mesh surface. More...
 
struct  MeshValidation
 A collection of methods for testing mesh quality and correctness. More...
 
struct  RadiusEdgeRatioMetric
 
class  RasterToMesh
 Converts raster data into an OGS mesh. More...
 
struct  SizeDifferenceMetric
 

Functions

template<typename ShapeFunction >
double computeElementVolumeNumerically (MeshLib::Element const &e)
 
double computeElementVolumeNumerically (MeshLib::Element const &e)
 
bool convertMeshToGeo (const MeshLib::Mesh &mesh, GeoLib::GEOObjects &geo_objects, double const eps)
 
MeshLib::MeshconvertSurfaceToMesh (const GeoLib::Surface &sfc, const std::string &mesh_name, double eps)
 
template<typename ElementType >
int getNumberOfElementIntegrationPointsGeneral (MeshLib::IntegrationPointMetaData const &ip_meta_data)
 
int getNumberOfElementIntegrationPoints (MeshLib::IntegrationPointMetaData const &ip_meta_data, MeshLib::Element const &e)
 
std::vector< std::size_t > getIntegrationPointDataOffsetsOfMeshElements (std::vector< MeshLib::Element * > const &mesh_elements, MeshLib::PropertyVectorBase const &pv, MeshLib::Properties const &properties)
 
MeshLib::ElementextrudeElement (std::vector< MeshLib::Node * > const &subsfc_nodes, MeshLib::Element const &sfc_elem, MeshLib::PropertyVector< std::size_t > const &sfc_to_subsfc_id_map, std::map< std::size_t, std::size_t > const &subsfc_sfc_id_map)
 
MeshLib::MeshaddLayerToMesh (MeshLib::Mesh const &mesh, double thickness, std::string const &name, bool on_top, bool copy_material_ids)
 
std::unique_ptr< MeshLib::MeshconvertToLinearMesh (MeshLib::Mesh const &mesh, std::string const &new_mesh_name)
 
std::unique_ptr< MeshLib::ElementcreateFlippedElement (MeshLib::Element const &elem, std::vector< MeshLib::Node * > const &nodes)
 
std::unique_ptr< MeshLib::MeshcreateFlippedMesh (MeshLib::Mesh const &mesh)
 
unsigned lutPrismThirdNode (unsigned const id1, unsigned const id2)
 
template<typename Iterator >
void moveMeshNodes (Iterator begin, Iterator end, Eigen::Vector3d const &displacement)
 
MeshLib::MeshremoveElements (const MeshLib::Mesh &mesh, const std::vector< std::size_t > &removed_element_ids, const std::string &new_mesh_name)
 
std::vector< bool > markUnusedNodes (std::vector< MeshLib::Element * > const &elements, std::vector< MeshLib::Node * > const &nodes)
 Marks nodes not used by any of the elements.
 
void removeMarkedNodes (std::vector< bool > const &nodes_to_delete, std::vector< MeshLib::Node * > &nodes)
 Deallocates and removes nodes marked true.
 
MeshLib::MeshremoveNodes (const MeshLib::Mesh &mesh, const std::vector< std::size_t > &del_nodes_idx, const std::string &new_mesh_name)
 
std::unique_ptr< MeshLib::MeshcreateQuadraticOrderMesh (MeshLib::Mesh const &linear_mesh, bool const add_centre_node)
 
static void trackSurface (MeshLib::Element const *const element, std::vector< unsigned > &sfc_idx, unsigned const current_index)
 
template<typename T >
void processPropertyVector (MeshLib::PropertyVector< T > const &property, std::vector< std::size_t > const &id_map, MeshLib::Mesh &sfc_mesh)
 
bool createSfcMeshProperties (MeshLib::Mesh &sfc_mesh, MeshLib::Properties const &properties, std::vector< std::size_t > const &node_ids_map, std::vector< std::size_t > const &element_ids_map)
 
std::tuple< std::vector< MeshLib::Node * >, std::vector< std::size_t > > createNodesAndIDMapFromElements (std::vector< MeshLib::Element * > const &elements, std::size_t const n_all_nodes)
 
void createSurfaceElementsFromElement (MeshLib::Element const &surface_element, std::vector< MeshLib::Element * > &surface_elements, std::vector< std::size_t > &element_to_bulk_element_id_map, std::vector< std::size_t > &element_to_bulk_face_id_map)
 
std::tuple< std::vector< MeshLib::Element * >, std::vector< std::size_t >, std::vector< std::size_t > > createBoundaryElements (MeshLib::Mesh const &bulk_mesh)
 
void addBulkIDPropertiesToMesh (MeshLib::Mesh &surface_mesh, std::string_view node_to_bulk_node_id_map_name, std::vector< std::size_t > const &node_to_bulk_node_id_map, std::string_view element_to_bulk_element_id_map_name, std::vector< std::size_t > const &element_to_bulk_element_id_map, std::string_view element_to_bulk_face_id_map_name, std::vector< std::size_t > const &element_to_bulk_face_id_map)
 
void zeroMeshFieldDataByMaterialIDs (MeshLib::Mesh &mesh, std::vector< int > const &selected_material_ids)
 

Function Documentation

◆ addBulkIDPropertiesToMesh()

void MeshToolsLib::addBulkIDPropertiesToMesh ( MeshLib::Mesh & surface_mesh,
std::string_view node_to_bulk_node_id_map_name,
std::vector< std::size_t > const & node_to_bulk_node_id_map,
std::string_view element_to_bulk_element_id_map_name,
std::vector< std::size_t > const & element_to_bulk_element_id_map,
std::string_view element_to_bulk_face_id_map_name,
std::vector< std::size_t > const & element_to_bulk_face_id_map )

Definition at line 529 of file MeshSurfaceExtraction.cpp.

537{
538 // transmit the original node ids of the bulk mesh as a property
539 if (!node_to_bulk_node_id_map_name.empty())
540 {
541 MeshLib::addPropertyToMesh(surface_mesh, node_to_bulk_node_id_map_name,
543 node_to_bulk_node_id_map);
544 }
545
546 // transmit the original bulk element ids as a property
547 if (!element_to_bulk_element_id_map_name.empty())
548 {
550 surface_mesh, element_to_bulk_element_id_map_name,
551 MeshLib::MeshItemType::Cell, 1, element_to_bulk_element_id_map);
552 }
553
554 // transmit the face id of the original bulk element as a property
555 if (!element_to_bulk_face_id_map_name.empty())
556 {
558 surface_mesh, element_to_bulk_face_id_map_name,
559 MeshLib::MeshItemType::Cell, 1, element_to_bulk_face_id_map);
560 }
561}
void addPropertyToMesh(Mesh &mesh, std::string_view name, MeshItemType item_type, std::size_t number_of_components, std::vector< T > const &values)

References MeshLib::addPropertyToMesh(), MeshLib::Cell, and MeshLib::Node.

Referenced by MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh(), and MeshToolsLib::MeshSurfaceExtraction::getMeshSurface().

◆ addLayerToMesh()

MeshLib::Mesh * MeshToolsLib::addLayerToMesh ( MeshLib::Mesh const & mesh,
double thickness,
std::string const & name,
bool on_top,
bool copy_material_ids )

Adds a layer to the mesh. If on_top is true, the layer is added on top, if it is false, the layer is added at the bottom.

Definition at line 84 of file AddLayerToMesh.cpp.

87{
88 INFO("Extracting top surface of mesh '{:s}' ... ", mesh.getName());
89 double const flag = on_top ? -1.0 : 1.0;
90 Eigen::Vector3d const dir({0, 0, flag});
91 double const angle(90);
92 std::unique_ptr<MeshLib::Mesh> sfc_mesh(nullptr);
93
94 auto const prop_name =
96
97 if (mesh.getDimension() == 3)
98 {
100 mesh, dir, angle, prop_name));
101 }
102 else
103 {
104 sfc_mesh = (on_top) ? std::make_unique<MeshLib::Mesh>(mesh)
105 : std::unique_ptr<MeshLib::Mesh>(
107 // add property storing node ids
108 auto* const pv =
109 sfc_mesh->getProperties().createNewPropertyVector<std::size_t>(
110 prop_name, MeshLib::MeshItemType::Node, 1);
111 if (pv)
112 {
113 pv->resize(sfc_mesh->getNumberOfNodes());
114 std::iota(pv->begin(), pv->end(), 0);
115 }
116 else
117 {
118 ERR("Could not create and initialize property '{}'.", prop_name);
119 return nullptr;
120 }
121 }
122 INFO("done.");
123
124 // *** add new surface nodes
125 std::vector<MeshLib::Node*> subsfc_nodes =
126 MeshLib::copyNodeVector(mesh.getNodes());
127 std::vector<MeshLib::Element*> subsfc_elements =
128 MeshLib::copyElementVector(mesh.getElements(), subsfc_nodes);
129
130 std::size_t const n_subsfc_nodes(subsfc_nodes.size());
131
132 std::vector<MeshLib::Node*> const& sfc_nodes(sfc_mesh->getNodes());
133 std::size_t const n_sfc_nodes(sfc_nodes.size());
134
135 if (!sfc_mesh->getProperties().existsPropertyVector<std::size_t>(prop_name))
136 {
137 ERR("Need subsurface node ids, but the property '{:s}' is not "
138 "available.",
139 prop_name);
140 return nullptr;
141 }
142 // fetch subsurface node ids PropertyVector
143 auto const* const node_id_pv =
144 sfc_mesh->getProperties().getPropertyVector<std::size_t>(prop_name);
145
146 // *** copy sfc nodes to subsfc mesh node
147 std::map<std::size_t, std::size_t> subsfc_sfc_id_map;
148 for (std::size_t k(0); k < n_sfc_nodes; ++k)
149 {
150 std::size_t const subsfc_id((*node_id_pv)[k]);
151 std::size_t const sfc_id(k + n_subsfc_nodes);
152 subsfc_sfc_id_map.insert(std::make_pair(subsfc_id, sfc_id));
153 MeshLib::Node const& node(*sfc_nodes[k]);
154 subsfc_nodes.push_back(new MeshLib::Node(
155 node[0], node[1], node[2] - (flag * thickness), sfc_id));
156 }
157
158 // *** insert new layer elements into subsfc_mesh
159 std::vector<MeshLib::Element*> const& sfc_elements(sfc_mesh->getElements());
160 std::size_t const n_sfc_elements(sfc_elements.size());
161 for (std::size_t k(0); k < n_sfc_elements; ++k)
162 {
163 subsfc_elements.push_back(extrudeElement(
164 subsfc_nodes, *sfc_elements[k], *node_id_pv, subsfc_sfc_id_map));
165 }
166
167 auto new_mesh = new MeshLib::Mesh(name, subsfc_nodes, subsfc_elements,
168 true /* compute_element_neighbors */);
169
170 if (!mesh.getProperties().existsPropertyVector<int>("MaterialIDs"))
171 {
172 ERR("Could not copy the property 'MaterialIDs' since the original mesh "
173 "does not contain such a property.");
174 return new_mesh;
175 }
176 auto const* const materials =
177 mesh.getProperties().getPropertyVector<int>("MaterialIDs");
178
179 auto* const new_materials =
180 new_mesh->getProperties().createNewPropertyVector<int>(
181 "MaterialIDs", MeshLib::MeshItemType::Cell, 1);
182 if (!new_materials)
183 {
184 ERR("Can not set material properties for new layer");
185 return new_mesh;
186 }
187
188 new_materials->reserve(subsfc_elements.size());
189 std::copy(materials->cbegin(), materials->cend(),
190 std::back_inserter(*new_materials));
191
192 if (copy_material_ids &&
193 sfc_mesh->getProperties().existsPropertyVector<int>("MaterialIDs"))
194 {
195 auto const& sfc_material_ids =
196 *sfc_mesh->getProperties().getPropertyVector<int>("MaterialIDs");
197 std::copy(sfc_material_ids.cbegin(), sfc_material_ids.cend(),
198 std::back_inserter(*new_materials));
199 }
200 else
201 {
202 int const new_layer_id(
203 *(std::max_element(materials->cbegin(), materials->cend())) + 1);
204 auto const n_new_props(subsfc_elements.size() -
205 mesh.getNumberOfElements());
206 std::fill_n(std::back_inserter(*new_materials), n_new_props,
207 new_layer_id);
208 }
209
210 return new_mesh;
211}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
static MeshLib::Mesh * getMeshSurface(const MeshLib::Mesh &subsfc_mesh, Eigen::Vector3d const &dir, double angle, std::string_view subsfc_node_id_prop_name="", std::string_view subsfc_element_id_prop_name="", std::string_view face_id_prop_name="")
std::vector< Node * > copyNodeVector(const std::vector< Node * > &nodes)
Creates a deep copy of a Node vector.
constexpr std::string_view getBulkIDString(MeshItemType mesh_item_type)
Definition Properties.h:188
std::vector< Element * > copyElementVector(std::vector< Element * > const &elements, std::vector< Node * > const &new_nodes, std::vector< std::size_t > const *const node_id_map)
std::unique_ptr< MeshLib::Mesh > createFlippedMesh(MeshLib::Mesh const &mesh)
MeshLib::Element * extrudeElement(std::vector< MeshLib::Node * > const &subsfc_nodes, MeshLib::Element const &sfc_elem, MeshLib::PropertyVector< std::size_t > const &sfc_to_subsfc_id_map, std::map< std::size_t, std::size_t > const &subsfc_sfc_id_map)

References MeshLib::Cell, MeshLib::copyElementVector(), MeshLib::copyNodeVector(), createFlippedMesh(), ERR(), MeshLib::Properties::existsPropertyVector(), extrudeElement(), MeshLib::getBulkIDString(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), MeshToolsLib::MeshSurfaceExtraction::getMeshSurface(), MeshLib::Mesh::getName(), MeshLib::Mesh::getNodes(), MeshLib::Mesh::getNumberOfElements(), MeshLib::Mesh::getProperties(), MeshLib::Properties::getPropertyVector(), INFO(), and MeshLib::Node.

Referenced by MeshToolsLib::MeshGenerator::generateRegularPrismMesh(), main(), and MeshView::openAddLayerDialog().

◆ computeElementVolumeNumerically() [1/2]

template<typename ShapeFunction >
double MeshToolsLib::computeElementVolumeNumerically ( MeshLib::Element const & e)

Definition at line 47 of file ComputeElementVolumeNumerically.cpp.

48{
49 // Space dimension is set to 3 in case that 1D or 2D element is inclined.
50 constexpr int space_dim = 3;
52
53 // Integration order is set to 3:
54 auto const& integration_method =
55 NumLib::IntegrationMethodRegistry::template getIntegrationMethod<
56 typename ShapeFunction::MeshElement>(NumLib::IntegrationOrder{3});
57
58 auto const shape_function_data =
59 NumLib::initShapeMatrices<ShapeFunction, ShapeMatricesType, space_dim>(
60 e, false /*is_axially_symmetric*/, integration_method);
61
62 auto const n_integration_points = integration_method.getNumberOfPoints();
63 double volume = 0.0;
64 for (unsigned ip = 0; ip < n_integration_points; ++ip)
65 {
66 auto const weight = integration_method.getWeightedPoint(ip).getWeight();
67 volume += shape_function_data[ip].detJ * weight;
68 }
69
70 return volume;
71}

Referenced by MeshToolsLib::ElementSizeMetric::calc1dQuality(), and MeshToolsLib::ElementSizeMetric::calc2dOr3dQuality().

◆ computeElementVolumeNumerically() [2/2]

double MeshToolsLib::computeElementVolumeNumerically ( MeshLib::Element const & e)

Definition at line 73 of file ComputeElementVolumeNumerically.cpp.

74{
75 switch (e.getCellType())
76 {
78 return computeElementVolumeNumerically<NumLib::ShapeLine2>(e);
80 return computeElementVolumeNumerically<NumLib::ShapeLine3>(e);
82 return computeElementVolumeNumerically<NumLib::ShapeTri3>(e);
84 return computeElementVolumeNumerically<NumLib::ShapeTri6>(e);
86 return computeElementVolumeNumerically<NumLib::ShapeQuad4>(e);
88 return computeElementVolumeNumerically<NumLib::ShapeQuad8>(e);
90 return computeElementVolumeNumerically<NumLib::ShapeQuad9>(e);
92 return computeElementVolumeNumerically<NumLib::ShapeTet4>(e);
94 return computeElementVolumeNumerically<NumLib::ShapeHex8>(e);
96 return computeElementVolumeNumerically<NumLib::ShapeHex20>(e);
98 return computeElementVolumeNumerically<NumLib::ShapeTet10>(e);
100 return computeElementVolumeNumerically<NumLib::ShapePrism6>(e);
102 return computeElementVolumeNumerically<NumLib::ShapePrism15>(e);
104 return computeElementVolumeNumerically<NumLib::ShapePyra5>(e);
106 return computeElementVolumeNumerically<NumLib::ShapePyra13>(e);
107 default:
108 OGS_FATAL(
109 "Numerical volume calculation is not available for element "
110 "with type {}. ",
111 MeshLib::CellType2String(e.getCellType()));
112 }
113}
#define OGS_FATAL(...)
Definition Error.h:26
std::string CellType2String(const CellType t)
Given a MeshElemType this returns the appropriate string.

References MeshLib::CellType2String(), MeshLib::Element::getCellType(), MeshLib::HEX20, MeshLib::HEX8, MeshLib::LINE2, MeshLib::LINE3, OGS_FATAL, MeshLib::PRISM15, MeshLib::PRISM6, MeshLib::PYRAMID13, MeshLib::PYRAMID5, MeshLib::QUAD4, MeshLib::QUAD8, MeshLib::QUAD9, MeshLib::TET10, MeshLib::TET4, MeshLib::TRI3, and MeshLib::TRI6.

◆ convertMeshToGeo()

bool MeshToolsLib::convertMeshToGeo ( const MeshLib::Mesh & mesh,
GeoLib::GEOObjects & geo_objects,
double eps = std::numeric_limits< double >::epsilon() )

Converts a 2D mesh into a geometry. A new geometry with the name of the mesh will be inserted into geo_objects, consisting of points identical with mesh nodes and one surface representing the mesh. Triangles are converted to geometric triangles, quads are split into two triangles, all other elements are ignored.

Definition at line 77 of file convertMeshToGeo.cpp.

80{
81 if (mesh.getDimension() != 2)
82 {
83 ERR("Mesh to geometry conversion is only working for 2D meshes.");
84 return false;
85 }
86
87 // Special handling of the bounds in case there are no materialIDs present.
88 auto get_material_ids_and_bounds = [&]()
89 -> std::tuple<MeshLib::PropertyVector<int> const*, std::pair<int, int>>
90 {
91 auto const materialIds = materialIDs(mesh);
92 if (!materialIds)
93 {
94 return std::make_tuple(nullptr, std::make_pair(0, 0));
95 }
96
97 auto const bounds =
99 if (!bounds)
100 {
101 OGS_FATAL(
102 "Could not get minimum/maximum ranges values for the "
103 "MaterialIDs property in the mesh '{:s}'.",
104 mesh.getName());
105 }
106 return std::make_tuple(materialIds, *bounds);
107 };
108
109 auto const [materialIds, bounds] = get_material_ids_and_bounds();
110 // elements to surface triangles conversion
111 const unsigned nMatGroups(bounds.second - bounds.first + 1);
112 std::vector<GeoLib::Surface*> sfcs;
113 sfcs.reserve(nMatGroups);
114 std::string const geoobject_name =
115 convertMeshNodesToGeoPoints(mesh, eps, geo_objects);
116 auto const& points = *geo_objects.getPointVec(geoobject_name);
117 for (unsigned i = 0; i < nMatGroups; ++i)
118 {
119 sfcs.push_back(new GeoLib::Surface(points));
120 }
121
122 const std::vector<std::size_t>& id_map(
123 geo_objects.getPointVecObj(geoobject_name)->getIDMap());
124 const std::vector<MeshLib::Element*>& elements = mesh.getElements();
125 const std::size_t nElems(mesh.getNumberOfElements());
126
127 for (unsigned i = 0; i < nElems; ++i)
128 {
129 auto surfaceId = !materialIds ? 0 : ((*materialIds)[i] - bounds.first);
130 addElementToSurface(*elements[i], id_map, *sfcs[surfaceId]);
131 }
132
133 std::for_each(sfcs.begin(), sfcs.end(),
134 [](GeoLib::Surface*& sfc)
135 {
136 if (sfc->getNumberOfTriangles() == 0)
137 {
138 delete sfc;
139 sfc = nullptr;
140 }
141 });
142 auto sfcs_end = std::remove(sfcs.begin(), sfcs.end(), nullptr);
143 sfcs.erase(sfcs_end, sfcs.end());
144
145 geo_objects.addSurfaceVec(std::move(sfcs), geoobject_name,
147 return true;
148}
const std::vector< Point * > * getPointVec(const std::string &name) const
const PointVec * getPointVecObj(const std::string &name) const
void addSurfaceVec(std::vector< Surface * > &&sfc, const std::string &name, SurfaceVec::NameIdMap &&sfc_names)
const std::vector< std::size_t > & getIDMap() const
Definition PointVec.h:94
A Surface is represented by Triangles. It consists of a reference to a vector of (pointers to) points...
Definition Surface.h:33
std::map< std::string, std::size_t > NameIdMap
Definition TemplateVec.h:41
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:103
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:97
static std::optional< std::pair< T, T > > const getValueBounds(MeshLib::PropertyVector< T > const &property)
void addElementToSurface(MeshLib::Element const &e, std::vector< std::size_t > const &id_map, GeoLib::Surface &surface)
std::string convertMeshNodesToGeoPoints(MeshLib::Mesh const &mesh, double const eps, GeoLib::GEOObjects &geo_objects)

References ERR(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), GeoLib::PointVec::getIDMap(), MeshLib::Mesh::getName(), MeshLib::Mesh::getNumberOfElements(), GeoLib::GEOObjects::getPointVec(), GeoLib::GEOObjects::getPointVecObj(), MeshToolsLib::MeshInformation::getValueBounds(), and OGS_FATAL.

Referenced by MainWindow::convertMeshToGeometry(), FileIO::createSurface(), and main().

◆ convertSurfaceToMesh()

MeshLib::Mesh * MeshToolsLib::convertSurfaceToMesh ( const GeoLib::Surface & sfc,
const std::string & mesh_name,
double eps = std::numeric_limits< double >::epsilon() )

Converts a surface into a triangular mesh

Parameters
sfcSurface object
mesh_nameNew mesh name
epsMinimum distance for nodes not to be collapsed
Returns
a pointer to a converted mesh object. nullptr is returned if the conversion fails.

Definition at line 150 of file convertMeshToGeo.cpp.

153{
154 // convert to a mesh including duplicated nodes
155 std::vector<MeshLib::Node*> nodes;
156 std::vector<MeshLib::Element*> elements;
157 std::size_t nodeId = 0;
158 for (std::size_t i = 0; i < sfc.getNumberOfTriangles(); i++)
159 {
160 auto* tri = sfc[i];
161 auto** tri_nodes = new MeshLib::Node*[3];
162 for (unsigned j = 0; j < 3; j++)
163 {
164 tri_nodes[j] =
165 new MeshLib::Node(tri->getPoint(j)->data(), nodeId++);
166 }
167 elements.push_back(new MeshLib::Tri(tri_nodes, i));
168 for (unsigned j = 0; j < 3; j++)
169 {
170 nodes.push_back(tri_nodes[j]);
171 }
172 }
173 MeshLib::Mesh mesh_with_duplicated_nodes(
174 mesh_name, nodes, elements, true /* compute_element_neighbors */);
175
176 // remove duplicated nodes
177 MeshToolsLib::MeshRevision rev(mesh_with_duplicated_nodes);
178 return rev.simplifyMesh(mesh_with_duplicated_nodes.getName(), eps);
179}
std::size_t getNumberOfTriangles() const
Definition Surface.cpp:82

References convertSurfaceToMesh(), MeshLib::Mesh::getName(), GeoLib::Surface::getNumberOfTriangles(), MeshLib::Node, and MeshToolsLib::MeshRevision::simplifyMesh().

Referenced by convertSurfaceToMesh(), and main().

◆ convertToLinearMesh()

std::unique_ptr< MeshLib::Mesh > MeshToolsLib::convertToLinearMesh ( const MeshLib::Mesh & mesh,
const std::string & new_mesh_name )

Converts a non-linear mesh to a linear mesh. All the mesh properties will be copied except for entries for non-linear nodes.

Definition at line 47 of file ConvertToLinearMesh.cpp.

49{
50 auto const& org_elements = mesh.getElements();
51
52 // mark base nodes
53 std::vector<bool> marked_base_nodes(mesh.getNodes().size(), false);
54 for (auto const org_element : org_elements)
55 {
56 for (std::size_t k = 0; k < org_element->getNumberOfBaseNodes(); ++k)
57 {
58 auto const& base_node = *org_element->getNode(k);
59 marked_base_nodes[base_node.getID()] = true;
60 }
61 }
62
63 // construct map and fill new_mesh_nodes
64 std::vector<MeshLib::Node*> new_mesh_nodes{static_cast<std::size_t>(
65 std::count(begin(marked_base_nodes), end(marked_base_nodes), true))};
66 std::size_t base_node_cnt = 0;
67 auto const& org_nodes = mesh.getNodes();
68 std::vector<std::size_t> base_node_map(org_nodes.size(), -1);
69 for (std::size_t k = 0; k < org_nodes.size(); ++k)
70 {
71 if (marked_base_nodes[k])
72 {
73 new_mesh_nodes[base_node_cnt] =
74 new MeshLib::Node(org_nodes[k]->data(), base_node_cnt);
75 base_node_map[k] = base_node_cnt;
76 base_node_cnt++;
77 }
78 }
79
80 // create new elements with the quadratic nodes
81 std::vector<MeshLib::Element*> vec_new_eles;
82 for (MeshLib::Element const* e : mesh.getElements())
83 {
84 if (e->getCellType() == MeshLib::CellType::LINE3)
85 {
86 vec_new_eles.push_back(createLinearElement<MeshLib::Line>(
87 e, new_mesh_nodes, base_node_map));
88 }
89 else if (e->getCellType() == MeshLib::CellType::QUAD8)
90 {
91 vec_new_eles.push_back(createLinearElement<MeshLib::Quad>(
92 e, new_mesh_nodes, base_node_map));
93 }
94 else if (e->getCellType() == MeshLib::CellType::TRI6)
95 {
96 vec_new_eles.push_back(createLinearElement<MeshLib::Tri>(
97 e, new_mesh_nodes, base_node_map));
98 }
99 else if (e->getCellType() == MeshLib::CellType::HEX20)
100 {
101 vec_new_eles.push_back(createLinearElement<MeshLib::Hex>(
102 e, new_mesh_nodes, base_node_map));
103 }
104 else if (e->getCellType() == MeshLib::CellType::TET10)
105 {
106 vec_new_eles.push_back(createLinearElement<MeshLib::Tet>(
107 e, new_mesh_nodes, base_node_map));
108 }
109 else
110 {
111 OGS_FATAL("Mesh element type {:s} is not supported",
112 MeshLib::CellType2String(e->getCellType()));
113 }
114 }
115
116 auto new_mesh = std::make_unique<MeshLib::Mesh>(
117 new_mesh_name, new_mesh_nodes, vec_new_eles,
118 true /* compute_element_neighbors */,
119 mesh.getProperties().excludeCopyProperties(
120 std::vector<MeshLib::MeshItemType>(1,
122
123 // copy property vectors for nodes
124 for (auto [name, property] : mesh.getProperties())
125 {
126 if (property->getMeshItemType() != MeshLib::MeshItemType::Node)
127 {
128 continue;
129 }
130 auto double_property =
131 dynamic_cast<MeshLib::PropertyVector<double>*>(property);
132 if (double_property == nullptr)
133 {
134 continue;
135 }
136 auto const n_src_comp = double_property->getNumberOfGlobalComponents();
137 auto new_prop =
138 new_mesh->getProperties().createNewPropertyVector<double>(
139 name, MeshLib::MeshItemType::Node, n_src_comp);
140 new_prop->resize(new_mesh->getNumberOfNodes() * n_src_comp);
141
142 for (std::size_t k = 0; k < org_nodes.size(); ++k)
143 {
144 if (!marked_base_nodes[k])
145 {
146 continue;
147 }
148 std::copy_n(double_property->begin() + k * n_src_comp,
149 n_src_comp,
150 new_prop->begin() + base_node_map[k] * n_src_comp);
151 }
152 }
153 return new_mesh;
154}
int getNumberOfGlobalComponents() const

References MeshLib::CellType2String(), MeshLib::Properties::excludeCopyProperties(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::PropertyVectorBase::getNumberOfGlobalComponents(), MeshLib::Mesh::getProperties(), MeshLib::HEX20, MeshLib::LINE3, MeshLib::Node, OGS_FATAL, MeshLib::QUAD8, MeshLib::TET10, and MeshLib::TRI6.

Referenced by main(), and anonymous_namespace{postLIE.cpp}::postVTU().

◆ createBoundaryElements()

std::tuple< std::vector< MeshLib::Element * >, std::vector< std::size_t >, std::vector< std::size_t > > MeshToolsLib::createBoundaryElements ( MeshLib::Mesh const & bulk_mesh)

Definition at line 444 of file MeshSurfaceExtraction.cpp.

445{
446 std::vector<std::size_t> element_to_bulk_element_id_map;
447 std::vector<std::size_t> element_to_bulk_face_id_map;
448 std::vector<MeshLib::Element*> surface_elements;
449
450 auto const& bulk_elements = bulk_mesh.getElements();
451 auto const mesh_dimension = bulk_mesh.getDimension();
452
453 for (auto const* elem : bulk_elements)
454 {
455 const unsigned element_dimension(elem->getDimension());
456 if (element_dimension < mesh_dimension)
457 {
458 continue;
459 }
460
461 if (!elem->isBoundaryElement())
462 {
463 continue;
464 }
465 createSurfaceElementsFromElement(*elem, surface_elements,
466 element_to_bulk_element_id_map,
467 element_to_bulk_face_id_map);
468 }
469 return {surface_elements, element_to_bulk_element_id_map,
470 element_to_bulk_face_id_map};
471}
void createSurfaceElementsFromElement(MeshLib::Element const &surface_element, std::vector< MeshLib::Element * > &surface_elements, std::vector< std::size_t > &element_to_bulk_element_id_map, std::vector< std::size_t > &element_to_bulk_face_id_map)

References createSurfaceElementsFromElement(), MeshLib::Mesh::getDimension(), and MeshLib::Mesh::getElements().

Referenced by MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh().

◆ createFlippedElement()

std::unique_ptr< MeshLib::Element > MeshToolsLib::createFlippedElement ( MeshLib::Element const & elem,
std::vector< MeshLib::Node * > const & nodes )

Creates a copy of an 1d / 2d element where the node order is reversed, such that the direction of a line changes and normals of 2D elements changes its sign.

Parameters
elemoriginal element
nodesnode vector used for the copy of the element
Returns
a flipped copy of the element

Definition at line 22 of file FlipElements.cpp.

24{
25 if (elem.getDimension() > 2)
26 {
27 return nullptr;
28 }
29
30 unsigned const n_nodes(elem.getNumberOfNodes());
31 auto elem_nodes = std::make_unique<MeshLib::Node*[]>(n_nodes);
32 for (unsigned i = 0; i < n_nodes; ++i)
33 {
34 elem_nodes[i] = nodes[elem.getNode(i)->getID()];
35 }
36 std::swap(elem_nodes[0], elem_nodes[1]);
37
38 if (elem.getGeomType() == MeshLib::MeshElemType::LINE)
39 {
40 return std::make_unique<MeshLib::Line>(elem_nodes.release(),
41 elem.getID());
42 }
43 if (elem.getGeomType() == MeshLib::MeshElemType::TRIANGLE)
44 {
45 return std::make_unique<MeshLib::Tri>(elem_nodes.release(),
46 elem.getID());
47 }
48 if (elem.getGeomType() == MeshLib::MeshElemType::QUAD)
49 {
50 std::swap(elem_nodes[2], elem_nodes[3]);
51 return std::make_unique<MeshLib::Quad>(elem_nodes.release(),
52 elem.getID());
53 }
54 return nullptr;
55}

References MeshLib::Element::getDimension(), MeshLib::Element::getGeomType(), MathLib::Point3dWithID::getID(), MeshLib::Element::getID(), MeshLib::Element::getNode(), MeshLib::Element::getNumberOfNodes(), MeshLib::LINE, MeshLib::QUAD, and MeshLib::TRIANGLE.

Referenced by createFlippedMesh().

◆ createFlippedMesh()

std::unique_ptr< MeshLib::Mesh > MeshToolsLib::createFlippedMesh ( MeshLib::Mesh const & mesh)

Creates a copy of a 1d / 2d mesh where the node order of all elements is reversed such that the direction of lines changes and normals of 2D elements changes their sign.

Parameters
meshinput mesh
Returns
a flipped copy of the input mesh

Definition at line 57 of file FlipElements.cpp.

58{
59 if (mesh.getDimension() > 2)
60 {
61 return nullptr;
62 }
63
64 std::vector<MeshLib::Node*> new_nodes(copyNodeVector(mesh.getNodes()));
65 std::vector<MeshLib::Element*> const& elems(mesh.getElements());
66 std::vector<MeshLib::Element*> new_elems;
67 std::size_t n_elems(mesh.getNumberOfElements());
68 new_elems.reserve(n_elems);
69
70 for (std::size_t i = 0; i < n_elems; ++i)
71 {
72 new_elems.push_back(
73 createFlippedElement(*elems[i], new_nodes).release());
74 }
75
76 return std::make_unique<MeshLib::Mesh>(
77 "FlippedElementMesh", new_nodes, new_elems,
78 true /* compute_element_neighbors */, mesh.getProperties());
79}

References createFlippedElement(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::Mesh::getNumberOfElements(), and MeshLib::Mesh::getProperties().

Referenced by addLayerToMesh().

◆ createNodesAndIDMapFromElements()

std::tuple< std::vector< MeshLib::Node * >, std::vector< std::size_t > > MeshToolsLib::createNodesAndIDMapFromElements ( std::vector< MeshLib::Element * > const & elements,
std::size_t const n_all_nodes )

Definition at line 175 of file MeshSurfaceExtraction.cpp.

177{
178 std::vector<MeshLib::Node const*> tmp_nodes(n_all_nodes, nullptr);
179 for (auto const* elem : elements)
180 {
181 auto const n_nodes = elem->getNumberOfNodes();
182 for (unsigned j = 0; j < n_nodes; ++j)
183 {
184 const MeshLib::Node* node(elem->getNode(j));
185 tmp_nodes[node->getID()] = node;
186 }
187 }
188
189 std::vector<MeshLib::Node*> nodes;
190 std::vector<std::size_t> node_id_map(n_all_nodes);
191 for (unsigned i = 0; i < n_all_nodes; ++i)
192 {
193 if (tmp_nodes[i])
194 {
195 node_id_map[i] = nodes.size();
196 nodes.push_back(new MeshLib::Node(*tmp_nodes[i]));
197 }
198 }
199 return {nodes, node_id_map};
200}

References MathLib::Point3dWithID::getID().

Referenced by MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh(), MeshToolsLib::MeshSurfaceExtraction::getMeshSurface(), and MeshToolsLib::MeshSurfaceExtraction::getSurfaceNodes().

◆ createQuadraticOrderMesh()

std::unique_ptr< MeshLib::Mesh > MeshToolsLib::createQuadraticOrderMesh ( MeshLib::Mesh const & linear_mesh,
bool const add_centre_node )

Create a quadratic order mesh from the linear order mesh. For some element types like Quad-4, a centre node might be added if the add_centre_node flag is set, yielding a Quad-9.

Definition at line 195 of file QuadraticMeshGenerator.cpp.

197{
198 // Clone the linear mesh nodes.
199 auto quadratic_mesh_nodes = MeshLib::copyNodeVector(linear_mesh.getNodes());
200
201 // Temporary container for unique quadratic nodes with O(log(n)) search.
202 std::set<MeshLib::Node*, nodeByCoordinatesComparator> unique_nodes;
203
204 // Create new elements with the quadratic nodes
205 std::vector<MeshLib::Element*> quadratic_elements;
206 auto const& linear_mesh_elements = linear_mesh.getElements();
207 for (MeshLib::Element const* e : linear_mesh_elements)
208 {
209 auto quadratic_element = createQuadraticElement(*e, add_centre_node);
210
211 // Replace the base nodes with cloned linear nodes.
212 int const number_base_nodes = quadratic_element->getNumberOfBaseNodes();
213 for (int i = 0; i < number_base_nodes; ++i)
214 {
215 quadratic_element->setNode(
216 i, quadratic_mesh_nodes[getNodeIndex(*quadratic_element, i)]);
217 }
218
219 // Make the new (middle-edge) nodes unique.
220 int const number_all_nodes = quadratic_element->getNumberOfNodes();
221 for (int i = number_base_nodes; i < number_all_nodes; ++i)
222 {
223 MeshLib::Node* original_node =
224 const_cast<MeshLib::Node*>(quadratic_element->getNode(i));
225
226 auto it = unique_nodes.insert(original_node);
227 if (!it.second) // same node was already inserted before, no
228 // insertion
229 {
230 // Replace the element's node with the unique node.
231 quadratic_element->setNode(i, *it.first);
232 // And delete the original node
233 delete original_node;
234 }
235 }
236
237 quadratic_elements.push_back(quadratic_element.release());
238 }
239
240 // Add the unique quadratic nodes to the cloned linear nodes.
241 quadratic_mesh_nodes.reserve(linear_mesh.getNodes().size() +
242 unique_nodes.size());
243 std::copy(unique_nodes.begin(), unique_nodes.end(),
244 std::back_inserter(quadratic_mesh_nodes));
245
246 return std::make_unique<MeshLib::Mesh>(
247 linear_mesh.getName(), quadratic_mesh_nodes, quadratic_elements,
248 true /* compute_element_neighbors */,
249 linear_mesh.getProperties().excludeCopyProperties(
250 std::vector<MeshLib::MeshItemType>(1,
252}
std::unique_ptr< MeshLib::Element > createQuadraticElement(MeshLib::Element const &e, bool const add_centre_node)
Return a new quadratic element corresponding to the linear element's type.

References MeshLib::copyNodeVector(), createQuadraticElement(), MeshLib::Properties::excludeCopyProperties(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getName(), MeshLib::Mesh::getNodes(), MeshLib::Mesh::getProperties(), and MeshLib::Node.

Referenced by main().

◆ createSfcMeshProperties()

bool MeshToolsLib::createSfcMeshProperties ( MeshLib::Mesh & sfc_mesh,
MeshLib::Properties const & properties,
std::vector< std::size_t > const & node_ids_map,
std::vector< std::size_t > const & element_ids_map )

Definition at line 54 of file MeshSurfaceExtraction.cpp.

58{
59 std::size_t const n_elems(sfc_mesh.getNumberOfElements());
60 std::size_t const n_nodes(sfc_mesh.getNumberOfNodes());
61 if (node_ids_map.size() != n_nodes)
62 {
63 ERR("createSfcMeshProperties() - Incorrect number of node IDs ({:d}) "
64 "compared to actual number of surface nodes ({:d}).",
65 node_ids_map.size(), n_nodes);
66 return false;
67 }
68
69 if (element_ids_map.size() != n_elems)
70 {
71 ERR("createSfcMeshProperties() - Incorrect number of element IDs "
72 "({:d}) compared to actual number of surface elements ({:d}).",
73 element_ids_map.size(), n_elems);
74 return false;
75 }
76 std::map<MeshLib::MeshItemType, std::vector<std::size_t> const*> const
77 id_maps = {{MeshLib::MeshItemType::Cell, &element_ids_map},
78 {MeshLib::MeshItemType::Node, &node_ids_map}};
79
80 std::size_t vectors_copied(0);
81 std::size_t vectors_skipped(0);
82 for (auto [name, property] : properties)
83 {
84 if (property->getMeshItemType() != MeshLib::MeshItemType::Cell &&
85 property->getMeshItemType() != MeshLib::MeshItemType::Node)
86 {
87 WARN(
88 "Skipping property vector '{:s}' - not defined on cells or "
89 "nodes.",
90 name);
91 vectors_skipped++;
92 continue;
93 }
94
95 auto const& id_map = *id_maps.at(property->getMeshItemType());
96 if (auto const* p =
97 dynamic_cast<MeshLib::PropertyVector<double>*>(property))
98 {
99 processPropertyVector(*p, id_map, sfc_mesh);
100 vectors_copied++;
101 }
102 else if (auto const* p =
103 dynamic_cast<MeshLib::PropertyVector<float>*>(property))
104 {
105 processPropertyVector(*p, id_map, sfc_mesh);
106 vectors_copied++;
107 }
108 else if (auto const* p =
109 dynamic_cast<MeshLib::PropertyVector<int>*>(property))
110 {
111 processPropertyVector(*p, id_map, sfc_mesh);
112 vectors_copied++;
113 }
114 else if (auto const* p =
115 dynamic_cast<MeshLib::PropertyVector<unsigned>*>(property))
116 {
117 processPropertyVector(*p, id_map, sfc_mesh);
118 vectors_copied++;
119 }
120 else if (auto const* p =
121 dynamic_cast<MeshLib::PropertyVector<long>*>(property))
122 {
123 processPropertyVector(*p, id_map, sfc_mesh);
124 vectors_copied++;
125 }
126 else if (auto const* p =
128 property))
129 {
130 processPropertyVector(*p, id_map, sfc_mesh);
131 vectors_copied++;
132 }
133 else if (auto const* p =
135 property))
136 {
137 processPropertyVector(*p, id_map, sfc_mesh);
138 vectors_copied++;
139 }
140 else if (auto const* p =
142 property))
143 {
144 processPropertyVector(*p, id_map, sfc_mesh);
145 vectors_copied++;
146 }
147 else if (auto const* p =
149 property))
150 {
151 processPropertyVector(*p, id_map, sfc_mesh);
152 vectors_copied++;
153 }
154 else if (auto const* p =
155 dynamic_cast<MeshLib::PropertyVector<char>*>(property))
156 {
157 processPropertyVector(*p, id_map, sfc_mesh);
158 vectors_copied++;
159 }
160 else
161 {
162 WARN(
163 "Skipping property vector '{:s}' - no matching data type "
164 "'{:s}' found.",
165 name, typeid(*property).name());
166 vectors_skipped++;
167 }
168 }
169 INFO("{:d} property vectors copied, {:d} vectors skipped.", vectors_copied,
170 vectors_skipped);
171 return true;
172}
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100
void processPropertyVector(MeshLib::PropertyVector< T > const &property, std::vector< std::size_t > const &id_map, MeshLib::Mesh &sfc_mesh)

References MeshLib::Cell, ERR(), MeshLib::Mesh::getNumberOfElements(), MeshLib::Mesh::getNumberOfNodes(), INFO(), MeshLib::Node, processPropertyVector(), and WARN().

Referenced by MeshToolsLib::BoundaryExtraction::getBoundaryElementsAsMesh(), and MeshToolsLib::MeshSurfaceExtraction::getMeshSurface().

◆ createSurfaceElementsFromElement()

void MeshToolsLib::createSurfaceElementsFromElement ( MeshLib::Element const & surface_element,
std::vector< MeshLib::Element * > & surface_elements,
std::vector< std::size_t > & element_to_bulk_element_id_map,
std::vector< std::size_t > & element_to_bulk_face_id_map )

Definition at line 421 of file MeshSurfaceExtraction.cpp.

426{
427 const unsigned n_faces(surface_element.getNumberOfBoundaries());
428 for (unsigned j = 0; j < n_faces; ++j)
429 {
430 if (surface_element.getNeighbor(j) != nullptr)
431 {
432 continue;
433 }
434
435 surface_elements.push_back(
436 const_cast<MeshLib::Element*>(surface_element.getBoundary(j)));
437 element_to_bulk_face_id_map.push_back(j);
438 element_to_bulk_element_id_map.push_back(surface_element.getID());
439 }
440}

References MeshLib::Element::getBoundary(), MeshLib::Element::getID(), MeshLib::Element::getNeighbor(), and MeshLib::Element::getNumberOfBoundaries().

Referenced by createBoundaryElements().

◆ extrudeElement()

MeshLib::Element * MeshToolsLib::extrudeElement ( std::vector< MeshLib::Node * > const & subsfc_nodes,
MeshLib::Element const & sfc_elem,
MeshLib::PropertyVector< std::size_t > const & sfc_to_subsfc_id_map,
std::map< std::size_t, std::size_t > const & subsfc_sfc_id_map )

Extrudes point, line, triangle or quad elements to its higher dimensional versions, i.e. line, quad, prism, hexahedron.

Parameters
subsfc_nodesthe nodes the elements are based on
sfc_elemthe element of the surface that will be extruded
sfc_to_subsfc_id_maprelation between the surface nodes of the surface element and the ids of the nodes of the subsurface mesh
subsfc_sfc_id_mapmapping of the surface nodes of the current mesh to the surface nodes of the extruded mesh
Returns
extruded element (point -> line, line -> quad, tri -> prism, quad -> hexahedron)

Definition at line 44 of file AddLayerToMesh.cpp.

49{
50 if (sfc_elem.getDimension() > 2)
51 {
52 return nullptr;
53 }
54
55 const unsigned nElemNodes(sfc_elem.getNumberOfBaseNodes());
56 auto new_nodes =
57 std::unique_ptr<MeshLib::Node*[]>{new MeshLib::Node*[2 * nElemNodes]};
58
59 for (unsigned j = 0; j < nElemNodes; ++j)
60 {
61 std::size_t const subsfc_id(
62 sfc_to_subsfc_id_map[sfc_elem.getNode(j)->getID()]);
63 new_nodes[j] = subsfc_nodes[subsfc_id];
64 std::size_t new_idx = (nElemNodes == 2) ? (3 - j) : (nElemNodes + j);
65 new_nodes[new_idx] = subsfc_nodes[subsfc_sfc_id_map.at(subsfc_id)];
66 }
67
68 if (sfc_elem.getGeomType() == MeshLib::MeshElemType::LINE)
69 {
70 return new MeshLib::Quad(new_nodes.release());
71 }
72 if (sfc_elem.getGeomType() == MeshLib::MeshElemType::TRIANGLE)
73 {
74 return new MeshLib::Prism(new_nodes.release());
75 }
76 if (sfc_elem.getGeomType() == MeshLib::MeshElemType::QUAD)
77 {
78 return new MeshLib::Hex(new_nodes.release());
79 }
80
81 return nullptr;
82}
TemplateElement< MeshLib::QuadRule4 > Quad
Definition Quad.h:28
TemplateElement< MeshLib::PrismRule6 > Prism
Definition Prism.h:25
TemplateElement< MeshLib::HexRule8 > Hex
Definition Hex.h:25

References MeshLib::Element::getDimension(), MeshLib::Element::getGeomType(), MathLib::Point3dWithID::getID(), MeshLib::Element::getNode(), MeshLib::Element::getNumberOfBaseNodes(), MeshLib::LINE, MeshLib::QUAD, and MeshLib::TRIANGLE.

Referenced by addLayerToMesh().

◆ getIntegrationPointDataOffsetsOfMeshElements()

std::vector< std::size_t > MeshToolsLib::getIntegrationPointDataOffsetsOfMeshElements ( std::vector< MeshLib::Element * > const & mesh_elements,
MeshLib::PropertyVectorBase const & pv,
MeshLib::Properties const & properties )

Definition at line 98 of file IntegrationPointDataTools.cpp.

102{
103 // For special field data such as OGS_VERSION, IntegrationPointMetaData,
104 // etc., which are not "real" integration points:
105 if (pv.getPropertyName().find("_ip") == std::string::npos)
106 {
107 return {};
108 }
109
110 auto const n_components = pv.getNumberOfGlobalComponents();
111
112 std::vector<std::size_t> element_ip_data_offsets(mesh_elements.size() + 1);
113 std::size_t counter = 0;
114 auto const ip_meta_data =
115 MeshLib::getIntegrationPointMetaData(properties, pv.getPropertyName());
116 for (std::size_t i = 0; i < mesh_elements.size(); i++)
117 {
118 auto const* const element = mesh_elements[i];
119
120 // Assuming that the order of elements in mesh_elements is not touched.
121 element_ip_data_offsets[i] = counter;
122 counter += getNumberOfElementIntegrationPoints(ip_meta_data, *element) *
123 n_components;
124 }
125 element_ip_data_offsets[mesh_elements.size()] = counter;
126
127 return element_ip_data_offsets;
128}
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Properties const &properties, std::string const &name)
int getNumberOfElementIntegrationPoints(MeshLib::IntegrationPointMetaData const &ip_meta_data, MeshLib::Element const &e)

References MeshLib::getIntegrationPointMetaData(), getNumberOfElementIntegrationPoints(), MeshLib::PropertyVectorBase::getNumberOfGlobalComponents(), and MeshLib::PropertyVectorBase::getPropertyName().

Referenced by ApplicationUtils::copyPropertyVector(), createPropertyVector(), and zeroMeshFieldDataByMaterialIDs().

◆ getNumberOfElementIntegrationPoints()

int MeshToolsLib::getNumberOfElementIntegrationPoints ( MeshLib::IntegrationPointMetaData const & ip_meta_data,
MeshLib::Element const & e )

Definition at line 36 of file IntegrationPointDataTools.cpp.

39{
40 switch (e.getCellType())
41 {
43 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Line>(
44 ip_meta_data);
46 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Line3>(
47 ip_meta_data);
49 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Tri>(
50 ip_meta_data);
52 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Tri6>(
53 ip_meta_data);
55 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Quad>(
56 ip_meta_data);
58 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Quad8>(
59 ip_meta_data);
61 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Quad9>(
62 ip_meta_data);
64 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Tet>(
65 ip_meta_data);
67 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Tet10>(
68 ip_meta_data);
70 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Hex>(
71 ip_meta_data);
73 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Hex20>(
74 ip_meta_data);
76 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Prism>(
77 ip_meta_data);
79 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Prism15>(
80 ip_meta_data);
82 return getNumberOfElementIntegrationPointsGeneral<MeshLib::Pyramid>(
83 ip_meta_data);
86 MeshLib::Pyramid13>(ip_meta_data);
89 OGS_FATAL("Mesh element type {:s} is not supported",
90 MeshLib::CellType2String(e.getCellType()));
92 default:
93 OGS_FATAL("Invalid Element Type");
94 }
95 return 0;
96}
int getNumberOfElementIntegrationPointsGeneral(MeshLib::IntegrationPointMetaData const &ip_meta_data)

References MeshLib::CellType2String(), MeshLib::Element::getCellType(), getNumberOfElementIntegrationPointsGeneral(), MeshLib::HEX20, MeshLib::HEX27, MeshLib::HEX8, MeshLib::INVALID, MeshLib::LINE2, MeshLib::LINE3, OGS_FATAL, MeshLib::PRISM15, MeshLib::PRISM18, MeshLib::PRISM6, MeshLib::PYRAMID13, MeshLib::PYRAMID5, MeshLib::QUAD4, MeshLib::QUAD8, MeshLib::QUAD9, MeshLib::TET10, MeshLib::TET4, MeshLib::TRI3, and MeshLib::TRI6.

Referenced by ApplicationUtils::checkFieldPropertyVectorSize(), ApplicationUtils::copyFieldPropertyDataToPartitions(), createPropertyVector(), getIntegrationPointDataOffsetsOfMeshElements(), and ApplicationUtils::setIntegrationPointNumberOfPartition().

◆ getNumberOfElementIntegrationPointsGeneral()

template<typename ElementType >
int MeshToolsLib::getNumberOfElementIntegrationPointsGeneral ( MeshLib::IntegrationPointMetaData const & ip_meta_data)

Definition at line 24 of file IntegrationPointDataTools.cpp.

26{
27 using IntegrationPolicy =
29 using NonGenericIntegrationMethod =
31 NonGenericIntegrationMethod int_met{
32 (unsigned int)ip_meta_data.integration_order};
33 return int_met.getNumberOfPoints();
34}
NumLib::IntegrationGaussLegendreRegular< MeshElement::dimension > IntegrationMethod

References MeshLib::IntegrationPointMetaData::integration_order.

Referenced by getNumberOfElementIntegrationPoints().

◆ lutPrismThirdNode()

unsigned MeshToolsLib::lutPrismThirdNode ( unsigned const id1,
unsigned const id2 )

Lookup-table for returning the third node of bottom or top triangle given the other two.

Definition at line 35 of file MeshRevision.cpp.

36{
37 if ((id1 == 0 && id2 == 1) || (id1 == 1 && id2 == 0))
38 {
39 return 2;
40 }
41 if ((id1 == 1 && id2 == 2) || (id1 == 2 && id2 == 1))
42 {
43 return 0;
44 }
45 if ((id1 == 0 && id2 == 2) || (id1 == 2 && id2 == 0))
46 {
47 return 1;
48 }
49 if ((id1 == 3 && id2 == 4) || (id1 == 4 && id2 == 3))
50 {
51 return 5;
52 }
53 if ((id1 == 4 && id2 == 5) || (id1 == 5 && id2 == 4))
54 {
55 return 3;
56 }
57 if ((id1 == 3 && id2 == 5) || (id1 == 5 && id2 == 3))
58 {
59 return 4;
60 }
61 return std::numeric_limits<unsigned>::max();
62}

Referenced by anonymous_namespace{MeshRevision.cpp}::reducePrism().

◆ markUnusedNodes()

std::vector< bool > MeshToolsLib::markUnusedNodes ( std::vector< MeshLib::Element * > const & elements,
std::vector< MeshLib::Node * > const & nodes )

Marks nodes not used by any of the elements.

Definition at line 101 of file RemoveMeshComponents.cpp.

104{
105 std::vector<bool> unused_nodes(nodes.size(), true);
106 for (auto e : elements)
107 {
108 for (unsigned i = 0; i < e->getNumberOfNodes(); i++)
109 {
110 unused_nodes[getNodeIndex(*e, i)] = false;
111 }
112 }
113
114 return unused_nodes;
115}

Referenced by createMeshFromElements(), and removeNodes().

◆ moveMeshNodes()

template<typename Iterator >
void MeshToolsLib::moveMeshNodes ( Iterator begin,
Iterator end,
Eigen::Vector3d const & displacement )

Function that moves mesh nodes.

The function iterates over all mesh nodes between the begin and the end iterator and moves them using the given displacement.

Parameters
beginbegin iterator
endend iterator
displacementthe displacement to use

Definition at line 33 of file moveMeshNodes.h.

36{
37 std::for_each(begin, end,
38 [&displacement](MathLib::Point3d* node)
39 { node->asEigenVector3d() += displacement; });
40};
Eigen::Vector3d const & asEigenVector3d() const
Definition Point3d.h:63

References MathLib::Point3d::asEigenVector3d().

Referenced by main(), TranslateDataDialog::moveGeometry(), and TranslateDataDialog::moveMesh().

◆ processPropertyVector()

template<typename T >
void MeshToolsLib::processPropertyVector ( MeshLib::PropertyVector< T > const & property,
std::vector< std::size_t > const & id_map,
MeshLib::Mesh & sfc_mesh )

Definition at line 35 of file MeshSurfaceExtraction.cpp.

38{
39 auto const number_of_components = property.getNumberOfGlobalComponents();
40
41 auto sfc_prop = MeshLib::getOrCreateMeshProperty<T>(
42 sfc_mesh, property.getPropertyName(), property.getMeshItemType(),
43 number_of_components);
44 sfc_prop->clear();
45 sfc_prop->reserve(id_map.size());
46
47 for (auto bulk_id : id_map)
48 {
49 std::copy_n(&property.getComponent(bulk_id, 0 /*component_id*/),
50 number_of_components, back_inserter(*sfc_prop));
51 }
52}
MeshItemType getMeshItemType() const
std::string const & getPropertyName() const
PROP_VAL_TYPE & getComponent(std::size_t tuple_index, int component)
Returns the value for the given component stored in the given tuple.

References MeshLib::PropertyVector< PROP_VAL_TYPE >::getComponent(), MeshLib::PropertyVectorBase::getMeshItemType(), and MeshLib::PropertyVectorBase::getPropertyName().

Referenced by createSfcMeshProperties().

◆ removeElements()

MeshLib::Mesh * MeshToolsLib::removeElements ( const MeshLib::Mesh & mesh,
const std::vector< std::size_t > & removed_element_ids,
const std::string & new_mesh_name )

Removes mesh elements and returns a new mesh object. The original mesh is kept unchanged.

Parameters
meshan original mesh whose elements are removed
removed_element_idsa vector of element indices to be removed
new_mesh_namea new mesh name
Returns
a new mesh object

Definition at line 52 of file RemoveMeshComponents.cpp.

56{
57 if (removed_element_ids.empty())
58 {
59 INFO("No elements to remove");
60 return nullptr;
61 }
62
63 INFO("Removing total {:d} elements...", removed_element_ids.size());
64 std::vector<MeshLib::Element*> tmp_elems =
65 details::excludeElementCopy(mesh.getElements(), removed_element_ids);
66 INFO("{:d} elements remain in mesh.", tmp_elems.size());
67
68 // copy node and element objects
69 std::vector<MeshLib::Node*> new_nodes =
71 std::vector<MeshLib::Element*> new_elems =
72 MeshLib::copyElementVector(tmp_elems, new_nodes);
73
74 // delete unused nodes
75 MeshLib::NodeSearch ns(mesh);
76 ns.searchNodesConnectedToOnlyGivenElements(removed_element_ids);
77 auto& removed_node_ids(ns.getSearchedNodeIDs());
78 INFO("Removing total {:d} nodes...", removed_node_ids.size());
79 for (auto nodeid : removed_node_ids)
80 {
81 delete new_nodes[nodeid];
82 new_nodes[nodeid] = nullptr;
83 }
84 new_nodes.erase(std::remove(new_nodes.begin(), new_nodes.end(), nullptr),
85 new_nodes.end());
86
87 if (!new_elems.empty())
88 {
89 MeshLib::Mesh* new_mesh =
90 new MeshLib::Mesh(new_mesh_name, new_nodes, new_elems,
91 true /* compute_element_neighbors */,
93 removed_element_ids, removed_node_ids));
94 return new_mesh;
95 }
96
97 INFO("Current selection removes all elements.");
98 return nullptr;
99}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
Properties & getProperties()
Definition Mesh.h:134
Node search class.
Definition NodeSearch.h:25
Properties excludeCopyProperties(std::vector< std::size_t > const &exclude_elem_ids, std::vector< std::size_t > const &exclude_node_ids) const

References MeshLib::copyElementVector(), MeshLib::copyNodeVector(), MeshLib::Properties::excludeCopyProperties(), MeshToolsLib::details::excludeElementCopy(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::Mesh::getProperties(), MeshLib::NodeSearch::getSearchedNodeIDs(), INFO(), and MeshLib::NodeSearch::searchNodesConnectedToOnlyGivenElements().

Referenced by MeshElementRemovalDialog::accept(), MeshToolsLib::RasterToMesh::convert(), LayeredVolume::createRasterLayers(), MeshToolsLib::MeshGenerators::VoxelFromLayeredMeshes::createVoxelFromLayeredMesh(), extractMatGroup(), MeshToolsLib::MeshGenerator::generateRegularPrismMesh(), main(), MeshToolsLib::MeshGenerator::VoxelGridFromMesh::removeUnusedGridCells(), and writeBoundary().

◆ removeMarkedNodes()

void MeshToolsLib::removeMarkedNodes ( std::vector< bool > const & nodes_to_delete,
std::vector< MeshLib::Node * > & nodes )

Deallocates and removes nodes marked true.

Definition at line 117 of file RemoveMeshComponents.cpp.

119{
120 assert(nodes_to_delete.size() == nodes.size());
121
122 for (std::size_t i = 0; i < nodes.size(); i++)
123 {
124 if (nodes_to_delete[i])
125 {
126 delete nodes[i];
127 nodes[i] = nullptr;
128 }
129 }
130 nodes.erase(remove(begin(nodes), end(nodes), nullptr), end(nodes));
131}

Referenced by createMeshFromElements(), and removeNodes().

◆ removeNodes()

MeshLib::Mesh * MeshToolsLib::removeNodes ( const MeshLib::Mesh & mesh,
const std::vector< std::size_t > & del_nodes_idx,
const std::string & new_mesh_name )

Removes the mesh nodes (and connected elements) given in the nodes-list from the mesh.

Parameters
meshan original mesh whose elements are removed
del_nodes_idxa vector of node indices to be removed
new_mesh_namea new mesh name
Returns
a new mesh object

Definition at line 133 of file RemoveMeshComponents.cpp.

136{
137 if (del_nodes_idx.empty())
138 {
139 return nullptr;
140 }
141
142 // copy node and element objects
143 std::vector<MeshLib::Node*> new_nodes =
145 std::vector<MeshLib::Element*> new_elems =
146 MeshLib::copyElementVector(mesh.getElements(), new_nodes);
147
148 // delete elements
149 MeshLib::ElementSearch es(mesh);
150 es.searchByNodeIDs(del_nodes_idx);
151 auto& removed_element_ids = es.getSearchedElementIDs();
152 for (auto eid : removed_element_ids)
153 {
154 delete new_elems[eid];
155 new_elems[eid] = nullptr;
156 }
157 new_elems.erase(std::remove(new_elems.begin(), new_elems.end(), nullptr),
158 new_elems.end());
159
160 // check unused nodes due to element deletion
161 std::vector<bool> const node_delete_flag =
162 markUnusedNodes(new_elems, new_nodes);
163
164 // delete unused nodes
165 removeMarkedNodes(node_delete_flag, new_nodes);
166
167 if (!new_elems.empty())
168 {
169 MeshLib::Mesh* new_mesh =
170 new MeshLib::Mesh(new_mesh_name, new_nodes, new_elems,
171 true /* compute_element_neighbors */,
173 removed_element_ids, del_nodes_idx));
174 return new_mesh;
175 }
176
177 return nullptr;
178}
Element search class.
std::vector< bool > markUnusedNodes(std::vector< MeshLib::Element * > const &elements, std::vector< MeshLib::Node * > const &nodes)
Marks nodes not used by any of the elements.
void removeMarkedNodes(std::vector< bool > const &nodes_to_delete, std::vector< MeshLib::Node * > &nodes)
Deallocates and removes nodes marked true.

References MeshLib::copyElementVector(), MeshLib::copyNodeVector(), MeshLib::Properties::excludeCopyProperties(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::Mesh::getProperties(), MeshLib::ElementSearch::getSearchedElementIDs(), markUnusedNodes(), removeMarkedNodes(), and MeshLib::ElementSearch::searchByNodeIDs().

Referenced by LayeredMeshGenerator::getMesh().

◆ trackSurface()

static void MeshToolsLib::trackSurface ( MeshLib::Element const *const element,
std::vector< unsigned > & sfc_idx,
unsigned const current_index )
static

Finds all surface elements that can be reached from element. All elements that are found in this way are marked in the global sfc_idx vector using the current_index.

Parameters
elementThe mesh element from which the search is started
sfc_idxThe global index vector notifying to which surface elements belong
current_indexThe index that all elements reachable from element will be assigned in sfc_idx

Definition at line 40 of file MeshValidation.cpp.

43{
44 std::stack<MeshLib::Element const*> elem_stack;
45 elem_stack.push(element);
46 while (!elem_stack.empty())
47 {
48 MeshLib::Element const* const elem = elem_stack.top();
49 elem_stack.pop();
50 sfc_idx[elem->getID()] = current_index;
51 std::size_t const n_neighbors(elem->getNumberOfNeighbors());
52 for (std::size_t i = 0; i < n_neighbors; ++i)
53 {
54 MeshLib::Element const* neighbor(elem->getNeighbor(i));
55 if (neighbor != nullptr && sfc_idx[neighbor->getID()] ==
56 std::numeric_limits<unsigned>::max())
57 {
58 elem_stack.push(neighbor);
59 }
60 }
61 }
62}
virtual unsigned getNumberOfNeighbors() const =0
Get the number of neighbors for this element.
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
virtual const Element * getNeighbor(unsigned i) const =0
Get the specified neighbor.

References MeshLib::Element::getID(), MeshLib::Element::getNeighbor(), and MeshLib::Element::getNumberOfNeighbors().

Referenced by MeshToolsLib::MeshValidation::detectHoles().

◆ zeroMeshFieldDataByMaterialIDs()

void MeshToolsLib::zeroMeshFieldDataByMaterialIDs ( MeshLib::Mesh & mesh,
std::vector< int > const & selected_material_ids )

Definition at line 24 of file ZeroMeshFieldDataByMaterialIDs.cpp.

26{
27 auto const materialIds = materialIDs(mesh);
28 if (!materialIds)
29 {
31 "Mesh contains no int-property vector named 'MaterialIDs' needed "
32 "by the 'zero_mesh_field_data_by_materialIDs'");
33 }
34
35 std::vector<std::size_t> element_ids_for_selected_materials;
36 for (std::size_t i = 0; i < materialIds->size(); ++i)
37 {
38 if (std::find(selected_material_ids.begin(),
39 selected_material_ids.end(),
40 (*materialIds)[i]) != selected_material_ids.end())
41 {
42 element_ids_for_selected_materials.push_back(i);
43 }
44 }
45
46 MeshLib::Properties& properties = mesh.getProperties();
47
48 std::vector<std::size_t> element_ip_data_offsets;
49
50 for (auto const& [name, property] : properties)
51 {
52 if (auto const item_type = property->getMeshItemType();
54 {
55 continue;
56 }
57
58 // For special field data such as OGS_VERSION,
59 // IntegrationPointMetaData,
60 // etc., which are not "real" integration points:
61 if (!property->getPropertyName().ends_with("_ip"))
62 {
63 continue;
64 }
65
66 if (properties.template hasPropertyVector<double>(
68 {
69 auto& pv = *properties.template getPropertyVector<double>(name);
70 const int n_components = pv.getNumberOfGlobalComponents();
71
72 if (element_ip_data_offsets.empty())
73 {
74 // The returned values has already been multiplied with
75 // pv.getNumberOfGlobalComponents()
76 element_ip_data_offsets =
78 mesh.getElements(), pv, properties);
79
80 // element_ip_data_offsets / pv.getNumberOfGlobalComponents()
81 std::transform(element_ip_data_offsets.begin(),
82 element_ip_data_offsets.end(),
83 element_ip_data_offsets.begin(),
84 [n = n_components](double const v)
85 { return v / n; });
86 }
87
88 for (auto const element_id : element_ids_for_selected_materials)
89 {
90 std::fill(
91 pv.begin() +
92 n_components * element_ip_data_offsets[element_id],
93 pv.begin() +
94 n_components * element_ip_data_offsets[element_id + 1],
95 0.0);
96 }
97 }
98 }
99}
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
Definition Properties.h:36
std::vector< std::size_t > getIntegrationPointDataOffsetsOfMeshElements(std::vector< MeshLib::Element * > const &mesh_elements, MeshLib::PropertyVectorBase const &pv, MeshLib::Properties const &properties)

References MeshLib::Mesh::getElements(), getIntegrationPointDataOffsetsOfMeshElements(), MeshLib::Mesh::getProperties(), MeshLib::IntegrationPoint, and OGS_FATAL.

Referenced by anonymous_namespace{ProjectData.cpp}::readMeshes().