libtrixi
trixi.h
Go to the documentation of this file.
1 #ifndef TRIXI_H_
2 #define TRIXI_H_
3 
9 // Setup
10 void trixi_initialize(const char * project_directory, const char * depot_path);
11 void trixi_finalize();
12 
13 // Version information
17 const char* trixi_version_library();
18 const char* trixi_version_julia();
19 const char* trixi_version_julia_extended();
20 
21 // Simulation control
22 int trixi_initialize_simulation(const char * libelixir);
23 void trixi_finalize_simulation(int handle);
24 int trixi_is_finished(int handle);
25 void trixi_step(int handle);
26 
27 // Simulation data
28 int trixi_ndims(int handle);
29 int trixi_nelements(int handle);
30 int trixi_nelements_global(int handle);
31 int trixi_nvariables(int handle);
32 double trixi_calculate_dt(int handle);
33 void trixi_load_cell_averages(double * data, int handle);
34 
35 // T8code
36 #if !defined(T8_H) && !defined(T8_FOREST_GENERAL_H)
37 typedef struct t8_forest *t8_forest_t;
38 #endif
40 
41 // Misc
42 void trixi_eval_julia(const char * code);
43 
48 #endif // ifndef LIBTRIXI_H_
void trixi_finalize()
Finalize Julia runtime environment.
Definition: api.c:143
int trixi_nelements_global(int handle)
Return number of global elements (cells).
Definition: api.c:481
const char * trixi_version_julia()
Return name and version of loaded Julia packages LibTrixi directly depends on.
Definition: api.c:276
int trixi_version_library_patch()
Return patch version number of libtrixi.
Definition: api.c:222
double trixi_calculate_dt(int handle)
Get time step length.
Definition: api.c:422
const char * trixi_version_julia_extended()
Return name and version of all loaded Julia packages.
Definition: api.c:303
int trixi_nvariables(int handle)
Return number of (conservative) variables.
Definition: api.c:498
void trixi_load_cell_averages(double *data, int handle)
Return cell averaged values.
Definition: api.c:522
void trixi_initialize(const char *project_directory, const char *depot_path)
Initialize Julia runtime environment.
Definition: api.c:87
void trixi_eval_julia(const char *code)
Execute Julia code.
Definition: api.c:574
int trixi_is_finished(int handle)
Check if simulation is finished.
Definition: api.c:357
int trixi_initialize_simulation(const char *libelixir)
Set up Trixi simulation.
Definition: api.c:335
t8_forest_t trixi_get_t8code_forest(int handle)
Definition: api.c:549
struct t8_forest * t8_forest_t
Definition: trixi.h:37
const char * trixi_version_library()
Return full version string of libtrixi.
Definition: api.c:250
void trixi_finalize_simulation(int handle)
Finalize simulation.
Definition: api.c:395
int trixi_version_library_major()
Return major version number of libtrixi.
Definition: api.c:182
int trixi_nelements(int handle)
Return number of local elements (cells).
Definition: api.c:460
int trixi_version_library_minor()
Return minor version number of libtrixi.
Definition: api.c:202
void trixi_step(int handle)
Perform next simulation step.
Definition: api.c:376
int trixi_ndims(int handle)
Return number of spatial dimensions.
Definition: api.c:439