preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
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
259 void initialize();
260
299 void advance(double computedTimeStepSize);
300
325 void finalize();
326
328
386 bool requiresWritingCheckpoint();
387
400 bool requiresReadingCheckpoint();
401
403
424 int getMeshDimensions(::precice::string_view meshName) const;
425
438 int getDataDimensions(::precice::string_view meshName, ::precice::string_view dataName) const;
439
456 bool isCouplingOngoing() const;
457
477 bool isTimeWindowComplete() const;
478
489 double getMaxTimeStepSize() const;
490
492
526 bool requiresMeshConnectivityFor(::precice::string_view meshName) const;
527
546 void resetMesh(::precice::string_view meshName);
547
560 VertexID setMeshVertex(
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
744 void setMeshTetrahedron(
745 ::precice::string_view meshName,
746 VertexID first,
747 VertexID second,
748 VertexID third,
749 VertexID fourth);
750
768 void setMeshTetrahedra(
769 ::precice::string_view meshName,
771
773
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
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
1067 void setMeshAccessRegion(
1068 ::precice::string_view meshName,
1069 ::precice::span<const double> boundingBox) const;
1070
1092 void getMeshVertexIDsAndCoordinates(
1093 ::precice::string_view meshName,
1095 ::precice::span<double> coordinates) const;
1096
1098
1103
1118 bool requiresGradientDataFor(::precice::string_view meshName,
1119 ::precice::string_view dataName) const;
1120
1163 void writeGradientData(
1164 ::precice::string_view meshName,
1165 ::precice::string_view dataName,
1168
1170
1185
1196 void startProfilingSection(::precice::string_view sectionName);
1197
1202 void stopLastProfilingSection();
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
Main Application Programming Interface of preCICE. Include using #include <precice/precice....
Participant(const Participant &copy)=delete
Disable copy construction.
std::unique_ptr< impl::ParticipantImpl > _impl
Pointer to implementation of Participant.
Participant & operator=(const Participant &assign)=delete
Disable assignment construction.
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
Main namespace of the precice library.
int VertexID
Definition Types.hpp:13
static std::unique_ptr< precice::Participant > impl
Definition preciceC.cpp:21
struct giving access to the impl of a befriended class or struct
Definition Testing.hpp:48