preCICE v3.3.0
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_vertex_(
264 const char *meshName,
265 const double *coordinates,
266 int *id,
267 int meshNameLength);
268
281PRECICE_API void precicef_get_mesh_vertex_size_(
282 const char *meshName,
283 int *meshSize,
284 int meshNameLength);
285
300PRECICE_API void precicef_set_vertices_(
301 const char *meshName,
302 const int *size,
303 double *coordinates,
304 int *ids,
305 int meshNameLength);
306
320PRECICE_API void precicef_set_edge_(
321 const char *meshName,
322 const int *firstVertexID,
323 const int *secondVertexID,
324 int meshNameLength);
325
339PRECICE_API void precicef_set_mesh_edges_(
340 const char *meshName,
341 const int *size,
342 const int *ids,
343 int meshNameLength);
344
359PRECICE_API void precicef_set_triangle_(
360 const char *meshName,
361 const int *firstVertexID,
362 const int *secondVertexID,
363 const int *thirdVertexID,
364 int meshNameLength);
365
379PRECICE_API void precicef_set_mesh_triangles_(
380 const char *meshName,
381 const int *size,
382 const int *ids,
383 int meshNameLength);
384
400PRECICE_API void precicef_set_quad_(
401 const char *meshName,
402 const int *firstVertexID,
403 const int *secondVertexID,
404 const int *thirdVertexID,
405 const int *fourthVertexID,
406 int meshNameLength);
407
421PRECICE_API void precicef_set_mesh_quads_(
422 const char *meshName,
423 const int *size,
424 const int *ids,
425 int meshNameLength);
426
442PRECICE_API void precicef_set_tetrahedron(
443 const char *meshName,
444 const int *firstVertexID,
445 const int *secondVertexID,
446 const int *thirdVertexID,
447 const int *fourthVertexID,
448 int meshNameLength);
449
463PRECICE_API void precicef_set_mesh_tetrahedra_(
464 const char *meshName,
465 const int *size,
466 const int *ids,
467 int meshNameLength);
468
470
473
484PRECICE_API void precicef_requires_initial_data_(
485 int *isRequired);
486
502PRECICE_API void precicef_write_data_(
503 const char *meshName,
504 const char *dataName,
505 const int *size,
506 int *ids,
507 double *values,
508 int meshNameLength,
509 int dataNameLength);
510
527PRECICE_API void precicef_read_data_(
528 const char *meshName,
529 const char *dataName,
530 const int *size,
531 int *ids,
532 const double *relativeReadTime,
533 double *values,
534 int meshNameLength,
535 int dataNameLength);
536
538
541
561PRECICE_API void precicef_write_and_map_data_(
562 const char *meshName,
563 const char *dataName,
564 const int *size,
565 double *coordinates,
566 double *values,
567 int meshNameLength,
568 int dataNameLength);
569
591PRECICE_API void precicef_map_and_read_data_(
592 const char *meshName,
593 const char *dataName,
594 const int *size,
595 double *coordinates,
596 const double *relativeReadTime,
597 double *values,
598 int meshNameLength,
599 int dataNameLength);
601
604
616PRECICE_API void precicef_set_mesh_access_region_(
617 const char *meshName,
618 const double *boundingBox,
619 int meshNameLength);
620
635 const char *meshName,
636 const int *size,
637 int *ids,
638 double *coordinates,
639 int meshNameLength);
640
642
647
661 const char *meshName,
662 const char *dataName, int *required,
663 int meshNameLength,
664 int dataNameLength);
665
680PRECICE_API void precicef_write_gradient_data_(
681 const char *meshName,
682 const char *dataName,
683 const int *size,
684 const int *ids,
685 const double *gradients,
686 int meshNameLength,
687 int dataNameLength);
688
690
693
695PRECICE_API void precicef_start_profiling_section_(const char *sectionName, int sectionNameLength);
696
699
701
702PRECICE_API void precicef_get_version_information_(
703 char *versionInfo,
704 int lengthVersionInfo);
705
706#ifdef __cplusplus
707}
708#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_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_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_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_start_profiling_section_(const char *sectionName, int sectionNameLength)