tblite
Light-weight tight-binding framework
Loading...
Searching...
No Matches
calculator.h
Go to the documentation of this file.
1/* This file is part of tblite.
2 * SPDX-Identifier: LGPL-3.0-or-later
3 *
4 * tblite is free software: you can redistribute it and/or modify it under
5 * the terms of the GNU Lesser General Public License as published by
6 * the Free Software Foundation, either version 3 of the License, or
7 * (at your option) any later version.
8 *
9 * tblite is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU Lesser General Public License for more details.
13 *
14 * You should have received a copy of the GNU Lesser General Public License
15 * along with tblite. If not, see <https://www.gnu.org/licenses/>.
16**/
17
29#pragma once
30
31#include "tblite/macros.h"
32#include "tblite/context.h"
33#include "tblite/structure.h"
34#include "tblite/result.h"
35#include "tblite/param.h"
36
44
46typedef struct _tblite_calculator* tblite_calculator;
47
56
65
74
84 tblite_param param);
85
91
100 double acc);
101
110 int max_iter);
111
120 double damping);
121
130 tblite_guess guess);
131
140 double etemp);
141
150 int save_integrals);
151
160 int* nsh);
161
170 int* sh2at);
171
180 int* am);
181
190 int* nao);
191
200 int* ao2sh);
201
212 tblite_result res);
213
220 tblite_calculator calc,
221 char* charptr);
222
230 tblite_param param);
tblite_guess
Possible initial guess of the wavefunction.
Definition calculator.h:38
@ TBLITE_GUESS_SAD
Use superposition of atomic densities as guess.
Definition calculator.h:40
@ TBLITE_GUESS_EEQ
Use partial charges obtained by electronegativity equilibration as guess.
Definition calculator.h:42
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_angular_momenta(tblite_context ctx, tblite_calculator calc, int *am)
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_set_calculator_guess(tblite_context ctx, tblite_calculator calc, tblite_guess guess)
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_max_iter(tblite_context ctx, tblite_calculator calc, int max_iter)
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_push_back_post_processing_param(tblite_context ctx, tblite_calculator calc, tblite_param param)
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_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)
struct _tblite_calculator * tblite_calculator
Single point calculator.
Definition calculator.h:46
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_delete_calculator(tblite_calculator *calc)
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_ipea1_calculator(tblite_context ctx, tblite_structure mol)
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_set_calculator_temperature(tblite_context ctx, tblite_calculator calc, double etemp)
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 void TBLITE_API_CALL tblite_get_calculator_shell_map(tblite_context ctx, tblite_calculator calc, int *sh2at)
Provides a persistent configuration object to modify the behaviour of a calculation....
General macro definitions for the tblite C API bindings.
#define TBLITE_API_ENTRY
Definition macros.h:49
#define TBLITE_API_CALL
Definition macros.h:57
Proxy module for the calculation context related data types.
Definition context.f90:59
Defininition of the parametrization data format.
Definition param.f90:25
Provides a representation of a parametrization of an xTB Hamiltonian which can be used to instantiate...
Provides a storage container tblite_result for all persistent quantities produced in a calculation by...
struct _tblite_result * tblite_result
Container to for storing and handling calculation results.
Definition result.h:37
The structure data tblite_structure is used to represent the system of interest in the library.
struct _tblite_structure * tblite_structure
Definition structure.h:47