preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Mesh.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <deque>
5#include <iosfwd>
6#include <list>
7#include <map>
8#include <string>
9#include <string_view>
10#include <vector>
11
12#include "logging/Logger.hpp"
13#include "mesh/BoundingBox.hpp"
14#include "mesh/Data.hpp"
15#include "mesh/Edge.hpp"
17#include "mesh/Tetrahedron.hpp"
18#include "mesh/Triangle.hpp"
19#include "mesh/Vertex.hpp"
21#include "query/Index.hpp"
23#include "utils/assertion.hpp"
24
25namespace precice::mesh {
26
38class Mesh {
39public:
46
49
52 using ConnectionMap = CommunicationMap; // until we decide on a name
53
55
57 static constexpr MeshID MESH_ID_UNDEFINED{-1};
58
66 Mesh(
68 int dimensions,
69 MeshID id,
70 bool isJustInTime = false);
71
74
76 const Vertex &vertex(VertexID id) const;
77
80
82 const VertexContainer &vertices() const;
83
85 std::size_t nVertices() const;
86
88 bool empty() const
89 {
90 return _vertices.empty();
91 }
92
95
97 const EdgeContainer &edges() const;
98
99 bool hasEdges() const
100 {
101 return !_edges.empty();
102 }
103
106
108 const TriangleContainer &triangles() const;
109
110 bool hasTriangles() const
111 {
112 return !_triangles.empty();
113 }
114
117
119 const TetraContainer &tetrahedra() const;
120
121 bool hasTetrahedra() const
122 {
123 return !_tetrahedra.empty();
124 }
125
126 bool hasConnectivity() const
127 {
128 return hasEdges() || hasTriangles() || hasTetrahedra();
129 }
130
131 int getDimensions() const;
132
134 Vertex &createVertex(const Eigen::Ref<const Eigen::VectorXd> &coords);
135
143 Vertex &vertexOne,
144 Vertex &vertexTwo);
145
154 Edge &edgeOne,
155 Edge &edgeTwo,
156 Edge &edgeThree);
157
166 Vertex &vertexOne,
167 Vertex &vertexTwo,
168 Vertex &vertexThree);
169
179 Vertex &vertexOne,
180 Vertex &vertexTwo,
181 Vertex &vertexThree,
182 Vertex &vertexFour);
183
186 int dimension,
187 DataID id,
188 int waveformDegree = time::Time::DEFAULT_WAVEFORM_DEGREE);
189
191 const DataContainer &data() const;
192
194 bool hasDataID(DataID dataID) const;
195
197 const PtrData &data(DataID dataID) const;
198
200 bool hasDataName(std::string_view dataName) const;
201
204
206 const PtrData &data(std::string_view dataName) const;
207
209 const std::string &getName() const;
210
212 MeshID getID() const;
213
215 bool isValidVertexID(VertexID vertexID) const;
216
218 void allocateDataValues(); //@todo Redesign mapping and remove this function. See https://github.com/precice/precice/issues/1651.
219
221 void computeBoundingBox();
222
231 void clear();
232
234 void clearPartitioning();
235
237 void clearDataStamples();
238
240 {
241 PRECICE_ASSERT(std::all_of(vd.begin(), vd.end(), [](const auto &p) { return std::is_sorted(p.second.begin(), p.second.end()); }));
242 _vertexDistribution = std::move(vd);
243 }
244
247 {
248 return _vertexDistribution;
249 }
250
252 {
253 return _vertexOffsets;
254 }
255
257 bool isPartitionEmpty(Rank rank) const;
258
260 void setVertexOffsets(VertexOffsets vertexOffsets)
261 {
262 _vertexOffsets = std::move(vertexOffsets);
263 }
264
266 {
268 }
269
271 {
273 }
274
275 // Get the data of owned vertices for given data ID
276 Eigen::VectorXd getOwnedVertexData(const Eigen::VectorXd &values);
277
278 // Tag all the vertices
279 void tagAll();
280
283 {
284 return _connectedRanks;
285 }
286
289 {
290 _connectedRanks = std::move(ranks);
291 }
292
298
299 void addMesh(Mesh &deltaMesh);
300
312 const BoundingBox &getBoundingBox() const;
313
314 void expandBoundingBox(const BoundingBox &bounding_box);
315
316 bool operator==(const Mesh &other) const;
317
318 bool operator!=(const Mesh &other) const;
319
321 const query::Index &index() const
322 {
323 return _index;
324 }
325
328 {
329 return _index;
330 }
331
332 bool isJustInTime() const
333 {
334 return _isJustInTime;
335 }
336
345 void preprocess();
346
347private:
348 mutable logging::Logger _log{"mesh::Mesh"};
349
352
355
358
364
367
374
376
381
388
394
400
402 bool _isJustInTime = false;
403
405
407
409 void removeDuplicates();
410
416};
417
419
420} // namespace precice::mesh
std::string name
T all_of(T... args)
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T begin(T... args)
This class provides a lightweight logger.
Definition Logger.hpp:17
An axis-aligned bounding box around a (partition of a) mesh.
Linear edge of a mesh, defined by two Vertex objects.
Definition Edge.hpp:15
Container and creator for meshes.
Definition Mesh.hpp:38
void expandBoundingBox(const BoundingBox &bounding_box)
Definition Mesh.cpp:393
Triangle & createTriangle(Edge &edgeOne, Edge &edgeTwo, Edge &edgeThree)
Creates and initializes a Triangle object.
Definition Mesh.cpp:121
void setGlobalNumberOfVertices(int num)
Definition Mesh.hpp:270
const VertexDistribution & getVertexDistribution() const
Returns a mapping from rank to used (not necessarily owned) vertex IDs.
Definition Mesh.hpp:246
MeshID _id
The ID of this mesh.
Definition Mesh.hpp:357
MeshID getID() const
Returns the base ID of the mesh.
Definition Mesh.cpp:225
int getGlobalNumberOfVertices() const
Definition Mesh.hpp:265
std::string _name
Name of the mesh.
Definition Mesh.hpp:351
std::map< Rank, std::vector< VertexID > > CommunicationMap
A mapping from remote local ranks to the IDs that must be communicated.
Definition Mesh.hpp:51
BoundingBox _boundingBox
Definition Mesh.hpp:404
std::deque< Triangle > TriangleContainer
Definition Mesh.hpp:42
bool hasTetrahedra() const
Definition Mesh.hpp:121
void setVertexDistribution(VertexDistribution vd)
Definition Mesh.hpp:239
int _globalNumberOfVertices
Number of unique vertices for complete distributed mesh.
Definition Mesh.hpp:387
int getDimensions() const
Definition Mesh.cpp:100
VertexContainer & vertices()
Returns modifieable container holding all vertices.
Definition Mesh.cpp:55
void clearDataStamples()
Clears all data stamples.
Definition Mesh.cpp:282
std::vector< PtrData > DataContainer
Definition Mesh.hpp:44
bool hasDataID(DataID dataID) const
Returns whether Mesh has Data with the matchingID.
Definition Mesh.cpp:177
Eigen::VectorXd getOwnedVertexData(const Eigen::VectorXd &values)
Definition Mesh.cpp:302
void clear()
Removes all mesh elements and data values (does not remove data or the bounding boxes).
Definition Mesh.cpp:258
DataContainer _data
Data hold by the vertices of the mesh.
Definition Mesh.hpp:366
bool hasEdges() const
Definition Mesh.hpp:99
void addMesh(Mesh &deltaMesh)
Definition Mesh.cpp:332
bool operator!=(const Mesh &other) const
Definition Mesh.cpp:532
std::vector< std::string > availableData() const
Returns the names of all available data.
Definition Mesh.cpp:202
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Definition Mesh.cpp:220
void removeDuplicates()
Removes all duplicate connectivity.
Definition Mesh.cpp:404
static constexpr MeshID MESH_ID_UNDEFINED
Use if the id of the mesh is not necessary.
Definition Mesh.hpp:57
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:65
bool hasConnectivity() const
Definition Mesh.hpp:126
VertexDistribution _vertexDistribution
Vertex distribution for the primary rank, holding for each secondary rank all vertex IDs it owns.
Definition Mesh.hpp:373
CommunicationMap & getCommunicationMap()
Returns a mapping from remote local connected ranks to the corresponding vertex IDs.
Definition Mesh.hpp:294
TetraContainer & tetrahedra()
Returns modifiable container holding all tetrahedra.
Definition Mesh.cpp:95
std::deque< Tetrahedron > TetraContainer
Definition Mesh.hpp:43
bool operator==(const Mesh &other) const
Definition Mesh.cpp:520
logging::Logger _log
Definition Mesh.hpp:348
std::vector< Rank > _connectedRanks
each rank stores list of connected remote ranks. In the m2n package, this is used to create the initi...
Definition Mesh.hpp:393
bool _isJustInTime
for just-in-time mapping, we need an artificial mesh, which we can use
Definition Mesh.hpp:402
TetraContainer _tetrahedra
Definition Mesh.hpp:363
Vertex & vertex(VertexID id)
Mutable access to a vertex by VertexID.
Definition Mesh.cpp:43
bool isJustInTime() const
Definition Mesh.hpp:332
Mesh(std::string name, int dimensions, MeshID id, bool isJustInTime=false)
Constructor.
Definition Mesh.cpp:27
const query::Index & index() const
Call preprocess() before index() to ensure correct projection handling.
Definition Mesh.hpp:321
query::Index _index
Definition Mesh.hpp:406
VertexContainer _vertices
Holds vertices, edges, triangles and tetrahedra.
Definition Mesh.hpp:360
bool hasDataName(std::string_view dataName) const
Returns whether Mesh has Data with the dataName.
Definition Mesh.cpp:194
std::deque< Edge > EdgeContainer
Definition Mesh.hpp:41
PtrData & createData(const std::string &name, int dimension, DataID id, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
Create only data for vertex.
Definition Mesh.cpp:153
CommunicationMap _communicationMap
each rank stores list of connected ranks and corresponding vertex IDs here. In the m2n package,...
Definition Mesh.hpp:399
bool empty() const
Does the mesh contain any vertices?
Definition Mesh.hpp:88
void generateImplictPrimitives()
Definition Mesh.cpp:455
const std::vector< Rank > & getConnectedRanks() const
Returns a vector of connected ranks.
Definition Mesh.hpp:282
const VertexOffsets & getVertexOffsets() const
Definition Mesh.hpp:251
void clearPartitioning()
Clears the partitioning information.
Definition Mesh.cpp:273
void computeBoundingBox()
Computes the boundingBox for the vertices.
Definition Mesh.cpp:244
TriangleContainer _triangles
Definition Mesh.hpp:362
bool isPartitionEmpty(Rank rank) const
checks if the given ranks partition is empty
Definition Mesh.cpp:289
Tetrahedron & createTetrahedron(Vertex &vertexOne, Vertex &vertexTwo, Vertex &vertexThree, Vertex &vertexFour)
Creates and initializes a Tetrahedron object.
Definition Mesh.cpp:143
VertexOffsets _vertexOffsets
Holds the index of the last vertex for each rank.
Definition Mesh.hpp:380
void setVertexOffsets(VertexOffsets vertexOffsets)
Only used for tests.
Definition Mesh.hpp:260
bool isValidVertexID(VertexID vertexID) const
Returns true if the given vertexID is valid.
Definition Mesh.cpp:230
EdgeContainer _edges
Definition Mesh.hpp:361
void setConnectedRanks(std::vector< Rank > ranks)
Returns a vector of connected ranks.
Definition Mesh.hpp:288
const DataContainer & data() const
Allows access to all data.
Definition Mesh.cpp:172
TriangleContainer & triangles()
Returns modifiable container holding all triangles.
Definition Mesh.cpp:80
const BoundingBox & getBoundingBox() const
Returns the bounding box of the mesh.
Definition Mesh.cpp:388
query::Index & index()
Call preprocess() before index() to ensure correct projection handling.
Definition Mesh.hpp:327
std::deque< Vertex > VertexContainer
Definition Mesh.hpp:40
Edge & createEdge(Vertex &vertexOne, Vertex &vertexTwo)
Creates and initializes an Edge object.
Definition Mesh.cpp:113
Vertex & createVertex(const Eigen::Ref< const Eigen::VectorXd > &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:105
EdgeContainer & edges()
Returns modifiable container holding all edges.
Definition Mesh.cpp:70
void allocateDataValues()
Allocates memory for the vertex data values and corresponding gradient values.
Definition Mesh.cpp:235
bool hasTriangles() const
Definition Mesh.hpp:110
int _dimensions
Dimension of mesh.
Definition Mesh.hpp:354
Tetrahedron of a mesh, defined by 4 vertices.
Triangle of a mesh, defined by three vertices.
Definition Triangle.hpp:24
Vertex of a mesh.
Definition Vertex.hpp:16
Class to query the index trees of the mesh.
Definition Index.hpp:64
static const int DEFAULT_WAVEFORM_DEGREE
To be used, when the interpolation degree is not defined.
Definition Time.hpp:8
T empty(T... args)
T end(T... args)
provides Mesh, Data and primitives.
std::ostream & operator<<(std::ostream &os, const BoundingBox &bb)
int MeshID
Definition Types.hpp:30
int VertexID
Definition Types.hpp:13
int Rank
Definition Types.hpp:37
int DataID
Definition Types.hpp:25