Code documentation
- class asf.ASFDMRG(mol: pyscf.gto.Mole, mo_coeff: ndarray, *args, nroots: int = 1, maxM: int = 500, tol: float = 1e-06, fcisolver_kwargs: dict | None = None, **kwargs)
Sets up and performs DMRG-CASCI calculations for entropy-based active space selection.
- calculate() None
Run the DMRG-(CAS)CI calculation for subsequent orbital selection.
- Raises:
ValueError – bad arguments
- cumulant_4idx(root: int = 0) ndarray
Calculates the complete two-body cumulant in the active space.
- Parameters:
root – the root for which the cumulant is computed
- Returns:
Four-index tensor with the cumulant. Note: 0 labels the first active orbital.
- Raises:
Exception – input error
- diagonal_cumulant(root: int = 0, full_space: bool = True) ndarray
The ‘diagonal’ elements of the two-body cumulant, which derive from <p^+ q^+ q p>.
Optionally, the matrix is filled up with zeros to cover the entire set of orbitals.
- Parameters:
root – the root for which the density is computed
full_space – if true, expand the cumulant to cover the full MO space (filled up with zeros)
- Returns:
Matrix with the relevant entries of the cumulant. Note: 0 labels the first active orbital.
- Raises:
Exception – input error
- entropy_selection(root: int = 0, threshold: float = 0.13862943611198905, plateau_threshold: float = 0.1, useCumulant: bool = True, verbose: bool = True) ActiveSpace
Deprecated: performs orbital selection via a simple entropy cutoff.
Note that this function is considered inferior to ‘find_one_entropy’. It is mainly retained to reproduce results of older ASF versions.
- Parameters:
root – Electronic state to perform the orbital selection for.
threshold – Threshold for the single-site entropies
plateau_threshold – Threshold to identify plateaus
useCumulant – Set whether to incorporate cumulant information in the selection
verbose – Print information or not
- Returns:
number of electrons list of selected orbitals
- find_many(root: int = 0, minimal_entropy: float = 0.07, mo_weight_type: Literal['min_entropy', 'min_edge_sum', 'none'] = 'min_entropy', comparison_tolerance: float = 0.05, keep_minimal: bool = True, filters: list[FilterFunction] | None = None) dict[tuple[int, int], list[MOListInfo]]
Selects a collection of active spaces based on pair information.
By default find_many applies sensible selection criteria to a collection of unfiltered spaces. For a detailed description of the criteria, see asf.filters.find_sensible. Active spaces may also be filtered further following custom criteria via the filters argument. Note that the default filters can only be partially disabled by supplying mo_weight_type=’none’ and minimal_entropy=0.0.
- Parameters:
root – Root to select the active spaces for.
minimal_entropy – Threshold for the minimal single-orbital entropy. Set to zero to disable.
mo_weight_type – Metric by which to compare orbital weights. Set to ‘none’ to disable.
comparison_tolerance – Relative tolerance to consider similar values being approximately equal.
keep_minimal – Whether to retain the minimal active space even if it does not satisfy filter criteria.
filters – A list of additional filter functions. Note that custom filter functions must follow the function signature f(x: MOListInfo) -> bool.
- find_one_entropy(root: int = 0, entropy_threshold: float = 0.13862943611198905, fallback_mode: Literal['strict', 'below_threshold', 'minimal', 'unfiltered'] = 'minimal', filters: list[FilterFunction] | None = None, **find_many_kwargs) ActiveSpace
Select a sensible active space based on a threshold for the minimal entropy.
Among the find_one variants, this is the recommended function to select a single active space. To determine a single active space, find_one_entropy first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. To select a single choice, only spaces with a minimal entropy of any orbital above the specified threshold are considered, and the space that maximizes the pair information sum is returned.
Since it may occur that no active space from the primary pool matches the search criteria several fallback modes are provided (see argument description), where the original (i.e. “strict”) search parameters are slightly altered to find a reasonable alternative. The active space searches modes of find_one_entropy are carried out following the order:
“strict” -> “below_threshold” -> “minimal” -> “unfiltered”.
By default, at most the “minimal” fallback is completed before an exception is thrown.
- Parameters:
root – Root to select the active space for.
entropy_threshold – Entropy threshold to select a single active space; all spaces need to have a minimum value above this threshold.
fallback_mode –
Sets the fallback behavior in case no active space could be found: - strict is equivalent to no fallback and raises an exception if
no space matching the criteria could be found.
below_threshold weakens the entropy criterion and searches for sensible spaces below the given entropy threshold.
minimal is the default and returns the minimal active space if the previous modes did not yield a result. This can, for instance, be useful for open-shell systems.
unfiltered extends the search for an active space to the collection of unfiltered spaces. This option can be useful, if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
filters – A list of additional filter functions. Note that custom filter must be of type FilterFunction as defined in asf.filters.
find_many_kwargs – Arguments for find_many excluding filters. Use these arguments to tune the number of active spaces in the primary pool of candidate spaces.
- Returns:
Single active space with minimal entropy above the specified threshold and with maximal pair information sum.
- find_one_generic(selector: SelectorFunction | Sequence[SelectorFunction], fallback_unfiltered: bool, fallback_minimal: bool, filters: list[FilterFunction] | None = None, root: int = 0, minimal_entropy: float = 0.07, mo_weight_type: Literal['min_entropy', 'min_edge_sum', 'none'] = 'min_entropy', comparison_tolerance: float = 0.05) ActiveSpace | None
Find single active space based on pair information and user-defined criteria.
This function does not make assumptions about filters/selectors (other than the criteria used in find_many) and is intended for expert users. To determine a single active space, find_one_generic first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. In case no active space could be selected from the primary pool of candidate spaces and fallback_unfiltered=True, the same search is applied to the collection of unfiltered active spaces (see unfiltered_pairinfo_spaces). If no active space can be found with either method, an error will be raised.
- Parameters:
root – Root to select the active space for.
fallback_unfiltered – Extends the search for an active space to the collection of unfiltered spaces. This option can be useful if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
fallback_minimal – Whether to retain the minimal space among the list of sensible spaces even if does not satisfy the filtering criteria otherwise.
filters – A list of filter functions. The custom filter functions must be of type FilterFunction defined in asf.filters.
selector – One or multiple selector functions. A selector is a custom function to pick one space from several viable choices. Note that a selector must be of type SelectorFunction as defined in asf.filters. If multiple selectors are provided, then all but the first selector act as fallback options. If the first selector fails to select one option, the second selector is applied, etc., until one selection succeeds.
minimal_entropy – Threshold for the minimal single-orbital entropy. Set to zero to disable.
mo_weight_type – Metric by which to compare orbital weights. Set to ‘none’ to disable.
comparison_tolerance – Relative tolerance to consider similar values being approximately equal.
- Returns:
Single active space matching given criteria.
- find_one_sized(norb: int, nel: int | None = None, root: int = 0, fallback_unfiltered: bool = False, filters: list[FilterFunction] | None = None, selector: SelectorFunction | Sequence[SelectorFunction] | None = None, **find_many_kwargs) ActiveSpace
Find optimal active space of given size based on pair information and other criteria.
To determine a single active space, find_one_sized first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. In case no active space could be selected from the primary pool of candidate spaces and the fallback is enabled (fallback_unfiltered=True), the same search is applied to the collections of unfiltered active spaces (see unfiltered_pairinfo_spaces). If no space could be determined with either method, an error will be raised.
- Parameters:
root – Root to select the active space for.
norb – Requested number of active orbitals.
nel – Requested number of active electrons.
fallback_unfiltered – Extends the search for an active space to the collection of unfiltered spaces. This option can be useful, if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
filters – A list of additional filter functions. Note that custom filter must be of type FilterFunction as defined in asf.filters.
selector – A custom function or a list of custom functions to pick one space from potentially several cases. By default, the space with maximum pair information sum is returned. Note that a selector must be of type SelectorFunction as defined in asf.filters, or a list thereof. If a list is provided, the selectors are applied successively until the first one does not return None.
find_many_kwargs – Arguments for find_many excluding filters. Use these arguments to tune the number of active spaces in the primary pool of candidate spaces.
- Raises:
ActiveSpaceSelectionError – no space found for given criteria
- Returns:
Single active space matching given size criteria.
- classmethod from_active_space(mol: pyscf.gto.Mole, space: ActiveSpace, **kwargs) Self
Initialize an ASF object from an ActiveSpace instance.
- Parameters:
mol – molecule (as pyscf.gto.Mole instance)
space – active space
**kwargs – Settings to pass on to the class initializer.
- Returns:
ASFBase instance (or one derived from ASFBase)
- classmethod from_preselection(strategy: Preselection, **kwargs) Self
Determine an initial space and initialize an ASFBase class with it.
This alternative initializer first carries out a step to determine an initial space following a preselection approach and then creates an ASFBase instance (e.g. ASFDMRG) with the MO coefficients, indices of active MOs, and number of active electrons arising from the initial space.
- Parameters:
strategy – Approach for selecting initial active space. The ASF package provides a small collection of such approaches in the asf.preselection module. Custom preselection procedures are also supported as long as they derive from the Preselection protocol (asf.asfbase.Preselection).
**kwargs – Settings to pass on to the class initializer.
- Returns:
A correlated calculation object derived from ASFBase instantiated with the determined initial active space.
- one_orbital_density(root: int = 0) ndarray
Determines the one-orbital density for the 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
- one_orbital_entropy(orbdens: ndarray | None = None, root: int = 0) ndarray
Calculates the one-orbital entropy based on a one-orbital density.
- Parameters:
orbdens – Orbital density to calculate the entropy with.
root – electronic state to obtain the density for
- Returns:
Vector of single-site entropies for the orbitals.
- rdm1s(root: int = 0) ndarray
Returns the active one-electron spin density matrix as an 2 x N x N array.
- rdm2(root: int = 0) ndarray
Returns the active two-electron spin-free density matrix as a four-dimensional array.
- unfiltered_pairinfo_spaces(root: int = 0) dict[tuple[int, int], list[MOListInfo]]
Determine a large number of candidate active spaces with very light screening.
For each number of active electrons and active MOs, (nel, nmo), multiple active spaces are determined. Most of those may not be suitable for converging CASSCF calculations.
- Parameters:
root – The electronic state for which to identify active spaces.
- Returns:
Dictionary with active space suggestions and various information as defined in ActiveSpacesDict. The MO lists are mapped to the full MO space.
- class asf.ASFCI(mol: pyscf.gto.Mole, mo_coeff: ndarray, *args, nroots: int = 1, spin_shift: float | None = None, fcisolver_kwargs: dict | None = None, **kwargs)
Sets up and performs CASCI calculations for entropy-based active space selection.
- calculate() None
Run the (CAS-)CI calculation for subsequent orbital selection.
- cumulant_4idx(root: int = 0) ndarray
Calculates the complete two-body cumulant in the active space.
- Parameters:
root – the root for which the cumulant is computed
- Returns:
Four-index tensor with the cumulant. Note: 0 labels the first active orbital.
- Raises:
Exception – input error
- diagonal_cumulant(root: int = 0, full_space: bool = True) ndarray
The ‘diagonal’ elements of the two-body cumulant, which derive from <p^+ q^+ q p>.
Optionally, the matrix is filled up with zeros to cover the entire set of orbitals.
- Parameters:
root – the root for which the density is computed
full_space – if true, expand the cumulant to cover the full MO space (filled up with zeros)
- Returns:
Matrix with the relevant entries of the cumulant. Note: 0 labels the first active orbital.
- Raises:
Exception – input error
- entropy_selection(root: int = 0, threshold: float = 0.13862943611198905, plateau_threshold: float = 0.1, useCumulant: bool = True, verbose: bool = True) ActiveSpace
Deprecated: performs orbital selection via a simple entropy cutoff.
Note that this function is considered inferior to ‘find_one_entropy’. It is mainly retained to reproduce results of older ASF versions.
- Parameters:
root – Electronic state to perform the orbital selection for.
threshold – Threshold for the single-site entropies
plateau_threshold – Threshold to identify plateaus
useCumulant – Set whether to incorporate cumulant information in the selection
verbose – Print information or not
- Returns:
number of electrons list of selected orbitals
- find_many(root: int = 0, minimal_entropy: float = 0.07, mo_weight_type: Literal['min_entropy', 'min_edge_sum', 'none'] = 'min_entropy', comparison_tolerance: float = 0.05, keep_minimal: bool = True, filters: list[FilterFunction] | None = None) dict[tuple[int, int], list[MOListInfo]]
Selects a collection of active spaces based on pair information.
By default find_many applies sensible selection criteria to a collection of unfiltered spaces. For a detailed description of the criteria, see asf.filters.find_sensible. Active spaces may also be filtered further following custom criteria via the filters argument. Note that the default filters can only be partially disabled by supplying mo_weight_type=’none’ and minimal_entropy=0.0.
- Parameters:
root – Root to select the active spaces for.
minimal_entropy – Threshold for the minimal single-orbital entropy. Set to zero to disable.
mo_weight_type – Metric by which to compare orbital weights. Set to ‘none’ to disable.
comparison_tolerance – Relative tolerance to consider similar values being approximately equal.
keep_minimal – Whether to retain the minimal active space even if it does not satisfy filter criteria.
filters – A list of additional filter functions. Note that custom filter functions must follow the function signature f(x: MOListInfo) -> bool.
- find_one_entropy(root: int = 0, entropy_threshold: float = 0.13862943611198905, fallback_mode: Literal['strict', 'below_threshold', 'minimal', 'unfiltered'] = 'minimal', filters: list[FilterFunction] | None = None, **find_many_kwargs) ActiveSpace
Select a sensible active space based on a threshold for the minimal entropy.
Among the find_one variants, this is the recommended function to select a single active space. To determine a single active space, find_one_entropy first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. To select a single choice, only spaces with a minimal entropy of any orbital above the specified threshold are considered, and the space that maximizes the pair information sum is returned.
Since it may occur that no active space from the primary pool matches the search criteria several fallback modes are provided (see argument description), where the original (i.e. “strict”) search parameters are slightly altered to find a reasonable alternative. The active space searches modes of find_one_entropy are carried out following the order:
“strict” -> “below_threshold” -> “minimal” -> “unfiltered”.
By default, at most the “minimal” fallback is completed before an exception is thrown.
- Parameters:
root – Root to select the active space for.
entropy_threshold – Entropy threshold to select a single active space; all spaces need to have a minimum value above this threshold.
fallback_mode –
Sets the fallback behavior in case no active space could be found: - strict is equivalent to no fallback and raises an exception if
no space matching the criteria could be found.
below_threshold weakens the entropy criterion and searches for sensible spaces below the given entropy threshold.
minimal is the default and returns the minimal active space if the previous modes did not yield a result. This can, for instance, be useful for open-shell systems.
unfiltered extends the search for an active space to the collection of unfiltered spaces. This option can be useful, if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
filters – A list of additional filter functions. Note that custom filter must be of type FilterFunction as defined in asf.filters.
find_many_kwargs – Arguments for find_many excluding filters. Use these arguments to tune the number of active spaces in the primary pool of candidate spaces.
- Returns:
Single active space with minimal entropy above the specified threshold and with maximal pair information sum.
- find_one_generic(selector: SelectorFunction | Sequence[SelectorFunction], fallback_unfiltered: bool, fallback_minimal: bool, filters: list[FilterFunction] | None = None, root: int = 0, minimal_entropy: float = 0.07, mo_weight_type: Literal['min_entropy', 'min_edge_sum', 'none'] = 'min_entropy', comparison_tolerance: float = 0.05) ActiveSpace | None
Find single active space based on pair information and user-defined criteria.
This function does not make assumptions about filters/selectors (other than the criteria used in find_many) and is intended for expert users. To determine a single active space, find_one_generic first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. In case no active space could be selected from the primary pool of candidate spaces and fallback_unfiltered=True, the same search is applied to the collection of unfiltered active spaces (see unfiltered_pairinfo_spaces). If no active space can be found with either method, an error will be raised.
- Parameters:
root – Root to select the active space for.
fallback_unfiltered – Extends the search for an active space to the collection of unfiltered spaces. This option can be useful if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
fallback_minimal – Whether to retain the minimal space among the list of sensible spaces even if does not satisfy the filtering criteria otherwise.
filters – A list of filter functions. The custom filter functions must be of type FilterFunction defined in asf.filters.
selector – One or multiple selector functions. A selector is a custom function to pick one space from several viable choices. Note that a selector must be of type SelectorFunction as defined in asf.filters. If multiple selectors are provided, then all but the first selector act as fallback options. If the first selector fails to select one option, the second selector is applied, etc., until one selection succeeds.
minimal_entropy – Threshold for the minimal single-orbital entropy. Set to zero to disable.
mo_weight_type – Metric by which to compare orbital weights. Set to ‘none’ to disable.
comparison_tolerance – Relative tolerance to consider similar values being approximately equal.
- Returns:
Single active space matching given criteria.
- find_one_sized(norb: int, nel: int | None = None, root: int = 0, fallback_unfiltered: bool = False, filters: list[FilterFunction] | None = None, selector: SelectorFunction | Sequence[SelectorFunction] | None = None, **find_many_kwargs) ActiveSpace
Find optimal active space of given size based on pair information and other criteria.
To determine a single active space, find_one_sized first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. In case no active space could be selected from the primary pool of candidate spaces and the fallback is enabled (fallback_unfiltered=True), the same search is applied to the collections of unfiltered active spaces (see unfiltered_pairinfo_spaces). If no space could be determined with either method, an error will be raised.
- Parameters:
root – Root to select the active space for.
norb – Requested number of active orbitals.
nel – Requested number of active electrons.
fallback_unfiltered – Extends the search for an active space to the collection of unfiltered spaces. This option can be useful, if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
filters – A list of additional filter functions. Note that custom filter must be of type FilterFunction as defined in asf.filters.
selector – A custom function or a list of custom functions to pick one space from potentially several cases. By default, the space with maximum pair information sum is returned. Note that a selector must be of type SelectorFunction as defined in asf.filters, or a list thereof. If a list is provided, the selectors are applied successively until the first one does not return None.
find_many_kwargs – Arguments for find_many excluding filters. Use these arguments to tune the number of active spaces in the primary pool of candidate spaces.
- Raises:
ActiveSpaceSelectionError – no space found for given criteria
- Returns:
Single active space matching given size criteria.
- classmethod from_active_space(mol: pyscf.gto.Mole, space: ActiveSpace, **kwargs) Self
Initialize an ASF object from an ActiveSpace instance.
- Parameters:
mol – molecule (as pyscf.gto.Mole instance)
space – active space
**kwargs – Settings to pass on to the class initializer.
- Returns:
ASFBase instance (or one derived from ASFBase)
- classmethod from_preselection(strategy: Preselection, **kwargs) Self
Determine an initial space and initialize an ASFBase class with it.
This alternative initializer first carries out a step to determine an initial space following a preselection approach and then creates an ASFBase instance (e.g. ASFDMRG) with the MO coefficients, indices of active MOs, and number of active electrons arising from the initial space.
- Parameters:
strategy – Approach for selecting initial active space. The ASF package provides a small collection of such approaches in the asf.preselection module. Custom preselection procedures are also supported as long as they derive from the Preselection protocol (asf.asfbase.Preselection).
**kwargs – Settings to pass on to the class initializer.
- Returns:
A correlated calculation object derived from ASFBase instantiated with the determined initial active space.
- one_orbital_density(root: int = 0) 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
- one_orbital_entropy(orbdens: ndarray | None = None, root: int = 0) ndarray
Calculates the one-orbital entropy based on a one-orbital density.
- Parameters:
orbdens – Orbital density to calculate the entropy with.
root – electronic state to obtain the density for
- Returns:
Vector of single-site entropies for the orbitals.
- rdm1s(root: int = 0) ndarray
Returns the active one-electron spin density matrix as an 2 x N x N array.
- rdm2(root: int = 0) ndarray
Returns the active two-electron spin-free density matrix as a four-dimensional array.
- unfiltered_pairinfo_spaces(root: int = 0) dict[tuple[int, int], list[MOListInfo]]
Determine a large number of candidate active spaces with very light screening.
For each number of active electrons and active MOs, (nel, nmo), multiple active spaces are determined. Most of those may not be suitable for converging CASSCF calculations.
- Parameters:
root – The electronic state for which to identify active spaces.
- Returns:
Dictionary with active space suggestions and various information as defined in ActiveSpacesDict. The MO lists are mapped to the full MO space.
asf.asfbase
The classes and functions contained in the asfbase
module fulfill a more fundamental role. For
instance, the ActiveSpace
class is a representation of an active space which is often used as
a return type. Another notable function is print_mo_table
which shows one-orbital densities and
the one-orbital entropy for selected molecular orbitals.
Core functions to calculate entropies and select orbitals.
- class asf.asfbase.ASFBase(mol: pyscf.gto.Mole, mo_coeff: ndarray, nel: int | None = None, norb: int | None = None, mo_list: list[int] | ndarray | None = None)
Abstract base class for correlated methods to perform orbital selection.
The general idea is to select orbitals based on criteria including single-site entropies, which are computed from an approximate correlated method, for example DMRG. Implementations for specific quantum chemical methods should be derived from this base class. A user can perform automatic active space selection by following the steps below:
# Initialize the object. ASFMethod is a subclass of ASFBase, implementing a specific quantum # chemical method. An initial orbital subspace to perform the calculation may be specified. sf = ASFMethod(…)
# The correlated calculation (e.g. DMRG) is performed, storing information for the subsequent # orbital selection in this object. sf.calculate(…)
# Finally, an active space is selected as a subset of the initially provided orbital space. # The method returns the number of selected active electrons and a list of active MOs. nel, mo_list = sf.entropy_selection(…)
- abstract calculate() None
Perform a correlated calculation, which generates the information to select orbitals.
- abstract diagonal_cumulant(root: int = 0, full_space: bool = True) ndarray
The ‘diagonal’ elements of the two-body cumulant, which derive from <p^+ q^+ q p>.
Optionally, the matrix is filled up with zeros to cover the entire set of orbitals.
- Parameters:
root – the root for which the density is computed
full_space – if true, expand the cumulant to cover the full MO space (filled up with zeros)
- Returns:
Matrix with the relevant entries of the cumulant. Note: 0 labels the first active orbital.
- Raises:
Exception – input error
- entropy_selection(root: int = 0, threshold: float = 0.13862943611198905, plateau_threshold: float = 0.1, useCumulant: bool = True, verbose: bool = True) ActiveSpace
Deprecated: performs orbital selection via a simple entropy cutoff.
Note that this function is considered inferior to ‘find_one_entropy’. It is mainly retained to reproduce results of older ASF versions.
- Parameters:
root – Electronic state to perform the orbital selection for.
threshold – Threshold for the single-site entropies
plateau_threshold – Threshold to identify plateaus
useCumulant – Set whether to incorporate cumulant information in the selection
verbose – Print information or not
- Returns:
number of electrons list of selected orbitals
- find_many(root: int = 0, minimal_entropy: float = 0.07, mo_weight_type: Literal['min_entropy', 'min_edge_sum', 'none'] = 'min_entropy', comparison_tolerance: float = 0.05, keep_minimal: bool = True, filters: list[FilterFunction] | None = None) dict[tuple[int, int], list[MOListInfo]]
Selects a collection of active spaces based on pair information.
By default find_many applies sensible selection criteria to a collection of unfiltered spaces. For a detailed description of the criteria, see asf.filters.find_sensible. Active spaces may also be filtered further following custom criteria via the filters argument. Note that the default filters can only be partially disabled by supplying mo_weight_type=’none’ and minimal_entropy=0.0.
- Parameters:
root – Root to select the active spaces for.
minimal_entropy – Threshold for the minimal single-orbital entropy. Set to zero to disable.
mo_weight_type – Metric by which to compare orbital weights. Set to ‘none’ to disable.
comparison_tolerance – Relative tolerance to consider similar values being approximately equal.
keep_minimal – Whether to retain the minimal active space even if it does not satisfy filter criteria.
filters – A list of additional filter functions. Note that custom filter functions must follow the function signature f(x: MOListInfo) -> bool.
- find_one_entropy(root: int = 0, entropy_threshold: float = 0.13862943611198905, fallback_mode: Literal['strict', 'below_threshold', 'minimal', 'unfiltered'] = 'minimal', filters: list[FilterFunction] | None = None, **find_many_kwargs) ActiveSpace
Select a sensible active space based on a threshold for the minimal entropy.
Among the find_one variants, this is the recommended function to select a single active space. To determine a single active space, find_one_entropy first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. To select a single choice, only spaces with a minimal entropy of any orbital above the specified threshold are considered, and the space that maximizes the pair information sum is returned.
Since it may occur that no active space from the primary pool matches the search criteria several fallback modes are provided (see argument description), where the original (i.e. “strict”) search parameters are slightly altered to find a reasonable alternative. The active space searches modes of find_one_entropy are carried out following the order:
“strict” -> “below_threshold” -> “minimal” -> “unfiltered”.
By default, at most the “minimal” fallback is completed before an exception is thrown.
- Parameters:
root – Root to select the active space for.
entropy_threshold – Entropy threshold to select a single active space; all spaces need to have a minimum value above this threshold.
fallback_mode –
Sets the fallback behavior in case no active space could be found: - strict is equivalent to no fallback and raises an exception if
no space matching the criteria could be found.
below_threshold weakens the entropy criterion and searches for sensible spaces below the given entropy threshold.
minimal is the default and returns the minimal active space if the previous modes did not yield a result. This can, for instance, be useful for open-shell systems.
unfiltered extends the search for an active space to the collection of unfiltered spaces. This option can be useful, if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
filters – A list of additional filter functions. Note that custom filter must be of type FilterFunction as defined in asf.filters.
find_many_kwargs – Arguments for find_many excluding filters. Use these arguments to tune the number of active spaces in the primary pool of candidate spaces.
- Returns:
Single active space with minimal entropy above the specified threshold and with maximal pair information sum.
- find_one_generic(selector: SelectorFunction | Sequence[SelectorFunction], fallback_unfiltered: bool, fallback_minimal: bool, filters: list[FilterFunction] | None = None, root: int = 0, minimal_entropy: float = 0.07, mo_weight_type: Literal['min_entropy', 'min_edge_sum', 'none'] = 'min_entropy', comparison_tolerance: float = 0.05) ActiveSpace | None
Find single active space based on pair information and user-defined criteria.
This function does not make assumptions about filters/selectors (other than the criteria used in find_many) and is intended for expert users. To determine a single active space, find_one_generic first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. In case no active space could be selected from the primary pool of candidate spaces and fallback_unfiltered=True, the same search is applied to the collection of unfiltered active spaces (see unfiltered_pairinfo_spaces). If no active space can be found with either method, an error will be raised.
- Parameters:
root – Root to select the active space for.
fallback_unfiltered – Extends the search for an active space to the collection of unfiltered spaces. This option can be useful if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
fallback_minimal – Whether to retain the minimal space among the list of sensible spaces even if does not satisfy the filtering criteria otherwise.
filters – A list of filter functions. The custom filter functions must be of type FilterFunction defined in asf.filters.
selector – One or multiple selector functions. A selector is a custom function to pick one space from several viable choices. Note that a selector must be of type SelectorFunction as defined in asf.filters. If multiple selectors are provided, then all but the first selector act as fallback options. If the first selector fails to select one option, the second selector is applied, etc., until one selection succeeds.
minimal_entropy – Threshold for the minimal single-orbital entropy. Set to zero to disable.
mo_weight_type – Metric by which to compare orbital weights. Set to ‘none’ to disable.
comparison_tolerance – Relative tolerance to consider similar values being approximately equal.
- Returns:
Single active space matching given criteria.
- find_one_sized(norb: int, nel: int | None = None, root: int = 0, fallback_unfiltered: bool = False, filters: list[FilterFunction] | None = None, selector: SelectorFunction | Sequence[SelectorFunction] | None = None, **find_many_kwargs) ActiveSpace
Find optimal active space of given size based on pair information and other criteria.
To determine a single active space, find_one_sized first considers a pool of reasonable active spaces (using the criteria set out in find_many) that will likely converge in a CASSCF or DMRG-SCF calculation. In case no active space could be selected from the primary pool of candidate spaces and the fallback is enabled (fallback_unfiltered=True), the same search is applied to the collections of unfiltered active spaces (see unfiltered_pairinfo_spaces). If no space could be determined with either method, an error will be raised.
- Parameters:
root – Root to select the active space for.
norb – Requested number of active orbitals.
nel – Requested number of active electrons.
fallback_unfiltered – Extends the search for an active space to the collection of unfiltered spaces. This option can be useful, if the desired space is not contained in the primary set of sensible spaces. Note that spaces selected from the fallback pool may not converge in a CASSCF calculation.
filters – A list of additional filter functions. Note that custom filter must be of type FilterFunction as defined in asf.filters.
selector – A custom function or a list of custom functions to pick one space from potentially several cases. By default, the space with maximum pair information sum is returned. Note that a selector must be of type SelectorFunction as defined in asf.filters, or a list thereof. If a list is provided, the selectors are applied successively until the first one does not return None.
find_many_kwargs – Arguments for find_many excluding filters. Use these arguments to tune the number of active spaces in the primary pool of candidate spaces.
- Raises:
ActiveSpaceSelectionError – no space found for given criteria
- Returns:
Single active space matching given size criteria.
- classmethod from_active_space(mol: pyscf.gto.Mole, space: ActiveSpace, **kwargs) Self
Initialize an ASF object from an ActiveSpace instance.
- Parameters:
mol – molecule (as pyscf.gto.Mole instance)
space – active space
**kwargs – Settings to pass on to the class initializer.
- Returns:
ASFBase instance (or one derived from ASFBase)
- classmethod from_preselection(strategy: Preselection, **kwargs) Self
Determine an initial space and initialize an ASFBase class with it.
This alternative initializer first carries out a step to determine an initial space following a preselection approach and then creates an ASFBase instance (e.g. ASFDMRG) with the MO coefficients, indices of active MOs, and number of active electrons arising from the initial space.
- Parameters:
strategy – Approach for selecting initial active space. The ASF package provides a small collection of such approaches in the asf.preselection module. Custom preselection procedures are also supported as long as they derive from the Preselection protocol (asf.asfbase.Preselection).
**kwargs – Settings to pass on to the class initializer.
- Returns:
A correlated calculation object derived from ASFBase instantiated with the determined initial active space.
- abstract one_orbital_density(root: int = 0) ndarray
Calculate the one-orbital density.
- Parameters:
root – electronic state to obtain the density for
- Returns:
The one-orbital density as a N x 4 matrix for N orbitals.
- one_orbital_entropy(orbdens: ndarray | None = None, root: int = 0) ndarray
Calculates the one-orbital entropy based on a one-orbital density.
- Parameters:
orbdens – Orbital density to calculate the entropy with.
root – electronic state to obtain the density for
- Returns:
Vector of single-site entropies for the orbitals.
- abstract rdm1s(root: int = 0) ndarray
Returns the active one-electron spin density matrix as an 2 x N x N array.
- abstract rdm2(root: int = 0) ndarray
Returns the active two-electron spin-free density matrix as a four-dimensional array.
- unfiltered_pairinfo_spaces(root: int = 0) dict[tuple[int, int], list[MOListInfo]]
Determine a large number of candidate active spaces with very light screening.
For each number of active electrons and active MOs, (nel, nmo), multiple active spaces are determined. Most of those may not be suitable for converging CASSCF calculations.
- Parameters:
root – The electronic state for which to identify active spaces.
- Returns:
Dictionary with active space suggestions and various information as defined in ActiveSpacesDict. The MO lists are mapped to the full MO space.
- class asf.asfbase.ActiveSpace(nel: int, mo_list: list[int] | ndarray, mo_coeff: ndarray)
Class to store information about an active space.
- nel
Number of active electrons.
- Type:
int
- norb
Number of active orbitals.
- Type:
int
- mo_list
List of active MO indices (referring to columns of mo_coeff).
- Type:
list[int]
- mo_coeff
Matrix of all MO coefficients (active and inactive).
- Type:
numpy.ndarray
- class DataDict
Serialized data of the attributes in an ActiveSpace object.
- from_active_indices(mo_list: list[int]) list[int]
Map relative MO indices within the active space onto the full set of MOs.
- Parameters:
mo_list – List of relative indices within the active space (relative to self.mo_list).
- Returns:
mo_list mapped to the full set of orbitals (= relative to the columns of mo_coeff).
- classmethod from_dict(data: DataDict) ActiveSpace
Create a new ActiveSpace object from a dictionary of serialized data.
- Parameters:
data – Data as returned by self.to_dict().
- Raises:
Exception – cannot instantiate with given data
- Returns:
A new instance of this class.
- merge_with(mol: pyscf.gto.Mole, other: ActiveSpace | None = None, mo_coeff_rtol: float = 1e-05, mo_coeff_atol: float = 1e-08) ActiveSpace
Merge two active spaces based on the same MOs.
- This function is based on the following assumptions:
A sensible active space can be obtained simply by merging mo_list1 and mo_list2.
Orbitals that are not already active in both active space definitions “bring” either two electrons (occup. MO) or no electrons (virtual MO) into the merged active space.
Inactive MOs are always ordered: doubly occupied MOs come before virtual orbitals.
For further details, see merge_active_spaces.
- Raises:
ValueError – MO coefficient matrices do not match
- Parameters:
mol – Mole object with molecular data.
other – Other ActiveSpace instance. If not provided, the current instance is returned.
mo_coeff_rtol – Relative tolerance for comparing a pair of MO coefficient matrices.
mo_coeff_atol – Absolute tolerance for comparing a pair of MO coefficient matrices.
- to_active_indices(mo_list: list[int]) list[int]
Map MO indices onto relative indices within the active space.
- Parameters:
mo_list – List of indices referring to the full set of MOs (= columns of mo_coeff).
- Returns:
mo_list mapped relative to self.mo_list (= relative to the active space indices).
- class asf.asfbase.Preselection(*args, **kwargs)
A general class interface for preselection strategies.
- property mol: pyscf.gto.Mole
Molecule associated with the preselection and active space search.
- select() ActiveSpace
Perform the preselection.
- Returns:
An active space as an ActiveSpace instance
- asf.asfbase.all_orbital_selections(orbdens: ndarray, **kwargs) dict[tuple[int, int], list[int]]
Generate a list of all active spaces that can be generated by varying the entropy threshold.
- Parameters:
orbdens – one-orbital density in an N x 4 array
**kwargs – keyword arguments for the entropy_selection function
- Returns:
Dictionary containing all the different active space suggestions. The keys are tuples containing (no. of electrons, no. of orbitals). The values are the lists of active orbitals.
- Raises:
ValueError – invalid arguments supplied
- asf.asfbase.calc_ncore(mol: pyscf.gto.Mole, nel_act: int) int
Calculates the number of inactive / core orbitals from the number of active electrons.
- Parameters:
mol – molecular information
nel_act – number of active electrons
- Returns:
the number of core orbitals
- Raises:
ValueError – error
- asf.asfbase.entropy_selection(orbdens: ndarray, pair_cumulant: ndarray | None = None, threshold: float = 0.13862943611198905, plateau_threshold: float = 0.1, zero_cutoff: float = 1e-12, cumulant_minimum_threshold: float = 1e-06) tuple[int, list[int]]
Select orbitals with entropies above a given threshold.
In addition, the list of chosen orbitals is extended - until the number of unpaired electrons is correct, - until the separation between entropies is sufficiently large (‘plateaus’), - and with orbitals which are predominantly singly occupied.
- Parameters:
orbdens – one-orbital density with columns representing empty, spin-up, spin-down, doubly occupied
pair_cumulant – If provided, incorporate cumulant information in the orbital selection.
threshold – Normal threshold for the entropy (unless there is a plateau)
plateau_threshold – Relative threshold to identify ‘plateaus’
zero_cutoff – Threshold to neglect zeroes when calculating the entropy
cumulant_minimum_threshold – Do not consider cumulant information below a ‘noise’ cutoff
- Returns:
number of electrons list of orbitals
- Raises:
Exception – various errors
- asf.asfbase.inactive_orbital_lists(ncore: int, nmo: int, active_list: list[int] | ndarray) tuple[ndarray, ndarray]
Obtain lists of inactive core and virtual orbitals from an active space specification.
- Parameters:
ncore – number of core orbitals
nmo – total number of orbitals
active_list – list of active orbitals
- Returns:
list of core orbitals, list of virtual orbitals
- Raises:
ValueError – erratic input encountered
- asf.asfbase.merge_active_spaces(mol: pyscf.gto.Mole, nel1: int, mo_list1: list[int], nel2: int, mo_list2: list[int]) tuple[int, list[int]]
Perform a simple merging of active spaces.
This function attempts to construct an active space for state-averaged CASSCF calculations, by forming the union of the lists of active orbitals for individual electronic states.
- This function is based on the following assumptions:
A sensible active space can be obtained simply by merging mo_list1 and mo_list2.
Orbitals that are not already active in both active space definitions “bring” either two electrons (occupied MO) or no electrons (virtual MO) into the merged active space.
Inactive MOs are always ordered: doubly occupied orbitals come before virtual orbitals.
- Parameters:
mol – Mole object with molecular data.
nel1 – Number of electrons in the first active space.
mo_list1 – List of orbitals in the first active space.
nel2 – number of electrons in the second active space.
mo_list2 – List of orbitals in the second active space.
- Returns:
(number of electrons, list of orbitals)
- Return type:
Tuple representing the merged active space
- Raises:
ValueError – Insonsistent input provided.
- asf.asfbase.print_mo_table(orbital_densities: ndarray, mo_list: list[int] | None = None, selections: list[str] | dict[str, Collection[int]] | None = None) None
Print an MO table with one-orbital densities, entropies and whether they are selected.
Note, that if selections are passed as a dictionary, the collections are not checked for overlapping indices.
- Parameters:
orbital_densities – single-orbital densities as N x 4 matrix.
mo_list – List of indices to print. If no list is provided, all MOs will be printed
selections – Indices of selected MO and their corresponding status labels OR a list of status labels.
asf.filters
The filters
module implements functions to filter or select from collections of active spaces.
Apart from find_sensible
which determines a set of reasonable active spaces, FilterFunction
and SelectorFunction
are useful objects, e.g. for defining custom filters or selectors.
Functions to filter and select from active spaces.
- exception asf.filters.ActiveSpaceSelectionError
ASF Exception related to selecting active spaces.
- class asf.filters.FilterFunction(*args, **kwargs)
General function interface for active space filter functions.
Given a set of active spaces (here represented as MOListInfo objects) it is often desired to create a subset following certain constraints. These constraints can be expressed in the form of filter functions as defined below. A subset can then be easily created via Python’s built-in filter(). Note: due to a more complex function signature and limitations of Callable in Python 3.9, the function type is defined as a callback protocol instead.
- class asf.filters.SelectorFunction(*args, **kwargs)
General function interface for active space selector functions.
Conceptually, selector functions ensure that only a single active space (or none) is picked from a set of active spaces (here represented as M̀OListInfo` objects). Hence, selectors often entail an optimization problem, i.e. finding the active space that optimizes an objective function. Moreover, they can involve filtering functionality prior to the selection. Note: due to a more complex function signature and limitations of Callable in Python 3.9, the function type is defined as a callback protocol instead.
We may encounter a situation where a selection fails, either due to an empty list of active spaces, or due to an unsatisfied selection criterion. In that case, the selector function is permitted to return None.
- asf.filters.apply_filters(spaces: list[MOListInfo] | dict[tuple[int, int], list[MOListInfo]], filters: list[FilterFunction]) list[MOListInfo]
Apply several filter functions to each element of an active spaces collection.
Note that if an ActiveSpacesDict is passed, the filtering will not occur per CAS tuple group but on the entire collection of active spaces as a flat list.
- Parameters:
spaces – List of active spaces (as MOListInfo).
filters – List of filter functions that follow the function signature f(x: MOListInfo) -> bool.
- Raises:
TypeError – invalid active space collection type
- Returns:
A list of active spaces matching the given criteria
- asf.filters.apply_selector(spaces: list[MOListInfo] | dict[tuple[int, int], list[MOListInfo]], selector: SelectorFunction | Sequence[SelectorFunction]) MOListInfo | None
Apply a selector to select at most a single element from an active spaces collection.
If multiple selectors are provided, then the first selector is applied initially. In case it fails to select an active space and returns None, the second selector is applied. This is continued until either a single active space has been selected, or the sequence of selectors has been exhausted.
Note that if an ActiveSpacesDict is passed, it will be converted to a flat list prior to the selection step.
- Parameters:
spaces – Active spaces as a list of MOListInfo objects, or as an ActiveSpacesDict.
selector – Individual selector function or a sequence of selector functions with the signature f(x: list[MOListInfo]) -> Optional[MOListInfo]. If multiple selectors are provided, their priority is determined by their ordering.
- Returns:
Object with the selected active space, or None if nothing could be selected.
- asf.filters.filter_entropy_above_thresh(space: MOListInfo, threshold: float, keep_minimal: bool = False) bool
Filter active spaces that have a minimal entropy above a given threshold.
Note that since the entropy is non-negative, also given the threshold must be larger than zero. If that is not the case, the filter is essentially disabled by returning True for each space. This function is supposed to be used in a list(filter(…)) statement.
- Parameters:
space – An element from a list of active spaces (as MOListInfo).
threshold – Threshold for the smallest single-orbital entropies.
keep_minimal – Whether to retain the minimal active space in the filtered set.
- asf.filters.filter_increments(space: MOListInfo, keep_minimal: bool, rel_tol: float = 1e-09, abs_tol: float = 1e-06) bool
Filter active spaces by increment and decrement information.
Active spaces are only retained if the largest increment is smaller than or approximately equal to the smallest decrement. This function is supposed to be used in a list(filter(…)) statement.
- Parameters:
space – An element from a list of active spaces (as MOListInfo).
keep_minimal – Whether to keep a space flagged as minimal.
rel_tol – Relative comparison tolerance between increments and decrements.
abs_tol – Absolute comparison tolerance between increments and decrements.
- asf.filters.filter_nel_strict(space: MOListInfo, nel: int) bool
Filter active spaces that strictly have a given number of MOs.
This function is supposed to be used in a list(filter(…)) statement.
- Parameters:
space – An element from a list of active spaces (as MOListInfo).
nel – Number of active orbitals the desired space should have.
- asf.filters.filter_node_weight_isclose(space: MOListInfo, weight_type: Literal['min_entropy', 'min_edge_sum'], target_value: float, rel_tol: float = 1e-09, abs_tol: float = 1e-06) bool
Filter active spaces that have a ‘node weight’ close to a given target value.
The ‘node weights’, i.e. orbital-assignable metrics, can be either the minimal single-orbital entropy (min_entropy) or the minimum sum over all edges (min_edge_sum).
- Parameters:
space – An element from a list of active spaces (as MOListInfo).
weight_type – Use entropies or edge sums for comparison.
target_value – Target value to compare with node weight.
rel_tol – Relative comparison tolerance between node weight values.
abs_tol – Absolute comparison tolerance between node weight values.
- asf.filters.filter_norb_strict(space: MOListInfo, norb: int) bool
Filter active spaces that strictly have a given number of MOs.
This function is supposed to be used in a list(filter(…)) statement.
- Parameters:
space – An element from a list of active spaces (as MOListInfo).
norb – Number of active orbitals the desired space should have.
- asf.filters.find_sensible(spaces: dict[tuple[int, int], list[MOListInfo]], minimal_entropy: float, mo_weight_type: Literal['min_entropy', 'min_edge_sum', 'none'], comparison_tolerance: float, keep_minimal: bool = False, filters: list[FilterFunction] | None = None) list[MOListInfo]
Selects a collection of reasonable active spaces based on pair information.
- By default the following criteria are used:
For each group of active spaces with same number of active electrons and active orbitals, (nel, nmo), the active space must include the orbitals with the largest node weights (to within a tolerance).
The smallest entropy of any orbital in the active space must be above the threshold.
The smallest decrement must exceed the largest increment.
The spaces can additionally be filtered according to custom criteria by passing filter functions via the filters argument. Note that when using custom filters, criteria 1), 2), and 3) are still applied. If a custom filtering without criteria 1-3) is desired, use apply_filters instead.
- Parameters:
spaces – Dictionary containing active spaces with associated information.
minimal_entropy – Threshold for criterion (2). Set to zero to disable.
mo_weight_type – Metric by which to compare orbital weights.
comparison_tolerance – Relative tolerance to establish equality.
keep_minimal – Whether to retain the minimal active space even if the filter criteria do not yield any space.
filters – A list of custom filter functions. Note that filter functions must follow the function signature f(x: MOListInfo) -> bool.
- Returns:
List of active space suggestions.
- asf.filters.screen_node_weights(spaces: dict[tuple[int, int], list[MOListInfo]], weight_type: Literal['min_entropy', 'min_edge_sum'], rel_tol: float, abs_tol: float = 1e-06) dict[tuple[int, int], list[MOListInfo]]
Screen active spaces for their minimal ‘node weight’ in each active space.
The ‘node weights’ can be either entropies or sums over all edges (to active and inactive MOs). Active spaces are not inspected individually, but in comparison with all other active spaces with the same number of active electrons and orbitals. The active space with the largest minimal node weight is kept, and in addition also active spaces if their minimum node weight is sufficiently close.
- Parameters:
spaces – Dictionary containing active spaces with associated information.
weight_type – Use entropies or edge sums for comparison.
rel_tol – Relative comparison tolerance between node weight values.
abs_tol – Absolute comparison tolerance between node weight values.
- Returns:
Dictionary of active spaces that has been filtered according to the criteria.
- asf.filters.select_max_pairinfo_sum(spaces: list[MOListInfo]) MOListInfo | None
Select the active space with maximum pair information sum.
- Parameters:
spaces – List of active spaces (as MOListInfo).
- Returns:
The selected space, or None.
- asf.filters.select_min_entropy_diff(spaces: list[MOListInfo], target_entropy: float, rel_tol: float) MOListInfo | None
Select active space that has minimum entropy closest to a target entropy.
The “distance” between target entropy and minimum entropy of a space is calculated as the absolute logarithm of their ratio. Because it is quite likely that the spaces contain two or more active spaces with identical values of min_entropy, the pair information sum is used as a secondary criterion to pick one active space that has a value of min_entropy closest to the target entropy and the maximum pair information sum. Active spaces with nearly identical values of min_entropy are identified using math.isclose.
- Parameters:
spaces – List of active spaces (as MOListInfo).
target_entropy – Target single-orbital entropy to find active space with closest minimum single-orbital entropy (‘min_entropy’).
rel_tol – Relative tolerance to determine whether an entropy is close to the target entropy.
- Returns:
Representation of the selected active space.
- Raises:
ValueError – invalid entropy value
- asf.filters.select_pairinfo_with_entropy(spaces: list[MOListInfo], entropy_threshold: float, keep_minimal: bool) MOListInfo | None
This selector combines entropy threshold filtering with pairinfo maximization.
First, it filters the spaces such that only those with a minimal entropy above a threshold are left. Among the surviving spaces, it selects the one that contains the largest pair information sum.
- Parameters:
spaces – List of active spaces (as MOListInfo) to select from.
entropy_threshold – Minimum entropy threshold used to filter out spaces.
keep_minimal – If set to True, the minimal space will always be retained in filtering regardless of its orbital entropies. If set to False, it will be filtered out if entropies are too low.
asf.molmath
Reusable functions involving molecular descriptors such as orbitals and densities.
- asf.molmath.corresponding_orbitals(S12: ndarray, mo_coeff1: ndarray, mo_coeff2: ndarray) tuple[ndarray, ndarray, ndarray]
Calculates a corresponding orbital transformation between different MO 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:
S12 – overlap matrix between MOs 1 and MOs 2, can be rectangular
mo_coeff1 – first set of MO coefficients
mo_coeff2 – second set of MO coefficients
- Returns:
corresponding orbital set 1, corresponding orbital set 2, singular values
- asf.molmath.overlap_square_roots(S: ndarray, Sthresh: float = 0.0) tuple[ndarray, ndarray]
Calculates square root and inverse square root of a matrix S.
- Parameters:
S – (overlap) matrix
Sthresh – eigenvalue cutoff to remove linear dependencies
- Returns:
S^1/2, S^-1/2
asf.mp2density_conventional
Implementation of modified Møller-Plesset reduced density matrices with conventional MP2.
- asf.mp2density_conventional.make_rdm1_rmp2(pt: pyscf.mp.mp2.MP2, mo_energy: ndarray) ndarray
Calculate MP(2, 1) one-particle reduced density matrix (RHF reference).
- Parameters:
pt – Restricted MP2 object.
mo_energy – MO orbital energies as a vector.
- Returns:
1-RDM in MO basis
- asf.mp2density_conventional.make_rdm1_ump2(pt: pyscf.mp.ump2.UMP2, mo_energy: tuple[ndarray, ndarray]) ndarray
Calculate MP(2, 1) one-particle reduced density matrix (UHF reference).
- Parameters:
pt – Unrestricted MP2 object.
mo_energy – MO orbital energies in tuple of vectors (alpha MO energies, beta MO energies).
- Returns:
1-RDM in MO basis as an array of dimensions 2 x N(MO) x N(MO), with alpha and beta parts.
- asf.mp2density_conventional.make_rdm2_0th_ump2(nocca: int, noccb: int, rdm2s: ndarray) None
Zeroth-order contribution to the UHF 2-RDM: Hartree-Fock terms.
- Parameters:
nocca – Number of occupied spin-up orbitals.
noccb – Number of occupied spin-down orbitals.
rdm2s – The 2-RDM contribution is added to this array.
- asf.mp2density_conventional.make_rdm2_1st_ump2(t2aa: ndarray, t2ab: ndarray, t2bb: ndarray, rdm2s: ndarray) None
First-order contribution to the UHF-MP(2, 1) 2-RDM.
- Parameters:
t2aa – All-alpha MP2 amplitudes.
t2ab – Alpha-beta MP2 amplitudes.
t2bb – All-beta MP2 amplitudes.
rdm2s – The 2-RDM contribution is added to this array.
- asf.mp2density_conventional.make_rdm2_2nd_separable_ump2(nocca: int, noccb: int, rdm1s_mp2: ndarray, rdm2s: ndarray) None
2nd order contributions to the UHF-MP(2, 1) 2-RDM that separate into products of 1-RDMs.
- Parameters:
nocca – Number of occupied alpha orbitals.
noccb – Number of occupied beta orbitals.
rdm1s_mp2 – Full UHF-MP(2, 1) 1-RDM, dimensions 2 x N(MO) x N(MO).
rdm2s – The 2-RDM contribution is added to this array.
- asf.mp2density_conventional.make_rdm2_contractions_ump2(mol: pyscf.gto.Mole, nocca: int, noccb: int, mo_coeff: ndarray | tuple[ndarray, ndarray], mo_energy: ndarray | tuple[ndarray, ndarray], t2aa: ndarray, t2ab: ndarray, t2bb: ndarray, rdm2s: ndarray) None
Second-order contributions to UHF-MP(2, 1) 2-RDM that are amplitude-integral contractions.
- Parameters:
mol – PySCF Mole object containing molecular data.
nocca – Number of occupied alpha orbitals.
noccb – Number of occupied beta orbitals.
mo_coeff – Molecular orbital coefficients.
mo_energy – Molecular orbital energies.
t2aa – Pure alpha MP2 amplitudes.
t2ab – Mixed alpha-beta MP2 amplitudes.
t2bb – Pure beta MP2 amplitudes.
rdm2s – The 2-RDM contribution is added to this array.
- Raises:
ValueError – Invalid input.
- asf.mp2density_conventional.make_rdm2_rmp2(pt: pyscf.mp.mp2.MP2, mo_energy: ndarray) ndarray
Calculate MP(2, 1) two-particle reduced density matrix (RHF reference).
- Parameters:
pt – Restricted MP2 object.
mo_energy – MO orbital energies as a vector.
- Returns:
2-RDM in MO basis
- asf.mp2density_conventional.make_rdm2_tprod_ump2(nocca: int, noccb: int, t2aa: ndarray, t2ab: ndarray, t2bb: ndarray, rdm2s: ndarray) None
Second-order contributions to UHF-MP(2, 1) 2-RDM that are tensor products of MP2 amplitudes.
- Parameters:
nocca – Number of occupied alpha orbitals.
noccb – Number of occupied beta orbitals.
t2aa – MP2 amplitudes with alpha orbitals.
t2ab – MP2 amplitudes with mixed alpha and beta orbitals (abab).
t2bb – MP2 amplitudes with beta orbitals.
rdm2s – The 2-RDM contribution is added to this array.
- asf.mp2density_conventional.make_rdm2_ump2(pt: pyscf.mp.ump2.UMP2, mo_energy: tuple[ndarray, ndarray]) ndarray
Calculate MP(2, 1) two-particle reduced density matrix (UHF reference).
- Parameters:
pt – Unrestricted MP2 object.
mo_energy – MO orbital energies in tuple of vectors (alpha MO energies, beta MO energies).
- Returns:
- 2-RDM in MO basis as an array of dimensions 3 x N(MO) x N(MO) x N(MO) x N(MO).
rdm2s[0, p, q, r, s] = <p(alpha)+ r(alpha)+ s(alpha) q(alpha)> rdm2s[1, p, q, r, s] = <p(alpha)+ r(beta)+ s(beta) q(alpha)> rdm2s[2, p, q, r, s] = <p(beta)+ r(beta)+ s(beta) q(beta)>
This follows the convention of PySCF for reduced density matrices.
- Raises:
ValueError – Invalid input.
asf.natorbs
Functions to calculate different types of natural orbitals.
- asf.natorbs.ao_natural_orbitals(rdm1ao: ndarray, S: ndarray, Sthresh: float = 1e-08) tuple[ndarray, 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.ccsd_natural_orbitals(cc: pyscf.cc.ccsd.CCSD | pyscf.cc.uccsd.UCCSD) tuple[ndarray, 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.count_active_electrons(mo_occ: Sequence[float] | ndarray, mo_list: Sequence[int] | 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.extend_orbital_space(mo_list: Sequence[int], mo_occ: Sequence[float], Kmat: ndarray, docc_thresh: float = 1.5, unocc_thresh: float = 0.5) tuple[int, list[int]]
Extend a set of active orbitals with correlation partners based on exchange integrals.
- Parameters:
mo_list – initial list of active orbitals
mo_occ – list of molecular orbital occupation numbers
Kmat – matrix of exchange integrals, K_pq = (pq|pq)
docc_thresh – occupation numbers above this threshold are considered doubly occupied
unocc_thresh – occupation numbers below this threshold are considered unoccupied
- Returns:
new number of electrons, extended list of orbitals
- Raises:
Exception – input error
- asf.natorbs.mp2_natural_orbitals(pt: pyscf.mp.mp2.RMP2 | pyscf.mp.ump2.UMP2) tuple[ndarray, 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.natural_spin_orbitals(rdm1mo: ndarray | tuple[ndarray, ndarray], mo_coeff: ndarray | tuple[ndarray, ndarray]) tuple[ndarray, ndarray]
Calculates natural spin orbitals from a density matrix in a basis of spin-unrestricted MOs.
- Parameters:
rdm1mo – one-particle reduced density matrix in MO basis, alpha and beta parts
mo_coeff – molecular orbital coefficients for alpha and beta orbitals
- Returns:
Tuple with (natural spin occupation numbers (a, b), natural spin orbitals (a, b))
- Raises:
ValueError – Invalid input.
- asf.natorbs.rccsd_natural_orbitals(cc: pyscf.cc.ccsd.CCSD) tuple[ndarray, 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.restricted_natural_orbitals(rdm1mo: ndarray, mo_coeff: ndarray, mol: pyscf.gto.Mole | None = None) tuple[ndarray, ndarray]
Calculates natural orbitals from a density matrix 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.rmp2_natural_orbitals(pt: pyscf.mp.mp2.RMP2) tuple[ndarray, 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.select_natural_occupations(natocc: ndarray, lower: float = 0.02, upper: float = 1.98, max_orb: int | None = None, min_orb: int | None = None) tuple[int, list[int]]
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.uccsd_natural_orbitals(cc: pyscf.cc.uccsd.UCCSD) tuple[ndarray, 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.uhf_natural_orbitals(mf: pyscf.scf.uhf.UHF) tuple[ndarray, 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.ump2_natural_orbitals(pt: pyscf.mp.ump2.UMP2) tuple[ndarray, 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.unrestricted_natural_orbitals(rdm1mo: tuple[ndarray, ndarray], mo_coeff: tuple[ndarray, ndarray], S: ndarray, mol: pyscf.gto.Mole | None = None) tuple[ndarray, ndarray]
Calculates natural orbitals from a density matrix 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
asf.pairinfo
The pairinfo
module implements the PairInfoAnalyzer
which computes the metrics which form
the basis of selection process. These metrics are collected in the MOListInfo
dataclass and
the ActiveSpacesDict
type, both contained in this module.
Functions to calculate pair information between orbitals.
- class asf.pairinfo.MOListInfo(mo_list: list[int], nel: int, minimal: bool, pairinfo_sum: float, max_increment: float, min_decrement: float, min_entropy: float, min_edge_sum: float)
Information on a molecular orbital list for an active space.
- mo_list
List of molecular orbital indices.
- Type:
list[int]
- nel
Number of electrons.
- Type:
int
- minimal
Flag indicating if this is a minimal active space.
- Type:
bool
- pairinfo_sum
Sum of pair information over the MOs in mo_list and their edges.
- Type:
float
- max_increment
Maximum amount of pair information that can be added with any orbital outside of the current active space.
- Type:
float
- min_decrement
Minimum amount of pair information that can be removed with any removable orbital (e.g. within the active space, but not part of the minimum space).
- Type:
float
- min_entropy
Minimum entropy among removable MOs.
- Type:
float
- min_edge_sum
Minimum sum over all edges (active and inactive) among removable MOs.
- Type:
float
- copy_with_indices(mo_list: list[int]) MOListInfo
Create copy of instance with different set of MO indices.
- Parameters:
mo_list – List of MO indices.
- Returns:
MOListInfo with new mo_list.
- class asf.pairinfo.PairInfoAnalyzer(pairinfo: ndarray, orbdens: ndarray, spin: int | None = None)
Functionality for the analysis of pair information to determine active space suggestions.
- generate_candidate_spaces(minimal_space: list[int] | None = None, mo_list_full: list[int] | None = None) dict[tuple[int, int], list[MOListInfo]]
Construct suggestions for active spaces based on orbital pair information.
Starting from a minimal (or empty) active space, active spaces are grown by extending them with one or two orbitals at a time. Multiple lists of orbitals are created for each number of active electrons in active orbitals, (nel, norb).
A heuristic based on pair information is used to truncate the sets of spaces: for a given (N, M) combination, all MO lists must have a pair information sum exceeding a threshold. This threshold is the maximum pair information sum for an (nel, norb) space, minus the associated decrement.
- Parameters:
minimal_space – Minimal active space with orbitals that must always be included. If None is specified, a minimal space is constructed from the orbital density to include MOs that substantially affect the spin.
mo_list_full – List of indices of initial MOs in the full MO space. The indices are used to map the relative indices to the full MO space.
- Returns:
Dictionary with active space suggestions. Keys: (nel, norb) where nel and norb are the numbers of active electrons and orbitals Values: List of different active space suggestions for each active space (nel, norb).
Each list item is an MOListInfo typed dictionary.
- info_for_mospace(mo_list: Sequence[int], minimal: Sequence[int]) MOListInfo
Calculate various information for a list of active orbitals.
- Parameters:
mo_list – List of active orbital indices.
minimal – Minimal MO space with orbitals that always need to be active. These MOs are excluded when calculating minimal decrements, entropies or edge sums.
- Returns:
Dictionary with information for the MO list provided.
- minimal_space(mo_list_full: list[int] | None = None) list[int]
Construct a minimal active space based on the spin of orbitals.
The active space is chosen such that it is guaranteed to round to the (integer) number of unpaired electrons from the orbital density; provided that the full set of MOs in orbdens adds to an integer number of unpaired electrons. This is achieved by ensuring that the spin of all orbitals outside the active space sums to no more than 0.5 for positive spin and no less than -0.5 for negative spin.
- Parameters:
map_to_full_space – Whether to map active MO indices onto the full MO space. Note that ‘initial_mos’ is required for the mapping.
mo_list_full – List of indices of initial MOs in the full MO space. The indices are used to map the relative indices to the full MO space.
- Returns:
List of MO indices for a minimal active space.
- pairinfo_per_orbital() ndarray
Calculate the contribution of each orbital to the total pair information sum.
- Returns:
Array containing the pair information sum per orbital.
- pairinfo_sum(mo_list: Iterable[int]) float
Calculate the pair information sum for the active space specified.
Calculated sums are stored in this object. Results that have been calculated previously are retrieved from the buffer instead of recalculating them.
- Parameters:
mo_list – List of orbital indices for which to calculate the pair information sum.
- Returns:
The pair information sum for the specified orbital subspace.
- sanitize_active_spaces(spaces: dict[tuple[int, int], list[MOListInfo]]) dict[tuple[int, int], list[MOListInfo]]
Filter active spaces for simple sanity criteria.
The active space must not be completely filled or empty. Moreover, number of unpaired electrons must be correct.
- Parameters:
spaces – Dictionary with active space suggestions.
- Returns:
Dictionary from the input with unreasonable choices removed.
- asf.pairinfo.group_active_spaces(spaces: list[MOListInfo]) dict[tuple[int, int], list[MOListInfo]]
Group active spaces with equal number of active electrons and orbitals.
The active spaces are grouped in a dictionary where the keys are CAS tuples (N, M), i.e. N active electrons in M active orbitals, and the values are lists of MOListInfo objects sorted in descending order by the pair information sum.
- Parameters:
spaces – List of active spaces as MOListInfo objects.
- Returns:
Active spaces grouped by CAS tuple (N, M).
- asf.pairinfo.less_or_equal(val1: float, val2: float, rel_tol: float, abs_tol: float = 1e-08) bool
Performs a less-than-or-equal comparison to within a numerical tolerance.
- Parameters:
val1 – First (smaller) value.
val2 – Second (larger) value.
rel_tol – Relative comparison tolerance as defined in math.isclose.
abs_tol – Absolute comparison tolerance as defined in math.isclose.
- Returns:
True if val1 < val2 or val1 is close to val2 to within the tolerances, False otherwise.
- asf.pairinfo.map_spaces_to_full(all_spaces: dict[tuple[int, int], list[MOListInfo]], mo_list_full: list[int]) dict[tuple[int, int], list[MOListInfo]]
Copy an active space collection with MO indices mapped onto the full MO space.
In some cases ActiveSpacesDict contains MOListInfo objects that refer to relative MO indices. Given MO indices in the full MO space, this function maps the relative MO indices onto the full space.
- Parameters:
all_spaces – active space collection.
mo_list_full – list of indices of initial MOs in the full MO space. The indices are used to map the relative indices to the full MO space.
- Returns:
active space collection with MO indices mapped to full space
asf.preselection
Procedures for determining an initial active space.
- class asf.preselection.MP2NatorbPreselection(scf: pyscf.scf.hf.SCF, lower: float = 0.01, upper: float = 1.99, min_orb: int | None = 10, max_orb: int | None = 30)
Preselection strategy based on natural occupations between a lower and an upper boundary.
Given a converged Hartree-Fock solution, this procedure calculates MP2 natural orbitals (from the orbital-unrelaxed density matrix). An active space is selected by including all natural orbitals with associated eigenvalues between a lower and an upper threshold. In addition, boundaries are imposed to select no less than a minimum and no more than a maximum number of orbitals for the active space.
- lower
lower boundary for the natural occupation numbers
- upper
upper boundary for the natural occupation numbers
- min_orb
minimal number of orbitals
- max_orb
maximal number of orbitals
- scf
SCF object (RHF or UHF)
- scf_type
the SCF object’s SCF type
- compute_natorbs() tuple[ndarray, ndarray]
Compute MP2 natural orbitals.
The density-fitting MP2 implementation (RHF-DF-MP2 or UHF-DF-MP2) will be used to compute the natural orbitals.
- Returns:
natural occupation numbers, natural orbitals
- Raises:
ValueError – bad natural occupation numbers
- property mol: pyscf.gto.Mole
Molecule associated with the preselection and active space search.
- select() ActiveSpace
Perform preselection based on MP2 natural orbital occupation numbers.
- Returns:
An active space based on MP2 natural orbitals.
- class asf.preselection.MP2PairinfoPreselection(reference: pyscf.scf.uhf.UHF | pyscf.mp.dfump2_native.DFUMP2, pair_thresh: float = 0.002, svd_thresh: float = 0.98, active_mp2no: bool = True, inactive_mp2no: bool = True)
Preselection strategy based on pair information from perturbation theory.
The procedure generates an initial active space of orbitals similar to MP2 natural orbitals.
- reference
Either a UHF (a DFUMP2 object is constructed) or an existing DFUMP2 object.
- pair_thresh
Cutoff to truncate pair information.
- svd_thresh
Cutoff to construct a minimal active space from UHF corresponding orbitals.
- active_mp2no
Orbitals in minimal active space are MP2 natural orbital-like.
- inactive_mp2no
All other orbitals are MP2 natural orbital-like.
- property mol: pyscf.gto.Mole
Molecule associated with the preselection and active space search.
- select() ActiveSpace
Perform preselection based on RI-MP(2, 1) pair information.
- Returns:
An initial active space of orbitals similar to MP2 natural orbitals.
asf.rimp2_pairinfo
Calculating pair information using Moller-Plesset theory with the RI approximation.
- asf.rimp2_pairinfo.UHF_corresponding_orbitals(mo_coeff: tuple[ndarray, ndarray] | ndarray, nocc: ndarray | tuple[int, int], S: ndarray) tuple[ndarray, ndarray]
Calculate corresponding orbitals between spin-unrestricted MOs.
Occupied corresponding orbitals are defined via the SVD of the rectangular matrix of overlaps between occupied alpha and beta orbitals, S[i, j] = <i_alpha | j_beta>. Likewise, virtual corresponding orbitals are obtained through the SVD of the overlap matrix between virtual alpha and beta orbitals, S[a, b] = <a_alpha | b_beta>. Corresponding occupied and virtual orbitals of the same spin are combined into one matrix. The ordering of the MO coefficient columns are such that:
Occupied orbitals come first / left, with descending singular values.
Singular values of zero, with associated orbitals, come in the center.
Virtual orbitals come last / right, with ascending singular values.
The are at least as many singular values of zero as there are unpaired electrons.
- Parameters:
mo_coeff – Unrestricted molecular orbital coefficients.
nocc – Number of occupied orbitals (alpha, beta).
S – Overlap matrix between atomic basis functions.
- Returns:
Tuple containing (singular values, corresponding orbital coefficients)
- Raises:
ValueError – Invalid input.
- asf.rimp2_pairinfo.corresponding_orbital_subspaces(pt: pyscf.mp.dfump2_native.DFUMP2, threshold: float = 0.98, active_mp2no: bool = True, inactive_mp2no: bool = True) tuple[ndarray, ndarray, ndarray]
Performs a transformation and partitioning of subspaces for RI-MP(2, 1) pair information.
The orbital transformation is performed in two steps: 1) Unrestricted corresponding orbitals are calculated for the unrestricted occupied and virtual
spin orbitals, respectively.
Pairs of occupied corresponding spin orbitals and pairs of virtual corresponding spin orbitals are selected for an initial “active space” if the associate singular value is below the specified threshold. Therefore, the resulting spin orbitals can be divided into four subspaces: occupied “inactive” MOs, occupied “active” MOs, virtual “active” MOs and virtual “inactive” MOs.
- Optionally, the spin orbitals are transformed further (switched on by default):
“Active” spin orbitals diagonalize the unrelaxed MP2 density matrix in their subspace. Separation of occupied and virtual orbitals is preserved.
“Inactive” spin orbitals are determined as a counterpart of corresponding orbitals, but with the spin-free unrelaxed MP2 density matrix as the metric instead of the overlap. Thus, they are similar to MP2 natural orbitals, but come in pairs of spin orbitals.
- Construct a set of restricted orbitals that the unrestricted orbitals can be mapped onto:
For “inactive” orbitals, a single restricted orbital is determined to be similar to each pair of alpha and beta orbitals.
For “active” orbitals, a corresponding orbital transformation of the entire active space is performed with overlap or with the MP2 density metric (mixing occupied and virtual MOs, but not active with inactive). These are used to determine restricted orbitals.
- Parameters:
pt – DF-UMP2 instance to perform MP2 calculations.
threshold – Threshold for singular values to distinguish identify active spin orbitals.
active_mp2no – Flag to diagonalize the MP2 density matrix for active spin orbitals.
inactive_mp2no – Flag to perform SVD of the MP2 density matrix for inactive spin orbitals.
- Returns:
Tuple containing (active indices, spin orbitals, restricted orbitals) - Active indices are in a single list, which applies to both alpha and beta orbitals. - Spin orbital coefficients stored as 2 x N(AO) x N(MO) array. The columns are ordered from
left to right: occupied inactive, occupied active, virtual active, virtual inactive.
Restricted orbital coefficients stored as N(AO) x N(MO) matrix. Columns ordered from left to right: occupied inactive, active, virtual inactive.
- asf.rimp2_pairinfo.diag_cumulant_ump2(pt: pyscf.mp.dfump2_native.DFUMP2, mo_coeff: ndarray, overlap_thresh: float = 1e-10) ndarray
Calculate diagonal elements of the two-particle-cumulant with UHF-RI-MP(2, 1).
Diagonal cumulant elements are calculated in a rotated target MO basis, not necessarily in the UHF MO basis. The only requirement is that occupied and virtual spin orbitals are not mixed (unrelaxed natural spin orbitals are fine, spin-free natural orbitals in general are not).
- Parameters:
pt – DFUMP2 object containing data to perform RI-MP2 calculations.
mo_coeff – MO coefficients of the target spin orbitals (preserving occupied-virtual separation).
overlap_thresh – Threshold to determine occupied-virtual MO mixing.
- Returns:
Contributions to the “diagonal” elements [p, p, q, q] of the two-body cumulant in a target orbital basis defined by the MO coefficients in mo_coeff. The result is an array of dimensions 3 x N(MO) x N(MO), containing alpha-alpha, alpha-beta and beta-beta contributions.
- Raises:
ValueError – Invalid input.
- asf.rimp2_pairinfo.make_rdm1_ump2(pt: pyscf.mp.dfump2_native.DFUMP2) ndarray
Calculate modified RI-MP2 one-particle reduced density matrix (UHF reference).
- Parameters:
pt – Unrestricted MP2 object as a native DF-MP2 instance from PySCF.
- Returns:
1-RDM in MO basis containing the alpha and beta components.
- Raises:
ValueError – Invalid input.
- asf.rimp2_pairinfo.make_rdm2diag_0th_ump2(pt: pyscf.mp.dfump2_native.DFUMP2, mo_coeff: ndarray) ndarray
SCF (0th order) contribution to diagonal 2-RDM elements from UHF in a rotated basis.
- Parameters:
pt – DFUMP2 object
mo_coeff – MO coefficients of the target spin orbitals (not the SCF orbitals).
- Returns:
Contributions to the “diagonal” elements [p, p, q, q] of the 2-RDM in a target orbital basis defined by the MO coefficients in mo_coeff. The result is an array of dimensions 3 x N(MO) x N(MO), containing alpha-alpha, alpha-beta and beta-beta contributions.
- asf.rimp2_pairinfo.make_rdm2diag_oooo_ump2(pt: pyscf.mp.dfump2_native.DFUMP2, mo_coeff: ndarray, overlap_thresh: float = 1e-10) ndarray
All-occupied contributions to diagonal elements of the UHF-MP(2, 1) 2-RDM.
This function calculates the contributions originating from products of amplitudes. It requires that the occupied and virtual spin orbitals do not mix.
- Parameters:
pt – DFUMP2 object containing data to perform RI-MP2 calculations.
mo_coeff – MO coefficients of the target spin orbitals (preserving occupied-virtual separation).
overlap_thresh – Threshold to determine occupied-virtual MO mixing.
- Returns:
Contributions to the “diagonal” elements [p, p, q, q] of the 2-RDM in a target orbital basis defined by the MO coefficients in mo_coeff. The result is an array of dimensions 3 x N(MO) x N(MO), containing alpha-alpha, alpha-beta and beta-beta contributions.
- Raises:
MemoryError – Insufficient memory to buffer RI integrals in memory.
- asf.rimp2_pairinfo.make_rdm2diag_oovv_ump2(pt: pyscf.mp.dfump2_native.DFUMP2, mo_coeff: ndarray, overlap_thresh: float = 1e-10) ndarray
Mixed occupied-virtual contributions to diagonal elements of the UHF-MP(2, 1) 2-RDM.
This function calculates the contributions originating from products of amplitudes. It requires that the occupied and virtual spin orbitals do not mix.
- Parameters:
pt – DFUMP2 object containing data to perform RI-MP2 calculations.
mo_coeff – MO coefficients of the target spin orbitals (preserving occupied-virtual separation).
overlap_thresh – Threshold to determine occupied-virtual MO mixing.
- Returns:
Contributions to the “diagonal” elements [p, p, q, q] of the 2-RDM in a target orbital basis defined by the MO coefficients in mo_coeff. The result is an array of dimensions 3 x N(MO) x N(MO), containing alpha-alpha, alpha-beta and beta-beta contributions.
- Raises:
MemoryError – Insufficient memory to buffer RI integrals in memory.
- asf.rimp2_pairinfo.make_rdm2diag_separable_ump2(pt: pyscf.mp.dfump2_native.DFUMP2, rdm1s_mp2: ndarray, mo_coeff: ndarray) ndarray
Separable contributions to diagonal elements of the UHF-MP2(2, 1) 2-RDM.
This function calculates contributions that can be separated into products of 1-RDMs.
- Parameters:
pt – DFUMP2 object containing data to perform RI-MP2 calculations.
rdm1s_mp2 – The complete MP2 1-RDM in UHF orbital basis: can be the unrelaxed 1-RDM if there is no mixing of occupied and virtual spin orbitals. With mixing, it must be the MP(2, 1) 1-RDM, as calculated with make_rdm1_ump2(…).
mo_coeff – MO coefficients of the target spin orbitals (not the SCF orbitals).
- Returns:
Contributions to the “diagonal” elements [p, p, q, q] of the 2-RDM in a target orbital basis defined by the MO coefficients in mo_coeff. The result is an array of dimensions 3 x N(MO) x N(MO), containing alpha-alpha, alpha-beta and beta-beta contributions.
- Raises:
ValueError – Invalid input.
- asf.rimp2_pairinfo.make_rdm2diag_vvvv_ump2(pt: pyscf.mp.dfump2_native.DFUMP2, mo_coeff: ndarray, overlap_thresh: float = 1e-10) ndarray
All-virtual contributions to diagonal elements of the UHF-MP(2, 1) 2-RDM.
The function requires that the occupied and virtual spin orbitals do not mix.
- Parameters:
pt – DFUMP2 object containing data to perform RI-MP2 calculations.
mo_coeff – MO coefficients of the target spin orbitals (preserving occupied-virtual separation).
overlap_thresh – Threshold to determine occupied-virtual MO mixing.
- Returns:
Contributions to the “diagonal” elements [p, p, q, q] of the 2-RDM in a target orbital basis defined by the MO coefficients in mo_coeff. The result is an array of dimensions 3 x N(MO) x N(MO), containing alpha-alpha, alpha-beta and beta-beta contributions.
- asf.rimp2_pairinfo.overlap_blocks(pt: pyscf.mp.dfump2_native.DFUMP2, mo_coeff: ndarray, overlap_thresh: float = 1e-10) tuple[tuple[ndarray, ndarray], tuple[ndarray, ndarray]]
Obtain transformation matrices of occupied and virtual spin orbital blocks separately.
The original MOs in pt are represented in the rows of the transformation matrices, the transformed basis mo_coeff is represented in the respective columns.
- Parameters:
pt – DFUMP2 object containing the original MO coefficients.
mo_coeff – The transformed MO coefficients (must not mix occupied and virtual MOs).
overlap_thresh – Threshold to determine that MO blocks do not get mixed.
- Returns:
Transformation matrices U between the old and new MO coefficients as nested tuples. ((U for occupied alpha, U for occupied beta), (U for virtual alpha, U for virtual beta))
- Raises:
ValueError – Mixing of occupied and virtual MOs.
- asf.rimp2_pairinfo.space_with_threshold(reference: pyscf.scf.uhf.UHF | pyscf.mp.dfump2_native.DFUMP2, pair_thresh: float = 0.002, svd_thresh: float = 0.98, active_mp2no: bool = True, inactive_mp2no: bool = True) tuple[int, list[int], ndarray]
Construct a guess active space using pair information from perturbation theory.
- Parameters:
reference – Either a UHF (a DFUMP2 object is constructed) or an existing DFUMP2 object.
pair_thresh – Cutoff to truncate pair information.
svd_thresh – Cutoff to construct a minimal active space from UHF corresponding orbitals.
active_mp2no – Orbitals in minimal active space are MP2 natural orbital-like.
inactive_mp2no – All other orbitals are MP2 natural orbital-like.
- Returns:
An initial active space of orbitals similar to MP2 natural orbitals. It is provided as a tuple containing:
The number of electrons in the initial active space.
The list of MO indices for the initial active space.
Restricted MO coefficients that the MO list refers to.
- asf.rimp2_pairinfo.ur_to_re_pairwise(mo_coeff: ndarray, S: ndarray) ndarray
Constructs restricted orbitals from pairwise combination of unrestricted orbitals.
The input MO coefficients are unrestricted orbitals with a one-to-one correspondence of their spatial parts; for example, if the orbitals were constructed via a corresponding orbital transformation of their respective occupied and virtual parts. Thus, it is expected that mo_coeff[0, :, i] and mo_coeff[1, :, i] have similar spatial parts. In the returned matrix, mo_orth[:, i] will contain a spin-restricted approximation of the two spin orbitals.
- Parameters:
mo_coeff – MO coefficients of the unrestricted orbitals (alpha, beta).
S – Overlap matrix of atomic basis functions.
- Returns:
Orthonormal MO coefficients of the restricted orbitals.
- Raises:
ValueError – Bad input.
asf.scf
Definitions and functions centered around the Self-Consistent Field Method.
- asf.scf.DEFAULT_MAX_RESTARTS = 5
Number of times to restart an SCF calculation if the previous attempt failed.
- asf.scf.DEFAULT_SCF_SETTINGS = {'conv_tol': 1e-07, 'damp': 0.2, 'level_shift': 0.2, 'max_cycle': 100}
Default settings for PySCF’s mean-field classes (e.g. RHF).
- exception asf.scf.SCFError
Failed to obtain a converged and stable SCF solution.
- class asf.scf.SCFtype(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
- asf.scf.stable_scf(mol: pyscf.gto.Mole, with_uhf: bool = True, newton_backup: bool = True, max_restarts: int = 5, scf_kwargs: dict[str, Any] | None = None, basic_print: bool = True) pyscf.scf.hf.RHF | pyscf.scf.uhf.UHF
Perform Hartree-Fock calculations with stability analysis on a Mole object.
If an SCF calculation did not converge or if the solution is internally unstable, the SCF calculation is restarted.
Many settings of the SCF object can be controlled via the scf_kwargs dictionary.
Example 1: PySCF does not print output for individual SCF iterations with by default. The level of output can be controlled either via the Mole object, or via the ‘verbose’ attribute of an SCF object:
scf_kwargs = {“verbose”: pyscf.lig.logger.INFO}
Example 2: the number of SCF iterations can be changed from its default value.
scf_kwargs = {“max_cycle”: 100}
Example 3: the guess type for the SCF calculation can be modified, for example by reading the starting orbitals from a checkpoint file.
IMPORTANT: if you use the keywords below, then the checkpoint file provided will be overwritten or deleted! Only do this with a disposable copy of the file.
- scf_kwargs = {
“init_guess”: “chkfile”, “chkfile”: “my_file.chk”,
}
User documentation for SCF: https://pyscf.org/user/scf.html API documentation for RHF: https://pyscf.org/pyscf_api_docs/pyscf.scf.html#module-pyscf.scf.hf API documentation for UHF: https://pyscf.org/pyscf_api_docs/pyscf.scf.html#module-pyscf.scf.uhf
- Parameters:
mol – Instance of the PySCF Mole class.
with_uhf – Perform a UHF calculation if True, an RHF calculation if False.
newton_backup – Switch to Newton solver (more expensive) if a previous calculation failed to converge if set to True. Set to False to disable.
max_restarts – Maximal number of times to restart after unconverged or unstable solutions.
scf_kwargs – Dictionary of attributes to set on the SCF object prior to performing the calculation.
basic_print – Print information on steps initiated by this function if set to True.
- Returns:
Either a UHF or an RHF object, depending on the input.
- Raises:
SCFError – No converged and stable result was achieved.
asf.utility
Besides the general functions contained in the utility
module, pictures_Jmol
and
write_jmol_script
are especially relevant for visualizing an active space.
Utility functions.
- asf.utility.ENTROPY_ZERO_CUTOFF = 1e-12
Default value to avoid calculating the logarithm of negative numbers.
- asf.utility.calculate_entropy_s1(orbdens: ndarray, zero_cutoff: float = 1e-12) ndarray
Calculates the single-site entropies from the one-orbital density.
- Parameters:
orbdens – one-orbital densities in an N x 4 array
zero_cutoff – Threshold to determine values that are effectively zero
- Returns:
single-site entropies
- Raises:
Exception – various errors
- asf.utility.compare_active_spaces(cas: pyscf.mcscf.casci.CASCI, mo_guess: ndarray, mo_list: Sequence[int] | ndarray, full_result: bool = False) float | tuple[ndarray, ndarray, 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: pyscf.gto.Mole, mo_coeff1: ndarray, mo_coeff2: ndarray) tuple[ndarray, ndarray, 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.iround(x: float) int
A safe variant of rounding that always returns an int.
Python’s round can return float, e.g. if its argument is a numpy.float.
- Parameters:
x – Floating point number to be rounded
- Returns:
rounded number
- asf.utility.orbdens_from_rdm(rdm1a: ndarray, rdm1b: ndarray, rdm2ab: ndarray) ndarray
Calculate the one-orbital density from the reduced density matrices.
- Parameters:
rdm1a – Spin-up one-particle reduced density matrix.
rdm1b – Spin-down one-particle reduced density matrix.
rdm2ab – Mixed-spin component of the two-particle density matrix. OR: 0.5 * the spin-free two-particle density matrix.
- Returns:
empty, spin-up, spin-down, doubly occupied.
- Return type:
One-orbital density. Columns
- Raises:
ValueError – Invalid input
- asf.utility.pictures_Jmol(mol: pyscf.gto.Mole, mo_coeff: ndarray, mo_list: list[int] | ndarray, name: str = 'mo', cutoff: float = 0.04, mo_index_title: str = 'MO %I of %N', ene: ndarray | list[float] | None = None, occ: ndarray | list[float] | None = None, rotate: tuple[float, float, float] = (0.0, 0.0, 0.0), zoom: int = 60, script_options: str = '', jmol_exec: str = 'jmol', jmol_opts: str = '--nodisplay', xvfb_cmd: str | None = None, suppress_out: bool = True) 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
mo_index_title – Title format for showing orbital numbers.
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.
zoom – Zooming percentage for Jmol. Factor 100 is suitable to show all atoms, but isosurfaces may appear cropped.
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.
suppress_out – Suppress stdout and stderr output of Jmol if set to True.
- asf.utility.rdm1s_from_rdm12(N: int, S: float, rdm1: ndarray, rdm2: ndarray) tuple[ndarray, 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.show_mos_grid(images: list[Path], mo_list: list[int] | None = None, columns: int = 5, figsize: tuple[int, int] = (14, 14), image_regex: str = 'mo_(\\d+)', **fig_kwargs) None
Show molecular orbitals images in a grid.
The function assumes that the MO images follow a systematic filename convention containing the molecular orbital index. By default, the same format as in pictures_Jmol is used, i.e. mo_[idx].png, where [idx] is the respective MO index (indexing starting from 0).
- Parameters:
images – List of image paths.
mo_list – MO indices to show.
columns – Number of grid columns.
figsize – Matplotlib figure size.
image_regex – Regular expression matching the index in the orbital image filename. Custom expressions should be written, so that the index is matched in a group.
fig_kwargs – Extra arguments that are forwared to matplotlib.pyplot.subplots().
- asf.utility.tempname(directory: str = pyscf.lib.parameters.TMPDIR, suffix: str | None = 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:
directory – Directory to store the temporary file.
suffix – Suffix for the temporary file name.
- Yields:
The name of the temporary file.
- asf.utility.transform_dm(dm: ndarray, U: ndarray) ndarray
Perform an orthogonal basis transformation of a density matrix.
- Parameters:
dm – The density matrix to be transformed.
U – An orthogonal transformation matrix. First index / rows: the old basis Second index / columns: the new basis
- Returns:
The transformed density matrix.
- asf.utility.write_jmol_script(stream: TextIO, mo_filename: Path | str, mo_list: list[int] | ndarray, mo_titleformat: str = 'MO %I of %N', name: str = 'mo', cutoff: float = 0.04, rotate: tuple[float, float, float] = (0.0, 0.0, 0.0), zoom: int = 60, script_options: str = '') None
Write a Jmol script for creating pictures of orbitals to an IO stream.
Note that the indices of orbitals in mo_list start from zero, whereas Jmol counts orbitals starting from one.
- Parameters:
stream – IO stream to which the contents of the script is written.
mo_filename – Path-like filename for file containing MOs to plot
mo_list – Indices of MOs to be plotted (indexing based on zero).
mo_titleformat – A Jmol template string for MO titles.
name – File names of pictures will be: [name]_[MO index].png
cutoff – Isosurface cutoff for orbital rendering.
rotate – Rotation angle about the (x-axis, y-axis, z-axis) in degrees.
zoom – Zooming percentage for Jmol. Factor 100 is suitable to show all atoms, but isosurfaces may appear cropped.
script_options – A string with arbitrary input to be included in the Jmol script.
asf.visualization
Analyze and draw data as networks.
- asf.visualization.draw_diagonal_cumulant(rdm1a: ndarray, rdm1b: ndarray, rdm2: ndarray, show_entropy: bool = True, **kwargs) None
Draw a network representing cumulant information calculated from reduced density matrices.
- Parameters:
rdm1a – Spin-up part of the one-particle reduced density matrix.
rdm1b – Spin-down part of the one-particle reduced density matrix.
rdm2 – Spin-free two-particle reduced density matrix.
show_entropy – Use entropy to color the nodes if set to true. Otherwise, color nodes based on the cumulant information.
**kwargs – Keyword arguments for draw_pair_information(…).
- asf.visualization.draw_overlap_active(S: ndarray, mo_lists: tuple[list[int] | ndarray, list[int] | ndarray] | None = None, labels: tuple[list | ndarray, list | ndarray] | None = None, lower_thresh: float = 0.15, upper_thresh: float = 0.85, active_colors: tuple[str, str] = ('royalblue', 'crimson'), inact_colors: tuple[str, str] = ('forestgreen', 'gold'), edge_color: str = 'black', edge_width: float = 5.0, directed: bool = True, seed: int = 0, check_overlap: bool = True, ax: Axes | None = None) None
Set up the overlap of orbitals as a network, and draw it.
This function is intended to visualize if two active spaces are consistent, or whether components of MOs outside of the active space have been picked up.
- Parameters:
S – Overlap matrix in MO basis.
mo_lists – Lists of indices of MO coefficients in the active space. The first element of the tuple contains the list of MOs in the first active space (rows of S). The second element of the tuple contains the list of MOs in the second active space (columns of S).
labels – The labels for both sets of orbitals. By default, orbitals are labelled as integers starting from 0.
lower_thresh – Do not draw a graph edge of the square of the overlap is below this cutoff.
upper_thresh – Omit orbitals (nodes) if the largest square of an overlap involving this orbital with another active partner orbital is above this threshold.
active_colors – Colors to use for active orbitals.
inact_colors – Colors to use for inactive orbitals.
edge_color – Color used to draw the edges.
edge_width – Maximum width to draw the edges.
directed – Draw as directed or as undirected graph.
seed – Random seed used to determine the node positions.
check_overlap – If true, ensure that the overlap projection magnitude does not exceed 1.
ax – Draw the overlap graph in the specified Matplotlib axes.
- Raises:
ValueError – invalid input
- asf.visualization.draw_overlap_with_unrestricted(S: ndarray | tuple[ndarray, ndarray], mo_list: list[int] | ndarray | None = None, mo_offset: int = 0, lower_thresh: float = 0.1, upper_thresh: float | None = None, colors: tuple[str, str] = ('crimson', 'royalblue'), **kwargs) None
Represent the overlap between a restricted and an unrestricted orbital set graphically.
- Parameters:
S – Overlap matrices between the restricted orbital set and each set of unrestricted orbitals. Rows represent the restricted orbitals, which are identical in both matrices. Columns represent the unrestricted orbitals.
mo_list – Optionally, subset of restricted orbitals to show (S contains full set).
mo_offset – Integer offset to add to all orbital indices in the network graph.
lower_thresh – Edge only drawn if the square of the overlap exceeds this threshold.
upper_thresh – Omit edges with the square of the overlap above this threshold.
colors – Colors to use for the restricted (1st) and unrestricted nodes (2nd).
**kwargs – Further keyword arguments for the drawing functions.
- Raises:
ValueError – Bad input.
- asf.visualization.draw_pair_information(pair_info: ndarray, pair_half_val: float = 0.03, pair_exponent: float = 1.5, orbital_info: ndarray | None = None, orbital_half_val: float = 0.05, orbital_exponent: float = 1.5, orbital_labels: ndarray | list | None = None, edge_width: float = 5.0, edge_cmap: str = 'seismic', node_cmap: str = 'RdPu', seed: int = 0, verbose: bool = True, ax: Axes | None = None) None
Represents a network of pair information graphically.
- Parameters:
pair_info – The pair information, typically minus the cumulant.
pair_half_val – Magnitude at which an edge representing pair information is shown with half its maximum width and intensity.
pair_exponent – The larger the exponent, the quicker the change from no to full color saturation and width of each edge.
orbital_info – Information to color the nodes (representing orbitals). If None is provided, use the pair information sum for each orbital.
orbital_half_val – Value for which to color the nodes with half their maximal saturation.
orbital_exponent – The larger the exponent, the quicker the change from no to full color saturation of each node.
orbital_labels – Labels to use for the nodes / orbitals. By default, use their indices starting with zero.
edge_width – Maximum width to draw the edges.
edge_cmap – Color map for the edges. Assuming that the pair information can be positive and negative, a “diverging” color map should be used.
node_cmap – Color map for the nodes. Assuming that the values are non-negative, a sequential color map should be used.
seed – Random seed to determine the positions of the nodes.
verbose – Information printing.
ax – Draw the pair information graph in the specified Matplotlib axes.
- Raises:
ValueError – Invalid input
- asf.visualization.draw_pair_information_uhf(pair_info: ndarray | tuple[ndarray, ndarray, ndarray], orbital_info: tuple[ndarray, ndarray] | ndarray | None = None, orbital_indices: tuple[ndarray, ndarray] | ndarray | None = None, **kwargs) None
Draw pair information for spin-unrestricted orbitals in one graph.
- Parameters:
pair_info – Matrices with pair information for spin orbitals: (aa, ab, bb)
orbital_info – Single-orbital information for alpha and beta MOs: (alpha, beta)
orbital_indices – Integer indices of alpha and beta MOs: (alpha indices, beta indices)
**kwargs – Further keyword arguments to pass to draw_pair_information.
- Raises:
ValueError – Invalid input.
asf.wrapper
Simple wrapper functions to find active spaces (semi-)automatically.
- asf.wrapper.create_asf_switched_cisolver(mol: pyscf.gto.Mole, initial_space: ActiveSpace, spin: int = 0, nroots: int = 1, switch_dmrg: int = 12, dmrg_kwargs: dict[str, Any] | None = None, fci_kwargs: dict[str, Any] | None = None) ASFBase
Create an ASFCI or ASFDMRG object depending on the number of active orbitals.
This function sets up an ASF object that manages either a CASCI or DMRGCI calculation. By default an ASFCI object is created, but if the number of active orbitals surpasses the value of switch_dmrg, an ASFDMRG object is created instead.
- Parameters:
mol – Mole instance containing the molecular data.
initial_space – Initial space to perform search for final active space.
spin – Spin of the requested CASCI/DMRGCI calculation.
nroots – Number of CI roots (excited states, 1=groundstate) to calculate.
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.
- Returns:
ASF object
- asf.wrapper.find_from_mol(mol: pyscf.gto.Mole, max_norb: int | None = None, min_norb: int | None = None, entropy_threshold: float = 0.13862943611198905, max_scf_restarts: int = 5, scf_kwargs: dict[str, Any] | None = None, verbose: bool = True, **kwargs) ActiveSpace
Select an active space with an entropy threshold, starting from a Molecule.
This function attempts to calculate a stable UHF solution for the Mole object provided. Subsequently, find_from_scf(…) is called to calculate MP2 natural orbitals, and to select an active space from those.
For more details, please refer to find_from_scf(…).
- Parameters:
mol – Mole instance containing the molecular data.
max_norb – Maximum number of active orbitals per root.
min_norb – Minimum number of active orbitals per root.
entropy_threshold – Entropy threshold to select a single active space. Among all reasonable active spaces, a space with the lowest entropy above this threshold and with maximal pair information sum is selected as the final choice.
max_scf_restarts – Maximal number of times to restart when an unconverged or unstable solution is found.
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 find_from_scf(…).
- Returns:
Active space suggestion as provided by find_from_scf(…).
- Raises:
SCFError – SCF calculation failed to converge or is unstable.
- asf.wrapper.find_from_scf(mf: pyscf.scf.hf.SCF, max_norb: int | None = None, min_norb: int | None = None, entropy_threshold: float = 0.13862943611198905, states: int | list[tuple[int, int]] | None = None, sort_mos: bool = False, mp2_kwargs: dict[str, Any] | None = None, switch_dmrg: int = 12, dmrg_kwargs: dict[str, Any] | None = None, fci_kwargs: dict[str, Any] | None = None, verbose: bool = True, ci_nroots: int = 1) ActiveSpace
Select an active space with an entropy threshold, starting from an SCF calculation.
- 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 DMRGCI or CASCI calculation (depending on the number of active orbitals in the initial space) using the previously determined orbital subset.
Among all sensible active spaces, select a solution such that the entropies of all active orbitals are above the entropy threshold, and simultaneously the pair information sum is maximized.
- The output of the function is an ActiveSpace instance that 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 sort_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.
max_norb – Maximum number of active orbitals per root.
min_norb – Minimum number of active orbitals per root.
entropy_threshold – Entropy threshold to select a single active space. Among all reasonable active spaces, a space with the lowest entropy above this threshold and with maximal pair information sum is selected as the final choice.
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.
mp2_kwargs – Keyword arguments for MP2 orbital selection (see ‘MP2NatorbPreselection’).
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.
ci_nroots – Number of CI roots to use in CI/DMRG calculation. When several roots are requested the CI solver may occasionally converge on the wrong eigenvector. Use this option to solve for more roots to improve convergence. This option has no direct consequence for the active space selection.
- Returns:
Active space including MO indices of active orbitals and MO coefficients
- asf.wrapper.reorder_mos(mol: pyscf.gto.Mole, nel: int, mo_list: list[int], mo_coeff: ndarray) tuple[list[int], ndarray]
Arrange MO columns in the order occupied/core, active, virtual from left to right.
This function is analogous to sort_mo(…) in pyscf.mcscf.
- Parameters:
mol – Mole object.
nel – Number of electrons in the active space.
mo_list – List of active orbitals.
mo_coeff – MO coefficients. Orbitals that are not in mo_list must be ordered such that occupied orbitals come first and virtual orbitals second (from left to right).
- Returns:
Tuple (new list of active orbitals, reordered MO coefficients)
- asf.wrapper.sized_space_from_mol(mol: pyscf.gto.Mole, size: int | tuple[int, int], max_scf_restarts: int = 5, scf_kwargs: dict[str, Any] | None = None, verbose: bool = True, **kwargs) ActiveSpace
Select an active space of certain size starting from a molecule.
This function attempts to calculate a stable UHF solution for the Mole object provided. Subsequently, sized_space_from_scf(…) is called to calculate MP2 natural orbitals, and to select an active space from those.
For more details, please refer to sized_space_from_scf(…).
- Parameters:
mol – Mole instance containing the molecular data.
size – Requested number of active orbitals M, or number of active electrons N and number of active orbitals as a CAS tuple (N, M).
max_scf_restarts – Maximal number of times to restart when an unconverged or unstable solution is found.
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 sized_space_from_scf(…).
- Returns:
Active space suggestion as provided by sized_space_from_scf(…).
- Raises:
SCFError – SCF calculation failed to converge or is unstable.
- asf.wrapper.sized_space_from_scf(mf: pyscf.scf.hf.SCF, size: int | tuple[int, int], state: int | tuple[int, int] | None = None, sort_mos: bool = False, mp2_kwargs: dict[str, Any] | None = None, switch_dmrg: int = 12, dmrg_kwargs: dict[str, Any] | None = None, fci_kwargs: dict[str, Any] | None = None, verbose: bool = True, ci_nroots: int = 1) ActiveSpace
Select an active space of certain size starting from an SCF calculation.
- 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 DMRGCI or CASCI calculation (depending on the number of active orbitals in the initial space) using the previously determined orbital subset.
Select the final active space based on the pair information.
- The output of the function is an ActiveSpace instance that 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 suggestions for state-averaged calculations, but for individual (excited) electronic states only. By default an active space is suggested for the ground state in the same spin configuration as the given SCF calculation. If a value is provided for ‘state’, the function performs an active space selection for the specific spin and electronic state.
- Parameters:
mf – RHF or UHF object containing converged orbitals.
size – Requested number of active orbitals M, or number of active electrons N and number of active orbitals as a CAS tuple (N, M).
state – Electronic state to construct the active space for. None (default): Ground state for the same spin as in mf. Integer n: electronic state n (indexed from 0), same spin as in mf. Tuple (s, n): electronic state n (indexed from 0) for spin s.
sort_mos – If True, reorder the MO columns as occupied, active, virtual.
mp2_kwargs – Keyword arguments for MP2 orbital selection (see ‘MP2NatorbPreselection’).
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.
ci_nroots – Number of CI roots to use in CI/DMRG calculation. When several roots are requested the CI solver may occasionally converge on the wrong eigenvector. Use this option to solve for more roots to improve convergence. This option has no direct consequence for the active space selection.
- Returns:
Active space including MO indices of active orbitals and MO coefficients