preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
preciceFortran.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <precice/export.h>
4
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
44
47
58PRECICE_API void precicef_initialize_();
59
70PRECICE_API void precicef_advance_(const double *timeStepSize);
71
79PRECICE_API void precicef_finalize_();
80
82
85
97 int *isRequired);
98
110 int *isRequired);
111
113
116
129PRECICE_API void precicef_get_mesh_dimensions_(
130 const char *meshName,
131 int *dimensions,
132 int meshNameLength);
133
147PRECICE_API void precicef_get_data_dimensions_(
148 const char *meshName,
149 const char *dataName,
150 int *dimensions,
151 int meshNameLength,
152 int dataNameLength);
153
164PRECICE_API void precicef_is_coupling_ongoing_(int *isOngoing);
165
176PRECICE_API void precicef_is_time_window_complete_(int *isComplete);
177
188PRECICE_API void precicef_get_max_time_step_size_(double *maxTimeStepSize);
189
191
195
208 const char *meshName,
209 int *required,
210 int meshNameLength);
211
221PRECICE_API void precicef_reset_mesh_(
222 const char *meshName,
223 int meshNameLength);
224
238PRECICE_API void precicef_set_vertex_(
239 const char *meshName,
240 const double *coordinates,
241 int *id,
242 int meshNameLength);
243
256PRECICE_API void precicef_get_mesh_vertex_size_(
257 const char *meshName,
258 int *meshSize,
259 int meshNameLength);
260
275PRECICE_API void precicef_set_vertices_(
276 const char *meshName,
277 const int *size,
278 double *coordinates,
279 int *ids,
280 int meshNameLength);
281
295PRECICE_API void precicef_set_edge_(
296 const char *meshName,
297 const int *firstVertexID,
298 const int *secondVertexID,
299 int meshNameLength);
300
314PRECICE_API void precicef_set_mesh_edges_(
315 const char *meshName,
316 const int *size,
317 const int *ids,
318 int meshNameLength);
319
334PRECICE_API void precicef_set_triangle_(
335 const char *meshName,
336 const int *firstVertexID,
337 const int *secondVertexID,
338 const int *thirdVertexID,
339 int meshNameLength);
340
354PRECICE_API void precicef_set_mesh_edges_(
355 const char *meshName,
356 const int *size,
357 const int *ids,
358 int meshNameLength);
359
375PRECICE_API void precicef_set_quad_(
376 const char *meshName,
377 const int *firstVertexID,
378 const int *secondVertexID,
379 const int *thirdVertexID,
380 const int *fourthVertexID,
381 int meshNameLength);
382
396PRECICE_API void precicef_set_mesh_quads_(
397 const char *meshName,
398 const int *size,
399 const int *ids,
400 int meshNameLength);
401
417PRECICE_API void precicef_set_tetrahedron(
418 const char *meshName,
419 const int *firstVertexID,
420 const int *secondVertexID,
421 const int *thirdVertexID,
422 const int *fourthVertexID,
423 int meshNameLength);
424
438PRECICE_API void precicef_set_mesh_tetrahedra_(
439 const char *meshName,
440 const int *size,
441 const int *ids,
442 int meshNameLength);
443
445
448
459PRECICE_API void precicef_requires_initial_data_(
460 int *isRequired);
461
477PRECICE_API void precicef_write_data_(
478 const char *meshName,
479 const char *dataName,
480 const int *size,
481 int *ids,
482 double *values,
483 int meshNameLength,
484 int dataNameLength);
485
502PRECICE_API void precicef_read_data_(
503 const char *meshName,
504 const char *dataName,
505 const int *size,
506 int *ids,
507 const double *relativeReadTime,
508 double *values,
509 int meshNameLength,
510 int dataNameLength);
511
513
516
536PRECICE_API void precicef_write_and_map_data_(
537 const char *meshName,
538 const char *dataName,
539 const int *size,
540 double *coordinates,
541 double *values,
542 int meshNameLength,
543 int dataNameLength);
544
566PRECICE_API void precicef_map_and_read_data_(
567 const char *meshName,
568 const char *dataName,
569 const int *size,
570 double *coordinates,
571 const double *relativeReadTime,
572 double *values,
573 int meshNameLength,
574 int dataNameLength);
576
579
591PRECICE_API void precicef_set_mesh_access_region_(
592 const char *meshName,
593 const double *boundingBox,
594 int meshNameLength);
595
610 const char *meshName,
611 const int size,
612 int *ids,
613 double *coordinates,
614 int meshNameLength);
615
617
622
636 const char *meshName,
637 const char *dataName, int *required,
638 int meshNameLength,
639 int dataNameLength);
640
655PRECICE_API void precicef_write_gradient_data_(
656 const char *meshName,
657 const char *dataName,
658 const int *size,
659 const int *ids,
660 const double *gradients,
661 int meshNameLength,
662 int dataNameLength);
663
665
668
670PRECICE_API void precicef_start_profiling_section_(const char *sectionName, int sectionNameLength);
671
674
676
677PRECICE_API void precicef_get_version_information_(
678 char *versionInfo,
679 int lengthVersionInfo);
680
681#ifdef __cplusplus
682}
683#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_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_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_get_mesh_vertex_ids_and_coordinates_(const char *meshName, const int size, int *ids, double *coordinates, 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_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)