preCICE v3.1.2
Loading...
Searching...
No Matches
SerializedMesh.cpp
Go to the documentation of this file.
1#include <algorithm>
2#include <map>
3
6#include "mesh/Mesh.hpp"
7#include "mesh/Vertex.hpp"
8#include "utils/assertion.hpp"
9
11
12void SerializedMesh ::assertValid() const
13{
14 PRECICE_ASSERT(sizes.size() == 5);
15 auto dim = sizes[0];
16 PRECICE_ASSERT(0 < dim && dim <= 3);
17 auto nVertices = sizes[1];
18 PRECICE_ASSERT(0 <= nVertices);
19
20 auto nEdges = sizes[2];
21 auto nTriangles = sizes[3];
22 auto nTetrahedra = sizes[4];
23
24 if (nVertices == 0) {
25 PRECICE_ASSERT(nEdges == 0);
26 PRECICE_ASSERT(nTriangles == 0);
27 PRECICE_ASSERT(nTetrahedra == 0);
28 PRECICE_ASSERT(ids.empty());
29 PRECICE_ASSERT(coords.empty());
30 return;
31 }
32 PRECICE_ASSERT(static_cast<std::size_t>(nVertices * dim) == coords.size());
33 bool hasConnectivity = (nEdges + nTriangles + nTetrahedra) > 0;
34 // Global IDs are allowed to have duplicates as they may not be initialized
35 if (hasConnectivity) {
36 PRECICE_ASSERT(ids.size() == static_cast<std::size_t>(2 * nVertices + 2 * nEdges + 3 * nTriangles + 4 * nTetrahedra));
37 std::set<int> validIDs;
38 for (int vertex = 0; vertex < nVertices; ++vertex) {
39 PRECICE_ASSERT(validIDs.count(ids[2 * vertex + 1]) == 0, "Duplicate IDs");
40 validIDs.insert(ids[2 * vertex + 1]);
41 }
42 for (std::size_t idx = 2 * nVertices; idx < ids.size(); ++idx) {
43 PRECICE_ASSERT(validIDs.count(ids[idx]) == 1, "Unknown ID");
44 }
45 } else {
46 PRECICE_ASSERT(ids.size() == static_cast<std::size_t>(nVertices));
47 }
48}
49
50void SerializedMesh::send(Communication &communication, int rankReceiver)
51{
52 communication.sendRange(sizes, rankReceiver);
53 if (sizes[1] > 0) {
54 communication.sendRange(coords, rankReceiver);
55 communication.sendRange(ids, rankReceiver);
56 }
57}
58
60{
62 sm.sizes = communication.receiveRange(rankSender, AsVectorTag<int>{});
63 PRECICE_ASSERT(sm.sizes.size() == 5);
64 auto nVertices = sm.sizes[1];
65 if (nVertices > 0) {
66 sm.coords = communication.receiveRange(rankSender, AsVectorTag<double>{});
67 sm.ids = communication.receiveRange(rankSender, AsVectorTag<int>{});
68 }
69 sm.assertValid();
70 return sm;
71}
72
74{
75 communication.broadcast(sizes);
76 if (sizes[1] > 0) {
77 communication.broadcast(coords);
78 communication.broadcast(ids);
79 }
80}
81
83{
84 constexpr int broadcasterRank{0};
86 communication.broadcast(sm.sizes, broadcasterRank);
87 PRECICE_ASSERT(sm.sizes.size() == 5);
88 auto nVertices = sm.sizes[1];
89 if (nVertices > 0) {
90 communication.broadcast(sm.coords, broadcasterRank);
91 communication.broadcast(sm.ids, broadcasterRank);
92 }
93 sm.assertValid();
94 return sm;
95}
96
98{
100
101 const auto numberOfVertices = sizes[1];
102 if (numberOfVertices == 0) {
103 return;
104 }
105 const auto numberOfEdges = sizes[2];
106 const auto numberOfTriangles = sizes[3];
107 const auto numberOfTetrahedra = sizes[4];
108
109 const bool hasConnectivity = (numberOfEdges + numberOfTriangles + numberOfTetrahedra) > 0;
110
111 const auto dim = mesh.getDimensions();
112
114 {
115 Eigen::VectorXd coord(dim);
116 for (std::size_t i = 0; i < static_cast<std::size_t>(numberOfVertices); ++i) {
117 std::copy_n(&coords[i * dim], dim, coord.data());
118 mesh::Vertex &v = mesh.createVertex(coord);
119
120 if (hasConnectivity) {
121 v.setGlobalIndex(ids[i * 2]);
122 vertices.emplace(ids[i * 2 + 1], &v);
123 } else {
124 v.setGlobalIndex(ids[i]);
125 }
126 }
127 }
128
129 if (!hasConnectivity) {
130 return;
131 }
132
133 const size_t offsetEdge = numberOfVertices * 2;
134 const size_t offsetTriangle = offsetEdge + 2 * numberOfEdges;
135 const size_t offsetTetrahedron = offsetTriangle + 3 * numberOfTriangles;
136
137 for (size_t idx = offsetEdge; idx != offsetTriangle; idx += 2) {
138 mesh.createEdge(
139 *vertices.at(ids[idx]),
140 *vertices.at(ids[idx + 1]));
141 }
142 for (size_t idx = offsetTriangle; idx != offsetTetrahedron; idx += 3) {
143 mesh.createTriangle(
144 *vertices.at(ids[idx]),
145 *vertices.at(ids[idx + 1]),
146 *vertices.at(ids[idx + 2]));
147 }
148 for (size_t idx = offsetTetrahedron; idx != ids.size(); idx += 4) {
150 *vertices.at(ids[idx]),
151 *vertices.at(ids[idx + 1]),
152 *vertices.at(ids[idx + 2]),
153 *vertices.at(ids[idx + 3]));
154 }
155}
156
158{
159 const auto &meshVertices = mesh.vertices();
160 const auto &meshEdges = mesh.edges();
161 const auto &meshTriangles = mesh.triangles();
162 const auto &meshTetrahedra = mesh.tetrahedra();
163
164 const auto numberOfVertices = meshVertices.size();
165 const auto numberOfEdges = meshEdges.size();
166 const auto numberOfTriangles = meshTriangles.size();
167 const auto numberOfTetrahedra = meshTetrahedra.size();
168
169 SerializedMesh result;
170
171 result.sizes = {
172 mesh.getDimensions(),
173 static_cast<int>(numberOfVertices),
174 static_cast<int>(numberOfEdges),
175 static_cast<int>(numberOfTriangles),
176 static_cast<int>(numberOfTetrahedra)};
177
178 // Empty mesh
179 if (numberOfVertices == 0) {
180 return result;
181 }
182
183 auto dim = static_cast<size_t>(mesh.getDimensions());
184
185 // we always need to send globalIDs
186 auto totalIDs = numberOfVertices;
187 const bool hasConnectivity = mesh.hasConnectivity();
188 if (hasConnectivity) {
189 totalIDs += numberOfVertices // ids for reconstruction, then connectivity
190 + numberOfEdges * 2 // 2 vertices per edge
191 + numberOfTriangles * 3 // 3 vertices per triangle
192 + numberOfTetrahedra * 4; // 4 vertices per tetrahedron
193 }
194 result.ids.reserve(totalIDs);
195
196 result.coords.resize(numberOfVertices * dim);
197 for (size_t i = 0; i < numberOfVertices; ++i) {
198 auto &v = meshVertices[i];
199 std::copy_n(v.rawCoords().begin(), dim, &result.coords[i * dim]);
200 result.ids.push_back(v.getGlobalIndex());
201 // local ids are only interleaved if required
202 if (hasConnectivity) {
203 result.ids.push_back(v.getID());
204 }
205 }
206
207 // Mesh without connectivity information
208 if (!hasConnectivity) {
209 PRECICE_ASSERT(result.ids.size() == numberOfVertices);
210 return result;
211 }
212
213 auto pushVID = [&result](const auto &element, auto... id) {
214 (result.ids.push_back(element.vertex(id).getID()), ...);
215 };
216
217 for (const auto &e : meshEdges) {
218 pushVID(e, 0, 1);
219 }
220
221 for (const auto &e : meshTriangles) {
222 pushVID(e, 0, 1, 2);
223 }
224
225 for (const auto &e : meshTetrahedra) {
226 pushVID(e, 0, 1, 2, 3);
227 }
228
229 result.assertValid();
230
231 // Mesh with connectivity information
232 return result;
233}
234
235} // namespace precice::com::serialize
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
T at(T... args)
Interface for all interprocess communication classes.
std::vector< int > receiveRange(Rank rankSender, AsVectorTag< int >)
Receives a range of ints as a vector<int>
void sendRange(precice::span< const double > itemsToSend, Rank rankReceiver)
Sends a range of doubles (size + content)
virtual void broadcast(precice::span< const int > itemsToSend)
serialized representation of mesh::Mesh
static SerializedMesh serialize(const mesh::Mesh &mesh)
std::vector< double > coords
sizes[0] * dimension coordinates for vertices
void broadcastSend(Communication &communication)
static SerializedMesh receive(Communication &communication, int rankSender)
receives a SerializedMesh and calls assertValid before returning
static SerializedMesh broadcastReceive(Communication &communication)
receives a SerializedMesh and calls assertValid before returning
void send(Communication &communication, int rankReceiver)
void addToMesh(mesh::Mesh &mesh) const
std::vector< int > sizes
contains the dimension, followed by the numbers of vertices, edges, triangles, and tetrahedra
void assertValid() const
asserts the content for correctness
Container and creator for meshes.
Definition Mesh.hpp:39
Triangle & createTriangle(Edge &edgeOne, Edge &edgeTwo, Edge &edgeThree)
Creates and initializes a Triangle object.
Definition Mesh.cpp:119
int getDimensions() const
Definition Mesh.cpp:98
VertexContainer & vertices()
Returns modifieable container holding all vertices.
Definition Mesh.cpp:53
bool hasConnectivity() const
Definition Mesh.hpp:126
TetraContainer & tetrahedra()
Returns modifiable container holding all tetrahedra.
Definition Mesh.cpp:93
Tetrahedron & createTetrahedron(Vertex &vertexOne, Vertex &vertexTwo, Vertex &vertexThree, Vertex &vertexFour)
Creates and initializes a Tetrahedron object.
Definition Mesh.cpp:141
TriangleContainer & triangles()
Returns modifiable container holding all triangles.
Definition Mesh.cpp:78
Edge & createEdge(Vertex &vertexOne, Vertex &vertexTwo)
Creates and initializes an Edge object.
Definition Mesh.cpp:111
EdgeContainer & edges()
Returns modifiable container holding all edges.
Definition Mesh.cpp:68
Vertex & createVertex(const Eigen::VectorXd &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:103
Vertex of a mesh.
Definition Vertex.hpp:16
void setGlobalIndex(int globalIndex)
Definition Vertex.cpp:17
T copy_n(T... args)
T count(T... args)
T data(T... args)
T emplace(T... args)
T insert(T... args)
contains serialization logic
T push_back(T... args)
T reserve(T... args)
T resize(T... args)
T size(T... args)