LOSC Python library¶
Introduction¶
py_losc
is the Python module that provides the Python interface to do
calculations of LOSC. It wraps the
LOSC C++ library with the help of
pybind11 library.
The functionalities provided in py_losc
library are similar to the
LOSC C++ library.
Basic Guide¶
To use Python interface for the LOSC library, all you need to do is
import py_losc
.
To construct LOSC curvature matrix, see this section.
To perform LOSC localization, see this section.
To calculate LOSC corrections, see this section.
To implement post-SCF-LOSC and SCF-LOSC calculations in Python with
the py_losc
module, please refer to this section, which demonstrates a real
example of using py_losc
in psi4 package.
Detailed References for the API¶
Python interface for localized orbital scaling correction (LOSC) library.
Notes
All the matrices are represented and stored in the 2-dimensional numpy.array object. The storage order of matrices is required to be C-style (row major). If the input matrix is not C-style, it will be converted interanally.
Following notations are used in the docstrings:
LO: the localized orbital.
AO: the atomic orbitals.
CO: the canonical orbital.
nlo: the number of LOs.
nbasis: the number of AOs.
nfitbasis: the number of fitting basis for density fitting.
npts: the number of grids.
DFA Representation in LOSC Library¶
-
class
py_losc.
DFAInfo
(gga_x, hf_x, name='')[source]¶ The information of a density functional approximation.
-
__init__
(gga_x, hf_x, name='')[source]¶ Constructor of DFAInfo class that represents a DFA.
- Parameters
name (str) – The description of the DFA. Default to an empty string.
hf_x (float) – The weight of HF exchange in the DFA.
gga_x (float) – The total weights of ALL GGA and LDA type exchange in the DFA.
Examples
B3LYP functional is
\[E^{\rm{B3LYP}}_{\rm{xc}} = E^{\rm{LDA}}_{\rm{x}} + a_0 (E^{\rm{HF}}_{\rm{x}} - E^{\rm{LDA}}_{\rm{x}}) + a_x (E^{\rm{GGA}}_{\rm{x}} - E^{\rm{LDA}}_{\rm{x}}) + E^{\rm{LDA}}_{\rm{c}} + a_c (E^{\rm{GGA}}_{\rm{c}} - E^{\rm{LDA}}_{\rm{c}}),\]in which exchanges end with suffix
_x
and correlations end with suffix_c
, \(a_0=0.20\), \(a_x=0.72\) and \(a_c=0.81\). Only the exchanges are considered and the correlations are ignored. The GGA and LDA exchanges are viewed the same type. Therefore, the total weights of GGA and LDA exchanges are \(1 - a_0 + a_x \times (1 - 1) = 1 - a_0 = 0.80\), and the total weights of HF exchanges is \(a_0 = 0.20\). To construct aDFAInfo
object for B3LYP functional, one should do the following>>> b3lyp = DFAInfo(0.80, 0.20, "B3LYP")
-
gga_x
(self: py_losc.py_losc_core.DFAInfo) → float¶ - Returns
The weight of the total GGA type exchange of the DFA.
- Return type
float
-
hf_x
(self: py_losc.py_losc_core.DFAInfo) → float¶ - Returns
The weight of HF exchange of the DFA.
- Return type
float
-
name
(self: py_losc.py_losc_core.DFAInfo) → str¶ - Returns
The description of the DFA.
- Return type
str
-
LOSC Curvature¶
-
class
py_losc.
CurvatureV1
(dfa_info: py_losc.py_losc.DFAInfo, df_pii: numpy.ndarray, df_vpq_inv: numpy.ndarray, grid_lo: numpy.ndarray, grid_weight: numpy.ndarray)[source]¶ LOSC curvature matrix version 1.
References
Eq. 10 in the original LOSC paper https://doi.org/10.1093/nsr/nwx11.
-
__init__
(dfa_info: py_losc.py_losc.DFAInfo, df_pii: numpy.ndarray, df_vpq_inv: numpy.ndarray, grid_lo: numpy.ndarray, grid_weight: numpy.ndarray)[source]¶ Constructor of LOSC curvature version1.
- Parameters
dfa_info (DFAInfo) – The information for the DFA.
df_pii (numpy.array) – The three-center integral \(\langle p | ii \rangle\) used in density fitting, in which index p refers to the = fitting basis and index i refers to the LO. The dimension of df_pii is [nfitbasis, nlo].
df_Vpq_inv (numpy.array) – The inverse of integral matrix \(\langle p | 1/\mathbf{r} | q \rangle\) used in density fitting, in which indices p and q refer to the fitting basis. The dimension of df_Vpq_inv is [nfitbasis, nfitbasis].
grid_lo (numpy.array) – The LOs values on grid points. The dimension of grid_lo is [npts, nlo].
grid_weight (numpy.array) – The weights for numerical integral over grid points. The dimension of grid_lo is [npts, ].
See also
-
kappa
(self: py_losc.py_losc_core.CurvatureBase) → numpy.ndarray[numpy.float64[m, n]]¶ - Returns
The LOSC curvature matrix with dimension [nlo, nlo].
- Return type
numpy.array
-
nfitbasis
(self: py_losc.py_losc_core.CurvatureBase) → int¶ - Returns
the number of fitting basis in density fitting.
- Return type
int
-
nlo
(self: py_losc.py_losc_core.CurvatureBase) → int¶ - Returns
the number of LOs.
- Return type
int
-
npts
(self: py_losc.py_losc_core.CurvatureBase) → int¶ - Returns
the number of grid points.
- Return type
int
-
set_tau
(self: py_losc.py_losc_core.CurvatureV1, arg0: float) → None¶ - Parameters
tau (float) – The curvature V1 parameter tau.
- Returns
- Return type
None
-
-
class
py_losc.
CurvatureV2
(dfa_info: py_losc.py_losc_core.DFAInfo, df_pii: numpy.ndarray, df_vpq_inv: numpy.ndarray, grid_lo: numpy.ndarray, grid_weight: numpy.ndarray)[source]¶ LOSC curvature matrix version 2.
References
Eq. 8 in the LOSC2 paper (J. Phys. Chem. Lett. 2020, 11, 4, 1528-1535).
-
__init__
(dfa_info: py_losc.py_losc_core.DFAInfo, df_pii: numpy.ndarray, df_vpq_inv: numpy.ndarray, grid_lo: numpy.ndarray, grid_weight: numpy.ndarray)[source]¶ Curvature version 2 class constructor.
See also
Notes
Same interface for input parameters as CurvatureV1.__init__.
-
kappa
(self: py_losc.py_losc_core.CurvatureBase) → numpy.ndarray[numpy.float64[m, n]]¶ - Returns
The LOSC curvature matrix with dimension [nlo, nlo].
- Return type
numpy.array
-
nfitbasis
(self: py_losc.py_losc_core.CurvatureBase) → int¶ - Returns
the number of fitting basis in density fitting.
- Return type
int
-
nlo
(self: py_losc.py_losc_core.CurvatureBase) → int¶ - Returns
the number of LOs.
- Return type
int
-
npts
(self: py_losc.py_losc_core.CurvatureBase) → int¶ - Returns
the number of grid points.
- Return type
int
-
set_tau
(self: py_losc.py_losc_core.CurvatureV2, arg0: float) → None¶ - Parameters
tau (float) – The curvature V2 parameter tau.
- Returns
- Return type
None
-
set_zeta
(self: py_losc.py_losc_core.CurvatureV2, arg0: float) → None¶ - Parameters
zeta (float) – The curvature V2 parameter zeta.
- Returns
- Return type
None
-
LOSC Localization¶
-
class
py_losc.
LocalizerV2
(C_lo_basis: numpy.ndarray, H_ao: numpy.ndarray, D_ao: List[numpy.ndarray])[source]¶ LOSC localization version 2.
References
Eq. 7 in the LOSC2 paper (J. Phys. Chem. Lett. 2020, 11, 4, 1528-1535).
-
__init__
(C_lo_basis: numpy.ndarray, H_ao: numpy.ndarray, D_ao: List[numpy.ndarray])[source]¶ Constructor of LOSC localizerV2.
- Parameters
C_lo_basis (numpy.array) – The coefficient matrix of LO basis that is expanded under AO with dimension [nbasis, nlo].
C_lo_basis[:, i]
is the coefficient of the i-th LO basis. The LO basis refers to a set of MOs that will be unitary transformed to be the LOs. The basis is usually being the COs from the associated DFA.H_ao (numpy.array) – The Hamiltonian matrix under AO representation for the associated DFA. The dimension is [nbasis, nbasis].
D_ao (list of numpy.array) – The dipole matrix under AO representation for x, y and z directions.
-
cost_func
(self: py_losc.py_losc_core.LocalizerBase, lo: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous]) → float¶ Calculate the localization cost function value for the given LOs.
- Parameters
lo (numpy.array) – The LOs coefficient matrix under AOs with dimension [nbasis, nlo].
- Returns
The cost function value in hartree.
- Return type
float
-
is_converged
(self: py_losc.py_losc_core.LocalizerBase) → bool¶ Return the convergence of the most recent localization performed by calling the lo_U() and lo() functions.
- Returns
The convergence of localization.
- Return type
bool
-
lo
(*args, **kwargs)¶ Overloaded function.
lo(self: py_losc.py_losc_core.LocalizerBase, guess: str = ‘identity’) -> numpy.ndarray[numpy.float64[m, n]]
Calculate the LOs’ coefficient matrix under AO.
- guess{‘identity’, ‘random’, ‘random_fixed_seed’}, default=’identity’
Initial guess of the unitary matrix to do localization.
‘identity’: initial U matrix is set as an identity matrix.
‘random’: initial U matrix is set as a random unitary matrix.
‘random_fixed_seed’: initial U matrix is set as a random unitary matrix with fixed random seed.
- C_lonumpy.array
The LO coefficient matrix.
C_lo[:, i]
is the coefficients of the i-th LO represented on AO.
lo_U
lo(self: py_losc.py_losc_core.LocalizerBase, U_guess: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], threshold: float = 1e-08) -> numpy.ndarray[numpy.float64[m, n]]
Calculate the LOs’ coefficient matrix under AO with a given U matrix as the initial guess.
- U_guessnumpy.array
The initial guess of U matrix. Its data will be copied for localization. Its unitarity will be verified and throw an exception if the validation fails.
- thresholdfloat, default=1.0e-8
The threshold used to check the unitarity.
- C_lonumpy.array
The LO coefficient matrix.
C_lo[:, i]
is the coefficients of the i-th LO represented on AO.
lo_U
-
lo_U
(*args, **kwargs)¶ Overloaded function.
lo_U(self: py_losc.py_losc_core.LocalizerBase, guess: str = ‘identity’) -> List[numpy.ndarray[numpy.float64[m, n]]]
Calculate the LOs and the unitary transformation matrix.
- guess{‘identity’, ‘random’, ‘random_fixed_seed’}, default=’identity’
Initial guess of the unitary matrix to do localization.
‘identity’: initial U matrix is set as an identity matrix.
‘random’: initial U matrix is set as a random unitary matrix.
‘random_fixed_seed’: initial U matrix is set as a random unitary matrix with fixed random seed.
- C_lonumpy.array
The LO coefficient matrix.
C_lo[:, i]
is the coefficients of the i-th LO represented on AO.- Unumpy.array
The corresponding unitary transformation matrix for C_lo.
U[:, i]
corresponds to the transformation of i-th LO.
lo_U(self: py_losc.py_losc_core.LocalizerBase, U_guess: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], threshold: float = 1e-08) -> List[numpy.ndarray[numpy.float64[m, n]]]
Calculate the LOs and the unitary transformation matrix with a given U matrix as the initial guess.
- U_guessnumpy.array
The initial guess of U matrix. Its data will be copied for localization. Its unitarity will be verified and throw an exception if the validation fails.
- thresholdfloat, default=1.0e-8
The threshold used to check the unitarity.
- C_lonumpy.array
The LO coefficient matrix.
C_lo[:, i]
is the coefficients of the i-th LO represented on AO.- Unumpy.array
The corresponding unitary transformation matrix for C_lo.
U[:, i]
corresponds to the transformation of i-th LO.
-
set_c
(self: py_losc.py_losc_core.LocalizerV2, c: float) → None¶ Set the localization parameter c.
- Returns
- Return type
None
-
set_convergence
(self: py_losc.py_losc_core.LocalizerBase, arg0: float) → None¶ Set the convergence tolerance for localization.
- Parameters
tol (float) – The convergence tolerance.
- Returns
- Return type
None
-
set_gamma
(self: py_losc.py_losc_core.LocalizerV2, gamma: float) → None¶ Set the localization parameter gamma.
- Returns
- Return type
None
-
set_max_iter
(self: py_losc.py_losc_core.LocalizerBase, arg0: int) → None¶ Set the maximum iteration number for localization.
- Parameters
max_iter (int) – The maximum number of iterations.
- Returns
- Return type
None
-
set_random_permutation
(self: py_losc.py_losc_core.LocalizerBase, arg0: bool) → None¶ Set the flag to perform random permutation or not in the localization with Jacobi-Sweep algorithm.
- Parameters
flag (bool) –
- Returns
- Return type
None
-
steps
(self: py_losc.py_losc_core.LocalizerBase) → int¶ Return the number of iteration steps for the most recent localization performed by calling the lo_U and lo functions.
- Returns
The number of iterations.
- Return type
int
-
LOSC Local Occupation¶
-
py_losc.
local_occupation
(arg0: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg1: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg2: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous]) → numpy.ndarray[numpy.float64[m, n]]¶ Calculate the LOSC local occupation matrix.
- Parameters
C_lo (numpy.array) – The LO coefficient matrix.
C_lo[:, i]
is the coefficients of the i-th LO represented on AO.S (numpy.ndarray) – The AO overlap matrix with dimension [nbasis, nbasis].
D (numpy.ndarray) – The spin density matrix under AO with dimension [nbasis, nbasis].
- Returns
The local occupation matrix with dimension [nlo, nlo].
- Return type
numpy.ndarray
See also
Notes
The local occupation matrix is defined as
\[\lambda_{ij} = \langle \phi_i|\rho|\phi_j \rangle,\]where \(\lambda_{ij}\) is the local occupation matrix element, \(\phi_i\) and \(\phi_j\) are the i-th and j-th localized orbital and \(\rho\) is the spin density operator.
LOSC Corrections¶
-
py_losc.
ao_hamiltonian_correction
(arg0: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg1: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg2: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg3: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous]) → numpy.ndarray[numpy.float64[m, n]]¶ Calculate LOSC effective Hamiltonian under AO basis.
- Parameters
S (numpy.ndarray) – AO overlap matrix with dimension [nbasis, nbasis].
C_lo (numpy.array) – LO coefficient matrix under AO basis with dimension of [nbasis, nbasis]. The coefficients of i-th LO is the i-th column of C_lo.
Curvature (numpy.array) – LOSC curvature matrix with dimension [nlo, nlo].
LocalOcc (numpy.array) – LOSC local occupation matrix with dimension [nlo, nlo].
- Returns
The LOSC effective Hamiltonian under AOs with dimension [nbasis, nbasis].
- Return type
numpy.array
Notes
The LOSC effective Hamiltonian is constructed with LOs fixed. The expression of the effective Hamiltonian is shown as Eq. S25 in the supporting information of the original LOSC paper. This effective Hamiltonian is exact in the developed version of SCF-LOSC
-
py_losc.
energy_correction
(arg0: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg1: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous]) → float¶ Calculate the total energy correction from LOSC.
- Parameters
Curvature (numpy.ndarray) – The LOSC curvature matrix with dimension [nlo, nlo].
LocalOcc (numpy.ndarray) – The LOSC local occupation matrix with dimension [nlo, nlo].
- Returns
The correction from LOSC to the total energy.
- Return type
float
Notes
This is just the energy correction from LOSC, NOT the total energy of LOSC-DFA. Total energy of LOSC-DFA is
E_losc_dfa = E_dfa + E_losc
.
-
py_losc.
orbital_energy_post_scf
(arg0: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg1: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous], arg2: numpy.ndarray[numpy.float64[m, n], flags.c_contiguous]) → List[float]¶ Calculate corrected orbital energy from LOSC in a post-SCF approach.
- Parameters
H_dfa (np.ndarray [nbasis, nbasis]) – The DFA Hamiltonian under AOs.
H_losc (np.ndarray [nbasis, nbasis]) – The LOSC effective Hamiltonian under AOs.
C_co (np.ndarray [nbasis, n], n <= basis, which is the number of COs.) – The coefficient matrix of converged COs under AOs from DFA.
- Returns
out – The corrected orbital energies from LOSC with size of n. The order of orbital energies match the order of input COs (order of columns in C_co matrix).
- Return type
list
See also
Notes
This function gives the final corrected orbital energies from LOSC. Note the difference to the function energy_correction.
The corrected orbital energies are the expectation values of converged DFA’s COs on the LOSC-DFA Hamiltonian, that is,
\[\epsilon_i = \langle \psi_i | H_{\rm{dfa}} + H_{\rm{losc}} | \psi_i \rangle.\]This is just one of the ways to calculate the LOSC corrected orbital energy in a post-SCF LOSC calculation. It is the way usually used to produce results in the published paper for the post-SCF LOSC calculations. Besides this way, there are another two ways to calculate corrected orbital energies: (1) diagonalize the corrected LOSC-DFA Hamiltonian; (2) Follow Eq. 11 to calculate the corrections to orbital energies. These three ways usually produce similar results.