C API#

The C API bindings are provided by using the iso_c_binding intrinsic module. Generally, objects are exported as opaque pointers and can only be manipulated within the library. The API user is required delete all objects created in the library by using the provided deconstructor functions to avoid memory leaks.

Overall four classes of objects are provided by the library

  • error handlers (tblite_error and tblite_context), used to communicate exceptional conditions and errors from the library to the user, also usable to communicate general setting to the library in case of a context handle

  • structure containers (tblite_structure), used to represent the system specific information and geometry data, only the latter are mutable for the user

  • calculator objects (tblite_calculator), used to store actual method parametrisation

  • result and restart data (tblite_result) used to hold calculated properties and persistent data for restarting between calculations

Note

Generally, all quantities provided to the library are assumed to be in atomic units.

Main header#

#include "tblite.h"

Provides convenience functionality for working with the C API bindings for tblite by including all relevant headers and defining general type generic macros for working with the object handles.

Defines

tblite_delete(ptr)#

Generic macro to free an object handle.

Parameters:
  • ptr – Opaque object handle

Version query#

#include "tblite/version.h"

Provides access to the version, compatibility and features exported by this API.

Functions

TBLITE_API_ENTRY int TBLITE_API_CALL tblite_get_version(void)#

Retrieve version of library used

Returns:

Compact version number in the format 10000 * major + 100 * minor + patch

Macro definitions#

#include "tblite/macros.h"

General macro definitions for the tblite C API bindings.

Defines

TBLITE_API_ENTRY#

Defines an external function exported by the Fortran library.

TBLITE_API_CALL#

Macro to define calling convention.

Error handler#

#include "tblite/error.h"

Provides a light-weight error handler for communicating with the library.

The library provides two different kinds handlers, a light error handle type tblite_error is defined here.

The error handle is used in the context of simple tasks and requires only small overhead to construct and use. It is mainly used in the context of retrieving data or building structures.

Typedefs

typedef struct _tblite_error *tblite_error#

Error instance.

Functions

TBLITE_API_ENTRY tblite_error TBLITE_API_CALL tblite_new_error(void)#

Create new error handle object.

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_error(tblite_error*)#

Delete an error handle object.

TBLITE_API_ENTRY int TBLITE_API_CALL tblite_check_error(tblite_error)#

Check error handle status.

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_clear_error(tblite_error)#

Clear error handle status.

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_error(tblite_error, char*, const int*)#

Get error message from error handle.

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_error(tblite_error, char*, const int*)#

Set error message to error handle.

Context handler#

#include "tblite/context.h"

Provides a persistent configuration object to modify the behaviour of a calculation. Acts as an error handler.

The environment context tblite_context can be considered a persistent setup for all calculations performed with the library, it is usually used together with calculator objects tblite_calculator. While the error handle can only contain a single error, multiple errors can be accumulated in a context object, which allows storing more complex error information like they can occur in an actual calculation.

Typedefs

typedef struct _tblite_context *tblite_context#

Context manager for the library usage.

typedef void (*tblite_logger_callback)(tblite_error, char*, int, void*)#

Define callback function for use in custom logger.

Functions

TBLITE_API_ENTRY tblite_context TBLITE_API_CALL tblite_new_context(void)#

Create new calculation environment object

Returns:

New context handle

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_context(tblite_context *ctx)#

Delete a calculation environment object

Parameters:
  • ctx – Context handle

TBLITE_API_ENTRY int TBLITE_API_CALL tblite_check_context(tblite_context ctx)#

Check calculation environment status

Parameters:
  • ctx – Context handle

Returns:

Current status, 0: okay, 1: error

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_context_error(tblite_context ctx, char *buffer, const int *buffersize)#

Get error message from calculation environment

Parameters:
  • ctx – Context handle

  • buffer – Character buffer for writing error message to

  • buffersize – Maximum length of buffer (optional)

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_context_logger(tblite_context ctx, tblite_logger_callback callback, void *userdata)#

Set custom logger function

Parameters:
  • ctx – Context handle

  • callback – Procedure pointer implementing logger

  • userdata – Passthrough data pointer

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_context_color(tblite_context ctx, int color)#

Enable colorful output

Parameters:
  • ctx – Context handle

  • color – Set color support, 0: disabled, 1: enabled

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_context_verbosity(tblite_context ctx, int verbosity)#

Set verbosity level of printout

Parameters:
  • ctx – Context handle

  • verbosity – Printout verbosity

Molecular structure data#

#include "tblite/structure.h"

The structure data tblite_structure is used to represent the system of interest in the library.

It contains immutable system specific information like the number of atoms, the unique atom groups and the boundary conditions as well as mutable geometry data like cartesian coordinates and lattice parameters.

To change immutable parameters of the tblite_structure data the object should be reinstantiated as well as all dependent objects, like tblite_calculator or tblite_result instances.

Typedefs

typedef struct _tblite_structure *tblite_structure#

Molecular structure data class.

The structure data is used to represent the system of interest in the library. It contains immutable system specific information like the number of atoms, the unique atom groups and the boundary conditions as well as mutable geometry data like cartesian coordinates and lattice parameters.

Functions

TBLITE_API_ENTRY tblite_structure TBLITE_API_CALL tblite_new_structure(tblite_error error, const int natoms, const int *numbers, const double *positions, const double *charge, const int *uhf, const double *lattice, const bool *periodic)#

Create new molecular structure data

Parameters:
  • error – Handle for error messages

  • natoms – Number of atoms

  • numbers – Atomic numbers for each atom, shape [natoms]

  • positions – Cartesian coordinates in Bohr for each atom, shape [natoms][3]

  • charge – Total charge of the system, (optional)

  • uhf – Number of unpaired electrons, (optional)

  • lattice – Lattice parameters in Bohr, shape [3][3], (optional)

  • periodic – Periodic dimensions, shape [3], (optional)

Returns:

New molecular structure data

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_structure(tblite_structure *mol)#

Delete molecular structure data

Parameters:
  • mol – Molecular structure data

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_update_structure_geometry(tblite_error error, tblite_structure mol, const double *positions, const double *lattice)#

Update coordinates and lattice parameters (quantities in Bohr)

Parameters:
  • error – Handle for error messages

  • mol – Molecular structure data

  • positions – Cartesian coordinates in Bohr for each atom, shape [natoms][3]

  • lattice – Lattice parameters in Bohr, shape [3][3], (optional)

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_update_structure_charge(tblite_error error, tblite_structure mol, const double *charge)#

Update total charge in structure object

Parameters:
  • error – Handle for error messages

  • mol – Molecular structure data

  • charge – Total charge of the system

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_update_structure_uhf(tblite_error error, tblite_structure mol, const int *uhf)#

Update number of unpaired electrons in structure object

Parameters:
  • error – Handle for error messages

  • mol – Molecular structure data

  • uhf – Number of unpaired electrons

Calculator instance#

#include "tblite/calculator.h"

Provides a single point calculator for performing extended tight-binding computations.

Provides a parametrized type tblite_calculator storing the method parameters required for the evaluation of single point calculations. The calculator is parametrized for a method as well as for the molecular structure data tblite_structure it was created with. Changes in the tblite_structure required to reinstantiate the structure data also require to reinstantiate the tblite_calculator object.

Typedefs

typedef struct _tblite_calculator *tblite_calculator#

Single point calculator.

Enums

enum tblite_guess#

Possible initial guess of the wavefunction.

Values:

enumerator TBLITE_GUESS_SAD#

Use superposition of atomic densities as guess.

enumerator TBLITE_GUESS_EEQ#

Use partial charges obtained by electronegativity equilibration as guess.

Functions

TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_gfn2_calculator(tblite_context ctx, tblite_structure mol)#

Construct calculator with GFN2-xTB parametrisation loaded

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

Returns:

New calculator instance

TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_gfn1_calculator(tblite_context ctx, tblite_structure mol)#

Construct calculator with GFN1-xTB parametrisation loaded

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

Returns:

New calculator instance

TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_ipea1_calculator(tblite_context ctx, tblite_structure mol)#

Construct calculator with IPEA1-xTB parametrisation loaded

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

Returns:

New calculator instance

TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_xtb_calculator(tblite_context ctx, tblite_structure mol, tblite_param param)#

Construct calculator from parametrization records

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

  • param – Parametrization records

Returns:

New calculator instance

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_calculator(tblite_calculator *calc)#

Delete calculator

Parameters:
  • calc – Calculator instance

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_accuracy(tblite_context ctx, tblite_calculator calc, double acc)#

Set calculation accuracy for the calculator object

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • acc – Accuracy value for numerical thresholds

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_max_iter(tblite_context ctx, tblite_calculator calc, int max_iter)#

Set maximum number of allowed iterations in calculator object

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • max_iter – Maximum number of allowed iterations

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_mixer_damping(tblite_context ctx, tblite_calculator calc, double damping)#

Set damping parameter for mixer in calculator object

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • damping – Value for mixer damping

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_guess(tblite_context ctx, tblite_calculator calc, tblite_guess guess)#

Set initial guess for creating new wavefunction objects

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • guess – Guess for initializing the wavefunction

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_temperature(tblite_context ctx, tblite_calculator calc, double etemp)#

Set electronic temperature for the calculator object (in Hartree)

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • etemp – Electronic temperature in Hartree

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_save_integrals(tblite_context ctx, tblite_calculator calc, int save_integrals)#

Set the flag in the calculator to retain the integral matrices

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • save_integrals – Flag to enable storing of integrals

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_shell_count(tblite_context ctx, tblite_calculator calc, int *nsh)#

Query calculator for the number of shells

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • nsh – Number of shells

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_shell_map(tblite_context ctx, tblite_calculator calc, int *sh2at)#

Query calculator for index mapping from shells to atomic centers

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • sh2at – Index mapping from shells to atomic centers

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_angular_momenta(tblite_context ctx, tblite_calculator calc, int *am)#

Query calculator for angular momenta of shells

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • am – Angular momenta of each shell

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_orbital_count(tblite_context ctx, tblite_calculator calc, int *nao)#

Query calculator for the number of orbitals

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • nao – Number of atomic orbitals

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_orbital_map(tblite_context ctx, tblite_calculator calc, int *ao2sh)#

Query calculator for index mapping from atomic orbitals to shells

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • ao2sh – Index mapping from atomic orbitals to shells

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_singlepoint(tblite_context ctx, tblite_structure mol, tblite_calculator calc, tblite_result res)#

Perform single point calculation

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

  • calc – Calculator instance

  • res – Result container

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_push_back_post_processing_str(tblite_context ctx, tblite_calculator calc, char *charptr)#

Push Back new conatiner to post processing construct

Parameters:
  • post_proc – Post Processing instance

  • charptr – String of the post processing desired

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_push_back_post_processing_param(tblite_context ctx, tblite_calculator calc, tblite_param param)#

Push Back new conatiner to post processing construct

Parameters:
  • post_proc – Post Processing instance

  • param – Param instance containing post processing information

Interaction container#

#include "tblite/container.h"

Provides an interaction container which can be added to a tblite_calculator.

Defines

tblite_new_cpcm_solvation(ctx, mol, calc, x)#
tblite_new_alpb_solvation(ctx, mol, calc, x)#

Typedefs

typedef struct _tblite_container *tblite_container#

Interaction container.

Functions

TBLITE_API_ENTRY tblite_container TBLITE_API_CALL tblite_new_electric_field(double *efield)#

Create new electric field container

Parameters:
  • efield – Electric field in atomic units (Hartree/(Bohr*e)), shape: [3]

Returns:

New interaction container

TBLITE_API_ENTRY tblite_container TBLITE_API_CALL tblite_new_spin_polarization(tblite_context ctx, tblite_structure mol, tblite_calculator calc, double wscale)#

Create new spin polarization container using internal parameters

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

  • calc – Calculator instance

  • wscale – Scaling factor for spin polarization (default: 1)

Returns:

New interaction container

TBLITE_API_ENTRY tblite_container TBLITE_API_CALL tblite_new_cpcm_solvation_epsilon(tblite_context ctx, tblite_structure mol, tblite_calculator calc, double eps)#

Create new CPCM implicit solvation container using internal parameters

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

  • calc – Calculator instance

  • eps – epsilon value for solvent

Returns:

New interaction container

TBLITE_API_ENTRY tblite_container TBLITE_API_CALL tblite_new_cpcm_solvation_solvent(tblite_context ctx, tblite_structure mol, tblite_calculator calc, char *solvent)#

Create new CPCM implicit solvation container using internal parameters

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

  • calc – Calculator instance

  • solvent – Solvent to be modelled, can be given as name of solvent or epsilon value

Returns:

New interaction container

TBLITE_API_ENTRY tblite_container TBLITE_API_CALL tblite_new_alpb_solvation_solvent(tblite_context ctx, tblite_structure mol, tblite_calculator calc, char *solvent)#

Create new ALPB implicit solvation container using internal parameters

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

  • calc – Calculator instance

  • solvent – Solvent to be modelled, can be given as name of solvent or epsilon value

Returns:

New interaction container

TBLITE_API_ENTRY tblite_container TBLITE_API_CALL tblite_new_alpb_solvation_epsilon(tblite_context ctx, tblite_structure mol, tblite_calculator calc, double eps)#

Create new ALPB implicit solvation container using internal parameters

Parameters:
  • ctx – Context handle

  • mol – Molecular structure data

  • calc – Calculator instance

  • eps – epsilon value of solvent

Returns:

New interaction container

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_calculator_push_back(tblite_context ctx, tblite_calculator calc, tblite_container *cont)#

Add container to calculator object.

Note: Ownership is transferred and container handle is destroyed after function call

Parameters:
  • ctx – Context handle

  • calc – Calculator instance

  • cont – Interaction container

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_container(tblite_container *cont)#

Delete container handle

Parameters:
  • cont – Container handle

Result container#

#include "tblite/result.h"

Provides a storage container tblite_result for all persistent quantities produced in a calculation by tblite_get_singlepoint. The data stored in tblite_result can be used as restart data for subsequent calculations.

The individual entries of the tblite_result container can be queried and retrieved. For some quantities, like the count of the shells and orbitals duplicated queries are implemented to allow the usage of the tblite_result container without requiring access to the tblite_calculator which produced it.

Typedefs

typedef struct _tblite_result *tblite_result#

Container to for storing and handling calculation results.

Functions

TBLITE_API_ENTRY tblite_result TBLITE_API_CALL tblite_new_result(void)#

Create new result container

Returns:

New result container

TBLITE_API_ENTRY tblite_result TBLITE_API_CALL tblite_copy_result(tblite_result res)#

Create new result container from existing container

Parameters:
  • res – Existing result container

Returns:

Copy of the results in a new container

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_result(tblite_result *res)#

Delete a calculation environment object

Parameters:
  • res – Result container

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_number_of_atoms(tblite_error error, tblite_result res, int *nat)#

Retrieve number of atoms from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • nat – Number of atoms

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_number_of_shells(tblite_error error, tblite_result res, int *nsh)#

Retrieve number of shells from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • nsh – Number of shells

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_number_of_orbitals(tblite_error error, tblite_result res, int *nao)#

Retrieve number of orbitals from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • nao – Number of orbitals

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_energy(tblite_error error, tblite_result res, double *energy)#

Retrieve energy from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • energy – Total energy

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_energies(tblite_error error, tblite_result res, double *energies)#

Retrieve atom-resolved energies from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • energies – Atom-resolved energies, shape [nat]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_gradient(tblite_error error, tblite_result res, double *gradient)#

Retrieve gradient from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • gradient – Cartesian gradient, shape [nat][3]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_virial(tblite_error error, tblite_result res, double *sigma)#

Retrieve virial from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • sigma – Strain derivatives, shape [3][3]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_charges(tblite_error error, tblite_result res, double *charges)#

Retrieve atomic charges from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • charges – Atomic partial charges, shape [nat]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_bond_orders(tblite_error error, tblite_result res, double *mbo)#

Retrieve bond orders from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • mbo – Bond orders, shape [nat][nat]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_dipole(tblite_error error, tblite_result res, double dipole[3])#

Retrieve dipole moment from result container (order x, y, z)

Parameters:
  • error – Handle for error messages

  • res – Result container

  • dipole – Total dipole moment, shape [3]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_quadrupole(tblite_error error, tblite_result res, double quadrupole[6])#

Retrieve traceless quadrupole moment from result container (packed xx, xy, yy, xz, yz, zz)

Parameters:
  • error – Handle for error messages

  • res – Result container

  • quadrupole – Total traceless quadrupole moment, shape [6]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_orbital_energies(tblite_error error, tblite_result res, double *emo)#

Retrieve orbital energies from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • emo – Eigenvalues for each orbital

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_orbital_occupations(tblite_error error, tblite_result res, double *occ)#

Retrieve orbital occupations from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • occ – Occupation numbers for each orbital

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_orbital_coefficients(tblite_error error, tblite_result res, double *cmat)#

Retrieve orbital coefficients from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • cmat – Orbital coefficient matrix, shape [nspin][nao][nao]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_density_matrix(tblite_error error, tblite_result res, double *pmat)#

Retrieve density matrix from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • pmat – Density matrix, shape [nspin][nao][nao]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_overlap_matrix(tblite_error error, tblite_result res, double *smat)#

Retrieve overlap matrix from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • smat – Overlap matrix, shape [nao][nao]

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_hamiltonian_matrix(tblite_error error, tblite_result res, double *hmat)#

Retrieve Hamiltonian matrix from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • hmat – Hamiltonian matrix, shape [nao][nao]

TBLITE_API_ENTRY tblite_double_dictionary TBLITE_API_CALL tblite_get_post_processing_dict(tblite_error error, tblite_result res)#

Retrieve Hamiltonian matrix from result container

Parameters:
  • error – Handle for error messages

  • res – Result container

  • dict – Pointer to dictionary object

Parametrization records#

#include "tblite/param.h"

Provides a representation of a parametrization of an xTB Hamiltonian which can be used to instantiate a tblite_calculator object.

The parametrization data itself can be represented as a tblite_table data structure, which provides the possibility to customize the parametrization programmatically.

Typedefs

typedef struct _tblite_param *tblite_param#

Parametrization records.

Functions

TBLITE_API_ENTRY tblite_param TBLITE_API_CALL tblite_new_param(void)#

Create new parametrization records object

Returns:

Parametrization records

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_param(tblite_param *param)#

Delete a parametrization records object

Parameters:
  • param – Parametrization records

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_load_param(tblite_error error, tblite_param param, tblite_table table)#

Load parametrization records from data table

Parameters:
  • error – Handle for error messages

  • param – Parametrization records

  • table – Table data structure

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_dump_param(tblite_error error, tblite_param param, tblite_table table)#

Dump parametrization records to data table

Parameters:
  • error – Handle for error messages

  • param – Parametrization records

  • table – Table data structure

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_export_gfn2_param(tblite_error error, tblite_param param)#

Export GFN2-xTB parametrization records

Parameters:
  • error – Handle for error messages

  • param – Parametrization records

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_export_gfn1_param(tblite_error error, tblite_param param)#

Export GFN1-xTB parametrization records

Parameters:
  • error – Handle for error messages

  • param – Parametrization records

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_export_ipea1_param(tblite_error error, tblite_param param)#

Export IPEA1-xTB parametrization records

Parameters:
  • error – Handle for error messages

  • param – Parametrization records

Table data structure#

#include "tblite/table.h"

Provides a representation of a generic table data structure.

Used to mirror the data available in the tblite_param object. It aims to provide a programmatic accessible representation of the parametrization records.

Defines

tblite_table_set_value(error, table, key, value, ...)#

Generic setter based on the type of the value

Parameters:
  • error – Error handle

  • table – Table data structure

  • key – Key to set value at

  • value – Value to set at key

Typedefs

typedef struct _tblite_table *tblite_table#

Handle for holding a table data structure

The table can either own its data or reference another table. Table references can be created from existing table data structures or by adding new table entries into an existing table, which returns a reference to the newly created table.

Functions

TBLITE_API_ENTRY tblite_table TBLITE_API_CALL tblite_new_table(tblite_table *table)#

Create new data table object

Parameters:
  • table – Table object to reference in new table (optional)

Returns:

: New table data structure

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_table(tblite_table *table)#

Delete a data table object

Parameters:
  • table – Table object to be deleted

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_double(tblite_error error, tblite_table table, char key[], double *value, int n)#

Set floating point number to data table

Parameters:
  • error – Error handle

  • table – Table data structure

  • key – Key to set value at

  • value – Double value or array to add to table

  • n – Number of entries to add, 0 for adding scalars

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_int64_t(tblite_error error, tblite_table table, char key[], int64_t *value, int n)#

Set integer number to data table (use int64_t rather than long)

Parameters:
  • error – Error handle

  • table – Table data structure

  • key – Key to set value at

  • value – Integer value or array to add to table

  • n – Number of entries to add, 0 for adding scalars

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_bool(tblite_error error, tblite_table table, char key[], bool *value, int n)#

Set boolean value to data table

Parameters:
  • error – Error handle

  • table – Table data structure

  • key – Key to set value at

  • value – Boolean value or array to add to table

  • n – Number of entries to add, 0 for adding scalars

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_char(tblite_error error, tblite_table table, char key[], char (*value)[], int n)#

Set character string to data table

Parameters:
  • error – Error handle

  • table – Table data structure

  • key – Key to set value at

  • value – Character value or array to add to table

  • n – Number of entries to add, 0 for adding scalars

TBLITE_API_ENTRY tblite_table TBLITE_API_CALL tblite_table_add_table(tblite_error error, tblite_table table, char key[])#

Create new subtable in existing data table

Parameters:
  • error – Error handle

  • table – Table data structure

  • key – Key to add new subtable at

Returns:

New table data structure