Ground-state structure optimization#
Molecular structure optimizations always start from an initial guess of the molecular geometry (constructed using SMILES \(-\) see tutorial \(-\), a chemistry software, e.g. Avogadro, or downloaded from a database like ChemSpider or MolView). Using the initial guess, we will see in the following how to perform ground-state geometry optimization using VeloxChem and GeomeTRIC.
SCF geometry optimization#
import veloxchem as vlx
import py3Dmol as p3d
from matplotlib import pyplot as plt
import numpy as np
import geometric
import sys
# Set up the molecule and basis set
water_xyz = """3
water molecule, bond angle 150 deg
O 0.0000000 0.0000000 0.0000000
H 0.2760264 0.8822741 0.6495351
H 0.1036576 -0.7768012 -0.5501827
"""
basis_set_label = 'sto-3g'
molecule = vlx.Molecule.read_xyz_string(water_xyz)
basis = vlx.MolecularBasis.read(molecule, basis_set_label)
Show code cell source
# View initial structure of the molecule
view = p3d.view(viewergrid=(1,1),width=300, height=300)
view.addModel(water_xyz, 'xyz', viewer=(0,0))
view.setViewStyle({"style": "outline", "width": 0.05})
view.setStyle({'stick': {}, 'sphere': {'scale':0.25}})
view.zoomTo()
view.show()
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol
# Run the initial SCF and set up the gradient and optimization drivers
scf_drv = vlx.ScfRestrictedDriver()
# Settings for DFT:
# scf_drv.xcfun = "b3lyp"
scf_results = scf_drv.compute(molecule, basis)
Show code cell output
Self Consistent Field Driver Setup
====================================
Wave Function Model : Spin-Restricted Hartree-Fock
Initial Guess Model : Superposition of Atomic Densities
Convergence Accelerator : Two Level Direct Inversion of Iterative Subspace
Max. Number of Iterations : 50
Max. Number of Error Vectors : 10
Convergence Threshold : 1.0e-06
ERI Screening Scheme : Cauchy Schwarz + Density
ERI Screening Mode : Dynamic
ERI Screening Threshold : 1.0e-12
Linear Dependence Threshold : 1.0e-06
* Info * Nuclear repulsion energy: 8.4257146308 a.u.
* Info * Overlap matrix computed in 0.00 sec.
* Info * Kinetic energy matrix computed in 0.00 sec.
* Info * Nuclear potential matrix computed in 0.00 sec.
* Info * Orthogonalization matrix computed in 0.01 sec.
* Info * SAD initial guess computed in 0.00 sec.
* Info * Starting Reduced Basis SCF calculation...
* Info * ...done. SCF energy in reduced basis set: -74.838008565136 a.u. Time: 0.04 sec.
* Info * Overlap matrix computed in 0.00 sec.
* Info * Kinetic energy matrix computed in 0.00 sec.
* Info * Nuclear potential matrix computed in 0.00 sec.
* Info * Orthogonalization matrix computed in 0.00 sec.
Iter. | Hartree-Fock Energy | Energy Change | Gradient Norm | Max. Gradient | Density Change
--------------------------------------------------------------------------------------------
1 -74.838008656380 0.0000000000 0.00020607 0.00005271 0.00000000
2 -74.838008679393 -0.0000000230 0.00008668 0.00002305 0.00023046
3 -74.838008685191 -0.0000000058 0.00001417 0.00000312 0.00018085
4 -74.838008685300 -0.0000000001 0.00000241 0.00000057 0.00002422
5 -74.838008685302 -0.0000000000 0.00000049 0.00000009 0.00000358
*** SCF converged in 5 iterations. Time: 0.04 sec.
Spin-Restricted Hartree-Fock:
-----------------------------
Total Energy : -74.8380086853 a.u.
Electronic Energy : -83.2637233161 a.u.
Nuclear Repulsion Energy : 8.4257146308 a.u.
------------------------------------
Gradient Norm : 0.0000004862 a.u.
Ground State Information
------------------------
Charge of Molecule : 0.0
Multiplicity (2S+1) : 1.0
Magnetic Quantum Number (M_S) : 0.0
Spin Restricted Orbitals
------------------------
Molecular Orbital No. 1:
--------------------------
Occupation: 2.000 Energy: -20.09946 a.u.
( 1 O 1s : -0.99)
Molecular Orbital No. 2:
--------------------------
Occupation: 2.000 Energy: -1.14596 a.u.
( 1 O 1s : -0.24) ( 1 O 2s : 0.85) ( 3 H 1s : 0.20)
Molecular Orbital No. 3:
--------------------------
Occupation: 2.000 Energy: -0.58710 a.u.
( 1 O 1p-1: 0.45) ( 1 O 1p0 : 0.32) ( 2 H 1s : 0.46)
( 3 H 1s : -0.42)
Molecular Orbital No. 4:
--------------------------
Occupation: 2.000 Energy: -0.30532 a.u.
( 1 O 1p-1: -0.58) ( 1 O 1p0 : 0.81)
Molecular Orbital No. 5:
--------------------------
Occupation: 2.000 Energy: -0.30273 a.u.
( 1 O 2s : 0.25) ( 1 O 1p+1: -0.95) ( 2 H 1s : -0.19)
Molecular Orbital No. 6:
--------------------------
Occupation: 0.000 Energy: 0.37318 a.u.
( 1 O 2s : 0.85) ( 1 O 1p+1: 0.32) ( 2 H 1s : -0.87)
( 3 H 1s : -0.66)
Molecular Orbital No. 7:
--------------------------
Occupation: 0.000 Energy: 0.84399 a.u.
( 1 O 2s : -0.29) ( 1 O 1p-1: 0.85) ( 1 O 1p0 : 0.61)
( 2 H 1s : -0.56) ( 3 H 1s : 0.98)
# Set up the gradient and optimization dirvers:
grad_drv = vlx.ScfGradientDriver()
opt_drv = vlx.OptimizationDriver(grad_drv)
opt_drv.filename = "../../data/structopt/h2o_gsstructopt"
# Run the optimization
opt_molecule = opt_drv.compute(molecule, basis, scf_drv)
Let’s plot the energies calculated during the optimization and check how the structure changes.
# The data we need are stored in the output file generated during the geometry optimization.
# Let's write a routine which reads these data points from file.
def read_optimization_file(file_name):
output_file = open(file_name, "r")
structures = output_file.read()
output_file.close()
output_file = open(file_name, "r")
steps = []
energies = []
i = 1
for line in output_file:
if 'Energy' in line:
parts = line.split()
if parts:
energy = float(parts[-1])
energies.append(energy)
steps.append(i)
i += 1
output_file.close()
return steps, energies, structures
steps, energies, structures = read_optimization_file('../../data/structopt/h2o_gsstructopt_optim.xyz')
# Plot the energies
plt.figure(figsize=(5,4))
plt.plot(steps, energies,'o--')
plt.xlabel('Iteration')
plt.ylabel('Energy (H)')
plt.title("H2O GS optimization")
plt.tight_layout(); plt.show()
# Animate optimization
viewer = p3d.view(width=300, height=200)
viewer.addModelsAsFrames(structures)
viewer.animate({"loop": "forward"})
viewer.setViewStyle({"style": "outline", "width": 0.05})
viewer.setStyle({'stick': {}, 'sphere': {'scale':0.25}})
# # rotate for a better initial view
# viewer.rotate(-90, "y")
# viewer.rotate(-90, "x")
viewer.show()

You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol
Internal coordinate scan#
Let’s say that we want to see how the total energy varies when we change one internal coordinate and we optimize the rest. For this, we can set up a scan using geomeTRIC
# Set up the gradient and optimization dirvers:
grad_drv = vlx.ScfGradientDriver()
opt_drv = vlx.OptimizationDriver(grad_drv)
opt_drv.filename = "../../data/structopt/h2o_scan"
angle_scan = "angle 2 1 3 70 140.0 20"
opt_drv.constraints = ["$scan", angle_scan]
# Run the scan
opt_molecule = opt_drv.compute(molecule, basis, scf_drv)
Show code cell source
angle_scan = "angle 2 1 3 70 140.0 20"
def read_scan(file_name, natms, constraint):
# start, stop and how many steps
parts = constraint.split()
start = float(parts[-3])
stop = float(parts[-2])
no_steps = int(parts[-1])
step = (stop - start) / no_steps
scan_coords = []
energies = []
structures = []
for i in range(no_steps):
energy = 500
if i+1 < 10:
output_name = file_name + "_scan-00%d.xyz" % (i+1)
elif i+1 < 100:
output_name = file_name + "_scan-0%d.xyz" % (i+1)
else:
output_name = file_name + "_scan-%d.xyz" % (i+1)
output_file = open(output_name, "r")
# The optimized geometry is the last geometry in the file
# in xyz format.
nlines = natms + 2
lines = output_file.readlines()
opt_struct = ""
for line in lines[-nlines:]:
opt_struct += line
if "Energy" in line:
parts = line.split()
energy = float(parts[-1])
scan_coords.append(start + i*step)
energies.append(energy)
structures.append(opt_struct)
return scan_coords, energies, structures
natms = molecule.number_of_atoms()
scan_coords, energies, structures = read_scan("../../data/structopt/h2o_scan", natms, angle_scan)
# Plot the energies
plt.figure(figsize=(5,4))
plt.plot(scan_coords, energies,'o--')
plt.xlabel('Angle (deg)')
plt.ylabel('Energy (H)')
plt.title("H2O bond angle scan")
plt.tight_layout(); plt.show()

Constrained optimization#
Now let’s see how to run a constrained optimization. For this example, we will constrain the OH bond lengths to specific values and we will optimize the bond angle of the water molecule to construct a “relaxed” potential energy surface. For this, we will need several water structures with different combinations of OH bond lengths. These will be read from file.
# Since we are considering many combinations, running this cell may take a while...
bond_lengths = np.arange(0.8, 1.9, 0.1)
folder = "../../data/structopt/h2o_relaxed_pes/"
i = 0
for r1 in bond_lengths:
# We will constrain the first OH bond length to bond_1.
bond_1 = "distance 1 2 %.2f" % r1
for r2 in bond_lengths[i:]:
file_name = folder + "h2o_%.2f_%.2f" % (r1, r2)
file = open(file_name + ".xyz" , "r")
h2o_xyz = file.read()
file.close()
new_molecule = vlx.Molecule.read_xyz_string(h2o_xyz)
scf_results = scf_drv.compute(new_molecule, basis)
# Set up the gradient and optimization dirvers:
new_grad_drv = vlx.ScfGradientDriver()
new_opt_drv = vlx.OptimizationDriver(new_grad_drv)
new_opt_drv.filename = file_name
# We will constrain the second OH bond length to value bond_2.
bond_2 = "distance 1 3 %.2f" % r2
# Possible options: set, scan, freeze.
new_opt_drv.constraints = ["$set", bond_1, bond_2]
# We now optimize th bond angle:
h2o_opt = new_opt_drv.compute(new_molecule, basis, scf_drv)
i += 1
def read_optimized_structure(file_name, natms):
output_file = open(file_name, "r")
# The optimized geometry is the last geometry in the file
# in xyz format.
nlines = natms + 2
lines = output_file.readlines()
opt_struct = ""
for line in lines[-nlines:]:
opt_struct += line
if "Energy" in line:
parts = line.split()
energy = float(parts[-1])
output_file.close()
return energy, opt_struct
# Now let's plot the 2D PES and the 2D map of relaxed bond angles.
npts = bond_lengths.shape[0]
natm = molecule.number_of_atoms()
pes = np.zeros((npts, npts))
theta = np.zeros((npts, npts))
# Read energies
i = 0
for r1 in bond_lengths:
j = i
for r2 in bond_lengths[i:]:
file_name = folder + "h2o_%.2f_%.2f_optim.xyz" % (r1, r2)
energy, xyz_opt = read_optimized_structure(file_name, natm)
pes[i,j] = energy
pes[j,i] = energy
# calculate the relaxed bond angle
opt_molecule = vlx.Molecule.read_xyz_string(xyz_opt)
angle = geometric.internal.Angle(1, 0, 2)
coords = opt_molecule.get_coordinates().reshape((3 * natm))
# transform from radians to degrees
angle_value = (angle.value(coords))*180.0/np.pi
theta[i,j] = angle_value
theta[j,i] = angle_value
j += 1
i += 1
fig = plt.figure(figsize=(13, 4))
plt.subplot(1, 2, 1)
plt.xlabel('r1 (Å)')
plt.ylabel('r2 (Å)')
plt.title("H2O relaxed PES")
conts = plt.contour(bond_lengths, bond_lengths, pes, 10, colors='black')
plt.clabel(conts, inline=True, fontsize=8)
CS=plt.contourf(bond_lengths, bond_lengths, pes, cmap="RdGy",
alpha=0.75)
cbar=plt.colorbar(CS)
cbar.set_label("Energy (H)")
plt.subplot(1, 2, 2)
plt.xlabel('r1 (Å)')
plt.ylabel('r2 (Å)')
plt.title("H2O relaxed bond angle")
conts = plt.contour(bond_lengths, bond_lengths, theta, 10, colors='black')
plt.clabel(conts, inline=True, fontsize=8)
CS = plt.contourf(bond_lengths, bond_lengths, theta,
cmap="Blues_r", alpha=0.75)
cbar=plt.colorbar(CS)
cbar.set_label("Angle (rad)")
plt.show()

XTB geometry optimization#
In the following we will perform a ground state structure optimization of cafestol and kahweol. Cafestol is a molecule present in Robusta coffee beans, while kahweol is a very similar molecule present in Arabica coffee beans. Kahweol has a distinctive Raman peak which allows to identify Arabica coffee beans by Raman spectroscopy. Here, we determine the relaxed geometries of these two molecules. Since these molecules are large, we will use the semiempirical tight binding method from the Xtb code. Veloxchem has dedicated drivers which connect to Xtb.
cafestol_xyz = """51
cafestol initial structure
O -4.5215 -1.3455 0.4926
O -5.2905 0.0812 -1.7033
O 5.1630 0.4345 -0.1101
C -1.2807 -0.7511 -0.0465
C -0.5244 0.5810 -0.4102
C 0.9909 0.6610 0.0407
C -2.9440 0.2916 1.2866
C -1.8143 -0.7311 1.3945
C -2.6188 -0.8159 -0.8642
C -3.7116 -0.2387 0.0606
C -0.4332 -2.0024 -0.3283
C -1.3742 1.8303 -0.0654
C 1.6844 -0.6873 -0.4019
C -2.3206 1.6833 1.1309
C 0.9997 -1.9155 0.1913
C 1.7285 1.8344 -0.7277
C 1.1437 0.9261 1.5588
C 3.1644 -0.6155 -0.2085
C -4.6490 0.7442 -0.6282
C 3.2525 1.9378 -0.4590
C 3.8298 0.5879 -0.2856
C 4.1654 -1.5888 0.0140
C 5.3604 -0.9015 0.0689
H -0.4604 0.5817 -1.5108
H -3.5646 0.2626 2.1897
H -1.0861 -0.4449 2.1498
H -2.1902 -1.7145 1.7003
H -2.5286 -0.2495 -1.7981
H -2.8534 -1.8505 -1.1460
H -0.3943 -2.1694 -1.4136
H -0.9135 -2.8938 0.0946
H -1.9810 2.0650 -0.9487
H -0.7587 2.7222 0.0818
H 1.5605 -0.7646 -1.4951
H -1.7642 1.9154 2.0477
H -3.1082 2.4431 1.0635
H 1.0176 -1.9046 1.2859
H 1.5276 -2.8282 -0.1118
H 1.2808 2.8061 -0.4944
H 1.5851 1.6842 -1.8068
H 0.5898 1.8173 1.8670
H 0.8147 0.0892 2.1727
H 2.1832 1.0992 1.8532
H -4.1337 1.6158 -1.0351
H -5.4265 1.0934 0.0595
H 3.4273 2.5264 0.4487
H 3.7385 2.4646 -1.2871
H -3.9519 -2.0832 0.7635
H 4.0442 -2.6573 0.1180
H -5.7603 -0.6874 -1.3370
H 6.3869 -1.2031 0.2186
"""
cafestol = vlx.Molecule.read_xyz_string(cafestol_xyz)
kahweol_xyz = """49
kahweol initial structure
O -4.4710 -1.5361 0.0905
O -5.6188 0.9753 -0.0811
O 5.1625 0.5695 -0.0997
C -1.2209 -0.8317 -0.2017
C -0.4821 0.5357 -0.4703
C -3.0038 0.1167 1.0462
C -1.8514 -0.8737 1.1997
C -2.4967 -0.9058 -1.1119
C 0.9884 0.6333 0.0944
C -3.6676 -0.3924 -0.2474
C -0.3220 -2.0551 -0.4569
C -1.3919 1.7479 -0.1363
C -2.4167 1.5304 0.9816
C 1.7435 -0.6786 -0.3542
C 1.0741 -1.9483 0.1552
C 1.0367 0.8330 1.6293
C -4.5845 0.5870 -0.9677
C 1.7610 1.8511 -0.4907
C 3.2146 -0.5444 -0.1101
C 3.1013 1.8722 -0.6493
C 3.8167 0.6946 -0.3059
C 4.2309 -1.4546 0.2464
C 5.3980 -0.7259 0.2399
H -0.3415 0.5842 -1.5631
H -3.6779 0.0397 1.9072
H -1.1861 -0.6021 2.0159
H -2.2356 -1.8757 1.4322
H -2.3648 -0.3090 -2.0216
H -2.6844 -1.9384 -1.4328
H -0.2098 -2.1924 -1.5413
H -0.8043 -2.9694 -0.0890
H -1.9406 2.0033 -1.0514
H -0.8134 2.6482 0.0921
H -1.9341 1.7452 1.9434
H -3.2193 2.2709 0.8849
H 1.6699 -0.7251 -1.4553
H 1.0259 -1.9792 1.2481
H 1.6508 -2.8307 -0.1490
H 2.0514 1.0518 1.9842
H 0.4166 1.6762 1.9486
H 0.7233 -0.0489 2.1864
H -5.0524 0.1195 -1.8415
H -4.0679 1.4861 -1.3088
H 1.2088 2.7568 -0.7218
H -5.0925 -1.2746 0.7905
H 3.6240 2.7536 -0.9997
H 4.1386 -2.5085 0.4650
H -6.2132 0.2149 0.0339
H 6.4308 -0.9731 0.4374
"""
kahweol = vlx.Molecule.read_xyz_string(kahweol_xyz)
The difference between cafestol and kahweol is very subtle, essentially only one double bond and two H atoms.
print("\n (a) Cafestol (b) Kahweol")
viewer = p3d.view(viewergrid=(1,2),width=700,height=250, linked=True)
viewer.addModel(cafestol_xyz, 'xyz', viewer=(0,0))
viewer.addModel(kahweol_xyz, 'xyz', viewer=(0,1))
viewer.setViewStyle({"style": "outline", "width": 0.05})
viewer.setStyle({'stick': {}, 'sphere': {'scale':0.25}})
viewer.zoomTo()
viewer.show()
(a) Cafestol (b) Kahweol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol
# Set up the xtb driver
xtbdrv = vlx.XtbDriver()
xtbdrv.set_method('gfn2')
xtbdrv.compute(cafestol)
Show code cell output
...................................................
: SETUP :
:.................................................:
: # basis functions 120 :
: # atomic orbitals 120 :
: # shells 74 :
: # electrons 126 :
: max. iterations 280 :
: Hamiltonian GFN2-xTB :
: restarted? true :
: GBSA solvation false :
: PC potential false :
: electronic temp. 300.0000000 K :
: accuracy 1.0000000 :
: -> integral cutoff 0.2500000E+02 :
: -> integral neglect 0.1000000E-07 :
: -> SCF convergence 0.1000000E-05 Eh :
: -> wf. convergence 0.1000000E-03 e :
: Broyden damping 0.4000000 :
...................................................
iter E dE RMSdq gap omega full diag
1 -70.1256729 -0.701257E+02 0.308E+00 4.51 0.0 T
2 -70.1969288 -0.712559E-01 0.183E+00 4.43 1.0 T
3 -70.1926047 0.432407E-02 0.777E-01 4.43 1.0 T
4 -70.2001181 -0.751331E-02 0.340E-01 4.33 1.0 T
5 -70.2039272 -0.380913E-02 0.779E-02 4.36 1.0 T
6 -70.2039738 -0.466179E-04 0.373E-02 4.36 1.0 T
7 -70.2039865 -0.126975E-04 0.138E-02 4.36 1.0 T
8 -70.2039876 -0.106686E-05 0.338E-03 4.36 4.1 T
9 -70.2039877 -0.118719E-06 0.156E-03 4.36 9.0 T
10 -70.2039877 -0.125795E-07 0.589E-04 4.36 23.8 T
11 -70.2039877 -0.255976E-08 0.214E-04 4.36 65.4 T
*** convergence criteria satisfied after 11 iterations ***
# Occupation Energy/Eh Energy/eV
-------------------------------------------------------------
1 2.0000 -0.7268208 -19.7778
... ... ... ...
57 2.0000 -0.4145936 -11.2817
58 2.0000 -0.4096877 -11.1482
59 2.0000 -0.4060970 -11.0505
60 2.0000 -0.4022968 -10.9471
61 2.0000 -0.4009046 -10.9092
62 2.0000 -0.3984565 -10.8426
63 2.0000 -0.3664995 -9.9730 (HOMO)
64 -0.2062144 -5.6114 (LUMO)
65 -0.1513512 -4.1185
66 -0.0221661 -0.6032
67 0.0014760 0.0402
68 0.0144309 0.3927
... ... ...
120 0.6705777 18.2473
-------------------------------------------------------------
HL-Gap 0.1602851 Eh 4.3616 eV
Fermi-level -0.2863569 Eh -7.7922 eV
SCC (total) 0 d, 0 h, 0 min, 0.067 sec
SCC setup ... 0 min, 0.002 sec ( 2.814%)
Dispersion ... 0 min, 0.001 sec ( 1.547%)
classical contributions ... 0 min, 0.000 sec ( 0.737%)
integral evaluation ... 0 min, 0.008 sec ( 11.807%)
iterations ... 0 min, 0.034 sec ( 50.278%)
molecular gradient ... 0 min, 0.021 sec ( 31.757%)
printout ... 0 min, 0.001 sec ( 1.000%)
:::::::::::::::::::::::::::::::::::::::::::::::::::::
:: SUMMARY ::
:::::::::::::::::::::::::::::::::::::::::::::::::::::
:: total energy -69.361512202626 Eh ::
:: gradient norm 0.059453607714 Eh/a0 ::
:: HOMO-LUMO gap 4.361580250884 eV ::
::.................................................::
:: SCC energy -70.203987702583 Eh ::
:: -> isotropic ES 0.064491640234 Eh ::
:: -> anisotropic ES 0.010230341319 Eh ::
:: -> anisotropic XC 0.020200140987 Eh ::
:: -> dispersion -0.059211487310 Eh ::
:: repulsion energy 0.841319621733 Eh ::
:: add. restraining 0.000000000000 Eh ::
:: total charge 0.000000000000 e ::
:::::::::::::::::::::::::::::::::::::::::::::::::::::
Set up the xtb gradient and optimization drivers and run the geometry optimization#
xtb_grad_drv = vlx.XtbGradientDriver()
xtb_opt_drv = vlx.OptimizationDriver(xtb_grad_drv)
xtb_opt_cafestol = xtb_opt_drv.compute(cafestol, xtbdrv)
Show code cell output
Optimization Driver Setup
===========================
Coordinate System : TRIC
Constraints : No
Max. Number of Steps : 300
Transition State : No
Hessian : never
* Info * Computing energy and gradient...
* Info * Energy : -69.3615122035 a.u.
* Info * Gradient : 8.325322e-03 a.u. (RMS)
* Info * 1.928454e-02 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3710058836 a.u.
* Info * Gradient : 3.923981e-03 a.u. (RMS)
* Info * 9.301846e-03 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3722458724 a.u.
* Info * Gradient : 2.159347e-03 a.u. (RMS)
* Info * 7.974968e-03 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3726427963 a.u.
* Info * Gradient : 1.107748e-03 a.u. (RMS)
* Info * 4.002520e-03 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3728625818 a.u.
* Info * Gradient : 5.724391e-04 a.u. (RMS)
* Info * 1.901848e-03 a.u. (Max)
* Info * Time : 0.11 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3729577548 a.u.
* Info * Gradient : 3.816103e-04 a.u. (RMS)
* Info * 9.881626e-04 a.u. (Max)
* Info * Time : 0.13 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730321297 a.u.
* Info * Gradient : 3.155273e-04 a.u. (RMS)
* Info * 8.318735e-04 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730679685 a.u.
* Info * Gradient : 2.001910e-04 a.u. (RMS)
* Info * 5.533687e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730753850 a.u.
* Info * Gradient : 1.166858e-04 a.u. (RMS)
* Info * 2.372876e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730804351 a.u.
* Info * Gradient : 6.970578e-05 a.u. (RMS)
* Info * 1.798257e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730830365 a.u.
* Info * Gradient : 6.708893e-05 a.u. (RMS)
* Info * 1.927371e-04 a.u. (Max)
* Info * Time : 0.08 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730839403 a.u.
* Info * Gradient : 4.328344e-05 a.u. (RMS)
* Info * 1.301321e-04 a.u. (Max)
* Info * Time : 0.14 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730847891 a.u.
* Info * Gradient : 3.104179e-05 a.u. (RMS)
* Info * 8.742235e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730852597 a.u.
* Info * Gradient : 3.926610e-05 a.u. (RMS)
* Info * 1.113513e-04 a.u. (Max)
* Info * Time : 0.11 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730857417 a.u.
* Info * Gradient : 3.673542e-05 a.u. (RMS)
* Info * 1.127451e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730861287 a.u.
* Info * Gradient : 2.352441e-05 a.u. (RMS)
* Info * 6.344279e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730864008 a.u.
* Info * Gradient : 2.201488e-05 a.u. (RMS)
* Info * 6.665192e-05 a.u. (Max)
* Info * Time : 0.08 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730866063 a.u.
* Info * Gradient : 2.456196e-05 a.u. (RMS)
* Info * 6.707557e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730867353 a.u.
* Info * Gradient : 1.788210e-05 a.u. (RMS)
* Info * 5.657216e-05 a.u. (Max)
* Info * Time : 0.08 sec
* Info * Computing energy and gradient...
* Info * Energy : -69.3730867912 a.u.
* Info * Gradient : 8.285115e-06 a.u. (RMS)
* Info * 2.364064e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Geometry optimization completed.
Molecular Geometry (Angstroms)
================================
Atom Coordinate X Coordinate Y Coordinate Z
O -4.586927894298 -1.279900580081 0.480283797793
O -5.319303242952 0.246562207320 -1.629047843300
O 5.096372964454 0.479142524541 0.047370527290
C -1.255190295744 -0.766697135858 -0.090694519716
C -0.509228219094 0.539152695692 -0.458496889438
C 0.978372037194 0.632243366820 -0.011990155161
C -2.880159755562 0.238836828861 1.312497729308
C -1.740436945819 -0.771913903412 1.364288546093
C -2.598892938446 -0.799816206160 -0.863851939900
C -3.682490893188 -0.248170686858 0.088636802757
C -0.416164506854 -1.999830242086 -0.434968241726
C -1.378410399217 1.760376004456 -0.090411193298
C 1.674856502875 -0.681263507644 -0.445474818245
C -2.250016721329 1.621564634949 1.161318676793
C 0.997770067356 -1.919050968336 0.128039543449
C 1.666376726005 1.766537760643 -0.811768240558
C 1.168837532826 0.910140988102 1.479849080800
C 3.135520466851 -0.569549508958 -0.157646718636
C -4.609457056370 0.795536374054 -0.551570433613
C 3.155683415864 1.954159444441 -0.478173900174
C 3.782186986193 0.630199905389 -0.219299021035
C 4.130802164950 -1.530150373797 0.161772362564
C 5.294632956634 -0.835297374349 0.274041555784
H -0.446342259653 0.529017873163 -1.555044787035
H -3.512301962729 0.221149517725 2.205105205326
H -0.987805192361 -0.489646713079 2.090189238482
H -2.099466069880 -1.770583425745 1.629844780980
H -2.561113930095 -0.227435609543 -1.788326352210
H -2.854655290972 -1.827144938935 -1.131680096908
H -0.352920373720 -2.090381541625 -1.522829914942
H -0.917392216442 -2.895613171681 -0.058907062896
H -2.030547320612 1.964782326180 -0.939715067955
H -0.751054032804 2.643276309973 0.031133878813
H 1.568213879475 -0.751868233951 -1.538315305172
H -1.648526534108 1.799632361004 2.055430463497
H -3.024936836985 2.390210621298 1.135414632642
H 0.980001478774 -1.887152726443 1.216951591943
H 1.555474968493 -2.809090968393 -0.168723212483
H 1.159647477215 2.715758648017 -0.638991630812
H 1.574182496031 1.537824463169 -1.875605731830
H 0.649837420195 1.820869381957 1.764759763362
H 0.808399132338 0.095444063030 2.096526843634
H 2.222626553193 1.033900090352 1.709265075873
H -4.054650314922 1.641880113076 -0.952434821920
H -5.302113708996 1.156095586962 0.223473746127
H 3.278594509806 2.596762656843 0.399147629118
H 3.657824474030 2.457385896016 -1.308818819672
H -4.087273901127 -1.977207925165 0.916416235667
H 4.002114599241 -2.584672649985 0.293398315288
H -5.695005406483 -0.585384959153 -1.313403501685
H 6.294955410769 -1.145419292801 0.501034196932
Summary of Geometry Optimization
==================================
Opt.Step Energy (a.u.) Energy Change (a.u.) Displacement (RMS, Max)
-------------------------------------------------------------------------------------
0 -69.361512203479 0.000000000000 0.000e+00 0.000e+00
1 -69.371005883633 -0.009493680153 9.371e-02 1.764e-01
2 -69.372245872446 -0.001239988813 3.111e-02 6.262e-02
3 -69.372642796341 -0.000396923895 2.277e-02 5.494e-02
4 -69.372862581760 -0.000219785419 1.642e-02 3.635e-02
5 -69.372957754801 -0.000095173041 1.237e-02 2.656e-02
6 -69.373032129673 -0.000074374871 1.314e-02 3.166e-02
7 -69.373067968527 -0.000035838854 1.261e-02 3.285e-02
8 -69.373075384965 -0.000007416438 3.850e-03 1.127e-02
9 -69.373080435105 -0.000005050141 2.618e-03 6.103e-03
10 -69.373083036530 -0.000002601425 3.036e-03 7.661e-03
11 -69.373083940347 -0.000000903817 1.273e-03 2.686e-03
12 -69.373084789081 -0.000000848734 1.658e-03 4.487e-03
13 -69.373085259683 -0.000000470602 9.193e-04 3.155e-03
14 -69.373085741656 -0.000000481973 1.153e-03 3.977e-03
15 -69.373086128723 -0.000000387067 1.139e-03 4.015e-03
16 -69.373086400847 -0.000000272125 9.826e-04 3.297e-03
17 -69.373086606303 -0.000000205455 8.375e-04 2.820e-03
18 -69.373086735312 -0.000000129009 6.423e-04 1.851e-03
19 -69.373086791229 -0.000000055918 3.349e-04 7.964e-04
Statistical Deviation between
Optimized Geometry and Initial Geometry
=========================================
Internal Coord. RMS deviation Max. deviation
-----------------------------------------------------------
Bonds 0.011 Angstrom 0.036 Angstrom
Angles 1.350 degree 3.923 degree
Dihedrals 3.195 degree 13.585 degree
*** Time spent in Optimization Driver: 10.74 sec
Run the XTB structure optimization for kahweol:#
xtbdrv = vlx.XtbDriver()
xtbdrv.set_method('gfn2')
xtbdrv.compute(kahweol)
Show code cell output
...................................................
: SETUP :
:.................................................:
: # basis functions 118 :
: # atomic orbitals 118 :
: # shells 72 :
: # electrons 124 :
: max. iterations 280 :
: Hamiltonian GFN2-xTB :
: restarted? true :
: GBSA solvation false :
: PC potential false :
: electronic temp. 300.0000000 K :
: accuracy 1.0000000 :
: -> integral cutoff 0.2500000E+02 :
: -> integral neglect 0.1000000E-07 :
: -> SCF convergence 0.1000000E-05 Eh :
: -> wf. convergence 0.1000000E-03 e :
: Broyden damping 0.4000000 :
...................................................
iter E dE RMSdq gap omega full diag
1 -69.0945256 -0.690945E+02 0.335E+00 3.05 0.0 T
2 -69.1552391 -0.607135E-01 0.197E+00 3.02 1.0 T
3 -69.1496721 0.556704E-02 0.776E-01 3.01 1.0 T
4 -69.1587260 -0.905395E-02 0.308E-01 2.99 1.0 T
5 -69.1617692 -0.304311E-02 0.807E-02 3.00 1.0 T
6 -69.1618210 -0.518211E-04 0.439E-02 3.00 1.0 T
7 -69.1618464 -0.254646E-04 0.237E-02 3.00 1.0 T
8 -69.1618577 -0.112592E-04 0.101E-02 3.00 1.4 T
9 -69.1618597 -0.199642E-05 0.261E-03 3.00 5.5 T
10 -69.1618598 -0.677972E-07 0.144E-03 3.00 9.9 T
11 -69.1618598 -0.38
1184E-07 0.545E-04 3.00 26.2 T
12 -69.1618598 -0.377007E-08 0.158E-04 3.00 90.1 T
*** convergence criteria satisfied after 12 iterations ***
# Occupation Energy/Eh Energy/eV
-------------------------------------------------------------
1 2.0000 -0.7255326 -19.7427
... ... ... ...
56 2.0000 -0.4203811 -11.4392
57 2.0000 -0.4084609 -11.1148
58 2.0000 -0.4076528 -11.0928
59 2.0000 -0.4020316 -10.9398
60 2.0000 -0.3988397 -10.8530
61 2.0000 -0.3963877 -10.7863
62 2.0000 -0.3564174 -9.6986 (HOMO)
63 -0.2461602 -6.6984 (LUMO)
64
-0.1788396 -4.8665
65 -0.1267271 -3.4484
66 -0.0099203 -0.2699
67 0.0051269 0.1395
... ... ...
118 0.6709176 18.2566
-------------------------------------------------------------
HL-Gap 0.1102572 Eh 3.0003 eV
Fermi-level -0.3012888 Eh -8.1985 eV
SCC (total) 0 d, 0 h, 0 min, 0.061 sec
SCC setup ... 0 min, 0.000 sec ( 0.440%)
Dispersion ... 0 min, 0.001 sec ( 0.891%)
classical contributions ... 0 min, 0.000 sec ( 0.175%)
integral evaluation ... 0 min, 0.005 sec ( 8.277%)
iterations ... 0 min, 0.034 sec ( 55.436%)
molecular gradient ... 0 min, 0.021 sec ( 34.266%)
printout ... 0 min, 0.000 sec ( 0.465%)
:::::::::::::::::::::::::::::::::::::::::::::::::::::
:: SUMMARY ::
:::::::::::::::::::::::::::::::::::::::::::::::::::::
:: total energy -68.310234156957 Eh ::
:: gradient norm 0.059015972263 Eh/a0 ::
:: HOMO-LUMO gap 3.000251409403 eV ::
::.................................................::
:: SCC energy -69.161859807846 Eh ::
:: -> isotropic ES 0.064055799707 Eh ::
:: -> anisotropic ES 0.010957972238 Eh ::
:: -> anisotropic XC 0.020993529454 Eh ::
:: -> dispersion -0.057983897157 Eh ::
:: repulsion energy 0.850629559592 Eh ::
:: add. restraining 0.000000000000 Eh ::
:: total charge -0.000000000000 e ::
:::::::::::::::::::::::::::::::::::::::::::::::::::::
xtb_grad_drv = vlx.XtbGradientDriver()
xtb_opt_drv = vlx.OptimizationDriver(xtb_grad_drv)
xtb_opt_kahweol = xtb_opt_drv.compute(kahweol, xtbdrv)
Show code cell output
Optimization Driver Setup
===========================
Coordinate System : TRIC
Constraints : No
Max. Number of Steps : 300
Transition State : No
Hessian : never
* Info * Computing energy and gradient...
* Info * Energy : -68.3102341572 a.u.
* Info * Gradient : 8.430888e-03 a.u. (RMS)
* Info * 3.152910e-02 a.u. (Max)
* Info * Time : 0.11 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3194377305 a.u.
* Info * Gradient : 3.788613e-03 a.u. (RMS)
* Info * 1.026022e-02 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3207078420 a.u.
* Info * Gradient : 1.696187e-03 a.u. (RMS)
* Info * 5.384314e-03 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3211017988 a.u.
* Info * Gradient : 8.702338e-04 a.u. (RMS)
* Info * 2.453755e-03 a.u. (Max)
* Info * Time : 0.11 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3212623136 a.u.
* Info * Gradient : 7.160991e-04 a.u. (RMS)
* Info * 1.889567e-03 a.u. (Max)
* Info * Time : 0.13 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3213472035 a.u.
* Info * Gradient : 6.526878e-04 a.u. (RMS)
* Info * 1.484448e-03 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3214034870 a.u.
* Info * Gradient : 3.020675e-04 a.u. (RMS)
* Info * 7.637624e-04 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3214290429 a.u.
* Info * Gradient : 1.803308e-04 a.u. (RMS)
* Info * 4.186287e-04 a.u. (Max)
* Info * Time : 0.10 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3214477989 a.u.
* Info * Gradient : 2.557944e-04 a.u. (RMS)
* Info * 6.645090e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3214677700 a.u.
* Info * Gradient : 2.352957e-04 a.u. (RMS)
* Info * 5.885042e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3214890583 a.u.
* Info * Gradient : 1.667861e-04 a.u. (RMS)
* Info * 3.682147e-04 a.u. (Max)
* Info * Time : 0.11 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3214960720 a.u.
* Info * Gradient : 1.055778e-04 a.u. (RMS)
* Info * 2.130604e-04 a.u. (Max)
* Info * Time : 0.08 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215007939 a.u.
* Info * Gradient : 8.680498e-05 a.u. (RMS)
* Info * 2.098592e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215030976 a.u.
* Info * Gradient : 8.548065e-05 a.u. (RMS)
* Info * 1.860242e-04 a.u. (Max)
* Info * Time : 0.12 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215045601 a.u.
* Info * Gradient : 5.725651e-05 a.u. (RMS)
* Info * 1.466322e-04 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215052972 a.u.
* Info * Gradient : 2.969516e-05 a.u. (RMS)
* Info * 5.772642e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215056017 a.u.
* Info * Gradient : 2.610916e-05 a.u. (RMS)
* Info * 5.841689e-05 a.u. (Max)
* Info * Time : 0.08 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215058954 a.u.
* Info * Gradient : 2.939108e-05 a.u. (RMS)
* Info * 8.178844e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215061317 a.u.
* Info * Gradient : 2.819444e-05 a.u. (RMS)
* Info * 8.585832e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215063064 a.u.
* Info * Gradient : 1.915187e-05 a.u. (RMS)
* Info * 4.266796e-05 a.u. (Max)
* Info * Time : 0.09 sec
* Info * Computing energy and gradient...
* Info * Energy : -68.3215063782 a.u.
* Info * Gradient : 1.020411e-05 a.u. (RMS)
* Info * 2.548941e-05 a.u. (Max)
* Info * Time : 0.08 sec
* Info * Geometry optimization completed.
Molecular Geometry (Angstroms)
================================
Atom Coordinate X Coordinate Y Coordinate Z
O -4.444246777041 -1.519500184760 0.133687531696
O -5.717823367908 0.941920247192 -0.067035378737
O 5.074906471664 0.622193623988 0.095473632113
C -1.218602072799 -0.830635685957 -0.247269720612
C -0.477652350991 0.497206859079 -0.529645873242
C -2.930032992048 0.109193794340 1.091620378561
C -1.787380235852 -0.895060417667 1.174536822172
C -2.515633689930 -0.856030106348 -1.095108430247
C 0.968427401326 0.594486819654 0.023105350455
C -3.654876617179 -0.372392602218 -0.180458472359
C -0.338761173834 -2.036945429660 -0.583487647170
C -1.380922979880 1.698249048613 -0.180015868996
C -2.322218578481 1.506829862993 1.012213415312
C 1.723651132776 -0.686675316567 -0.399653676283
C 1.033746130303 -1.956638037485 0.075497470977
C 1.068775515460 0.833466747697 1.537760243536
C -4.585076305760 0.655631548432 -0.853891425080
C 1.691792477249 1.770178546276 -0.620456373727
C 3.162384674665 -0.539524931275 -0.024308835986
C 3.026453417345 1.832993169416 -0.661213458299
C 3.761356117396 0.686850100864 -0.200988462347
C 4.186105772490 -1.411077541156 0.402955911206
C 5.318402057353 -0.648904554007 0.464791088959
H -0.335044283025 0.523188240595 -1.619121935637
H -3.599210157051 0.065664226708 1.957648973711
H -1.080844988373 -0.645952746940 1.956895670602
H -2.177196951506 -1.898085777230 1.358665298029
H -2.433748208795 -0.246586339554 -1.993187408162
H -2.750555988910 -1.876641394325 -1.397583155030
H -0.202960809259 -2.084227408498 -1.667524196228
H -0.851689802471 -2.949995114334 -0.273727664073
H -1.983529167451 1.923955505205 -1.061256403684
H -0.768118578497 2.580194152344 0.010452737792
H -1.779681306265 1.670316886081 1.945888113628
H -3.112544172753 2.258627623935 0.965687394824
H 1.702617585291 -0.706788044839 -1.502641974234
H 0.941020366065 -1.968972909219 1.160844552883
H 1.629998855522 -2.823282549407 -0.216933715679
H 2.103646971685 1.037195455404 1.803383367123
H 0.473056818750 1.693597069299 1.827356187713
H 0.749260115453 -0.028531398807 2.112178950134
H -4.890278846740 0.251727445199 -1.829217855767
H -4.081465688582 1.606980030046 -1.012092786443
H 1.101025022226 2.598128915550 -0.978016602208
H -4.939041486337 -1.335806703585 0.940565644080
H 3.560854657454 2.688860507163 -1.042623805224
H 4.105797619339 -2.451938367526 0.640845591965
H -6.244852481130 0.135937465410 -0.009556359705
H 6.323310879036 -0.889980330118 0.749463157684
Summary of Geometry Optimization
==================================
Opt.Step Energy (a.u.) Energy Change (a.u.) Displacement (RMS, Max)
-------------------------------------------------------------------------------------
0 -68.310234157198 0.000000000000 0.000e+00 0.000e+00
1 -68.319437730465 -0.009203573267 9.448e-02 2.163e-01
2 -68.320707842027 -0.001270111562 3.422e-02 7.342e-02
3 -68.321101798833 -0.000393956806 3.393e-02 7.123e-02
4 -68.321262313555 -0.000160514722 1.908e-02 4.657e-02
5 -68.321347203517 -0.000084889962 2.115e-02 5.549e-02
6 -68.321403486976 -0.000056283459 6.872e-03 1.656e-02
7 -68.321429042901 -0.000025555925 7.176e-03 2.097e-02
8 -68.321447798926 -0.000018756025 7.326e-03 2.508e-02
9 -68.321467769997 -0.000019971071 7.674e-03 3.013e-02
10 -68.321489058336 -0.000021288339 1.244e-02 5.091e-02
11 -68.321496071988 -0.000007013652 5.918e-03 2.461e-02
12 -68.321500793911 -0.000004721923 5.777e-03 2.325e-02
13 -68.321503097569 -0.000002303658 2.600e-03 9.960e-03
14 -68.321504560138 -0.000001462569 1.516e-03 4.717e-03
15 -68.321505297248 -0.000000737109 1.157e-03 3.325e-03
16 -68.321505601721 -0.000000304473 8.243e-04 3.052e-03
17 -68.321505895430 -0.000000293709 1.093e-03 4.356e-03
18 -68.321506131708 -0.000000236278 9.455e-04 3.618e-03
19 -68.321506306422 -0.000000174713 7.641e-04 2.729e-03
20 -68.321506378185 -0.000000071763 3.670e-04 1.459e-03
Statistical Deviation between
Optimized Geometry and Initial Geometry
=========================================
Internal Coord. RMS deviation Max. deviation
-----------------------------------------------------------
Bonds 0.012 Angstrom 0.033 Angstrom
Angles 1.359 degree 3.480 degree
Dihedrals 3.908 degree 13.107 degree
*** Time spent in Optimization Driver: 11.39 sec
Compare initial and final geometries#
xyz_xtb_opt_cafestol = xtb_opt_cafestol.get_xyz_string()
xyz_xtb_opt_kahweol = xtb_opt_kahweol.get_xyz_string()
print('\n (a) Cafestol, initial vs. XTB opt. geom. (b) Kahweol, initial vs. XTB opt. geom.:')
viewer = p3d.view(viewergrid=(1,2),width=700,height=250, linked=True)
viewer.addModel(cafestol_xyz, 'xyz', viewer=(0,0))
viewer.addModel(xyz_xtb_opt_cafestol, 'xyz', viewer=(0,0))
viewer.addModel(kahweol_xyz, 'xyz', viewer=(0,1))
viewer.addModel(xyz_xtb_opt_kahweol, 'xyz', viewer=(0,1))
viewer.setViewStyle({"style": "outline", "width": 0.05})
viewer.setStyle({'stick': {}, 'sphere': {'scale':0.25}})
viewer.zoomTo()
viewer.show()
(a) Cafestol, initial vs. XTB opt. geom. (b) Kahweol, initial vs. XTB opt. geom.:
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol