coul_gauss

mdhelper.openmm.pair.coul_gauss(cutoff: float | Quantity, tol: float = 0.0001, *, g_ewald: float | Quantity = None, dims: ndarray[float] | Quantity = None, mix: str = 'default', per_params: list = None, global_params: dict[str, float | Quantity] = None, tab_funcs: dict[str, ndarray[int | float] | Quantity | Discrete2DFunction] = None) tuple[CustomNonbondedForce, NonbondedForce][source]

Implements the smeared Coulomb pair potential.

The charges have spherical Gaussian-distributed charge distributions [1] [2]

\[u_\mathrm{Coul}(r_{12})=\frac{q_1q_2}{4\pi\varepsilon_0r_{12}} \mathrm{erf}(\alpha_{12}r_{12})\]

where \(q_1\) and \(q_2\) are the particle charges in \(\mathrm{e}\), \(\varepsilon_0\) is the vacuum permittivity, and \(\alpha_{12}=\sqrt{\alpha_1^2+\alpha_2^2}\) is an inverse screening length or damping parameter in \(\mathrm{nm}^{-1}\). Effectively, electrostatic interactions are dampened for charges at small separation distances but remain unchanged at large separation distances.

Additionally, an implementation of the electrostatic pair potential in the Gaussian core model

\[u_\mathrm{Coul}(r_{12})=\frac{q_1q_2}{4\pi\varepsilon_0r_{12}} \mathrm{erf}\left(\frac{\pi^{1/2}}{2^{1/2}a_{12}}r_{12}\right)\]

is available, where \(a_{12}=\sqrt{a_1^2+a_2^2}\) with \(a_1\) and \(a_2\) being the electrostatic smearing radii in \(\mathrm{nm}\).

To account for the solvent polarization implicitly using its relative permittivity \(\varepsilon_\mathrm{r}\), scale the particle charges \(q_1\) and \(q_2\) by \(1/\sqrt{\varepsilon_\mathrm{r}}\).

After creating the pair potentials, particles should be registered using openmm.openmm.NonbondedForce.addParticle() and openmm.openmm.CustomNonbondedForce.addParticle().

Parameters:
cutofffloat or openmm.unit.Quantity

Shared cutoff distance for all nonbonded interactions in the simulation sytem. Must be less than half the minimum periodic simulation box dimension.

Reference unit: \(\mathrm{nm}\).

tolfloat, default: 1e-4

Error tolerance \(\delta\) for Ewald summation.

g_ewaldfloat or openmm.unit.Quantity, keyword-only, optional

Ewald splitting parameter \(g_\mathrm{Ewald}\). If not provided, it is computed using

\[g_\mathrm{Ewald}=\frac{\sqrt{-\log(2\delta)}}{r_\mathrm{cut}}\]

Reference unit: \(\mathrm{nm}^{-1}\).

dimsarray_like or openmm.unit.Quantity, keyword-only, optional

Simulation system dimensions. Must be provided with g_ewald to calculate the number of mesh nodes \(n_\mathrm{mesh}\). Both dims and g_ewald must either have units or be unitless.

Shape: \((3,)\).

Reference unit: \(\mathrm{nm}\).

mixstr, keyword-only, default: "default"

Mixing rule for \(\alpha_{12}\).

Valid values:

  • "default": Default mixing rule.

    \[\alpha_{12}=\frac{\alpha_1\alpha_2} {\sqrt{\alpha_1^2+\alpha_2^2}}\]

    Per-particle parameters:

  • "core": Gaussian core model.

    \[\begin{split}\begin{gather*} a_{12}^2=a_1^2+a_2^2\\ \alpha_{12}=\sqrt{\frac{\pi}{2a_{12}^2}} \end{gather*}\end{split}\]

    This is equivalent to setting \(\alpha_i=\sqrt{\pi/(2a_i^2)}\) in the default mixing rule.

  • "alpha12 = ...;": Custom mixing rule. The string containing the expression for \(\alpha_{12}\) must be written in valid C++ syntax, with any custom global and per-particle parameters and tabulated functions defined in global_params, per_params, and tab_funcs, respectively.

    Per-particle parameters:

Tip

To disable the Lennard-Jones potential, set \(\sigma_i=0\,\mathrm{nm}\) and \(\epsilon_i=0\,\mathrm{kJ/mol}\) for all particles.

global_paramsdict, keyword-only, optional

Additional global parameters for use in the definition of \(\alpha_{12}\).

per_paramslist, keyword-only, optional

Additional per-particle parameters for use in the definition of \(\alpha_{12}\).

tab_funcsdict, keyword-only, optional

Optional tabulated functions for use in the definition of \(\alpha_{12}\).

Returns:
pair_coul_gauss_diropenmm.CustomNonbondedForce

Short-range electrostatic contribution evaluated in real space.

pair_coul_gauss_recopenmm.NonbondedForce

Long-range electrostatic contribution evaluated in Fourier (reciprocal) space.

References

[1]

Kiss, P. T.; Sega, M.; Baranyai, A. Efficient Handling of Gaussian Charge Distributions: An Application to Polarizable Molecular Models. J. Chem. Theory Comput. 2014, 10 (12), 5513–5519. https://doi.org/10.1021/ct5009069.

[2]

Eslami, H.; Khani, M.; Müller-Plathe, F. Gaussian Charge Distributions for Incorporation of Electrostatic Interactions in Dissipative Particle Dynamics: Application to Self-Assembly of Surfactants. J. Chem. Theory Comput. 2019, 15 (7), 4197–4207. https://doi.org/10.1021/acs.jctc.9b00174.