Solver Routines

Orbital Localization

Molecular orbital localization

molbe.lo.localize(self, lo_method, mol=None, valence_basis='sto-3g', hstack=False, pop_method=None, init_guess=None, valence_only=False, nosave=False)

Molecular orbital localization

Performs molecular orbital localization computations. For large basis, IAO is recommended augmented with PAO orbitals.

Parameters:
  • lo_method (str) – Localization method in quantum chemistry. ‘lowdin’, ‘boys’, and ‘iao’ are supported.

  • mol (pyscf.gto.Molecule) – pyscf.gto.Molecule object.

  • valence_basis (str) – Name of minimal basis set for IAO scheme. ‘sto-3g’ suffice for most cases.

  • valence_only (bool) – If this option is set to True, all calculation will be performed in the valence basis in the IAO partitioning. This is an experimental feature.

Crystalline orbital localization

kbe.lo.localize(self, lo_method, mol=None, valence_basis='sto-3g', iao_wannier=True, valence_only=False, iao_val_core=True)

Orbital localization

Performs orbital localization computations for periodic systems. For large basis, IAO is recommended augmented with PAO orbitals.

Parameters:
  • lo_method (str) – Localization method in quantum chemistry. ‘lowdin’, ‘boys’,’iao’, and ‘wannier’ are supported.

  • mol (pyscf.gto.Molecule) – pyscf.gto.Molecule object.

  • valence_basis (str) – Name of minimal basis set for IAO scheme. ‘sto-3g’ suffice for most cases.

  • valence_only (bool) – If this option is set to True, all calculation will be performed in the valence basis in the IAO partitioning. This is an experimental feature.

  • iao_wannier (bool) – Whether to perform Wannier localization in the IAO space

Density Matching Error

molbe.solver.solve_error(Fobjs, Nocc, only_chem=False)

Compute the error for self-consistent fragment density matrix matching.

This function calculates the error in the one-particle density matrix for a given fragment, matching the density matrix elements of the edges and centers. It returns the norm of the error vector and the error vector itself.

Parameters:
  • Fobjs (list of MolBE.fragpart) – List of fragment objects.

  • Nocc (int) – Number of occupied orbitals.

Returns:

  • float – Norm of the error vector.

  • numpy.ndarray – Error vector.

Interface to Quantum Chemistry Methods

molbe.solver.solve_mp2(mf, frozen=None, mo_coeff=None, mo_occ=None, mo_energy=None)

Perform an MP2 (2nd order Moller-Plesset perturbation theory) calculation.

This function sets up and runs an MP2 calculation using the provided mean-field object. It returns the MP2 object after the calculation.

Parameters:
  • mf (pyscf.scf.hf.RHF) – Mean-field object from PySCF.

  • frozen (list or int, optional) – List of frozen orbitals or number of frozen core orbitals. Defaults to None.

  • mo_coeff (numpy.ndarray, optional) – Molecular orbital coefficients. Defaults to None.

  • mo_occ (numpy.ndarray, optional) – Molecular orbital occupations. Defaults to None.

  • mo_energy (numpy.ndarray, optional) – Molecular orbital energies. Defaults to None.

Returns:

The MP2 object after running the calculation.

Return type:

pyscf.mp.mp2.MP2

molbe.solver.solve_ccsd(mf, frozen=None, mo_coeff=None, relax=False, use_cumulant=False, with_dm1=True, rdm2_return=False, mo_occ=None, mo_energy=None, rdm_return=False, verbose=0)

Solve the CCSD (Coupled Cluster with Single and Double excitations) equations.

This function sets up and solves the CCSD equations using the provided mean-field object. It can return the CCSD amplitudes (t1, t2), the one- and two-particle density matrices, and the CCSD object.

Parameters:
  • mf (pyscf.scf.hf.RHF) – Mean-field object from PySCF.

  • frozen (list or int, optional) – List of frozen orbitals or number of frozen core orbitals. Defaults to None.

  • mo_coeff (numpy.ndarray, optional) – Molecular orbital coefficients. Defaults to None.

  • relax (bool, optional) – Whether to use relaxed density matrices. Defaults to False.

  • use_cumulant (bool, optional) – Whether to use cumulant-based energy expression. Defaults to False.

  • with_dm1 (bool, optional) – Whether to include one-particle density matrix in the two-particle density matrix calculation. Defaults to True.

  • rdm2_return (bool, optional) – Whether to return the two-particle density matrix. Defaults to False.

  • mo_occ (numpy.ndarray, optional) – Molecular orbital occupations. Defaults to None.

  • mo_energy (numpy.ndarray, optional) – Molecular orbital energies. Defaults to None.

  • rdm_return (bool, optional) – Whether to return the one-particle density matrix. Defaults to False.

  • verbose (int, optional) – Verbosity level. Defaults to 0.

Returns:

  • t1 (numpy.ndarray): Single excitation amplitudes.

  • t2 (numpy.ndarray): Double excitation amplitudes.

  • rdm1a (numpy.ndarray, optional): One-particle density matrix (if rdm_return is True).

  • rdm2s (numpy.ndarray, optional): Two-particle density matrix (if rdm2_return is True and rdm_return is True).

  • cc__ (pyscf.cc.ccsd.CCSD, optional): CCSD object (if rdm_return is True and rdm2_return is False).

Return type:

tuple

molbe.helper.get_scfObj(h1, Eri, nocc, dm0=None, enuc=0.0, pert_h=False, pert_list=None, save_chkfile=False, fname='f0')

Initialize and run a restricted Hartree-Fock (RHF) calculation.

This function sets up an SCF (Self-Consistent Field) object using the provided one-electron Hamiltonian, electron repulsion integrals, and number of occupied orbitals. It then runs the SCF procedure, optionally using an initial density matrix.

Parameters:
  • h1 (numpy.ndarray) – One-electron Hamiltonian matrix.

  • Eri (numpy.ndarray) – Electron repulsion integrals.

  • nocc (int) – Number of occupied orbitals.

  • dm0 (numpy.ndarray, optional) – Initial density matrix. If not provided, the SCF calculation will start from scratch. Defaults to None.

  • enuc (float, optional) – Nuclear repulsion energy. Defaults to 0.0.

Returns:

mf_ – The SCF object after running the Hartree-Fock calculation.

Return type:

pyscf.scf.hf.RHF

Schmidt Decomposition

Molecular Schmidt decomposition

molbe.solver.schmidt_decomposition(mo_coeff, nocc, Frag_sites, cinv=None, rdm=None, norb=None, return_orb_count=False)

Perform Schmidt decomposition on the molecular orbital coefficients.

This function decomposes the molecular orbitals into fragment and environment parts using the Schmidt decomposition method. It computes the transformation matrix (TA) which includes both the fragment orbitals and the entangled bath.

Parameters:
  • mo_coeff (numpy.ndarray) – Molecular orbital coefficients.

  • nocc (int) – Number of occupied orbitals.

  • Frag_sites (list of int) – List of fragment sites (indices).

  • cinv (numpy.ndarray, optional) – Inverse of the transformation matrix. Defaults to None.

  • rdm (numpy.ndarray, optional) – Reduced density matrix. If not provided, it will be computed from the molecular orbitals. Defaults to None.

  • norb (int, optional) – Specifies number of bath orbitals. Used for UBE to make alpha and beta spaces the same size. Defaults to None

  • return_orb_count (bool, optional) – Return more information about the number of orbitals. Used in UBE. Defaults to False

Returns:

  • numpy.ndarray – Transformation matrix (TA) including both fragment and entangled bath orbitals.

  • if return_orb_count – numpy.ndarray, int, int returns TA (above), number of orbitals in the fragment space, and number of orbitals in bath space

Periodic Schmidt decomposition

kbe.solver.schmidt_decomp_svd(rdm, Frag_sites)

Perform Schmidt decomposition on the orbital coefficients in the real space.

This function decomposes the molecular orbitals into fragment and environment parts using the Schmidt decomposition method. It computes the transformation matrix (TA) which includes both the fragment orbitals and the entangled bath.

Parameters:
  • rdm (numpy.ndarray) – Density matrix (HF) in the real space.

  • Frag_sites (list of int) – List of fragment sites (indices).

Returns:

Transformation matrix (TA) including both fragment and entangled bath orbitals.

Return type:

numpy.ndarray

Handling Hamiltonian

molbe.helper.get_eri(i_frag, Nao, symm=8, ignore_symm=False, eri_file='eri_file.h5')

Retrieve and optionally restore electron repulsion integrals (ERI) from an HDF5 file.

This function reads the ERI for a given fragment from an HDF5 file, and optionally restores the symmetry of the ERI.

Parameters:
  • i_frag (str) – Fragment identifier used to locate the ERI data in the HDF5 file.

  • Nao (int) – Number of atomic orbitals.

  • symm (int, optional) – Symmetry of the ERI. Defaults to 8.

  • ignore_symm (bool, optional) – If True, the symmetry step is skipped. Defaults to False.

  • eri_file (str, optional) – Filename of the HDF5 file containing the electron repulsion integrals. Defaults to ‘eri_file.h5’.

Returns:

Electron repulsion integrals, possibly restored with symmetry.

Return type:

numpy.ndarray

molbe.helper.get_core(mol)

Calculate the number of cores for each atom in the molecule.

Parameters:

mol (pyscf.gto.Mole) – Molecule object from PySCF.

Returns:

(Ncore, idx, corelist)

Return type:

tuple

Build molecular HF potential

molbe.helper.get_veff(eri_, dm, S, TA, hf_veff)

Calculate the effective HF potential (Veff) for a given density matrix and electron repulsion integrals.

This function computes the effective potential by transforming the density matrix, computing the Coulomb (J) and exchange (K) integrals.

Parameters:
  • eri (numpy.ndarray) – Electron repulsion integrals.

  • dm (numpy.ndarray) – Density matrix. 2D array.

  • S (numpy.ndarray) – Overlap matrix.

  • TA (numpy.ndarray) – Transformation matrix.

  • hf_veff (numpy.ndarray) – Hartree-Fock effective potential for the full system.

Returns:

Effective HF potential in the embedding basis.

Return type:

numpy.ndarray

Build perioidic HF potential

kbe.helper.get_veff(eri_, dm, S, TA, hf_veff, return_veff0=False)

Handling Energies

molbe.helper.get_frag_energy(mo_coeffs, nsocc, nfsites, efac, TA, h1, hf_veff, rdm1, rdm2s, dname, eri_file='eri_file.h5', veff0=None)

Compute the fragment energy.

This function calculates the energy contribution of a fragment within a larger molecular system using the provided molecular orbital coefficients, density matrices, and effective potentials.

Parameters:
  • mo_coeffs (numpy.ndarray) – Molecular orbital coefficients.

  • nsocc (int) – Number of occupied orbitals.

  • nfsites (int) – Number of fragment sites.

  • efac (list) – List containing energy scaling factors and indices.

  • TA (numpy.ndarray) – Transformation matrix.

  • h1 (numpy.ndarray) – One-electron Hamiltonian.

  • hf_veff (numpy.ndarray) – Hartree-Fock effective potential.

  • rdm1 (numpy.ndarray) – One-particle density matrix.

  • rdm2s (numpy.ndarray) – Two-particle density matrix.

  • dname (str) – Dataset name in the HDF5 file.

  • eri_file (str, optional) – Filename of the HDF5 file containing the electron repulsion integrals. Defaults to ‘eri_file.h5’.

Returns:

List containing the energy contributions: [e1_tmp, e2_tmp, ec_tmp].

Return type:

list

molbe.rdm.compute_energy_full(self, approx_cumulant=False, use_full_rdm=False, return_rdm=True)

Compute the total energy using rdms in the full basis.

Parameters:
  • approx_cumulant (bool, optional) – If True, use an approximate cumulant for the energy computation. Default is False.

  • use_full_rdm (bool, optional) – If True, use the full two-particle RDM for energy computation. Default is False.

  • return_rdm (bool, optional) – If True, return the computed reduced density matrices (RDMs). Default is True.

Returns:

If return_rdm is True, returns a tuple containing the one-particle and two-particle reduced density matrices (RDM1 and RDM2). Otherwise, returns None.

Return type:

tuple of numpy.ndarray or None

Notes

This function computes the total energy in the full basis, with options to use approximate or true cumulants, and to return the reduced density matrices (RDMs). The energy components are printed as part of the function’s output.

Handling Densities

molbe.rdm.rdm1_fullbasis(self, return_ao=True, only_rdm1=False, only_rdm2=False, return_lo=False, return_RDM2=True, print_energy=False)

Compute the one-particle and two-particle reduced density matrices (RDM1 and RDM2).

Parameters:

return_aobool, optional

Whether to return the RDMs in the AO basis. Default is True.

only_rdm1bool, optional

Whether to compute only the RDM1. Default is False.

only_rdm2bool, optional

Whether to compute only the RDM2. Default is False.

return_lobool, optional

Whether to return the RDMs in the localized orbital (LO) basis. Default is False.

return_RDM2bool, optional

Whether to return the two-particle RDM (RDM2). Default is True.

print_energybool, optional

Whether to print the energy contributions. Default is False.

Returns:

rdm1AOnumpy.ndarray

The one-particle RDM in the AO basis.

rdm2AOnumpy.ndarray

The two-particle RDM in the AO basis (if return_RDM2 is True).

rdm1LOnumpy.ndarray

The one-particle RDM in the LO basis (if return_lo is True).

rdm2LOnumpy.ndarray

The two-particle RDM in the LO basis (if return_lo and return_RDM2 are True).

rdm1MOnumpy.ndarray

The one-particle RDM in the molecular orbital (MO) basis (if return_ao is False).

rdm2MOnumpy.ndarray

The two-particle RDM in the MO basis (if return_ao is False and return_RDM2 is True).