preCICE v3.1.2
Loading...
Searching...
No Matches
Participant.cpp
Go to the documentation of this file.
1#include <optional>
2#include <string_view>
3
7#include "precice/impl/versions.hpp"
8
9namespace precice {
10
11namespace {
13{
14 // The given string_view may contain null chars.
15 // We trim to the first null char here.
16 std::string_view s{sv.data(), sv.size()};
17 auto trim_pos = s.find('\0');
18 if (trim_pos != s.npos) {
19 s.remove_suffix(s.size() - trim_pos);
20 }
21 return s;
22}
23} // namespace
24
26 ::precice::string_view participantName,
27 ::precice::string_view configurationFileName,
28 int solverProcessIndex,
29 int solverProcessSize)
30 : _impl(new impl::ParticipantImpl(toSV(participantName), toSV(configurationFileName), solverProcessIndex, solverProcessSize, std::nullopt))
31{
32}
33
35 ::precice::string_view participantName,
36 ::precice::string_view configurationFileName,
37 int solverProcessIndex,
38 int solverProcessSize,
39 void * communicator)
40 : _impl(new impl::ParticipantImpl(toSV(participantName), toSV(configurationFileName), solverProcessIndex, solverProcessSize, {communicator}))
41{
42}
43
45
47{
48 _impl->initialize();
49}
50
52 double computedTimeStepSize)
53{
54 _impl->advance(computedTimeStepSize);
55}
56
58{
59 return _impl->finalize();
60}
61
63{
64 return _impl->getMeshDimensions(toSV(meshName));
65}
66
68{
69 return _impl->getDataDimensions(toSV(meshName), toSV(dataName));
70}
71
73{
74 return _impl->isCouplingOngoing();
75}
76
78{
79 return _impl->isTimeWindowComplete();
80}
81
83{
84 return _impl->getMaxTimeStepSize();
85}
86
88{
89 return _impl->requiresInitialData();
90}
91
93{
94 return _impl->requiresReadingCheckpoint();
95}
96
98{
99 return _impl->requiresWritingCheckpoint();
100}
101
103{
104 return _impl->requiresMeshConnectivityFor(toSV(meshName));
105}
106
108 ::precice::string_view dataName) const
109{
110 return _impl->requiresGradientDataFor(toSV(meshName), toSV(dataName));
111}
112
114 ::precice::string_view meshName,
116{
117 return _impl->setMeshVertex(toSV(meshName), coordinates);
118}
119
121 ::precice::string_view meshName) const
122{
123 return _impl->getMeshVertexSize(toSV(meshName));
124}
125
127 ::precice::string_view meshName,
130{
131 _impl->setMeshVertices(toSV(meshName), coordinates, ids);
132}
133
135 ::precice::string_view meshName,
136 VertexID first,
137 VertexID second)
138{
139 _impl->setMeshEdge(toSV(meshName), first, second);
140}
141
143 ::precice::string_view meshName,
145{
146 _impl->setMeshEdges(toSV(meshName), ids);
147}
148
150 ::precice::string_view meshName,
151 VertexID first,
152 VertexID second,
153 VertexID third)
154{
155 _impl->setMeshTriangle(toSV(meshName), first, second, third);
156}
157
159 ::precice::string_view meshName,
161{
162 _impl->setMeshTriangles(toSV(meshName), ids);
163}
164
166 ::precice::string_view meshName,
167 VertexID first,
168 VertexID second,
169 VertexID third,
170 VertexID fourth)
171{
172 _impl->setMeshQuad(toSV(meshName), first, second, third, fourth);
173}
174
176 ::precice::string_view meshName,
178{
179 _impl->setMeshQuads(toSV(meshName), ids);
180}
181
183 ::precice::string_view meshName,
184 VertexID first,
185 VertexID second,
186 VertexID third,
187 VertexID fourth)
188{
189 _impl->setMeshTetrahedron(toSV(meshName), first, second, third, fourth);
190}
191
193 ::precice::string_view meshName,
195{
196 _impl->setMeshTetrahedra(toSV(meshName), ids);
197}
198
200 ::precice::string_view meshName,
201 ::precice::string_view dataName,
204{
205 _impl->writeData(toSV(meshName), toSV(dataName), ids, values);
206}
207
209 ::precice::string_view meshName,
210 ::precice::string_view dataName,
212 double relativeReadTime,
213 ::precice::span<double> values) const
214{
215 _impl->readData(toSV(meshName), toSV(dataName), ids, relativeReadTime, values);
216}
217
219 ::precice::span<const double> boundingBox) const
220{
221 _impl->setMeshAccessRegion(toSV(meshName), boundingBox);
222}
223
226 ::precice::span<double> coordinates) const
227{
228 _impl->getMeshVertexIDsAndCoordinates(toSV(meshName), ids, coordinates);
229}
230
232 ::precice::string_view meshName,
233 ::precice::string_view dataName,
236{
237 _impl->writeGradientData(toSV(meshName), toSV(dataName), ids, gradients);
238}
239
240} // namespace precice
int getMeshVertexSize(::precice::string_view meshName) const
Returns the number of vertices of a mesh.
void setMeshAccessRegion(::precice::string_view meshName, ::precice::span< const double > boundingBox) const
setMeshAccessRegion Define a region of interest on a received mesh (<receive-mesh ....
void setMeshQuads(::precice::string_view meshName, ::precice::span< const VertexID > ids)
Sets multiple mesh quads from vertex IDs.
int getDataDimensions(::precice::string_view meshName, ::precice::string_view dataName) const
Returns the spatial dimensionality of the given data on the given mesh.
Participant(::precice::string_view participantName, ::precice::string_view configurationFileName, int solverProcessIndex, int solverProcessSize)
Constructs a Participant for the given participant.
bool isCouplingOngoing() const
Checks if the coupled simulation is still ongoing.
void getMeshVertexIDsAndCoordinates(::precice::string_view meshName, ::precice::span< VertexID > ids, ::precice::span< double > coordinates) const
getMeshVertexIDsAndCoordinates Iterates over the region of interest defined by bounding boxes and rea...
VertexID setMeshVertex(::precice::string_view meshName, ::precice::span< const double > position)
Creates a mesh vertex.
void finalize()
Finalizes preCICE.
void setMeshTetrahedra(::precice::string_view meshName, ::precice::span< const VertexID > ids)
Sets multiple mesh tetrahedra from vertex IDs.
void initialize()
Fully initializes preCICE and coupling data.
void writeData(::precice::string_view meshName, ::precice::string_view dataName, ::precice::span< const VertexID > ids, ::precice::span< const double > values)
Writes data to a mesh.
void setMeshTriangle(::precice::string_view meshName, VertexID first, VertexID second, VertexID third)
Sets mesh triangle from vertex IDs.
std::unique_ptr< impl::ParticipantImpl > _impl
Pointer to implementation of Participant.
void setMeshVertices(::precice::string_view meshName, ::precice::span< const double > coordinates, ::precice::span< VertexID > ids)
Creates multiple mesh vertices.
void writeGradientData(::precice::string_view meshName, ::precice::string_view dataName, ::precice::span< const VertexID > ids, ::precice::span< const double > gradients)
Writes vector gradient data to a mesh.
double getMaxTimeStepSize() const
Get the maximum allowed time step size of the current window.
bool requiresGradientDataFor(::precice::string_view meshName, ::precice::string_view dataName) const
Checks if the given data set requires gradient data. We check if the data object has been initialized...
bool requiresMeshConnectivityFor(::precice::string_view meshName) const
Checks if the given mesh requires connectivity.
void advance(double computedTimeStepSize)
Advances preCICE after the solver has computed one time step.
int getMeshDimensions(::precice::string_view meshName) const
Returns the spatial dimensionality of the given mesh.
void setMeshTriangles(::precice::string_view meshName, ::precice::span< const VertexID > ids)
Sets multiple mesh triangles from vertex IDs.
void setMeshEdges(::precice::string_view meshName, ::precice::span< const VertexID > ids)
Sets multiple mesh edges from vertex IDs.
void readData(::precice::string_view meshName, ::precice::string_view dataName, ::precice::span< const VertexID > ids, double relativeReadTime, ::precice::span< double > values) const
Reads data values from a mesh. Values correspond to a given point in time relative to the beginning o...
void setMeshEdge(::precice::string_view meshName, VertexID first, VertexID second)
Sets a mesh edge from vertex IDs.
bool isTimeWindowComplete() const
Checks if the current coupling window is completed.
void setMeshTetrahedron(::precice::string_view meshName, VertexID first, VertexID second, VertexID third, VertexID fourth)
Set tetrahedron in 3D mesh from vertex ID.
void setMeshQuad(::precice::string_view meshName, VertexID first, VertexID second, VertexID third, VertexID fourth)
Sets a planar surface mesh quadrangle from vertex IDs.
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
constexpr pointer data() const noexcept
Definition span.hpp:500
constexpr size_type size() const noexcept
Definition span.hpp:469
Main namespace of the precice library.
int VertexID
Definition Types.hpp:13
STL namespace.