Function Reference¶
The following section contains a detailed description of all modules present in ORBKIT.
orbkit.qcinfo¶
Module for processing the data read from the output files of quantum chemical software.
-
class
orbkit.qcinfo.
CIinfo
(method='ci')[source]¶ Class managing all information from the from the output files of quantum chemical software for CI calculations.
The CI related features are in ongoing development.
-
class
orbkit.qcinfo.
QCinfo
(data=None)[source]¶ Class managing all information from the from the output files of quantum chemical software.
See Central Variables in the manual for details.
-
atoms
(bbox=None, **kwargs)¶ Create an ASE atoms object. (cf. https://wiki.fysik.dtu.dk/ase/ase/atoms.html )
Parameters:
- bbox : list of floats (bbox=[xmin,xmax,ymin,ymax,zmin,zmax]), optional
- If not None, sets the unit cell to the grid boundaries and moves the molecule in its center.
Returns:
- atoms : Atoms object
- See https://wiki.fysik.dtu.dk/ase/ase/atoms.html for details
Note
ASE has to be in the PYTHONPATH
-
format_geo
(is_angstrom=False)[source]¶ Converts geo_info and geo_spec to a universal format. Parameters:
- is_angstrom : bool, optional
- If True, input is assumed to be in Angstrom and positions are converted to Bohr radii.
-
get_ase_atoms
(bbox=None, **kwargs)[source]¶ Create an ASE atoms object. (cf. https://wiki.fysik.dtu.dk/ase/ase/atoms.html )
Parameters:
- bbox : list of floats (bbox=[xmin,xmax,ymin,ymax,zmin,zmax]), optional
- If not None, sets the unit cell to the grid boundaries and moves the molecule in its center.
Returns:
- atoms : Atoms object
- See https://wiki.fysik.dtu.dk/ase/ase/atoms.html for details
Note
ASE has to be in the PYTHONPATH
-
nuclear_repulsion
¶ Calculates nuclear repulsion energy.
-
select_spin
(restricted, spin=None)[source]¶ For an unrestricted calculation, the name of the MO (‘sym’ keyword in
qc.mo_spec
) is modified, e.g., 3.1_b for MO 3.1 with beta spin and 3.1_a for MO 3.1 with alpha spin. For restricted calculation, the ‘spin’ keyword fromqc.mo_spec
is removed.Parameters:
- restricted : bool
- If True, removes the ‘spin’ keyword from
qc.mo_spec
. - spin : {None, ‘alpha’, or ‘beta’}, optional
- If not None, returns exclusively ‘alpha’ or ‘beta’ molecular orbitals.
-
view
(select=slice(None, None, None), bbox=None, **kwargs)[source]¶ Opens ase-gui with the atoms of the QCinfo class. (cf. https://wiki.fysik.dtu.dk/ase/ase/visualize/visualize.html )
Parameters:
- select : slice or (array of int), default: all atoms
- Specifies the atoms to be shown.
- bbox : list of floats (bbox=[xmin,xmax,ymin,ymax,zmin,zmax]), optional
- If not None, sets the unit cell to the grid boundaries and moves the molecule in its center.
Note
ASE has to be in the PYTHONPATH
-
orbkit.orbitals¶
-
class
orbkit.orbitals.
AOClass
(seq=[], restart=None)[source] AO base class which contains all information on atomic orbitals. Two types of dataformats are available:
Numpy-style data:
- cont2atoms : numpy.ndarray, dtype=numpy.intc, shape = (NAO)
Transformation matrix between contracted GTO’s and atoms.
- _assign_prim_to_cont : numpy.ndarray, dtype=numpy.intc, shape = (NPAO)
Transformation matrix between contracted GTO’s and primitive GTO’s.
- pg_expcont : numpy.ndarray, dtype=float64, shape = (NPAO, 2)
- Information on primitive GTO’s:
1st element exponent and 2nd element contraction
- contspher : numpy.ndarray, dtype=numpy.intc, shape = (NAO, 2)
Same information as ao_spherical as numpy.ndarray
- lxlylz : numpy.ndarray, dtype=numpy.intc, shape = (NAO, 3)
Contains the expontents lx, ly, lz for the Cartesian Gaussians.
Lists of dictionaries / list of tuples:
Member of dict Content ——————– ————————————- ‘atom’ Index of atom ‘pnum’ Number of primitives ‘type’ Type of AO ‘coeffs’ AO coefficients ‘lxlylz’ Exponents of Cartesian Gaussians ‘lm’ (if spherical) Quantum number of Spherical Gaussians
See Central Variables in the manual for details.
-
append
(item)[source] S.append(value) – append value to the end of the sequence
-
extend
(item)[source] S.extend(iterable) – extend sequence by appending elements from the iterable
-
get_old_ao_spherical
()[source] Compatability funtction to allow access to old-style ao_spherical.
-
internal_to_dict
()[source] Transforms Numpy-style data to lists of dictionary style data for compatability.
-
is_normlized
(force=False)[source] Check if orbitals in AOClass are normalized.
-
remove
(item)[source] S.remove(value) – remove first occurrence of value. Raise ValueError if the value is not present.
-
set_lm_dict
(p=[1, 0])[source] Sets the l,m quantum numbers for the contracted spherical harmonic Gaussian basis set.
-
update
()[source] Transfers UserList data dictionary to internal variables
-
update_lxlylz
()[source] Extracts the exponents lx, ly, lz for the Cartesian Gaussians.
Parameters:
- get_assign : bool, optional
- Specifies, if the index of the atomic orbital shall be returned as well.
Returns:
- _lxlylz : numpy.ndarray, dtype=numpy.intc, shape = (NAO,3)
- Contains the expontents lx, ly, lz for the Cartesian Gaussians.
- _assign_lxlylz_to_cont : list of int, optional
- Contains the index of the atomic orbital in ao_spec.
-
class
orbkit.orbitals.
MOClass
(seq=[], restart=None)[source] MO base class which contains all information on atomic orbitals. Two types of dataformats are available:
Numpy-style data:
- coeffs : numpy.ndarray, dtype=float64, shape = (NMO, NAO)
Molecular orbital coefficients.
- occ : numpy.ndarray, dtype=float64, shape = (NMO)
Occupation numbers for molecular orbitals.
- eig : numpy.ndarray, dtype=float64, shape = (NMO)
Eigenvalues for molecular orbitals.
- sym : numpy.ndarray, dtype=str, shape = (NMO)
MOLPRO-like symmetry label of molecular orbitals.
Lists of dictionaries:
mo_spec
See Central Variables in the manual for details.
-
append
(item)[source] S.append(value) – append value to the end of the sequence
-
combine_restricted_orbitals
()[source] Combines unrestricted orbitals with same MO coefficients to restricted orbitals. (For example Psi4 creates a molden file with unrestricted orbitals after a ROHF calculation.)
-
extend
(item)[source] S.extend(iterable) – extend sequence by appending elements from the iterable
-
get_coeffs
()[source] Get function for numpy array version of molecular orbital coefficients.
Returns:
coeffs: numpy.ndarray, dtype=numpy.float64, shape = (NMO, NAO)
-
get_eig
()[source] Get function for numpy array version of molecular orbital energies.
Returns:
eig: numpy.ndarray, dtype=numpy.float64, shape = (NMO)
-
get_homo
(tol=1e-05, sort=True)[source] Returns index of highest occupied MO.
-
get_lastbound
(sort=True)[source] Returns index of highest bound MO.
-
get_lumo
(tol=1e-05, sort=True)[source] Returns index of lowest unoccupied MO.
-
get_occ
(return_alpha_beta=False, return_int=False, tol_int=1e-05, sum_occ=False)[source] Get function for numpy array version of molecular orbital occupancies.
Parameters:
return_alpha_beta: bool, optional, determines whether occupancies should be returned as a shape = (NMO) array (False) or a shape = (2,NMO/2) array (True) return_int: bool, optional, determines whether occupancies should be returned as a numpy.intc array. Raises an error if there are partially occupied orbitals. tol_int: float, optional, defines the tolerance for when an occupation number is considered an integer sum_occ: bool, optional, determines whether occupancies should be summed overNote:
return_alpha_beta=True is only supported for spin-polarized calculations.Returns:
occ: numpy.ndarray, dtype=numpy.float64, shape = (NMO)
-
get_spin_index
(spin)[source] - Function used to select MO’s by spin. A numpy.ndarray is returned
- which contains the indexes of MO’s of the selected spin.
Parameters:
spin : ‘str’, can be either ‘alpha’ or ‘beta’
Returns:
alpha : numpy.ndarray, dtype=numpy.intc
-
get_spinstate
()[source] Determines whether the MOClass has alpha and beta spins and removes the _a/_b spin labels. For spin-paired calculations all spins are set to ‘alpha’.
-
get_sym
()[source] Get function for numpy array version of molecular orbital symmetries.
Returns:
sym: numpy.ndarray, dtype=numpy.str, shape = (NMO)
-
remove
(item)[source] S.remove(value) – remove first occurrence of value. Raise ValueError if the value is not present.
-
select
(fid_mo_list, flatten_input=True, sort_indices=True)[source] - Selects molecular orbitals from an external file or a list of molecular
- orbital labels.
Parameters:
- mo_spec :
- See Central Variables for details.
- fid_mo_list : str, ‘all_mo’, or list
- If fid_mo_list is a str, specifies the filename of the molecular orbitals list.If fid_mo_list is ‘all_mo’, creates a list containing all molecular orbitals.If fid_mo_list is a list, provides a list (or a list of lists) of molecular orbital labels.
- flatten_input : boolean, optional
- Specifies wheter lists of lists should be flattened so a single MOClass instance can be returned rather than a list of MOClass instances.
- sort_indices : boolean, optional
- Specifies wheter list of indexes should be sorted before it is returned. This is only supported if flatten_input is set to
True
.
Supported Formats:
Integer List (Counting from Zero!):
1 2 3 5 4 homo lumo+2:lumo+4
List with Symmetry Labels:
1.1 2.1 1.3 1.1 4.1 4.1 2.3 2.1
alpha
andbeta
can be used together with symmetry labels to restrict the selection to orbitals of that symmetry. This option is not supported for integer lists. Note also thatalpha
andbeta
only restrict selection within one string of inputs. If you which to spin-restrict orbitlas given as a list of strings please useall_alpha
orall_beta
. You may also select all MO’s of a given symmetry by specifying only the symmetry label with a preceeding asterisc, e.g.*.A1
to get all orbitals of A1 symmetry.Returns:
List of MOClass instances containing the selected orbitals as well as further information on the selection criteria used If a sinlge list is used as in input and/or flatten_input=True, an MOClass instance is returned instead.
-
set_coeffs
(item)[source] Set function for numpy array version of molecular orbital coefficients.
Returns:
coeff: numpy.ndarray, dtype=numpy.float64, shape = (NMO, NAO)
-
set_eig
(item)[source] Set function for numpy array version of molecular orbital energies.
Parameters:
eig: numpy.ndarray, dtype=numpy.float64, shape = (NMO)
-
set_occ
(item)[source] Set function for numpy array version of molecular orbital occupancies.
Parameters:
occ: numpy.ndarray, dtype=numpy.float64, shape = (NMO)
-
set_sym
(item)[source] Set function for numpy array version of molecular orbital symmetries.
Parameters:
sym: numpy.ndarray, dtype=numpy.str, shape = (NMO)
-
sort_by_sym
()[source] Sorts mo_spec by symmetry.
orbkit.options¶
Module containing and processing all orbkit options.
Functions for developers¶
-
orbkit.options.
check_options
(error=<function raise_error>, display=<function print_message>, interactive=False, info=True, check_io=True)[source]¶ Checks options for errors.
Parameters:
- error : function, optional
- Handles the errors.
- display : function, optional
- Handles the print commands.
- interactive : bool, optional
- If True and a file does not exist, asks the user to insert name of existing file.
- info : bool, optional
- If True, some additional information is printed.
Default Error and Exception Handling: Prints the errors and continues.
-
orbkit.options.
check_if_exists
(fid, what='', error=<class 'OSError'>, display=<built-in method write of _io.TextIOWrapper object>, interactive=False)[source]¶ Checks the existence of a file.
Parameters:
- fid : string
- Specifies filename of the requested file.
- what : string, optional
- Describes the file.
- error : function, optional
- Handles the errors.
- display : function, optional
- Handles the print commands.
- interactive : bool, optional
- If True and a file does not exist, asks the user to insert name of existing file.
Returns:
- fid : string
- Specifies filename of the requested file.
Parameters¶
-
orbkit.options.
itypes
= ['auto', 'molden', 'aomix', 'gamess', 'gaussian.log', 'gaussian.fchk', 'wfn', 'wfx', 'cclib', 'native', 'orbkit.dump']¶ Specifies possible input types.
-
orbkit.options.
otypes
= ['h5', 'hdf5', 'npz', 'cube', 'cb', 'am', 'hx', 'vmd', 'mayavi', 'native', 'auto']¶ Specifies possible output types.
-
orbkit.options.
drv_options
= ['None', 'x', 'y', 'z', 'xx', 'yy', 'zz', 'x2', 'y2', 'z2', 'xy', 'yx', 'xz', 'zx', 'yz', 'zy']¶ Specifies possible derivative variables.
orbkit.main¶
Module for controlling all computational tasks.
-
orbkit.main.
init
(reset_display=True)[source]¶ Resets all
orbkit.options
andorbkit.display
.
-
orbkit.main.
run_orbkit
(use_qc=None, check_options=True, standalone=False)[source]¶ Controls the execution of all computational tasks.
Parameters:
- use_qc : QCinfo, optional
- If not None, the reading of a quantum chemistry output is omitted and the given QCinfo class is used for all computational tasks. (See Central Variables in the manual for details on QCinfo.)
- check_options : bool, optional
- If True, the specified options will be validated.
Returns:
- data : type and shape depend on the options.
- Contains orbkit’s output. See ORBKIT’s High-Level Interface in the manual for details.
orbkit.read¶
Build new ORBKIT io interface
orbkit.grid¶
Module for creating and manipulating the grid on which all computations are performed.
-
orbkit.grid.
N_
= [101, 101, 101]¶ Specifies the number of grid points (regular grid).
-
orbkit.grid.
adjust_to_geo
(qc, extend=5.0, step=0.1)[source]¶ Adjusts the grid boundaries to the molecular geometry.
Parameters:
- qc : QCinfo class
- See Central Variables for details.
- extend : float
- Specifies the value by which the grid boundaries are extended in each direction.
- step : float
- Specifies the grid spacing.
-
orbkit.grid.
center_grid
(ac, display=<built-in method write of _io.TextIOWrapper object>)[source]¶ Centers the grid to the point ac and to the origin (0,0,0).
-
orbkit.grid.
cyl2cart_vector
(r, phi, zed)[source]¶ Converts a cylindrical regular grid matrix (r, phi, zed) to a Cartesian grid matrix (vector grid) with the shape (3, (Nr*Nphi*Nzed)).
Parameters:
- r : numpy.ndarray, shape=(Nr,)
- Specifies radial distance.
- phi : numpy.ndarray, shape=(Nphi,)
- Specifies azimuth angle.
- zed : numpy.ndarray, shape=(Nz,)
- Specifies z distance.
-
orbkit.grid.
d3r
= 0.0¶ A volume element
-
orbkit.grid.
delta_
= array([[0.], [0.], [0.]])¶ Contains the grid spacing.
-
orbkit.grid.
grid2vector
()[source]¶ Converts the regular grid characterized by x-, y-, z-vectors to a (3, (Nx*Ny*Nz)) grid matrix (vector grid). Reverse operation:
orbkit.grid.vector2grid
-
orbkit.grid.
grid_init
(is_vector=False, force=False)[source]¶ Sets up the regular x-, y-, z-grid specified by the global lists:
min_: List of minimum grid values max_: List of maximum grid values N_: List of number of grid points Parameters:
- is_vector : bool, optional
- If True, converts the regular grid to a vector grid.
-
orbkit.grid.
init
(is_vector=False, force=False)¶ Sets up the regular x-, y-, z-grid specified by the global lists:
min_: List of minimum grid values max_: List of maximum grid values N_: List of number of grid points Parameters:
- is_vector : bool, optional
- If True, converts the regular grid to a vector grid.
-
orbkit.grid.
init_grid
(is_vector=False, force=False)¶ Sets up the regular x-, y-, z-grid specified by the global lists:
min_: List of minimum grid values max_: List of maximum grid values N_: List of number of grid points Parameters:
- is_vector : bool, optional
- If True, converts the regular grid to a vector grid.
-
orbkit.grid.
is_initialized
= False¶ If True, the grid is assumed to be initialized.
-
orbkit.grid.
is_regular
= False¶ If True, the grid is assumed to be regular, i.e., a conversion of a vector grid to a regular grid is possible, if
N_
is set.
-
orbkit.grid.
is_vector
= True¶ If True, the grid is assumed to be vector grid.
-
orbkit.grid.
matrix_grid2vector
(matrix)[source]¶ Converts the (Nx,Ny,Nz) data matrix back to the regular grid (Nx,Nz,Ny)
-
orbkit.grid.
matrix_vector2grid
(matrix, Nx=None, Ny=None, Nz=None)[source]¶ Converts the (Nx*Ny*Nz) data matrix back to the (Nx,Nz,Ny)
-
orbkit.grid.
max_
= [8.0, 8.0, 8.0]¶ Specifies maximum grid values (regular grid).
-
orbkit.grid.
min_
= [-8.0, -8.0, -8.0]¶ Specifies minimum grid values (regular grid).
-
orbkit.grid.
mv2g
(**kwargs)[source]¶ Converts all numpy.ndarrays given as the keyword arguments (**kwargs) from a vector grid of shape=(…, Nx*Ny*Nz, …,) to a regular grid of shape=(…, Nx, Ny, Nz, …,), and, if more than one **kwargs is given, returns it as a dictionary.
Hint: The global values for the grid dimensionality, i.e.,
grid.N_
, are used for reshaping.
-
orbkit.grid.
random_grid
(geo_spec, N=1000000.0, scale=0.5)[source]¶ Creates a normally distributed grid around the atom postions (geo_spec).
Parameters:
- geo_spec :
- See Central Variables for details.
- N : int
- Number of points distributed around each atom
- scale : float
- Width of normal distribution
-
orbkit.grid.
read
(filename, comment='#')[source]¶ Reads a grid from a plain text file.
Parameters:
- fid : str
- Specifies the filename of the grid file.
Returns:
- is_vector : bool
- If True, a vector grid is used for the computations.
Supported Formats:
Regular Grid:
# Regular Grid Example File # Format: # x xmin xmax Nx # y ymin ymax Ny # z zmin zmax Nz x -5 5 101 y -2 2 51 z 0 0 1
Vector Grid:
# Vectorized Grid Example File # The header 'x y z' is mandatory! x y z -5 -5 0 -4 -5 0 -3 -5 0 0 0 0 2 -1e-1 9.78
Hint: If a line starts with ‘#’, it will be skipped. Please, do not use ‘#’ at the end of a line!
-
orbkit.grid.
reflect
(plane)[source]¶ Creates matrix representation for reflection Plane has to be specified as follows: xy-plane -> plane= numpy.array([0,1]) xz-plane -> plane= numpy.array([0,2]) yz-plane -> plane= numpy.array([1,2])
-
orbkit.grid.
rot
(ang, axis)[source]¶ Creates matrix representation for arbitrary rotations Angle has to be defined in radians, e.g., numpy.pi/2.0 Axis has to be specified as follows: x-axis -> axis=0, y-axis -> axis=1, z-axis -> axis=2,
-
orbkit.grid.
sph2cart_vector
(r, theta, phi)[source]¶ Converts a spherical regular grid matrix (r, theta, phi) to a Cartesian grid matrix (vector grid) with the shape (3, (Nr*Ntheta*Nphi)).
Parameters:
- r : numpy.ndarray, shape=(Nr,)
- Specifies radial distance.
- theta : numpy.ndarray, shape=(Ntheta,)
- Specifies polar angle.
- phi : numpy.ndarray, shape=(Nphi,)
- Specifies azimuth angle.
-
orbkit.grid.
vector2grid
(Nx, Ny, Nz)[source]¶ Converts the (3, (Nx*Ny*Nz)) grid matrix (vector grid) back to the regular grid characterized by the x-, y-, z-vectors. Reverse operation:
orbkit.grid.grid2vector
-
orbkit.grid.
x
= array([0.])¶ Contains the x-coordinates.
-
orbkit.grid.
y
= array([0.])¶ Contains the y-coordinates.
-
orbkit.grid.
z
= array([0.])¶ Contains the z-coordinates.
orbkit.core¶
Module performing all computational tasks.
-
orbkit.core.
ao_creator
(geo_spec, ao_spec, drv=None, x=None, y=None, z=None, is_vector=None)[source]¶ Calculates all contracted atomic orbitals or its derivatives with respect to a specific variable (e.g. drv = ‘x’ or drv = 0).
Parameters:
- geo_spec,ao_spec :
- See Central Variables in the manual for details.
- sel_ao : int
- Index of the requested atomic orbital
- drv : int or string, {None, ‘x’, ‘y’, ‘z’, 0, 1, 2}, optional
- If not None, an analytical calculation of the derivatives for the atomic orbitals with respect to DRV is requested.
- x,y,z : None or list of floats, optional
- If not None, provides a list of Cartesian coordinates, else the respective coordinates of grid. will be used
- is_vector : bool, optional
- If True, a vector grid will be applied
Returns:
- ao_list : numpy.ndarray, shape=((NAO,) + N)
- Contains the computed NAO atomic orbitals on a grid.
-
orbkit.core.
calc_mo_matrix
(qc_a, qc_b=None, drv=None, numproc=1, slice_length=10000.0, save_hdf5=False, **kwargs)[source]¶ Calculates and saves the selected molecular orbitals or the derivatives thereof. Parameters:
- qc.geo_spec, qc.geo_info, qc.ao_spec, qc.mo_spec :
- See Central Variables for details.
- fid_mo_list : str
- Specifies the filename of the molecular orbitals list or list of molecular
orbital labels (cf.
orbkit.read.mo_select
for details). If fid_mo_list is ‘all_mo’, creates a list containing all molecular orbitals. - drv : string or list of strings {None,’x’,’y’, ‘z’, ‘xx’, ‘xy’, …}, optional
- If not None, a derivative calculation of the molecular orbitals is requested.
- otype : str or list of str, optional
- Specifies output file type. See
otypes
for details. - ofid : str, optional
- Specifies output file name. If None, the filename will be based on
orbkit.options.outputname
. - numproc : int, optional
- Specifies number of subprocesses for multiprocessing.
If None, uses the value from
options.numproc
. - slice_length : int, optional
- Specifies the number of points per subprocess.
If None, uses the value from
options.slice_length
. - full_matrix : bool
- if False, numpy.tril_indices(NMO)
- Returns:
- mo_list : numpy.ndarray, shape=((NMO,) + N)
- Contains the NMO=len(qc.mo_spec) molecular orbitals on a grid.
-
orbkit.core.
cartesian2spherical
(ao_list, ao_spec)[source]¶ Transforms the atomic orbitals from a Cartesian Gaussian basis to a (real) pure spherical harmonic Gaussian basis set.
Adapted from H.B. Schlegel and M.J. Frisch, International Journal of Quantum Chemistry, Vol. 54, 83-87 (1995).
Parameters:
- ao_list : numpy.ndarray, shape=((NAO,) + N)
- Contains the NAO atomic orbitals on a grid.
Returns:
- ao_list_sph : numpy.ndarray, shape=((NAO,) + N)
- Contains the NAO spherical atomic orbitals on a grid.
..hint:
The conversion is currently only supported up to g atomic orbitals.
-
orbkit.core.
mo_creator
(ao_list, mo_spec)[source]¶ Calculates the molecular orbitals.
Parameters:
- ao_list : numpy.ndarray, shape=((NAO,) + N)
- Contains the NAO atomic orbitals on a grid.
- mo_spec : List of dictionaries
- See Central Variables for details.
- mo_coeff : numpy.ndarray, shape = (NMO,NAO)
- Contains the molecular orbital coefficients of all orbitals.
Returns:
- mo_list : numpy.ndarray, shape=((NMO,) + N)
- Contains the NMO=len(mo_spec) molecular orbitals on a grid.
-
orbkit.core.
rho_compute
(qc, calc_ao=False, calc_mo=False, drv=None, laplacian=False, numproc=1, slice_length=10000.0, vector=None, save_hdf5=False, **kwargs)[source]¶ Calculate the density, the molecular orbitals, or the derivatives thereof.
orbkit divides 3-dimensional regular grids into 2-dimensional slices and 1-dimensional vector grids into 1-dimensional slices of equal length. By default, 3-dimensional grids are used (
vector=None
). The computational tasks are distributed to the worker processes.Parameters:
- qc : class or dict
- QCinfo class or dictionary containing the following attributes/keys. See Central Variables for details.
- qc.geo_spec : numpy.ndarray, shape=(3,NATOMS)
- See Central Variables for details.
- qc.ao_spec : List of dictionaries
- See Central Variables for details.
- qc.mo_spec : List of dictionaries
- See Central Variables for details.
- calc_ao : bool, optional
- If True, the computation of the atomic orbitals is only carried out.
- calc_mo : bool, optional
- If True, the computation of the molecular orbitals requested is only carried out.
- slice_length : int, optional
- Specifies the number of points per subprocess.
- drv : string or list of strings {None,’x’,’y’, ‘z’, ‘xx’, ‘xy’, …}, optional
- If not None, computes the analytical derivative of the requested quantities with respect to DRV.
- laplacian : bool, optional
- If True, computes the laplacian of the density.
- numproc : int
- Specifies number of subprocesses for multiprocessing.
- grid : module or class, global
- Contains the grid, i.e., grid.x, grid.y, and grid.z. If grid.is_initialized is not True, functions runs grid.grid_init().
Returns:
If calc_mo and drv is None: - mo_list
If calc_mo and drv is not None: - delta_mo_list
If not calc_mo and drv is None: - rho
If not calc_mo and drv is not None: - rho, delta_rho
If not calc_mo and laplacian: - rho, delta_rho, laplacian_rho
- mo_list : numpy.ndarray, shape=((NMO,) + N)
- Contains the NMO=len(qc.mo_spec) molecular orbitals on a grid.
- delta_mo_list : numpy.ndarray, shape=((NDRV,NMO) + N)
- Contains the derivatives with respect to drv (NDRV=len(drv)) of the NMO=len(qc.mo_spec) molecular orbitals on a grid.
- mo_norm : numpy.ndarray, shape=(NMO,)
- Contains the numerical norms of the molecular orbitals.
- rho : numpy.ndarray, shape=(N)
- Contains the density on a grid.
- delta_rho : numpy.ndarray, shape=((NDRV,) + N)
- Contains derivatives with respect to drv (NDRV=len(drv)) of the density on a grid.
- laplacian_rho : numpy.ndarray, shape=(N)
- Contains the laplacian of the density on a grid, i.e. \(\nabla^2 \rho = \nabla^2_x \rho + \nabla^2_y \rho + \nabla^2_z \rho\).
-
orbkit.core.
rho_compute_no_slice
(qc, calc_ao=False, calc_mo=False, drv=None, laplacian=False, return_components=False, x=None, y=None, z=None, is_vector=None, **kwargs)[source]¶ Calculates the density, the molecular orbitals, or the derivatives thereof without slicing the grid.
Parameters:
- qc : class or dict
- QCinfo class or dictionary containing the following attributes/keys. See Central Variables for details.
- qc.geo_spec : numpy.ndarray, shape=(3,NATOMS)
- See Central Variables for details.
- qc.ao_spec : List of dictionaries
- See Central Variables for details.
- qc.mo_spec : List of dictionaries
- See Central Variables for details.
- calc_ao : bool, optional
- If True, the computation of the atomic orbitals is only carried out.
- calc_mo : bool, optional
- If True, the computation of the molecular orbitals requested is only carried out.
- is_vector : bool, optional
- If True, performs the computations for a vector grid, i.e., with x, y, and z as vectors.
- drv : string or list of strings {None,’x’,’y’, ‘z’, ‘xx’, ‘xy’, …}, optional
- If not None, computes the analytical derivative of the requested quantities with respect to DRV.
- laplacian : bool, optional
- If True, computes the laplacian of the density.
- return_components : bool, optional
- If True, returns the atomic and molecular orbitals, and the density, and if requested, the derivatives thereof as well.
- x,y,z : numpy.ndarray, optional
- If not None, provides a list of Cartesian coordinates,
else the respective coordinates of the module
orbkit.grid
will be used.
Returns:
If not return_components: if calc_mo and drv is None: - mo_list
if calc_mo and drv is not None: - delta_mo_list
if not calc_mo and drv is None: - rho
if not calc_mo and drv is not None: - rho, delta_rho
if not calc_mo and laplacian: - rho, delta_rho, laplacian_rho
Else: if calc_mo and drv is None: - ao_list,mo_list
if calc_mo and drv is not None: - delta_ao_list,delta_mo_list
if not calc_mo and drv is None: - ao_list,mo_list,rho
if not calc_mo and drv is not None: - ao_list, mo_list, rho, delta_ao_list, delta_mo_list, delta_rho
if not calc_mo and laplacian: - ao_list, mo_list, rho, delta_ao_list, delta_mo_list, delta_rho, laplacian_rho
- ao_list : numpy.ndarray, shape=((NAO,) + N)
- Contains the NAO=len(ao_spec) atomic orbitals on a grid.
- delta_ao_list : numpy.ndarray, shape=((NDRV,NAO) + N)
- Contains the derivatives with respect to drv (NDRV=len(drv)) of the NAO=len(ao_spec) atomic orbitals on a grid.
- mo_list : numpy.ndarray, shape=((NMO,) + N)
- Contains the NMO=len(qc.mo_spec) molecular orbitals on a grid.
- delta_mo_list : numpy.ndarray, shape=((NDRV,NMO) + N)
- Contains the derivatives with respect to drv (NDRV=len(drv)) of the NMO=len(qc.mo_spec) molecular orbitals on a grid.
- mo_norm : numpy.ndarray, shape=(NMO,)
- Contains the numerical norms of the molecular orbitals.
- rho : numpy.ndarray, shape=(N)
- Contains the density on a grid.
- delta_rho : numpy.ndarray, shape=((NDRV,) + N)
- Contains the derivatives with respect to drv (NDRV=len(drv)) of the density on a grid.
- laplacian_rho : numpy.ndarray, shape=(N)
- Contains the laplacian of the density on a grid, i.e. \(\nabla^2 \rho = \nabla^2_x \rho + \nabla^2_y \rho + \nabla^2_z \rho\).
-
orbkit.core.
slice_rho
(xx)[source]¶ Calculates the density, the molecular orbitals, or the derivatives thereof with respect to Spec[‘Derivative’] for one slice (xx)
This function is called by the multiprocessing module in the
orbkit.core.rho_compute
.Parameters:
- xx : [float] or [int, int]
Specifies which slice in x-direction shall be computed.
If not is_vector: One slice at x=xx will be computed.Else: One slice from index xx[0] to xx[1] will be calculated.- Spec : dict, global
- Dictionary containing all required varibles:
geo_spec: List of floats, shape=(NATOMS, 3) (see Central Variables for details). ao_spec: List of dictionaries (see Central Variables for details). mo_spec: List of dictionaries (see Central Variables for details). calc_mo: Bool if only the molecular orbitals are requested. is_vector: Bool if a vector grid is used. Derivative: List of strings, choices={‘x’,’y’, or ‘z’}. If not None, derivative calculation will be carried out.
- grid : module or class, global
- Contains the grid, i.e., grid.x, grid.y, and grid.z.
Returns:
If calc_mo and drv is None: - mo_list
If calc_mo and drv is not None: - delta_mo_list
If not calc_mo and drv is None: - rho, mo_norm
If not calc_mo and drv is not None: - rho, mo_norm, delta_rho
- mo_list : numpy.ndarray, shape=((NMO,) + N)
- Contains the NMO=len(mo_spec) molecular orbitals on a grid.
- delta_mo_list : numpy.ndarray, shape=((NDRV,NMO) + N)
- Contains the derivatives with respect to drv (NDRV=len(drv)) of the NMO=len(mo_spec) molecular orbitals on a grid.
- rho : numpy.ndarray, shape=(N)
- Contains the density on a grid.
- delta_rho : numpy.ndarray, shape=((NDRV,) + N)
- Contains the derivatives with respect to drv (NDRV=len(drv)) of the density on a grid.
orbkit.analytical_integrals¶
Module performing analytical integrals between atomic and molecular orbitals.
Code for the computation of the overlap between primitive atomic basis functions adapted from
- Hô, J. M. Hernandez-Perez: “Evaluation of Gaussian Molecular Integrals”, DOI:10.3888/tmj.14-3
-
orbkit.analytical_integrals.
cartesian2spherical_aoom
(ao_overlap_matrix, ao_spec)[source]¶ Transforms the atomic orbitals overlap matrix from a Cartesian Gaussian basis set to a (real) pure spherical harmonic Gaussian basis set.
Adapted from H.B. Schlegel and M.J. Frisch, International Journal of Quantum Chemistry, Vol. 54, 83-87 (1995).
Parameters:
- ao_overlap_matrix : numpy.ndarray, shape = (NAO,NAO)
- Contains the overlap matrix of the Cartesian basis set.
- ao_spec,ao_spherical :
- See Central Variables in the manual for details.
Returns:
- aoom_sph : numpy.ndarray, shape = (NAO,NAO)
- Contains the overlap matrix of the spherical harmonic Gaussian basis set.
..hint:
Only supported up to g atomic orbitals and only for contracted atomic orbitals.
-
orbkit.analytical_integrals.
get_ao_dipole_matrix
(qc, component='x')[source]¶ Computes the expectation value of the dipole moment operator between all atomic orbitals.
Parameters:
- qc : class
- QCinfo class. (See Central Variables for details.)
- component : int or string, {‘x’,’y’, ‘z’}
- Specifies the compontent of the dipole moment operator which shall be applied.
Returns:
- ao_dipole_matrix : numpy.ndarray, shape=(NAO,NAO)
- Contains the expectation value matrix.
-
orbkit.analytical_integrals.
get_ao_overlap
(coord_a, coord_b, ao_spec, lxlylz_b=None, drv=None)[source]¶ Computes the overlap matrix of a basis set, where Bra basis set corresponds to the geometry
coord_a
and Ket basis set corresponds to the geometrycoord_b
.In order to enable the computation of analytical expectation values, the exponents lx, ly, lz for the primitive Cartesian Gaussians of the Ket basis set can be set manually with
lxlylz_b
. Please note that for the normalization of the primitive Cartesian Gaussians the exponents fromao_spec
are used.Parameters:
- coord_a : geo_spec
- Specifies the geometry of the Bra basis set. See Central Variables in the manual for details.
- coord_b : geo_spec
- Specifies the geometry of the Ket basis set. See Central Variables in the manual for details.
- ao_spec :
- See Central Variables in the manual for details.
- lxlylz_b : numpy.ndarray, dtype=numpy.int64, shape = (NAO,3), optional
- Contains the expontents lx, ly, lz for the primitive Cartesian Gaussians of the Ket basis set.
Returns:
- ao_overlap_matrix : numpy.ndarray, shape = (NAO,NAO)
- Contains the overlap matrix.
-
orbkit.analytical_integrals.
get_atom2mo
(qc)[source]¶ Assigns atom indices to molecular orbital coefficients.
Parameters:
- qc.ao_spec :
- See Central Variables for details.
Returns:
- atom2mo : numpy.ndarray, shape = (NAO,)
- Contains indices of atoms assigned to the molecular orbital coefficients.
-
orbkit.analytical_integrals.
get_dipole_moment
(qc, component=['x', 'y', 'z'])[source]¶ Computes the dipole moment analytically.
Parameters:
- qc : class
- QCinfo class. (See Central Variables for details.)
- component : string or list of strings, {‘x’,’y’, or ‘z’}
- Specifies the compontent(s) of the dipole moment which shall be computed.
Returns:
- dipole_moment : 1D numpy.array, shape[0]=len(component)
- Contains the dipole moment.
-
orbkit.analytical_integrals.
get_lc
(atoms, atom2mo, strict=False)[source]¶ Returns indices of molecular orbital coefficients corresponding to the selected atoms.
Parameters:
- atoms : int or list of int
- Contains the indices of the selected atoms.
- atom2mo : numpy.ndarray, shape = (NAO,)
- Contains indices of atoms assigned to the molecular orbital coefficients. >> atom2mo = get_atom2mo(qc)
Returns:
- lc : numpy.ndarray, shape = (NAO_atom,)
- Contains the NAO_atom indices molecular orbital coefficients corresponding to the selected atoms.
Example:
>> atom2mo = get_atom2mo(qc) >> get_lc([0,3],atom2mo)
-
orbkit.analytical_integrals.
get_mo_overlap
(mo_a, mo_b, ao_overlap_matrix)[source]¶ Computes the overlap of two molecular orbitals.
Parameters:
- mo_a : numpy.ndarray with shape = (,NAO)
- Contains the molecular orbital coefficients of all Bra orbitals.
- mo_b : numpy.ndarray with shape = (,NAO)
- Contains the molecular orbital coefficients of all Ket orbitals.
- ao_overlap_matrix : numpy.ndarray, shape = (NAO,NAO)
- Contains the overlap matrix of the basis set.
Returns:
- mo_overlap : float
- Contains the overlap of the two input molecular orbitals.
-
orbkit.analytical_integrals.
get_mo_overlap_matrix
(mo_a, mo_b, ao_overlap_matrix, numproc=1)[source]¶ Computes the overlap of two sets of molecular orbitals.
Parameters:
- mo_a : numpy.ndarray with shape = (NMO,NAO), dict, or MOClass instance
- Contains the molecular orbital coefficients of all Bra orbitals.
- mo_b : numpy.ndarray with shape = (NMO,NAO), dict, or MOClass instance
- Contains the molecular orbital coefficients of all Ket orbitals.
- ao_overlap_matrix : numpy.ndarray, shape = (NAO,NAO)
- Contains the overlap matrix of the basis set.
- numproc : int
- Specifies number of subprocesses for multiprocessing.
Returns:
- mo_overlap_matrix : numpy.ndarray, shape = (NMO,NMO)
- Contains the overlap matrix between the two sets of input molecular orbitals.
-
orbkit.analytical_integrals.
get_moom_atoms
(atoms, qc, mo_a, mo_b, ao_overlap_matrix, numproc=1)[source]¶ Computes the molecular orbital overlap matrix for selected atoms.
Parameters:
- atoms : int or list of int
- Contains the indices of the selected atoms.
- mo_a : numpy.ndarray with shape = (NMO,NAO) or mo_spec (cf. Central Variables)
- Contains the molecular orbital coefficients of all Bra orbitals.
- mo_b : numpy.ndarray with shape = (NMO,NAO) or mo_spec (cf. Central Variables)
- Contains the molecular orbital coefficients of all Ket orbitals.
- ao_overlap_matrix : numpy.ndarray, shape = (NAO,NAO)
- Contains the overlap matrix of the basis set.
- numproc : int
- Specifies number of subprocesses for multiprocessing.
Returns:
- mo_overlap_matrix : numpy.ndarray, shape = (NMO,NMO)
- Contains the overlap matrix between the two sets of input molecular orbitals.
-
orbkit.analytical_integrals.
get_nuclear_dipole_moment
(qc, component='x')[source]¶ Computes the nuclear part of the dipole moment.
Parameters:
- qc : class
- QCinfo class. (See Central Variables for details.)
- component : string, {‘x’,’y’, ‘z’}
- Specifies the compontent of the dipole moment operator which shall be applied.
Returns:
- nuclear_dipole_moment : float
- Contains the nuclear dipole moment.
orbkit.multiple_files¶
Module for processing multiple output files of quantum chemical software.
Capabilities
- Read a list of files
- Order the molecular orbital coefficients by using analytical integrals, by extrapolation, or by hand (e.g. for the preparation of an interpolation)
- Save and read the information obtained to and from an HDF5-file
- Depiction of molecular orbitals
Example for the Execution:
# Create a list of input filenames
fid_list = []
for i in range(101):
fid_list.append('~/molden_files/%03d.molden' % i)
init_display(name = 'MO_ordering') # Specify a filename for the oklog file
# (We want a oklog file but we have no
# options.outputname here.)
# Start reading and ordering routine.
order_using_analytical_overlap(fid_list)
orbkit.display¶
Module for writing the .oklog files and printing the terminal output.
-
orbkit.display.
display
(string)[source]¶ Prints
string
to the terminal output and to the .oklog file.
-
orbkit.display.
init_display
(name=None)[source]¶ Sets the name of the .oklog file and removes the old .oklog file.
-
orbkit.display.
is_initiated
= False¶ If True, logfile is initialized.
-
orbkit.display.
log_fid
= None¶ Specifies the filename of the oklog file
orbkit.output¶
Build new ORBKIT io interface
-
orbkit.output.
main_output
(data, qc=None, outputname='new', otype='auto', gname='', drv=None, omit=[], datalabels='', geo_info=None, geo_spec=None, **kwargs)[source]¶ Creates the requested output.
Parameters:
- data : numpy.ndarray, shape=N, shape=((NDRV,) + N), shape=(n, (NDRV,) + N) or list of numpy.ndarrays
- Contains the output data. The shape (N) depends on the grid and the data, i.e., 3d for regular grid, 1d for vector grid.
- geo_info, geo_spec :
- See Central Variables for details.
- outputname : str or list of str
- Contains the base name of the output file.
- otype : str or list of str
- Contains the output file type. Possible options: ‘auto’, ‘h5’, ‘cb’, ‘am’, ‘hx’, ‘vmd’, ‘mayavi’
- drv : None, list of str or list of list of str optional
- If not None, a 4d(regular)/2d(vector) input data array will be expected with NDRV = len(drv).
- omit : list of str, optional
- If not empty, the output file types specified here are omitted.
Note:
All additional keyword arguments are forwarded to the output functions.
-
orbkit.output.
amira_creator
(data, filename)[source]¶ Creates a ZIBAmira mesh file. (plain text)
Parameters:
- data : numpy.ndarray, shape=N or (3,)+N
- Contains the output data.
- filename : str
- Contains the base name of the output file.
-
orbkit.output.
hx_network_creator
(rho, filename)[source]¶ Creates a ZIBAmira hx-network file including a colormap file (.cmap) adjusted to the density for the easy depiction of the density.
-
orbkit.output.
cube_creator
(data, filename, geo_info, geo_spec, comments='', labels=None, **kwargs)[source]¶ Creates a plain text Gaussian cube file.
Parameters:
- data : numpy.ndarray, shape=N
- Contains the output data.
- filename : str
- Contains the base name of the output file.
- geo_info, geo_spec :
- See Central Variables for details.
- comments : str, optional
- Specifies the second (comment) line of the cube file.
orbkit.extras¶
Module for all additional features of ORBKIT.
-
orbkit.extras.
atom2index
(atom, geo_info=None)[source]¶ Converts a list of atom numbers to indices of
geo_info
.Parameters:
- atom : int or list of int
- Specifies the indices of the atoms (counting from one).
- geo_info : optional
- See Central Variables for details.
Returns:
- atom : list of int
- Contains the indices of the atoms (counting from one).
- index : numpy.ndarray, dtype=int
- Contains the indices of the atoms as occuring in qc.geo_info (counting from zero).
-
orbkit.extras.
calc_ao
(qc, drv=None, otype=None, ofid=None, numproc=None, slice_length=None)[source]¶ Computes and saves all atomic orbital or a derivative thereof.
Parameters:
- qc.geo_spec, qc.geo_info, qc.ao_spec, qc.mo_spec :
- See Central Variables for details.
- drv : string or list of strings {None,’x’,’y’, ‘z’, ‘xx’, ‘xy’, …}, optional
- If not None, a derivative calculation of the atomic orbitals is requested.
- otype : str or list of str, optional
- Specifies output file type. See
otypes
for details. - ofid : str, optional
- Specifies output file name. If None, the filename will be based on
orbkit.options.outputname
. - numproc : int, optional
- Specifies number of subprocesses for multiprocessing.
If None, uses the value from
options.numproc
. - slice_length : int, optional
- Specifies the number of points per subprocess.
If None, uses the value from
options.slice_length
.
- Returns:
- ao_list : numpy.ndarray, shape=((NAO,) + N)
- Contains the computed NAO atomic orbitals on a grid.
-
orbkit.extras.
calc_jmo
(qc, ij, drv=['x', 'y', 'z'], numproc=1, otype=None, ofid='', **kwargs)[source]¶ Calculate one component (e.g. drv=’x’) of the transition electoronic flux density between the molecular orbitals i and j.
\[moTEFD_{i,j}(r) = <mo_i|\delta(r-r')\nabla_x|mo_j>\]Parameters:
- i: int
- index of the first mo (Bra)
- j: int
- index of the second mo (Ket)
- drv: {0,1,2,’x’,’y’,’z’}
- The desired component of the vector field of the transition electoronic flux density
Returns:
mo_tefd : numpy.ndarray
-
orbkit.extras.
calc_mo
(qc, fid_mo_list, drv=None, otype=None, ofid=None, numproc=None, slice_length=None)[source]¶ Calculates and saves the selected molecular orbitals or the derivatives thereof. Parameters:
- qc.geo_spec, qc.geo_info, qc.ao_spec, qc.mo_spec :
- See Central Variables for details.
- fid_mo_list : str
- Specifies the filename of the molecular orbitals list or list of molecular
orbital labels (cf.
orbkit.read.mo_select
for details). If fid_mo_list is ‘all_mo’, creates a list containing all molecular orbitals. - drv : string or list of strings {None,’x’,’y’, ‘z’, ‘xx’, ‘xy’, …}, optional
- If not None, a derivative calculation of the molecular orbitals is requested.
- otype : str or list of str, optional
- Specifies output file type. See
otypes
for details. - ofid : str, optional
- Specifies output file name. If None, the filename will be based on
orbkit.options.outputname
. - numproc : int, optional
- Specifies number of subprocesses for multiprocessing.
If None, uses the value from
options.numproc
. - slice_length : int, optional
- Specifies the number of points per subprocess.
If None, uses the value from
options.slice_length
.
- Returns:
- mo_list : numpy.ndarray, shape=((NMO,) + N)
- Contains the NMO=len(qc.mo_spec) molecular orbitals on a grid.
-
orbkit.extras.
gross_atomic_density
(atom, qc, bReturnmo=False, ao_list=None, mo_list=None, drv=None)[source]¶ Computes the gross atomic density with respect to the selected atoms.
\[\rho^a = \sum_i^{N_{\rm MO}} {\rm occ}_i \cdot \varphi_i^a \cdot \varphi_i\]Parameters:
- atom : ‘all’ or int or list of int
- Specifies the atoms (counting from one) for which the gross atomic density will be computed. If (atom == ‘all’) or (atom == -1), computes the gross atomic density for all atoms.
- qc.geo_spec, qc.geo_info, qc.ao_spec, qc.mo_spec :
- See Central Variables for details.
- bReturnmo : bool, optional
- If True, the gross atomic molecular orbitals are additionally returned.
Returns:
- rho_atom : list of numpy.ndarrays, shape=(len(atoms,) + N)
- Contains the atom gross atomic density on a grid.
- mo_atom : list of numpy.ndarrays, shape=(len(atoms,NMO) + N)
- Contains the NMO=len(mo_spec) gross atomic molecular orbitals on a grid.
-
orbkit.extras.
mo_set
(qc, fid_mo_list, drv=None, laplacian=None, otype=None, ofid=None, return_all=True, numproc=None, slice_length=None)[source]¶ Calculates and saves the density or the derivative thereof using selected molecular orbitals.
Parameters:
- qc.geo_spec, qc.geo_info, qc.ao_spec, qc.mo_spec :
- See Central Variables for details.
- fid_mo_list : str
- Specifies the filename of the molecular orbitals list or list of molecular
orbital labels (cf.
orbkit.orbitals.MOClass.select
for details). - otype : str or list of str, optional
- Specifies output file type. See
otypes
for details. - ofid : str, optional
- Specifies output file name. If None, the filename will be based on
orbkit.options.outputname
. - drv : string or list of strings {None,’x’,’y’, ‘z’, ‘xx’, ‘xy’, …}, optional
- If not None, a derivative calculation is requested.
- return_all : bool
- If False, no data will be returned.
- numproc : int, optional
- Specifies number of subprocesses for multiprocessing.
If None, uses the value from
options.numproc
. - slice_length : int, optional
- Specifies the number of points per subprocess.
If None, uses the value from
options.slice_length
.
- Returns:
- datasets : numpy.ndarray, shape=((NSET,) + N)
- Contains the NSET molecular orbital sets on a grid.
- delta_datasets : numpy.ndarray, shape=((NSET,NDRV) + N)
- Contains the NDRV NSET molecular orbital set on a grid. This is only present if derivatives are requested.
-
orbkit.extras.
numerical_mulliken_charges
(atom, qc, ao_list=None, mo_list=None, rho_atom=None)[source]¶ Compute the Mulliken charges and gross atomic populations of the selected atoms numerically using the respective gross atomic densities.
\[P^a = \int \rho^a(x,y,z) {\rm d}^3r\]Parameters:
- atom : ‘all’ or int or list of int
- Specifies the atom (numbering starting from one) to which the numerical Mulliken charges and gross populations will be computed. If (atom == ‘all’) or (atom == -1), computes the charges and populations for all atoms.
Returns:
- rho_atom : list of numpy.ndarrays, shape=(len(atoms,) + N)
- Contains the atom gross atomic density on a grid.
- mulliken_num : dict, (present if not is_vector)
- Contains information of Mulliken charge analysis and has following members:
population: - Mulliken population for each atom.
charges: - Mulliken charges for each atom.