OGS
Mesh.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <cstdlib>
7#include <memory>
8#include <range/v3/view/iota.hpp>
9#include <range/v3/view/transform.hpp>
10#include <string>
11#include <vector>
12
13#include "BaseLib/Algorithm.h"
14#include "BaseLib/Error.h"
15#include "Location.h"
16#include "MathLib/Point3d.h"
17#include "MeshEnums.h"
18#include "Properties.h"
19
20namespace ApplicationUtils
21{
23}
24
25namespace MeshLib
26{
27class Node;
28class Element;
29
33class Mesh
34{
35 /* friend functions: */
36 friend void removeMeshNodes(Mesh& mesh,
37 const std::vector<std::size_t>& nodes);
38
40
41public:
49 Mesh(std::string name,
50 std::vector<Node*>
51 nodes,
52 std::vector<Element*>
53 elements,
54 bool const compute_element_neighbors = false,
55 Properties const& properties = Properties());
56
58 Mesh(const Mesh& mesh);
59
60 Mesh(Mesh&& mesh);
61
62 Mesh& operator=(const Mesh& mesh) = delete;
63 Mesh& operator=(Mesh&& mesh) = delete;
64
66 virtual ~Mesh();
67
73 void shallowClean();
74
76 void addElement(Element* elem);
77
80 unsigned getDimension() const { return _mesh_dimension; }
81
83 const Node* getNode(std::size_t idx) const { return _nodes[idx]; }
84
86 const Element* getElement(std::size_t idx) const { return _elements[idx]; }
87
89 std::size_t getNumberOfElements() const { return _elements.size(); }
90
92 std::size_t getNumberOfNodes() const { return _nodes.size(); }
93
95 const std::string getName() const { return _name; }
96
98 std::vector<Node*> const& getNodes() const { return _nodes; }
99
101 std::vector<Element*> const& getElements() const { return _elements; }
102
105 void resetElementIDs();
106
108 void resetNodeIDs();
109
111 void setName(const std::string& name) { this->_name = name; }
112
114 std::size_t getID() const { return _id; }
115
117 std::size_t computeNumberOfBaseNodes() const;
118
120 bool hasNonlinearElement() const;
121
122 std::vector<Element const*> const& getElementsConnectedToNode(
123 std::size_t node_id) const;
124 std::vector<Element const*> const& getElementsConnectedToNode(
125 Node const& node) const;
126
128 Properties const& getProperties() const { return _properties; }
129
131 void setAxiallySymmetric(bool is_axial_symmetric)
132 {
133 _is_axially_symmetric = is_axial_symmetric;
134 }
135
136protected:
138 void setDimension();
139
143 void setElementNeighbors();
144
145 std::size_t const _id;
149 std::pair<double, double> _node_distance;
150 std::string _name;
151 std::vector<Node*> _nodes;
152 std::vector<Element*> _elements;
154
155 std::vector<std::vector<Element const*>> _elements_connected_to_nodes;
156
159}; /* class */
160
163std::vector<std::vector<Node*>> calculateNodesConnectedByElements(
164 Mesh const& mesh);
165
167inline bool operator==(Mesh const& a, Mesh const& b)
168{
169 return a.getID() == b.getID();
170}
171
172inline bool operator!=(Mesh const& a, Mesh const& b)
173{
174 return !(a == b);
175}
176
182PropertyVector<int> const* materialIDs(Mesh const& mesh);
183PropertyVector<int>* materialIDs(Mesh& mesh);
184PropertyVector<std::size_t> const* bulkNodeIDs(Mesh const& mesh);
185PropertyVector<std::size_t> const* bulkElementIDs(Mesh const& mesh);
186
189bool isBaseNode(Node const& node,
190 std::vector<Element const*> const& elements_connected_to_node);
191
193std::pair<double, double> minMaxEdgeLength(
194 std::vector<Element*> const& elements);
195
198template <typename T>
199bool idsComparator(T const a, T const b)
200{
201 if constexpr (std::is_pointer_v<T>)
202 {
203 return a->getID() < b->getID();
204 }
205 else
206 {
207 return a.getID() < b.getID();
208 }
209}
210
211Mesh& findMeshByName(std::vector<std::unique_ptr<Mesh>> const& meshes,
212 std::string_view const name);
213
215namespace views
216{
218inline constexpr ranges::views::view_closure ids =
219 ranges::views::transform([](auto const& a) { return a->getID(); });
220
222inline constexpr ranges::views::view_closure names =
223 ranges::views::transform([](auto const& a) { return a->getName(); });
224
225inline constexpr ranges::views::view_closure coords =
226 ranges::views::transform([](MathLib::Point3d const* n)
227 { return std::span(n->data(), n->data() + 3); });
228
229inline auto meshLocations(Mesh const& mesh, MeshItemType const item_type)
230{
231 return ranges::views::iota(std::size_t{0}, mesh.getNumberOfNodes()) |
232 ranges::views::transform(
233 [mesh_id = mesh.getID(), item_type](std::size_t const n)
234 { return Location{mesh_id, item_type, n}; });
235}
236} // namespace views
237} // namespace MeshLib
const double * data() const
Definition Point3d.h:51
bool isAxiallySymmetric() const
Definition Mesh.h:130
Properties _properties
Definition Mesh.h:153
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:98
std::size_t const _id
Definition Mesh.h:145
bool const _compute_element_neighbors
Definition Mesh.h:158
Mesh(Mesh &&mesh)
std::vector< std::vector< Element const * > > _elements_connected_to_nodes
Definition Mesh.h:155
void setName(const std::string &name)
Changes the name of the mesh.
Definition Mesh.h:111
void addElement(Element *elem)
Add an element to the mesh.
Definition Mesh.cpp:144
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:101
unsigned _mesh_dimension
Definition Mesh.h:146
std::size_t computeNumberOfBaseNodes() const
Get the number of base nodes.
Definition Mesh.cpp:228
friend void removeMeshNodes(Mesh &mesh, const std::vector< std::size_t > &nodes)
std::string _name
Definition Mesh.h:150
Mesh & operator=(Mesh &&mesh)=delete
unsigned getDimension() const
Definition Mesh.h:80
std::size_t getID() const
Get id of the mesh.
Definition Mesh.h:114
void resetNodeIDs()
Resets the IDs of all mesh-nodes to their position in the node vector.
Definition Mesh.cpp:149
const Node * getNode(std::size_t idx) const
Get the node with the given index.
Definition Mesh.h:83
Properties & getProperties()
Definition Mesh.h:127
const Element * getElement(std::size_t idx) const
Get the element with the given index.
Definition Mesh.h:86
std::vector< Element * > _elements
Definition Mesh.h:152
void setDimension()
Sets the dimension of the mesh.
Definition Mesh.cpp:167
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:95
void resetElementIDs()
Definition Mesh.cpp:158
void setElementNeighbors()
Definition Mesh.cpp:194
virtual ~Mesh()
Destructor.
Definition Mesh.cpp:129
Mesh(std::string name, std::vector< Node * > nodes, std::vector< Element * > elements, bool const compute_element_neighbors=false, Properties const &properties=Properties())
Definition Mesh.cpp:52
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:92
std::vector< Element const * > const & getElementsConnectedToNode(std::size_t node_id) const
Definition Mesh.cpp:246
void setAxiallySymmetric(bool is_axial_symmetric)
Definition Mesh.h:131
std::pair< double, double > _node_distance
Definition Mesh.h:149
Mesh & operator=(const Mesh &mesh)=delete
void shallowClean()
Definition Mesh.cpp:123
Properties const & getProperties() const
Definition Mesh.h:128
bool hasNonlinearElement() const
Check if the mesh contains any nonlinear element.
Definition Mesh.cpp:238
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:89
std::vector< Node * > _nodes
Definition Mesh.h:151
bool _is_axially_symmetric
Definition Mesh.h:157
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
MeshLib specific, lazy, non-owning, non-mutating, composable range views.
Definition Mesh.h:216
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:218
auto meshLocations(Mesh const &mesh, MeshItemType const item_type)
Definition Mesh.h:229
constexpr ranges::views::view_closure coords
Definition Mesh.h:225
constexpr ranges::views::view_closure names
For an element of a range view return its name.
Definition Mesh.h:222
std::vector< std::vector< Node * > > calculateNodesConnectedByElements(Mesh const &mesh)
Definition Mesh.cpp:298
Mesh & findMeshByName(std::vector< std::unique_ptr< Mesh > > const &meshes, std::string_view const name)
Definition Mesh.cpp:354
bool idsComparator(T const a, T const b)
Definition Mesh.h:199
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition Mesh.cpp:258
PropertyVector< std::size_t > const * bulkElementIDs(Mesh const &mesh)
Definition Mesh.cpp:290
std::pair< double, double > minMaxEdgeLength(std::vector< Element * > const &elements)
Returns the minimum and maximum edge length for given elements.
Definition Mesh.cpp:179
bool operator==(Mesh const &a, Mesh const &b)
Meshes are equal if their id's are equal.
Definition Mesh.h:167
bool operator!=(Mesh const &a, Mesh const &b)
Definition Mesh.h:172
bool isBaseNode(Node const &node, std::vector< Element const * > const &elements_connected_to_node)
Definition Mesh.cpp:336
PropertyVector< std::size_t > const * bulkNodeIDs(Mesh const &mesh)
Definition Mesh.cpp:282