preCICE v3.2.0
|
Implementation of Participant. See also pimpl ideom (https://en.cppreference.com/w/cpp/language/pimpl). More...
#include <ParticipantImpl.hpp>
Classes | |
struct | MappedSamples |
Public Member Functions | |
const mesh::Mesh & | mesh (const std::string &meshName) const |
Allows to access a registered mesh. | |
MappedSamples | mappedSamples () const |
Returns the amount of mapped read and write samples in the last call to advance. | |
ParticipantImpl (ParticipantImpl const &)=delete | |
Disable copy construction. | |
ParticipantImpl & | operator= (ParticipantImpl const &)=delete |
Disable assignment construction. | |
ParticipantImpl (ParticipantImpl &&)=delete | |
Disable move construction. | |
ParticipantImpl & | operator= (ParticipantImpl &&)=delete |
Disable move assignment. | |
Construction and Configuration | |
ParticipantImpl (std::string_view participantName, std::string_view configurationFileName, int solverProcessIndex, int solverProcessSize, std::optional< void * > communicator) | |
Generic constructor for ParticipantImpl. | |
~ParticipantImpl () | |
Destructor. | |
Steering Methods | |
void | initialize () |
Fully initializes preCICE and coupling data. | |
void | advance (double computedTimeStepSize) |
Advances preCICE after the solver has computed one time step. | |
void | finalize () |
Finalizes preCICE. | |
Status Queries | |
int | getMeshDimensions (std::string_view meshName) const |
Returns the spatial dimensionality of the given mesh. | |
int | getDataDimensions (std::string_view meshName, std::string_view dataName) const |
Returns the spatial dimensionality of the given data on the given mesh. | |
bool | isCouplingOngoing () const |
Checks if the coupled simulation is still ongoing. | |
bool | isTimeWindowComplete () const |
Checks if the current coupling window is completed. | |
double | getMaxTimeStepSize () const |
Get the maximum allowed time step size of the current window. | |
Requirements | |
bool | requiresInitialData () |
bool | requiresReadingCheckpoint () |
bool | requiresWritingCheckpoint () |
Mesh Access | |
void | resetMesh (std::string_view meshName) |
Removes all vertices and connectivity information from the mesh. | |
bool | requiresMeshConnectivityFor (std::string_view meshName) const |
Checks if the given mesh requires connectivity. | |
bool | requiresGradientDataFor (std::string_view meshName, std::string_view dataName) const |
Checks if the given data set requires gradient data. We check if the data object has been initialized with the gradient flag. | |
VertexID | setMeshVertex (std::string_view meshName, ::precice::span< const double > position) |
Creates a mesh vertex. | |
int | getMeshVertexSize (std::string_view meshName) const |
Returns the number of vertices of a mesh. | |
void | setMeshVertices (std::string_view meshName, ::precice::span< const double > positions, ::precice::span< VertexID > ids) |
Creates multiple mesh vertices. | |
void | setMeshEdge (std::string_view meshName, VertexID first, VertexID second) |
Sets a mesh edge from vertex IDs. | |
void | setMeshEdges (std::string_view meshName, ::precice::span< const VertexID > vertices) |
Sets multiple mesh edges from vertex IDs. | |
void | setMeshTriangle (std::string_view meshName, VertexID first, VertexID second, VertexID third) |
Sets mesh triangle from vertex IDs. | |
void | setMeshTriangles (std::string_view meshName, ::precice::span< const VertexID > vertices) |
Sets multiple mesh triangles from vertex IDs. | |
void | setMeshQuad (std::string_view meshName, VertexID first, VertexID second, VertexID third, VertexID fourth) |
Sets a planar surface mesh quadrangle from vertex IDs. | |
void | setMeshQuads (std::string_view meshName, ::precice::span< const VertexID > vertices) |
Sets multiple mesh quads from vertex IDs. | |
void | setMeshTetrahedron (std::string_view meshName, VertexID first, VertexID second, VertexID third, VertexID fourth) |
Set tetrahedron in 3D mesh from vertex ID. | |
void | setMeshTetrahedra (std::string_view meshName, ::precice::span< const VertexID > vertices) |
Sets multiple mesh tetrahedra from vertex IDs. | |
Data Access | |
void | readData (std::string_view meshName, std::string_view dataName, ::precice::span< const VertexID > vertices, 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 of the current timestep. | |
void | mapAndReadData (std::string_view meshName, std::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 in time relative to the beginning of the current timestep (experimental). | |
void | writeAndMapData (std::string_view meshName, std::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). | |
void | writeData (std::string_view meshName, std::string_view dataName, ::precice::span< const VertexID > vertices, ::precice::span< const double > values) |
Writes data to a mesh. | |
void | writeGradientData (std::string_view meshName, std::string_view dataName, ::precice::span< const VertexID > vertices, ::precice::span< const double > gradients) |
Writes vector gradient data to a mesh. | |
Direct Access | |
void | setMeshAccessRegion (std::string_view meshName, ::precice::span< const double > boundingBox) const |
setMeshAccessRegion Define a region of interest on a received mesh (<receive-mesh ... from="otherParticipant />") in order to receive only a certain mesh region. Have a look at the website under https://precice.org/couple-your-code-direct-access.html or navigate manually to the page Docs->Couple your code -> Advanced topics -> Accessing received meshes directly for a comprehensive documentation | |
void | getMeshVertexIDsAndCoordinates (std::string_view meshName, ::precice::span< VertexID > ids, ::precice::span< double > coordinates) const |
getMeshVertexIDsAndCoordinates Iterates over the region of interest defined by bounding boxes and reads the corresponding coordinates omitting the mapping. | |
User profiling | |
void | startProfilingSection (std::string_view eventName) |
void | stopLastProfilingSection () |
Private Types | |
enum struct | State { Constructed , Initialized , Finalized } |
Represents the various states a Participant can be in. More... | |
enum struct | ExportTiming : bool { Advance = false , Initial = true } |
enum class | CloseChannels : bool { All = false , Distributed = true } |
Which channels to close in closeCommunicationChannels() More... | |
using | MeshChanges = std::vector<int> |
How many ranks have changed each used mesh. | |
Private Member Functions | |
void | configure (std::string_view configurationFileName) |
Configures the coupling interface from the given xml file. | |
void | configure (const config::Configuration &configuration) |
Configures the coupling interface with a prepared configuration. | |
void | configureM2Ns (const m2n::M2NConfiguration::SharedPointer &config) |
void | handleExports (ExportTiming timing) |
void | configurePartitions (const m2n::M2NConfiguration::SharedPointer &m2nConfig) |
Determines participants providing meshes to other participants. | |
void | compareBoundingBoxes () |
Communicate bounding boxes and look for overlaps. | |
void | computePartitions () |
Communicate meshes and create partitions. | |
void | computeMappings (std::vector< MappingContext > &contexts, const std::string &mappingType) |
Helper for mapWrittenData and mapReadData. | |
void | mapInitialWrittenData () |
Computes, and performs write mappings of the initial data in initialize. | |
void | mapWrittenData (std::optional< double > after=std::nullopt) |
Computes, and performs suitable write mappings either entirely or after given time. | |
void | mapInitialReadData () |
void | mapReadData () |
void | trimReadMappedData (double timeAfterAdvance, bool isTimeWindowComplete, const cplscheme::ImplicitData &fromData) |
Removes samples in mapped to data connected to received data via a mapping. | |
void | performDataActions (const std::set< action::Action::Timing > &timings) |
Performs all data actions with given timing. | |
void | resetWrittenData () |
Resets written data. | |
impl::PtrParticipant | determineAccessingParticipant (const config::Configuration &config) |
Determines participant accessing this interface from the configuration. | |
void | initializeIntraCommunication () |
Initializes intra-participant communication. | |
void | advanceCouplingScheme () |
Advances the coupling schemes. | |
void | syncTimestep (double computedTimeStepSize) |
Syncs the time step size between all ranks (all time steps sizes should be the same!) | |
void | closeCommunicationChannels (CloseChannels cc) |
Syncs the primary ranks of all connected participants. | |
void | handleDataBeforeAdvance (bool reachedTimeWindowEnd, double timeSteppedTo) |
Completes everything data-related between adding time to and advancing the coupling scheme. | |
void | handleDataAfterAdvance (bool reachedTimeWindowEnd, bool isTimeWindowComplete, double timeSteppedTo, double timeAfterAdvance, const cplscheme::ImplicitData &receivedData) |
Completes everything data-related after advancing the coupling scheme. | |
void | samplizeWriteData (double time) |
Creates a Stample at the given time for each write Data and zeros the buffers. | |
void | trimOldDataBefore (double time) |
Discards data before the given time for all meshes and data known by this participant. | |
void | trimSendDataAfter (double time) |
Discards send (currently write) data of a participant after a given time when another iteration is required. | |
MeshChanges | getTotalMeshChanges () const |
void | clearStamplesOfChangedMeshes (MeshChanges totalMeshChanges) |
Clears stample of changed meshes to make them consistent after the reinitialization. | |
bool | reinitHandshake (bool requestReinit) const |
void | reinitialize () |
Reinitializes preCICE. | |
void | setupCommunication () |
Connect participants including repartitioning. | |
void | setupWatcher () |
Setup mesh watcher such as WatchPoints. | |
bool | requiresUserDefinedAccessRegion (std::string_view meshName) const |
Private Attributes | |
logging::Logger | _log {"impl::ParticipantImpl"} |
std::string | _accessorName |
int | _accessorProcessRank |
int | _accessorCommunicatorSize |
impl::PtrParticipant | _accessor |
int | _dimensions = 0 |
Spatial dimensions of problem. | |
utils::MultiLock< std::string > | _meshLock |
std::map< std::string, int > | _meshIDs |
mesh name to mesh ID mapping. | |
std::map< std::string, m2n::BoundM2N > | _m2ns |
std::vector< impl::PtrParticipant > | _participants |
Holds information about solvers participating in the coupled simulation. | |
cplscheme::PtrCouplingScheme | _couplingScheme |
bool | _allowsExperimental = false |
Are experimental API calls allowed? | |
bool | _allowsRemeshing = false |
Are experimental remeshing API calls allowed? | |
bool | _waitInFinalize = false |
Are participants waiting for each other in finalize? | |
State | _state {State::Constructed} |
The current State of the Participant. | |
long int | _numberAdvanceCalls = 0 |
Counts calls to advance for plotting. | |
int | _executedWriteMappings = 0 |
Counts the amount of samples mapped in write mappings executed in the latest advance. | |
int | _executedReadMappings = 0 |
Counts the amount of samples mapped in read mappings executed in the latest advance. | |
std::string | _configHash |
The hash of the configuration file used to configure this participant. | |
std::unique_ptr< profiling::Event > | _solverInitEvent |
std::unique_ptr< profiling::Event > | _solverAdvanceEvent |
std::vector< profiling::Event > | _userEvents |
Friends | |
struct | Integration::Serial::Whitebox::TestConfigurationPeano |
To allow white box tests. | |
struct | Integration::Serial::Whitebox::TestConfigurationComsol |
Implementation of Participant. See also pimpl ideom (https://en.cppreference.com/w/cpp/language/pimpl).
Definition at line 53 of file ParticipantImpl.hpp.
|
private |
How many ranks have changed each used mesh.
Definition at line 484 of file ParticipantImpl.hpp.
|
strongprivate |
Which channels to close in closeCommunicationChannels()
Enumerator | |
---|---|
All | |
Distributed |
Definition at line 460 of file ParticipantImpl.hpp.
|
strongprivate |
Enumerator | |
---|---|
Advance | |
Initial |
Definition at line 390 of file ParticipantImpl.hpp.
|
strongprivate |
Represents the various states a Participant can be in.
Enumerator | |
---|---|
Constructed | |
Initialized | |
Finalized |
Definition at line 340 of file ParticipantImpl.hpp.
precice::impl::ParticipantImpl::ParticipantImpl | ( | std::string_view | participantName, |
std::string_view | configurationFileName, | ||
int | solverProcessIndex, | ||
int | solverProcessSize, | ||
std::optional< void * > | communicator ) |
Generic constructor for ParticipantImpl.
Use the parameter communicator to specify a custom global MPI communicator. Pass std::nullopt to signal preCICE to use MPI_COMM_WORLD.
[in] | participantName | Name of the participant using the interface. Has to match the name given for a participant in the xml configuration file. |
[in] | configurationFileName | Name (with path) of the xml configuration file. |
[in] | solverProcessIndex | If the solver code runs with several processes, each process using preCICE has to specify its index, which has to start from 0 and end with solverProcessSize - 1. |
[in] | solverProcessSize | The number of solver processes using preCICE. |
[in] | communicator | An optional pointer to an MPI_Comm to use as communicator. |
Definition at line 82 of file ParticipantImpl.cpp.
precice::impl::ParticipantImpl::~ParticipantImpl | ( | ) |
Destructor.
Ensures that finalize() has been called.
Definition at line 161 of file ParticipantImpl.cpp.
|
delete |
Disable copy construction.
|
delete |
Disable move construction.
void precice::impl::ParticipantImpl::advance | ( | double | computedTimeStepSize | ) |
Advances preCICE after the solver has computed one time step.
There are two cases to distinguish at this point. If computedTimeStepSize
== getMaxTimeStepSize(), then the solver has reached the end of a time window and proceeds the coupling. This call is a computationally expensive process as it involves among other tasks:
If computedTimeStepSize
< getMaxTimeStepSize(), then the solver hasn't reached the end of a time window and it is subcycling. Depending on the configuration, written data can be used by preCICE to generate additional samples allowing for time interpolation using readData(). This call is computationally inexpensive.
[in] | computedTimeStepSize | Size of time step used by the solver. |
computedTimeStep
to advance().computedTimeStepSize
.Definition at line 399 of file ParticipantImpl.cpp.
|
private |
Advances the coupling schemes.
Definition at line 1890 of file ParticipantImpl.cpp.
|
private |
Clears stample of changed meshes to make them consistent after the reinitialization.
Definition at line 1988 of file ParticipantImpl.cpp.
|
private |
Syncs the primary ranks of all connected participants.
Definition at line 1904 of file ParticipantImpl.cpp.
|
private |
Communicate bounding boxes and look for overlaps.
Definition at line 1624 of file ParticipantImpl.cpp.
|
private |
Helper for mapWrittenData and mapReadData.
Definition at line 1697 of file ParticipantImpl.cpp.
|
private |
Communicate meshes and create partitions.
Definition at line 1644 of file ParticipantImpl.cpp.
|
private |
Configures the coupling interface with a prepared configuration.
Can be used to configure the ParticipantImpl without xml file. Requires to manually setup the configuration object.
|
private |
Configures the coupling interface from the given xml file.
Only after the configuration a reasonable state of a ParticipantImpl object is achieved.
[in] | configurationFileName | Name (with path) of the xml config. file. |
Definition at line 169 of file ParticipantImpl.cpp.
|
private |
|
private |
Determines participants providing meshes to other participants.
Definition at line 1574 of file ParticipantImpl.cpp.
|
private |
Determines participant accessing this interface from the configuration.
Definition at line 1846 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::finalize | ( | ) |
Finalizes preCICE.
If initialize() has been called:
Always:
Definition at line 589 of file ParticipantImpl.cpp.
int precice::impl::ParticipantImpl::getDataDimensions | ( | std::string_view | meshName, |
std::string_view | dataName ) const |
Returns the spatial dimensionality of the given data on the given mesh.
Note that vectorial data dimensionality directly depends on the spatial dimensionality of the mesh.
[in] | meshName | the name of the associated mesh |
[in] | dataName | the name of the data to get the dimensions for |
Definition at line 648 of file ParticipantImpl.cpp.
double precice::impl::ParticipantImpl::getMaxTimeStepSize | ( | ) | const |
Get the maximum allowed time step size of the current window.
Allows the user to query the maximum allowed time step size in the current window. This should be used to compute the actual time step that the solver uses.
Definition at line 672 of file ParticipantImpl.cpp.
int precice::impl::ParticipantImpl::getMeshDimensions | ( | std::string_view | meshName | ) | const |
Returns the spatial dimensionality of the given mesh.
[in] | meshName | the name of the associated mesh |
Definition at line 641 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::getMeshVertexIDsAndCoordinates | ( | std::string_view | meshName, |
::precice::span< VertexID > | ids, | ||
::precice::span< double > | coordinates ) const |
getMeshVertexIDsAndCoordinates Iterates over the region of interest defined by bounding boxes and reads the corresponding coordinates omitting the mapping.
[in] | meshName | corresponding mesh name |
[out] | ids | ids corresponding to the coordinates |
[out] | coordinates | the coordinates associated to the ids and corresponding data values |
initialize()
, if the meshName
corresponds to a received mesh, since the relevant mesh data is exchanged during the initialize()
call.Definition at line 1471 of file ParticipantImpl.cpp.
int precice::impl::ParticipantImpl::getMeshVertexSize | ( | std::string_view | meshName | ) | const |
Returns the number of vertices of a mesh.
[in] | meshName | the name of the mesh |
setMeshAccessRegion()
).initialize()
, if the meshName
corresponds to a received mesh, since the relevant mesh data is exchanged during the initialize()
call. Definition at line 741 of file ParticipantImpl.cpp.
|
private |
Allreduce of the amount of changed meshes on each rank.
Definition at line 1965 of file ParticipantImpl.cpp.
|
private |
Completes everything data-related after advancing the coupling scheme.
Definition at line 499 of file ParticipantImpl.cpp.
|
private |
Completes everything data-related between adding time to and advancing the coupling scheme.
Definition at line 471 of file ParticipantImpl.cpp.
|
private |
Exports meshes with data and watch point data.
[in] | timing | when the exports are requested |
Definition at line 1814 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::initialize | ( | ) |
Fully initializes preCICE and coupling data.
- Sets up a connection to the other participants of the coupled simulation.
Definition at line 249 of file ParticipantImpl.cpp.
|
private |
Initializes intra-participant communication.
Definition at line 1861 of file ParticipantImpl.cpp.
bool precice::impl::ParticipantImpl::isCouplingOngoing | ( | ) | const |
Checks if the coupled simulation is still ongoing.
A coupling is ongoing as long as
Definition at line 656 of file ParticipantImpl.cpp.
bool precice::impl::ParticipantImpl::isTimeWindowComplete | ( | ) | const |
Checks if the current coupling window is completed.
The following reasons require several solver time steps per time window:
Hence, a time window is complete if we reach the end of the time window and the implicit coupling has converged.
For implicit coupling this condition is equivalent with the requirement to write an iteration checkpoint. This is, however, not the case for explicit coupling.
Definition at line 664 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::mapAndReadData | ( | std::string_view | meshName, |
std::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 in time relative to the beginning of the current timestep (experimental).
This function reads values at temporary locations from data of a mesh. As opposed to the readData function using VertexIDs, this function allows reading data via coordinates, which don't have to be specified during the initialization. This is particularly useful for meshes, which vary over time. Note that using this function comes at a performance cost, since the specified mapping needs to be computed locally for the given locations, whereas the other variant (readData) can typically exploit the static interface mesh and pre-compute data structures more efficient.
Values are read into a block of continuous memory defined by values in the order specified by vertices.
The 1D/Scalar-format of values is (d0, d1, ..., dn) The 2D-format of values is (d0x, d0y, d1x, d1y, ..., dnx, dny) The 3D-format of values is (d0x, d0y, d0z, d1x, d1y, d1z, ..., dnx, dny, dnz)
The data is read at relativeReadTime, which indicates the point in time measured from the beginning of the current time step. relativeReadTime = 0 corresponds to data at the beginning of the time step. Assuming that the user will call advance(dt) at the end of the time step, dt indicates the size of the current time step. Then relativeReadTime = dt corresponds to the data at the end of the time step.
[in] | meshName | The name of the mesh that holds the data, needs to be a mesh received from another participant. |
[in] | dataName | The name of the data to read from. |
[in] | coordinates | a span to the coordinates of the vertices The 2D-format is (d0x, d0y, d1x, d1y, ..., dnx, dny) The 3D-format is (d0x, d0y, d0z, d1x, d1y, d1z, ..., dnx, dny, dnz) |
[in] | relativeReadTime | Point in time where data is read relative to the beginning of the current time step. |
[out] | values | The destination memory to read the data from. |
values
contain the read data as specified in the above format.relativeReadTime
on each individual rank. As a result, calling the function for different coordinates
while keeping the relativeReadTime
constant is more efficient than calling it for varying relativeReadTime
values at each coordinate
. In practice, this means that iterating over time (i.e., varying relativeReadTime
) in the outer loop and over space (i.e., iterating over coordinates
) in the inner loop gives a significantly better performance than vice versa.Definition at line 1262 of file ParticipantImpl.cpp.
|
private |
|
private |
Computes, and performs write mappings of the initial data in initialize.
Definition at line 1717 of file ParticipantImpl.cpp.
ParticipantImpl::MappedSamples precice::impl::ParticipantImpl::mappedSamples | ( | ) | const |
Returns the amount of mapped read and write samples in the last call to advance.
Definition at line 1955 of file ParticipantImpl.cpp.
|
private |
|
private |
Computes, and performs suitable write mappings either entirely or after given time.
Definition at line 1734 of file ParticipantImpl.cpp.
const mesh::Mesh & precice::impl::ParticipantImpl::mesh | ( | const std::string & | meshName | ) | const |
Allows to access a registered mesh.
Definition at line 1949 of file ParticipantImpl.cpp.
|
delete |
Disable move assignment.
|
delete |
Disable assignment construction.
|
private |
Performs all data actions with given timing.
[in] | timings | the timings of the action. |
Definition at line 1804 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::readData | ( | std::string_view | meshName, |
std::string_view | dataName, | ||
::precice::span< const VertexID > | vertices, | ||
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 of the current timestep.
This function reads values of specified vertices from data of a mesh. Values are read into a block of continuous memory defined by values in the order specified by vertices.
The 1D/Scalar-format of values is (d0, d1, ..., dn) The 2D-format of values is (d0x, d0y, d1x, d1y, ..., dnx, dny) The 3D-format of values is (d0x, d0y, d0z, d1x, d1y, d1z, ..., dnx, dny, dnz)
The data is read at relativeReadTime, which indicates the point in time measured from the beginning of the current time step. relativeReadTime = 0 corresponds to data at the beginning of the time step. Assuming that the user will call advance(dt) at the end of the time step, dt indicates the size of the current time step. Then relativeReadTime = dt corresponds to the data at the end of the time step.
[in] | meshName | the name of mesh that hold the data. |
[in] | dataName | the name of the data to read from. |
[in] | ids | the vertex ids of the vertices to read data from. |
[in] | relativeReadTime | Point in time where data is read relative to the beginning of the current time step. |
[out] | values | the destination memory to read the data from. |
Definition at line 1211 of file ParticipantImpl.cpp.
|
private |
Exchanges request to remesh with all connecting participants.
[in] | requestReinit | does this participant request to remesh? |
Definition at line 1998 of file ParticipantImpl.cpp.
|
private |
Reinitializes preCICE.
Definition at line 305 of file ParticipantImpl.cpp.
bool precice::impl::ParticipantImpl::requiresGradientDataFor | ( | std::string_view | meshName, |
std::string_view | dataName ) const |
Checks if the given data set requires gradient data. We check if the data object has been initialized with the gradient flag.
preCICE may require gradient data information from the solver and ignores any API calls regarding gradient data if it is not required. (When applying a nearest-neighbor-gradient mapping)
[in] | meshName | the name of mesh that hold the data. |
[in] | dataName | the name of the data. |
Definition at line 729 of file ParticipantImpl.cpp.
bool precice::impl::ParticipantImpl::requiresInitialData | ( | ) |
Checks if the participant is required to provide initial data.
If true, then the participant needs to write initial data to defined vertices prior to calling initialize().
Definition at line 689 of file ParticipantImpl.cpp.
bool precice::impl::ParticipantImpl::requiresMeshConnectivityFor | ( | std::string_view | meshName | ) | const |
Checks if the given mesh requires connectivity.
preCICE may require connectivity information from the solver and ignores any API calls regarding connectivity if it is not required. Use this function to conditionally generate this connectivity.
[in] | meshName | the name of the mesh |
Definition at line 722 of file ParticipantImpl.cpp.
bool precice::impl::ParticipantImpl::requiresReadingCheckpoint | ( | ) |
Checks if the participant is required to read an iteration checkpoint.
If true, the participant is required to read an iteration checkpoint before calling advance().
Definition at line 711 of file ParticipantImpl.cpp.
|
private |
Returns if a user has to define an access region for direct mesh access and just-in-time mapping or not Right now, that's required in parallel runs on received meshes
Definition at line 1944 of file ParticipantImpl.cpp.
bool precice::impl::ParticipantImpl::requiresWritingCheckpoint | ( | ) |
Checks if the participant is required to write an iteration checkpoint.
If true, the participant is required to write an iteration checkpoint before calling advance().
Definition at line 700 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::resetMesh | ( | std::string_view | meshName | ) |
Removes all vertices and connectivity information from the mesh.
Allows redefining a mesh during runtime. After the call to resetMesh(), the mesh vertices need to be set with setMeshVertex() and setMeshVertices() again. Connectivity information may be set as well.
Reading data from this mesh using readData() is not possible until the next call to advance().
[in] | meshName | the name of the mesh to reset |
Definition at line 777 of file ParticipantImpl.cpp.
|
private |
Resets written data.
isAtWindowEnd | set true, if function is called at end of window to also trim the time sample storage |
isTimeWindowComplete | set true, if function is called at end of converged window to trim and move the sample storage. |
Definition at line 1837 of file ParticipantImpl.cpp.
|
private |
Creates a Stample at the given time for each write Data and zeros the buffers.
Definition at line 548 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshAccessRegion | ( | std::string_view | meshName, |
::precice::span< const double > | boundingBox ) const |
setMeshAccessRegion Define a region of interest on a received mesh (<receive-mesh ... from="otherParticipant />") in order to receive only a certain mesh region. Have a look at the website under https://precice.org/couple-your-code-direct-access.html or navigate manually to the page Docs->Couple your code -> Advanced topics -> Accessing received meshes directly for a comprehensive documentation
This function is required if you don't want to use the mapping schemes in preCICE, but rather want to use your own solver for data mapping. As opposed to the usual preCICE mapping, only a single mesh (from the other participant) is now involved in this situation since an 'own' mesh defined by the participant itself is not required any more. In order to re-partition the received mesh, the participant needs to define the mesh region it wants read data from and write data to. The mesh region is specified through an axis-aligned bounding box given by the lower and upper [min and max] bounding-box limits in each space dimension [x, y, z].
[in] | meshName | name of the mesh you want to access through the bounding box |
[in] | boundingBox | Axis aligned bounding boxes which has in 3D the format [x_min, x_max, y_min, y_max, z_min, z_max] |
initialize()
has not yet been called. Definition at line 1430 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshEdge | ( | std::string_view | meshName, |
VertexID | first, | ||
VertexID | second ) |
Sets a mesh edge from vertex IDs.
Ignored if preCICE doesn't require connectivity for the mesh.
[in] | meshName | name of the mesh to add the edge to |
[in] | first | ID of the first vertex of the edge |
[in] | second | ID of the second vertex of the edge |
Definition at line 851 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshEdges | ( | std::string_view | meshName, |
::precice::span< const VertexID > | vertices ) |
Sets multiple mesh edges from vertex IDs.
vertices contain pairs of vertex indices for each edge to define. The format follows: e1a, e1b, e2a, e2b, ... Ignored if preCICE doesn't require connectivity for the mesh.
[in] | meshName | the name of the mesh to add the n edges to |
[in] | ids | an array containing 2n vertex IDs for n edges |
ids
were added to the mesh with the name meshName ids.size()
is multiple of 2Definition at line 871 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshQuad | ( | std::string_view | meshName, |
VertexID | first, | ||
VertexID | second, | ||
VertexID | third, | ||
VertexID | fourth ) |
Sets a planar surface mesh quadrangle from vertex IDs.
The planar quad will be triangulated, maximizing area-to-circumference. Ignored if preCICE doesn't require connectivity for the mesh.
[in] | meshName | name of the mesh to add the Quad to |
[in] | first | ID of the first vertex of the Quad |
[in] | second | ID of the second vertex of the Quad |
[in] | third | ID of the third vertex of the Quad |
[in] | fourth | ID of the fourth vertex of the Quad |
Definition at line 984 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshQuads | ( | std::string_view | meshName, |
::precice::span< const VertexID > | vertices ) |
Sets multiple mesh quads from vertex IDs.
vertices contain quadruples of vertex indices for each quad to define. The format follows: q1a, q1b, q1c, q1d, q2a, q2b, q2c, q2d, ...
Each planar quad will be triangulated, maximizing area-to-circumference. Ignored if preCICE doesn't require connectivity for the mesh.
[in] | meshName | name of the mesh to add the n quads to |
[in] | ids | an array containing 4n vertex IDs for n quads |
ids
were added to the mesh with the name meshName ids.size()
is multiple of 4Definition at line 1038 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshTetrahedra | ( | std::string_view | meshName, |
::precice::span< const VertexID > | vertices ) |
Sets multiple mesh tetrahedra from vertex IDs.
vertices contain quadruples of vertex indices for each tetrahedron to define. The format follows: t1a, t1b, t1c, t1d, t2a, t2b, t2c, t2d, ... Ignored if preCICE doesn't require connectivity for the mesh.
[in] | meshName | name of the mesh to add the n tetrahedra to |
[in] | ids | an array containing 4n vertex IDs for n tetrahedra |
ids
were added to the mesh with the name meshName Definition at line 1133 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshTetrahedron | ( | std::string_view | meshName, |
VertexID | first, | ||
VertexID | second, | ||
VertexID | third, | ||
VertexID | fourth ) |
Set tetrahedron in 3D mesh from vertex ID.
[in] | meshName | name of the mesh to add the Tetrahedron to |
[in] | first | ID of the first vertex of the Tetrahedron |
[in] | second | ID of the second vertex of the Tetrahedron |
[in] | third | ID of the third vertex of the Tetrahedron |
[in] | fourth | ID of the fourth vertex of the Tetrahedron |
Definition at line 1103 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshTriangle | ( | std::string_view | meshName, |
VertexID | first, | ||
VertexID | second, | ||
VertexID | third ) |
Sets mesh triangle from vertex IDs.
[in] | meshName | name of the mesh to add the triangle to |
[in] | first | ID of the first vertex of the triangle |
[in] | second | ID of the second vertex of the triangle |
[in] | third | ID of the third vertex of the triangle |
Definition at line 907 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshTriangles | ( | std::string_view | meshName, |
::precice::span< const VertexID > | vertices ) |
Sets multiple mesh triangles from vertex IDs.
vertices contain triples of vertex indices for each triangle to define. The format follows: t1a, t1b, t1c, t2a, t2b, t2c, ... Ignored if preCICE doesn't require connectivity for the mesh.
[in] | meshName | name of the mesh to add the n triangles to |
[in] | ids | an array containing 3n vertex IDs for n triangles |
ids
were added to the mesh with the name meshName ids.size()
is multiple of 3Definition at line 945 of file ParticipantImpl.cpp.
VertexID precice::impl::ParticipantImpl::setMeshVertex | ( | std::string_view | meshName, |
::precice::span< const double > | position ) |
Creates a mesh vertex.
[in] | meshName | the name of the mesh to add the vertex to. |
[in] | position | the coordinates of the vertex. |
Definition at line 794 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::setMeshVertices | ( | std::string_view | meshName, |
::precice::span< const double > | positions, | ||
::precice::span< VertexID > | ids ) |
Creates multiple mesh vertices.
[in] | meshName | the name of the mesh to add the vertices to. |
[in] | coordinates | a span to the coordinates of the vertices The 2D-format is (d0x, d0y, d1x, d1y, ..., dnx, dny) The 3D-format is (d0x, d0y, d0z, d1x, d1y, d1z, ..., dnx, dny, dnz) |
[out] | ids | The ids of the created vertices |
coordinates.size()
== getMeshDimensions(meshName) * ids.size()Definition at line 818 of file ParticipantImpl.cpp.
|
private |
Connect participants including repartitioning.
Definition at line 326 of file ParticipantImpl.cpp.
|
private |
Setup mesh watcher such as WatchPoints.
Definition at line 387 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::startProfilingSection | ( | std::string_view | eventName | ) |
Start a named user-defined profiling section
Starts a profiling section to record with the given name. This action is active until it is stopped by stopLastProfilingSection().
All active sections must be stopped before calling initialize() or advance().
[in] | sectionName | the name of the profiling section to start. The name may not contain forward slashes / . |
Definition at line 2032 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::stopLastProfilingSection | ( | ) |
Stop the last profiling section
Definition at line 2040 of file ParticipantImpl.cpp.
|
private |
Syncs the time step size between all ranks (all time steps sizes should be the same!)
Definition at line 1872 of file ParticipantImpl.cpp.
|
private |
Discards data before the given time for all meshes and data known by this participant.
Definition at line 573 of file ParticipantImpl.cpp.
|
private |
Removes samples in mapped to data connected to received data via a mapping.
This prevents old samples from blocking remappings.
Definition at line 1751 of file ParticipantImpl.cpp.
|
private |
Discards send (currently write) data of a participant after a given time when another iteration is required.
Definition at line 582 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::writeAndMapData | ( | std::string_view | meshName, |
std::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).
This function writes values at temporary locations to data of a mesh. As opposed to the writeData function using VertexIDs, this function allows to write data via coordinates, which don't have to be specified during the initialization. This is particularly useful for meshes, which vary over time. Note that using this function comes at a performance cost, since the specified mapping needs to be computed locally for the given locations, whereas the other variant (writeData) can typically exploit the static interface mesh and pre-compute data structures more efficiently.
Values are passed via a block of continuous memory defined by values in the order specified by vertices.
The 1D/Scalar-format of values is (d0, d1, ..., dn) The 2D-format of values is (d0x, d0y, d1x, d1y, ..., dnx, dny) The 3D-format of values is (d0x, d0y, d0z, d1x, d1y, d1z, ..., dnx, dny, dnz)
[in] | meshName | The name of the mesh that holds the data, needs to be a mesh received from another participant. |
[in] | dataName | The name of the data to write. |
[in] | coordinates | A span to the coordinates of the vertices The 2D-format is (d0x, d0y, d1x, d1y, ..., dnx, dny) The 3D-format is (d0x, d0y, d0z, d1x, d1y, d1z, ..., dnx, dny, dnz) |
[in] | values | The values containing the write data. |
Definition at line 1325 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::writeData | ( | std::string_view | meshName, |
std::string_view | dataName, | ||
::precice::span< const VertexID > | vertices, | ||
::precice::span< const double > | values ) |
Writes data to a mesh.
This function writes values of specified vertices to data of a mesh. Values are provided as a block of continuous memory defined by values. The order of the provided data follows the order specified by vertices.
The 1D/Scalar-format of values is (d0, d1, ..., dn) The 2D-format of values is (d0x, d0y, d1x, d1y, ..., dnx, dny) The 3D-format of values is (d0x, d0y, d0z, d1x, d1y, d1z, ..., dnx, dny, dnz)
[in] | meshName | the name of mesh that hold the data. |
[in] | dataName | the name of the data to write to. |
[in] | ids | the vertex ids of the vertices to write data to. |
[in] | values | the values to write to preCICE. |
ids
is a return value of setMeshVertex or setMeshVertices Definition at line 1174 of file ParticipantImpl.cpp.
void precice::impl::ParticipantImpl::writeGradientData | ( | std::string_view | meshName, |
std::string_view | dataName, | ||
::precice::span< const VertexID > | vertices, | ||
::precice::span< const double > | gradients ) |
Writes vector gradient data to a mesh.
This function writes gradient values of specified vertices to data of a mesh. Values are provided as a block of continuous memory defined by gradients. The order of the provided gradient data follows the order specified by vertices.
Each gradient or Jacobian depends on the dimensionality of the mesh and data. Each gradient has a total of getMeshDimensions(meshName) * getDataDimensions(meshName, dataName) components and is stored in a linearised format as follows:
Spatial Dimensions | Scalar Data | Vectorial Data |
---|---|---|
2D | s dx, s dy | x dx, y dx, x dy, y dy |
3D | s dy, s dy, s dz | x dx, y dx, z dx, x dy, y dy, z dy, x dz, y dz, z dz |
The gradients/Jacobian for all vertices are then contiguously saved in memory.
Example for 2D Vectorial:
Index | 0 | 1 | 2 | 3 | ... | 4n | 4n+1 | 4n+2 | 4n+3 |
---|---|---|---|---|---|---|---|---|---|
Component | x0 dx | y0 dx | x0 dy | y0 dy | ... | xn dx | yn dx | xn dy | yn dy |
[in] | meshName | the name of mesh that hold the data. |
[in] | dataName | the name of the data to write to. |
[in] | ids | the vertex ids of the vertices to write gradient data to. |
[in] | gradients | the linearised gradient data to write to preCICE. |
Definition at line 1381 of file ParticipantImpl.cpp.
|
friend |
Definition at line 518 of file ParticipantImpl.hpp.
|
friend |
To allow white box tests.
Definition at line 517 of file ParticipantImpl.hpp.
|
private |
Definition at line 322 of file ParticipantImpl.hpp.
|
private |
Definition at line 320 of file ParticipantImpl.hpp.
|
private |
Definition at line 316 of file ParticipantImpl.hpp.
|
private |
Definition at line 318 of file ParticipantImpl.hpp.
|
private |
Are experimental API calls allowed?
Definition at line 347 of file ParticipantImpl.hpp.
|
private |
Are experimental remeshing API calls allowed?
Definition at line 350 of file ParticipantImpl.hpp.
|
private |
The hash of the configuration file used to configure this participant.
Definition at line 368 of file ParticipantImpl.hpp.
|
private |
Definition at line 337 of file ParticipantImpl.hpp.
|
private |
Spatial dimensions of problem.
Definition at line 325 of file ParticipantImpl.hpp.
|
private |
Counts the amount of samples mapped in read mappings executed in the latest advance.
Definition at line 365 of file ParticipantImpl.hpp.
|
private |
Counts the amount of samples mapped in write mappings executed in the latest advance.
Definition at line 362 of file ParticipantImpl.hpp.
|
mutableprivate |
Definition at line 314 of file ParticipantImpl.hpp.
|
private |
Definition at line 332 of file ParticipantImpl.hpp.
|
private |
mesh name to mesh ID mapping.
Definition at line 330 of file ParticipantImpl.hpp.
|
private |
Definition at line 327 of file ParticipantImpl.hpp.
|
private |
Counts calls to advance for plotting.
Definition at line 359 of file ParticipantImpl.hpp.
|
private |
Holds information about solvers participating in the coupled simulation.
Definition at line 335 of file ParticipantImpl.hpp.
|
private |
Definition at line 521 of file ParticipantImpl.hpp.
|
private |
Definition at line 520 of file ParticipantImpl.hpp.
|
private |
The current State of the Participant.
Definition at line 356 of file ParticipantImpl.hpp.
|
private |
Definition at line 523 of file ParticipantImpl.hpp.
|
private |
Are participants waiting for each other in finalize?
Definition at line 353 of file ParticipantImpl.hpp.