tangelo.toolboxes.molecular_computation package

Subpackages

Submodules

tangelo.toolboxes.molecular_computation.coefficients module

Module containing functions to manipulate molecular coefficient arrays.

tangelo.toolboxes.molecular_computation.coefficients.spatial_from_spinorb(one_body_coefficients, two_body_coefficients)

Function to reverse openfermion.chem.molecular_data.spinorb_from_spatial.

Parameters:
  • one_body_coefficients – One-body coefficients (array of 2N*2N, where N is the number of molecular orbitals).

  • two_body_coefficients – Two-body coefficients (array of 2N*2N*2N*2N, where N is the number of molecular orbitals).

Returns:

(array of floats, array of floats) – One- and two-body integrals (arrays of N*N and N*N*N*N elements, where N is the number of molecular orbitals.

tangelo.toolboxes.molecular_computation.fno module

Module containing datastructures for interfacing with the Frozen Natural Orbitals (FNOs), to automatically truncate the virtual orbital space.

Reference: arXiv:2002.07901

class tangelo.toolboxes.molecular_computation.fno.FNO(sqmol, threshold)

Bases: object

Class to interface with the Frozen Natural Orbitals protocol, that aims at reducing the computational cost of a molecular problem by truncating the virtual space. In general, the virtual orbitals are ranked according to their MP2 occupancies, and selected with a given threshold. They are also transformed to the FNO basis, using the eigenvectors of the MP2 virtual-virtual density matrix.

sqmol

Self-explanatory.

Type:

SecondQuantizedMolecule

uhf

Flag indicating the type of mean field used.

Type:

bool

n_mos

Number of molecular orbitals.

Type:

int

fock_ao

Fock matrix in atomic orbital form.

Type:

np.array

frozen_occupied

List of indices of frozen occupied orbitals.

Type:

list

threshold

Threshold(s) for FNO occupancy.

Type:

float or list

compute_fno(threshold)

Method to compute and truncate the FNO orbitals. It calls the _compute_rfno or the `_compute_ufno`method, whichever is appropriate.

static diagonalize_and_reorder(m, reorder=True)

Method to diagonalize a matrix and reorder the eigenvalues and eigenvectors based on occupations.

Parameters:
  • m (np.array) – The matrix to be diagonalized.

  • reorder (bool) – Flag indicating whether to reorder the eigenvalues and eigenvectors based on occupations. Defaults to True.

Returns:

tuple – A tuple containing the reordered eigenvalues and the corresponding eigenvectors. The eigenvalues represent occupations, and the eigenvectors represent rotation operators.

property fermionic_hamiltonian

Property that returns a FNO fermionic hamiltonian object, with the truncated active space and updated MO coefficients.

Returns:

FermionOperator – Self-explanatory.

get_frozen_indices()

Method to determine the indices of the frozen orbitals, and it calls the _get_restricted_frozen_indices or the _get_unrestricted_frozen_indices method, whether is appropriate.

static get_number_of_fnos_from_frac_occupancies(fno_occ, threshold_frac_occ)

Method to calculate the number of Frozen Natural Orbitals (FNOs) to consider, based on fractional occupancies.

Parameters:
  • fno_occ (np.array) – Array containing fractional occupancies of FNOs.

  • threshold_frac_occ (float) – Threshold value for fractional occupancy.

Returns:

int – Number of FNOs determined by the cumulative sum of fractional occupancies satisfying the provided threshold.

tangelo.toolboxes.molecular_computation.frozen_orbitals module

This module defines functions to get suggestions for freezing orbitals. Those functions take a molecule and return an integer or a list of orbital indices for freezing orbitals. Depending on the function, a Molecule or a SecondQuantizedMolecule object can be used.

tangelo.toolboxes.molecular_computation.frozen_orbitals.convert_frozen_orbitals(sec_mol, frozen_orbitals)

This function converts an int or a list of frozen_orbitals into four categories: - Active and occupied MOs; - Active and virtual MOs; - Frozen and occupied MOs; - Frozen and virtual MOs. Each of them are list with MOs indices (first one is 0). Note that they are MOs labelled, not spin-orbitals (MOs * 2) indices.

Parameters:
  • sec_mol (SecondQuantizedMolecule) – Self-explanatory

  • frozen_orbitals (int or list of int) – Number of MOs or MOs indices to freeze.

Returns:

tuple of list – Active occupied, frozen occupied, active virtual and frozen virtual orbital indices.

tangelo.toolboxes.molecular_computation.frozen_orbitals.get_frozen_core(molecule)

Function to compute the number of frozen orbitals. This function is only for the core (occupied orbitals).

Parameters:

molecule (SecondQuantizedMolecule) – Molecule to be evaluated.

Returns:

int – First N molecular orbitals to freeze.

tangelo.toolboxes.molecular_computation.frozen_orbitals.get_orbitals_excluding_homo_lumo(molecule, homo_minus_n=0, lumo_plus_n=0)

Function that returns a list of orbitals to freeze if the user wants to consider only a subset from HOMO(-homo_min_n) to LUMO(+lumo_plus_n) orbitals. Users should be aware of degeneracies, as this function does not take this property into account.

Parameters:
  • molecule (SecondQuantizedMolecule) – Molecule to be evaluated.

  • homo_minus_n (int) – Starting point at HOMO - homo_minus_n.

  • lumo_plus_n (int) – Ending point at LUMO + lumo_plus_n.

Returns:

list of int – Frozen orbitals not detected in the active space.

tangelo.toolboxes.molecular_computation.integral_solver module

class tangelo.toolboxes.molecular_computation.integral_solver.IntegralSolver

Bases: ABC

Instantiate electronic integral solver

abstract compute_mean_field(sqmol)

Run a unrestricted/restricted (openshell-)Hartree-Fock calculation and modify/add the following variables to sqmol

Modify sqmol variables.

sqmol.mf_energy (float): Mean-field energy (RHF or ROHF energy depending on the spin). sqmol.mo_energies (list of float): Molecular orbital energies. sqmol.mo_occ (list of float): Molecular orbital occupancies (between 0. and 2.). sqmol.n_mos (int): Number of molecular orbitals with a given basis set. sqmol.n_sos (int): Number of spin-orbitals with a given basis set.

Add to sqmol:
self.mo_coeff (ndarray or List[ndarray]): array of molecular orbital coefficients (MO coeffs) if RHF ROHF

list of arrays [alpha MO coeffs, beta MO coeffs] if UHF

Parameters:

sqmol (SecondQuantizedMolecule) –

Populated variables of Molecule plus sqmol.basis (string): Basis set. sqmol.ecp (dict): The effective core potential (ecp) for any atoms in the molecule.

e.g. {“C”: “crenbl”} use CRENBL ecp for Carbon atoms.

sqmol.symmetry (bool or str): Whether to use symmetry in RHF or ROHF calculation.

Can also specify point group using string. e.g. “Dooh”, “D2h”, “C2v”, …

sqmol.uhf (bool): If True, Use UHF instead of RHF or ROHF reference. Default False

abstract get_integrals(sqmol, mo_coeff=None)

Computes core constant, one_body, and two-body integrals for all orbitals

one-body integrals should be in the form h[p,q]= int phi_p(x)* (T + V_{ext}) phi_q(x) dx

two-body integrals should be in the form h[p,q,r,s] = int phi_p(x) * phi_q(y) * V_{elec-elec} phi_r(y) phi_s(x) dxdy

Using molecular orbitals phi_j(x) = sum_{ij} A_i(x) mo_coeff_{i,j} where A_i(x) are the atomic orbitals.

For UHF (if sqmol.uhf is True) one_body coefficients are [alpha one_body, beta one_body] two_body coefficients are [alpha-alpha two_body, alpha-beta two_body, beta-beta two_body]

where one_body and two_body are appropriately sized arrays for each spin sector.

Parameters:
  • sqmol (SecondQuantizedMolecule) – SecondQuantizedMolecule populated with all variables defined above

  • mo_coeff – Molecular orbital coefficients to use for calculating the integrals, instead of self.mo_coeff

Returns:

(float, array or List[array], array or List[array]) – (core_constant, one_body coefficients, two_body coefficients)

abstract set_physical_data(mol)

Set molecular data that is independant of basis set in mol

Modify mol variable:
mol.xyz to (list): Nested array-like structure with elements and coordinates

(ex:[ [“H”, (0., 0., 0.)], …]) in angstrom

Add to mol:

mol.n_electrons (int): Self-explanatory. mol.n_atoms (int): Self-explanatory.

Parameters:

mol (Molecule or SecondQuantizedMolecule) – Class to add the other variables given populated. mol.xyz (in appropriate format for solver): Definition of molecular geometry. mol.q (float): Total charge. mol.spin (int): Absolute difference between alpha and beta electron number.

class tangelo.toolboxes.molecular_computation.integral_solver.IntegralSolverEmpty

Bases: IntegralSolver

compute_mean_field(sqmol)

Run a unrestricted/restricted (openshell-)Hartree-Fock calculation and modify/add the following variables to sqmol

Modify sqmol variables.

sqmol.mf_energy (float): Mean-field energy (RHF or ROHF energy depending on the spin). sqmol.mo_energies (list of float): Molecular orbital energies. sqmol.mo_occ (list of float): Molecular orbital occupancies (between 0. and 2.). sqmol.n_mos (int): Number of molecular orbitals with a given basis set. sqmol.n_sos (int): Number of spin-orbitals with a given basis set.

Add to sqmol:
self.mo_coeff (ndarray or List[ndarray]): array of molecular orbital coefficients (MO coeffs) if RHF ROHF

list of arrays [alpha MO coeffs, beta MO coeffs] if UHF

Parameters:

sqmol (SecondQuantizedMolecule) –

Populated variables of Molecule plus sqmol.basis (string): Basis set. sqmol.ecp (dict): The effective core potential (ecp) for any atoms in the molecule.

e.g. {“C”: “crenbl”} use CRENBL ecp for Carbon atoms.

sqmol.symmetry (bool or str): Whether to use symmetry in RHF or ROHF calculation.

Can also specify point group using string. e.g. “Dooh”, “D2h”, “C2v”, …

sqmol.uhf (bool): If True, Use UHF instead of RHF or ROHF reference. Default False

get_integrals(sqmol, mo_coeff=None)

Computes core constant, one_body, and two-body integrals for all orbitals

one-body integrals should be in the form h[p,q]= int phi_p(x)* (T + V_{ext}) phi_q(x) dx

two-body integrals should be in the form h[p,q,r,s] = int phi_p(x) * phi_q(y) * V_{elec-elec} phi_r(y) phi_s(x) dxdy

Using molecular orbitals phi_j(x) = sum_{ij} A_i(x) mo_coeff_{i,j} where A_i(x) are the atomic orbitals.

For UHF (if sqmol.uhf is True) one_body coefficients are [alpha one_body, beta one_body] two_body coefficients are [alpha-alpha two_body, alpha-beta two_body, beta-beta two_body]

where one_body and two_body are appropriately sized arrays for each spin sector.

Parameters:
  • sqmol (SecondQuantizedMolecule) – SecondQuantizedMolecule populated with all variables defined above

  • mo_coeff – Molecular orbital coefficients to use for calculating the integrals, instead of self.mo_coeff

Returns:

(float, array or List[array], array or List[array]) – (core_constant, one_body coefficients, two_body coefficients)

set_physical_data(mol)

Set molecular data that is independant of basis set in mol

Modify mol variable:
mol.xyz to (list): Nested array-like structure with elements and coordinates

(ex:[ [“H”, (0., 0., 0.)], …]) in angstrom

Add to mol:

mol.n_electrons (int): Self-explanatory. mol.n_atoms (int): Self-explanatory.

Parameters:

mol (Molecule or SecondQuantizedMolecule) – Class to add the other variables given populated. mol.xyz (in appropriate format for solver): Definition of molecular geometry. mol.q (float): Total charge. mol.spin (int): Absolute difference between alpha and beta electron number.

tangelo.toolboxes.molecular_computation.integral_solver_psi4 module

class tangelo.toolboxes.molecular_computation.integral_solver_psi4.IntegralSolverPsi4

Bases: IntegralSolver

psi4 IntegralSolver class

compute_mean_field(sqmol)

Run a unrestricted/restricted (openshell-)Hartree-Fock calculation and modify/add the following variables to sqmol

Modify sqmol variables.

sqmol.mf_energy (float): Mean-field energy (RHF or ROHF energy depending on the spin). sqmol.mo_energies (list of float): Molecular orbital energies. sqmol.mo_occ (list of float): Molecular orbital occupancies (between 0. and 2.). sqmol.n_mos (int): Number of molecular orbitals with a given basis set. sqmol.n_sos (int): Number of spin-orbitals with a given basis set.

Add to sqmol:
self.mo_coeff (ndarray or List[ndarray]): array of molecular orbital coefficients (MO coeffs) if RHF ROHF

list of arrays [alpha MO coeffs, beta MO coeffs] if UHF

Parameters:

sqmol (SecondQuantizedMolecule) –

Populated variables of Molecule plus sqmol.basis (string): Basis set. sqmol.ecp (dict): The effective core potential (ecp) for any atoms in the molecule.

e.g. {“C”: “crenbl”} use CRENBL ecp for Carbon atoms.

sqmol.symmetry (bool or str): Whether to use symmetry in RHF or ROHF calculation.

Can also specify point group using string. e.g. “Dooh”, “D2h”, “C2v”, …

sqmol.uhf (bool): If True, Use UHF instead of RHF or ROHF reference. Default False

compute_uhf_integrals(mo_coeff)

Compute 1-electron and 2-electron integrals The return is formatted as [numpy.ndarray]*2 numpy array h_{pq} for alpha and beta blocks [numpy.ndarray]*3 numpy array storing h_{pqrs} for alpha-alpha, alpha-beta, beta-beta blocks

Parameters:

mo_coeff (List[array]) – The molecular orbital coefficients for both spins [alpha, beta]

Returns:

List[array], List[array] – One and two body integrals

get_integrals(sqmol, mo_coeff=None)

Computes core constant, one_body, and two-body integrals for all orbitals

one-body integrals should be in the form h[p,q]= int phi_p(x)* (T + V_{ext}) phi_q(x) dx

two-body integrals should be in the form h[p,q,r,s] = int phi_p(x) * phi_q(y) * V_{elec-elec} phi_r(y) phi_s(x) dxdy

Using molecular orbitals phi_j(x) = sum_{ij} A_i(x) mo_coeff_{i,j} where A_i(x) are the atomic orbitals.

For UHF (if sqmol.uhf is True) one_body coefficients are [alpha one_body, beta one_body] two_body coefficients are [alpha-alpha two_body, alpha-beta two_body, beta-beta two_body]

where one_body and two_body are appropriately sized arrays for each spin sector.

Parameters:
  • sqmol (SecondQuantizedMolecule) – SecondQuantizedMolecule populated with all variables defined above

  • mo_coeff – Molecular orbital coefficients to use for calculating the integrals, instead of self.mo_coeff

Returns:

(float, array or List[array], array or List[array]) – (core_constant, one_body coefficients, two_body coefficients)

modify_c(wfn, mo_coeff, a=True)

Modify Psi4 WaveFunction coefficient Ca or Cb to be equal to mo_coeff :param wfn: The Psi4Wavefunction to be modified :type wfn: Psi4 WaveFunction :param mo_coeff: The values to place in the wavefunction :type mo_coeff: array :param a: If True, replace the alpha coefficients Ca, If False replace the beta coefficients Cb :type a: bool

modify_solver_mo_coeff(sqmol)

Change the molecular coefficients in the self.wfn object to match self.mo_coeff :param sqmol: The SecondQuantizedMolecule object with mo_coeff :type sqmol: SecondQuantizedMolecule

set_physical_data(mol)

Set molecular data that is independant of basis set in mol

Modify mol variable:
mol.xyz to (list): Nested array-like structure with elements and coordinates

(ex:[ [“H”, (0., 0., 0.)], …]) in angstrom

Add to mol:

mol.n_electrons (int): Self-explanatory. mol.n_atoms (int): Self-explanatory.

Parameters:

mol (Molecule or SecondQuantizedMolecule) – Class to add the other variables given populated. mol.xyz (in appropriate format for solver): Definition of molecular geometry. mol.q (float): Total charge. mol.spin (int): Absolute difference between alpha and beta electron number.

class tangelo.toolboxes.molecular_computation.integral_solver_psi4.IntegralSolverPsi4QMMM(charges)

Bases: IntegralSolverPsi4

Psi4 IntegralSolver class with charges supplied for electrostatic embedding.

Args: charges (List[Tuple[float, Tuple[float, float, float]]]): The partial charges for the electrostatic embedding as

a list with elements (charge, (x, y, z))

compute_mean_field(sqmol)

Run a unrestricted/restricted (openshell-)Hartree-Fock calculation and modify/add the following variables to sqmol

Modify sqmol variables.

sqmol.mf_energy (float): Mean-field energy (RHF or ROHF energy depending on the spin). sqmol.mo_energies (list of float): Molecular orbital energies. sqmol.mo_occ (list of float): Molecular orbital occupancies (between 0. and 2.). sqmol.n_mos (int): Number of molecular orbitals with a given basis set. sqmol.n_sos (int): Number of spin-orbitals with a given basis set.

Add to sqmol:
self.mo_coeff (ndarray or List[ndarray]): array of molecular orbital coefficients (MO coeffs) if RHF ROHF

list of arrays [alpha MO coeffs, beta MO coeffs] if UHF

Parameters:

sqmol (SecondQuantizedMolecule) –

Populated variables of Molecule plus sqmol.basis (string): Basis set. sqmol.ecp (dict): The effective core potential (ecp) for any atoms in the molecule.

e.g. {“C”: “crenbl”} use CRENBL ecp for Carbon atoms.

sqmol.symmetry (bool or str): Whether to use symmetry in RHF or ROHF calculation.

Can also specify point group using string. e.g. “Dooh”, “D2h”, “C2v”, …

sqmol.uhf (bool): If True, Use UHF instead of RHF or ROHF reference. Default False

set_physical_data(mol)

Set molecular data that is independant of basis set in mol

Modify mol variable:
mol.xyz to (list): Nested array-like structure with elements and coordinates

(ex:[ [“H”, (0., 0., 0.)], …]) in angstrom

Add to mol:

mol.n_electrons (int): Self-explanatory. mol.n_atoms (int): Self-explanatory.

Parameters:

mol (Molecule or SecondQuantizedMolecule) –

Class to add the other variables given populated. mol.xyz (in appropriate format for solver): Definition of molecular geometry. If supplying an input file,

nocom and noreorient should be provided.

mol.q (float): Total charge. mol.spin (int): Absolute difference between alpha and beta electron number.

tangelo.toolboxes.molecular_computation.integral_solver_pyscf module

class tangelo.toolboxes.molecular_computation.integral_solver_pyscf.IntegralSolverPySCF(chkfile=None, use_newton=False)

Bases: IntegralSolver

Electronic Structure integration for pyscf

assign_mo_coeff_symmetries(sqmol)

Assigns the symmetry labels for the current molecular coefficients.

Modify sqmol variables mo_symm_ids and mo_symm_labels to match current molecular coefficients.

Parameters:

sqmol (SecondQuantizedMolecule) – The SecondQuantizedMolecule with the symmetrized molecular coefficients to be assigned.

compute_mean_field(sqmol)

Run a unrestricted/restricted (openshell-)Hartree-Fock calculation and modify/add the following variables to sqmol

Modify sqmol variables.

sqmol.mf_energy (float): Mean-field energy (RHF or ROHF energy depending on the spin). sqmol.mo_energies (list of float): Molecular orbital energies. sqmol.mo_occ (list of float): Molecular orbital occupancies (between 0. and 2.). sqmol.n_mos (int): Number of molecular orbitals with a given basis set. sqmol.n_sos (int): Number of spin-orbitals with a given basis set.

Add to sqmol:
self.mo_coeff (ndarray or List[ndarray]): array of molecular orbital coefficients (MO coeffs) if RHF ROHF

list of arrays [alpha MO coeffs, beta MO coeffs] if UHF

Parameters:

sqmol (SecondQuantizedMolecule) –

Populated variables of Molecule plus sqmol.basis (string): Basis set. sqmol.ecp (dict): The effective core potential (ecp) for any atoms in the molecule.

e.g. {“C”: “crenbl”} use CRENBL ecp for Carbon atoms.

sqmol.symmetry (bool or str): Whether to use symmetry in RHF or ROHF calculation.

Can also specify point group using string. e.g. “Dooh”, “D2h”, “C2v”, …

sqmol.uhf (bool): If True, Use UHF instead of RHF or ROHF reference. Default False

compute_uhf_integrals(sqmol, mo_coeff)

Compute 1-electron and 2-electron integrals The return is formatted as [numpy.ndarray]*2 numpy array h_{pq} for alpha and beta blocks [numpy.ndarray]*3 numpy array storing h_{pqrs} for alpha-alpha, alpha-beta, beta-beta blocks

Parameters:
  • sqmol (SecondQuantizedMolecule) – The SecondQuantizedMolecule object to calculated UHF integrals for.

  • mo_coeff (List[array]) – The molecular orbital coefficients for both spins [alpha, beta]

Returns:

List[array], List[array] – One and two body integrals

get_integrals(sqmol, mo_coeff=None)

Computes core constant, one_body, and two-body integrals for all orbitals

one-body integrals should be in the form h[p,q]= int phi_p(x)* (T + V_{ext}) phi_q(x) dx

two-body integrals should be in the form h[p,q,r,s] = int phi_p(x) * phi_q(y) * V_{elec-elec} phi_r(y) phi_s(x) dxdy

Using molecular orbitals phi_j(x) = sum_{ij} A_i(x) mo_coeff_{i,j} where A_i(x) are the atomic orbitals.

For UHF (if sqmol.uhf is True) one_body coefficients are [alpha one_body, beta one_body] two_body coefficients are [alpha-alpha two_body, alpha-beta two_body, beta-beta two_body]

where one_body and two_body are appropriately sized arrays for each spin sector.

Parameters:
  • sqmol (SecondQuantizedMolecule) – SecondQuantizedMolecule populated with all variables defined above

  • mo_coeff – Molecular orbital coefficients to use for calculating the integrals, instead of self.mo_coeff

Returns:

(float, array or List[array], array or List[array]) – (core_constant, one_body coefficients, two_body coefficients)

modify_solver_mo_coeff(sqmol)

Change the molecular coefficients in the SecondQuantizedMolecule mean_field object using self.mo_coeff :param sqmol: The SecondQuantizedMolecule object with mo_coeff :type sqmol: SecondQuantizedMolecule

set_physical_data(mol)

Set molecular data that is independant of basis set in mol

Modify mol variable:
mol.xyz to (list): Nested array-like structure with elements and coordinates

(ex:[ [“H”, (0., 0., 0.)], …]) in angstrom

Add to mol:

mol.n_electrons (int): Self-explanatory. mol.n_atoms (int): Self-explanatory.

Parameters:

mol (Molecule or SecondQuantizedMolecule) – Class to add the other variables given populated. mol.xyz (in appropriate format for solver): Definition of molecular geometry. mol.q (float): Total charge. mol.spin (int): Absolute difference between alpha and beta electron number.

class tangelo.toolboxes.molecular_computation.integral_solver_pyscf.IntegralSolverPySCFQMMM(charges, chkfile=None, use_newton=False)

Bases: IntegralSolverPySCF

Instantiate Electronic Structure integration with PySCF using electrostatic embedding

compute_mean_field(sqmol)

Run a unrestricted/restricted (openshell-)Hartree-Fock calculation and modify/add the following variables to sqmol

Modify sqmol variables.

sqmol.mf_energy (float): Mean-field energy (RHF or ROHF energy depending on the spin). sqmol.mo_energies (list of float): Molecular orbital energies. sqmol.mo_occ (list of float): Molecular orbital occupancies (between 0. and 2.). sqmol.n_mos (int): Number of molecular orbitals with a given basis set. sqmol.n_sos (int): Number of spin-orbitals with a given basis set.

Add to sqmol:
self.mo_coeff (ndarray or List[ndarray]): array of molecular orbital coefficients (MO coeffs) if RHF ROHF

list of arrays [alpha MO coeffs, beta MO coeffs] if UHF

Parameters:

sqmol (SecondQuantizedMolecule) –

Populated variables of Molecule plus sqmol.basis (string): Basis set. sqmol.ecp (dict): The effective core potential (ecp) for any atoms in the molecule.

e.g. {“C”: “crenbl”} use CRENBL ecp for Carbon atoms.

sqmol.symmetry (bool or str): Whether to use symmetry in RHF or ROHF calculation.

Can also specify point group using string. e.g. “Dooh”, “D2h”, “C2v”, …

sqmol.uhf (bool): If True, Use UHF instead of RHF or ROHF reference. Default False

tangelo.toolboxes.molecular_computation.integral_solver_pyscf.mol_to_pyscf(mol, basis='CRENBL', symmetry=False, ecp=None)

Method to return a pyscf.gto.Mole object.

Parameters:
  • sqmol (SecondQuantizedMolecule or Molecule) – The molecule to export to a pyscf molecule.

  • basis (string) – Basis set.

  • symmetry (bool) – Flag to turn symmetry on

  • ecp (dict) – Dictionary with ecp definition for each atom e.g. {“Cu”: “crenbl”}

Returns:

pyscf.gto.Mole – PySCF compatible object.

tangelo.toolboxes.molecular_computation.mm_charges_solver module

class tangelo.toolboxes.molecular_computation.mm_charges_solver.MMChargesSolver

Bases: ABC

Instantiate electrostatic charges solver

abstract get_charges() List[Tuple[float, Tuple[float, float, float]]]

Obtain the charges for the given low parameters of the input fragment (e.g. file names, list(s) of coordinates)

Returns:
  • List[float] – The list of charges

  • List[Tuple (float, Tuple(float, float, float)) – The atoms and geometries of the charges

class tangelo.toolboxes.molecular_computation.mm_charges_solver.MMChargesSolverOpenMM(force_field_params=('amber14-all.xml', 'amber14/tip3pfb.xml'))

Bases: MMChargesSolver

get_charges(files) List[Tuple[float, Tuple[float, float, float]]]

Obtain the charges for the given low parameters of the input fragment (e.g. file names, list(s) of coordinates)

Returns:
  • List[float] – The list of charges

  • List[Tuple (float, Tuple(float, float, float)) – The atoms and geometries of the charges

class tangelo.toolboxes.molecular_computation.mm_charges_solver.MMChargesSolverRDKit(mmffVariant='MMFF94')

Bases: MMChargesSolver

get_charges(files) List[Tuple[float, Tuple[float, float, float]]]

Obtain the charges for the given low parameters of the input fragment (e.g. file names, list(s) of coordinates)

Returns:
  • List[float] – The list of charges

  • List[Tuple (float, Tuple(float, float, float)) – The atoms and geometries of the charges

tangelo.toolboxes.molecular_computation.mm_charges_solver.convert_files_to_pdbs(input_files: str | List[str])

Convert file or list of files to pdb files using openbabel.

Parameters:

input_files (Union[str, List[str]]) – The file name(s) that describe the MM region

Returns:

List[str] – The list of pdb equivalent file name(s).

tangelo.toolboxes.molecular_computation.mm_charges_solver.get_default_mm_package() Type[MMChargesSolver]
tangelo.toolboxes.molecular_computation.mm_charges_solver.get_mm_package(package: str = 'default')

Return the MMCharges instance

Parameters:

package (str) – The MM backend to use

Returns:

Type[MMChargesSolver] – The MMChargesSolver that computes the partial charges for given geometries

tangelo.toolboxes.molecular_computation.molecule module

Module containing datastructures for interfacing with this package functionalities.

class tangelo.toolboxes.molecular_computation.molecule.Molecule(xyz: list, q: int = 0, spin: int = 0, solver: ~tangelo.toolboxes.molecular_computation.integral_solver.IntegralSolver = <factory>)

Bases: object

Custom datastructure to store information about a Molecule. This contains only physical information.

xyz

Nested array-like structure with elements and coordinates (ex:[ [“H”, (0., 0., 0.)], …]). Can also be a multi-line string.

Type:

array-like or string

q

Total charge.

Type:

float

spin

Absolute difference between alpha and beta electron number.

Type:

int

solver

The class that performs the mean field and integral computation.

Type:

IntegralSolver

n_atom

Self-explanatory.

Type:

int

n_electrons

Self-explanatory.

Type:

int

n_min_orbitals

Number of orbitals with a minimal basis set.

Type:

int

Properties:

elements (list): List of all elements in the molecule. coords (array of float): N x 3 coordinates matrix.

property coords

N x 3 coordinates matrix.

Type:

(array of float)

property elements

List of all elements in the molecule.

Type:

(list)

n_atoms: int
n_electrons: int
q: int = 0
solver: IntegralSolver
spin: int = 0
to_file(filename, format=None)

Write molecule geometry to filename in specified format

Parameters:
  • filename (str) – The name of the file to output the geometry.

  • format (str) – The output type of “raw”, “xyz”, or “zmat”. If None, will be inferred by the filename Unless using IntegralSolverPySCF, only “xyz” is supported.

to_openfermion(basis='sto-3g')

Method to return a openfermion.MolecularData object.

Parameters:

basis (string) – Basis set.

Returns:

openfermion.MolecularData – Openfermion compatible object.

xyz: list
class tangelo.toolboxes.molecular_computation.molecule.SecondQuantizedMolecule(xyz: list, q: int = 0, spin: int = 0, solver: ~tangelo.toolboxes.molecular_computation.integral_solver.IntegralSolver = <factory>, basis: str = 'sto-3g', ecp: dict = <factory>, symmetry: bool = False, uhf: bool = False, frozen_orbitals: list = 'frozen_core')

Bases: Molecule

Custom datastructure to store information about a mean field derived from a molecule. This class inherits from Molecule and add a number of attributes defined by the second quantization.

basis

Basis set.

Type:

string

ecp

The effective core potential (ecp) for any atoms in the molecule. e.g. {“C”: “crenbl”} use CRENBL ecp for Carbon atoms.

Type:

dict

symmetry

Whether to use symmetry in RHF or ROHF calculation. Can also specify point group using pyscf allowed string. e.g. “Dooh”, “D2h”, “C2v”, …

Type:

bool or str

uhf

If True, Use UHF instead of RHF or ROHF reference. Default False

Type:

bool

mf_energy

Mean-field energy (RHF or ROHF energy depending on the spin).

Type:

float

mo_energies

Molecular orbital energies.

Type:

list of float

mo_occ

Molecular orbital occupancies (between 0. and 2.).

Type:

list of float

mean_field

Mean-field object (used by CCSD and FCI).

Type:

pyscf.scf

n_mos

Number of molecular orbitals with a given basis set.

Type:

int

n_sos

Number of spin-orbitals with a given basis set.

Type:

int

active_occupied

Occupied molecular orbital indexes that are considered active.

Type:

list of int

frozen_occupied

Occupied molecular orbital indexes that are considered frozen.

Type:

list of int

active_virtual

Virtual molecular orbital indexes that are considered active.

Type:

list of int

frozen_virtual

Virtual molecular orbital indexes that are considered frozen.

Type:

list of int

fermionic_hamiltonian

Self-explanatory.

Type:

FermionOperator

freeze_mos()

Change frozen orbitals attributes. It can be done inplace or not.

Properties:
n_active_electrons (int): Difference between number of total

electrons and number of electrons in frozen orbitals.

n_active_sos (int): Number of active spin-orbitals. n_active_mos (int): Number of active molecular orbitals. frozen_mos (list or None): Frozen MOs indexes. actives_mos (list): Active MOs indexes.

property active_mos

This property returns MOs indexes for the active orbitals.

Returns:

list – MOs indexes that are active (occupied + virtual).

active_occupied: list
property active_spin

This property returns the spin of the active space.

Returns:

int – n_alpha - n_beta electrons of the active occupied orbital space.

active_virtual: list
basis: str = 'sto-3g'
ecp: dict
energy_from_rdms(one_rdm, two_rdm)

Computes the molecular energy from one- and two-particle reduced density matrices (RDMs). Coefficients (integrals) are computed on-the-fly using a pyscf object and the mean-field. Frozen orbitals are supported with this method.

Parameters:
  • one_rdm (array or List[array]) – One-particle density matrix in MO basis. If UHF [alpha one_rdm, beta one_rdm]

  • two_rdm (array or List[array]) – Two-particle density matrix in MO basis. If UHF [alpha-alpha two_rdm, alpha-beta two_rdm, beta-beta two_rdm]

Returns:

float – Molecular energy.

property fermionic_hamiltonian
freeze_mos(frozen_orbitals, inplace=True)

This method recomputes frozen orbitals with the provided input.

property frozen_mos

This property returns MOs indexes for the frozen orbitals. It was written to take into account if one of the two possibilities (occ or virt) is None. In fact, list + None, None + list or None + None return an error. An empty list cannot be sent because PySCF mp2 returns “IndexError: list index out of range”.

Returns:

list – MOs indexes frozen (occupied + virtual).

frozen_occupied: list
frozen_orbitals: list = 'frozen_core'
frozen_virtual: list
get_active_space_integrals(mo_coeff=None)

Computes core constant, one_body, and two-body coefficients with frozen orbitals folded into one-body coefficients and core constant For UHF one_body coefficients are [alpha one_body, beta one_body] two_body coefficients are [alpha-alpha two_body, alpha-beta two_body, beta-beta two_body]

Parameters:

mo_coeff (array) – The molecular orbital coefficients to use to generate the integrals

Returns:

(float, array or List[array], array or List[array]) – (core_constant, one_body coefficients, two_body coefficients)

get_full_space_integrals(mo_coeff=None)

Computes core constant, one_body, and two-body integrals for all orbitals For UHF one_body coefficients are [alpha one_body, beta one_body] two_body coefficients are [alpha-alpha two_body, alpha-beta two_body, beta-beta two_body]

Parameters:

mo_coeff (array) – The molecular orbital coefficients to use to generate the integrals.

Returns:

(float, array or List[array], array or List[array]) – (core_constant, one_body coefficients, two_body coefficients)

get_integrals(mo_coeff=None, fold_frozen=True)

Computes core constant, one_body, and two-body coefficients with frozen orbitals folded into one-body coefficients and core constant for mo_coeff if fold_frozen is True

For UHF one_body coefficients are [alpha one_body, beta one_body] two_body coefficients are [alpha-alpha two_body, alpha-beta two_body, beta-beta two_body]

Parameters:
  • mo_coeff (array) – The molecular orbital coefficients to use to generate the integrals

  • fold_frozen (bool) – If False, the full integral matrices and core constant are returned. If True, the core constant, one_body, and two-body coefficients are calculated with frozen orbitals folded into one-body coefficients and core constant. Default True

Returns:

(float, array or List[array], array or List[array]) – (core_constant, one_body coefficients, two_body coefficients)

mf_energy: float
property mo_coeff

This property returns the molecular orbital coefficients.

Returns:

np.array – MO coefficient as a numpy array.

mo_energies: list
mo_occ: list
property n_active_ab_electrons
property n_active_electrons
property n_active_mos
property n_active_sos
n_mos: int
n_sos: int
symmetry: bool = False
uhf: bool = False
tangelo.toolboxes.molecular_computation.molecule.atom_string_to_list(atom_string)

Convert atom coordinate string (typically stored in text files) into a list/tuple representation suitable for openfermion.MolecularData.

tangelo.toolboxes.molecular_computation.molecule.get_default_integral_solver(qmmm=False)
tangelo.toolboxes.molecular_computation.molecule.get_integral_solver(name: str, qmmm=False)
tangelo.toolboxes.molecular_computation.molecule.molecule_to_secondquantizedmolecule(mol, basis_set='sto-3g', frozen_orbitals=None)

Function to convert a Molecule into a SecondQuantizedMolecule.

Parameters:
  • mol (Molecule) – Self-explanatory.

  • basis_set (string) – String representing the basis set.

  • frozen_orbitals (int or list of int) – Number of MOs or MOs indexes to freeze.

Returns:

SecondQuantizedMolecule – Mean-field data structure for a molecule.

tangelo.toolboxes.molecular_computation.rdms module

Module containing functions to manipulate 1- and 2-RDMs.

tangelo.toolboxes.molecular_computation.rdms.compute_rdms(ferm_ham, mapping, up_then_down, exp_vals=None, exp_data=None, shadow=None, return_spinsum=True, **eval_args)

Compute the 1- and 2-RDM and their spin-summed versions using a FermionOperator and experimental frequency data in the form of a classical shadow or a dictionary of frequency histograms. Exactly one of the following must be provided by the user: exp_vals, exp_data or shadow

Parameters:
  • ferm_ham (FermionicOperator) – Fermionic operator with n_spinorbitals, n_electrons, and spin defined

  • mapping (str) – Qubit mapping

  • up_then_down (bool) – Spin ordering for the mapping

  • exp_vals (dict) – Optional, dictionary of Pauli word expectation values

  • exp_data (dict) – Optional, dictionary of {basis: histogram} where basis is the measurement basis and histogram is a {bitstring: frequency} dictionary

  • shadow (ClassicalShadow) – Optional, a classical shadow object

  • return_spinsum (bool) – Optional, if True, return also the spin-summed RDMs

  • eval_args – Optional arguments to pass to the ClassicalShadow object

Returns:
  • complex array – 1-RDM

  • complex array – 2-RDM

  • complex array – spin-summed 1-RDM

  • complex array – spin-summed 2-RDM

tangelo.toolboxes.molecular_computation.rdms.energy_from_rdms(ferm_op, one_rdm, two_rdm)

Computes the molecular energy from one- and two-particle reduced density matrices (RDMs). Coefficients (integrals) are computed from the fermionic Hamiltonian provided.

Parameters:
  • ferm_op (FermionOperator) – Self-explanatory.

  • one_rdm (numpy.array) – One-particle density matrix in MO basis.

  • two_rdm (numpy.array) – Two-particle density matrix in MO basis.

Returns:

float – Molecular energy.

tangelo.toolboxes.molecular_computation.rdms.matricize_2rdm(two_rdm, n_orbitals)

Turns the two_rdm tensor into a matrix for test purposes.

tangelo.toolboxes.molecular_computation.rdms.pad_rdms_with_frozen_orbitals_restricted(sec_mol, onerdm, twordm)

Function to pad the RDMs with the frozen orbitals data. It is based on the pyscf.cccsd_rdm code, where we can set with_frozen=True.

Source:

https://github.com/pyscf/pyscf/blob/master/pyscf/cc/ccsd_rdm.py

Parameters:
  • sec_mol (SecondQuantizedMolecule) – Self-explanatory.

  • onerdm (numpy.array) – One-particle reduced density matrix (shape of (N_active_mos,)*2).

  • twordm (numpy.array) – Two-particle reduced density matrix (shape of (N_active_mos,)*4).

Returns:
  • numpy.array – One-particle reduced density matrix (shape of (N_total_mos,)*2).

  • numpy.array – Two-particle reduced density matrix (shape of (N_total_mos,)*4).

tangelo.toolboxes.molecular_computation.rdms.pad_rdms_with_frozen_orbitals_unrestricted(sec_mol, onerdm, twordm)

Function to pad the RDMs with the frozen orbitals data. It is based on the pyscf.ucccsd_rdm code, where we can set with_frozen=True.

Source:

https://github.com/pyscf/pyscf/blob/master/pyscf/cc/uccsd_rdm.py

Parameters:
  • sec_mol (SecondQuantizedMolecule) – Self-explanatory.

  • onerdm (tuple of arrays) – Tuple of alpha and beta one-particle reduced density matrix (shape of (N_active_mos,)*2).

  • twordm (tuple of arrays) – Tuple of alpha-alpha, alpha-beta and beta-beta two-particle reduced density matrix (shape of (N_active_mos,)*4).

Returns:
  • tuple of arrays – Tuple of alpha and beta one-particle reduced density matrix (shape of (N_total_mos,)*2).

  • tuple of arrays – Tuple of alpha-alpha, alpha-beta and beta-beta two-particle reduced density matrix (shape of (N_total_mos,)*4).

Module contents