ufjc.potential
A module for handling potentials.
This module contains several different classes representing potentials,
each having methods to compute relevant nondimensional quantities as
functions of nondimensional force or stretch.
This module also contains the parent class Potential
that is used to
assign a potential to a given model using keyword arguments.
Examples
Create a Lennard-Jones potential model with a nondimensional potential energy scale of 8 and evaluate the nondimensional potential energy at a stretch of 1.23:
>>> from ufjc.potential import LennardJonesPotential
>>> model = LennardJonesPotential(varepsilon=8)
>>> model.beta_u(1.23)
4.046654314368616
Do the same with the Lenard-Jones-FENE potential:
>>> from ufjc.potential import LJFENEPotential
>>> model = LJFENEPotential(varepsilon=(8, 8))
>>> model.beta_u(1.23)
8.510502022381505
Create a single-link model in one dimension, instantiate it with the Morse potential, and compute the incremental link stretch under a nondimensional force of 8:
>>> from ufjc.potential import Potential
>>> class Link1D(Potential):
... def __init__(self, **kwargs):
... Potential.__init__(self, **kwargs)
>>> _ = Link1D(potential='morse').delta_lambda(8)
>>> Link1D(potential='lj-fene').eta_link(1)
184.0
- class CustomPotential(**kwargs)[source]
Bases:
object
A custom user-defined potential.
- potential
The potential name.
- Type:
str
- varepsilon
The nondimensional energy scale.
- Type:
float
- phi
The scaled nondimensional energy function.
- Type:
function
- eta_link
The nondimensional force as a function of stretch.
- Type:
function
- delta_lambda
The incremental stretch as a function of the nondimensional force (optional).
- Type:
function
- kappa
The nondimensional stiffness.
- Type:
float
- c
The correction parameter.
- Type:
float
- class HarmonicPotential(**kwargs)[source]
Bases:
object
The harmonic potential.
- varepsilon
The nondimensional energy scale.
- Type:
float
- kappa
The nondimensional stiffness \(\kappa=\varepsilon\).
- Type:
float
- c
The correction parameter \(c=1\).
- Type:
float
- beta_u(lambda_)[source]
The nondimensional potential energy function,
\[\beta u(\lambda) = \varepsilon\phi(\lambda) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional potential energy(s).
- Return type:
numpy.ndarray
- delta_lambda(eta)[source]
The incremental stretch as a function of nondimensional force,
\[\Delta\lambda(\eta) = \frac{\eta}{\varepsilon} .\]- Parameters:
eta (array_like) – The nondimensional force(s).
- Returns:
The incremental stretch(s).
- Return type:
numpy.ndarray
- eta_link(lambda_)[source]
The nondimensional force as a function of stretch,
\[\eta(\lambda) = \varepsilon(\lambda - 1) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional force(s).
- Return type:
numpy.ndarray
Example
Compute the nondimensional force at a sample stretch:
>>> from ufjc.potential import HarmonicPotential >>> HarmonicPotential().eta_link(1.8) 70.4
- class LJFENEPotential(**kwargs)[source]
Bases:
object
The Lennard-Jones-FENE potential [2].
- varepsilon
The nondimensional energy scale.
- Type:
float
- kappa
The nondimensional stiffness.
- Type:
float
- c
The correction parameter.
- Type:
float
- eta_max
The maximum nondimensional force \(\eta_\mathrm{max} = \eta(\lambda_\mathrm{max})\).
- Type:
float
- lambda_max
The stretch at the maximum nondimensional force.
- Type:
float
- beta_u(lambda_)[source]
The nondimensional potential energy function,
\[\beta u(\lambda) = \varepsilon\phi(\lambda) = \varepsilon_1\left( \frac{1}{\lambda^{12}} - \frac{2}{\lambda^6} + 1 \right) - \frac{\varepsilon_2}{2}\,\ln\left[ 1 - \left(\frac{\lambda}{\lambda_\mathrm{max}}\right)^2 \right] .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional potential energy(s).
- Return type:
numpy.ndarray
- eta_link(lambda_)[source]
The nondimensional force as a function of stretch,
\[\eta(\lambda) = 12\varepsilon_1\left( \frac{1}{\lambda^7} - \frac{1}{\lambda^{13}} \right) + \frac{\varepsilon_2\lambda} {\lambda_\mathrm{max}^2 - \lambda^2} .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional force(s).
- Return type:
numpy.ndarray
- phi(lambda_)[source]
The scaled nondimensional potential energy function,
\[\phi(\lambda) = \frac{\varepsilon_1}{\varepsilon_2}\left( \frac{1}{\lambda^{12}} - \frac{2}{\lambda^6} + 1 \right) - \frac{1}{2}\,\ln\left[ 1 - \left(\frac{\lambda}{\lambda_\mathrm{max}}\right)^2 \right] .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The scaled nondimensional potential energy(s).
- Return type:
numpy.ndarray
- class LennardJonesPotential(**kwargs)[source]
Bases:
object
The Lennard-Jones potential [1].
- varepsilon
The nondimensional energy scale.
- Type:
float
- kappa
The nondimensional stiffness \(\kappa=72\varepsilon\).
- Type:
float
- c
The correction parameter \(c=2/23\).
- Type:
float
- eta_max
The maximum nondimensional force \(\eta_\mathrm{max} = \eta(\lambda_\mathrm{max})\).
- Type:
float
- lambda_max
The stretch at the maximum nondimensional force, \(\lambda_\mathrm{max} = (13/7)^{1/6}\).
- Type:
float
- beta_u(lambda_)[source]
The nondimensional potential energy function,
\[\beta u(\lambda) = \varepsilon\phi(\lambda) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional potential energy(s).
- Return type:
numpy.ndarray
- class LogSquaredPotential(**kwargs)[source]
Bases:
object
The log-squared potential [3].
- varepsilon
The nondimensional energy scale.
- Type:
float
- kappa
The nondimensional stiffness \(\kappa=\varepsilon\).
- Type:
float
- c
The correction parameter \(c=2/5\).
- Type:
float
- eta_max
The maximum nondimensional force \(\eta_\mathrm{max} = e^{-1}\varepsilon\).
- Type:
float
- lambda_max
The stretch at the maximum nondimensional force, \(\lambda_\mathrm{max} = e^{1}\).
- Type:
float
- beta_u(lambda_)[source]
The nondimensional potential energy function,
\[\beta u(\lambda) = \varepsilon\phi(\lambda) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional potential energy(s).
- Return type:
numpy.ndarray
- delta_lambda(eta)[source]
The incremental stretch as a function of nondimensional force,
\[\Delta\lambda(\eta) = e^{-\mathcal{W}(-\eta/\varepsilon)} ,\qquad \eta\in[0,\eta_\mathrm{max}] .\]- Parameters:
eta (array_like) – The nondimensional force(s).
- Returns:
The incremental stretch(s).
- Return type:
numpy.ndarray
- eta_link(lambda_)[source]
The nondimensional force as a function of stretch,
\[\eta(\lambda) = \varepsilon\,\frac{\ln(\lambda)}{\lambda} .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional force(s).
- Return type:
numpy.ndarray
Example
Compute the nondimensional force at a sample stretch:
>>> from ufjc.potential import LogSquaredPotential >>> LogSquaredPotential().eta_link(1.8) 28.736236950770266
- class MiePotential(**kwargs)[source]
Bases:
object
The Mie potential [4].
- varepsilon
The nondimensional energy scale.
- Type:
float
- n
The repulsive exponent.
- Type:
float
- m
The attractive exponent.
- Type:
float
- kappa
The nondimensional stiffness \(\kappa=nm\varepsilon\).
- Type:
float
- c
The correction parameter \(c=\frac{4m(m+1)-2n(n+1)}{2m(m^2+5m+4)-n(n^2+5n+4)}\).
- Type:
float
- eta_max
The maximum nondimensional force \(\eta_\mathrm{max} = \eta(\lambda_\mathrm{max})\).
- Type:
float
- lambda_max
The stretch at the maximum nondimensional force, \(\lambda_\mathrm{max}=\left[\frac{n+1}{m+1}\right]^{1/(n-m)}\) .
- Type:
float
- beta_u(lambda_)[source]
The nondimensional potential energy function,
\[\beta u(\lambda) = \varepsilon\phi(\lambda) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional potential energy(s).
- Return type:
numpy.ndarray
- eta_link(lambda_)[source]
The nondimensional force as a function of stretch,
\[\eta(\lambda) = \varepsilon\, \frac{nm}{(n-m)} \left(\frac{1}{\lambda^{m+1}} - \frac{1}{\lambda^{n+1}}\right) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional force(s).
- Return type:
numpy.ndarray
- phi(lambda_)[source]
The scaled nondimensional potential energy function,
\[\phi(\lambda) = \frac{1}{(n-m)} \left(\frac{m}{\lambda^n} - \frac{n}{\lambda^m}\right) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The scaled nondimensional potential energy(s).
- Return type:
numpy.ndarray
- class MorsePotential(**kwargs)[source]
Bases:
object
The Morse potential [5].
- varepsilon
The nondimensional energy scale.
- Type:
float
- alpha
The Morse parameter.
- Type:
float
- kappa
The nondimensional stiffness \(\kappa=2\varepsilon\alpha^2\).
- Type:
float
- c
The correction parameter \(c=1/(1+3\alpha/2)\).
- Type:
float
- eta_max
The maximum nondimensional force \(\eta_\mathrm{max} = \sqrt{\kappa\varepsilon/8}\).
- Type:
float
- lambda_max
The stretch at the maximum nondimensional force, \(\lambda_\mathrm{max} = 1+\ln(2)/\alpha\).
- Type:
float
- beta_u(lambda_)[source]
The nondimensional potential energy function,
\[\beta u(\lambda) = \varepsilon\phi(\lambda) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional potential energy(s).
- Return type:
numpy.ndarray
- delta_lambda(eta)[source]
The incremental stretch as a function of nondimensional force,
\[\Delta\lambda(\eta) = \ln\left(\frac{2}{1 + \sqrt{1 - \eta/\eta_\mathrm{max}}}\right)^{1/\alpha} ,\qquad \eta\in[0,\eta_\mathrm{max}] .\]- Parameters:
eta (array_like) – The nondimensional force(s).
- Returns:
The incremental stretch(s).
- Return type:
numpy.ndarray
- eta_link(lambda_)[source]
The nondimensional force as a function of stretch,
\[\eta(\lambda) = 2\alpha\varepsilon e^{-\alpha(\lambda - 1)} \left[1 - e^{-\alpha(\lambda - 1)}\right]^2 .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional force(s).
- Return type:
numpy.ndarray
Example
Compute the nondimensional force at a sample stretch:
>>> from ufjc.potential import MorsePotential >>> _ = MorsePotential().eta_link(1.23)
- class PolynomialPotential(**kwargs)[source]
Bases:
object
A polynomial potential.
- varepsilon
The nondimensional energy scale.
- Type:
float
- coefficients
The coefficients \(a_k\).
- Type:
array_like
- kappa
The nondimensional stiffness \(\kappa=a_1\varepsilon\).
- Type:
float
- c
The correction parameter \(c=(1 - \frac{a_2}{2a_1})^{-1}\).
- Type:
float
- beta_u(lambda_)[source]
The nondimensional potential energy function,
\[\beta u(\lambda) = \varepsilon\phi(\lambda) .\]- Parameters:
lambda (array_like) – The stretch(s).
- Returns:
The nondimensional potential energy(s).
- Return type:
numpy.ndarray
- class Potential(**kwargs)[source]
Bases:
object
A class to assign a potential to a given model through inheritance.
- potential
The potential type.
- Type:
str
- pot
The potential model instance.
- Type:
object
References
John Edward Jones. On the determination of molecular fields. II. From the equation of state of a gas. Proceedings of the Royal Society of London. Series A, Containing Papers of a Mathematical and Physical Character, 106(738):463–477, 1924. doi:10.1098/rspa.1924.0082.
Kurt Kremer and Gary S Grest. Dynamics of entangled linear polymer melts: A molecular-dynamics simulation. The Journal of Chemical Physics, 92(8):5057–5086, 1990. doi:10.1063/1.458541.
Yunwei Mao, Brandon Talamini, and Lallit Anand. Rupture of polymers by chain scission. Extreme Mechanics Letters, 13:17–24, 2017. doi:10.1016/j.eml.2017.01.003.
Gustav Mie. Zur kinetischen theorie der einatomigen körper. Annalen der Physik, 316(8):657–697, 1903. doi:10.1002/andp.19033160802.
Philip M Morse. Diatomic molecules according to the wave mechanics. II. Vibrational levels. Physical Review, 34(1):57, 1929. doi:10.1103/PhysRev.34.57.