preCICE v3.3.0
Loading...
Searching...
No Matches
preciceFortran.cpp
Go to the documentation of this file.
2#include <cstddef>
3#include <iostream>
4#include <memory>
5#ifndef PRECICE_NO_MPI
6#include <mpi.h>
7#endif
8#include <stdexcept>
9#include <string>
10#include <string_view>
11#include "logging/LogMacros.hpp"
12#include "logging/Logger.hpp"
13#include "precice/impl/versions.hpp"
14#include "precice/precice.hpp"
15#include "utils/assertion.hpp"
16
17#ifdef __GNUC__
18#pragma GCC diagnostic push
19#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
20#endif
21#ifdef __clang__
22#pragma clang diagnostic push
23#pragma clang diagnostic ignored "-Wdeprecated-declarations"
24#endif
25
26using namespace std;
27
29
30static precice::logging::Logger _log("preciceFortran");
31
32static std::string errormsg = "preCICE has not been created properly. Be sure to call \"precicef_create\" before any other call to preCICE.";
33
34namespace precice::impl {
38int strippedLength(const char *string, int length);
39
40std::string_view strippedStringView(const char *string, int length);
41
42} // namespace precice::impl
43
45 const char *participantName,
46 const char *configFileName,
47 const int *solverProcessIndex,
48 const int *solverProcessSize,
49 int participantNameLength,
50 int configFileNameLength)
51try {
53 precice::impl::strippedStringView(participantName, participantNameLength),
54 precice::impl::strippedStringView(configFileName, configFileNameLength),
55 *solverProcessIndex,
56 *solverProcessSize);
57} catch (::precice::Error &e) {
58 std::abort();
59}
60
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)
69try {
70#ifndef PRECICE_NO_MPI
71 // MPI_Comm_f2c requires an MPI_Fint argument, which may differ from int.
72 // However, the Fortran standard guarantees interoperability of the int type,
73 // so that is what is passed. Do the conversion in case int is not identical
74 // to MPI_Fint.
75 MPI_Fint f_communicator = static_cast<MPI_Fint>(*communicator);
76 auto c_communicator = MPI_Comm_f2c(f_communicator);
78 precice::impl::strippedStringView(participantName, participantNameLength),
79 precice::impl::strippedStringView(configFileName, configFileNameLength),
80 *solverProcessIndex,
81 *solverProcessSize,
82 &c_communicator);
83#else
84 PRECICE_WARN("preCICE was configured without MPI but you passed an MPI communicator. preCICE ignores the communicator and continues.");
86 precice::impl::strippedStringView(participantName, participantNameLength),
87 precice::impl::strippedStringView(configFileName, configFileNameLength),
88 *solverProcessIndex,
89 *solverProcessSize);
90#endif
91} catch (::precice::Error &e) {
92 std::abort();
93}
94
96try {
97 PRECICE_CHECK(impl != nullptr, errormsg);
98 impl->initialize();
99} catch (::precice::Error &e) {
100 std::abort();
101}
102
104 const double *timeStepSize)
105try {
106 PRECICE_CHECK(impl != nullptr, errormsg);
107 impl->advance(*timeStepSize);
108} catch (::precice::Error &e) {
109 std::abort();
110}
111
113try {
114 PRECICE_CHECK(impl != nullptr, errormsg);
115 impl->finalize();
116 impl.reset();
117} catch (::precice::Error &e) {
118 std::abort();
119}
120
122 const char *meshName,
123 int *dimensions,
124 int meshNameLength)
125try {
126 PRECICE_CHECK(impl != nullptr, errormsg);
127 *dimensions = impl->getMeshDimensions(precice::impl::strippedStringView(meshName, meshNameLength));
128} catch (::precice::Error &e) {
129 std::abort();
130}
131
133 const char *meshName,
134 const char *dataName,
135 int *dimensions,
136 int meshNameLength,
137 int dataNameLength)
138try {
139 PRECICE_CHECK(impl != nullptr, errormsg);
140 *dimensions = impl->getDataDimensions(precice::impl::strippedStringView(meshName, meshNameLength), precice::impl::strippedStringView(dataName, dataNameLength));
141} catch (::precice::Error &e) {
142 std::abort();
143}
144
146 int *isOngoing)
147try {
148 PRECICE_CHECK(impl != nullptr, errormsg);
149 if (impl->isCouplingOngoing()) {
150 *isOngoing = 1;
151 } else {
152 *isOngoing = 0;
153 }
154} catch (::precice::Error &e) {
155 std::abort();
156}
157
159 int *isComplete)
160try {
161 PRECICE_CHECK(impl != nullptr, errormsg);
162 if (impl->isTimeWindowComplete()) {
163 *isComplete = 1;
164 } else {
165 *isComplete = 0;
166 }
167} catch (::precice::Error &e) {
168 std::abort();
169}
170
172 double *maxTimeStepSize)
173try {
174 PRECICE_CHECK(impl != nullptr, errormsg);
175 *maxTimeStepSize = impl->getMaxTimeStepSize();
176} catch (::precice::Error &e) {
177 std::abort();
178}
179
181 int *isRequired)
182try {
183 PRECICE_CHECK(impl != nullptr, errormsg);
184 *isRequired = impl->requiresInitialData() ? 1 : 0;
185} catch (::precice::Error &e) {
186 std::abort();
187}
188
190 int *isRequired)
191try {
192 PRECICE_CHECK(impl != nullptr, errormsg);
193 *isRequired = impl->requiresWritingCheckpoint() ? 1 : 0;
194} catch (::precice::Error &e) {
195 std::abort();
196}
197
199 int *isRequired)
200try {
201 PRECICE_CHECK(impl != nullptr, errormsg);
202 *isRequired = impl->requiresReadingCheckpoint() ? 1 : 0;
203} catch (::precice::Error &e) {
204 std::abort();
205}
206
208 const char *meshName,
209 int *required,
210 int meshNameLength)
211try {
212 PRECICE_CHECK(impl != nullptr, errormsg);
213 if (impl->requiresMeshConnectivityFor(precice::impl::strippedStringView(meshName, meshNameLength))) {
214 *required = 1;
215 } else {
216 *required = 0;
217 }
218} catch (::precice::Error &e) {
219 std::abort();
220}
221
223 const char *meshName,
224 int meshNameLength)
225try {
226 PRECICE_CHECK(impl != nullptr, errormsg);
227 impl->resetMesh(precice::impl::strippedStringView(meshName, meshNameLength));
228} catch (::precice::Error &e) {
229 std::abort();
230}
231
233 const char *meshName,
234 const double *position,
235 int *vertexID,
236 int meshNameLength)
237try {
238 PRECICE_CHECK(impl != nullptr, errormsg);
239 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
240 auto positionSize = static_cast<unsigned long>(impl->getMeshDimensions(sv));
241 *vertexID = impl->setMeshVertex(sv, {position, positionSize});
242} catch (::precice::Error &e) {
243 std::abort();
244}
245
247 const char *meshName,
248 int *meshSize,
249 int meshNameLength)
250try {
251 PRECICE_CHECK(impl != nullptr, errormsg);
252 *meshSize = impl->getMeshVertexSize(precice::impl::strippedStringView(meshName, meshNameLength));
253} catch (::precice::Error &e) {
254 std::abort();
255}
256
258 const char *meshName,
259 const int *size,
260 double *coordinates,
261 int *ids,
262 int meshNameLength)
263try {
264 PRECICE_CHECK(impl != nullptr, errormsg);
265 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
266 auto positionSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * *size);
267 impl->setMeshVertices(sv, {coordinates, positionSize}, {ids, static_cast<unsigned long>(*size)});
268} catch (::precice::Error &e) {
269 std::abort();
270}
271
273 const char *meshName,
274 const int *firstVertexID,
275 const int *secondVertexID,
276 int meshNameLength)
277try {
278 PRECICE_CHECK(impl != nullptr, errormsg);
279 impl->setMeshEdge(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID);
280} catch (::precice::Error &e) {
281 std::abort();
282}
283
285 const char *meshName,
286 const int *size,
287 const int *ids,
288 int meshNameLength)
289try {
290 PRECICE_CHECK(impl != nullptr, errormsg);
291 auto idsSize = static_cast<unsigned long>(*size) * 2;
292 impl->setMeshEdges(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
293} catch (::precice::Error &e) {
294 std::abort();
295}
296
298 const char *meshName,
299 const int *firstVertexID,
300 const int *secondVertexID,
301 const int *thirdVertexID,
302 int meshNameLength)
303try {
304 PRECICE_CHECK(impl != nullptr, errormsg);
305 impl->setMeshTriangle(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID, *thirdVertexID);
306} catch (::precice::Error &e) {
307 std::abort();
308}
309
311 const char *meshName,
312 const int *size,
313 const int *ids,
314 int meshNameLength)
315try {
316 PRECICE_CHECK(impl != nullptr, errormsg);
317 auto idsSize = static_cast<unsigned long>(*size) * 3;
318 impl->setMeshTriangles(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
319} catch (::precice::Error &e) {
320 std::abort();
321}
322
324 const char *meshName,
325 const int *firstVertexID,
326 const int *secondVertexID,
327 const int *thirdVertexID,
328 const int *fourthVertexID,
329 int meshNameLength)
330try {
331 PRECICE_CHECK(impl != nullptr, errormsg);
332 impl->setMeshQuad(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID, *thirdVertexID, *fourthVertexID);
333} catch (::precice::Error &e) {
334 std::abort();
335}
336
338 const char *meshName,
339 const int *size,
340 const int *ids,
341 int meshNameLength)
342try {
343 PRECICE_CHECK(impl != nullptr, errormsg);
344 auto idsSize = static_cast<unsigned long>(*size) * 4;
345 impl->setMeshQuads(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
346} catch (::precice::Error &e) {
347 std::abort();
348}
349
351 const char *meshName,
352 const int *firstVertexID,
353 const int *secondVertexID,
354 const int *thirdVertexID,
355 const int *fourthVertexID,
356 int meshNameLength)
357try {
358 PRECICE_CHECK(impl != nullptr, errormsg);
359 impl->setMeshTetrahedron(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID, *thirdVertexID, *fourthVertexID);
360} catch (::precice::Error &e) {
361 std::abort();
362}
363
365 const char *meshName,
366 const int *size,
367 const int *ids,
368 int meshNameLength)
369try {
370 PRECICE_CHECK(impl != nullptr, errormsg);
371 auto idsSize = static_cast<unsigned long>(*size) * 4;
372 impl->setMeshTetrahedra(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
373} catch (::precice::Error &e) {
374 std::abort();
375}
376
378 const char *meshName,
379 const char *dataName,
380 const int *size,
381 int *ids,
382 double *values,
383 int meshNameLength,
384 int dataNameLength)
385try {
386 PRECICE_CHECK(impl != nullptr, errormsg);
387 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
388 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
389 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
390 impl->writeData(strippedMeshName,
391 strippedDataName,
392 {ids, static_cast<unsigned long>(*size)},
393 {values, static_cast<unsigned long>(dataSize)});
394} catch (::precice::Error &e) {
395 std::abort();
396}
397
399 const char *meshName,
400 const char *dataName,
401 const int *size,
402 int *ids,
403 const double *relativeReadTime,
404 double *values,
405 int meshNameLength,
406 int dataNameLength)
407try {
408 PRECICE_CHECK(impl != nullptr, errormsg);
409 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
410 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
411 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
412 impl->readData(
413 strippedMeshName,
414 strippedDataName,
415 {ids, static_cast<unsigned long>(*size)},
416 *relativeReadTime,
417 {values, static_cast<unsigned long>(dataSize)});
418} catch (::precice::Error &e) {
419 std::abort();
420}
421
423 const char *string,
424 int length)
425{
426 int i = length - 1;
427 while (((string[i] == ' ') || (string[i] == 0)) && (i >= 0)) {
428 i--;
429 }
430 return i + 1;
431}
432
434{
435 return {string, static_cast<std::string_view::size_type>(strippedLength(string, length))};
436}
437
439 char *versionInfo,
440 int lengthVersionInfo)
441{
442 std::strncpy(versionInfo, precice::versionInformation, lengthVersionInfo);
443}
444
446 const char *meshName,
447 const char *dataName, int *required,
448 int meshNameLength,
449 int dataNameLength)
450try {
451 PRECICE_CHECK(impl != nullptr, errormsg);
452 if (impl->requiresGradientDataFor(precice::impl::strippedStringView(meshName, meshNameLength), precice::impl::strippedStringView(dataName, dataNameLength))) {
453 *required = 1;
454 } else {
455 *required = 0;
456 }
457} catch (::precice::Error &e) {
458 std::abort();
459}
460
462 const char *meshName,
463 const char *dataName,
464 const int *size,
465 const int *ids,
466 const double *gradients,
467 int meshNameLength,
468 int dataNameLength)
469try {
470 PRECICE_CHECK(impl != nullptr, errormsg);
471 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
472 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
473 auto gradientComponents = impl->getMeshDimensions(strippedMeshName) * impl->getDataDimensions(strippedMeshName, strippedDataName);
474 auto gradientSize = *size * gradientComponents;
475
476 impl->writeGradientData(strippedMeshName,
477 strippedDataName,
478 {ids, static_cast<unsigned long>(*size)},
479 {gradients, static_cast<unsigned long>(gradientSize)});
480} catch (::precice::Error &e) {
481 std::abort();
482}
483
485 const char *meshName,
486 const char *dataName,
487 const int *size,
488 double *coordinates,
489 double *values,
490 int meshNameLength,
491 int dataNameLength)
492try {
493 PRECICE_CHECK(impl != nullptr, errormsg);
494 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
495 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
496 auto coordinatesSize = *size * impl->getMeshDimensions(strippedMeshName);
497 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
498 impl->writeAndMapData(strippedMeshName,
499 strippedDataName,
500 {coordinates, static_cast<unsigned long>(coordinatesSize)},
501 {values, static_cast<unsigned long>(dataSize)});
502} catch (::precice::Error &e) {
503 std::abort();
504}
505
507 const char *meshName,
508 const char *dataName,
509 const int *size,
510 double *coordinates,
511 const double *relativeReadTime,
512 double *values,
513 int meshNameLength,
514 int dataNameLength)
515try {
516 PRECICE_CHECK(impl != nullptr, errormsg);
517 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
518 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
519 auto coordinatesSize = *size * impl->getMeshDimensions(strippedMeshName);
520 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
521 impl->mapAndReadData(
522 strippedMeshName,
523 strippedDataName,
524 {coordinates, static_cast<unsigned long>(coordinatesSize)},
525 *relativeReadTime,
526 {values, static_cast<unsigned long>(dataSize)});
527} catch (::precice::Error &e) {
528 std::abort();
529}
530
532 const char *meshName,
533 const double *boundingBox,
534 int meshNameLength)
535try {
536 PRECICE_CHECK(impl != nullptr, errormsg);
537 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
538 auto bbSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * 2);
539 impl->setMeshAccessRegion(sv, {boundingBox, bbSize});
540} catch (::precice::Error &e) {
541 std::abort();
542}
543
545 const char *meshName,
546 const int *size,
547 int *ids,
548 double *coordinates,
549 int meshNameLength)
550try {
551 PRECICE_CHECK(impl != nullptr, errormsg);
552 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
553 auto coordinatesSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * *size);
554 impl->getMeshVertexIDsAndCoordinates(sv, {ids, static_cast<unsigned long>(*size)}, {coordinates, coordinatesSize});
555} catch (::precice::Error &e) {
556 std::abort();
557}
558
560 const char *sectionName,
561 int sectionNameLength)
562try {
563 PRECICE_CHECK(impl != nullptr, errormsg);
564 auto sv = precice::impl::strippedStringView(sectionName, sectionNameLength);
565 impl->startProfilingSection(sv);
566} catch (::precice::Error &e) {
567 std::abort();
568}
569
571try {
572 PRECICE_CHECK(impl != nullptr, errormsg);
573 impl->stopLastProfilingSection();
574} catch (::precice::Error &e) {
575 std::abort();
576}
577
578#ifdef __GNUC__
579#pragma GCC diagnostic pop
580#endif
581#ifdef __clang__
582#pragma clang diagnostic pop
583#endif
#define PRECICE_WARN(...)
Definition LogMacros.hpp:12
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
T abort(T... args)
This class provides a lightweight logger.
Definition Logger.hpp:17
T make_unique(T... args)
int strippedLength(const char *string, int length)
Returns length of string without trailing whitespace.
std::string_view strippedStringView(const char *string, int length)
STL namespace.
static std::unique_ptr< precice::Participant > impl
Definition preciceC.cpp:21
static std::string errormsg
Definition preciceC.cpp:25
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)
T size(T... args)
T strncpy(T... args)