10#include "precice/impl/versions.hpp"
15#pragma GCC diagnostic push
16#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
19#pragma clang diagnostic push
20#pragma clang diagnostic ignored "-Wdeprecated-declarations"
29static std::string errormsg =
"preCICE has not been created properly. Be sure to call \"precicef_create\" before any other call to preCICE.";
42 const char *participantName,
43 const char *configFileName,
44 const int * solverProcessIndex,
45 const int * solverProcessSize,
46 int participantNameLength,
47 int configFileNameLength)
63 const double *timeStepSize)
66 impl->advance(*timeStepSize);
100 if (
impl->isCouplingOngoing()) {
111 if (
impl->isTimeWindowComplete()) {
119 double *maxTimeStepSize)
122 *maxTimeStepSize =
impl->getMaxTimeStepSize();
129 *isRequired =
impl->requiresInitialData() ? 1 : 0;
136 *isRequired =
impl->requiresWritingCheckpoint() ? 1 : 0;
143 *isRequired =
impl->requiresReadingCheckpoint() ? 1 : 0;
147 const char *meshName,
160 const char * meshName,
161 const double *position,
167 auto positionSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv));
168 *vertexID =
impl->setMeshVertex(sv, {position, positionSize});
172 const char *meshName,
181 const char *meshName,
183 double * coordinates,
189 auto positionSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv) * *
size);
190 impl->setMeshVertices(sv, {coordinates, positionSize}, {ids,
static_cast<unsigned long>(*size)});
194 const char *meshName,
195 const int * firstVertexID,
196 const int * secondVertexID,
204 const char *meshName,
210 auto idsSize =
static_cast<unsigned long>(*size) * 2;
215 const char *meshName,
216 const int * firstVertexID,
217 const int * secondVertexID,
218 const int * thirdVertexID,
226 const char *meshName,
232 auto idsSize =
static_cast<unsigned long>(*size) * 3;
237 const char *meshName,
238 const int * firstVertexID,
239 const int * secondVertexID,
240 const int * thirdVertexID,
241 const int * fourthVertexID,
249 const char *meshName,
255 auto idsSize =
static_cast<unsigned long>(*size) * 4;
260 const char *meshName,
261 const int * firstVertexID,
262 const int * secondVertexID,
263 const int * thirdVertexID,
264 const int * fourthVertexID,
272 const char *meshName,
278 auto idsSize =
static_cast<unsigned long>(*size) * 4;
283 const char *meshName,
284 const char *dataName,
294 auto dataSize = *
size *
impl->getDataDimensions(strippedMeshName, strippedDataName);
295 impl->writeData(strippedMeshName,
297 {ids,
static_cast<unsigned long>(*size)},
298 {values,
static_cast<unsigned long>(dataSize)});
302 const char * meshName,
303 const char * dataName,
306 const double *relativeReadTime,
314 auto dataSize = *
size *
impl->getDataDimensions(strippedMeshName, strippedDataName);
318 {ids,
static_cast<unsigned long>(*size)},
320 {values,
static_cast<unsigned long>(dataSize)});
328 while (((
string[i] ==
' ') || (
string[i] == 0)) && (i >= 0)) {
341 int lengthVersionInfo)
343 const std::string &versionInformation = precice::versionInformation;
344 PRECICE_ASSERT(versionInformation.
size() < (
size_t) lengthVersionInfo, versionInformation.
size(), lengthVersionInfo);
345 for (
size_t i = 0; i < versionInformation.
size(); i++) {
346 versionInfo[i] = versionInformation[i];
351 const char *meshName,
352 const char *dataName,
int *required,
365 const char * meshName,
366 const char * dataName,
369 const double *gradients,
376 auto gradientComponents =
impl->getMeshDimensions(strippedMeshName) *
impl->getDataDimensions(strippedMeshName, strippedDataName);
377 auto gradientSize = *
size * gradientComponents;
379 impl->writeGradientData(strippedMeshName,
381 {ids,
static_cast<unsigned long>(*size)},
382 {gradients,
static_cast<unsigned long>(gradientSize)});
386 const char * meshName,
387 const double *boundingBox,
392 auto bbSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv) * 2);
393 impl->setMeshAccessRegion(sv, {boundingBox, bbSize});
397 const char *meshName,
400 double * coordinates,
405 auto coordinatesSize =
static_cast<unsigned long>(
impl->getMeshDimensions(sv) *
size);
406 impl->getMeshVertexIDsAndCoordinates(sv, {ids,
static_cast<unsigned long>(
size)}, {coordinates, coordinatesSize});
410#pragma GCC diagnostic pop
413#pragma clang diagnostic pop
#define PRECICE_CHECK(check,...)
#define PRECICE_ASSERT(...)
Main Application Programming Interface of preCICE. Include using #include <precice/precice....
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)
void precicef_get_mesh_vertex_ids_and_coordinates_(const char *meshName, const int size, int *ids, double *coordinates, int meshNameLength)
void precicef_advance_(const double *timeStepSize)
void precicef_requires_initial_data_(int *isRequired)
void precicef_get_version_information_(char *versionInfo, int lengthVersionInfo)
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_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)
static std::unique_ptr< precice::Participant > impl
void precicef_set_tetrahedron(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
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)
static std::string errormsg
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_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_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)