3.2.3. horton.correlatedwfn.restricted_ap1rog – Correlated wavefunction implementations

This module contains restricted AP1roG

Variables used in this module:
nocc:number of (active) occupied orbitals in the principle configuration
pairs:number of electron pairs
nvirt:number of (active) virtual orbitals in the principle configuration
nbasis:total number of basis functions
Indexing convention:
i,j,k,..:occupied orbitals of principle configuration
a,b,c,..:virtual orbitals of principle configuration
p,q,r,..:general indices (occupied, virtual)
Abbreviations used (if not mentioned in doc-strings):
L_pqrs:2<pq|rs>-<pq|sr>
g_pqrs:<pq|rs>

For more information see doc-strings.

class horton.correlatedwfn.restricted_ap1rog.RAp1rog(lf, occ_model, npairs=None, nvirt=None)

Bases: horton.correlatedwfn.geminal.Geminal

Arguments:

lf
A LinalgFactory instance.
occ_model
Occupation model

Optional arguments:

npairs
Number of electron pairs, if not specified, npairs = number of occupied orbitals
nvirt
Number of virtual orbitals, if not specified, nvirt = (nbasis-npairs)
__call__(one, two, core, orb, olp, scf, **kwargs)

Optimize geminal coefficients and—if required—find optimal set of orbitals.

Arguments:

one, two
One- and two-body integrals (some Hamiltonian matrix elements).
core
The core energy (not included in ‘one’ and ‘two’).
orb
An expansion instance. It contains the MO coefficients (orbitals).
olp
The AO overlap matrix. A TwoIndex instance.
scf
A boolean. If True: Initializes orbital optimization.
Keywords:
See RAp1rog.solve() and RAp1rog.solve_scf()
__init__(lf, occ_model, npairs=None, nvirt=None)

Arguments:

lf
A LinalgFactory instance.
occ_model
Occupation model

Optional arguments:

npairs
Number of electron pairs, if not specified, npairs = number of occupied orbitals
nvirt
Number of virtual orbitals, if not specified, nvirt = (nbasis-npairs)
check_convergence(e0, e1, gradient, thresh)

Check convergence.

Arguements:

e0, e1
Used to calculate energy difference e0-e1
gradient
The gradient, a OneIndex instance
thresh
Dictionary containing threshold parameters (‘energy’, ‘gradientmax’, ‘gradientnorm’)
Returns:
True if energy difference, norm of orbital gradient, largest element of orbital gradient are smaller than some threshold values.
check_keywords(guess, solver, maxiter, dumpci, thresh, printoptions)

Check dictionaries if they contain proper keys.

Arguments:

guess, solver, maxiter, dumpci, thresh, printoptions
See RAp1rog.solve()
check_keywords_scf(guess, solver, maxiter, dumpci, thresh, printoptions, stepsearch)

Check dictionaries if they contain proper keys.

Arguments:

guess, solver, maxiter, dumpci, thresh, printoptions, stepseach
See RAp1rog.solve_scf()
check_stepsearch(linesearch)

Check trustradius. Abort calculation if trustradius is smaller than 1e-8

clear()

Clear all wavefunction information

clear_auxmatrix()

Clear auxiliary matrices

clear_dm()

Clear RDM information

clear_geminal()

Clear geminal information

clear_lagrange()

Clear lagrange information

compute_1dm(dmout, mat1, mat2, factor=1.0, response=True)

Compute 1-RDM for AP1roG

Arguments:

mat1, mat2
A DenseTwoIndex instance used to calculated 1-RDM. For response RDM, mat1 is the geminal matrix, mat2 the Lagrange multipliers

Optional arguments:

factor
A scalar factor
select
Switch between response (True) and PS2 (False) 1-RDM
compute_2dm(dmout, dm1, mat1, mat2, select, response=True)

Compute response 2-RDM for AP1roG

** Arguments **

one_dm
A 1DM. A OneIndex instance.
mat1, mat2
TwoIndex instances used to calculate the 2-RDM. To get the response DM, mat1 is the geminal coefficient matrix, mat2 are the Lagrange multipliers
select
Either ‘ppqq’ or ‘pqpq’. Note that the elements (iiii), i.e., the 1DM, are stored in pqpq, while the elements (aaaa) are stored in ppqq.
response
If True, calculate response 2-RDM. Otherwise the PS2 2-RDM is calculated.
compute_correlation_energy(arg=None)

Get correlation energy of restricted AP1roG

Optional arguments:

arg
The geminal coefficient matrix (np.array or TwoIndex instance). If not provided, the correlation energy is calculated from self.geminal (default None)
compute_objective_function(coeff=None)

Objective function for line search optimization

Optional arguments:

coeff
See RAp1rog.compute_total_energy()
compute_orbital_gradient(*args, **kwargs)

Determine orbital gradient for all non-reduntant orbital rotations (oo,vo,vv).

Optional arguments:

ps2c
(boolean) If True, switches to PS2c orbital optimization (default False)
compute_orbital_hessian(lshift=1e-08, pos=False, ps2c=False)

Construct diagonal approximation to Hessian for orbital optimization of AP1roG.

Optional arguments:

lshift
Level shift (float) (default 1e-8)
pos
Set all elements positive (boolean) (default False)
ps2c
(boolean) If True, switches to PS2c orbital optimization (default False)
compute_reference_energy()

Get energy of reference determinant for restricted AP1roG

compute_rotation_matrix(coeff)

Determine orbital rotation matrix for (oo), (vo), and (vv) blocks

Arguments:

coeff
The nonreduntant orbital rotations k_pq (1-dim np.array). The elements are sorted w.r.t. lower triangular indices (p>q).
compute_total_energy(coeff=None)

Get total energy (reference+correlation) including nuclear-repulsion/ core energy for restricted AP1roG.

Optional arguments:

arg
The geminal coefficient matrix (np.array or TwoIndex instance). If not provided, the correlation energy is calculated from self.geminal (default None)
do_checkpoint(orb, olp, checkpoint_fn)

Dump orbitals for restart

Arguments:

orb
An expansion instance. AO/MO coefficients stored to disk
olp
A TwoIndex instance. The AO overlap matrix. Required for restart at different molecular geometry.
checkpoint_fn
The filename of the checkpoint file.
dump_final(orb, olp, printoptions, dumpci, checkpoint, checkpoint_fn=’checkpoint.h5’)

Dump final solution (orbitals, wavefunction amplitudes, geminal coefficients)

Arguments:

orb
An expansion instance. AO/MO coefficients stored to disk
olp
A TwoIndex instance. The AO overlap matrix. Required for restart
printoptions, dumpci
A dictionary. See RAp1rog.solve_scf()
checkpoint:
An integer. See RAp1rog.solve_scf()
checkpoint_fn
The filename of the checkpoint file.
generate_guess(guess, dim=None)

Generate a guess of type ‘guess’.

Arguments:

guess
A dictionary, containing the type of guess.

Optional arguments:

dim
Length of guess.
get_auxmatrix(select)

Get an auxiliary matrix from cache.

Arguments:

select
‘t’, ‘gpqpq’, ‘gpqqp’, ‘lpqpq’, ‘pqrq’, ‘gpqrr’, ‘fock’.
get_exact_hessian(*args, **kwargs)

Construct exact Hessian for orbital optimization of restricted OAP1roG.

The usage of the exact orbital Hessian for the orbital optimization is currently not supported. The exact Hessian can only be evaluated a posteriori.

Arguements

mo1, mo2
The 1- and 2-el integrals in the MO basis (TwoIndex and FourIndex instances)
get_four_dm(select)

Get a density matrix (4-RDM). If not available, it will be created (if possible)

Arguments:

select

get_one_dm(select)

Get a density matrix (1-RDM). If not available, it will be created (if possible)

Arguments:

select
‘ps2’, or ‘response’.
get_slater_determinant(*indices)

Return excited Slater Determinant.

Arguments:

indices
A (list of) multi-index. First element contains occupied index w.r.t. reference determinant, second element contains virtual index w.r.t. reference determinant.
get_three_dm(select)

Get a density matrix (3-RDM). If not available, it will be created (if possible)

Arguments:

select

get_two_dm(select)

Get a density matrix (2-RDM). If not available, it will be created (if possible)

Arguments:

select
‘(r(esponse))ppqq’, or ‘(r(esponse))pqpq’.
init_auxmatrix(select)

Initialize auxiliary matrices

Arguments:

select
One of t, fock, gppqq, gpqpq, gpqqp, lpqpq, lpqrq, gpqrr
init_four_dm(select)

Initialize 4-RDM

Arguments

select

init_one_dm(select)

Initialize 1-RDM as OneIndex object

The 1-RDM expressed in the natural orbital basis is diagonal and only the diagonal elements are stored.

Arguments

select
‘ps2’ or ‘response’.
init_three_dm(select)

Initialize 3-RDM

Arguments

select

init_two_dm(select)

Initialize 2-RDM as TwoIndex object

Only the symmetry-unique elements of the (response) 2-RDM are stored. These are matrix elements of type

\[Gamma_{p\bar{q}p\bar{q}}\]

(spin-up and spin-down (bar-sign)) or

\[Gamma_{p\bar{p}q\bar{q}}\]

and are stored as elements \({pq}\) of two_dm_pqpq, and two_dm_ppqq.

Arguments

select
‘(r(esponse))ppqq’, or ‘(r(esponse))pqpq’.
jacobian_ap1rog(coeff, miiaa, miaia, one, fock)

Construct Jacobian for optimization of geminal coefficients.

jacobian_lambda(lagrange, gmat, miiaa, miaia, one, diagfock)

Construct Jacobian for optimization of Lagrange multipliers for restricted AP1roG.

orbital_rotation_step(lshift=1e-08, pos=True, optimizer=’variational’)

Get orbital rotation step (Newton–Raphson)

Arguments:

Optional arguments:

lshift
Level shift (float) for approximate Hessian added to small elements (default 1e-8)
pos
(boolean) Make all elements of Hessian positive if set to True (default True)
optimizer
Orbital otimization method (str) (default ‘variational’)
perm(a)

Calculate the permament of a matrix

Arguements

a
A np array
print_final(s=’Final’)

Print energies

Optional arguments:

s
A string.
print_geminal_coefficients()

Print geminal coefficients

print_options(guess, solver, maxiter, thresh, printoptions, indextrans)

Print optimization options.

Arguments:

guess, solver, maxiter, thresh, printoptions, indextrans
See RAp1rog.solve()
print_options_scf(guess, solver, maxiter, lshift, stepsearch, thresh, printoptions, checkpoint, checkpoint_fn, indextrans, orbitaloptimizer, sort)

Print optimization options.

Arguments:

See RAp1rog.solve_scf()

print_solution(cithresh=0.01, excitationlevel=2, amplitudestofile=False, filename=’./ap1rog_amplitudes’)

Print coefficients \({ci}\) and Slater Determinant if \(|ci|\) > threshold. Prints up to hextuply excited pairs.

Optional arguments:

cithresh
Upper threshold (a float) for CI coefficients to be reconstructed. (Default 0.01)
excitationslevel
The maximum number of substitutions w.r.t. the reference determinant for which the CI coefficients are reconstructed. (Default 2)
amplitudestofile
A boolean. If True, the CI coefficients are stored in a separate file. (Default False)
filename
The file name for the wfn amplitudes. (Default ap1rog_amplitudes)
prod(lst)
solve(*args, **kwargs)

Optimize AP1roG coefficients for some Hamiltonian. For restricted, closed-shell AP1roG.

Arguments:

one, two
One- (TwoIndex instance) and two-body integrals (FourIndex or Cholesky instance) (some Hamiltonian matrix elements)
core
The core energy (not included in ‘one’ and ‘two’); a float.
orb
An expansion instance. It contains the MO coefficients.
olp
The AO overlap matrix. A TwoIndex instance.

Keywords:

indextrans:

4-index Transformation (str). Choice between tensordot (default) and einsum

warning:

Print warnings (boolean) (default False)

guess:

initial guess (dictionary) containing:

  • type: guess type (str). One of random (random numbers, default), const (a constant scaled by a factor)
  • factor: a scaling factor (float) (default -0.1)
  • geminal: external guess for geminal coefficients (1-d np.array); if provided, ‘type’ is ignored (default None)
solver:

wfn solver (dictionary) containing:

  • wfn: wavefunction solver (str) (default krylov)
maxiter:

max number of iterations (dictionary) containing:

  • wfniter: maximum number of iterations (int) for wfn solver (default 200)
dumpci:

dump ci coefficients (dictionary):

  • amplitudestofile: write wfn amplitudes to file (boolean) (default False)
  • amplitudesfilename: (str) (default ap1rog_amplitudes.dat)
thresh:

optimization thresholds (dictionary) containing:

  • wfn: threshold for geminal coefficients (float) (default 1e-12)
printoptions:

print level (dictionary) containing:

  • geminal: (boolean), if True, geminal matrix is printed (default True)
  • ci: threshold for CI coefficients (float) (requires evaluation of a permanent), all coefficients (for a given excitation order) larger than ‘ci’ are printed (default 0.01)
  • excitationlevel: number of excited pairs w.r.t. reference determinant for which wfn amplitudes are reconstructed (int) (default 1)
swapa:

swap orbitals (numpy 2-dim array), each row contains 2 orbital indices to be swapped (default np.array([[]]))

solve_geminal(*args, **kwargs)

Solves for geminal matrix

Arguments:

guess
The initial guess. A 1-dim np array.
solver
The solver used. A dictionary with element ‘wfn’.
wfnthreshold
The optimization threshold. A float.
wfnmaxiter
The maximum number of iterations. An integer.
solve_lagrange(*args, **kwargs)

Solves for Lagrange multipliers

Arguments:

guess
The initial guess. A 1-dim np array.
solver
The solver used. A dictionary with element ‘lagrange’.
wfnthreshold
The optimization threshold. A float.
wfnmaxiter
The maximum number of iterations. An integer.
solve_model(one, two, orb, **kwargs)

Solve for geminal model.

Arguments:

one, two
One- and two-body integrals (some Hamiltonian matrix elements). A TwoIndex and FourIndex/Cholesky instance
orb
An expansion instance which contains the MO coefficients.
Keywords:

guess: initial guess for wfn (1-dim np.array) guesslm: initial guess Lagrange multipliers (1-dim np.array) solver: wfn/Lagrange solver (dictionary) indextrans: 4-index Transformation (str). maxiter: maximum number of iterations (dictionary) thresh: thresholds (dictionary) orbitaloptimizer: orbital optimization method (str)

For more details, see RAp1rog.solve_scf()

solve_scf(*args, **kwargs)

Find Geminal expansion coefficient for some Hamiltonian. For restricted, closed-shell AP1roG. Perform orbital optimization.

Arguments:

one, two
One- (TwoIndex instance) and two-body integrals (FourIndex or Cholesky instance) (some Hamiltonian matrix elements)
core
The core energy (not included in ‘one’ and ‘two’).
orb
An expansion instance. It contains the MO coefficients.
olp
The AO overlap matrix. A TwoIndex instance.

Keywords:

indextrans:

4-index Transformation (str). Choice between tensordot (default) and einsum

warning:

print warnings (boolean) (default False)

guess:

initial guess (dictionary) containing:

  • type: guess type (str). One of random (random numbers, default), const (a constant scaled by a factor)
  • factor: a scaling factor (float) (default -0.1)
  • geminal: external guess for geminal coefficients (1-dim np.array) (default None)
  • lagrange: external guess for Lagrange multipliers (1-dim np.array) (default None)
solver:

wfn/Lagrange solver (dictionary) containing:

  • wfn: wavefunction solver (str) (default krylov)
  • lagrange: Lagrange multiplier solver (str) (default krylov)
maxiter:

max number of iterations (dictionary) containing:

  • wfniter: maximum number of iterations for
    wfn/lagrange solver (int) (default 200)
  • orbiter: maximum number of orbital optimization
    steps (int) (default 100)
dumpci:

dump ci coefficient (dictionary) containing:

  • amplitudestofile: write wfn amplitudes to file (boolean) (default False)
  • amplitudesfilename: (str) (default ap1rog_amplitudes.dat)
thresh:

optimization thresholds (dictionary) containing:

  • wfn: threshold for geminal coefficients and Lagrange multipliers (float) (defaulf 1e-12)
  • energy: convergence threshold for energy (float) (default 1e-8)
  • gradientnorm: convergence threshold for norm of orbital gradient (float) (default 1e-4)
  • gradientmax: threshold for maximum absolute value of orbital gradient (float) (default 5e-5)
printoptions:

print level; dictionary containing:

  • geminal: (boolean) if True, geminal matrix is printed (default True)
  • ci: threshold for CI coefficients (requires evaluation of a permanent) (float). All coefficients (for a given excitation order) larger than ci are printed (default 0.01)
  • excitationlevel: number of excited pairs w.r.t. reference determinant for which wfn amplitudes are reconstructed (int) (default 1)
stepsearch:

step search options (dictionary) containing:

  • method: step search method used (str). One of trust-region (default), None, backtracking
  • alpha: scaling factor for Newton step (float), used in backtracking and None method (default 1.00)
  • c1: parameter used in backtracking (float) (default 1e-4)
  • minalpha: minimum scaling factor used in backracking (float) (default 1e-6)
  • maxiterouter: maximum number of search steps (int) (default 10)
  • maxiterinner: maximum number of optimization step in each search step (int) (used only in pcg, default 500)
  • maxeta: upper bound for estimated vs actual change in trust-region (float) (default 0.75)
  • mineta: lower bound for estimated vs actual change in trust-region (float) (default 0.25)
  • upscale: scaling factor to increase trustradius in trust-region (float) (default 2.0)
  • downscale: scaling factor to decrease trustradius in trust-region and step length in backtracking (float) (default 0.25)
  • trustradius: initial trustradius (float) (default 0.75)
  • maxtrustradius: maximum trustradius (float) (default 0.75)
  • threshold: trust-region optimization threshold, only used in pcg (float) (default 1e-8)
  • optimizer: optimizes step to boundary of trustradius (str). One of pcg, dogleg, ddl (default ddl)
checkpoint:

frequency of checkpointing (int). If > 0, writes orbitals and overlap to a checkpont file (defatul 1)

checkpoint_fn:

filename to use for the checkpoint file (default “checkpoint.h5”)

levelshift:

level shift of Hessian (float) (default 1e-8)

absolute:

(boolean), if True, take absolute values of Hessian (default False)

sort:

(boolean), if True, orbitals are sorted according to their natural occupation numbers. This requires us to solve for the wavefunction again. Works only if orbitaloptimizer is set to variational. (default True)

swapa:

swap orbitals (numpy 2-dim array) each row contains 2 orbital indices to be swapped (default np.array([[]]))

givensrot:

rotate two orbitals (numpy 2-dim array) each row contains 2 orbital indices and rotation angle (default np.array([[]]))

orbitaloptimizer:
 

(str) switch between variational orbital optimization (variational) and PS2c (ps2c) (default variational).

sort_natural_orbitals(orb)

Sort orbitals w.r.t. the natural occupation numbers

Arguments:

orb
The AO/MO coefficients (natural orbitals) and the natural occupation numbers (Expansion instance)
start_up(orb, swapa, givensrot=array([], shape=(1, 0), dtype=float64))

Modify orbitals prior to optimization

Arguments:

swapa
An integer numpy array with two columns where every row corresponds to one swap.

Optional arguments:

givensrot
Givens rotation of orbital pairs. A numpy array where every row corresponds to one Givens rotation with elements (index1, index2, angle in deg).
update_auxmatrix(*args, **kwargs)

Derive all auxiliary matrices. gppqq: <pp|qq>, gpqpq: <pq|pq>, gpqqp: <pq|qp>, lpqpq: 2<pq|pq>-<pq|qp>, lpqrq: 2<pq|rq>-<pq|qr>, gpqrr: <pq|rr>, fock: h_pp + sum_i(2<pi|pi>-<pi|ip>)

Arguments:

select
scf.
two_mo
two-electron integrals in the MO basis. A FourIndex instance.
one_mo
one-electron integrals in the MO basis. A TwoIndex instance.
update_ecore(new)

Update core energy

update_four_dm(select, four_dm=None)

Update 4-RDM

Optional arguments:

four_dm
When provided, this 4-RDM is stored.
update_geminal(geminal=None)

Update geminal matrix

Optional arguments:

geminal
When provided, this geminal matrix is stored.
update_lagrange(lagrange=None, dim1=None, dim2=None)

Update Lagragne multipliers

Optional arguments:

lagrange
When provided, this set of Lagrange multipliers is stored.
update_one_dm(select, one_dm=None)

Update 1-RDM

Arguments:

select
One of ps2, response

Optional arguments:

one_dm
When provided, this 1-RDM is stored. A OneIndex instance.
update_three_dm(select, three_dm=None)

Update 3-RDM

Optional arguments:

three_dm
When provided, this 3-RDM is stored.
update_two_dm(select, two_dm=None)

Update 2-RDM

Arguments:

select
One of ps2, response

Optional arguments:

two_dm
When provided, this 2-RDM is stored.
vector_function_geminal(*args, **kwargs)

Construct vector function for optimization of geminal coefficients.

Arguments:

coeff
The geminal coefficients (np.array)
miiaa, miaia
Auxiliary matrices (TwoIndex instances)
one
1-electron integrals (TwoIndex instance)
diagfock
Diagonal inactive Fock matrix (OneIndex instances)
vector_function_lagrange(*args, **kwargs)

Construct vector function for optimization of Lagrange multipliers.

Arguments:

lagrange
The lagrange multipliers (np array)
gmat
The geminal coefficients (TwoIndex instance)
miiaa, miaia
Auxiliary matrices (TwoIndex instances)
one
1-electron integrals (TwoIndex instance)
diagfock
Diagonal inactive Fock matrix (OneIndex instance)
dimension

The number of unknowns (i.e. the number of geminal coefficients)

ecore

The core energy

geminal

The geminal coefficients

lagrange

The Lagrange multipliers

lf

The LinalgFactory instance

nbasis

The number of basis functions

nocc

The number of occupied orbitals

npairs

The number of electron pairs

nvirt

The number of virtual orbitals

one_dm_ps2

Alpha 1-RDM

one_dm_response

Alpha 1-RDM

two_dm_ppqq

Alpha-beta PS2 (ppqq) 2-RDM

two_dm_pqpq

Alpha-beta PS2 (pqpq) 2-RDM

two_dm_rppqq

Alpha-beta (ppqq) 2-RDM

two_dm_rpqpq

Alpha-beta (pqpq) 2-RDM