13#include "precice/impl/versions.hpp"
18#pragma GCC diagnostic push
19#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
22#pragma clang diagnostic push
23#pragma clang diagnostic ignored "-Wdeprecated-declarations"
32static std::string errormsg =
"preCICE has not been created properly. Be sure to call \"precicef_create\" before any other call to preCICE.";
45 const char *participantName,
46 const char *configFileName,
47 const int *solverProcessIndex,
48 const int *solverProcessSize,
49 int participantNameLength,
50 int configFileNameLength)
62 const char *participantName,
63 const char *configFileName,
64 const int *solverProcessIndex,
65 const int *solverProcessSize,
66 const int *communicator,
67 int participantNameLength,
68 int configFileNameLength)
75 MPI_Fint f_communicator =
static_cast<MPI_Fint
>(*communicator);
76 auto c_communicator = MPI_Comm_f2c(f_communicator);
84 PRECICE_WARN(
"preCICE was configured without MPI but you passed an MPI communicator. preCICE ignores the communicator and continues.");
104 const double *timeStepSize)
107 impl->advance(*timeStepSize);
122 const char *meshName,
133 const char *meshName,
134 const char *dataName,
149 if (
impl->isCouplingOngoing()) {
162 if (
impl->isTimeWindowComplete()) {
172 double *maxTimeStepSize)
175 *maxTimeStepSize =
impl->getMaxTimeStepSize();
184 *isRequired =
impl->requiresInitialData() ? 1 : 0;
193 *isRequired =
impl->requiresWritingCheckpoint() ? 1 : 0;
202 *isRequired =
impl->requiresReadingCheckpoint() ? 1 : 0;
208 const char *meshName,
223 const char *meshName,
233 const char *meshName,
234 const double *position,
240 auto positionSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv));
241 *vertexID =
impl->setMeshVertex(sv, {position, positionSize});
247 const char *meshName,
258 const char *meshName,
266 auto positionSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv) * *
size);
267 impl->setMeshVertices(sv, {coordinates, positionSize}, {ids,
static_cast<unsigned long>(*size)});
273 const char *meshName,
274 const int *firstVertexID,
275 const int *secondVertexID,
285 const char *meshName,
291 auto idsSize =
static_cast<unsigned long>(*size) * 2;
298 const char *meshName,
299 const int *firstVertexID,
300 const int *secondVertexID,
301 const int *thirdVertexID,
311 const char *meshName,
317 auto idsSize =
static_cast<unsigned long>(*size) * 3;
324 const char *meshName,
325 const int *firstVertexID,
326 const int *secondVertexID,
327 const int *thirdVertexID,
328 const int *fourthVertexID,
338 const char *meshName,
344 auto idsSize =
static_cast<unsigned long>(*size) * 4;
351 const char *meshName,
352 const int *firstVertexID,
353 const int *secondVertexID,
354 const int *thirdVertexID,
355 const int *fourthVertexID,
365 const char *meshName,
371 auto idsSize =
static_cast<unsigned long>(*size) * 4;
378 const char *meshName,
379 const char *dataName,
389 auto dataSize = *
size *
impl->getDataDimensions(strippedMeshName, strippedDataName);
390 impl->writeData(strippedMeshName,
392 {ids,
static_cast<unsigned long>(*size)},
393 {values,
static_cast<unsigned long>(dataSize)});
399 const char *meshName,
400 const char *dataName,
403 const double *relativeReadTime,
411 auto dataSize = *
size *
impl->getDataDimensions(strippedMeshName, strippedDataName);
415 {ids,
static_cast<unsigned long>(*size)},
417 {values,
static_cast<unsigned long>(dataSize)});
427 while (((
string[i] ==
' ') || (
string[i] == 0)) && (i >= 0)) {
440 int lengthVersionInfo)
442 std::strncpy(versionInfo, precice::versionInformation, lengthVersionInfo);
446 const char *meshName,
447 const char *dataName,
int *required,
462 const char *meshName,
463 const char *dataName,
466 const double *gradients,
473 auto gradientComponents =
impl->getMeshDimensions(strippedMeshName) *
impl->getDataDimensions(strippedMeshName, strippedDataName);
474 auto gradientSize = *
size * gradientComponents;
476 impl->writeGradientData(strippedMeshName,
478 {ids,
static_cast<unsigned long>(*size)},
479 {gradients,
static_cast<unsigned long>(gradientSize)});
485 const char *meshName,
486 const char *dataName,
496 auto coordinatesSize = *
size *
impl->getMeshDimensions(strippedMeshName);
497 auto dataSize = *
size *
impl->getDataDimensions(strippedMeshName, strippedDataName);
498 impl->writeAndMapData(strippedMeshName,
500 {coordinates,
static_cast<unsigned long>(coordinatesSize)},
501 {values,
static_cast<unsigned long>(dataSize)});
507 const char *meshName,
508 const char *dataName,
511 const double *relativeReadTime,
519 auto coordinatesSize = *
size *
impl->getMeshDimensions(strippedMeshName);
520 auto dataSize = *
size *
impl->getDataDimensions(strippedMeshName, strippedDataName);
521 impl->mapAndReadData(
524 {coordinates,
static_cast<unsigned long>(coordinatesSize)},
526 {values,
static_cast<unsigned long>(dataSize)});
532 const char *meshName,
533 const double *boundingBox,
538 auto bbSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv) * 2);
539 impl->setMeshAccessRegion(sv, {boundingBox, bbSize});
545 const char *meshName,
553 auto coordinatesSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv) * *
size);
554 impl->getMeshVertexIDsAndCoordinates(sv, {ids,
static_cast<unsigned long>(*size)}, {coordinates, coordinatesSize});
560 const char *sectionName,
561 int sectionNameLength)
565 impl->startProfilingSection(sv);
573 impl->stopLastProfilingSection();
579#pragma GCC diagnostic pop
582#pragma clang diagnostic pop
#define PRECICE_WARN(...)
#define PRECICE_CHECK(check,...)
This class provides a lightweight logger.
int strippedLength(const char *string, int length)
Returns length of string without trailing whitespace.
std::string_view strippedStringView(const char *string, int length)
static std::unique_ptr< precice::Participant > impl
static std::string errormsg
void precicef_advance_(const double *timeStepSize)
void precicef_requires_initial_data_(int *isRequired)
void precicef_get_version_information_(char *versionInfo, int lengthVersionInfo)
void precicef_start_profiling_section_(const char *sectionName, int sectionNameLength)
static precice::logging::Logger _log("preciceFortran")
void precicef_requires_gradient_data_for_(const char *meshName, const char *dataName, int *required, int meshNameLength, int dataNameLength)
void precicef_read_data_(const char *meshName, const char *dataName, const int *size, int *ids, const double *relativeReadTime, double *values, int meshNameLength, int dataNameLength)
void precicef_get_mesh_vertex_size_(const char *meshName, int *meshSize, int meshNameLength)
void precicef_get_mesh_dimensions_(const char *meshName, int *dimensions, int meshNameLength)
void precicef_create_(const char *participantName, const char *configFileName, const int *solverProcessIndex, const int *solverProcessSize, int participantNameLength, int configFileNameLength)
void precicef_set_vertex_(const char *meshName, const double *position, int *vertexID, int meshNameLength)
void precicef_reset_mesh_(const char *meshName, int meshNameLength)
void precicef_set_edge_(const char *meshName, const int *firstVertexID, const int *secondVertexID, int meshNameLength)
void precicef_is_time_window_complete_(int *isComplete)
void precicef_requires_mesh_connectivity_for_(const char *meshName, int *required, int meshNameLength)
void precicef_map_and_read_data_(const char *meshName, const char *dataName, const int *size, double *coordinates, const double *relativeReadTime, double *values, int meshNameLength, int dataNameLength)
Reads data using just-in-time data mapping. See.
void precicef_set_tetrahedron(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
void precicef_create_with_communicator_(const char *participantName, const char *configFileName, const int *solverProcessIndex, const int *solverProcessSize, const int *communicator, int participantNameLength, int configFileNameLength)
void precicef_set_mesh_quads_(const char *meshName, const int *size, const int *ids, int meshNameLength)
void precicef_get_data_dimensions_(const char *meshName, const char *dataName, int *dimensions, int meshNameLength, int dataNameLength)
void precicef_write_gradient_data_(const char *meshName, const char *dataName, const int *size, const int *ids, const double *gradients, int meshNameLength, int dataNameLength)
void precicef_set_quad_(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
void precicef_write_and_map_data_(const char *meshName, const char *dataName, const int *size, double *coordinates, double *values, int meshNameLength, int dataNameLength)
Writes data using just-in-time data mapping See.
void precicef_requires_reading_checkpoint_(int *isRequired)
void precicef_initialize_()
void precicef_requires_writing_checkpoint_(int *isRequired)
void precicef_finalize_()
void precicef_set_mesh_tetrahedra_(const char *meshName, const int *size, const int *ids, int meshNameLength)
void precicef_set_vertices_(const char *meshName, const int *size, double *coordinates, int *ids, int meshNameLength)
void precicef_set_triangle_(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, int meshNameLength)
void precicef_stop_last_profiling_section_()
void precicef_get_mesh_vertex_ids_and_coordinates_(const char *meshName, const int *size, int *ids, double *coordinates, int meshNameLength)
void precicef_write_data_(const char *meshName, const char *dataName, const int *size, int *ids, double *values, int meshNameLength, int dataNameLength)
void precicef_set_mesh_edges_(const char *meshName, const int *size, const int *ids, int meshNameLength)
void precicef_set_mesh_access_region_(const char *meshName, const double *boundingBox, int meshNameLength)
void precicef_get_max_time_step_size_(double *maxTimeStepSize)
void precicef_set_mesh_triangles_(const char *meshName, const int *size, const int *ids, int meshNameLength)
void precicef_is_coupling_ongoing_(int *isOngoing)