tblite
Light-weight tight-binding framework
Loading...
Searching...
No Matches
calculator.h File Reference

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

#include "tblite/macros.h"
#include "tblite/context.h"
#include "tblite/structure.h"
#include "tblite/result.h"
#include "tblite/param.h"

Go to the source code of this file.

Typedefs

typedef struct _tblite_calculator * tblite_calculator
 Single point calculator.
 

Enumerations

enum  tblite_guess { TBLITE_GUESS_SAD = 0 , TBLITE_GUESS_EEQ = 1 }
 Possible initial guess of the wavefunction. More...
 

Functions

TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_gfn2_calculator (tblite_context ctx, tblite_structure mol)
 
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_gfn1_calculator (tblite_context ctx, tblite_structure mol)
 
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_ipea1_calculator (tblite_context ctx, tblite_structure mol)
 
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_xtb_calculator (tblite_context ctx, tblite_structure mol, tblite_param param)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_calculator (tblite_calculator *calc)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_accuracy (tblite_context ctx, tblite_calculator calc, double acc)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_max_iter (tblite_context ctx, tblite_calculator calc, int max_iter)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_mixer_damping (tblite_context ctx, tblite_calculator calc, double damping)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_guess (tblite_context ctx, tblite_calculator calc, tblite_guess guess)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_temperature (tblite_context ctx, tblite_calculator calc, double etemp)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_save_integrals (tblite_context ctx, tblite_calculator calc, int save_integrals)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_shell_count (tblite_context ctx, tblite_calculator calc, int *nsh)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_shell_map (tblite_context ctx, tblite_calculator calc, int *sh2at)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_angular_momenta (tblite_context ctx, tblite_calculator calc, int *am)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_orbital_count (tblite_context ctx, tblite_calculator calc, int *nao)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_orbital_map (tblite_context ctx, tblite_calculator calc, int *ao2sh)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_singlepoint (tblite_context ctx, tblite_structure mol, tblite_calculator calc, tblite_result res)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_push_back_post_processing_str (tblite_context ctx, tblite_calculator calc, char *charptr)
 
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_push_back_post_processing_param (tblite_context ctx, tblite_calculator calc, tblite_param param)
 

Detailed Description

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.

Enumeration Type Documentation

◆ tblite_guess

Possible initial guess of the wavefunction.

Enumerator
TBLITE_GUESS_SAD 

Use superposition of atomic densities as guess.

TBLITE_GUESS_EEQ 

Use partial charges obtained by electronegativity equilibration as guess.

Function Documentation

◆ tblite_delete_calculator()

TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_calculator ( tblite_calculator * calc)

Delete calculator

Parameters
calcCalculator instance

◆ tblite_get_calculator_angular_momenta()

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
ctxContext handle
calcCalculator instance
amAngular momenta of each shell

◆ tblite_get_calculator_orbital_count()

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
ctxContext handle
calcCalculator instance
naoNumber of atomic orbitals

◆ tblite_get_calculator_orbital_map()

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
ctxContext handle
calcCalculator instance
ao2shIndex mapping from atomic orbitals to shells

◆ tblite_get_calculator_shell_count()

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
ctxContext handle
calcCalculator instance
nshNumber of shells

◆ tblite_get_calculator_shell_map()

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
ctxContext handle
calcCalculator instance
sh2atIndex mapping from shells to atomic centers

◆ tblite_get_singlepoint()

Perform single point calculation

Parameters
ctxContext handle
molMolecular structure data
calcCalculator instance
resResult container

◆ tblite_new_gfn1_calculator()

Construct calculator with GFN1-xTB parametrisation loaded

Parameters
ctxContext handle
molMolecular structure data
Returns
New calculator instance

◆ tblite_new_gfn2_calculator()

Construct calculator with GFN2-xTB parametrisation loaded

Parameters
ctxContext handle
molMolecular structure data
Returns
New calculator instance

◆ tblite_new_ipea1_calculator()

Construct calculator with IPEA1-xTB parametrisation loaded

Parameters
ctxContext handle
molMolecular structure data
Returns
New calculator instance

◆ tblite_new_xtb_calculator()

Construct calculator from parametrization records

Parameters
ctxContext handle
molMolecular structure data
Returns
New calculator instance
Parameters
paramParametrization records

◆ tblite_push_back_post_processing_param()

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_procPost Processing instance
paramParam instance containing post processing information

◆ tblite_push_back_post_processing_str()

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_procPost Processing instance
charptrString of the post processing desired

◆ tblite_set_calculator_accuracy()

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
ctxContext handle
calcCalculator instance
accAccuracy value for numerical thresholds

◆ tblite_set_calculator_guess()

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
ctxContext handle
calcCalculator instance
guessGuess for initializing the wavefunction

◆ tblite_set_calculator_max_iter()

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
ctxContext handle
calcCalculator instance
max_iterMaximum number of allowed iterations

◆ tblite_set_calculator_mixer_damping()

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
ctxContext handle
calcCalculator instance
dampingValue for mixer damping

◆ tblite_set_calculator_save_integrals()

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
ctxContext handle
calcCalculator instance
save_integralsFlag to enable storing of integrals

◆ tblite_set_calculator_temperature()

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
ctxContext handle
calcCalculator instance
etempElectronic temperature in Hartree