Geometry Optimization¶
Code author: Rollin A. King and Alexander G. Heide
Section author: Rollin A. King, Alexander G. Heide, and Lori A. Burns
Note
As of version 1.10, the validation of OptKing’s options has been moved from PSI4 into
OptKing. Any options can still be set directly in PSI4 like in previous versions; however,
validation of any OptKing keywords will not take place until optimize is called.
If an option is not available in your version of PSI4 please make sure that PSI4
and OptKing are up-to-date. If an option is still not available, use the optimizer_keywords
argument to optimize()
.
PSI4 carries out molecular optimizations using a Python module called optking. The optking program takes as input nuclear gradients and, optionally, nuclear second derivatives — both in Cartesian coordinates. The default minimization algorithm employs an empirical model Hessian, redundant internal coordinates, an RFO step with trust radius scaling, and the BFGS Hessian update.
The principal literature references include the introduction of redundant internal coordinates by Peng et al. [Peng:1996:49]. The general approach employed in this code is similar to the “model Hessian plus RF method” described and tested by Bakken and Helgaker [Bakken:2002:9160]. However, for separated fragments, we have chosen not to employ their “extra-redundant” coordinates.
The internal coordinates are generated automatically based on an assumed bond
connectivity. The connectivity is determined by testing if the interatomic
distance is less than the sum of atomic radii times the value of
covalent_connect
. If the user finds that some
connectivity is lacking by default, then this value may be increased.
Warning
The selection of a Z-matrix input, and in particular the inclusion of dummy atoms, has no effect on the behavior of the optimizer, which begins from a Cartesian representation of the system.
Presently, by default, separate fragments are bonded by the
nearest atoms, and the whole system is treated as if it were part of one
molecule. However, with the option frag_mode
, fragments
may instead be related by a minimal set of interfragment coordinates defined by
reference points within each fragment. The reference points can be atomic
positions (current default) or linear combinations of atomic positions
(automatic use of principal axes is under development).
These dimer coordinates can be directly specified through interfrag_coords
)
See here <DimerSection_> for two examples of their use.
Basic Keywords¶
-
field OptParams.OPT_TYPE:
str
= 'MIN' One of
["MIN", "TS", or "IRC"]
.OPT_TYPE
will be changed ifOPT_TYPE
is not provided, butSTEP_TYPE
is provided, and the two are inconsistent. If both are provided but are inconsistent, an error will be raised.Allowed
opt_type
andstep_type
valuesopt_type
compatible
step_type
MIN
RFO
NR
SD
LINESEARCH
Conjugate
TS
RS_I_RFO
P_RFO
IRC
N/A
- Constraints:
pattern = (?i)^(?:MIN|TS|IRC)$
-
field OptParams.STEP_TYPE:
str
= 'RFO' One of
["RFO", "RS_I_RFO", "P_RFO", "NR", "SD", "LINESEARCH", "CONJUGATE"]
. IfOPT_TYPE
isTS
andSTEP_TYPE
is not specified thenSTEP_TYPE
will be set toRS_I_RFO
.- Constraints:
pattern = (?i)^(?:RFO|RS_I_RFO|P_RFO|NR|SD|LINESEARCH|CONJUGATE)$
-
field OptParams.GEOM_MAXITER:
int
= 50 The maximum number of geometry optimization steps allowed - Technically this is the maximum number of gradients that Optking is allowed to calculate.
- Constraints:
exclusiveMinimum = 0
-
field OptParams.G_CONVERGENCE:
str
= 'QCHEM' 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"]
. Specification of anyMAX_*_G_CONVERGENCE
orRMS_*_G_CONVERGENCE
options will overwrite the criteria set here. Ifflexible_g_convergence
is also on then the specified keyword will be appended. See Table Geometry Convergence for details.- Constraints:
pattern = (?i)^(?:QCHEM|MOLPRO|GAU|GAU_LOOSE|GAU_TIGHT|GAU_VERYTIGHT|TURBOMOLE|CFOUR|NWCHEM_LOOSE|INTERFRAG_TIGHT)$
-
field OptParams.FULL_HESS_EVERY:
int
= -1 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:
minimum = -1
Optimizing Minima¶
First, define the molecule and basis in the input.
molecule h2o {
O
H 1 1.0
H 1 1.0 2 105.0
}
set basis dz
Then the following are examples of various types of calculations that can be completed.
Optimize a geometry using default methods (RFO step):
optimize('scf')
Optimize using Newton-Raphson steps instead of RFO steps:
set step_type nr optimize('scf')
Optimize using finite differences of energies instead of gradients:
optimize('scf', dertype='energy')
Optimize while limiting the initial step size to 0.1 au:
set intrafrag_step_limit 0.1 optimize('scf')
Optimize while always limiting the step size to 0.1 au:
set {
intrafrag_step_limit 0.1
intrafrag_step_limit_min 0.1
intrafrag_step_limit_max 0.1
}
optimize('scf')
Optimize while calculating the Hessian at every step:
set full_hess_every 1
optimize('scf')
import optking
Hessian¶
If Cartesian second derivatives are available, optking can read them
and transform them into internal coordinates to make an initial Hessian in
internal coordinates. Otherwise, several empirical Hessians are available,
including those of Schlegel [Schlegel:1984:333] and Fischer and Almlof
[Fischer:1992:9770].
Either of these or a simple diagonal Hessian may be selected using the
intrafrag_hess
keyword.
All the common Hessian update schemes are available. For formulas, see Schlegel [Schlegel:1987:AIMQC] and Bofill [Bofill:1994:1].
The Hessian may be computed during an optimization using the
full_hess_every
keyword.
Transition States and Reaction Paths¶
Optking currently has two transition state algorithms. The current default is the
newer RS_I_RFO algorithm [Besalu:1998:265] . The old algorithm can be used by setting
STEP_TYPE P_RFO
for OPT_TYPE TS
Calculate a starting Hessian and optimize the “transition state” of linear water (note that without a reasonable starting geometry and Hessian, such a straightforward search often fails):
molecule h2o { O H 1 1.0 H 1 1.0 2 160.0 } set { basis dz full_hess_every 0 opt_type ts } optimize('scf')
At a transition state (planar HOOH), compute the second derivative, and then follow the intrinsic reaction path to the minimum:
molecule hooh { symmetry c1 H O 1 0.946347 O 2 1.397780 1 107.243777 H 3 0.946347 2 107.243777 1 0.0 } set { basis dzp opt_type irc geom_maxiter 50 } frequencies('scf') optimize('scf')
Constrained Optimizations¶
Optimize a geometry (HOOH) at a frozen dihedral angle of 90 degrees.
molecule { H O 1 0.90 O 2 1.40 1 100.0 H 3 0.90 2 100.0 1 90.0 } set optking { frozen_dihedral = (" 1 2 3 4 ") } optimize('scf')
To instead freeze the two O-H bond distances
set optking { frozen_distance = (" 1 2 3 4 ") }
For bends, the corresponding keyword is “frozen_bend”.
To freeze the cartesian coordinates of atom 2
freeze_list = """
2 xyz
"""
set optking frozen_cartesian $freeze_list
To freeze only the y coordinates of atoms 2 and 3
freeze_list = """
2 y
3 y
"""
set optking frozen_cartesian $freeze_list
To optimize toward a value of 0.95 Angstroms for the distance between atoms 1 and 3, as well as that between 2 and 4
set optking {
ranged_distance = ("
1 3 0.949 0.95
2 4 0.949 0.95
")
}
Note
The effect of the frozen and ranged keywords is generally independent of how the geometry of the molecule was input (whether Z-matrix or Cartesian, etc.).. At this time; however, enforcing Cartesian constraints when using a zmatrix for molecular input is not supported. Freezing or constraining Cartesian coordinates requires Cartesian molecule input. If numerical errors results in symmetry breaking, while Cartesian constraints are active, symmetrization cannot occur and an error will be raised, prompting you to restart the job.
As a shortcut, the entire set of dihedral angles can be frozen. A subset can then be unfrozen if desired.
set {
freeze_all_dihedrals true
unfreeze_dihedrals "1 2 3 4"
}
To scan the potential energy surface by optimizing at several fixed values of the dihedral angle of HOOH.
molecule hooh {
0 1
H 0.850718 0.772960 0.563468
O 0.120432 0.684669 -0.035503
O -0.120432 -0.684669 -0.035503
H -0.850718 -0.772960 0.563468
}
set {
basis cc-pvdz
intrafrag_step_limit 0.1
}
lower_bound = [99.99, 109.99, 119.99, 129.99, 149.99]
upper_bound = [100, 110, 120, 130, 140, 150]
PES = []
for lower, upper in zip(lower_bound, upper_bound):
my_string = f"1 2 3 4 {lower} {upper}"
set optking ranged_dihedral = $my_string
E = optimize('scf')
PES.append((upper, E))
print("\n\tcc-pVDZ SCF energy as a function of phi\n")
for point in PES:
print("\t%5.1f%20.10f" % (point[0], point[1]))
To scan the potential energy surface without the
ranged_dihedral
keyword, a zmatrix can be used.
Warning
Rotating dihedrals in large increments without allowing the molecule to relax in between increments can lead to unphysical geometries with overlapping functional groups in larger molecules, which may prevent successful constrained optimzations. Furthermore, such a relaxed scan of the PES does not always procude a result close to an IRC, or even a reaction path along which the energy changes in a continuous way.
molecule hooh {
0 1
H
O 1 0.95
O 2 1.39 1 103
H 3 0.95 2 103 1 D
D = 99
units ang
}
set {
basis cc-pvdz
intrafrag_step_limit 0.1
frozen_dihedral (" 1 2 3 4 ")
}
dihedrals = [100, 110, 120, 130, 140, 150]
PES = []
for phi in dihedrals:
hooh.D = phi
E = optimize('scf')
PES.append((phi, E))
print("\n\tcc-pVDZ SCF energy as a function of phi\n")
for point in PES:
print("\t%5.1f%20.10f" % (point[0], point[1]))
Multi-Fragment Optimizations¶
In previous versions of optking, the metric for connecting atoms was increased until all atoms
were connected. This is the current behavior for frag_mode
single`
.
Setting frag_mode
to multi
will now add a special
set of intermolecular coordinates between fragments - internally referred to as DimerFrag
coordinates (see here <DimerIntro_> for the brief description).
For each pair of molecular fragments, a set of up to 3 reference points
are chosen on each fragment. Each reference point will be either an atom or a linear combination
of positions of atoms within that fragment. Stretches, bends, and dihedral angles between the two
fragments will be created using these reference points. See
Dimer coordinate table for how reference points are created.
For a set of three dimers A, B, and C, sets of coordinates are created between each pair:
AB, AC, and BC. Each DimerFrag may use different reference points.
Creation of the intermolecular coordinates can be controlled through frag_ref_atoms
and interfrag_coords
. frag_ref_atoms
specifies which atoms
(or linear combination of atoms) to use for the reference points and interfrag_coords
,
which encompasses frag_ref_atoms
, allows for constraints and labels to be added to the
intermolecular coordinates.
Note
Manual specification of the interfragment coordinates is supported for power users,
and provides complete control of fragments’ relative orientations.
Setting interfrag_mode
to multi should suffice in almost all cases.
Dimer coordinate table. provides the name and ordering
convention for the coordinates.
Basic multi-fragment optimization. Use automatically generated reference points.
memory 4GB
molecule mol {
0 1
O -0.5026452583 -0.9681078610 -0.4772692868
H -2.3292990446 -1.1611084524 -0.4772692868
H -0.8887241813 0.8340933116 -0.4772692868
--
0 1
C 0.8853463281 -5.2235996493 0.5504918473
C 1.8139169342 -2.1992967152 3.8040686146
C 2.8624456357 -4.1143863257 0.5409035710
C -0.6240195463 -4.8153482424 2.1904642137
C -0.1646305764 -3.3031992532 3.8141619690
C 3.3271056135 -2.6064153737 2.1669340785
H 0.5244823836 -6.4459192939 -0.7478283184
H 4.0823309159 -4.4449979205 -0.7680411190
H -2.2074914566 -5.7109913627 2.2110247636
H -1.3768100495 -2.9846751653 5.1327625515
H 4.9209603634 -1.7288723155 2.1638694922
H 2.1923374156 -0.9964630692 5.1155773223
nocom
units au
}
set {
basis 6-31+G
frag_mode MULTI
}
optimize("mp2")
Warning
The molecule input for psi4 has no effect upon optking, expect to provide Cartesian
coordinates. Specifying independent fragments with the – seperator, will not trigger
optking to add specific interfragment coordinates. Use frag_mode
multi
.
Specify the reference points to use for coordinates via
frag_ref_atoms
. Each list corresponds to a fragment. A list of indices denotes a linear combination of the atoms. In this case, the first reference point for the second fragment is the center of the benzene ring. Indexing starts at 1, so the second fragment in this example starts at index 4.
memory 4GB
molecule mol {
0 1
O -0.5026452583 -0.9681078610 -0.4772692868
H -2.3292990446 -1.1611084524 -0.4772692868
H -0.8887241813 0.8340933116 -0.4772692868
--
0 1
C 0.8853463281 -5.2235996493 0.5504918473
C 1.8139169342 -2.1992967152 3.8040686146
C 2.8624456357 -4.1143863257 0.5409035710
C -0.6240195463 -4.8153482424 2.1904642137
C -0.1646305764 -3.3031992532 3.8141619690
C 3.3271056135 -2.6064153737 2.1669340785
H 0.5244823836 -6.4459192939 -0.7478283184
H 4.0823309159 -4.4449979205 -0.7680411190
H -2.2074914566 -5.7109913627 2.2110247636
H -1.3768100495 -2.9846751653 5.1327625515
H 4.9209603634 -1.7288723155 2.1638694922
H 2.1923374156 -0.9964630692 5.1155773223
nocom
units au
}
set {
basis 6-31+G
frag_mode MULTI
# The line below specifies the reference points that will be used to construct the
# interfragment coordinates between the two fragments (called A and B).
# The format is the following:
# [[A-1], [A-2], [A-3]], [[B-1], [B-2], [B-3]]
#
# In terms of atoms within each fragment, the line below chooses, for water:
# H3 of water for the first reference point, O1 of water for the second reference point, and
# H2 of water for the third reference point.
# For benzene: the mean of the positions of all the C atoms, C2, one of the Carbon atoms,
# and C6, another one of the carbon atoms.
frag_ref_atoms [
[[3], [1], [2]], [[4, 5, 6, 7, 8, 9], [5], [9]]
]
}
optimize("mp2")
For even greater control, a dictionary can be passed to interfrag_coords
The coordinates that are created between two dimers depend upon the number of atoms present The fragments A and B have up to 3 reference atoms each as shown in Dimer coordinate table. The interfragment coordinates are named and can be frozen according to their names as show in example below. For specifying reference points, use 1 based indexing.
name |
type |
atom-labels |
present, if |
---|---|---|---|
RAB |
distance |
A0-B0 |
always |
theta_A |
angle |
A1-A0-B0 |
A has > 1 atom |
theta_B |
angle |
A0-B0-B1 |
B has > 1 atom |
tau |
dihedral |
A1-A0-B0-B1 |
A and B have > 1 atom |
phi_A |
dihedral |
A2-A1-A0-B0 |
A has > 2 atoms. Is not linear |
phi_B |
dihedral |
A0-B0-B1-B2 |
B has > 2 atoms. Is not linear |
A constrained optimization is performed where the orientation of the two fragments is fixed but the distance between the fragments and all intrafragment coordinates are allowed to relax. In this example, the centers of the benzene and thiophene rings are selected for the first reference points. The methyl groups carbon and one hydrogen are selected for the other two reference points on the first fragments. For fragment two, two carbons of the benzene ring are chosen for the other reference points.
memory 4GB
molecule mol {
C -1.258686 0.546935 0.436840
H -0.683650 1.200389 1.102833
C -0.699036 -0.349093 -0.396608
C -2.693370 0.550414 0.355311
H -3.336987 1.206824 0.952052
C -3.159324 -0.343127 -0.536418
H -4.199699 -0.558111 -0.805894
S -1.883829 -1.212288 -1.301525
C 0.786082 -0.656530 -0.606057
H 1.387673 -0.016033 0.048976
H 1.054892 -0.465272 -1.651226
H 0.978834 -1.708370 -0.365860
--
C -6.955593 -0.119764 -1.395442
C -6.977905 -0.135060 1.376787
C -7.111625 1.067403 -0.697024
C -6.810717 -1.314577 -0.707746
C -6.821873 -1.322226 0.678369
C -7.122781 1.059754 0.689090
H -7.226173 2.012097 -1.240759
H -6.687348 -2.253224 -1.259958
H -6.707325 -2.266920 1.222105
H -7.246150 1.998400 1.241304
O -6.944245 -0.111984 -2.805375
H -7.058224 0.807436 -3.049180
C -6.990227 -0.143507 2.907714
H -8.018305 -0.274985 3.264065
H -6.592753 0.807024 3.281508
H -6.368443 -0.968607 3.273516
nocom
unit angstrom
}
# Create a python dictionary and convert to string for pass through to optking
MTdimer = """{
"Natoms per frag": [12, 16],
"A Frag": 1,
"A Ref Atoms": [[1, 3, 4, 6, 8], [8], [11]],
"A Label": "methylthiophene",
"B Frag": 2,
"B Ref Atoms": [[13, 14, 15, 16, 17, 18], [13], [15]],
"B Label": "tyrosine",
"Frozen": ["theta_A", "theta_B", "tau", "phi_A", "phi_B"],
}"""
set {
basis 6-31+G
frag_mode MULTI
interfrag_coords $MTdimer
}
optimize("mp2")
Dealing with problematic optimizations¶
Although optking is continuously improved with robustness in mind, some attempted optimizations will inevitably fail to converge to the desired minima. For difficult cases, the following suggestions are made.
As for any optimizer, computing the Hessian and limiting the step size will successfully converge a higher percentage of cases. The default settings have been chosen because they perform efficiently for common, representative test sets. More restrictive, cautious steps are sometimes necessary.
dynamic_level
allows optking to change the method of optimization toward algorithms that, while often less efficient, may help to converge difficult cases. If this is initially set to 1, then optking, as poor steps are detected, will increase the dynamic level through several forms of more robust and cautious algorithms. The changes will reduce the trust radius, allow backward steps (partial line searching), add cartesian coordinates, switch to cartesian coordinates, and take steepest-descent steps.The developers have found the
opt_coordinates
set to “BOTH” which includes both the redundant internal coordinate set, as well as cartesian coordinates, works well for systems with long ‘arms’ or floppy portions of a molecule poorly described by local internals.Optking does support the specification of ghost atoms. Certain internal coordinates such as torsions become poorly defined when they contain near-linear bends. An internal error AlgError may be raised in such cases. Optking will avoid such coordinates when choosing an initial coordinate system; however, they may arise in the course of an optimization. In such cases, try restarting from the most recent geometry. Alternatively, setting
opt_coordinates
to cartesian will avoid any internal coordinate difficulties altogether. These coordinate changes can be automatically performed by turningdynamic_level
to 1.
Warning
In some cases, such as the coordinate issues described above, optking will reset to maintain
a consistent history. If an error occurs in Psi4 due to geom_maxiter
being exceeded but
the final step report indicates that optking has not taken geom_maxiter
steps, such a
reset has occured. Inspection will show that the step counter was reset to 1 somewhere in the
optimization.
Convergence Criteria¶
Optking monitors five quantities to evaluate the progress of a geometry
optimization. These are (with their keywords) the change in energy
(max_energy_g_convergence
), the maximum element of
the gradient (max_force_g_convergence
), the root-mean-square
of the gradient (rms_force_g_convergence
), the maximum element
of displacement (max_disp_g_convergence
), and the
root-mean-square of displacement (rms_disp_g_convergence
),
all in internal coordinates and atomic units. Usually, these options will not
be set directly. Primary control for geometry convergence lies with the keyword
g_convergence
which sets the aforementioned in accordance
with Table Geometry Convergence.
|
Max Energy |
Max Force |
RMS Force |
Max Disp |
RMS Disp |
---|---|---|---|---|---|
NWCHEM_LOOSE [4] |
\(4.5 \times 10^{-3}\) |
\(3.0 \times 10^{-3}\) |
\(5.4 \times 10^{-3}\) |
\(3.6 \times 10^{-3}\) |
|
GAU_LOOSE [6] |
\(2.5 \times 10^{-3}\) |
\(1.7 \times 10^{-3}\) |
\(1.0 \times 10^{-2}\) |
\(6.7 \times 10^{-3}\) |
|
TURBOMOLE [4] |
\(1.0 \times 10^{-6}\) |
\(1.0 \times 10^{-3}\) |
\(5.0 \times 10^{-4}\) |
\(1.0 \times 10^{-3}\) |
\(5.0 \times 10^{-4}\) |
\(4.5 \times 10^{-4}\) |
\(3.0 \times 10^{-4}\) |
\(1.8 \times 10^{-3}\) |
\(1.2 \times 10^{-3}\) |
||
CFOUR [4] |
\(1.0 \times 10^{-4}\) |
||||
\(1.0 \times 10^{-6}\) |
\(3.0 \times 10^{-4}\) |
\(1.2 \times 10^{-3}\) |
|||
\(1.0 \times 10^{-6}\) |
\(3.0 \times 10^{-4}\) |
\(3.0 \times 10^{-4}\) |
|||
INTERFRAG_TIGHT [7] |
\(1.0 \times 10^{-6}\) |
\(1.5 \times 10^{-5}\) |
\(1.0 \times 10^{-5}\) |
\(6.0 \times 10^{-4}\) |
\(4.0 \times 10^{-4}\) |
\(1.5 \times 10^{-5}\) |
\(1.0 \times 10^{-5}\) |
\(6.0 \times 10^{-5}\) |
\(4.0 \times 10^{-5}\) |
||
GAU_VERYTIGHT [6] |
\(2.0 \times 10^{-6}\) |
\(1.0 \times 10^{-6}\) |
\(6.0 \times 10^{-6}\) |
\(4.0 \times 10^{-6}\) |
Footnotes
For ultimate control, specifying a value for any of the five monitored options activates that
criterium and overwrites/appends it to the criteria set by g_convergence
.
Note that this revokes the special convergence arrangements detailed in notes [5] and [6]
and instead requires all active criteria to be fulfilled to
achieve convergence. To avoid this revokation, turn on keyword flexible_g_convergence
.
Interface to GeomeTRIC¶
The GeomeTRIC optimizer developed by Wang and Song [Wang:2016:214108] may be used in place of Psi4’s native Optking optimizer. GeomeTRIC uses a translation-rotation-internal coordinate (TRIC) system that works well for optimizing geometries of systems containing noncovalent interactions.
Use of the GeomeTRIC optimizer is specified with the engine
argument to
optimize()
. The optimization will respect the keywords g_convergence
and geom_maxiter
. Any other GeomeTRIC-specific options (including constraints)
may be specified with the optimizer_keywords
argument to optimize()
.
Constraints may be placed on cartesian coordinates, bonds, angles, and dihedrals, and they can be
used to either freeze a coordinate or set it to a specific value. See the GeomeTRIC github
for more information on keywords and JSON specification of constraints.
Optimize the water molecule using GeomeTRIC:
molecule h2o { O H 1 1.0 H 1 1.0 2 160.0 } set { maxiter 100 g_convergence gau } optimize('hf/cc-pvdz', engine='geometric')
Optimize the water molecule using GeomeTRIC, with one of the two OH bonds constrained to 2.0 au and the HOH angle constrained to 104.5 degrees:
molecule h2o { O H 1 1.0 H 1 1.0 2 160.0 } set { maxiter 100 g_convergence gau } geometric_keywords = { 'coordsys' : 'tric', 'constraints' : { 'set' : [{'type' : 'distance', 'indices' : [0, 1], 'value' : 2.0 }, {'type' : 'angle', 'indices' : [1, 0, 2], 'value' : 104.5 }] } } optimize('hf/cc-pvdz', engine='geometric', optimizer_keywords=geometric_keywords)
Optimize the benzene/water dimer using GeomeTRIC, with the 6 carbon atoms of benzene frozen in place:
molecule h2o { C 0.833 1.221 -0.504 H 1.482 2.086 -0.518 C 1.379 -0.055 -0.486 H 2.453 -0.184 -0.483 C 0.546 -1.167 -0.474 H 0.971 -2.162 -0.466 C -0.833 -1.001 -0.475 H -1.482 -1.867 -0.468 C -1.379 0.275 -0.490 H -2.453 0.404 -0.491 C -0.546 1.386 -0.506 H -0.971 2.381 -0.524 -- O 0.000 0.147 3.265 H 0.000 -0.505 2.581 H 0.000 0.965 2.790 no_com no_reorient } set { maxiter 100 g_convergence gau } geometric_keywords = { 'coordsys' : 'tric', 'constraints' : { 'freeze' : [{'type' : 'xyz', 'indices' : [0, 2, 4, 6, 8, 10]}] } } optimize('hf/cc-pvdz', engine='geometric', optimizer_keywords=geometric_keywords)
Output¶
The progress of a geometry optimization can be monitored by grepping the output file for the
tilde character (~
). This produces a table like the one below that shows
for each iteration the value for each of the five quantities and whether the criterion
is active and fulfilled (*
), active and unfulfilled ( ), or inactive (o
).
--------------------------------------------------------------------------------------------- ~
Step Total Energy Delta E MAX Force RMS Force MAX Disp RMS Disp ~
--------------------------------------------------------------------------------------------- ~
Convergence Criteria 1.00e-06 * 3.00e-04 * o 1.20e-03 * o ~
--------------------------------------------------------------------------------------------- ~
1 -38.91591820 -3.89e+01 6.91e-02 5.72e-02 o 1.42e-01 1.19e-01 o ~
2 -38.92529543 -9.38e-03 6.21e-03 3.91e-03 o 2.00e-02 1.18e-02 o ~
3 -38.92540669 -1.11e-04 4.04e-03 2.46e-03 o 3.63e-02 2.12e-02 o ~
4 -38.92548668 -8.00e-05 2.30e-04 * 1.92e-04 o 1.99e-03 1.17e-03 o ~
5 -38.92548698 -2.98e-07 * 3.95e-05 * 3.35e-05 o 1.37e-04 * 1.05e-04 o ~
Information on the Psithon function that drives geometry optimizations is provided
at optimize()
.
Important User Changes from cpp-optking¶
FIXED_COORD keywords have been generalized to RANGED_COORD e.g.
ranged_distance
Detailed optimization is now printed through the python logging system. If more information about the optimization is needed. Please see <output_name>.log
Keywords¶
- pydantic model optking.v1.optparams.OptParams[source]¶
-
field ACCEPT_SYMMETRY_BREAKING:
bool
= False¶ 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¶ Add bond coordinates for atoms separated by less than \(2.5 \times\) their covalent radii
-
field AUXILIARY_BOND_FACTOR:
float
= 2.5¶ This factor times the standard covalent distance is used to add extra stretch coordinates.
- Constraints:
exclusiveMinimum = 1.0
-
field BT_DX_CONV:
float
= 1e-07¶ Threshold for the change in any given Cartesian coordinate during iterative back-transformation.
- Constraints:
exclusiveMinimum = 0.0
-
field BT_DX_RMS_CHANGE_CONV:
float
= 1e-12¶ Threshold for RMS change in Cartesian coordinates during iterative back-transformation.
- Constraints:
exclusiveMinimum = 0.0
-
field BT_MAX_ITER:
int
= 25¶ Maximum number of iterations allowed to converge back-transformation
- Constraints:
exclusiveMinimum = 0
-
field BT_PINV_RCOND:
float
= 1e-06¶ 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:
exclusiveMinimum = 0.0
-
field CART_HESS_READ:
bool
= False¶ Do read Cartesian Hessian? Recommended to use
full_hess_every
instead. cfour format or.json
file (AtomicResult) allowed. The filetype is determined by the presence of a.json
extension. 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.
-
field CONJUGATE_GRADIENT_TYPE:
str
= 'FLETCHER'¶ One of
["POLAK", "FLETCHER", "DESCENT"]
. Changes how the step direction is calculated.- Constraints:
pattern = (?i)^(?:FLETCHER|DESCENT|POLAK)$
-
field CONSECUTIVE_BACKSTEPS:
int
= 0¶ Sets the number of consecutive backward steps allowed in an optimization. This option can be updated by
Optking
ifdynamic_lvl
is > 0. Not recommended for general use.- Constraints:
minimum = 0
-
field MAX_ENERGY_G_CONVERGENCE:
float
= 1e-06¶ Convergence criterion for geometry optimization: maximum energy change
-
field MAX_DISP_G_CONVERGENCE:
float
= 0.0012¶ Convergence criterion for geometry optimization: maximum displacement (internal coordinates, au)
-
field MAX_FORCE_G_CONVERGENCE:
float
= 0.0003¶ Convergence criterion for geometry optimization: maximum force (internal coordinates, au)
-
field RMS_DISP_G_CONVERGENCE:
float
= 0.0012¶ Convergence criterion for geometry optimization: rms displacement (internal coordinates, au)
-
field RMS_FORCE_G_CONVERGENCE:
float
= 0.0003¶ Convergence criterion for geometry optimization: maximum force (internal coordinates, au)
-
field COVALENT_CONNECT:
float
= 1.3¶ When determining connectivity, a bond is assigned if the interatomic distance is less than
COVANLENT_CONNECT
:math:` times ` the sum of covalent radii. When connecting disparate fragments andFRAG_MODE
is SINGLE, a “bond” is assigned if interatomic distance is less than (COVALENT_CONNECT
) \(\times\) sum of covalent radii. The value is then increased until all the fragments are connected directly or indirectly.- Constraints:
exclusiveMinimum = 0.0
-
field DYNAMIC_LVL:
int
= 0¶ An integer between 0 and 6. Larger values reflect less aggressive optimization techniques If
dynamic_lvl
is not set,Optking
will not change thedynamic_lvl
. Thedynamic_lvl
must 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
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:
minimum = 0
maximum = 6
-
field DYNAMIC_LVL_MAX:
int
= 0¶ How large
dynamic_lvl
is allowed to grow. Ifdynamic_lvl
\(> 0\),dynamic_lvl_max
will default to 6- Constraints:
minimum = 0
maximum = 6
-
field ENSURE_BT_CONVERGENCE:
bool
= False¶ Reduces step size as necessary to ensure convergence of back-transformation of internal coordinate step to Cartesian coordinates
-
field EXT_FORCE_BEND:
str
= ''¶ 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
= ''¶ 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
= ''¶ 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
= ''¶ 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
= ''¶ 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¶ Normally, any specified
*_G_CONVERGENCE
keyword likeMAX_FORCE_G_CONVERGENCE
will be obeyed exclusively. If active,FLEXIBLE_G_CONVERGENCE
appends toG_CONVERGENCE
with the value from*_G_CONVERGENCE
instead of overriding.
-
field FRAG_MODE:
str
= 'SINGLE'¶ For multi-fragment molecules, treat as single bonded molecule or via interfragment coordinates. A primary difference is that in
MULTI
mode, the interfragment coordinates are not redundant.- Constraints:
pattern = (?i)^(?:SINGLE|MULTI)$
-
field FRAG_REF_ATOMS:
list
[list
[list
[int
]]] = []¶ 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. Multi-Fragment Optimizations
-
field FREEZE_ALL_DIHEDRALS:
bool
= False¶ A shortcut to request that all dihedrals should be frozen.
-
field FREEZE_INTRAFRAG:
bool
= False¶ Whether to freeze all intrafragment coordinates (rigid molecules). Only optimize the interfragment coordinates.
-
field FROZEN_BEND:
str
= ''¶ 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
= ''¶ 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
= ''¶ 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
= ''¶ 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
= ''¶ 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¶ 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:
minimum = -1
-
field G_CONVERGENCE:
str
= 'QCHEM'¶ 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"]
. Specification of anyMAX_*_G_CONVERGENCE
orRMS_*_G_CONVERGENCE
options will overwrite the criteria set here. Ifflexible_g_convergence
is also on then the specified keyword will be appended. See Table Geometry Convergence for details.- Constraints:
pattern = (?i)^(?:QCHEM|MOLPRO|GAU|GAU_LOOSE|GAU_TIGHT|GAU_VERYTIGHT|TURBOMOLE|CFOUR|NWCHEM_LOOSE|INTERFRAG_TIGHT)$
-
field GEOM_MAXITER:
int
= 50¶ The maximum number of geometry optimization steps allowed - Technically this is the maximum number of gradients that Optking is allowed to calculate.
- Constraints:
exclusiveMinimum = 0
-
field H_BOND_CONNECT:
float
= 4.3¶ General, maximum distance for the definition of H-bonds.
- Constraints:
exclusiveMinimum = 0.0
-
field H_GUESS_EVERY:
bool
= False¶ Re-estimate the Hessian at every step, i.e., ignore the currently stored Hessian. This is NOT recommended
-
field HESS_UPDATE:
str
= 'BFGS'¶ one of:
[NONE, "BFGS", "MS", "POWELL", "BOFILL"]
Update scheme for the hessian. Default depends onOPT_TYPE
- Constraints:
pattern = (?i)^(?:NONE|BFGS|MS|POWELL|BOFILL)$
-
field HESS_UPDATE_DEN_TOL:
float
= 1e-07¶ Denominator check for hessian update.
- Constraints:
exclusiveMinimum = 0.0
-
field HESS_UPDATE_DQ_TOL:
float
= 0.5¶ Hessian update is avoided if any internal coordinate has changed by more than this in radians/au
- Constraints:
minimum = 0.0
-
field HESS_UPDATE_LIMIT:
bool
= True¶ Do limit the magnitude of changes caused by the Hessian update? If
hess_update_limit
is 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¶ Absolute upper limit for how much any given Hessian value can be changed when updating
- Constraints:
minimum = 0.0
-
field HESS_UPDATE_LIMIT_SCALE:
float
= 0.5¶ Relative upper limit for how much any given Hessian value can be changed when updating
- Constraints:
minimum = 0.0
maximum = 1.0
-
field HESS_UPDATE_USE_LAST:
int
= 4¶ Number of previous steps to use in Hessian update, 0 uses all steps.
- Constraints:
minimum = 0
-
field HESSIAN_FILE:
Path
= PosixPath('.')¶ Accompanies
CART_HESS_READ
. path to file where hessian has been saved. WARNING: As of Psi4 v1.10~nightly psi4.optimize() overrides this variable. If you have written a hessian to disk, copy the file topsi4.core.write_file_prefix(psi4.core.get_active_molecule().name())
or useoptking.optimize_psi4()
-
field INTERFRAG_COLLINEAR_TOL:
float
= 0.01¶ 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:
exclusiveMinimum = 0.0
-
field INTERFRAG_COORDS:
list
[dict
] = []¶ Let the user submit a dictionary (or array of dictionaries) for the interfrag coordinates. The input may also be given as a string, but the string input must be “loadable” as a python dictionary. See input examples Multi-Fragment Optimzations.
-
field INTERFRAG_DIST_INV:
bool
= False¶ Do use 1/R for the interfragment stretching coordinate instead of R?
-
field INTERFRAG_HESS:
str
= 'DEFAULT'¶ Model Hessian to guess interfragment force constants. One of
["DEFAULT", "FISCHER_LIKE"]
- Constraints:
pattern = (?i)^(?:DEFAULT|FISCHER_LIKE)$
-
field INTERFRAG_MODE:
str
= 'FIXED'¶ 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 = (?i)^(?:FIXED|PRINCIPAL_AXES)$
-
field INTERFRAG_STEP_LIMIT:
float
= 0.5¶ Initial maximum step size in bohr or radian along an interfragment coordinate
- Constraints:
exclusiveMinimum = 0.0
-
field INTERFRAG_STEP_LIMIT_MAX:
float
= 1.0¶ Upper bound for dynamic trust radius [au] for interfragment coordinates
- Constraints:
exclusiveMinimum = 0.0
-
field INTERFRAG_STEP_LIMIT_MIN:
float
= 0.001¶ Lower bound for dynamic trust radius [au] for interfragment coordinates
- Constraints:
exclusiveMinimum = 0.0
-
field INTRAFRAG_HESS:
str
= 'SCHLEGEL'¶ Model Hessian to guess intrafragment force constants. One of
["SCHLEGEL", "FISCHER", "SIMPLE", "LINDH", "LINDH_SIMPLE"]
- Constraints:
pattern = (?i)^(?:SCHLEGEL|FISCHER|SIMPLE|LINDH|LINDH_SIMPLE)$
-
field INTRAFRAG_STEP_LIMIT:
float
= 0.5¶ Initial maximum step size in bohr or radian in internal coordinates for trust region methods (
RFO
andRS_I_RFO
). This value will be updated throughout optimization.- Constraints:
exclusiveMinimum = 0.0
-
field INTRAFRAG_STEP_LIMIT_MAX:
float
= 1.0¶ Upper bound for dynamic trust radius [au]
- Constraints:
exclusiveMinimum = 0.0
-
field INTRAFRAG_STEP_LIMIT_MIN:
float
= 0.001¶ Lower bound for dynamic trust radius [au]
- Constraints:
exclusiveMinimum = 0.0
-
field IRC_CONVERGENCE:
float
= -0.7¶ 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:
exclusiveMinimum = -1.0
exclusiveMaximum = -0.5
-
field IRC_DIRECTION:
str
= 'FORWARD'¶ One of
["FORWARD", "BACKWARD"]
. Whether to step in the forward (+) direction along the transition state mode (smallest mode of hessian) or backward (-)- Constraints:
pattern = (?i)^(?:FORWARD|BACKWARD)$
-
field IRC_MODE:
str
= 'NORMAL'¶ 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 original
covalent_connect
must be increased by more than 0.4 au.- Constraints:
pattern = (?i)^(?:NORMAL|CONFIRM)$
-
field IRC_POINTS:
int
= 20¶ 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:
exclusiveMinimum = 0
-
field IRC_STEP_SIZE:
float
= 0.2¶ Specifies the distance between each converged point along the IRC reaction path in \(bohr amu^{1/2}\)
- Constraints:
exclusiveMinimum = 0.0
-
field LINESEARCH_STEP:
float
= 0.1¶ Initial stepsize when performing a 1-D linesearch
- Constraints:
exclusiveMinimum = 0.0
-
field OPT_COORDINATES:
str
= 'INTERNAL'¶ 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 = (?i)^(?:REDUNDANT|INTERNAL|DELOCALIZED|NATURAL|CARTESIAN|BOTH)$
-
field OPT_TYPE:
str
= 'MIN'¶ One of
["MIN", "TS", or "IRC"]
.OPT_TYPE
will be changed ifOPT_TYPE
is not provided, butSTEP_TYPE
is provided, and the two are inconsistent. If both are provided but are inconsistent, an error will be raised.Allowed
opt_type
andstep_type
valuesopt_type
compatible
step_type
MIN
RFO
NR
SD
LINESEARCH
Conjugate
TS
RS_I_RFO
P_RFO
IRC
N/A
- Constraints:
pattern = (?i)^(?:MIN|TS|IRC)$
-
field PRINT:
int
= 1¶ An integer between 1 (least printing) and 5 (most printing). This has been largely, but not entirely, replaced by using the logging modules
DEBUG
andINFO
levels. Consider changing the logging handler inloggingconfig.py
or if using Psi4 change the logging level from the command line.psi4 --loglevel=10...
- Constraints:
minimum = 1
maximum = 5
-
field PRINT_TRAJECTORY_XYZ_FILE:
bool
= False¶ Deprecated use
write_trajectory
instead. Create an xyz file with geometries for each step of the Optimization. Ifopt_type
isIRC
, 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'¶ What program to use for running gradient and energy calculations through
qcengine
.
-
field RANGED_BEND:
str
= ''¶ 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
= ''¶ 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
= ''¶ 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
= ''¶ 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
= ''¶ 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¶ Whether or not to optimize along the previously chosen mode of the augmented hessian matrix
-
field RFO_NORMALIZATION_MAX:
float
= 100¶ Eigenvectors of RFO matrix with elements greater than this are ignored as candidates for the step direction.
-
field RFO_ROOT:
int
= 0¶ root for
RFO
orRS_I_RFO
to follow. Changing rfo_root for aTS
may lead to a higher-order stationary point.- Constraints:
minimum = 0
-
field RSRFO_ALPHA_MAX:
float
= 100000000.0¶ Absolute maximum value of step scaling parameter in
RFO
andRS_I_RFO
.
-
field SD_HESSIAN:
float
= 1.0¶ Guess at Hessian in steepest-descent direction (acts as a stepsize control).
- Constraints:
exclusiveMinimum = 0.0
-
field SIMPLE_STEP_SCALING:
bool
= False¶ Do simple, linear scaling of internal coordinates to limit step instead of restricted-step (dynamic trust radius) approaches like
RS_RFO
orRS_I_RFO
-
field STEEPEST_DESCENT_TYPE:
str
= 'OVERLAP'¶ One of
["OVERLAP", "BARZILAI_BORWEIN"]
. Change how theSD
step is calculated (scaled)- Constraints:
pattern = (?i)^(?:OVERLAP|BARZILAI_BORWEIN)$
-
field STEP_TYPE:
str
= 'RFO'¶ One of
["RFO", "RS_I_RFO", "P_RFO", "NR", "SD", "LINESEARCH", "CONJUGATE"]
. IfOPT_TYPE
isTS
andSTEP_TYPE
is not specified thenSTEP_TYPE
will be set toRS_I_RFO
.- Constraints:
pattern = (?i)^(?:RFO|RS_I_RFO|P_RFO|NR|SD|LINESEARCH|CONJUGATE)$
-
field UNFREEZE_DIHEDRALS:
str
= ''¶ 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¶ Create an xyz file with geometries for each step of the Optimization. If
opt_type
isIRC
, 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()[source]¶
Returns the currently active values for each convergence criteria. Not the original user input / presets
- Return type:
- 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# o dict It’s probably preferable to dump by alias and then recreate instead of using this
-
field ACCEPT_SYMMETRY_BREAKING: