preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
preciceFortran.cpp
Go to the documentation of this file.
2#include <cstddef>
3#include <iostream>
4#include <memory>
5#include <stdexcept>
6#include <string>
7#include <string_view>
9#include "logging/Logger.hpp"
10#include "precice/impl/versions.hpp"
11#include "precice/precice.hpp"
12#include "utils/assertion.hpp"
13
14#ifdef __GNUC__
15#pragma GCC diagnostic push
16#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
17#endif
18#ifdef __clang__
19#pragma clang diagnostic push
20#pragma clang diagnostic ignored "-Wdeprecated-declarations"
21#endif
22
23using namespace std;
24
26
27static precice::logging::Logger _log("preciceFortran");
28
29static std::string errormsg = "preCICE has not been created properly. Be sure to call \"precicef_create\" before any other call to preCICE.";
30
31namespace precice::impl {
35int strippedLength(const char *string, int length);
36
37std::string_view strippedStringView(const char *string, int length);
38
39} // namespace precice::impl
40
57
59try {
60 PRECICE_CHECK(impl != nullptr, errormsg);
61 impl->initialize();
62} catch (::precice::Error &e) {
63 std::abort();
64}
65
67 const double *timeStepSize)
68try {
69 PRECICE_CHECK(impl != nullptr, errormsg);
70 impl->advance(*timeStepSize);
71} catch (::precice::Error &e) {
72 std::abort();
73}
74
76try {
77 PRECICE_CHECK(impl != nullptr, errormsg);
78 impl->finalize();
79 impl.reset();
80} catch (::precice::Error &e) {
81 std::abort();
82}
83
85 const char *meshName,
86 int *dimensions,
88try {
89 PRECICE_CHECK(impl != nullptr, errormsg);
90 *dimensions = impl->getMeshDimensions(precice::impl::strippedStringView(meshName, meshNameLength));
91} catch (::precice::Error &e) {
92 std::abort();
93}
94
96 const char *meshName,
97 const char *dataName,
98 int *dimensions,
100 int dataNameLength)
101try {
102 PRECICE_CHECK(impl != nullptr, errormsg);
104} catch (::precice::Error &e) {
105 std::abort();
106}
107
109 int *isOngoing)
110try {
111 PRECICE_CHECK(impl != nullptr, errormsg);
112 if (impl->isCouplingOngoing()) {
113 *isOngoing = 1;
114 } else {
115 *isOngoing = 0;
116 }
117} catch (::precice::Error &e) {
118 std::abort();
119}
120
122 int *isComplete)
123try {
124 PRECICE_CHECK(impl != nullptr, errormsg);
125 if (impl->isTimeWindowComplete()) {
126 *isComplete = 1;
127 } else {
128 *isComplete = 0;
129 }
130} catch (::precice::Error &e) {
131 std::abort();
132}
133
135 double *maxTimeStepSize)
136try {
137 PRECICE_CHECK(impl != nullptr, errormsg);
138 *maxTimeStepSize = impl->getMaxTimeStepSize();
139} catch (::precice::Error &e) {
140 std::abort();
141}
142
144 int *isRequired)
145try {
146 PRECICE_CHECK(impl != nullptr, errormsg);
147 *isRequired = impl->requiresInitialData() ? 1 : 0;
148} catch (::precice::Error &e) {
149 std::abort();
150}
151
153 int *isRequired)
154try {
155 PRECICE_CHECK(impl != nullptr, errormsg);
156 *isRequired = impl->requiresWritingCheckpoint() ? 1 : 0;
157} catch (::precice::Error &e) {
158 std::abort();
159}
160
162 int *isRequired)
163try {
164 PRECICE_CHECK(impl != nullptr, errormsg);
165 *isRequired = impl->requiresReadingCheckpoint() ? 1 : 0;
166} catch (::precice::Error &e) {
167 std::abort();
168}
169
171 const char *meshName,
172 int *required,
173 int meshNameLength)
174try {
175 PRECICE_CHECK(impl != nullptr, errormsg);
176 if (impl->requiresMeshConnectivityFor(precice::impl::strippedStringView(meshName, meshNameLength))) {
177 *required = 1;
178 } else {
179 *required = 0;
180 }
181} catch (::precice::Error &e) {
182 std::abort();
183}
184
186 const char *meshName,
187 int meshNameLength)
188try {
189 PRECICE_CHECK(impl != nullptr, errormsg);
191} catch (::precice::Error &e) {
192 std::abort();
193}
194
196 const char *meshName,
197 const double *position,
198 int *vertexID,
199 int meshNameLength)
200try {
201 PRECICE_CHECK(impl != nullptr, errormsg);
203 auto positionSize = static_cast<unsigned long>(impl->getMeshDimensions(sv));
204 *vertexID = impl->setMeshVertex(sv, {position, positionSize});
205} catch (::precice::Error &e) {
206 std::abort();
207}
208
210 const char *meshName,
211 int *meshSize,
212 int meshNameLength)
213try {
214 PRECICE_CHECK(impl != nullptr, errormsg);
215 *meshSize = impl->getMeshVertexSize(precice::impl::strippedStringView(meshName, meshNameLength));
216} catch (::precice::Error &e) {
217 std::abort();
218}
219
221 const char *meshName,
222 const int *size,
223 double *coordinates,
224 int *ids,
225 int meshNameLength)
226try {
227 PRECICE_CHECK(impl != nullptr, errormsg);
229 auto positionSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * *size);
230 impl->setMeshVertices(sv, {coordinates, positionSize}, {ids, static_cast<unsigned long>(*size)});
231} catch (::precice::Error &e) {
232 std::abort();
233}
234
236 const char *meshName,
237 const int *firstVertexID,
238 const int *secondVertexID,
239 int meshNameLength)
240try {
241 PRECICE_CHECK(impl != nullptr, errormsg);
243} catch (::precice::Error &e) {
244 std::abort();
245}
246
248 const char *meshName,
249 const int *size,
250 const int *ids,
251 int meshNameLength)
252try {
253 PRECICE_CHECK(impl != nullptr, errormsg);
254 auto idsSize = static_cast<unsigned long>(*size) * 2;
255 impl->setMeshEdges(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
256} catch (::precice::Error &e) {
257 std::abort();
258}
259
261 const char *meshName,
262 const int *firstVertexID,
263 const int *secondVertexID,
264 const int *thirdVertexID,
265 int meshNameLength)
266try {
267 PRECICE_CHECK(impl != nullptr, errormsg);
269} catch (::precice::Error &e) {
270 std::abort();
271}
272
274 const char *meshName,
275 const int *size,
276 const int *ids,
277 int meshNameLength)
278try {
279 PRECICE_CHECK(impl != nullptr, errormsg);
280 auto idsSize = static_cast<unsigned long>(*size) * 3;
281 impl->setMeshTriangles(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
282} catch (::precice::Error &e) {
283 std::abort();
284}
285
287 const char *meshName,
288 const int *firstVertexID,
289 const int *secondVertexID,
290 const int *thirdVertexID,
291 const int *fourthVertexID,
292 int meshNameLength)
293try {
294 PRECICE_CHECK(impl != nullptr, errormsg);
296} catch (::precice::Error &e) {
297 std::abort();
298}
299
301 const char *meshName,
302 const int *size,
303 const int *ids,
304 int meshNameLength)
305try {
306 PRECICE_CHECK(impl != nullptr, errormsg);
307 auto idsSize = static_cast<unsigned long>(*size) * 4;
308 impl->setMeshQuads(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
309} catch (::precice::Error &e) {
310 std::abort();
311}
312
314 const char *meshName,
315 const int *firstVertexID,
316 const int *secondVertexID,
317 const int *thirdVertexID,
318 const int *fourthVertexID,
319 int meshNameLength)
320try {
321 PRECICE_CHECK(impl != nullptr, errormsg);
323} catch (::precice::Error &e) {
324 std::abort();
325}
326
328 const char *meshName,
329 const int *size,
330 const int *ids,
331 int meshNameLength)
332try {
333 PRECICE_CHECK(impl != nullptr, errormsg);
334 auto idsSize = static_cast<unsigned long>(*size) * 4;
335 impl->setMeshTetrahedra(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
336} catch (::precice::Error &e) {
337 std::abort();
338}
339
341 const char *meshName,
342 const char *dataName,
343 const int *size,
344 int *ids,
345 double *values,
346 int meshNameLength,
347 int dataNameLength)
348try {
349 PRECICE_CHECK(impl != nullptr, errormsg);
352 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
353 impl->writeData(strippedMeshName,
355 {ids, static_cast<unsigned long>(*size)},
356 {values, static_cast<unsigned long>(dataSize)});
357} catch (::precice::Error &e) {
358 std::abort();
359}
360
362 const char *meshName,
363 const char *dataName,
364 const int *size,
365 int *ids,
366 const double *relativeReadTime,
367 double *values,
368 int meshNameLength,
369 int dataNameLength)
370try {
371 PRECICE_CHECK(impl != nullptr, errormsg);
374 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
375 impl->readData(
378 {ids, static_cast<unsigned long>(*size)},
380 {values, static_cast<unsigned long>(dataSize)});
381} catch (::precice::Error &e) {
382 std::abort();
383}
384
386 const char *string,
387 int length)
388{
389 int i = length - 1;
390 while (((string[i] == ' ') || (string[i] == 0)) && (i >= 0)) {
391 i--;
392 }
393 return i + 1;
394}
395
397{
398 return {string, static_cast<std::string_view::size_type>(strippedLength(string, length))};
399}
400
402 char *versionInfo,
404{
405 std::strncpy(versionInfo, precice::versionInformation, lengthVersionInfo);
406}
407
409 const char *meshName,
410 const char *dataName, int *required,
411 int meshNameLength,
412 int dataNameLength)
413try {
414 PRECICE_CHECK(impl != nullptr, errormsg);
416 *required = 1;
417 } else {
418 *required = 0;
419 }
420} catch (::precice::Error &e) {
421 std::abort();
422}
423
425 const char *meshName,
426 const char *dataName,
427 const int *size,
428 const int *ids,
429 const double *gradients,
430 int meshNameLength,
431 int dataNameLength)
432try {
433 PRECICE_CHECK(impl != nullptr, errormsg);
436 auto gradientComponents = impl->getMeshDimensions(strippedMeshName) * impl->getDataDimensions(strippedMeshName, strippedDataName);
438
439 impl->writeGradientData(strippedMeshName,
441 {ids, static_cast<unsigned long>(*size)},
442 {gradients, static_cast<unsigned long>(gradientSize)});
443} catch (::precice::Error &e) {
444 std::abort();
445}
446
448 const char *meshName,
449 const char *dataName,
450 const int *size,
451 double *coordinates,
452 double *values,
453 int meshNameLength,
454 int dataNameLength)
455try {
456 PRECICE_CHECK(impl != nullptr, errormsg);
459 auto coordinatesSize = *size * impl->getMeshDimensions(strippedMeshName);
460 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
461 impl->writeAndMapData(strippedMeshName,
463 {coordinates, static_cast<unsigned long>(coordinatesSize)},
464 {values, static_cast<unsigned long>(dataSize)});
465} catch (::precice::Error &e) {
466 std::abort();
467}
468
470 const char *meshName,
471 const char *dataName,
472 const int *size,
473 double *coordinates,
474 const double *relativeReadTime,
475 double *values,
476 int meshNameLength,
477 int dataNameLength)
478try {
479 PRECICE_CHECK(impl != nullptr, errormsg);
482 auto coordinatesSize = *size * impl->getMeshDimensions(strippedMeshName);
483 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
484 impl->mapAndReadData(
487 {coordinates, static_cast<unsigned long>(coordinatesSize)},
489 {values, static_cast<unsigned long>(dataSize)});
490} catch (::precice::Error &e) {
491 std::abort();
492}
493
495 const char *meshName,
496 const double *boundingBox,
497 int meshNameLength)
498try {
499 PRECICE_CHECK(impl != nullptr, errormsg);
501 auto bbSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * 2);
502 impl->setMeshAccessRegion(sv, {boundingBox, bbSize});
503} catch (::precice::Error &e) {
504 std::abort();
505}
506
508 const char *meshName,
509 const int size,
510 int *ids,
511 double *coordinates,
512 int meshNameLength)
513try {
514 PRECICE_CHECK(impl != nullptr, errormsg);
516 auto coordinatesSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * size);
517 impl->getMeshVertexIDsAndCoordinates(sv, {ids, static_cast<unsigned long>(size)}, {coordinates, coordinatesSize});
518} catch (::precice::Error &e) {
519 std::abort();
520}
521
523 const char *sectionName,
525try {
526 PRECICE_CHECK(impl != nullptr, errormsg);
528 impl->startProfilingSection(sv);
529} catch (::precice::Error &e) {
530 std::abort();
531}
532
534try {
535 PRECICE_CHECK(impl != nullptr, errormsg);
536 impl->stopLastProfilingSection();
537} catch (::precice::Error &e) {
538 std::abort();
539}
540
541#ifdef __GNUC__
542#pragma GCC diagnostic pop
543#endif
544#ifdef __clang__
545#pragma clang diagnostic pop
546#endif
#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.
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)
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)
static std::unique_ptr< precice::Participant > impl
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_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_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_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)