tblite
Light-weight tight-binding framework
|
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... | |
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.
enum tblite_guess |
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_calculator | ( | tblite_calculator * | calc | ) |
Delete calculator
calc | Calculator instance |
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
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
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
ctx | Context handle |
calc | Calculator instance |
ao2sh | Index mapping from atomic orbitals to shells |
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
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
ctx | Context handle |
calc | Calculator instance |
sh2at | Index mapping from shells to atomic centers |
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
ctx | Context handle |
mol | Molecular structure data |
calc | Calculator instance |
res | Result container |
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
ctx | Context handle |
mol | Molecular structure data |
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
ctx | Context handle |
mol | Molecular structure data |
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
ctx | Context handle |
mol | Molecular structure data |
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
ctx | Context handle |
mol | Molecular structure data |
param | Parametrization records |
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
post_proc | Post Processing instance |
param | Param instance containing post processing information |
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
post_proc | Post Processing instance |
charptr | String of the post processing desired |
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
ctx | Context handle |
calc | Calculator instance |
acc | Accuracy value for numerical thresholds |
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
ctx | Context handle |
calc | Calculator instance |
guess | Guess for initializing the wavefunction |
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
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
ctx | Context handle |
calc | Calculator instance |
damping | Value for mixer damping |
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
ctx | Context handle |
calc | Calculator instance |
save_integrals | Flag to enable storing of integrals |
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)
ctx | Context handle |
calc | Calculator instance |
etemp | Electronic temperature in Hartree |