RESP charges#

The restrained electrostatic potential (RESP) charge model [BCCK93, CCBK95] represents an improvement to the Merz–Kollman (MK) scheme as the ESP figure-of-merit is rather insensitive to variations in charges of atoms buried inside the molecule.

../../_images/chi_square.svg

Fig. 9 Dependence of figure-of-merit, \(\chi^2_\mathrm{esp}\), with respect to variations in atomic charges. Four separate atoms are here considered.#

To avoid unphysically high magnitudes of the charges of interior atoms, a hyperbolic penalty function is added

\[ \chi_{\mathrm{rstr}}^2 = \alpha \sum_{A=1}^M \big[ (q_A^2 + \beta^2)^{1/2} - \beta \big] \]

so that the diagonal matrix elements of the \(A\)-matrix in MK scheme become equal to

\[ A_{JJ} = \frac{1}{4\pi\varepsilon_0} \sum_{a} \frac{1}{r_{aJ}^2} + \alpha \, (q_J^2+\beta^2)^{-1/2} \]

with a dependency on the partial charge. Consequently, RESP charges are obtained by solving the matrix equation iteratively until the charges and Lagrange multipliers become self-consistent. In addition to that, the RESP charge model allows for the introduction of constraints on charges of equivalent atoms due to symmetry operations or bond rotations.

Let us determine the RESP charges for methanol at the Hartree–Fock/6-31G(d) level of theory. We will require the partial charges of the hydrogen atoms in the methyl group to be identical.

Hide code cell source
import py3Dmol as p3d

methanol_xyz = """6

  H      1.2001      0.0363      0.8431
  C      0.7031      0.0083     -0.1305
  H      0.9877      0.8943     -0.7114
  H      1.0155     -0.8918     -0.6742
  O     -0.6582     -0.0067      0.1730
  H     -1.1326     -0.0311     -0.6482
"""

viewer = p3d.view(width=300, height=200)

viewer.addModel(methanol_xyz, "xyz")
viewer.setStyle({"stick": {}})

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

import veloxchem as vlx

molecule = vlx.Molecule.read_xyz_string(methanol_xyz)
basis = vlx.MolecularBasis.read(molecule, "6-31G*")

A calculation of RESP charges is performed with the compute method of the RespChargesDriver class.

resp_drv = vlx.RespChargesDriver()

resp_drv.update_settings({"equal_charges": "1 = 3, 1 = 4"})

resp_charges = resp_drv.compute(molecule, basis, "resp")
print("Atom   RESP charge")
print(20 * "-")

for label, resp_charge in zip(molecule.get_labels(), resp_charges):

    print(f"{label :s} {resp_charge : 18.6f}")

print(20 * "-")

print(f"Total: {resp_charges.sum() : 13.6f}")
Atom   RESP charge
--------------------
H           0.033747
C           0.118610
H           0.033747
H           0.033747
O          -0.639004
H           0.419154
--------------------
Total:      0.000000

It is noted that the partial charge of the buried carbon atom (0.12 a.u.) is significantly lower than the corresponding ESP charge (0.20 a.u.).

Note

The RESP charge model is based on the quantum mechanical electrostatic potential calculated at the Hartree–Fock/6-31G(d) level of theory because it overestimates the gas-phase dipole moments in a way that it imitates approximately polarization effects in aqueous protein environments.