cis
- adcc.cis(*args, **kwargs)
Run an ADC calculation.
Main entry point to run an ADC calculation. The reference to build the ADC calculation upon is supplied using the data_or_matrix argument. adcc is pretty flexible here. Possible options include:
Hartree-Fock data from a host program, e.g. a molsturm SCF state, a pyscf SCF object or any class implementing the
adcc.HartreeFockProvider
interface. From this data all objects mentioned in (b) to (d) will be implicitly created and will become available in the returned state.A
adcc.ReferenceState
objectA
adcc.LazyMp
objectA
adcc.AdcMatrix
object
- Parameters
data_or_matrix – Data containing the SCF reference
n_states (int, optional) –
kind (str, optional) –
n_singlets (int, optional) –
n_triplets (int, optional) –
n_spin_flip (int, optional) – Specify the number and kind of states to be computed. Possible values for kind are “singlet”, “triplet”, “spin_flip” and “any”, which is the default. For unrestricted references clamping spin-pure singlets/triplets is currently not possible and kind has to remain as “any”. For restricted references kind=”singlets” or kind=”triplets” may be employed to enforce a particular excited states manifold. Specifying n_singlets is equivalent to setting kind=”singlet” and n_states=5. Similarly for n_triplets and n_spin_flip. n_spin_flip is only valid for unrestricted references.
conv_tol (float, optional) – Convergence tolerance to employ in the iterative solver for obtaining the ADC vectors (default: 1e-6 or 10 * SCF tolerance, whatever is larger)
eigensolver (str, optional) – The eigensolver algorithm to use.
n_guesses (int, optional) – Total number of guesses to compute. By default only guesses derived from the singles block of the ADC matrix are employed. See n_guesses_doubles for alternatives. If no number is given here n_guesses = min(4, 2 * number of excited states to compute) or a smaller number if the number of excitation is estimated to be less than the outcome of above formula.
n_guesses_doubles (int, optional) – Number of guesses to derive from the doubles block. By default none unless n_guesses as explicitly given or automatically determined is larger than the number of singles guesses, which can be possibly found.
guesses (list, optional) – Provide the guess vectors to be employed for the ADC run. Takes preference over n_guesses and n_guesses_doubles, such that these parameters are ignored.
output (stream, optional) – Python stream to which output will be written. If None all output is disabled.
core_orbitals (int or list or tuple, optional) – The orbitals to be put into the core-occupied space. For ways to define the core orbitals see the description in
adcc.ReferenceState
. Required if core-valence separation is applied and the input data is given as data from the host program (i.e. option (a) discussed above)frozen_core (int or list or tuple, optional) – The orbitals to select as frozen core orbitals (i.e. inactive occupied orbitals for both the MP and ADC methods performed). For ways to define these see the description in
adcc.ReferenceState
.frozen_virtual (int or list or tuple, optional) – The orbitals to select as frozen virtual orbitals (i.e. inactive virtuals for both the MP and ADC methods performed). For ways to define these see the description in
adcc.ReferenceState
.environment (bool or list or dict, optional) – The keywords to specify how coupling to an environment model, e.g. PE, is treated. For details see environment.
max_subspace (int, optional) – Maximal subspace size
max_iter (int, optional) – Maximal number of iterations
- Returns
An
adcc.ExcitedStates
object containing theadcc.AdcMatrix
, theadcc.LazyMp
ground state and theadcc.ReferenceState
as well as computed eigenpairs.- Return type
Examples
Run an ADC(2) calculation on top of a pyscf RHF reference of hydrogen flouride.
>>> from pyscf import gto, scf ... mol = gto.mole.M(atom="H 0 0 0; F 0 0 1.1", basis="sto-3g") ... mf = scf.RHF(mol) ... mf.conv_tol_grad = 1e-8 ... mf.kernel() ... ... state = adcc.run_adc(mf, method="adc2", n_singlets=3)
The same thing can also be achieved using the adcc.adcN family of short-hands (see e.g.
adcc.adc2()
,adcc.cvs_adc2x()
):>>> state = adcc.adc2(mf, n_singlets=3)
Run a CVS-ADC(3) calculation of O2 with one core-occupied orbital
>>> from pyscf import gto, scf ... mol = gto.mole.M(atom="O 0 0 0; O 0 0 1.2", basis="sto-3g") ... mf = scf.RHF(mol) ... mf.conv_tol_grad = 1e-8 ... mf.kernel() ... ... state = adcc.cvs_adc3(mf, core_orbitals=1, n_singlets=3)