preCICE
Loading...
Searching...
No Matches
preciceFortran.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <precice/export.h>
4
13
14#ifdef __cplusplus
15extern "C" {
16#endif
17
20
35PRECICE_API void precicef_create_(
36 const char *participantName,
37 const char *configFileName,
38 const int *solverProcessIndex,
39 const int *solverProcessSize,
40 int participantNameLength,
41 int configFileNameLength);
42
60 const char *participantName,
61 const char *configFileName,
62 const int *solverProcessIndex,
63 const int *solverProcessSize,
64 const int *communicator,
65 int participantNameLength,
66 int configFileNameLength);
67
69
72
83PRECICE_API void precicef_initialize_();
84
95PRECICE_API void precicef_advance_(const double *timeStepSize);
96
104PRECICE_API void precicef_finalize_();
105
107
110
122 int *isRequired);
123
135 int *isRequired);
136
138
141
154PRECICE_API void precicef_get_mesh_dimensions_(
155 const char *meshName,
156 int *dimensions,
157 int meshNameLength);
158
172PRECICE_API void precicef_get_data_dimensions_(
173 const char *meshName,
174 const char *dataName,
175 int *dimensions,
176 int meshNameLength,
177 int dataNameLength);
178
189PRECICE_API void precicef_is_coupling_ongoing_(int *isOngoing);
190
201PRECICE_API void precicef_is_time_window_complete_(int *isComplete);
202
213PRECICE_API void precicef_get_max_time_step_size_(double *maxTimeStepSize);
214
216
220
233 const char *meshName,
234 int *required,
235 int meshNameLength);
236
246PRECICE_API void precicef_reset_mesh_(
247 const char *meshName,
248 int meshNameLength);
249
263PRECICE_API void precicef_set_mesh_vertex_(
264 const char *meshName,
265 const double *coordinates,
266 int *id,
267 int meshNameLength);
268
273PRECICE_API void precicef_set_vertex_(
274 const char *meshName,
275 const double *coordinates,
276 int *id,
277 int meshNameLength);
278
291PRECICE_API void precicef_get_mesh_vertex_size_(
292 const char *meshName,
293 int *meshSize,
294 int meshNameLength);
295
310PRECICE_API void precicef_set_mesh_vertices_(
311 const char *meshName,
312 const int *size,
313 double *coordinates,
314 int *ids,
315 int meshNameLength);
316
321PRECICE_API void precicef_set_vertices_(
322 const char *meshName,
323 const int *size,
324 double *coordinates,
325 int *ids,
326 int meshNameLength);
327
341PRECICE_API void precicef_set_mesh_edge_(
342 const char *meshName,
343 const int *firstVertexID,
344 const int *secondVertexID,
345 int meshNameLength);
346
351PRECICE_API void precicef_set_edge_(
352 const char *meshName,
353 const int *firstVertexID,
354 const int *secondVertexID,
355 int meshNameLength);
356
370PRECICE_API void precicef_set_mesh_edges_(
371 const char *meshName,
372 const int *size,
373 const int *ids,
374 int meshNameLength);
375
390PRECICE_API void precicef_set_mesh_triangle_(
391 const char *meshName,
392 const int *firstVertexID,
393 const int *secondVertexID,
394 const int *thirdVertexID,
395 int meshNameLength);
396
401PRECICE_API void precicef_set_triangle_(
402 const char *meshName,
403 const int *firstVertexID,
404 const int *secondVertexID,
405 const int *thirdVertexID,
406 int meshNameLength);
407
421PRECICE_API void precicef_set_mesh_triangles_(
422 const char *meshName,
423 const int *size,
424 const int *ids,
425 int meshNameLength);
426
442PRECICE_API void precicef_set_mesh_quad_(
443 const char *meshName,
444 const int *firstVertexID,
445 const int *secondVertexID,
446 const int *thirdVertexID,
447 const int *fourthVertexID,
448 int meshNameLength);
449
454PRECICE_API void precicef_set_quad_(
455 const char *meshName,
456 const int *firstVertexID,
457 const int *secondVertexID,
458 const int *thirdVertexID,
459 const int *fourthVertexID,
460 int meshNameLength);
461
475PRECICE_API void precicef_set_mesh_quads_(
476 const char *meshName,
477 const int *size,
478 const int *ids,
479 int meshNameLength);
480
496PRECICE_API void precicef_set_mesh_tetrahedron_(
497 const char *meshName,
498 const int *firstVertexID,
499 const int *secondVertexID,
500 const int *thirdVertexID,
501 const int *fourthVertexID,
502 int meshNameLength);
503
508PRECICE_API void precicef_set_tetrahedron(
509 const char *meshName,
510 const int *firstVertexID,
511 const int *secondVertexID,
512 const int *thirdVertexID,
513 const int *fourthVertexID,
514 int meshNameLength);
515
529PRECICE_API void precicef_set_mesh_tetrahedra_(
530 const char *meshName,
531 const int *size,
532 const int *ids,
533 int meshNameLength);
534
536
539
550PRECICE_API void precicef_requires_initial_data_(
551 int *isRequired);
552
568PRECICE_API void precicef_write_data_(
569 const char *meshName,
570 const char *dataName,
571 const int *size,
572 int *ids,
573 double *values,
574 int meshNameLength,
575 int dataNameLength);
576
593PRECICE_API void precicef_read_data_(
594 const char *meshName,
595 const char *dataName,
596 const int *size,
597 int *ids,
598 const double *relativeReadTime,
599 double *values,
600 int meshNameLength,
601 int dataNameLength);
602
604
607
627PRECICE_API void precicef_write_and_map_data_(
628 const char *meshName,
629 const char *dataName,
630 const int *size,
631 double *coordinates,
632 double *values,
633 int meshNameLength,
634 int dataNameLength);
635
657PRECICE_API void precicef_map_and_read_data_(
658 const char *meshName,
659 const char *dataName,
660 const int *size,
661 double *coordinates,
662 const double *relativeReadTime,
663 double *values,
664 int meshNameLength,
665 int dataNameLength);
667
670
682PRECICE_API void precicef_set_mesh_access_region_(
683 const char *meshName,
684 const double *boundingBox,
685 int meshNameLength);
686
701 const char *meshName,
702 const int *size,
703 int *ids,
704 double *coordinates,
705 int meshNameLength);
706
708
713
727 const char *meshName,
728 const char *dataName, int *required,
729 int meshNameLength,
730 int dataNameLength);
731
746PRECICE_API void precicef_write_gradient_data_(
747 const char *meshName,
748 const char *dataName,
749 const int *size,
750 const int *ids,
751 const double *gradients,
752 int meshNameLength,
753 int dataNameLength);
754
756
759
761PRECICE_API void precicef_start_profiling_section_(const char *sectionName, int sectionNameLength);
762
765
767
768PRECICE_API void precicef_get_version_information_(
769 char *versionInfo,
770 int lengthVersionInfo);
771
772#ifdef __cplusplus
773}
774#endif
PRECICE_API void precicef_get_version_information_(char *versionInfo, int lengthVersionInfo)
PRECICE_API void precicef_requires_writing_checkpoint_(int *isRequired)
PRECICE_API void precicef_set_mesh_access_region_(const char *meshName, const double *boundingBox, int meshNameLength)
PRECICE_API void precicef_set_vertices_(const char *meshName, const int *size, double *coordinates, int *ids, int meshNameLength)
PRECICE_API void precicef_set_quad_(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
PRECICE_API void precicef_requires_initial_data_(int *isRequired)
PRECICE_API void precicef_reset_mesh_(const char *meshName, int meshNameLength)
PRECICE_API void precicef_set_edge_(const char *meshName, const int *firstVertexID, const int *secondVertexID, int meshNameLength)
PRECICE_API void precicef_get_mesh_vertex_ids_and_coordinates_(const char *meshName, const int *size, int *ids, double *coordinates, int meshNameLength)
PRECICE_API void precicef_create_(const char *participantName, const char *configFileName, const int *solverProcessIndex, const int *solverProcessSize, int participantNameLength, int configFileNameLength)
PRECICE_API void precicef_set_vertex_(const char *meshName, const double *coordinates, int *id, int meshNameLength)
PRECICE_API void precicef_write_gradient_data_(const char *meshName, const char *dataName, const int *size, const int *ids, const double *gradients, int meshNameLength, int dataNameLength)
PRECICE_API void precicef_set_mesh_edges_(const char *meshName, const int *size, const int *ids, int meshNameLength)
PRECICE_API void precicef_set_mesh_triangles_(const char *meshName, const int *size, const int *ids, int meshNameLength)
PRECICE_API void precicef_requires_reading_checkpoint_(int *isRequired)
PRECICE_API void precicef_set_mesh_quads_(const char *meshName, const int *size, const int *ids, int meshNameLength)
PRECICE_API void precicef_set_tetrahedron(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
PRECICE_API void precicef_initialize_()
PRECICE_API 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.
PRECICE_API void precicef_advance_(const double *timeStepSize)
PRECICE_API void precicef_is_coupling_ongoing_(int *isOngoing)
PRECICE_API void precicef_set_mesh_triangle_(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, int meshNameLength)
PRECICE_API 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.
PRECICE_API void precicef_get_max_time_step_size_(double *maxTimeStepSize)
PRECICE_API void precicef_set_mesh_vertices_(const char *meshName, const int *size, double *coordinates, int *ids, int meshNameLength)
PRECICE_API void precicef_requires_mesh_connectivity_for_(const char *meshName, int *required, int meshNameLength)
PRECICE_API void precicef_is_time_window_complete_(int *isComplete)
PRECICE_API void precicef_get_mesh_dimensions_(const char *meshName, int *dimensions, int meshNameLength)
PRECICE_API void precicef_get_data_dimensions_(const char *meshName, const char *dataName, int *dimensions, int meshNameLength, int dataNameLength)
PRECICE_API void precicef_stop_last_profiling_section_()
PRECICE_API void precicef_write_data_(const char *meshName, const char *dataName, const int *size, int *ids, double *values, int meshNameLength, int dataNameLength)
PRECICE_API void precicef_create_with_communicator_(const char *participantName, const char *configFileName, const int *solverProcessIndex, const int *solverProcessSize, const int *communicator, int participantNameLength, int configFileNameLength)
PRECICE_API void precicef_set_mesh_tetrahedron_(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
PRECICE_API void precicef_requires_gradient_data_for_(const char *meshName, const char *dataName, int *required, int meshNameLength, int dataNameLength)
PRECICE_API void precicef_get_mesh_vertex_size_(const char *meshName, int *meshSize, int meshNameLength)
PRECICE_API void precicef_finalize_()
PRECICE_API void precicef_set_mesh_tetrahedra_(const char *meshName, const int *size, const int *ids, int meshNameLength)
PRECICE_API void precicef_set_triangle_(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, int meshNameLength)
PRECICE_API void precicef_read_data_(const char *meshName, const char *dataName, const int *size, int *ids, const double *relativeReadTime, double *values, int meshNameLength, int dataNameLength)
PRECICE_API void precicef_set_mesh_edge_(const char *meshName, const int *firstVertexID, const int *secondVertexID, int meshNameLength)
PRECICE_API void precicef_set_mesh_quad_(const char *meshName, const int *firstVertexID, const int *secondVertexID, const int *thirdVertexID, const int *fourthVertexID, int meshNameLength)
PRECICE_API void precicef_set_mesh_vertex_(const char *meshName, const double *coordinates, int *id, int meshNameLength)
PRECICE_API void precicef_start_profiling_section_(const char *sectionName, int sectionNameLength)