Python API#
Python API of the tblite project
Note
The documentation of the ASE integration can be found in Usage in ASE.
Library interface#
Definition of the basic interface to library for most further integration in other Python frameworks. The classes defined here allow a more Pythonic usage of the library in actual workflows than the low-level access provided in the CFFI generated wrappers.
Structure#
- class tblite.interface.Structure(numbers: numpy.ndarray | List[int], positions: numpy.ndarray, charge: float | None = None, uhf: int | None = None, lattice: numpy.ndarray | None = None, periodic: numpy.ndarray | None = None)[source]#
Represents a wrapped structure object in
tblite
. The molecular structure data object has a fixed number of atoms and immutable atomic identifiers.Example
>>> from tblite.interface import Structure >>> import numpy as np >>> mol = Structure( ... positions=np.array([ ... [+0.00000000000000, +0.00000000000000, -0.73578586109551], ... [+1.44183152868459, +0.00000000000000, +0.36789293054775], ... [-1.44183152868459, +0.00000000000000, +0.36789293054775], ... ]), ... numbers = np.array([8, 1, 1]), ... ) ... >>> len(mol) 3
- Raises:
TBLiteValueError – on invalid input, like incorrect shape / type of the passed arrays
- update(positions: numpy.ndarray | None = None, lattice: numpy.ndarray | None = None, charge: float | None = None, uhf: int | None = None) None [source]#
Update coordinates and lattice parameters, both provided in atomic units (Bohr). The lattice update is optional also for periodic structures.
Generally, only the cartesian coordinates and the lattice parameters can be updated, every other modification, regarding total charge, total spin, boundary condition, atomic types or number of atoms requires the complete reconstruction of the object.
- Raises:
TBLiteValueError – on invalid input, like incorrect shape / type of the passed arrays
Calculator#
- class tblite.interface.Calculator(method: str, numbers: numpy.ndarray | List[int], positions: numpy.ndarray, charge: float | None = None, uhf: int | None = None, lattice: numpy.ndarray | None = None, periodic: numpy.ndarray | None = None, **context_kwargs)[source]#
Represents a wrapped calculator object in
tblite
and the associated structure data. The calculator is instantiated for the respective structure and is immutable once created. The cartesian coordinates and the lattice parameters of the structure can be updated, while changing the boundary conditions, atomic types or number of atoms require the complete reconstruction of the calculator instance.Available methods to parametrization of the calculator are:
GFN2-xTB:
Self-consistent extended tight binding Hamiltonian with anisotropic second order electrostatic contributions, third order on-site contributions and self-consistent D4 dispersion.
Geometry, frequency and non-covalent interactions parametrisation for elements up to Z=86.
Cite as:
C. Bannwarth, S. Ehlert and S. Grimme., J. Chem. Theory Comput. (2019), 15, 1652-1671. DOI: 10.1021/acs.jctc.8b01176
GFN1-xTB:
Self-consistent extended tight binding Hamiltonian with isotropic second order electrostatic contributions and third order on-site contributions.
Geometry, frequency and non-covalent interactions parametrisation for elements up to Z=86.
Cite as:
S. Grimme, C. Bannwarth, P. Shushkov, J. Chem. Theory Comput. (2017), 13, 1989-2009. DOI: 10.1021/acs.jctc.7b00118
IPEA1-xTB:
Special parametrisation for the GFN1-xTB Hamiltonian to improve the description of vertical ionisation potentials and electron affinities. Uses additional diffuse s-functions on light main group elements. Parametrised up to Z=86.
Cite as:
V. Asgeirsson, C. Bauer and S. Grimme, Chem. Sci. (2017), 8, 4879. DOI: 10.1039/c7sc00601b
Example
>>> from tblite.interface import Calculator >>> import numpy as np >>> numbers = np.array([1, 1, 6, 5, 1, 15, 8, 17, 13, 15, 5, 1, 9, 15, 1, 15]) >>> positions = np.array([ # Coordinates in Bohr ... [+2.79274810283778, +3.82998228828316, -2.79287054959216], ... [-1.43447454186833, +0.43418729987882, +5.53854345129809], ... [-3.26268343665218, -2.50644032426151, -1.56631149351046], ... [+2.14548759959147, -0.88798018953965, -2.24592534506187], ... [-4.30233097423181, -3.93631518670031, -0.48930754109119], ... [+0.06107643564880, -3.82467931731366, -2.22333344469482], ... [+0.41168550401858, +0.58105573172764, +5.56854609916143], ... [+4.41363836635653, +3.92515871809283, +2.57961724984000], ... [+1.33707758998700, +1.40194471661647, +1.97530004949523], ... [+3.08342709834868, +1.72520024666801, -4.42666116106828], ... [-3.02346932078505, +0.04438199934191, -0.27636197425010], ... [+1.11508390868455, -0.97617412809198, +6.25462847718180], ... [+0.61938955433011, +2.17903547389232, -6.21279842416963], ... [-2.67491681346835, +3.00175899761859, +1.05038813614845], ... [-4.13181080289514, -2.34226739863660, -3.44356159392859], ... [+2.85007173009739, -2.64884892757600, +0.71010806424206], ... ]) >>> calc = Calculator("GFN2-xTB", numbers, positions) >>> res = calc.singlepoint() >>> res.get("energy") # Results in atomic units -31.716159156026254
- Raises:
TBLiteValueError – on invalid input, like incorrect shape / type of the passed arrays
- add(interaction, *args) None [source]#
Add an interaction to the calculator instance. Supported interactions are
name
description
Arguments
electric-field
Uniform electric field
Field vector (3,)
spin-polarization
Spin polarization
Scaling factor
alpb-solvation
ALPB implicit solvation
Solvent name, solution state (optional)
gbsa-solvation
GBSA implicit solvation
Solvent name, solution state (optional)
cpcm-solvation
CPCM implicit solvation
Epsilon
gbe-solvation
GBε implicit solvation
Epsilon, Born kernel
gb-solvation
GB implicit solvation
Epsilon, Born kernel
Note
For GSBA and ALPB: For named solvents, uses parametrized GBSA/ALPB with CDS and empirical shift. For unnamed solvents (dielectric constant), uses non-empirical GBSA/ALPB. Optional solution state correction: gsolv (default), bar1mol, reference.
- get(attribute: str) Any [source]#
Set an attribute from the calculator instance. Supported attributes are
name
description
angular-momenta
Angular momenta of shells
orbital-map
Mapping from orbitals to shells
shell-map
Mapping from shells to atoms
- Raises:
TBLiteValueError – on invalid attributes
- set(attribute: str, value) None [source]#
Set an attribute in the calculator instance. Supported attributes are
name
description
default
accuracy
Numerical thresholds for SCC
1.0
guess
Initial guess for wavefunction
0 (SAD)
max-iter
Maximum number of SCC iterations
250
mixer-damping
Parameter for the SCC mixer
0.4
save-integrals
Keep integral matrices in results
0 (False)
temperature
Electronic temperature for filling
9.500e-4
verbosity
Set verbosity of printout
1
Note
The electronic temperature is given in Hartree, rather than Kelvin. The conversion factor from Kelvin to Hartree is the Boltzmann constant in Hartree/Kelvin (3.166808578545117e-6).
- Raises:
TBLiteValueError – on invalid input, like incorrect shape / type of the passed arrays
- singlepoint(res: Result | None = None, copy: bool = False) Result [source]#
Perform actual single point calculation in the library backend. The output of the library will be forwarded to the standard output.
The routine returns an object containing the results, which can be used to restart the calculation. Unless specified the restart object will be updated inplace rather than copied.
- Raises:
TBLiteRuntimeError – in case the calculation fails
Result#
- class tblite.interface.Result(other=None)[source]#
Container for calculation results, can be passed to a single point calculation as restart data. Allows to retrieve individual results or export all results as dict.
Example
>>> from tblite.interface import Calculator >>> import numpy as np >>> calc = Calculator( ... method="GFN2-xTB", ... numbers=np.array([14, 1, 1, 1, 1]), ... positions=np.array([ ... [ 0.000000000000, 0.000000000000, 0.000000000000], ... [ 1.617683897558, 1.617683897558,-1.617683897558], ... [-1.617683897558,-1.617683897558,-1.617683897558], ... [ 1.617683897558,-1.617683897558, 1.617683897558], ... [-1.617683897558, 1.617683897558, 1.617683897558], ... ]), ... ) >>> res = calc.singlepoint() ------------------------------------------------------------ cycle total energy energy error density error ------------------------------------------------------------ 1 -3.710873811182 -3.7424104E+00 1.5535536E-01 2 -3.763115011046 -5.2241200E-02 4.8803250E-02 3 -3.763205560205 -9.0549159E-05 2.0456319E-02 4 -3.763213200846 -7.6406409E-06 2.7461272E-03 5 -3.763250843634 -3.7642788E-05 2.4071582E-03 6 -3.763236758520 1.4085114E-05 2.3635321E-04 7 -3.763237287151 -5.2863152E-07 3.5812281E-04 8 -3.763233829618 3.4575334E-06 8.5040576E-06 9 -3.763233750684 7.8934355E-08 1.1095285E-07 ------------------------------------------------------------ >>> res.get("energy") -3.7632337506836944 >>> res = calc.singlepoint(res) ------------------------------------------------------------ cycle total energy energy error density error ------------------------------------------------------------ 1 -3.763233752583 -3.7947704E+00 1.3823013E-07 2 -3.763233751557 1.0256169E-09 1.4215249E-08 ------------------------------------------------------------
- Raises:
TBliteValueError – on invalid input, like incorrect shape / type of the passed arrays
TBLiteRuntimeError – if a requested quantity is not available in the container
- dict() dict [source]#
Return all quantities inside the result container as dict. In case no results are present an empty dict is returned.
- get(attribute: str)[source]#
Get a quantity stored instade the result container. The following quantities are available
property
dimension [spin-polarized case]
unit
energy
scalar
Hartree
energies
nat
Hartree
gradient
nat, 3
Hartree/Bohr
virial
3, 3
Hartree
charges
nat
e
bond-orders
nat, nat
e
dipole
3
e·Bohr
quadrupole
6
e·Bohr²
orbital-energies
norb [2, norb]
Hartree
orbital-occupations
norb [2, norb]
e
orbital-coefficients
norb, norb [2, norb, norb]
unitless
overlap-matrix
norb, norb
unitless
hamiltonian-matrix
norb, norb
Hartree
density-matrix
norb, norb [2, norb, norb]
e
natoms
scalar
unitless
norbitals
scalar
unitless
Notes
The Hamiltonian matrix is the core Hamiltonian rather than the converged full Hamiltonian after selfconsistency. To reconstruct it transform the orbital energies from the MO to the AO basis using the orbital coefficients.
- Raises:
TBLiteValueError – on invalid input, like incorrect shape / type of the passed arrays
TBLiteRuntimeError – if a requested quantity is not available in the container
- set(attribute: str, value)[source]#
Get a quantity stored instade the result container. Currently, no quantities can be set in the result container.
- Raises:
TBLiteValueError – on invalid input, like incorrect shape / type of the passed arrays
TBLiteRuntimeError – if a requested quantity cannot be set in the container