preCICE v3.2.0
Loading...
Searching...
No Matches
Participant.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <memory>
4#include <precice/Version.h>
5#include <precice/export.h>
6#include <precice/span.hpp>
7#include "precice/Types.hpp"
8
12namespace precice {
13namespace impl {
14class ParticipantImpl;
15}
16namespace testing {
17struct WhiteboxAccessor;
18}
19} // namespace precice
20
21// ----------------------------------------------------------- CLASS DEFINITION
22
23namespace precice {
24
64
168class PRECICE_API Participant {
169public:
188
201 ::precice::string_view participantName,
202 ::precice::string_view configurationFileName,
203 int solverProcessIndex,
204 int solverProcessSize);
205
219 ::precice::string_view participantName,
220 ::precice::string_view configurationFileName,
221 int solverProcessIndex,
222 int solverProcessSize,
223 void *communicator);
224
226
228
241
259 void initialize();
260
299 void advance(double computedTimeStepSize);
300
325 void finalize();
326
328
375
387
401
403
416
424 int getMeshDimensions(::precice::string_view meshName) const;
425
439
456 bool isCouplingOngoing() const;
457
477 bool isTimeWindowComplete() const;
478
489 double getMaxTimeStepSize() const;
490
492
515
527
546 void resetMesh(::precice::string_view meshName);
547
561 ::precice::string_view meshName,
563
577 int getMeshVertexSize(::precice::string_view meshName) const;
578
594 void setMeshVertices(
595 ::precice::string_view meshName,
598
612 void setMeshEdge(
613 ::precice::string_view meshName,
614 VertexID first,
615 VertexID second);
616
634 void setMeshEdges(
635 ::precice::string_view meshName,
637
654 void setMeshTriangle(
655 ::precice::string_view meshName,
656 VertexID first,
657 VertexID second,
658 VertexID third);
659
677 void setMeshTriangles(
678 ::precice::string_view meshName,
680
699 void setMeshQuad(
700 ::precice::string_view meshName,
701 VertexID first,
702 VertexID second,
703 VertexID third,
704 VertexID fourth);
705
725 void setMeshQuads(
726 ::precice::string_view meshName,
728
745 ::precice::string_view meshName,
746 VertexID first,
747 VertexID second,
748 VertexID third,
749 VertexID fourth);
750
769 ::precice::string_view meshName,
771
773
795
805 bool requiresInitialData();
806
830 void writeData(
831 ::precice::string_view meshName,
832 ::precice::string_view dataName,
835
867 void readData(
868 ::precice::string_view meshName,
869 ::precice::string_view dataName,
871 double relativeReadTime,
872 ::precice::span<double> values) const;
874
899
935 void writeAndMapData(
936 ::precice::string_view meshName,
937 ::precice::string_view dataName,
940
989 void mapAndReadData(
990 ::precice::string_view meshName,
991 ::precice::string_view dataName,
993 double relativeReadTime,
994 ::precice::span<double> values) const;
995
997
1010
1068 ::precice::string_view meshName,
1069 ::precice::span<const double> boundingBox) const;
1070
1093 ::precice::string_view meshName,
1095 ::precice::span<double> coordinates) const;
1096
1098
1103
1119 ::precice::string_view dataName) const;
1120
1163 void writeGradientData(
1164 ::precice::string_view meshName,
1165 ::precice::string_view dataName,
1168
1170
1185
1197
1203
1205
1207 Participant(const Participant &copy) = delete;
1208
1210 Participant &operator=(const Participant &assign) = delete;
1211
1212private:
1215
1216 // @brief To allow white box tests.
1218};
1219
1220} // 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.
Participant(const Participant &copy)=delete
Disable copy construction.
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.
void startProfilingSection(::precice::string_view sectionName)
int getMeshDimensions(::precice::string_view meshName) const
Returns the spatial dimensionality of the given mesh.
void mapAndReadData(::precice::string_view meshName, ::precice::string_view dataName, ::precice::span< const double > coordinates, double relativeReadTime, ::precice::span< double > values) const
Reads data values from a mesh using a just-in-time data mapping. Values correspond to a given point i...
Participant & operator=(const Participant &assign)=delete
Disable assignment construction.
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.
void resetMesh(::precice::string_view meshName)
Removes all vertices and connectivity information from the mesh.
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.
void writeAndMapData(::precice::string_view meshName, ::precice::string_view dataName, ::precice::span< const double > coordinates, ::precice::span< const double > values)
Writes data values to a mesh using a just-in-time mapping (experimental).
Implementation of Participant. See also pimpl ideom (https://en.cppreference.com/w/cpp/language/pimpl...
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
contains the testing framework.
Definition helper.hpp:9
Main namespace of the precice library.
::precice::span< const char > string_view
forwards-compatible typedef for C++17 std::string_view.
int VertexID
Definition Types.hpp:13
struct giving access to the impl of a befriended class or struct
Definition Testing.hpp:48