Geometry Optimization — optimize()
and gradient()
¶
- Psi4 Native Gradient Methods
- Psi4 Native DFT Gradient Methods (excepting double-hybrids)
- CFOUR Interfaced Gradient Methods
For further discussion of geometry optimization, see Sec. Geometry Optimization.
optimize()
is the only command most users will ever
need to access directly to perform geometry optimizations. Behind
the scenes, optimize()
is a wrapper that repeatedly
calls gradient()
that computes the gradient then adds a
call to the geometry optimization module.
-
psi4.
optimize
(name[, molecule, return_wfn, func, mode, dertype, hessian_with])[source]¶ Function to perform a geometry optimization.
Aliases: opt() Returns: float – Total electronic energy of optimized structure in Hartrees. Returns: (float, Wavefunction
) – energy and wavefunction when return_wfn specified.Raises: psi4.OptimizationConvergenceError if GEOM_MAXITER exceeded without reaching geometry convergence. PSI variables: Parameters: - name (string) –
'scf'
||'mp2'
||'ci5'
|| etc.First argument, usually unlabeled. Indicates the computational method to be applied to the database. May be any valid argument to
energy()
. - molecule (molecule) –
h2o
|| etc.The target molecule, if not the last molecule defined.
- return_wfn (boolean) –
'on'
|| \(\Rightarrow\)'off'
\(\Leftarrow\)Indicate to additionally return the
Wavefunction
calculation result as the second element (after float energy) of a tuple. - return_history (boolean) –
'on'
|| \(\Rightarrow\)'off'
\(\Leftarrow\)Indicate to additionally return dictionary of lists of geometries, energies, and gradients at each step in the optimization.
- func (function) –
\(\Rightarrow\)
gradient
\(\Leftarrow\) ||energy
||cbs
Indicates the type of calculation to be performed on the molecule. The default dertype accesses
'gradient'
or'energy'
, while'cbs'
performs a multistage finite difference calculation. If a nested series of python functions is intended (see Function Intercalls), use keywordopt_func
instead offunc
. - dertype (dertype) –
'gradient'
||'energy'
Indicates whether analytic (if available) or finite difference optimization is to be performed.
- hessian_with (string) –
'scf'
||'mp2'
|| etc.Indicates the computational method with which to perform a hessian analysis to guide the geometry optimization.
Warning
Optimizations where the molecule is specified in Z-matrix format with dummy atoms will result in the geometry being converted to a Cartesian representation.
Note
Analytic gradients area available for all methods in the table below. Optimizations with other methods in the energy table proceed by finite differences.
name calls method efp efp-only optimizations scf Hartree–Fock (HF) or density functional theory (DFT) [manual] hf HF self consistent field (SCF) [manual] dcft density cumulant functional theory [manual] mp2 2nd-order Møller–Plesset perturbation theory (MP2) [manual] [details] mp3 3rd-order Møller–Plesset perturbation theory (MP3) [manual] [details] mp2.5 average of MP2 and MP3 [manual] [details] omp2 orbital-optimized second-order MP perturbation theory [manual] omp3 orbital-optimized third-order MP perturbation theory [manual] omp2.5 orbital-optimized MP2.5 [manual] lccd Linear CCD [manual] [details] olccd orbital optimized LCCD [manual] ccd coupled cluster doubles (CCD) [manual] ccsd coupled cluster singles and doubles (CCSD) [manual] [details] ccsd(t) CCSD with perturbative triples (CCSD(T)) [manual] [details] eom-ccsd equation of motion (EOM) CCSD [manual] Examples: >>> # [1] Analytic hf optimization >>> optimize('hf')
>>> # [2] Finite difference mp5 optimization with gradient >>> # printed to output file >>> e, wfn = opt('mp5', return_wfn='yes') >>> wfn.gradient().print_out()
>>> # [3] Can automatically perform complete basis set extrapolations >>> optimize('MP2/cc-pV([D,T]+d)Z')
>>> # [4] Can automatically perform delta corrections that include extrapolations >>> # even with a user-defined extrapolation formula. See sample inputs named >>> # cbs-xtpl* for more examples of this input style >>> optimize("MP2/aug-cc-pv([d,t]+d)z + d:ccsd(t)/cc-pvdz", corl_scheme=myxtplfn_2)
>>> # [5] Get info like geometry, gradient, energy back after an >>> # optimization fails. Note that the energy and gradient >>> # correspond to the last optimization cycle, whereas the >>> # geometry (by default) is the anticipated *next* optimization step. >>> try: >>> optimize('hf/cc-pvtz') >>> except psi4.OptimizationConvergenceError as ex: >>> next_geom_coords_as_numpy_array = np.asarray(ex.wfn.molecule().geometry())
- name (string) –
-
psi4.
gradient
(name[, molecule, return_wfn, func, dertype])[source]¶ Function complementary to :py:func:~driver.optimize(). Carries out one gradient pass, deciding analytic or finite difference.
Returns: Matrix
– Total electronic gradient in Hartrees/Bohr.Returns: ( Matrix
,Wavefunction
) – gradient and wavefunction when return_wfn specified.Examples: >>> # [1] Single-point dft gradient getting the gradient >>> # in file, core.Matrix, and np.array forms >>> set gradient_write on >>> G, wfn = gradient('b3lyp-d', return_wfn=True) >>> wfn.gradient().print_out() >>> np.array(G)