Code documentation¶
-
class
asf.ASFDMRG(*args, **kwargs)¶ Sets up and performs DMRG-CASCI calculations for entropy-based active space selection.
-
calculate(*args, **kwargs) → None¶ Wrapper for the runDMRG method.
Parameters: - *args – positional arguments for self.runDMRG(…)
- **kwargs – keyword arguments for self.runDMRG(…)
-
getOrbDens(root: int = 0) → numpy.ndarray¶ Determines the one-orbital density for the previously calculated DMRG-(CAS)CI wave function.
Parameters: root – state that the density is calculated for Returns: The one-orbital density for all orbitals as an N x 4 numpy array. Order of the columns is empty, spin-up, spin-down, doubly occupied. Raises: Exception– various errors
-
runDMRG(nroots: int = 1, maxM: int = 500, tol: float = 1e-06, **kwargs) → None¶ Run the DMRG-(CAS)CI calculation for subsequent orbital selection.
Parameters: - nroots – number of electronic state to calculate
- maxM – bond dimension
- tol – tolerance
- **kwargs – all other arguments used in the DMRGCI class
Raises: ValueError– bad arguments
-
-
class
asf.ASFCI(*args, **kwargs)¶ Sets up and performs CASCI calculations for entropy-based active space selection.
-
calculate(*args, **kwargs) → None¶ Wrapper for the runCI method.
Parameters: - *args – positional arguments for self.runCI(…)
- **kwargs – keyword arguments for self.runCI(…)
-
getOrbDens(root: int = 0) → numpy.ndarray¶ Determines the one-orbital density for the previously calculated (CAS-)CI wave function.
Parameters: root – state that the density is calculated for Returns: The one-orbital density for all orbitals as an N x 4 numpy array. Order of the columns is empty, spin-up, spin-down, doubly occupied. Raises: Exception– various errors
-
runCI(nroots: int = 1, spin_shift: float = None, **kwargs) → None¶ Run the (CAS-)CI calculation for subsequent orbital selection.
Parameters: - nroots – number of electronic states to calculate
- spin_shift – optionally, the shift for PySCF’s fix_spin function
- **kwargs – arguments that can be passed on the the FCI solver
-
Utility functions
-
asf.utility.compare_active_spaces(cas: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a1d50>, mo_guess: numpy.ndarray, mo_list: Union[Sequence[int], numpy.ndarray], full_result: bool = False) → Union[float, Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]]¶ Given a CASSCF (or CASCI) object, calculate its consistency with some guess orbitals.
This function employs the corresponding orbital transformation.
One way to use the function is just to examine the smallest singular value. If it is close to 1.0, the orbitals are likely consistent. With a value of 0.0 they are most probably inconsistent.
Another way to use the function is to obtain the full set of corresponding orbitals.
Parameters: - cas – Instance of CASSCF (or CASCI).
- mo_guess – Guess orbitals that were used for CASSCF.
- mo_list – List of guess orbitals that were included in the active space.
- full_result – Return the full set of corresponding orbitals if set to True. Otherwise, just provide the smallest singular value.
Returns: The smallest singular value if full_result is False. Tuple (singular values, corresponding guess MOs, corresponding CASSCF MOs) otherwise.
-
asf.utility.corresponding_orbitals(mol: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a1350>, mo_coeff1: numpy.ndarray, mo_coeff2: numpy.ndarray) → Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]¶ Perform corresponding orbital transformation for two orbital sets.
The original idea dates back to Amos and Hall (https://doi.org/10.1098/rspa.1961.0175). See also the work of Neese (https://doi.org/10.1016/j.jpcs.2003.11.015).
Parameters: - mol – Mole object with the basis set and molecular data.
- mo_coeff1 – First set of MO coefficients.
- mo_coeff2 – Second set of MO coefficients.
Returns: Tuple (singular values, first corresponding orbital set, second corresponding orbital set)
-
asf.utility.pictures_Jmol(mol: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a1350>, mo_coeff: numpy.ndarray, mo_list: Union[List[int], numpy.ndarray], name: str = 'mo', cutoff: float = 0.04, ene: Union[numpy.ndarray, List[float]] = None, occ: Union[numpy.ndarray, List[float]] = None, rotate: Tuple[float, float, float] = (0.0, 0.0, 0.0), script_options: str = '', jmol_exec: str = 'jmol', jmol_opts: str = '--nodisplay', xvfb_cmd: str = None) → None¶ Create pictures of orbitals using Jmol.
Note that the indices of orbitals in mo_list start from zero, whereas Jmol counts orbitals starting from one.
Parameters: - mol – The molecule as an instance of pyscf.gto.Mole.
- mo_coeff – Matrix of MO coefficients.
- mo_list – Indices of MOs to be plotted (indexing based on zero).
- name – File names of pictures will be: [name]_[MO index].png
- cutoff – Isosurface cutoff for orbital rendering.
- ene – If provided, orbital energies will included in the MO titles.
- occ – If provided, orbital occupancies will be included in the MO titles.
- rotate – Rotation angle about the (x-axis, y-axis, z-axis) in degrees.
- script_options – A string with arbitrary input to be included in the Jmol script.
- jmol_exec – Name of the Jmol executable.
- jmol_opts – String containing to be provided to Jmol.
- xvfb_cmd – X virtual framebuffer command. None (default): Use ‘xvfb-run’ if it is available. string: Enforces the supplied name for the xfvb-run command. ‘’: Call Jmol without xvfb-run.
-
asf.utility.rdm1s_from_rdm12(N: int, S: float, rdm1: numpy.ndarray, rdm2: numpy.ndarray) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculate the alpha and beta components of the 1-RDM from the spin-free 1-RDM and 2-RDM.
Implements equation 3 from Gidofalvi, Shepard, IJQC 109 (2009), 3552 DOI: 10.1002/qua.22320
Parameters: - N – Number of electrons.
- S – The spin (0.0, 0.5, 1.0, …), such that 2S+1 equals the multiplicity.
- rdm1 – Spin-free one-particle density matrix.
- rdm2 – Spin-free two-particle density matrix in PySCF convention.
Returns: alpha-1-rdm, beta-1-rdm
Raises: ValueError– invalid input
-
asf.utility.tempname(dir: str = <sphinx.ext.autodoc.importer._MockObject object>, suffix: str = None) → Generator[str, None, None]¶ Context manager that creates a temporary file name, and deletes the file at exit.
In contrast to tempfile.TemporaryFile, it only returns the name and not a file object. This permits the user to hand it over to functions which expect a file name instead of a file object, hand it over to other programs for opening and closing, etc.
Parameters: - dir – Directory to store the temporary file.
- suffix – Suffix for the temporary file name.
Yields: The name of the temporary file.
Functions to calculate different types of natural orbitals.
-
asf.natorbs.AOnaturalOrbitals(rdm1ao: numpy.ndarray, S: numpy.ndarray, Sthresh: float = 1e-08) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates natural orbitals for a density matrix in AO basis.
Parameters: - rdm1ao – one-particle reduced density matrix in AO basis
- S – overlap matrix
- Sthresh – eigenvalue cutoff to remove linear dependencies
Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.CCSDNaturalOrbitals(cc: Union[<sphinx.ext.autodoc.importer._MockObject object at 0x7f18d2ff2290>, <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d0d06fd0>]) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates CCSD natural orbitals.
If symmetry is enabled, the natural orbitals will be symmetry-adapted.
Parameters: cc – A CCSD object. Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.MP2NaturalOrbitals(pt: Union[<sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a7b50>, <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a7090>]) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates MP2 natural orbitals.
Attempts to identify restricted/unrestricted basis automatically. If symmetry is enabled, the natural orbitals will be symmetry-adapted.
Parameters: pt – An MP2 object. Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.RCCSDNaturalOrbitals(cc: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d2ff2290>) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates spin-restricted CCSD natural orbitals.
If symmetry is enabled, the natural orbitals will be symmetry-adapted.
Parameters: cc – A CCSD object. Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.RMP2NaturalOrbitals(pt: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a7b50>) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates restricted MP2 natural orbitals.
If symmetry is enabled, the natural orbitals will be symmetry-adapted.
Parameters: pt – An MP2 object. Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.UCCSDNaturalOrbitals(cc: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d0d06fd0>) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates spin-unrestricted CCSD natural orbitals.
If symmetry is enabled, the natural orbitals will be symmetry-adapted.
Parameters: cc – A CCSD object. Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.UHFNaturalOrbitals(mf: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a7a90>) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates the natural orbitals for a converged UHF-type object.
If symmetry is enabled in mf.mol, the natural orbitals will be symmetry-adapted.
Parameters: mf – UHF object with broken spin symmetry Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.UMP2NaturalOrbitals(pt: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a7090>) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates unrestricted MP2 natural orbitals.
If symmetry is enabled, the natural orbitals will be symmetry-adapted.
Parameters: pt – An MP2 object. Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.countActiveElectrons(mo_occ: Union[Sequence[float], numpy.ndarray], mo_list: Union[Sequence[int], numpy.ndarray]) → int¶ Counts the number of active electrons based on occupation numbers.
It is assumed that the occupation numbers outside the active space round to 2 or 0.
Parameters: - mo_occ – list of orbital occupation numbers
- mo_list – list of orbitals in the active space
Returns: number of active electrons
Raises: Exception– various errors
-
asf.natorbs.restrictedNaturalOrbitals(rdm1mo: numpy.ndarray, mo_coeff: numpy.ndarray, mol: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d0cf8650> = None) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates natural orbitals in a basis of spin-restricted MOs.
Parameters: - rdm1mo – one-particle reduced density matrix in MO basis
- mo_coeff – molecular orbital coefficients
- mol – Instance of pyscf.gto.Mole. Note that if mol.symmetry is enabled, the natural orbitals will be symmetry adapted.
Returns: natural occupation numbers, natural orbitals
-
asf.natorbs.selectNaturalOccupations(natocc: numpy.ndarray, lower: float = 0.02, upper: float = 1.98, max_orb: int = None, min_orb: int = None) → Tuple[int, numpy.ndarray]¶ Selects natural orbitals based on their eigenvalues between a lower and an upper boundary.
If a maximal number of orbitals is provided, this function will truncate the orbital list. Orbitals with the highest occupation numbers will be removed first if the space is more than half occupied. Likewise, orbitals with the lowest occupation numbers will be removed first if the space is less than half occupied.
In a similar way, providing a minimum number of orbitals will lead to an extension of the orbital list. The function will proceed such as to arrive closest to a half-occupied space.
Parameters: - natocc – list of natural occupation numbers
- lower – lower boundary for the natural occupation numbers
- upper – upper boundary for the natural occupation numbers
- max_orb – maximal number of orbitals
- min_orb – minimal number of orbitals
Returns: number of electrons, list of MO indices (counting from zero)
Raises: Exception– various errors
-
asf.natorbs.unrestrictedNaturalOrbitals(rdm1mo: Tuple[numpy.ndarray, numpy.ndarray], mo_coeff: Tuple[numpy.ndarray, numpy.ndarray], S: numpy.ndarray, mol: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d0cf8650> = None) → Tuple[numpy.ndarray, numpy.ndarray]¶ Calculates natural orbitals in a basis of spin-unrestricted MOs.
Parameters: - rdm1mo – one-particle reduced density matrix in MO basis
- mo_coeff – molecular orbital coefficients
- S – atomic orbital overlap matrix
- mol – Instance of pyscf.gto.Mole. Note that if mol.symmetry is enabled, the natural orbitals will be symmetry adapted.
Returns: natural occupation numbers, natural orbitals
Simple wrapper functions to find active spaces (semi-)automatically.
-
asf.wrapper.runasf_from_scf(mf: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a38d0>, entropy_threshold: float = 0.13862943611198905, plateau_threshold: float = 0.1, states: Union[int, List[Tuple[int, int]]] = None, sort_mos: bool = False, maxM: int = 250, tol: float = 1e-06, mp2_kwargs: Dict[str, Any] = None, switch_dmrg: int = 12, dmrg_kwargs: Dict[str, Any] = None, fci_kwargs: Dict[str, Any] = None, verbose: bool = True) → Tuple[int, List[int], numpy.ndarray]¶ Entropy-based active space suggestion with DMRG from a set of MP2 natural orbitals.
- The function executes the following steps:
- Calculate MP2 natural orbitals using the SCF object provided in the input.
- Select an initial set of MP2 natural orbitals based on their eigenvalues.
- Perform a DMRG-CASCI calculation using the previously determined orbital subset.
- Select the final active space based on the entropies from the DMRG calculation.
- The output of the function contains:
- the number of active electrons,
- the list of MP2 natural orbitals selected for the active space,
- and the MO coefficient matrix containing all (active and inactive) MP2 natural orbitals.
The columns representing the inactive orbitals in the MO coefficient matrix are ordered according to the convention of PySCF: occupied MOs to the left and virtual MOs to the right. By default, the original ordering of the MP2 natural orbitals is preserved. This can be changed by setting reorder_mos to True: then, the columns of the MO coefficients are reordered such that they can be passed directly to CASCI/CASSCF without any further sorting.
This function does not provide active spaces for individual excited electronic states. However, it can be used to obtain active space suggestions for state-averaged calculations. If a value is provided for ‘states’, the function performs active space selections for all spin and electronic states, and returns a representation of the combined active space.
Parameters: - mf – RHF or UHF object containing converged orbitals.
- entropy_threshold – Single-orbital entropy cutoff for active space selection.
- plateau_threshold – Cutoff for “plateaus” of orbitals with similar entropies.
- states – Electronic states to construct the active space for. None (default): Ground state for the same spin as in mf. Integer n: lowest n electronic states, same spin as in mf. List of Tuples [(s, n)]: lowest n electronic states for each spin s.
- sort_mos – If True, reorder the MO columns as occupied, active, virtual.
- maxM – Maximal DMRG bond dimension.
- tol – Convergence tolerance for DMRG.
- mp2_kwargs – Keyword arguments for MP2 orbital selection (see mp2_from_scf(…)).
- switch_dmrg – Number of orbitals to switch between the regular FCI solver and DMRG.
- dmrg_kwargs – Further keyword arguments to be passed to the DMRG fcisolver.
- fci_kwargs – Further keyword arguments to be passed to CASCI.fcisolver.
- verbose – Print information if set to True.
Returns: Tuple with (no. of act. electrons, list of act. orbitals, orbital coefficients)
-
asf.wrapper.runasf_from_mole(mol: <sphinx.ext.autodoc.importer._MockObject object at 0x7f18d30a3710>, scf_kwargs: Dict[str, Any] = None, stability_analysis: bool = True, verbose: bool = True, **kwargs) → Tuple[int, List[int], numpy.ndarray]¶ Perform SCF calculation and suggest an active space for a molecule.
This function attempts to calculate a stable UHF solution for the Mole object provided. Subsequently, runasf_from_scf(…) is called to calculate MP2 natural orbitals, and to select an active space from those.
For more details, please refer to runasf_from_scf(…).
Parameters: - mol – Mole instance containing the molecular data.
- scf_kwargs – Dictionary of options to be set for the SCF object. Note: this function deviates from PySCF defaults for a few settings.
- stability_analysis – Perform a stability analysis of the SCF solution if True. Re-start the SCF calculation once of the solution is unstable.
- verbose – Print information if set to True.
- **kwargs – Keyword arguments to be passed to runasf_from_scf(…).
Returns: Active space suggestion as provided by runasf_from_scf(…).
Raises: FailedSCF– SCF calculation failed to converge or is unstable.