preCICE v3.1.2
Loading...
Searching...
No Matches
Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
precice::impl::ParticipantImpl Class Reference

Implementation of Participant. See also pimpl ideom (https://en.cppreference.com/w/cpp/language/pimpl). More...

#include <ParticipantImpl.hpp>

Collaboration diagram for precice::impl::ParticipantImpl:
[legend]

Classes

struct  MappedSamples
 

Public Member Functions

const mesh::Meshmesh (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.
 
ParticipantImploperator= (ParticipantImpl const &)=delete
 Disable assignment construction.
 
 ParticipantImpl (ParticipantImpl &&)=delete
 Disable move construction.
 
ParticipantImploperator= (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)
 
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 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.
 

Private Types

enum struct  State { Constructed , Initialized , Finalized }
 Represents the various states a Participant can be in. More...
 
enum class  CloseChannels : bool { All = false , Distributed = true }
 Which channels to close in closeCommunicationChannels() More...
 

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 ()
 Exports meshes with data and watch point data.
 
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.
 

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 _waitInFinalize = false
 Are participants waiting for each other in finalize?
 
bool _accessRegionDefined = false
 setMeshAccessRegion may only be called once
 
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::unique_ptr< profiling::Event_solverInitEvent
 
std::unique_ptr< profiling::Event_solverAdvanceEvent
 

Friends

struct Integration::Serial::Whitebox::TestConfigurationPeano
 To allow white box tests.
 
struct Integration::Serial::Whitebox::TestConfigurationComsol
 

Detailed Description

Implementation of Participant. See also pimpl ideom (https://en.cppreference.com/w/cpp/language/pimpl).

Definition at line 56 of file ParticipantImpl.hpp.

Member Enumeration Documentation

◆ CloseChannels

enum class precice::impl::ParticipantImpl::CloseChannels : bool
strongprivate

Which channels to close in closeCommunicationChannels()

Enumerator
All 
Distributed 

Definition at line 427 of file ParticipantImpl.hpp.

◆ State

enum struct precice::impl::ParticipantImpl::State
strongprivate

Represents the various states a Participant can be in.

Enumerator
Constructed 
Initialized 
Finalized 

Definition at line 316 of file ParticipantImpl.hpp.

Constructor & Destructor Documentation

◆ ParticipantImpl() [1/3]

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.

Parameters
[in]participantNameName of the participant using the interface. Has to match the name given for a participant in the xml configuration file.
[in]configurationFileNameName (with path) of the xml configuration file.
[in]solverProcessIndexIf 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]solverProcessSizeThe number of solver processes using preCICE.
[in]communicatorAn optional pointer to an MPI_Comm to use as communicator.

Definition at line 79 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ ~ParticipantImpl()

precice::impl::ParticipantImpl::~ParticipantImpl ( )

Destructor.

Ensures that finalize() has been called.

See also
finalize

Definition at line 160 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ ParticipantImpl() [2/3]

precice::impl::ParticipantImpl::ParticipantImpl ( ParticipantImpl const & )
delete

Disable copy construction.

◆ ParticipantImpl() [3/3]

precice::impl::ParticipantImpl::ParticipantImpl ( ParticipantImpl && )
delete

Disable move construction.

Member Function Documentation

◆ advance()

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:

  • Sending and resetting coupling data written by solver to coupling partners.
  • Receiving coupling data read by solver.
  • Computing and applying data mappings.
  • Computing convergence measures in implicit coupling schemes.
  • Computing acceleration of coupling data.
  • Exchanging and computing information regarding the state of the coupled simulation.

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.

Parameters
[in]computedTimeStepSizeSize of time step used by the solver.
Attention
All ranks of participants running in parallel must pass the same computedTimeStep to advance().
See also
getMaxTimeStepSize to get the maximum allowed value for computedTimeStepSize.
Precondition
initialize() has been called successfully.
The solver has computed one time step.
The solver has written all coupling data.
isCouplingOngoing() returns true.
finalize() has not yet been called.
Postcondition
Coupling data values specified in the configuration are exchanged.
Coupling scheme state (computed time, computed time steps, ...) is updated.
The coupling state is logged.
Configured data mapping schemes are applied.
[Second Participant] Configured acceleration schemes are applied.
Meshes with data are exported to files if configured.

Definition at line 347 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ advanceCouplingScheme()

void precice::impl::ParticipantImpl::advanceCouplingScheme ( )
private

Advances the coupling schemes.

Definition at line 1544 of file ParticipantImpl.cpp.

◆ closeCommunicationChannels()

void precice::impl::ParticipantImpl::closeCommunicationChannels ( CloseChannels cc)
private

Syncs the primary ranks of all connected participants.

Definition at line 1557 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ compareBoundingBoxes()

void precice::impl::ParticipantImpl::compareBoundingBoxes ( )
private

Communicate bounding boxes and look for overlaps.

Definition at line 1315 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ computeMappings()

void precice::impl::ParticipantImpl::computeMappings ( std::vector< MappingContext > & contexts,
const std::string & mappingType )
private

Helper for mapWrittenData and mapReadData.

Definition at line 1388 of file ParticipantImpl.cpp.

◆ computePartitions()

void precice::impl::ParticipantImpl::computePartitions ( )
private

Communicate meshes and create partitions.

Definition at line 1335 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ configure() [1/2]

void precice::impl::ParticipantImpl::configure ( const config::Configuration & configuration)
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.

◆ configure() [2/2]

void precice::impl::ParticipantImpl::configure ( std::string_view configurationFileName)
private

Configures the coupling interface from the given xml file.

Only after the configuration a reasonable state of a ParticipantImpl object is achieved.

Parameters
[in]configurationFileNameName (with path) of the xml config. file.

Definition at line 168 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ configureM2Ns()

void precice::impl::ParticipantImpl::configureM2Ns ( const m2n::M2NConfiguration::SharedPointer & config)
private

Definition at line 1227 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ configurePartitions()

void precice::impl::ParticipantImpl::configurePartitions ( const m2n::M2NConfiguration::SharedPointer & m2nConfig)
private

Determines participants providing meshes to other participants.

Definition at line 1265 of file ParticipantImpl.cpp.

◆ determineAccessingParticipant()

PtrParticipant precice::impl::ParticipantImpl::determineAccessingParticipant ( const config::Configuration & config)
private

Determines participant accessing this interface from the configuration.

Definition at line 1501 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ finalize()

void precice::impl::ParticipantImpl::finalize ( )

Finalizes preCICE.

If initialize() has been called:

  • Synchronizes with remote participants
  • handles final exports
  • cleans up general state

Always:

  • flushes and finalizes Events
  • finalizes managed PETSc
  • finalizes managed MPI
Precondition
finalize() has not been called.
Postcondition
Communication channels are closed.
Meshes and data are deallocated
Finalized managed PETSc
Finalized managed MPI
See also
isCouplingOngoing()

Definition at line 483 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ getDataDimensions()

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 spacial dimensionality of the mesh.

Parameters
[in]meshNamethe name of the associated mesh
[in]dataNamethe name of the data to get the dimensions for
Returns
the dimensions of the given Data
See also
getMeshDimensions

Definition at line 535 of file ParticipantImpl.cpp.

◆ getMaxTimeStepSize()

double precice::impl::ParticipantImpl::getMaxTimeStepSize ( ) const

Get the maximum allowed time step size of the current window.

Returns
Maximum size of time step to be computed by solver.

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.

Precondition
initialize() has been called successfully.

Definition at line 559 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ getMeshDimensions()

int precice::impl::ParticipantImpl::getMeshDimensions ( std::string_view meshName) const

Returns the spatial dimensionality of the given mesh.

Parameters
[in]meshNamethe name of the associated mesh
Returns
the dimensions of the given mesh

Definition at line 528 of file ParticipantImpl.cpp.

◆ getMeshVertexIDsAndCoordinates()

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.

Parameters
[in]meshNamecorresponding mesh name
[out]idsids corresponding to the coordinates
[out]coordinatesthe coordinates associated to the ids and corresponding data values
Precondition
ids.size() == getMeshVertexSize(meshName)
coordinates.size() == getMeshVertexSize(meshName) * getMeshDimensions(meshName)
This function can be called on received meshes as well as provided meshes. However, you need to call this function after initialize(), if the meshName corresponds to a received mesh, since the relevant mesh data is exchanged during the initialize() call.
See also
getMeshVertexSize() to get the amount of vertices in the mesh
getMeshDimensions() to get the spacial dimensionality of the mesh

Definition at line 1181 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ getMeshVertexSize()

int precice::impl::ParticipantImpl::getMeshVertexSize ( std::string_view meshName) const

Returns the number of vertices of a mesh.

Parameters
[in]meshNamethe name of the mesh
Returns
the amount of the vertices of the mesh
Precondition
This function can be called on received meshes as well as provided meshes. However, you need to call this function after initialize(), if the meshName corresponds to a received mesh, since the relevant mesh data is exchanged during the initialize() call.

Definition at line 628 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ handleDataAfterAdvance()

void precice::impl::ParticipantImpl::handleDataAfterAdvance ( bool reachedTimeWindowEnd,
bool isTimeWindowComplete,
double timeSteppedTo,
double timeAfterAdvance,
const cplscheme::ImplicitData & receivedData )
private

Completes everything data-related after advancing the coupling scheme.

Definition at line 420 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ handleDataBeforeAdvance()

void precice::impl::ParticipantImpl::handleDataBeforeAdvance ( bool reachedTimeWindowEnd,
double timeSteppedTo )
private

Completes everything data-related between adding time to and advancing the coupling scheme.

Definition at line 403 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ handleExports()

void precice::impl::ParticipantImpl::handleExports ( )
private

Exports meshes with data and watch point data.

Definition at line 1481 of file ParticipantImpl.cpp.

◆ initialize()

void precice::impl::ParticipantImpl::initialize ( )

Fully initializes preCICE and coupling data.

- Sets up a connection to the other participants of the coupled simulation.

  • Pre-processes defined meshes and handles partitions in parallel.
  • Receives first coupling data. The starting values for coupling data are zero by default.
  • Determines maximum allowed size of the first time step to be computed.
Precondition
initialize() has not yet been called.
Postcondition
Parallel communication to the coupling partner(s) is setup.
Meshes are exchanged between coupling partners and the parallel partitions are created.
Initial coupling data was exchanged.
See also
getMaxTimeStepSize()
requiresInitialData()

Definition at line 242 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ initializeIntraCommunication()

void precice::impl::ParticipantImpl::initializeIntraCommunication ( )
private

Initializes intra-participant communication.

Definition at line 1516 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ isCouplingOngoing()

bool precice::impl::ParticipantImpl::isCouplingOngoing ( ) const

Checks if the coupled simulation is still ongoing.

Returns
whether the coupling is ongoing.

A coupling is ongoing as long as

  • the maximum number of time windows has not been reached, and
  • the final time has not been reached.
Precondition
initialize() has been called successfully.
See also
advance()
Note
The user should call finalize() after this function returns false.

Definition at line 543 of file ParticipantImpl.cpp.

◆ isTimeWindowComplete()

bool precice::impl::ParticipantImpl::isTimeWindowComplete ( ) const

Checks if the current coupling window is completed.

Returns
whether the current time window is complete.

The following reasons require several solver time steps per time window:

  • A solver chooses to perform subcycling, i.e. using a smaller timestep than the time window.
  • An implicit coupling iteration is not yet converged.

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.

Precondition
initialize() has been called successfully.

Definition at line 551 of file ParticipantImpl.cpp.

◆ mapInitialReadData()

void precice::impl::ParticipantImpl::mapInitialReadData ( )
private

Definition at line 1445 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ mapInitialWrittenData()

void precice::impl::ParticipantImpl::mapInitialWrittenData ( )
private

Computes, and performs write mappings of the initial data in initialize.

Definition at line 1404 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ mappedSamples()

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 1598 of file ParticipantImpl.cpp.

◆ mapReadData()

void precice::impl::ParticipantImpl::mapReadData ( )
private

Definition at line 1458 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ mapWrittenData()

void precice::impl::ParticipantImpl::mapWrittenData ( std::optional< double > after = std::nullopt)
private

Computes, and performs suitable write mappings either entirely or after given time.

Definition at line 1416 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ mesh()

const mesh::Mesh & precice::impl::ParticipantImpl::mesh ( const std::string & meshName) const

Allows to access a registered mesh.

Todo
try to remove or make private. See https://github.com/precice/precice/issues/1269

Definition at line 1592 of file ParticipantImpl.cpp.

◆ operator=() [1/2]

ParticipantImpl & precice::impl::ParticipantImpl::operator= ( ParticipantImpl && )
delete

Disable move assignment.

◆ operator=() [2/2]

ParticipantImpl & precice::impl::ParticipantImpl::operator= ( ParticipantImpl const & )
delete

Disable assignment construction.

◆ performDataActions()

void precice::impl::ParticipantImpl::performDataActions ( const std::set< action::Action::Timing > & timings)
private

Performs all data actions with given timing.

Parameters
[in]timingsthe timings of the action.

Definition at line 1471 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ readData()

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.

Parameters
[in]meshNamethe name of mesh that hold the data.
[in]dataNamethe name of the data to read from.
[in]idsthe vertex ids of the vertices to read data from.
[in]relativeReadTimePoint in time where data is read relative to the beginning of the current time step.
[out]valuesthe destination memory to read the data from.
Precondition
every VertexID in ids is a return value of setMeshVertex or setMeshVertices
values.size() == getDataDimensions(meshName, dataName) * ids.size()
Postcondition
values contain the read data as specified in the above format.
See also
Participant::setMeshVertex()
Participant::setMeshVertices()
Participant::getDataDimensions()

Definition at line 1057 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ requiresGradientDataFor()

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.

Attention
This API function is experimental and may change in the future!

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)

Parameters
[in]meshNamethe name of mesh that hold the data.
[in]dataNamethe name of the data.
Returns
whether gradient is required

Definition at line 616 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ requiresInitialData()

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().

Note
If initial data is configured, then this function needs to be called.
Precondition
initialize() has not yet been called

Definition at line 576 of file ParticipantImpl.cpp.

◆ requiresMeshConnectivityFor()

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.

Parameters
[in]meshNamethe name of the mesh
Returns
whether connectivity is required

Definition at line 609 of file ParticipantImpl.cpp.

◆ requiresReadingCheckpoint()

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().

Note
If implicit coupling is configured for this Participant, then this function needs to be called.
This function returns false before the first call to advance().
Precondition
initialize() has been called
See also
requiresWritingCheckpoint()

Definition at line 598 of file ParticipantImpl.cpp.

◆ requiresWritingCheckpoint()

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().

Note
If implicit coupling is configured for this Participant, then this function needs to be called.
Precondition
initialize() has been called
See also
requiresReadingCheckpoint()

Definition at line 587 of file ParticipantImpl.cpp.

◆ resetMesh()

void precice::impl::ParticipantImpl::resetMesh ( std::string_view meshName)

Todo
Currently not supported as we would need to re-compute the re-partition

Definition at line 644 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ resetWrittenData()

void precice::impl::ParticipantImpl::resetWrittenData ( )
private

Resets written data.

Parameters
isAtWindowEndset true, if function is called at end of window to also trim the time sample storage
isTimeWindowCompleteset true, if function is called at end of converged window to trim and move the sample storage.

Definition at line 1493 of file ParticipantImpl.cpp.

◆ samplizeWriteData()

void precice::impl::ParticipantImpl::samplizeWriteData ( double time)
private

Creates a Stample at the given time for each write Data and zeros the buffers.

Definition at line 459 of file ParticipantImpl.cpp.

◆ setMeshAccessRegion()

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].

Note
Defining a bounding box for serial runs of the solver (not to be confused with serial coupling mode) is valid. However, a warning is raised in case vertices are filtered out completely on the receiving side, since the associated data values of the filtered vertices are filled with zero data.
This function can only be called once per participant and rank and trying to call it more than once results in an error.
If you combine the direct access with a mpping (say you want to read data from a defined mesh, as usual, but you want to directly access and write data on a received mesh without a mapping) you may not need this function at all since the region of interest is already defined through the defined mesh used for data reading. This is the case if you define any mapping involving the directly accessed mesh on the receiving participant. (In parallel, only the cases read-consistent and write-conservative are relevant, as usual).
The safety factor scaling (see safety-factor in the configuration file) is not applied to the defined access region and a specified safety will be ignored in case there is no additional mapping involved. However, in case a mapping is in addition to the direct access involved, you will receive (and gain access to) vertices inside the defined access region plus vertices inside the safety factor region resulting from the mapping. The default value of the safety factor is 0.5,i.e., the defined access region as computed through the involved provided mesh is by 50% enlarged.
Parameters
[in]meshNamename of the mesh you want to access through the bounding box
[in]boundingBoxAxis aligned bounding boxes which has in 3D the format [x_min, x_max, y_min, y_max, z_min, z_max]
Precondition
initialize() has not yet been called.
boundingBox.size() == 2 * getMeshDimensions(meshName)

Definition at line 1144 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshEdge()

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.

Note
The order of vertices does not matter.
Parameters
[in]meshNamename of the mesh to add the edge to
[in]firstID of the first vertex of the edge
[in]secondID of the second vertex of the edge
Precondition
vertices with IDs first and second were added to the mesh with the name meshName

Definition at line 712 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshEdges()

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.

Note
The order of vertices per edge does not matter.
Parameters
[in]meshNamethe name of the mesh to add the n edges to
[in]idsan array containing 2n vertex IDs for n edges
Precondition
vertices in ids were added to the mesh with the name meshName
ids.size() is multiple of 2
See also
requiresMeshConnectivityFor()

Definition at line 731 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshQuad()

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.

Warning
The order of vertices does not matter, however, only planar quads are allowed.
Parameters
[in]meshNamename of the mesh to add the Quad to
[in]firstID of the first vertex of the Quad
[in]secondID of the second vertex of the Quad
[in]thirdID of the third vertex of the Quad
[in]fourthID of the fourth vertex of the Quad
Precondition
vertices with IDs first, second, third, and fourth were added to the mesh with the name meshName
See also
requiresMeshConnectivityFor()

Definition at line 839 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshQuads()

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.

Warning
The order of vertices per quad does not matter, however, only planar quads are allowed.
Parameters
[in]meshNamename of the mesh to add the n quads to
[in]idsan array containing 4n vertex IDs for n quads
Precondition
vertices in ids were added to the mesh with the name meshName
ids.size() is multiple of 4
See also
requiresMeshConnectivityFor()

Definition at line 891 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshTetrahedra()

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.

Note
The order of vertices per tetrahedron does not matter.
Parameters
[in]meshNamename of the mesh to add the n tetrahedra to
[in]idsan array containing 4n vertex IDs for n tetrahedra
Precondition
vertices in ids were added to the mesh with the name meshName
ids.size() is multiple of 4
See also
requiresMeshConnectivityFor()

Definition at line 982 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshTetrahedron()

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.

Note
The order of vertices does not matter.
Parameters
[in]meshNamename of the mesh to add the Tetrahedron to
[in]firstID of the first vertex of the Tetrahedron
[in]secondID of the second vertex of the Tetrahedron
[in]thirdID of the third vertex of the Tetrahedron
[in]fourthID of the fourth vertex of the Tetrahedron
Precondition
vertices with IDs first, second, third, and fourth were added to the mesh with the name meshName
See also
requiresMeshConnectivityFor()

Definition at line 954 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshTriangle()

void precice::impl::ParticipantImpl::setMeshTriangle ( std::string_view meshName,
VertexID first,
VertexID second,
VertexID third )

Sets mesh triangle from vertex IDs.

Note
The order of vertices does not matter.
Parameters
[in]meshNamename of the mesh to add the triangle to
[in]firstID of the first vertex of the triangle
[in]secondID of the second vertex of the triangle
[in]thirdID of the third vertex of the triangle
Precondition
vertices with IDs first, second, and third were added to the mesh with the name meshName
See also
requiresMeshConnectivityFor()

Definition at line 765 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshTriangles()

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.

Note
The order of vertices per triangle does not matter.
Parameters
[in]meshNamename of the mesh to add the n triangles to
[in]idsan array containing 3n vertex IDs for n triangles
Precondition
vertices in ids were added to the mesh with the name meshName
ids.size() is multiple of 3
See also
requiresMeshConnectivityFor()

Definition at line 802 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshVertex()

VertexID precice::impl::ParticipantImpl::setMeshVertex ( std::string_view meshName,
::precice::span< const double > position )

Creates a mesh vertex.

Parameters
[in]meshNamethe name of the mesh to add the vertex to.
[in]positionthe coordinates of the vertex.
Returns
the id of the created vertex
Precondition
initialize() has not yet been called
position.size() == getMeshDimensions(meshName)
See also
getMeshDimensions()

Definition at line 657 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ setMeshVertices()

void precice::impl::ParticipantImpl::setMeshVertices ( std::string_view meshName,
::precice::span< const double > positions,
::precice::span< VertexID > ids )

Creates multiple mesh vertices.

Parameters
[in]meshNamethe name of the mesh to add the vertices to.
[in]coordinatesa 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]idsThe ids of the created vertices
Precondition
initialize() has not yet been called
coordinates.size() == getMeshDimensions(meshName) * ids.size()
See also
getDimensions()

Definition at line 680 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ syncTimestep()

void precice::impl::ParticipantImpl::syncTimestep ( double computedTimeStepSize)
private

Syncs the time step size between all ranks (all time steps sizes should be the same!)

Definition at line 1527 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ trimOldDataBefore()

void precice::impl::ParticipantImpl::trimOldDataBefore ( double time)
private

Discards data before the given time for all meshes and data known by this participant.

Definition at line 467 of file ParticipantImpl.cpp.

◆ trimReadMappedData()

void precice::impl::ParticipantImpl::trimReadMappedData ( double timeAfterAdvance,
bool isTimeWindowComplete,
const cplscheme::ImplicitData & fromData )
private

Removes samples in mapped to data connected to received data via a mapping.

This prevents old samples from blocking remappings.

Definition at line 1428 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ trimSendDataAfter()

void precice::impl::ParticipantImpl::trimSendDataAfter ( double time)
private

Discards send (currently write) data of a participant after a given time when another iteration is required.

Definition at line 476 of file ParticipantImpl.cpp.

◆ writeData()

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)

Parameters
[in]meshNamethe name of mesh that hold the data.
[in]dataNamethe name of the data to write to.
[in]idsthe vertex ids of the vertices to write data to.
[in]valuesthe values to write to preCICE.
Precondition
every VertexID in ids is a return value of setMeshVertex or setMeshVertices
values.size() == getDataDimensions(meshName, dataName) * ids.size()
See also
Participant::setMeshVertex()
Participant::setMeshVertices()
Participant::getDataDimensions()

Definition at line 1021 of file ParticipantImpl.cpp.

Here is the call graph for this function:

◆ writeGradientData()

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.

Attention
This API function is experimental and may change in the future!

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
Parameters
[in]meshNamethe name of mesh that hold the data.
[in]dataNamethe name of the data to write to.
[in]idsthe vertex ids of the vertices to write gradient data to.
[in]gradientsthe linearised gradient data to write to preCICE.
Precondition
Data has attribute hasGradient = true
every VertexID in vertices is a return value of setMeshVertex or setMeshVertices
gradients.size() == ids.size() * getMeshDimensions(meshName) * getDataDimensions(meshName, dataName)
See also
Participant::setMeshVertex()
Participant::setMeshVertices()
Participant::getMeshDimensions()
Participant::getDataDimensions()

Definition at line 1097 of file ParticipantImpl.cpp.

Here is the call graph for this function:

Friends And Related Symbol Documentation

◆ Integration::Serial::Whitebox::TestConfigurationComsol

friend struct Integration::Serial::Whitebox::TestConfigurationComsol
friend

Definition at line 452 of file ParticipantImpl.hpp.

◆ Integration::Serial::Whitebox::TestConfigurationPeano

friend struct Integration::Serial::Whitebox::TestConfigurationPeano
friend

To allow white box tests.

Definition at line 451 of file ParticipantImpl.hpp.

Member Data Documentation

◆ _accessor

impl::PtrParticipant precice::impl::ParticipantImpl::_accessor
private

Definition at line 298 of file ParticipantImpl.hpp.

◆ _accessorCommunicatorSize

int precice::impl::ParticipantImpl::_accessorCommunicatorSize
private

Definition at line 296 of file ParticipantImpl.hpp.

◆ _accessorName

std::string precice::impl::ParticipantImpl::_accessorName
private

Definition at line 292 of file ParticipantImpl.hpp.

◆ _accessorProcessRank

int precice::impl::ParticipantImpl::_accessorProcessRank
private

Definition at line 294 of file ParticipantImpl.hpp.

◆ _accessRegionDefined

bool precice::impl::ParticipantImpl::_accessRegionDefined = false
mutableprivate

setMeshAccessRegion may only be called once

Definition at line 329 of file ParticipantImpl.hpp.

◆ _allowsExperimental

bool precice::impl::ParticipantImpl::_allowsExperimental = false
private

Are experimental API calls allowed?

Definition at line 323 of file ParticipantImpl.hpp.

◆ _couplingScheme

cplscheme::PtrCouplingScheme precice::impl::ParticipantImpl::_couplingScheme
private

Definition at line 313 of file ParticipantImpl.hpp.

◆ _dimensions

int precice::impl::ParticipantImpl::_dimensions = 0
private

Spatial dimensions of problem.

Definition at line 301 of file ParticipantImpl.hpp.

◆ _executedReadMappings

int precice::impl::ParticipantImpl::_executedReadMappings = 0
private

Counts the amount of samples mapped in read mappings executed in the latest advance.

Definition at line 341 of file ParticipantImpl.hpp.

◆ _executedWriteMappings

int precice::impl::ParticipantImpl::_executedWriteMappings = 0
private

Counts the amount of samples mapped in write mappings executed in the latest advance.

Definition at line 338 of file ParticipantImpl.hpp.

◆ _log

logging::Logger precice::impl::ParticipantImpl::_log {"impl::ParticipantImpl"}
mutableprivate

Definition at line 290 of file ParticipantImpl.hpp.

◆ _m2ns

std::map<std::string, m2n::BoundM2N> precice::impl::ParticipantImpl::_m2ns
private

Definition at line 308 of file ParticipantImpl.hpp.

◆ _meshIDs

std::map<std::string, int> precice::impl::ParticipantImpl::_meshIDs
private

mesh name to mesh ID mapping.

Definition at line 306 of file ParticipantImpl.hpp.

◆ _meshLock

utils::MultiLock<std::string> precice::impl::ParticipantImpl::_meshLock
private

Definition at line 303 of file ParticipantImpl.hpp.

◆ _numberAdvanceCalls

long int precice::impl::ParticipantImpl::_numberAdvanceCalls = 0
private

Counts calls to advance for plotting.

Definition at line 335 of file ParticipantImpl.hpp.

◆ _participants

std::vector<impl::PtrParticipant> precice::impl::ParticipantImpl::_participants
private

Holds information about solvers participating in the coupled simulation.

Definition at line 311 of file ParticipantImpl.hpp.

◆ _solverAdvanceEvent

std::unique_ptr<profiling::Event> precice::impl::ParticipantImpl::_solverAdvanceEvent
private

Definition at line 455 of file ParticipantImpl.hpp.

◆ _solverInitEvent

std::unique_ptr<profiling::Event> precice::impl::ParticipantImpl::_solverInitEvent
private

Definition at line 454 of file ParticipantImpl.hpp.

◆ _state

State precice::impl::ParticipantImpl::_state {State::Constructed}
private

The current State of the Participant.

Definition at line 332 of file ParticipantImpl.hpp.

◆ _waitInFinalize

bool precice::impl::ParticipantImpl::_waitInFinalize = false
private

Are participants waiting for each other in finalize?

Definition at line 326 of file ParticipantImpl.hpp.


The documentation for this class was generated from the following files: