3.7.1. horton.localization.localization – Orbital localization procedures

Currently supported localization flavours:
  • Pipek-Mezey

This is work in progress.

class horton.localization.localization.Localization(lf, occ_model, projector)

Bases: object

Localize canonical HF orbitals.

Arguments:

lf
A LinalgFactory instance.
occ_model
Occupation model.
Projector
Projectors for atomic basis function. A list of TwoIndex instances.

Optional arguments:

__call__(*args, **kwargs)

Localizes the orbitals using a unitary transformation to rotate the AO/MO coefficient matrix. The orbitals are optimized by minimizing an objective function.

This works only for restricted orbitals.

Arguments:

orb
The AO/MO coefficients. An Expansion instance.
select
The orbital block to be localised (str). Any of occ (occupied orbitals), virt (virtual orbitals)

Keywords:

Maxiter:

(int) maximum number of iterations for localization (default 2000)

Threshold:

(float) localization threshold for objective function (default 1e-6)

Levelshift:

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

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 0.75)
  • c1: parameter used in backtracking (float) (default 1e-4)
  • minalpha: minimum step length used in backracking (float) (default 1e-6)
  • maxiterouter: maximum number of search steps (int) (default 10)
  • maxiterinner: maximum number of optimization steps 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 (float) and scaling factor in backtracking (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)
__init__(lf, occ_model, projector)

Localize canonical HF orbitals.

Arguments:

lf
A LinalgFactory instance.
occ_model
Occupation model.
Projector
Projectors for atomic basis function. A list of TwoIndex instances.

Optional arguments:

clear()

Clear all wavefunction information

compute_population_matrix(exp)

Determine population matrix

Arguments:

exp
The current AO/MO coefficients. An Expansion instance
compute_rotation_matrix(coeff)

Determine orbital rotation matrix

Arguments:

coeff
The non-reduntant orbital rotations, we need only values for p<q
update_locblock(new)

Update localization block

lf

The LinalgFactory

locblock

The orbital block to be localized

nbasis

The number of basis functions

nocc

The number of occupied orbitals

nvirt

The number of virtual orbitals

popmatrix

The population matrix. A list of TwoIndex instances

proj

The Projectors. A list of TwoIndex instances

class horton.localization.localization.PipekMezey(lf, occ_model, projector)

Bases: horton.localization.localization.Localization

Localize canonical HF orbitals.

Arguments:

lf
A LinalgFactory instance.
occ_model
Occupation model.
Projector
Projectors for atomic basis function. A list of TwoIndex instances.

Optional arguments:

__call__(*args, **kwargs)

Localizes the orbitals using a unitary transformation to rotate the AO/MO coefficient matrix. The orbitals are optimized by minimizing an objective function.

This works only for restricted orbitals.

Arguments:

orb
The AO/MO coefficients. An Expansion instance.
select
The orbital block to be localised (str). Any of occ (occupied orbitals), virt (virtual orbitals)

Keywords:

Maxiter:

(int) maximum number of iterations for localization (default 2000)

Threshold:

(float) localization threshold for objective function (default 1e-6)

Levelshift:

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

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 0.75)
  • c1: parameter used in backtracking (float) (default 1e-4)
  • minalpha: minimum step length used in backracking (float) (default 1e-6)
  • maxiterouter: maximum number of search steps (int) (default 10)
  • maxiterinner: maximum number of optimization steps 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 (float) and scaling factor in backtracking (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)
__init__(lf, occ_model, projector)

Localize canonical HF orbitals.

Arguments:

lf
A LinalgFactory instance.
occ_model
Occupation model.
Projector
Projectors for atomic basis function. A list of TwoIndex instances.

Optional arguments:

assign_locblock()

Get localization block. A OneIndex instance

clear()

Clear all wavefunction information

compute_objective_function()

Objective function of PM localization to be maximized. The current implementation minimizes -(objective_function).

compute_population_matrix(exp)

Determine population matrix

Arguments:

exp
The current AO/MO coefficients. An Expansion instance
compute_rotation_matrix(coeff)

Determine orbital rotation matrix

Arguments:

coeff
The non-reduntant orbital rotations, we need only values for p<q
grad()

Gradient of objective function

Arguments:

popmatrix
The population matrix. A list of TwoIndex instances
hessian(lshift=1e-08)

Diagonal Hessian of objective function

Arguments:

popmatrix
The population matrix. A list of TwoIndex instances

Optinal arguments:

lshift
Shift elements of Hessian (float)
orbital_rotation_step(lshift)

Calculate gradient, hessian, and orbital rotation step

Arguments:

lshift
Shift elements of Hessian
solve_model(*args, **kwargs)

Update population matrix used to calculate objective function

update_locblock(new)

Update localization block

lf

The LinalgFactory

locblock

The orbital block to be localized

nbasis

The number of basis functions

nocc

The number of occupied orbitals

nvirt

The number of virtual orbitals

popmatrix

The population matrix. A list of TwoIndex instances

proj

The Projectors. A list of TwoIndex instances