preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Mesh.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <Eigen/src/Core/Matrix.h>
3#include <algorithm>
4#include <array>
5#include <boost/container/flat_map.hpp>
6#include <functional>
7#include <memory>
8#include <ostream>
9#include <type_traits>
10#include <utility>
11#include <vector>
12
13#include "Edge.hpp"
14#include "Mesh.hpp"
15#include "Tetrahedron.hpp"
16#include "Triangle.hpp"
17#include "logging/LogMacros.hpp"
18#include "math/geometry.hpp"
19#include "mesh/Data.hpp"
20#include "mesh/Vertex.hpp"
22#include "query/Index.hpp"
23#include "utils/assertion.hpp"
24
25namespace precice::mesh {
26
29 int dimensions,
30 MeshID id,
31 bool isJustInTime)
32 : _name(std::move(name)),
33 _dimensions(dimensions),
34 _id(id),
35 _isJustInTime(isJustInTime),
36 _boundingBox(dimensions),
37 _index(*this)
38{
41}
42
44{
46 return _vertices.at(id);
47}
48
50{
52 return _vertices.at(id);
53}
54
59
61{
62 return _vertices;
63}
64
66{
67 return _vertices.size();
68}
69
74
76{
77 return _edges;
78}
79
84
86{
87 return _triangles;
88}
89
91{
92 return _tetrahedra;
93}
94
99
101{
102 return _dimensions;
103}
104
105Vertex &Mesh::createVertex(const Eigen::Ref<const Eigen::VectorXd> &coords)
106{
107 PRECICE_ASSERT(coords.size() == _dimensions, coords.size(), _dimensions);
108 auto nextID = _vertices.size();
109 _vertices.emplace_back(coords, nextID);
110 return _vertices.back();
111}
112
114 Vertex &vertexOne,
115 Vertex &vertexTwo)
116{
117 _edges.emplace_back(vertexOne, vertexTwo);
118 return _edges.back();
119}
120
122 Edge &edgeOne,
123 Edge &edgeTwo,
124 Edge &edgeThree)
125{
127 edgeOne.connectedTo(edgeTwo) &&
128 edgeTwo.connectedTo(edgeThree) &&
129 edgeThree.connectedTo(edgeOne));
130 _triangles.emplace_back(edgeOne, edgeTwo, edgeThree);
131 return _triangles.back();
132}
133
135 Vertex &vertexOne,
136 Vertex &vertexTwo,
137 Vertex &vertexThree)
138{
139 _triangles.emplace_back(vertexOne, vertexTwo, vertexThree);
140 return _triangles.back();
141}
142
144 Vertex &vertexOne,
145 Vertex &vertexTwo,
146 Vertex &vertexThree,
147 Vertex &vertexFour)
148{
149 _tetrahedra.emplace_back(vertexOne, vertexTwo, vertexThree, vertexFour);
150 return _tetrahedra.back();
151}
152
154 const std::string &name,
155 int dimension,
156 DataID id,
157 int waveformDegree)
158{
159 PRECICE_TRACE(name, dimension);
160 for (const PtrData &data : _data) {
161 PRECICE_CHECK(data->getName() != name,
162 "Data \"{}\" cannot be created twice for mesh \"{}\". "
163 "Please rename or remove one of the use-data tags with name \"{}\".",
164 name, _name, name);
165 }
166 // #rows = dimensions of current mesh #columns = dimensions of corresponding data set
167 PtrData data(new Data(name, id, dimension, _dimensions, waveformDegree));
169 return _data.back();
170}
171
173{
174 return _data;
175}
176
177bool Mesh::hasDataID(DataID dataID) const
178{
179 auto iter = std::find_if(_data.begin(), _data.end(), [dataID](const auto &dptr) {
180 return dptr->getID() == dataID;
181 });
182 return iter != _data.end(); // if id was not found in mesh, iter == _data.end()
183}
184
185const PtrData &Mesh::data(DataID dataID) const
186{
187 auto iter = std::find_if(_data.begin(), _data.end(), [dataID](const auto &dptr) {
188 return dptr->getID() == dataID;
189 });
190 PRECICE_ASSERT(iter != _data.end(), "Data with id not found in mesh.", dataID, _name);
191 return *iter;
192}
193
195{
196 auto iter = std::find_if(_data.begin(), _data.end(), [&dataName](const auto &dptr) {
197 return dptr->getName() == dataName;
198 });
199 return iter != _data.end(); // if name was not found in mesh, iter == _data.end()
200}
201
203{
205 for (const auto &data : _data) {
206 names.push_back(data->getName());
207 }
208 return names;
209}
210
211const PtrData &Mesh::data(std::string_view dataName) const
212{
213 auto iter = std::find_if(_data.begin(), _data.end(), [&dataName](const auto &dptr) {
214 return dptr->getName() == dataName;
215 });
216 PRECICE_ASSERT(iter != _data.end(), "Data not found in mesh", dataName, _name);
217 return *iter;
218}
219
221{
222 return _name;
223}
224
226{
227 return _id;
228}
229
230bool Mesh::isValidVertexID(VertexID vertexID) const
231{
232 return (0 <= vertexID) && (static_cast<size_t>(vertexID) < nVertices());
233}
234
236{
238 const auto expectedCount = _vertices.size();
239 for (PtrData &data : _data) {
240 data->allocateValues(expectedCount);
241 }
242}
243
245{
247
248 // Keep the bounding box if set via the API function.
250
251 for (const Vertex &vertex : _vertices) {
252 bb.expandBy(vertex);
253 }
254 _boundingBox = std::move(bb);
255 PRECICE_DEBUG("Bounding Box, {}", _boundingBox);
256}
257
259{
261 _edges.clear();
264 _index.clear();
265
267 for (mesh::PtrData &data : _data) {
268 data->values().resize(0);
269 }
270}
271
281
283{
284 for (mesh::PtrData &data : _data) {
285 data->timeStepsStorage().clear();
286 }
287}
288
290{
291 // Without offset data, we assume non-empty partitions
292 if (_vertexOffsets.empty()) {
293 return false;
294 }
295 PRECICE_ASSERT(_vertexOffsets.size() >= static_cast<std::size_t>(rank));
296 if (rank == 0) {
297 return _vertexOffsets[0] == 0;
298 }
299 return _vertexOffsets[rank] - _vertexOffsets[rank - 1] == 0;
300}
301
302Eigen::VectorXd Mesh::getOwnedVertexData(const Eigen::VectorXd &values)
303{
304 std::vector<double> ownedDataVector;
305 PRECICE_ASSERT(static_cast<std::size_t>(values.size()) >= nVertices());
306 if (empty()) {
307 return {};
308 }
309 int valueDim = values.size() / nVertices();
310 int index = 0;
311
312 for (const auto &vertex : vertices()) {
313 if (vertex.isOwner()) {
314 for (int dim = 0; dim < valueDim; ++dim) {
315 ownedDataVector.push_back(values[index * valueDim + dim]);
316 }
317 }
318 ++index;
319 }
320 Eigen::Map<Eigen::VectorXd> ownedData(ownedDataVector.data(), ownedDataVector.size());
321
322 return ownedData;
323}
324
326{
327 for (auto &vertex : _vertices) {
328 vertex.tag();
329 }
330}
331
333 Mesh &deltaMesh)
334{
337
338 boost::container::flat_map<VertexID, Vertex *> vertexMap;
339 vertexMap.reserve(deltaMesh.nVertices());
340 Eigen::VectorXd coords(_dimensions);
341 for (const Vertex &vertex : deltaMesh.vertices()) {
342 coords = vertex.getCoords();
343 Vertex &v = createVertex(coords);
345 if (vertex.isTagged())
346 v.tag();
349 vertexMap[vertex.getID()] = &v;
350 }
351
352 // you cannot just take the vertices from the edge and add them,
353 // since you need the vertices from the new mesh
354 // (which may differ in IDs)
355 for (const Edge &edge : deltaMesh.edges()) {
356 VertexID vertexIndex1 = edge.vertex(0).getID();
357 VertexID vertexIndex2 = edge.vertex(1).getID();
358 PRECICE_ASSERT((vertexMap.count(vertexIndex1) == 1) &&
359 (vertexMap.count(vertexIndex2) == 1));
360 createEdge(*vertexMap[vertexIndex1], *vertexMap[vertexIndex2]);
361 }
362
363 for (const Triangle &triangle : deltaMesh.triangles()) {
364 VertexID vertexIndex1 = triangle.vertex(0).getID();
365 VertexID vertexIndex2 = triangle.vertex(1).getID();
366 VertexID vertexIndex3 = triangle.vertex(2).getID();
367 PRECICE_ASSERT((vertexMap.count(vertexIndex1) == 1) &&
368 (vertexMap.count(vertexIndex2) == 1) &&
369 (vertexMap.count(vertexIndex3) == 1));
370 createTriangle(*vertexMap[vertexIndex1], *vertexMap[vertexIndex2], *vertexMap[vertexIndex3]);
371 }
372
373 for (const Tetrahedron &tetra : deltaMesh.tetrahedra()) {
374 VertexID vertexIndex1 = tetra.vertex(0).getID();
375 VertexID vertexIndex2 = tetra.vertex(1).getID();
376 VertexID vertexIndex3 = tetra.vertex(2).getID();
377 VertexID vertexIndex4 = tetra.vertex(3).getID();
378
379 PRECICE_ASSERT((vertexMap.count(vertexIndex1) == 1) &&
380 (vertexMap.count(vertexIndex2) == 1) &&
381 (vertexMap.count(vertexIndex3) == 1) &&
382 (vertexMap.count(vertexIndex4) == 1));
383 createTetrahedron(*vertexMap[vertexIndex1], *vertexMap[vertexIndex2], *vertexMap[vertexIndex3], *vertexMap[vertexIndex4]);
384 }
385 _index.clear();
386}
387
389{
390 return _boundingBox;
391}
392
393void Mesh::expandBoundingBox(const BoundingBox &boundingBox)
394{
395 _boundingBox.expandBy(boundingBox);
396}
397
403
405{
406 // Remove duplicate tetrahedra
407 auto tetrahedraCnt = _tetrahedra.size();
409 auto lastTetrahedron = std::unique(_tetrahedra.begin(), _tetrahedra.end());
410 _tetrahedra = TetraContainer{_tetrahedra.begin(), lastTetrahedron};
411
412 // Remove duplicate triangles
413 auto triangleCnt = _triangles.size();
415 auto lastTriangle = std::unique(_triangles.begin(), _triangles.end());
416 _triangles = TriangleContainer{_triangles.begin(), lastTriangle};
417
418 // Remove duplicate edges
419 auto edgeCnt = _edges.size();
421 auto lastEdge = std::unique(_edges.begin(), _edges.end());
422 _edges = EdgeContainer{_edges.begin(), lastEdge};
423
424 PRECICE_DEBUG("Compression removed {} tetrahedra ({} to {}), {} triangles ({} to {}), and {} edges ({} to {})",
425 tetrahedraCnt - _tetrahedra.size(), tetrahedraCnt, _tetrahedra.size(),
426 triangleCnt - _triangles.size(), triangleCnt, _triangles.size(),
427 edgeCnt - _edges.size(), edgeCnt, _edges.size());
428}
429
430namespace {
431
432template <class Primitive, int... Indices>
433auto sortedVertexPtrsForImpl(Primitive &p, std::integer_sequence<int, Indices...>)
434{
435 std::array<Vertex *, Primitive::vertexCount> vs{&p.vertex(Indices)...};
436 std::sort(vs.begin(), vs.end());
437 return std::tuple_cat(vs);
438}
439
448template <class Primitive>
449auto sortedVertexPtrsFor(Primitive &p)
450{
451 return sortedVertexPtrsForImpl(p, std::make_integer_sequence<int, Primitive::vertexCount>{});
452}
453} // namespace
454
456{
457 if (_triangles.empty() && _tetrahedra.empty()) {
458 PRECICE_DEBUG("No implicit primitives required");
459 return; // no implicit primitives
460 }
461
462 // count explicit primitives for debug
463 const auto explTriangles = _triangles.size();
464 const auto explEdges = _triangles.size();
465
466 // First handle all explicit tetrahedra
467
468 // Build a set of all explicit triangles
469 using ExisitingTriangle = std::tuple<Vertex *, Vertex *, Vertex *>;
471 for (auto &t : _triangles) {
472 triangles.insert(sortedVertexPtrsFor(t));
473 }
474
475 // Generate all missing implicit triangles of explicit tetrahedra
476 // Update the triangles set used by the implicit edge generation
477 auto createTriangleIfMissing = [&](Vertex *a, Vertex *b, Vertex *c) {
478 if (triangles.count({a, b, c}) == 0) {
479 triangles.emplace(a, b, c);
480 createTriangle(*a, *b, *c);
481 };
482 };
483 for (auto &t : _tetrahedra) {
484 auto [a, b, c, d] = sortedVertexPtrsFor(t);
485 // Make sure these are in the same order as above
486 createTriangleIfMissing(a, b, c);
487 createTriangleIfMissing(a, b, d);
488 createTriangleIfMissing(a, c, d);
489 createTriangleIfMissing(b, c, d);
490 }
491
492 // Second handle all triangles, both explicit and implicit from the tetrahedron phase
493 // Build an set of all explicit triangles
494 using ExisitingEdge = std::tuple<Vertex *, Vertex *>;
496 for (auto &e : _edges) {
497 edges.emplace(sortedVertexPtrsFor(e));
498 }
499
500 // generate all missing implicit edges of implicit and explicit triangles
501 auto createEdgeIfMissing = [&](Vertex *a, Vertex *b) {
502 if (edges.count({a, b}) == 0) {
503 edges.emplace(a, b);
504 createEdge(*a, *b);
505 };
506 };
507 for (auto &t : _triangles) {
508 auto [a, b, c] = sortedVertexPtrsFor(t);
509 // Make sure these are in the same order as above
510 createEdgeIfMissing(a, b);
511 createEdgeIfMissing(a, c);
512 createEdgeIfMissing(b, c);
513 }
514
515 PRECICE_DEBUG("Generated {} implicit triangles and {} implicit edges",
516 _triangles.size() - explTriangles,
517 _edges.size() - explEdges);
518}
519
520bool Mesh::operator==(const Mesh &other) const
521{
522 bool equal = true;
523 equal &= _vertices.size() == other._vertices.size() &&
524 std::is_permutation(_vertices.begin(), _vertices.end(), other._vertices.begin());
525 equal &= _edges.size() == other._edges.size() &&
526 std::is_permutation(_edges.begin(), _edges.end(), other._edges.begin());
527 equal &= _triangles.size() == other._triangles.size() &&
528 std::is_permutation(_triangles.begin(), _triangles.end(), other._triangles.begin());
529 return equal;
530}
531
532bool Mesh::operator!=(const Mesh &other) const
533{
534 return !(*this == other);
535}
536
537std::ostream &operator<<(std::ostream &os, const Mesh &m)
538{
539 os << "Mesh \"" << m.getName() << "\", dimensionality = " << m.getDimensions() << ":\n";
540 os << "GEOMETRYCOLLECTION(\n";
541 const auto token = ", ";
542 const auto *sep = "";
543 for (auto &vertex : m.vertices()) {
544 os << sep << vertex;
545 sep = token;
546 }
547 sep = ",\n";
548 for (auto &edge : m.edges()) {
549 os << sep << edge;
550 sep = token;
551 }
552 sep = ",\n";
553 for (auto &triangle : m.triangles()) {
554 os << sep << triangle;
555 sep = token;
556 }
557 os << "\n)";
558 return os;
559}
560
561} // namespace precice::mesh
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
std::string name
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T at(T... args)
T back(T... args)
T begin(T... args)
An axis-aligned bounding box around a (partition of a) mesh.
void expandBy(const BoundingBox &otherBB)
Expand bounding box using another bounding box.
Describes a set of data values belonging to the vertices of a mesh.
Definition Data.hpp:27
Linear edge of a mesh, defined by two Vertex objects.
Definition Edge.hpp:15
bool connectedTo(const Edge &other) const
Checks whether both edges share a vertex.
Definition Edge.cpp:39
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
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
std::string _name
Name of the mesh.
Definition Mesh.hpp:351
BoundingBox _boundingBox
Definition Mesh.hpp:404
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
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
void addMesh(Mesh &deltaMesh)
Definition Mesh.cpp:332
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
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:65
VertexDistribution _vertexDistribution
Vertex distribution for the primary rank, holding for each secondary rank all vertex IDs it owns.
Definition Mesh.hpp:373
TetraContainer & tetrahedra()
Returns modifiable container holding all tetrahedra.
Definition Mesh.cpp:95
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
TetraContainer _tetrahedra
Definition Mesh.hpp:363
Vertex & vertex(VertexID id)
Mutable access to a vertex by VertexID.
Definition Mesh.cpp:43
Mesh(std::string name, int dimensions, MeshID id, bool isJustInTime=false)
Constructor.
Definition Mesh.cpp:27
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
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
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
bool isValidVertexID(VertexID vertexID) const
Returns true if the given vertexID is valid.
Definition Mesh.cpp:230
EdgeContainer _edges
Definition Mesh.hpp:361
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
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
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
VertexID getID() const
Returns the unique (among vertices of one mesh on one processor) ID of the vertex.
Definition Vertex.hpp:109
bool isTagged() const
Definition Vertex.cpp:32
void setGlobalIndex(int globalIndex)
Definition Vertex.cpp:17
void setOwner(bool owner)
Definition Vertex.cpp:27
bool isOwner() const
Definition Vertex.cpp:22
Eigen::VectorXd getCoords() const
Returns the coordinates of the vertex.
Definition Vertex.hpp:114
int getGlobalIndex() const
Globally unique index.
Definition Vertex.cpp:12
void clear()
Clear the index.
Definition Index.cpp:391
T clear(T... args)
T count(T... args)
T data(T... args)
T emplace_back(T... args)
T emplace(T... args)
T empty(T... args)
T end(T... args)
T find_if(T... args)
T insert(T... args)
T is_permutation(T... args)
provides Mesh, Data and primitives.
int MeshID
Definition Types.hpp:30
int VertexID
Definition Types.hpp:13
int Rank
Definition Types.hpp:37
int DataID
Definition Types.hpp:25
STL namespace.
T push_back(T... args)
T resize(T... args)
T size(T... args)
T sort(T... args)
T tuple_cat(T... args)
T unique(T... args)