libtrixi
Loading...
Searching...
No Matches
trixi.h
Go to the documentation of this file.
1#ifndef TRIXI_H_
2#define TRIXI_H_
3
9// Setup
10void trixi_initialize(const char * project_directory, const char * depot_path);
11void trixi_finalize();
12
13// Version information
17const char* trixi_version_library();
18const char* trixi_version_julia();
20
21// Simulation control
22int trixi_initialize_simulation(const char * libelixir);
23void trixi_finalize_simulation(int handle);
24int trixi_is_finished(int handle);
25void trixi_step(int handle);
26
27// Simulation data
28int trixi_ndims(int handle);
29int trixi_nelements(int handle);
30int trixi_nelementsglobal(int handle);
31int trixi_ndofs(int handle);
32int trixi_ndofsglobal(int handle);
33int trixi_ndofselement(int handle);
34int trixi_nvariables(int handle);
35int trixi_nnodes(int handle);
36double trixi_calculate_dt(int handle);
37double trixi_get_simulation_time(int handle);
38void trixi_load_node_reference_coordinates(int handle, double* node_coords);
39void trixi_load_node_weights(int handle, double* node_weights);
40void trixi_load_conservative_var(int handle, int variable_id, double * data);
41void trixi_load_primitive_var(int handle, int variable_id, double * data);
42void trixi_load_element_averaged_primitive_var(int handle, int variable_id, double * data);
43void trixi_store_conservative_var(int handle, int variable_id, double * data);
44void trixi_register_data(int handle, int index, int size, const double * data);
45double * trixi_get_conservative_vars_pointer(int handle);
46
47// T8code
48#if !defined(T8_H) && !defined(T8_FOREST_GENERAL_H)
49typedef struct t8_forest *t8_forest_t;
50#endif
52
53// Misc
54void trixi_eval_julia(const char * code);
55
60#endif // ifndef LIBTRIXI_H_
void trixi_finalize()
Finalize Julia runtime environment.
Definition api.c:167
int trixi_version_library_patch()
Return patch version number of libtrixi.
Definition api.c:246
double trixi_calculate_dt(int handle)
Get time step length.
Definition api.c:446
double trixi_get_simulation_time(int handle)
Return current physical time.
Definition api.c:819
int trixi_nvariables(int handle)
Return number of (conservative) variables.
Definition api.c:581
void trixi_initialize(const char *project_directory, const char *depot_path)
Initialize Julia runtime environment.
Definition api.c:111
int trixi_ndofs(int handle)
Return number of local degrees of freedom.
Definition api.c:526
void trixi_eval_julia(const char *code)
Execute Julia code.
Definition api.c:871
int trixi_nelementsglobal(int handle)
Return global number of elements.
Definition api.c:505
void trixi_load_conservative_var(int handle, int variable_id, double *data)
Load conservative variable.
Definition api.c:667
int trixi_is_finished(int handle)
Check if simulation is finished.
Definition api.c:381
int trixi_initialize_simulation(const char *libelixir)
Set up Trixi simulation.
Definition api.c:359
t8_forest_t trixi_get_t8code_forest(int handle)
Definition api.c:846
int trixi_ndofselement(int handle)
Return number of degrees of freedom per element.
Definition api.c:564
struct t8_forest * t8_forest_t
Definition trixi.h:49
int trixi_ndofsglobal(int handle)
Return global number of degrees of freedom.
Definition api.c:547
void trixi_finalize_simulation(int handle)
Finalize simulation.
Definition api.c:419
int trixi_version_library_major()
Return major version number of libtrixi.
Definition api.c:206
int trixi_nelements(int handle)
Return number of local elements.
Definition api.c:484
void trixi_load_node_weights(int handle, double *node_weights)
Get weights of 1D quadrature nodes.
Definition api.c:642
const char * trixi_version_julia_extended()
Return name and version of all loaded Julia packages.
Definition api.c:327
const char * trixi_version_library()
Return full version string of libtrixi.
Definition api.c:274
void trixi_register_data(int handle, int index, int size, const double *data)
Store data vector in current simulation's registry.
Definition api.c:775
int trixi_version_library_minor()
Return minor version number of libtrixi.
Definition api.c:226
void trixi_step(int handle)
Perform next simulation step.
Definition api.c:400
const char * trixi_version_julia()
Return name and version of loaded Julia packages LibTrixi directly depends on.
Definition api.c:300
void trixi_store_conservative_var(int handle, int variable_id, double *data)
Store conservative variable.
Definition api.c:744
double * trixi_get_conservative_vars_pointer(int handle)
Return pointer to internal data vector.
Definition api.c:799
void trixi_load_primitive_var(int handle, int variable_id, double *data)
Load primitive variable.
Definition api.c:693
int trixi_ndims(int handle)
Return number of spatial dimensions.
Definition api.c:463
int trixi_nnodes(int handle)
Return number of quadrature nodes per dimension.
Definition api.c:598
void trixi_load_element_averaged_primitive_var(int handle, int variable_id, double *data)
Load element averages for primitive variable.
Definition api.c:719
void trixi_load_node_reference_coordinates(int handle, double *node_coords)
Get reference coordinates of 1D quadrature nodes.
Definition api.c:620