17.1.1.12. pysisyphus.modefollow package

17.1.1.12.1. Submodules

17.1.1.12.2. pysisyphus.modefollow.NormalMode module

class pysisyphus.modefollow.NormalMode.NormalMode(l, masses)[source]

Bases: object

See http://gaussian.com/vib/

__init__(l, masses)[source]

NormalMode class.

Cartesian displacements are normalized to 1.

Parameters:
  • l (np.array) -- Cartesian, non-mass-weighted displacements.

  • masses (np.array) -- Atomic masses.

property l_mw
mw_norm_for_norm(norm=0.01)[source]
property red_mass

17.1.1.12.3. pysisyphus.modefollow.davidson module

class pysisyphus.modefollow.davidson.DavidsonResult(cur_cycle, converged, final_modes, qs, nus, mode_inds, res_rms)

Bases: tuple

converged

Alias for field number 1

cur_cycle

Alias for field number 0

final_modes

Alias for field number 2

mode_inds

Alias for field number 5

nus

Alias for field number 4

qs

Alias for field number 3

res_rms

Alias for field number 6

pysisyphus.modefollow.davidson.block_davidson(cart_coords, masses, forces_getter, guess_modes, lowest=None, trial_step_size=0.01, hessian_precon=None, max_cycles=25, res_rms_thresh=0.0001, start_precon=5, remove_trans_rot=True, print_level=1)[source]
pysisyphus.modefollow.davidson.forces_fin_diff(forces_getter, coords, b, step_size)[source]
pysisyphus.modefollow.davidson.geom_davidson(geom, *args, **kwargs)[source]

17.1.1.12.4. pysisyphus.modefollow.lanczos module

pysisyphus.modefollow.lanczos.geom_lanczos(geom, *args, **kwargs)[source]

Wraps Lanczos algorithm for use with Geometry objects.

pysisyphus.modefollow.lanczos.lanczos(ncoords, grad_getter, dx=0.005, dl=0.01, guess=None, max_cycles=25, reortho=True, logger=None, fix_sign=True)[source]

Lanczos method to determine smallest eigenvalue & -vector.

See [1] for description of algorithm.

Parameters:
  • ncoords (int) -- Dimensionality of the problem, e.g., number of rows/columns of the Hessian.

  • grad_getter (Callable) -- Function that calculates the gradient for a given displacement. It is NOT called with full coordinates, BUT with a displacement. Calculation of the full coordinates must be handled inside the function. See 'grad_getter()' in 'geom_lanczos()' below.

  • dx (float, default: 0.005) -- Step size for finite differences calculation.

  • dl (float, default: 0.01) -- Eigenvalue convergence threshold. See eq. (8) in [1]. When abs((w_min - w_min_prev) / w_min_prev) < dl, convergence is signalled.

  • guess (Optional[ndarray], default: None) -- Initial guess vector for the lowest eigenvector.

  • max_cycles (int, default: 25) -- Maximum number of Lanczos cycles. Every cycle requires one gradient calculation.

  • reortho (bool, default: True) -- Whether reorthogonalization is carried out.

  • logger (Optional[Logger], default: None) -- Logger object.

  • fix_sign (bool, default: True) -- If true the first item of every eigenvector will be >= 0.0, e.g.

Return type:

tuple[float, ndarray, bool]

Returns:

  • w_min -- Smallest eigenvector.

  • eigenvector -- Normalized eigenvector belonging to the smallest eigenvalue.

  • converged -- Whether the Lanczos iterations converged.

17.1.1.12.5. Module contents

class pysisyphus.modefollow.NormalMode(l, masses)[source]

Bases: object

See http://gaussian.com/vib/

__init__(l, masses)[source]

NormalMode class.

Cartesian displacements are normalized to 1.

Parameters:
  • l (np.array) -- Cartesian, non-mass-weighted displacements.

  • masses (np.array) -- Atomic masses.

property l_mw
mw_norm_for_norm(norm=0.01)[source]
property red_mass
pysisyphus.modefollow.geom_davidson(geom, *args, **kwargs)[source]
pysisyphus.modefollow.geom_lanczos(geom, *args, **kwargs)[source]

Wraps Lanczos algorithm for use with Geometry objects.