libtrixi
Typedefs | Functions
C API

Typedefs

typedef struct t8_forest * t8_forest_t
 

Functions

void trixi_initialize (const char *project_directory, const char *depot_path)
 Initialize Julia runtime environment. More...
 
void trixi_finalize ()
 Finalize Julia runtime environment. More...
 
int trixi_version_library_major ()
 Return major version number of libtrixi. More...
 
int trixi_version_library_minor ()
 Return minor version number of libtrixi. More...
 
int trixi_version_library_patch ()
 Return patch version number of libtrixi. More...
 
const char * trixi_version_library ()
 Return full version string of libtrixi. More...
 
const char * trixi_version_julia ()
 Return name and version of loaded Julia packages LibTrixi directly depends on. More...
 
const char * trixi_version_julia_extended ()
 Return name and version of all loaded Julia packages. More...
 
int trixi_initialize_simulation (const char *libelixir)
 Set up Trixi simulation. More...
 
void trixi_finalize_simulation (int handle)
 Finalize simulation. More...
 
int trixi_is_finished (int handle)
 Check if simulation is finished. More...
 
void trixi_step (int handle)
 Perform next simulation step. More...
 
int trixi_ndims (int handle)
 Return number of spatial dimensions. More...
 
int trixi_nelements (int handle)
 Return number of local elements. More...
 
int trixi_nelementsglobal (int handle)
 Return global number of elements. More...
 
int trixi_ndofs (int handle)
 Return number of local degrees of freedom. More...
 
int trixi_ndofsglobal (int handle)
 Return global number of degrees of freedom. More...
 
int trixi_ndofselement (int handle)
 Return number of degrees of freedom per element. More...
 
int trixi_nvariables (int handle)
 Return number of (conservative) variables. More...
 
int trixi_nnodes (int handle)
 Return number of quadrature nodes per dimension. More...
 
double trixi_calculate_dt (int handle)
 Get time step length. More...
 
void trixi_load_node_reference_coordinates (int handle, double *node_coords)
 Get reference coordinates of 1D quadrature nodes. More...
 
void trixi_load_node_weights (int handle, double *node_weights)
 Get weights of 1D quadrature nodes. More...
 
void trixi_load_primitive_vars (int handle, int variable_id, double *data)
 Load primitive variable. More...
 
void trixi_load_element_averaged_primitive_vars (int handle, int variable_id, double *data)
 Load element averages for primitive variable. More...
 
t8_forest_t trixi_get_t8code_forest (int handle)
 
void trixi_eval_julia (const char *code)
 Execute Julia code. More...
 

Detailed Description

Typedef Documentation

◆ t8_forest_t

typedef struct t8_forest* t8_forest_t

Function Documentation

◆ trixi_calculate_dt()

double trixi_calculate_dt ( int  handle)

Get time step length.

Get the current time step length of the simulation identified by handle.

Parameters
[in]handlesimulation handle
Returns
Time step length

◆ trixi_eval_julia()

void trixi_eval_julia ( const char *  code)

Execute Julia code.

Execute the provided code in the current Julia runtime environment.

Warning
Only for development. Code is not checked prior to execution.

◆ trixi_finalize()

void trixi_finalize ( )

Finalize Julia runtime environment.

Clean up internal states. This function should be executed near the end of the process' lifetime. After the call to trixi_finalize, no other libtrixi functions may be called anymore, including trixi_finalize itself.

◆ trixi_finalize_simulation()

void trixi_finalize_simulation ( int  handle)

Finalize simulation.

Finalize the simulation identified by handle. This will also release the handle.

Parameters
[in]handlesimulation handle

◆ trixi_get_t8code_forest()

t8_forest_t trixi_get_t8code_forest ( int  handle)

Get t8code forest

For Trixi simulations on t8code meshes, the t8code forest is returned.

Parameters
[in]handlesimulation handle
Warning
The interface to t8code is experimental and implementation details may change at any time without warning.
Returns
t8code forest

◆ trixi_initialize()

void trixi_initialize ( const char *  project_directory,
const char *  depot_path 
)

Initialize Julia runtime environment.

Initialize Julia and activate the project at project_directory. If depot_path is not a null pointer, forcefully set the environment variable JULIA_DEPOT_PATH to the value of depot_path. If depot_path is null, then proceed as follows: If JULIA_DEPOT_PATH is already set, do not touch it. Otherwise, set JULIA_DEPOT_PATH to project_directory + default_depot_path

This function must be called before most other libtrixi functions can be used. Libtrixi maybe only be initialized once; subsequent calls to trixi_initialize are erroneous.

Parameters
[in]project_directoryPath to project directory.
[in]depot_pathPath to Julia depot path (optional; can be null pointer).

◆ trixi_initialize_simulation()

int trixi_initialize_simulation ( const char *  libelixir)

Set up Trixi simulation.

Set up a Trixi simulation by reading the provided libelixir file. It resembles Trixi's typical elixir files with the following differences:

  • Everything (except using ...) has to be inside a function init_simstate()
  • OrdinaryDiffEq's integrator has to be created by calling init (instead of solve)
  • A SimulationState has to be created from the semidiscretization and the integrator See the examples in the LibTrixi.jl/examples folder
Parameters
[in]libelixirPath to libelexir file.
Returns
handle (integer) to Trixi simulation instance

◆ trixi_is_finished()

int trixi_is_finished ( int  handle)

Check if simulation is finished.

Checks if the simulation identified by handle has reached its final time.

Parameters
[in]handlesimulation handle
Returns
1 if finished, 0 if not

◆ trixi_load_element_averaged_primitive_vars()

void trixi_load_element_averaged_primitive_vars ( int  handle,
int  variable_id,
double *  data 
)

Load element averages for primitive variable.

Element averaged values for the primitive variable at position variable_id for each element are stored in the given array data.

The given array has to be of correct size (nelements) and memory has to be allocated beforehand.

Parameters
[in]handlesimulation handle
[in]variable_idindex of variable
[out]dataelement averaged values for all elements

◆ trixi_load_node_reference_coordinates()

void trixi_load_node_reference_coordinates ( int  handle,
double *  node_coords 
)

Get reference coordinates of 1D quadrature nodes.

The reference coordinates in [-1,1] of the quadrature nodes in the current DG scheme are stored in the provided array node_coords. The given array has to be of correct size, i.e. nnodes, and memory has to be allocated beforehand.

Parameters
[in]handlesimulation handle
[out]node_coordsnode reference coordinates

◆ trixi_load_node_weights()

void trixi_load_node_weights ( int  handle,
double *  node_weights 
)

Get weights of 1D quadrature nodes.

The weights of the quadrature nodes in the current DG scheme are stored in the provided array node_weights. The given array has to be of correct size, i.e. nnodes, and memory has to be allocated beforehand.

Parameters
[in]handlesimulation handle
[out]node_weightsnode weights

◆ trixi_load_primitive_vars()

void trixi_load_primitive_vars ( int  handle,
int  variable_id,
double *  data 
)

Load primitive variable.

The values for the primitive variable at position variable_id at every degree of freedom are stored in the given array data.

The given array has to be of correct size (ndofs) and memory has to be allocated beforehand.

Parameters
[in]handlesimulation handle
[in]variable_idindex of variable
[out]datavalues for all degrees of freedom

◆ trixi_ndims()

int trixi_ndims ( int  handle)

Return number of spatial dimensions.

Parameters
[in]handlesimulation handle

◆ trixi_ndofs()

int trixi_ndofs ( int  handle)

Return number of local degrees of freedom.

These usually differ from the global count when doing parallel computations.

Parameters
[in]handlesimulation handle
See also
trixi_ndofsglobal_api_c

◆ trixi_ndofselement()

int trixi_ndofselement ( int  handle)

Return number of degrees of freedom per element.

Parameters
[in]handlesimulation handle

◆ trixi_ndofsglobal()

int trixi_ndofsglobal ( int  handle)

Return global number of degrees of freedom.

These usually differ from the local count when doing parallel computations.

Parameters
[in]handlesimulation handle
See also
trixi_ndofs_api_c

◆ trixi_nelements()

int trixi_nelements ( int  handle)

Return number of local elements.

These usually differ from the global count when doing parallel computations.

Parameters
[in]handlesimulation handle
See also
trixi_nelementsglobal_api_c

◆ trixi_nelementsglobal()

int trixi_nelementsglobal ( int  handle)

Return global number of elements.

These usually differ from the local count when doing parallel computations.

Parameters
[in]handlesimulation handle
See also
trixi_nelements_api_c

◆ trixi_nnodes()

int trixi_nnodes ( int  handle)

Return number of quadrature nodes per dimension.

Parameters
[in]handlesimulation handle

◆ trixi_nvariables()

int trixi_nvariables ( int  handle)

Return number of (conservative) variables.

Parameters
[in]handlesimulation handle

◆ trixi_step()

void trixi_step ( int  handle)

Perform next simulation step.

Let the simulation identified by handle advance by one step.

Parameters
[in]handlesimulation handle

◆ trixi_version_julia()

const char* trixi_version_julia ( )

Return name and version of loaded Julia packages LibTrixi directly depends on.

The return value is a read-only pointer to a NULL-terminated string with the name and version information of the loaded Julia packages, separated by newlines.

The returned pointer is to static memory and must not be used to change the contents of the version string. Multiple calls to the function will return the same address.

This function is thread-safe. It must be run after trixi_initialize has been called.

Returns
Pointer to a read-only, NULL-terminated character array with the names and versions of loaded Julia packages.

◆ trixi_version_julia_extended()

const char* trixi_version_julia_extended ( )

Return name and version of all loaded Julia packages.

The return value is a read-only pointer to a NULL-terminated string with the name and version information of all loaded Julia packages, including implicit dependencies, separated by newlines.

The returned pointer is to static memory and must not be used to change the contents of the version string. Multiple calls to the function will return the same address.

This function is thread-safe. It must be run after trixi_initialize has been called.

Returns
Pointer to a read-only, NULL-terminated character array with the names and versions of all loaded Julia packages.

◆ trixi_version_library()

const char* trixi_version_library ( )

Return full version string of libtrixi.

The return value is a read-only pointer to a NULL-terminated string with the version information. This may include not just MAJOR.MINOR.PATCH but possibly also additional build or development version information.

The returned pointer is to static memory and must not be used to change the contents of the version string. Multiple calls to the function will return the same address.

This function is thread-safe and may be run before trixi_initialize has been called.

Returns
Pointer to a read-only, NULL-terminated character array with the full version of libtrixi.

◆ trixi_version_library_major()

int trixi_version_library_major ( )

Return major version number of libtrixi.

This function may be run before trixi_initialize has been called.

Returns
Major version of libtrixi.

◆ trixi_version_library_minor()

int trixi_version_library_minor ( )

Return minor version number of libtrixi.

This function may be run before trixi_initialize has been called.

Returns
Minor version of libtrixi.

◆ trixi_version_library_patch()

int trixi_version_library_patch ( )

Return patch version number of libtrixi.

This function may be run before trixi_initialize has been called.

Returns
Patch version of libtrixi.