preCICE v3.1.1
Loading...
Searching...
No Matches
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
42 const char *participantName,
43 const char *configFileName,
44 const int * solverProcessIndex,
45 const int * solverProcessSize,
46 int participantNameLength,
47 int configFileNameLength)
48{
49 impl.reset(new precice::Participant(
50 precice::impl::strippedStringView(participantName, participantNameLength),
51 precice::impl::strippedStringView(configFileName, configFileNameLength),
52 *solverProcessIndex,
53 *solverProcessSize));
54}
55
57{
58 PRECICE_CHECK(impl != nullptr, errormsg);
59 impl->initialize();
60}
61
63 const double *timeStepSize)
64{
65 PRECICE_CHECK(impl != nullptr, errormsg);
66 impl->advance(*timeStepSize);
67}
68
70{
71 PRECICE_CHECK(impl != nullptr, errormsg);
72 impl->finalize();
73 impl.reset();
74}
75
77 const char *meshName,
78 int * dimensions,
79 int meshNameLength)
80{
81 PRECICE_CHECK(impl != nullptr, errormsg);
82 *dimensions = impl->getMeshDimensions(precice::impl::strippedStringView(meshName, meshNameLength));
83}
84
86 const char *meshName,
87 const char *dataName,
88 int * dimensions,
89 int meshNameLength,
90 int dataNameLength)
91{
92 PRECICE_CHECK(impl != nullptr, errormsg);
93 *dimensions = impl->getDataDimensions(precice::impl::strippedStringView(meshName, meshNameLength), precice::impl::strippedStringView(dataName, dataNameLength));
94}
95
97 int *isOngoing)
98{
99 PRECICE_CHECK(impl != nullptr, errormsg);
100 if (impl->isCouplingOngoing()) {
101 *isOngoing = 1;
102 } else {
103 *isOngoing = 0;
104 }
105}
106
108 int *isComplete)
109{
110 PRECICE_CHECK(impl != nullptr, errormsg);
111 if (impl->isTimeWindowComplete()) {
112 *isComplete = 1;
113 } else {
114 *isComplete = 0;
115 }
116}
117
119 double *maxTimeStepSize)
120{
121 PRECICE_CHECK(impl != nullptr, errormsg);
122 *maxTimeStepSize = impl->getMaxTimeStepSize();
123}
124
126 int *isRequired)
127{
128 PRECICE_CHECK(impl != nullptr, errormsg);
129 *isRequired = impl->requiresInitialData() ? 1 : 0;
130}
131
133 int *isRequired)
134{
135 PRECICE_CHECK(impl != nullptr, errormsg);
136 *isRequired = impl->requiresWritingCheckpoint() ? 1 : 0;
137}
138
140 int *isRequired)
141{
142 PRECICE_CHECK(impl != nullptr, errormsg);
143 *isRequired = impl->requiresReadingCheckpoint() ? 1 : 0;
144}
145
147 const char *meshName,
148 int * required,
149 int meshNameLength)
150{
151 PRECICE_CHECK(impl != nullptr, errormsg);
152 if (impl->requiresMeshConnectivityFor(precice::impl::strippedStringView(meshName, meshNameLength))) {
153 *required = 1;
154 } else {
155 *required = 0;
156 }
157}
158
160 const char * meshName,
161 const double *position,
162 int * vertexID,
163 int meshNameLength)
164{
165 PRECICE_CHECK(impl != nullptr, errormsg);
166 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
167 auto positionSize = static_cast<unsigned long>(impl->getMeshDimensions(sv));
168 *vertexID = impl->setMeshVertex(sv, {position, positionSize});
169}
170
172 const char *meshName,
173 int * meshSize,
174 int meshNameLength)
175{
176 PRECICE_CHECK(impl != nullptr, errormsg);
177 *meshSize = impl->getMeshVertexSize(precice::impl::strippedStringView(meshName, meshNameLength));
178}
179
181 const char *meshName,
182 const int * size,
183 double * coordinates,
184 int * ids,
185 int meshNameLength)
186{
187 PRECICE_CHECK(impl != nullptr, errormsg);
188 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
189 auto positionSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * *size);
190 impl->setMeshVertices(sv, {coordinates, positionSize}, {ids, static_cast<unsigned long>(*size)});
191}
192
194 const char *meshName,
195 const int * firstVertexID,
196 const int * secondVertexID,
197 int meshNameLength)
198{
199 PRECICE_CHECK(impl != nullptr, errormsg);
200 impl->setMeshEdge(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID);
201}
202
204 const char *meshName,
205 const int * size,
206 const int * ids,
207 int meshNameLength)
208{
209 PRECICE_CHECK(impl != nullptr, errormsg);
210 auto idsSize = static_cast<unsigned long>(*size) * 2;
211 impl->setMeshEdges(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
212}
213
215 const char *meshName,
216 const int * firstVertexID,
217 const int * secondVertexID,
218 const int * thirdVertexID,
219 int meshNameLength)
220{
221 PRECICE_CHECK(impl != nullptr, errormsg);
222 impl->setMeshTriangle(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID, *thirdVertexID);
223}
224
226 const char *meshName,
227 const int * size,
228 const int * ids,
229 int meshNameLength)
230{
231 PRECICE_CHECK(impl != nullptr, errormsg);
232 auto idsSize = static_cast<unsigned long>(*size) * 3;
233 impl->setMeshTriangles(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
234}
235
237 const char *meshName,
238 const int * firstVertexID,
239 const int * secondVertexID,
240 const int * thirdVertexID,
241 const int * fourthVertexID,
242 int meshNameLength)
243{
244 PRECICE_CHECK(impl != nullptr, errormsg);
245 impl->setMeshQuad(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID, *thirdVertexID, *fourthVertexID);
246}
247
249 const char *meshName,
250 const int * size,
251 const int * ids,
252 int meshNameLength)
253{
254 PRECICE_CHECK(impl != nullptr, errormsg);
255 auto idsSize = static_cast<unsigned long>(*size) * 4;
256 impl->setMeshQuads(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
257}
258
260 const char *meshName,
261 const int * firstVertexID,
262 const int * secondVertexID,
263 const int * thirdVertexID,
264 const int * fourthVertexID,
265 int meshNameLength)
266{
267 PRECICE_CHECK(impl != nullptr, errormsg);
268 impl->setMeshTetrahedron(precice::impl::strippedStringView(meshName, meshNameLength), *firstVertexID, *secondVertexID, *thirdVertexID, *fourthVertexID);
269}
270
272 const char *meshName,
273 const int * size,
274 const int * ids,
275 int meshNameLength)
276{
277 PRECICE_CHECK(impl != nullptr, errormsg);
278 auto idsSize = static_cast<unsigned long>(*size) * 4;
279 impl->setMeshTetrahedra(precice::impl::strippedStringView(meshName, meshNameLength), {ids, idsSize});
280}
281
283 const char *meshName,
284 const char *dataName,
285 const int * size,
286 int * ids,
287 double * values,
288 int meshNameLength,
289 int dataNameLength)
290{
291 PRECICE_CHECK(impl != nullptr, errormsg);
292 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
293 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
294 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
295 impl->writeData(strippedMeshName,
296 strippedDataName,
297 {ids, static_cast<unsigned long>(*size)},
298 {values, static_cast<unsigned long>(dataSize)});
299}
300
302 const char * meshName,
303 const char * dataName,
304 const int * size,
305 int * ids,
306 const double *relativeReadTime,
307 double * values,
308 int meshNameLength,
309 int dataNameLength)
310{
311 PRECICE_CHECK(impl != nullptr, errormsg);
312 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
313 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
314 auto dataSize = *size * impl->getDataDimensions(strippedMeshName, strippedDataName);
315 impl->readData(
316 strippedMeshName,
317 strippedDataName,
318 {ids, static_cast<unsigned long>(*size)},
319 *relativeReadTime,
320 {values, static_cast<unsigned long>(dataSize)});
321}
322
324 const char *string,
325 int length)
326{
327 int i = length - 1;
328 while (((string[i] == ' ') || (string[i] == 0)) && (i >= 0)) {
329 i--;
330 }
331 return i + 1;
332}
333
335{
336 return {string, static_cast<std::string_view::size_type>(strippedLength(string, length))};
337}
338
340 char *versionInfo,
341 int lengthVersionInfo)
342{
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];
347 }
348}
349
351 const char *meshName,
352 const char *dataName, int *required,
353 int meshNameLength,
354 int dataNameLength)
355{
356 PRECICE_CHECK(impl != nullptr, errormsg);
357 if (impl->requiresGradientDataFor(precice::impl::strippedStringView(meshName, meshNameLength), precice::impl::strippedStringView(dataName, dataNameLength))) {
358 *required = 1;
359 } else {
360 *required = 0;
361 }
362}
363
365 const char * meshName,
366 const char * dataName,
367 const int * size,
368 const int * ids,
369 const double *gradients,
370 int meshNameLength,
371 int dataNameLength)
372{
373 PRECICE_CHECK(impl != nullptr, errormsg);
374 auto strippedMeshName = precice::impl::strippedStringView(meshName, meshNameLength);
375 auto strippedDataName = precice::impl::strippedStringView(dataName, dataNameLength);
376 auto gradientComponents = impl->getMeshDimensions(strippedMeshName) * impl->getDataDimensions(strippedMeshName, strippedDataName);
377 auto gradientSize = *size * gradientComponents;
378
379 impl->writeGradientData(strippedMeshName,
380 strippedDataName,
381 {ids, static_cast<unsigned long>(*size)},
382 {gradients, static_cast<unsigned long>(gradientSize)});
383}
384
386 const char * meshName,
387 const double *boundingBox,
388 int meshNameLength)
389{
390 PRECICE_CHECK(impl != nullptr, errormsg);
391 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
392 auto bbSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * 2);
393 impl->setMeshAccessRegion(sv, {boundingBox, bbSize});
394}
395
397 const char *meshName,
398 const int size,
399 int * ids,
400 double * coordinates,
401 int meshNameLength)
402{
403 PRECICE_CHECK(impl != nullptr, errormsg);
404 auto sv = precice::impl::strippedStringView(meshName, meshNameLength);
405 auto coordinatesSize = static_cast<unsigned long>(impl->getMeshDimensions(sv) * size);
406 impl->getMeshVertexIDsAndCoordinates(sv, {ids, static_cast<unsigned long>(size)}, {coordinates, coordinatesSize});
407}
408
409#ifdef __GNUC__
410#pragma GCC diagnostic pop
411#endif
412#ifdef __clang__
413#pragma clang diagnostic pop
414#endif
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:35
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
Main Application Programming Interface of preCICE. Include using #include <precice/precice....
This class provides a lightweight logger.
Definition Logger.hpp:16
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)
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)
T size(T... args)