Options#
Introduction#
Optking uses pydantic to validate user input for types and values. Brief descriptions of each option is presented below. All options have a reasonable default. In general tightening or loosening convergence, adding constraints, and changing optimization type should be the only options that need to be adjusted.
Option Short Cuts#
Some of the more important options for general use of OptKing. See below for the full list.
Global and Local Options#
For developers or users interacting with optking’s internals, Optking historically has utilized a single global options dictionary. Some parts of Optking still utilize this global dictionary Most of the current code’s classes and functions; however, accept local options by passing a dictionary, options object, or explicit parameters. Note - the options reported here are the names that the user should use when providing keywords to Optking. For developers, Optking may use different names internally. For starters, variables (with the exception of matrices) should adhere to PEP8 (lower snake-case).
OptParams#
Alphabetized Keywords#
- pydantic model optking.v2.optparams.OptParams[source]#
- Config:
alias_generator: function = <function OptParams.<lambda> at 0x747c6d1d8860>
str_to_upper: bool = True
- field ACCEPT_SYMMETRY_BREAKING: bool = False (name 'accept_symmetry_breaking')#
Whether to accept steps that lower the molecular point group. Within optking this check is not rigorous and if the only reasonable step is symmetry breaking it will be taken. This keyword affects optking’s symmetrization not Psi4’s
- field ADD_AUXILIARY_BONDS: bool = False (name 'add_auxiliary_bonds')#
Do add bond coordinates at nearby atoms for non-bonded systems?
- field AUXILIARY_BOND_FACTOR: float = 2.5 (name 'auxiliary_bond_factor')#
This factor times the standard covalent distance is used to add extra stretch coordinates.
- Constraints:
gt = 1.0
- field BT_DX_CONV: float = 1e-07 (name 'bt_dx_conv')#
Threshold for the change in any given Cartesian coordinate during iterative back-transformation.
- Constraints:
gt = 0.0
- field BT_DX_RMS_CHANGE_CONV: float = 1e-12 (name 'bt_dx_rms_change_conv')#
Threshold for RMS change in Cartesian coordinates during iterative back-transformation.
- Constraints:
gt = 0.0
- field BT_MAX_ITER: int = 25 (name 'bt_max_iter')#
Maximum number of iterations allowed to converge back-transformation
- Constraints:
gt = 0
- field BT_PINV_RCOND: float = 1e-06 (name 'bt_pinv_rcond')#
Threshold to remove redundancies from generalized inverse. Corresponds to the rcond from numpy The following should be used whenever redundancies in the coordinates are removed, in particular when forces and Hessian are projected and in back-transformation from \(\delta(q)\) to \(\delta(x)\).
- Constraints:
gt = 0.0
- field CART_HESS_READ: bool = False (name 'cart_hess_read')#
Do read Cartesian Hessian? Recommended to use
full_hess_everyinstead. cfour format or.jsonfile (AtomicOutput) allowed. The filetype is determined by the presence of a.jsonextension. The cfour hessian format specifies that the first line contains the number of atoms. Each subsequent line contains three hessian values provided in row-major order. see psi4 docs for details on [cfour format]
- field CONJUGATE_GRADIENT_TYPE: str = 'FLETCHER' (name 'conjugate_gradient_type')#
One of
["POLAK", "FLETCHER", "DESCENT"]. Changes how the step direction is calculated.- Constraints:
pattern = re.compile(‘^(?:FLETCHER|DESCENT|POLAK)$’, re.IGNORECASE)
- field CONSECUTIVE_BACKSTEPS: int = 0 (name 'consecutive_backsteps_allowed')#
Sets the number of consecutive backward steps allowed in an optimization. This option can be updated by
Optkingifdynamic_lvlis > 0. Not recommended for general use.- Constraints:
ge = 0
- field MAX_ENERGY_G_CONVERGENCE: float = 1e-06 (name 'conv_max_DE')#
Convergence criterion for geometry optimization: maximum energy change
- field MAX_DISP_G_CONVERGENCE: float = 0.0012 (name 'conv_max_disp')#
Convergence criterion for geometry optimization: maximum displacement (internal coordinates, au)
- field MAX_FORCE_G_CONVERGENCE: float = 0.0003 (name 'conv_max_force')#
Convergence criterion for geometry optimization: maximum force (internal coordinates, au)
- field RMS_DISP_G_CONVERGENCE: float = 0.0012 (name 'conv_rms_disp')#
Convergence criterion for geometry optimization: rms displacement (internal coordinates, au)
- field RMS_FORCE_G_CONVERGENCE: float = 0.0003 (name 'conv_rms_force')#
Convergence criterion for geometry optimization: maximum force (internal coordinates, au)
- field COVALENT_CONNECT: float = 1.3 (name 'covalent_connect')#
When determining connectivity, a bond is assigned if interatomic distance is less than (this number) * sum of covalent radii. When connecting disparate fragments when frag_mode = SINGLE, a “bond” is assigned if interatomic distance is less than (this number) * sum of covalent radii. The value is then increased until all the fragments are connected directly or indirectly.
- Constraints:
gt = 0.0
- field DYNAMIC_LVL: int = 0 (name 'dynamic_level')#
An integer between 0 and 6. Larger values reflect less aggressive optimization techniques If
dynamic_lvlis not set,Optkingwill not change thedynamic_lvl. Thedynamic_lvlmust be > 0 for alternative approaches to be tried. A backstep will be triggered (if allowed) by \(\Delta E > 0\) in a minimization. A step is considered “bad” if \(\Delta E > 0\) when no more backsteps are allowed and iterations \(> 5\), or there are badly defined internal coordinates or derivatives. Default = 0dynamic
step
coord
trust
backsteps
criteria to change dynamic_lvl
run_level
decrease
increase
0
RFO
RI
dynamic
no
none
none
1
RFO
RI
dynamic(d)
no
1 bad step
none
2
RFO
RI
smaller
yes (1)
1 bad step
none
3
RFO
BOTH
small
yes (1)
1 bad step
none
4
RFO
XYZ
large
yes (1)
1 bad step
none
5
RFO
XYZ
small
yes (1)
1 bad step
none
6
SD
XYZ
large
yes (1)
1 bad step
none
7
SD
XYZ
small
yes (1)
1 bad step
none
- Constraints:
ge = 0
le = 6
- field DYNAMIC_LVL_MAX: int = 0 (name 'dynamic_lvl_max')#
How large
dynamic_lvlis allowed to grow. Ifdynamic_lvl\(> 0\),dynamic_lvl_maxwill default to 6- Constraints:
ge = 0
le = 6
- field ENSURE_BT_CONVERGENCE: bool = False (name 'ensure_bt_convergence')#
Reduces step size as necessary to ensure convergence of back-transformation of internal coordinate step to Cartesian coordinates.
- field EXT_FORCE_BEND: str = '' (name 'ext_force_bend')#
A string of white-space separated atomic indices (3) followed by a single variable equation surrounded in either a single or double quotation mark. Example:
"1 2 3 'Sin(x)'"evaluates the force along the coordinate as a 1-D sinusoidal function where x is the “value” (angle [radians]) of the coordinate (bend)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+s*W*s*s*[‘”].*[‘”]W?)*$
- field EXT_FORCE_CARTESIAN: str = '' (name 'ext_force_cartesian')#
A string of whitecaps separated atomic indices (1) and Cartesian labels, followed by a single variable equation surrounded in either a single or double quotation mark. Example:
"1 X 'Sin(x)'"evaluates the force along the coordinate as 1 1-D sinusoidal function where x is the “value” (angle [bohr]) of the coordinate (bohr)- Constraints:
pattern = (?i)^s*(?:W?d+s+(?:xyz|xy|yz|x|y|z)W?s+[’"].*[’"]W?)*$
- field EXT_FORCE_DIHEDRAL: str = '' (name 'ext_force_dihedral')#
A string of white-space separated atomic indices (4) followed by a single variable equation surrounded in either a single or double quotation mark. Example:
"1 2 3 4 'Sin(x)'"evaluates the force along the coordinate as a 1-D sinusoidal function where x is the “value” (angle [radians]) of the coordinate (torsion)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+W?s+d+s*W*s*s*[‘”].*[‘”]W?)*$
- field EXT_FORCE_DISTANCE: str = '' (name 'ext_force_distance')#
A string of white-space separated, atomic indices (2) followed by a single variable equation surrounded in either a single or double quotation mark. Example:
"1 2 'Sin(x)'"or'1 2 "Sin(x)"'evaluates the force along the coordinate as a 1-dimensional sinusoidal function where x is the “value” (distance [bohr]) of the coordinate (stretch).- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+s*W*s*s*[‘”].*[‘”]W?)*$
- field EXT_FORCE_OOFP: str = '' (name 'ext_force_oofp')#
A string of white-space separated atomic indices (4) followed by a single variable equation surrounded in either a single or double quotation mark. Example:
"1 2 3 4 'Sin(x)'"evaluates the force along the coordinate as a 1-D sinusoidal function where x is the “value” (angle [radians]) of the coordinate (oofp)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+W?s+d+s*W*s*s*[‘”].*[‘”]W?)*$
- field FLEXIBLE_G_CONVERGENCE: bool = False (name 'flexible_g_convergence')#
Normally, any specified
<any>_G_CONVERGENCEkeyword likeMAX_FORCE_G_CONVERGENCEwill be obeyed exclusively. If active,FLEXIBLE_G_CONVERGENCEappends to G_CONVERGENCE with the value from<any>_G_CONVERGENCEinstead of overriding.
- field FRAG_MODE: str = 'SINGLE' (name 'frag_mode')#
For multi-fragment molecules, treat as single bonded molecule or via interfragment coordinates. A primary difference is that in
MULTImode, the interfragment coordinates are not redundant.- Constraints:
pattern = re.compile(‘^(?:SINGLE|MULTI)$’, re.IGNORECASE)
- field FRAG_REF_ATOMS: list[list[list[int]]] = [] (name 'frag_ref_atoms')#
Which atoms define the reference points for interfragment coordinates? Example for a simple diatomic dimer like \(Ne_2\)
[[[1]], [[2]]]. Please see the section on multi-fragment optimizations for more information.
- field FREEZE_ALL_DIHEDRALS: bool = False (name 'freeze_all_dihedrals')#
A shortcut to request that all dihedrals should be frozen.
- field FREEZE_INTRAFRAG: bool = False (name 'freeze_intrafrag')#
Whether to freeze all intrafragment coordinates (rigid molecules). Only optimize the interfragment coordinates.
- field FROZEN_BEND: str = '' (name 'frozen_bend')#
A string of white-space separated atomic indices to specify that the distances between the atoms should be frozen (unchanged).
OptParams({"frozen_bend": "1 2 3 2 3 4"})FreezesBend(1 2 3)andBend(2 3 4)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+s*W*s*)*$
- field FROZEN_CARTESIAN: str = '' (name 'frozen_cartesian')#
A string of white-space separated atomic indices and Cartesian labels to specify that the Cartesian coordinates for a given atom should be frozen (unchanged).
OptParams({"frozen_cartesian": "1 XYZ 2 XY 2 Z"})FreezesCART(1, X),CART(1, Y),CART(1, Z),CART(2, X), etc…- Constraints:
pattern = (?i)^s*(?:W?ds(?:xyz|xy|yz|x|y|z)W*s*)*$
- field FROZEN_DIHEDRAL: str = '' (name 'frozen_dihedral')#
A string of white-space separated atomic indices to specify that the corresponding dihedral angle should be frozen (unchanged).
OptParams({"frozen_tors": "1 2 3 4 2 3 4 5"})FreezesTors(1, 2, 3, 4)andTors(2, 3, 4, 5)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+W?s+d+s*W*s*)*$
- field FROZEN_DISTANCE: str = '' (name 'frozen_distance')#
A string of white-space separated atomic indices to specify that the distances between the atoms should be frozen (unchanged).
OptParams({"frozen_distance": "1 2 3 4"})FreezesStre(1, 2)andStre(3, 4)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+s*W*s*)*$
- field FROZEN_OOFP: str = '' (name 'frozen_oofp')#
A string of white-space separated atomic indices to specify that the corresponding out-of-plane angle should be frozen. atoms should be frozen (unchanged).
OptParams({"frozen_oofp": "1 2 3 4"})FreezesOOFP(1, 2, 3, 4)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+W?s+d+s*W*s*)*$
- field FULL_HESS_EVERY: int = -1 (name 'full_hess_every')#
Frequency with which to compute the full Hessian in the course of a geometry optimization. 0 means to compute the initial Hessian only, 1 means recompute every step, and N means recompute every N steps. -1 indicates that the full hessian should never be computed.
- Constraints:
ge = -1
- field G_CONVERGENCE: str = 'QCHEM' (name 'g_convergence')#
A set of optimization criteria covering the change in energy, magnitude of the forces and the step_size. One of
["QCHEM", "MOLPRO", "GAU", "GAU_LOOSE", "GAU_TIGHT", "GAU_VERYTIGHT", "TURBOMOLE", "CFOUR", "NWCHEM_LOOSE", "INTERFRAG_TIGHT"]. Set of optimization criteria. Specification of anyMAX_<any>_G_CONVERGENCEorRMS_<any>_G_CONVERGENCEoptions will append to overwrite the criteria set here ifflexible_g_convergenceis also on. See Table Geometry Convergence for details.- Constraints:
pattern = re.compile(‘^(?:QCHEM|MOLPRO|GAU|GAU_LOOSE|GAU_TIGHT|GAU_VERYTIGHT|TURBOMOLE|CFOUR|NWCHEM_LOOSE|INTERFRAG_TIGHT)$’, re.IGNORECASE)
- field GEOM_MAXITER: int = 50 (name 'geom_maxiter')#
The maximum number of geometry optimization steps allowed - Technically this is the maximum number of gradients that Optking is allowed to calculate.
- Constraints:
gt = 0
- field H_BOND_CONNECT: float = 4.3 (name 'h_bond_connect')#
General, maximum distance for the definition of H-bonds.
- Constraints:
gt = 0.0
- field H_GUESS_EVERY: bool = False (name 'h_guess_every')#
Re-estimate the Hessian at every step, i.e., ignore the currently stored Hessian. This is NOT recommended
- field HESS_UPDATE: str = 'BFGS' (name 'hess_update')#
one of:
[NONE, "BFGS", "MS", "POWELL", "BOFILL"]Update scheme for the hessian. Default depends on OPT_TYPE- Constraints:
pattern = ^(?:NONE|BFGS|MS|POWELL|BOFILL)$
- field HESS_UPDATE_DEN_TOL: float = 1e-07 (name 'hess_update_den_tol')#
Denominator check for hessian update.
- Constraints:
gt = 0.0
- field HESS_UPDATE_DQ_TOL: float = 0.5 (name 'hess_update_dq_tol')#
Hessian update is avoided if any internal coordinate has changed by more than this in radians/au
- Constraints:
ge = 0.0
- field HESS_UPDATE_LIMIT: bool = True (name 'hess_update_limit')#
Do limit the magnitude of changes caused by the Hessian update? If
hess_update_limitis True, changes to the Hessian from the update are limited to the larger ofhess_update_limit_scale* (current value) andhess_update_limit_max[au]. By default, a Hessian value cannot be changed by more than 50% and 1 au.
- field HESS_UPDATE_LIMIT_MAX: float = 1.0 (name 'hess_update_limit_max')#
Absolute upper limit for how much any given Hessian value can be changed when updating
- Constraints:
ge = 0.0
- field HESS_UPDATE_LIMIT_SCALE: float = 0.5 (name 'hess_update_limit_scale')#
Relative upper limit for how much any given Hessian value can be changed when updating
- Constraints:
ge = 0.0
le = 1.0
- field HESS_UPDATE_USE_LAST: int = 4 (name 'hess_update_use_last')#
Number of previous steps to use in Hessian update, 0 uses all steps.
- Constraints:
ge = 0
- field HESSIAN_FILE: Path = PosixPath('.') (name 'hessian_file')#
Accompanies
cart_hess_read. path to file where hessian has been saved.
- field INCLUDE_OOFP: bool = False (name 'include_oofp')#
Add out-of-plane angles (usually not needed)
- field INTERFRAG_COLLINEAR_TOL: float = 0.01 (name 'interfrag_collinear_tol')#
Used for determining which atoms in a system are too collinear to be chosen as default reference atoms. We avoid collinearity. Greater is more restrictive.
- Constraints:
gt = 0.0
- field INTERFRAG_COORDS: list[dict] = [] (name 'interfrag_coords')#
Let the user submit a dictionary (or array of dictionaries) for the interfrag coordinates. The string input must be “loadable” as a python dictionary. See input examples.
- field INTERFRAG_DIST_INV: bool = False (name 'interfrag_dist_inv')#
Do use 1/R for the interfragment stretching coordinate instead of R?
- field INTERFRAG_HESS: str = 'DEFAULT' (name 'interfrag_hess')#
Model Hessian to guess interfragment force constants. One of [“DEFAULT”, “FISCHER_LIKE”]
- Constraints:
pattern = re.compile(‘^(?:DEFAULT|FISCHER_LIKE)$’, re.IGNORECASE)
- field INTERFRAG_MODE: str = 'FIXED' (name 'interfrag_mode')#
One of [‘FIXED’, ‘PRINCIPAL_AXES’]. Use either principal axes or fixed linear combinations of atoms as reference points for generating the interfragment coordinates.
- Constraints:
pattern = re.compile(‘^(?:FIXED|PRINCIPAL_AXES)$’, re.IGNORECASE)
- field INTERFRAG_STEP_LIMIT: float = 0.5 (name 'interfrag_trust')#
Initial maximum step size in bohr or radian along an interfragment coordinate
- Constraints:
gt = 0.0
- field INTERFRAG_STEP_LIMIT_MAX: float = 1.0 (name 'interfrag_trust_max')#
Upper bound for dynamic trust radius [au] for interfragment coordinates
- Constraints:
gt = 0.0
- field INTERFRAG_STEP_LIMIT_MIN: float = 0.001 (name 'interfrag_trust_min')#
Lower bound for dynamic trust radius [au] for interfragment coordinates
- Constraints:
gt = 0.0
- field INTERFRAGMENT_CONNECT: float = 1.8 (name 'interfragment_connect')#
A string of whitecaps separated atomic indices (1) and Cartesian labels, followed by a single variable equation surrounded in either a single or double quotation mark. Example:
"1 X 'Sin(x)'"evaluates the force along the coordinate as 1 1-D sinusoidal function where x is the “value” (angle [bohr]) of the coordinate (bohr)- Constraints:
gt = 0.0
- field INTRAFRAG_HESS: str = 'SCHLEGEL' (name 'intrafrag_hess')#
Model Hessian to guess intrafragment force constants. One of
["SCHLEGEL", "FISCHER", "SIMPLE", "LINDH", "LINDH_SIMPLE"]- Constraints:
pattern = re.compile(‘^(?:SCHLEGEL|FISCHER|SIMPLE|LINDH|LINDH_SIMPLE)$’, re.IGNORECASE)
- field INTRAFRAG_STEP_LIMIT: float = 0.5 (name 'intrafrag_trust')#
Initial maximum step size in bohr or radian along an internal coordinate for trust region methods (
RFOandRS_I_RFO)- Constraints:
gt = 0.0
- field INTRAFRAG_STEP_LIMIT_MAX: float = 1.0 (name 'intrafrag_trust_max')#
Upper bound for dynamic trust radius [au]
- Constraints:
gt = 0.0
- field INTRAFRAG_STEP_LIMIT_MIN: float = 0.001 (name 'intrafrag_trust_min')#
Lower bound for dynamic trust radius [au]
- Constraints:
gt = 0.0
- field IRC_CONVERGENCE: float = -0.7 (name 'irc_convergence')#
Main criteria for declaring convergence for an IRC. The overlap between the unit forces at two points of the IRC is compared to this value to assess whether a minimum has been stepped over. If \(overlap < irc_convergence\), declare convergence. If an IRC terminates too early, this may be symptomatic of a highly curved reaction-path, decrease try
irc_converence = -0.9- Constraints:
gt = -1.0
lt = -0.5
- field IRC_DIRECTION: str = 'FORWARD' (name 'irc_direction')#
One of
["FORWARD", "BACKWARD"]. Whether to step in the forward (+) direction along the transition state mode (smallest mode of hessian) or backward (-)- Constraints:
pattern = re.compile(‘^(?:FORWARD|BACKWARD)$’, re.IGNORECASE)
- field IRC_MODE: str = 'NORMAL' (name 'irc_mode')#
Experimental - One of
['NORMAL', 'CONFIRM'].'CONFIRM'is meant to be used for dissociation reactions. The IRC is terminated once the molecule’s connectivity has changed. Convergence is declared once the originalcovalent_connectmust be increased by more than 0.4 au. to connect all atoms- Constraints:
pattern = re.compile(‘^(?:NORMAL|CONFIRM)$’, re.IGNORECASE)
- field IRC_POINTS: int = 20 (name 'irc_points')#
Maximum number of converged points along the IRC path to map out before quitting. For dissociation reactions, where the reaction path may not terminate in a minimum, this is needed to cap the number of step’s Optking is allowed to take
- Constraints:
gt = 0
- field IRC_STEP_SIZE: float = 0.2 (name 'irc_step_size')#
Specifies the distance between each converged point along the IRC reaction path in \(bohr amu^{1/2}\)
- Constraints:
gt = 0.0
- field LINESEARCH_STEP: float = 0.1 (name 'linesearch_step')#
stepsize to start with when displacing to perform linesearch
- Constraints:
gt = 0.0
- field OPT_COORDINATES: str = 'INTERNAL' (name 'opt_coordinates')#
One of
["REDUNDANT", "INTERNAL", "CARTESIAN", "BOTH"]."INTERNAL"is just a synonym for"REDUNDANT"."BOTH"utilizes a full set of redundant internal coordinates and cartesian \((3N - 6+) + (3N) = (6N - 6+)\) coordinates.- Constraints:
pattern = re.compile(‘^(?:REDUNDANT|INTERNAL|DELOCALIZED|NATURAL|CARTESIAN|BOTH)$’, re.IGNORECASE)
- field OPT_TYPE: str = 'MIN' (name 'opt_type')#
- One of
["MIN", "TS", or "IRC"].OPT_TYPEwill be changed ifOPT_TYPEis not provided, but STEP_TYPEis provided, and the two are inconsistent. If both are provided but are inconsistent, an error will be raised.
Allowed
opt_typeandstep_typevaluesopt_typecompatible
step_typeMIN
RFO
NR
SD
LINESEARCH
Conjugate
TS
RS_I_RFO
P_RFO
IRC
N/A
- Constraints:
pattern = re.compile(‘^(?:MIN|TS|IRC)$’, re.IGNORECASE)
- One of
- field PRINT: int = 1 (name 'print_lvl')#
An integer between 1 (least printing) and 5 (most printing). This has been largely, but not entirely, replaced by using the logging modules
DEBUGandINFOlevels. Consider changing the logging handler in theloggingconfig.pyor if using a program like as Psi4 change the logging level from the command line.psi4 --loglevel=10...- Constraints:
ge = 1
le = 5
- field PRINT_TRAJECTORY_XYZ_FILE: bool = False (name 'print_trajectory_xyz_file')#
Deprecated- use
write_trajectoryinstead. Create an xyz file with geometries for each step of the Optimization. Ifopt_typeisIRC, only the converged IRC points will be included and the file will be named irc_traj.<pid>.json. Otherwise the file will contain all points and be namedopt_traj.<pid>.json
- field PROGRAM: str = 'psi4' (name 'program')#
What program to use for running gradient and energy calculations through
qcengine.
- field RANGED_BEND: str = '' (name 'ranged_bend')#
A string of white-space separated atomic indices and bounds for the angle between three atoms.
OptParams({1 2 3 100 110})Forces \(100^{\circ} <\)Bend(1, 2, 3)\(< 110^{\circ}\)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+s*W*s*W?-?d+.d+W?s+-?d+.d+s*W*s*)*$
- field RANGED_CARTESIAN: str = '' (name 'ranged_cartesian')#
A string of white-space separated atomic indices, Cartesian labels, and bounds for the Cartesian coordinates of a given atom.
OptParams({"ranged_cart": "1 XYZ 2.0 2.1"})Forces \(2.0 <\)Cart(1, X), Cart(1, Y), Cart(1, Z)\(< 2.1\) (Angstroms)- Constraints:
pattern = (?i)^s*(?:W?d+s+(?:xyz|xy|yz|x|y|z)W*s*W?-?d+.d+W?s+-?d+.d+s*W*s*)*$
- field RANGED_DIHEDRAL: str = '' (name 'ranged_dihedral')#
A string of white-space separated atomic indices and bounds for the torsion angle of four atoms. The order of specification determines whether the dihedral is a proper or improper torsion/dihedral.
OptParams({"ranged_dihedral": "1 2 3 4 100 110"})Forces \(100^{\circ} <\)Tors(1, 2, 3, 4)\(< 110^{\circ}\)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+W?s+d+s*W*s*W?-?d+.d+W?s+-?d+.d+s*W*s*)*$
- field RANGED_DISTANCE: str = '' (name 'ranged_distance')#
A string of white-space separated atomic indices and bounds for the distance between two atoms.
OptParams({"ranged_distance": 1 2 2.3 2.4})Forces \(2.3 <\)Stre(1, 2)\(< 2.4\)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+s*W*s*W?-?d+.d+W?s+-?d+.d+s*W*s*)*$
- field RANGED_OOFP: str = '' (name 'ranged_oofp')#
A string of white-space separated atomic indices and bounds for the out of plane angle defined by four atoms where the second atom is the central atom.
OptParams({"ranged_oofp": "1 2 3 4 100 110"})Forces \(100^{\circ} <\)Oofp(1, 2, 3, 4)\(< 110^{\circ}\)- Constraints:
pattern = ^s*(?:W?s*d+W?s+d+W?s+d+W?s+d+s*W*s*W?-?d+.d+W?s+-?d+.d+s*W*s*)*$
- field RFO_FOLLOW_ROOT: bool = False (name 'rfo_follow_root')#
Whether or not to optimize along the previously chosen mode of the augmented hessian matrix
- field RFO_NORMALIZATION_MAX: float = 100 (name 'rfo_normalization_max')#
Eigenvectors of RFO matrix with elements greater than this are ignored as candidates for the step direction.
- field RFO_ROOT: int = 0 (name 'rfo_root')#
root for
RFOorRS_I_RFOto follow. Changing rfo_root for aTSmay lead to a higher-order stationary point.- Constraints:
ge = 0
- field RSRFO_ALPHA_MAX: float = 100000000.0 (name 'rsrfo_alpha_max')#
Absolute maximum value of step scaling parameter in
RFOandRS_I_RFO.
- field SD_HESSIAN: float = 1.0 (name 'sd_hessian')#
Guess at Hessian in steepest-descent direction (acts as a stepsize control).
- Constraints:
gt = 0.0
- field SIMPLE_STEP_SCALING: bool = False (name 'simple_step_scaling')#
Do simple, linear scaling of internal coordinates to step limit instead of restricted-step (dynamic trust radius) approaches like
RS_RFOorRS_I_RFO
- field STEEPEST_DESCENT_TYPE: str = 'OVERLAP' (name 'steepest_descent_type')#
One of
["OVERLAP", "BARZILAI_BORWEIN"]. Change how theSDstep is calculated (scaled)- Constraints:
pattern = re.compile(‘^(?:OVERLAP|BARZILAI_BORWEIN)$’, re.IGNORECASE)
- field STEP_TYPE: str = 'RFO' (name 'step_type')#
One of
["RFO", "RS_I_RFO", "P_RFO", "NR", "SD", "LINESEARCH", "CONJUGATE"]. IfOPT_TYPEis set toTSandSTEP_TYPEis not specified.STEP_TYPEwill be set toRS_I_RFO.- Constraints:
pattern = re.compile(‘^(?:RFO|RS_I_RFO|P_RFO|NR|SD|LINESEARCH|CONJUGATE)$’, re.IGNORECASE)
- field UNFREEZE_DIHEDRALS: str = '' (name 'unfreeze_dihedrals')#
A string of white-space separated atomic indices to specify that the corresponding dihedral angle should be unfrozen. This keyword is meant to be used in conjunction with
FREEZE_ALL_DIHEDRALS- Constraints:
pattern = s*(?:W?s*d+W?s+d+W?s+d+W?s+d+s*W*s*)*$
- field WRITE_TRAJECTORY: bool = False (name 'write_trajectory')#
Create an xyz file with geometries for each step of the Optimization. If
opt_typeisIRC, only the converged IRC points will be included and the file will be named irc_traj.<pid>.json. Otherwise the file will contain all points and be namedopt_traj.<pid>.json
- conv_criteria() dict[source]#
Returns the currently active values for each convergence criteria. Not the original user input / presets
- classmethod from_internal_dict(params)[source]#
Assumes that params does not use the input key and syntax, but uses the internal names and internal syntax. Meant to be used for recreating options object after dump to dict It’s probably preferable to dump by alias and then recreate instead of using this
- model_post_init(context: Any, /) None#
This function is meant to behave like a BaseModel method to initialise private attributes.
It takes context as an argument since that’s what pydantic-core passes when calling it.
- Parameters:
self – The BaseModel instance.
context – The context.