preCICE v3.1.1
Loading...
Searching...
No Matches
Functions
preciceFortran.hpp File Reference
#include <precice/export.h>
Include dependency graph for preciceFortran.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

PRECICE_API void precicef_get_version_information_ (char *versionInfo, int lengthVersionInfo)
 
Construction and Configuration
PRECICE_API void precicef_create_ (const char *participantName, const char *configFileName, const int *solverProcessIndex, const int *solverProcessSize, int participantNameLength, int configFileNameLength)
 
Steering Methods
PRECICE_API void precicef_initialize_ ()
 
PRECICE_API void precicef_advance_ (const double *timeStepSize)
 
PRECICE_API void precicef_finalize_ ()
 
Implicit coupling
PRECICE_API void precicef_requires_reading_checkpoint_ (int *isRequired)
 
PRECICE_API void precicef_requires_writing_checkpoint_ (int *isRequired)
 
Status Queries
PRECICE_API void precicef_get_mesh_dimensions_ (const char *meshName, int *dimensions, int meshNameLength)
 
PRECICE_API void precicef_get_data_dimensions_ (const char *meshName, const char *dataName, int *dimensions, int meshNameLength, int dataNameLength)
 
PRECICE_API void precicef_is_coupling_ongoing_ (int *isOngoing)
 
PRECICE_API void precicef_is_time_window_complete_ (int *isComplete)
 
PRECICE_API void precicef_get_max_time_step_size_ (double *maxTimeStepSize)
 
Mesh Access

PRECICE_API void precicef_requires_mesh_connectivity_for_ (const char *meshName, int *required, int meshNameLength)
 
PRECICE_API void precicef_set_vertex_ (const char *meshName, const double *coordinates, int *id, int meshNameLength)
 
PRECICE_API void precicef_get_mesh_vertex_size_ (const char *meshName, int *meshSize, int meshNameLength)
 
PRECICE_API void precicef_set_vertices_ (const char *meshName, const int *size, double *coordinates, int *ids, int meshNameLength)
 
PRECICE_API void precicef_set_edge_ (const char *meshName, const int *firstVertexID, const int *secondVertexID, int meshNameLength)
 
PRECICE_API void precicef_set_mesh_edges_ (const char *meshName, const int *size, const int *ids, int meshNameLength)
 
PRECICE_API void precicef_set_triangle_ (const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, int meshNameLength)
 
PRECICE_API void precicef_set_quad_ (const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
 
PRECICE_API void precicef_set_mesh_quads_ (const char *meshName, const int *size, const int *ids, int meshNameLength)
 
PRECICE_API void precicef_set_tetrahedron (const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
 
PRECICE_API void precicef_set_mesh_tetrahedra_ (const char *meshName, const int *size, const int *ids, int meshNameLength)
 
Data Access
PRECICE_API void precicef_requires_initial_data_ (int *isRequired)
 
PRECICE_API void precicef_write_data_ (const char *meshName, const char *dataName, const int *size, int *ids, double *values, int meshNameLength, int dataNameLength)
 
PRECICE_API void precicef_read_data_ (const char *meshName, const char *dataName, const int *size, int *ids, const double *relativeReadTime, double *values, int meshNameLength, int dataNameLength)
 
Direct mesh access
PRECICE_API void precicef_set_mesh_access_region_ (const char *meshName, const double *boundingBox, int meshNameLength)
 
PRECICE_API void precicef_get_mesh_vertex_ids_and_coordinates_ (const char *meshName, const int size, int *ids, double *coordinates, int meshNameLength)
 
Experimental Gradient Data

These API functions are experimental and may change in future versions.

PRECICE_API void precicef_requires_gradient_data_for_ (const char *meshName, const char *dataName, int *required, int meshNameLength, int dataNameLength)
 
PRECICE_API void precicef_write_gradient_data_ (const char *meshName, const char *dataName, const int *size, const int *ids, const double *gradients, int meshNameLength, int dataNameLength)
 

Detailed Description

This file contains a Fortran 77 compatible interface written in C/C++. The Fortran bindings are thin wrappers around the C++ API.

Every method has a Fortran syntax equivalent in the method comment, and a listing for input and output variables. A variable can be input and output at the same time.

Definition in file preciceFortran.hpp.

Function Documentation

◆ precicef_advance_()

PRECICE_API void precicef_advance_ ( const double * timeStepSize)

Fortran syntax: precicef_advance( DOUBLE PRECISION timeStepSize )

IN: timeStepSize OUT: -

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 62 of file preciceFortran.cpp.

◆ precicef_create_()

PRECICE_API void precicef_create_ ( const char * participantName,
const char * configFileName,
const int * solverProcessIndex,
const int * solverProcessSize,
int participantNameLength,
int configFileNameLength )

Fortran syntax: precicef_create( CHARACTER participantName(*), CHARACTER configFileName(*), INTEGER solverProcessIndex, INTEGER solverProcessSize )

IN: participantName, configFileName, solverProcessIndex, solverProcessSize OUT: -

Constructs a Participant for the given participant.

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.

Definition at line 41 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_finalize_()

PRECICE_API void precicef_finalize_ ( )

Fortran syntax: precicef_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 69 of file preciceFortran.cpp.

◆ precicef_get_data_dimensions_()

PRECICE_API void precicef_get_data_dimensions_ ( const char * meshName,
const char * dataName,
int * dimensions,
int meshNameLength,
int dataNameLength )

Fortran syntax: precicef_get_data_dimensions_( CHARACTER meshName(*), CHARACTER dataName(*), INTEGER dimensions)

IN: mesh, data, meshNameLength, dataNameLength OUT: dimensions

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 85 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_get_max_time_step_size_()

PRECICE_API void precicef_get_max_time_step_size_ ( double * maxTimeStepSize)

Fortran syntax: precicef_get_max_time_step_size( DOUBLE PRECISION maxTimeStepSize )

IN: - OUT: maxTimeStepSize

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 118 of file preciceFortran.cpp.

◆ precicef_get_mesh_dimensions_()

PRECICE_API void precicef_get_mesh_dimensions_ ( const char * meshName,
int * dimensions,
int meshNameLength )

Fortran syntax: precicef_get_mesh_dimensions_( CHARACTER meshName(*), INTEGER dimensions)

IN: mesh, meshNameLength OUT: dimensions

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 76 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_get_mesh_vertex_ids_and_coordinates_()

PRECICE_API void precicef_get_mesh_vertex_ids_and_coordinates_ ( const char * meshName,
const int size,
int * ids,
double * coordinates,
int meshNameLength )

Fortran syntax: precicef_get_mesh_vertex_ids_and_coordinates_( CHARACTER meshName(*), INTEGER size, INTEGER ids(size), DOUBLE PRECISION coordinates(dim*size))

IN: mesh, size, meshNameLength OUT: ids, coordinates

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 396 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_get_mesh_vertex_size_()

PRECICE_API void precicef_get_mesh_vertex_size_ ( const char * meshName,
int * meshSize,
int meshNameLength )

Fortran syntax: precicef_get_mesh_vertex_size( CHARACTER meshName(*), INTEGER meshSize )

IN: mesh, meshNameLength OUT: meshSize

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 171 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_get_version_information_()

PRECICE_API void precicef_get_version_information_ ( char * versionInfo,
int lengthVersionInfo )

Definition at line 339 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_initialize_()

PRECICE_API void precicef_initialize_ ( )

Fortran syntax: precicef_initialize()

IN: - OUT: -

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 56 of file preciceFortran.cpp.

◆ precicef_is_coupling_ongoing_()

PRECICE_API void precicef_is_coupling_ongoing_ ( int * isOngoing)

Fortran syntax: precicef_is_coupling_ongoing( INTEGER isOngoing )

IN: - OUT: isOngoing(1:true, 0:false)

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 96 of file preciceFortran.cpp.

◆ precicef_is_time_window_complete_()

PRECICE_API void precicef_is_time_window_complete_ ( int * isComplete)

Fortran syntax: precicef_is_time_window_complete( INTEGER isComplete );

IN: - OUT: isComplete(1:true, 0:false)

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 107 of file preciceFortran.cpp.

◆ precicef_read_data_()

PRECICE_API void precicef_read_data_ ( const char * meshName,
const char * dataName,
const int * size,
int * ids,
const double * relativeReadTime,
double * values,
int meshNameLength,
int dataNameLength )

Fortran syntax: precicef_read_data( CHARACTER meshName(*), CHARACTER dataName(*), INTEGER size, INTEGER ids, DOUBLE PRECISION relativeReadTime, DOUBLE PRECISION values(dim*size) )

IN: mesh, data, size, ids, meshNameLength, dataNameLength OUT: values

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 301 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_requires_gradient_data_for_()

PRECICE_API void precicef_requires_gradient_data_for_ ( const char * meshName,
const char * dataName,
int * required,
int meshNameLength,
int dataNameLength )

Fortran syntax: precicef_requires_gradient_data_for_( CHARACTER meshName(*), CHARACTER dataName(*), INTEGER required )

IN: dataID OUT: required(1:true, 0:false)

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 350 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_requires_initial_data_()

PRECICE_API void precicef_requires_initial_data_ ( int * isRequired)

Fortran syntax: precicef_requires_initial_data_( INTEGER isRequired )

IN: - OUT: isRequired(1:true, 0:false)

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 125 of file preciceFortran.cpp.

◆ precicef_requires_mesh_connectivity_for_()

PRECICE_API void precicef_requires_mesh_connectivity_for_ ( const char * meshName,
int * required,
int meshNameLength )

Fortran syntax: precicef_precicef_requires_mesh_connectivity_for_( CHARACTER meshName(*), INTEGER required)

IN: mesh, meshNameLength OUT: required(1:true, 0:false)

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 146 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_requires_reading_checkpoint_()

PRECICE_API void precicef_requires_reading_checkpoint_ ( int * isRequired)

Fortran syntax: precicef_requires_reading_checkpoint_( INTEGER isRequired )

IN: - OUT: isRequired(1:true, 0:false)

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 139 of file preciceFortran.cpp.

◆ precicef_requires_writing_checkpoint_()

PRECICE_API void precicef_requires_writing_checkpoint_ ( int * isRequired)

Fortran syntax: precicef_requires_writing_checkpoint_( INTEGER isRequired )

IN: - OUT: isRequired(1:true, 0:false)

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 132 of file preciceFortran.cpp.

◆ precicef_set_edge_()

PRECICE_API void precicef_set_edge_ ( const char * meshName,
const int * firstVertexID,
const int * secondVertexID,
int meshNameLength )

Fortran syntax: precicef_set_edge( CHARACTER meshName(*), INTEGER firstVertexID, INTEGER secondVertexID )

IN: mesh, firstVertexID, secondVertexID, meshNameLength OUT: -

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 193 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_mesh_access_region_()

PRECICE_API void precicef_set_mesh_access_region_ ( const char * meshName,
const double * boundingBox,
int meshNameLength )

Fortran syntax: precicef_set_mesh_access_region_( CHARACTER meshName(*), DOUBLE PRECISION bounding_box(dim*2))

IN: mesh, bounding_box, meshNameLength OUT: -

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 385 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_mesh_edges_()

PRECICE_API void precicef_set_mesh_edges_ ( const char * meshName,
const int * size,
const int * ids,
int meshNameLength )

Fortran syntax: precicef_set_mesh_edges_( CHARACTER meshName(*), INTEGER size, INTEGER ids(size*2) )

IN: mesh, size, ids, meshNameLength OUT: -

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

Fortran syntax: precicef_set_mesh_triangles_( CHARACTER meshName(*), INTEGER size, INTEGER ids(size*3) )

IN: mesh, size, ids, meshNameLength OUT: -

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 203 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_mesh_quads_()

PRECICE_API void precicef_set_mesh_quads_ ( const char * meshName,
const int * size,
const int * ids,
int meshNameLength )

Fortran syntax: precicef_set_mesh_quads( CHARACTER meshName(*), INTEGER size, INTEGER ids(size*4) )

IN: mesh, size, ids, meshNameLength OUT: -

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 248 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_mesh_tetrahedra_()

PRECICE_API void precicef_set_mesh_tetrahedra_ ( const char * meshName,
const int * size,
const int * ids,
int meshNameLength )

Fortran syntax: precicef_set_mesh_tetrahedra_( CHARACTER meshName(*), INTEGER size, INTEGER ids(size*4) )

IN: mesh, size, ids, meshNameLength OUT: -

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 271 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_quad_()

PRECICE_API void precicef_set_quad_ ( const char * meshName,
const int * firstVertexID,
const int * secondVertexID,
const int * thirdVertexID,
const int * fourthVertexID,
int meshNameLength )

Fortran syntax: precicef_set_quad_( CHARACTER meshName(*), INTEGER firstVertexID, INTEGER secondVertexID, INTEGER thirdVertexID, INTEGER fourthVertexID )

IN: mesh, firstVertexID, secondVertexID, thirdVertexID, fourthVertexID, meshNameLength OUT: -

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 236 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_tetrahedron()

PRECICE_API void precicef_set_tetrahedron ( const char * meshName,
const int * firstVertexID,
const int * secondVertexID,
const int * thirdVertexID,
const int * fourthVertexID,
int meshNameLength )

Fortran syntax: precicef_set_tetrahedron( CHARACTER meshName(*), INTEGER firstVertexID, INTEGER secondVertexID, INTEGER thirdVertexID, INTEGER fourthVertexID )

IN: mesh, firstVertexID, secondVertexID, thirdVertexID, fourthVertexID, meshNameLength OUT: -

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 259 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_triangle_()

PRECICE_API void precicef_set_triangle_ ( const char * meshName,
const int * firstVertexID,
const int * secondVertexID,
const int * thirdVertexID,
int meshNameLength )

Fortran syntax: precicef_set_triangle_( CHARACTER meshName(*), INTEGER firstVertexID, INTEGER secondVertexID, INTEGER thirdVertexID )

IN: mesh, firstVertexID, secondVertexID, thirdVertexID, meshNameLength OUT: -

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 214 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_vertex_()

PRECICE_API void precicef_set_vertex_ ( const char * meshName,
const double * coordinates,
int * id,
int meshNameLength )

Fortran syntax: precicef_set_vertex( CHARACTER meshName(*), DOUBLE PRECISION coordinates(dim), INTEGER id )

IN: mesh, coordinates, meshNameLength OUT: id

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 159 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_set_vertices_()

PRECICE_API void precicef_set_vertices_ ( const char * meshName,
const int * size,
double * coordinates,
int * ids,
int meshNameLength )

Fortran syntax: precicef_set_vertices( CHARACTER meshName(*), INTEGER size, DOUBLE PRECISION coordinates(dim*size), INTEGER ids(size) )

IN: mesh, size, coordinates, meshNameLength OUT: 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 180 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_write_data_()

PRECICE_API void precicef_write_data_ ( const char * meshName,
const char * dataName,
const int * size,
int * ids,
double * values,
int meshNameLength,
int dataNameLength )

Fortran syntax: precicef_write_data( CHARACTER meshName(*), CHARACTER dataName(*), INTEGER size, INTEGER ids, DOUBLE PRECISION values(dim*size) )

IN: mesh, data, size, ids, values, meshNameLength, dataNameLength OUT: -

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 282 of file preciceFortran.cpp.

Here is the call graph for this function:

◆ precicef_write_gradient_data_()

PRECICE_API void precicef_write_gradient_data_ ( const char * meshName,
const char * dataName,
const int * size,
const int * ids,
const double * gradients,
int meshNameLength,
int dataNameLength )

Fortran syntax: precicef_write_gradient_data_( CHARACTER meshName(*), CHARACTER dataName(*), INTEGER size, INTEGER ids, DOUBLE PRECISION gradients )

IN: mesh, data, size, ids, gradients, meshNameLength, dataNameLength OUT: -

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 364 of file preciceFortran.cpp.

Here is the call graph for this function: