Documantation of shin1koda/dmf

dmf.py is a Python module of the direct MaxFlux method. When you use this module in your research, please cite [1] Shin-ichi Koda and Shinji Saito, JCTC (2024). doi.

class dmf.DirectMaxFlux(ref_images, coefs=None, nsegs=4, dspl=3, remove_rotation_and_translation=True, mass_weighted=False, parallel=False, world=None, t_eval=None, w_eval=None, n_vel=None, n_trans=None, n_rot=None, eps_vel=0.01, eps_rot=0.01, beta=10.0, nmove=5, update_teval=False, params_t_update={'de': 0.15, 'dia': 1.0, 'dib': 0.2, 'epsb': 0.02, 'max_alpha0': 0.1, 'mua': 5.0, 'mub': 5.0})

Class for the direct MaxFlux method.

This class defines the variational problem of the direct MaxFlux method.

Parameters:
  • ref_images (list of ase.Atoms) – Reference images for generating an initial path. If coefs is not present, coefs of the initial path is generated by a linear interpolation with ref_images. If coefs is present, only ref_images[0] and ref_images[-1] are used as the endpoints of the path.

  • coefs (numpy.ndarray shape (nbasis,natoms,3),optional) – Coefficients of basis functions. If coefs is present, the interpolation with ref_images is skipped, and coefs is stored as it is. Default is None.

  • nsegs (int,optional) – Determines the number of B-spline functions. See Ref. 1 for details. Default is 4.

  • dspl (int,optional) – Degree of B-spline functions. Default is 3.

  • remove_rotation_and_translation (bool,optional) – Remove redundancy regarding rotaional and translational symmetry. Default is True.

  • mass_weighted (bool,optional) – Use mass-weighted coordinates. False is recommended for a rapid convergence. Default is False.

  • parallel (bool,optional) – Calculate forces of images in parallel. Both Threading and MPI are supported. Default is False. Usage is basically the same as the NEB method in ASE. See ASE’s document.

  • world (MPI world object,optional) – If not present, ase.parallel.world is used. Default is None.

  • t_eval (numpy.ndarray shape (nmove+2),optional) – Initial energy evaluation points. If not present or update_teval is True, an equally distributed t_eval is used. Default is None.

  • w_eval (numpy.ndarray shape (nmove+2),optional) – Weights of the numerical integuration of the objective function. If not present, w_eval is determined by the trapezoidal formula. Default is None.

  • n_vel (int,optional) – Determines the number of velocity constraints. See Ref. 1 for details. If not present, n_vel is 4*nsegs. Default is None.

  • n_trans (int,optional) – Determines the number of translational constraints. See Ref. 1 for details. If not present, n_trans is 2*nsegs. Default is None.

  • n_rot (int,optional) – Determines the number of rotational constraints. See Ref. 1 for details. If not present, n_rot is 2*nsegs. Default is None.

  • eps_vel (float,optional) – Parameter for loosening velocity constraints. See Ref. 1 for details. Default is 0.01.

  • eps_rot (float,optional) – Parameter for loosening rotational constraints. See Ref. 1 for details. Default is 0.01.

  • beta (float,optional) – Reciprocal temperature of the direct MaxFlux method. Default is 10 [/eV].

  • nmove (int,optional) – The number of movable energy evaluation points (images). Default is 5.

  • update_teval (bool,optional) – Update t_eval at each iteration. Default is False.

  • params_t_update (dict,optional) – Parameters for t_eval updating. See Ref. 1 for details.

images

nmove+2 images at t_eval.

Type:

list of ase.Atoms

coefs

Coefficients of B-spline functions.

Type:

numpy.ndarray shape (nbasis,natoms,3)

angs

Euler angles of images[-1], which is used when remove_rotation_and_translation is True.

Type:

numpy.ndarray shape (3)

energies

Energies of all images, which are shifted by e0, are stored after calling get_forces().

Type:

numpy.ndarray shape (nmove+2)

e0

min(E(images[0]),E(images[-1]))

Type:

float

forces

Forces of each image are stored after calling get_forces().

Type:

numpy.ndarray shape (nmove+2,natoms,3)

history

Object that stores histories of some properties. See History object below.

Type:

History object

natoms

The number of atoms in the system.

Type:

int

nbasis

The number of the B-spline functions. nbasis = nsegs + dspl.

Type:

int

ipopt_options

Non-default options of IPOPT. On initialization, {‘tol’:1.0, ‘dual_inf_tol’:0.04, ‘constr_viol_tol’:0.01, ‘compl_inf_tol’:0.01, ‘nlp_scaling_method’:’user-scaling’, ‘obj_scaling_factor’:0.1, ‘limited_memory_initialization’:’constant’, ‘limited_memory_init_val’:2.5, ‘accept_every_trial_step’:’yes’, ‘output_file’:’dmf.out’,} is set.

Type:

dict

remove_rotation_and_translation

Same as the above parameter.

Type:

bool

nsegs

Same as the above parameter.

Type:

int

dspl

Same as the above parameter.

Type:

int

t_eval

Same as the above parameter.

Type:

numpy.ndarray

w_eval

Same as the above parameter.

Type:

numpy.ndarray

n_vel

Same as the above parameter.

Type:

int

n_trans

Same as the above parameter.

Type:

int

n_rot

Same as the above parameter.

Type:

int

eps_vel

Same as the above parameter.

Type:

float

eps_rot

Same as the above parameter.

Type:

float

beta

Same as the above parameter.

Type:

float

nmove

Same as the above parameter.

Type:

int

update_teval

Same as the above parameter.

Type:

bool

params_t_update

Same as the above parameter.

Type:

dict

class History

Object storing histories of properties listed below.

forces
Type:

list of numpy.ndarray shape (nmove+2,natoms,3)

energies
Type:

list of numpy.ndarray shape (nmove+2)

coefs
Type:

list of numpy.ndarray shape (nbasis,natoms,3)

angs
Type:

list of numpy.ndarray shape (3)

t_eval
Type:

list of numpy.ndarray shape (nmove+2)

tmax

History of tmax. tmax is the highest energy point of the interpoated energy along the path. See Ref. 1 for details.

Type:

list of float

images_tmax

History of the image at t=tmax, an approximation of TS.

Type:

list of ase.Atoms

duals

History of the scaled dual infeasibility.

Type:

list of float

get_positions(t=None, P=None, nu=0)

Get all positions of images at t or their derivatives.

Parameters:
  • t (1D numpy.ndarray, optional) – All positions of images at t or their derivatives are returned. If not present, t_eval is used. Default is None.

  • P (numpy.ndarray, optional) – Return of _get_basis_values(t). See source code. Default is None. If both t and P are present, P is used in priority.

  • nu (int, optional) – Degree of derivative with respect to t. nu must be 0, 1, or 2. Default is 0.

Returns:

numpy.ndarray shape (nimages,natoms,3) – Positions or their nu-th derivative at t.

set_positions(coefs=None, angs=None)

Set positions of all images in self.images. These positions are determined by coefs and angs.

Parameters:
  • coefs (numpy.ndarray shape (nbasis,natoms,3), optional) – If not present, self.coefs is used. Default is None.

  • angs (numpy.ndarray shape (3), optional) – If not present, self.angs is used. Default is None.

get_forces()

Get forces of all images in self.images. After calling this method, forces and energies are stored in self.forces and self.energies, respectively.

Returns:

numpy.ndarray shape (nimages,natoms,3) – Forces of all images in self.images.

interpolate_energies(t_eval=None, energies=None, forces=None, coefs=None, delta_e=None)

Interpolate the enegy along the path. See Ref. 1 for details.

Parameters:
  • optional) (t_eval (numpy.ndarray shape (: ),) – Points where energies and forces were evaluated. If not present, self.t_eval is used. Default is None.

  • energies (numpy.ndarray shape (len(t_eval)), optional) – Energies evaluated at t_eval. If not present, self.energies is used. Default is None.

  • forces (numpy.ndarray shape (len(t_eval),natoms,3), optional) – Forces evaluated at t_eval. If not present, self.forces is used. Default is None.

  • coefs (numpy.ndarray shape (nbasis,natoms,3), optional) – Coefficients of B-spline functions. If not present, self.coefs is used. Default is None.

  • delta_e (list of float, optional) – If present, return points that satisfy \(\tilde{E}(t) = E_{\max} - \text{delta_e}\).

Returns:
  • polys (numpy.ndarray shape (len(t_eval)-1,4)) – Coefficients of piecewise cubic polynomial.

  • t_max (float) – Maximum point of \(\tilde{E}(t)\).

  • e_max (float) – Maximum value of \(\tilde{E}(t)\).

  • t_de (list of float) – Points that satisfy \(\tilde{E}(t) = E_{\max} - \text{delta_e}\).

get_x()

Get variables of the variational problem in a 1D array.

Returns:

x (1D numpy.ndarray) – self.coefs[1:-1].flatten(). If remove_rotation_and_translation is True, self.angs is appended.

set_x(x)

Set self.coefs and self.angs from x.

Parameters:

x (1D numpy.ndarray) – Variables of the variational problem in a 1D array.

solve(tol='tight')

Solve the variational problem.

Parameters:

tol (float or str, optional) – Change IOPOT option dual_inf_tol. If tol is float, dual_inf_tol is set to tol. If tol is either ‘tight’, ‘middle’, or ‘loose’ (keywords used in Ref. 1), dual_inf_tol is set to 0.04, 0.1, or 0.2, respectively. Default is ‘tight’.

objective(x)

Objective function.

Parameters:

x (1D numpy.ndarray) – Variables of the variational problem in a 1D array.

Returns:

float – Objective function.

gradient(x)

Gradient of the objective function.

Parameters:

x (1D numpy.ndarray) – Variables of the variational problem in a 1D array.

Returns:

numpy.ndarray shape (len(x)) – Gradient of the objective function.

constraints(x)

Constraints. :param x: Variables of the variational problem in a 1D array. :type x: 1D numpy.ndarray

Returns:

numpy.ndarray shape (# of constraints) – Constraints.

jacobian(x)

Jacobian of the constraints. :param x: Variables of the variational problem in a 1D array. :type x: 1D numpy.ndarray

Returns:

numpy.ndarray shape (# of constraints,len(x)) – Jacobian of the constraints.

intermediate(alg_mod, iter_count, obj_value, inf_pr, inf_du, mu, d_norm, regularization_size, alpha_du, alpha_pr, ls_trials)

Method called at the end of each iteration. Updating of t_eval is defined in this method.

Indices and tables