ufjc.core module

The core module for the uFJC single-chain model.

This module consist of the class uFJC which, upon instantiation, becomes a uFJC single-chain model instance with methods for computing single-chain quantities in either thermodynamic ensemble.

Example

Import and create an instance of the model:

>>> from ufjc import uFJC
>>> model = uFJC()
class uFJC(**kwargs)[source]

Bases: Potential, uFJCIsometric

The uFJC single-chain model class.

This class represents the uFJC single-chain model, meaning that an instance of this class is a uFJC single-chain model instance. It inherits all attributes and methods from the uFJCIsometric class, which inherits all attributes and methods from the uFJCIsotensional class, which inherits all attributes and methods from the BasicUtility class. It also inherits a potential from the Potential class as the attribute pot, a model instance itself. Keyword arguments are utilized during instantiation in order to specify model parameters; see the inherited classes and various examples through the documentation for more information.

N_b

The number of links in the chain.

Type:

int

init_config

The initial configuration for Monte Carlo.

Type:

np.ndarray

c_kappa

The constant used for ideal/Gaussian approximations.

Type:

float

P_eq_normalizations

The normalizations for the nondimensional equilibrium distribution P_eq for each approach.

Type:

dict

pot

The link potential model instance.

Type:

object

Note

The attributes of a instantiated model should not be changed. There are several quantities that are computed and stored during instantiation that depend on these attributes and are not re-computed.

P_eq(gamma, **kwargs)[source]

The nondimensional equilibrium distribution function.

This function computes the nondimensional probability density distribution of nondimensional end-to-end vectors at equilibrium,

\[\mathscr{P}_\mathrm{eq}(\boldsymbol{\gamma}) = \frac{e^{-N_b\,\vartheta(\boldsymbol{\gamma})}} {\iiint e^{-N_b\,\vartheta(\tilde{\boldsymbol{\gamma}})} \,d^3\tilde{\boldsymbol{\gamma}}} .\]

Since the uFJC model is spherically-symmetric in \(\boldsymbol{\gamma}\), P_eq is only a function of the scalar nondimensional end-to-end length \(\gamma\), and the normalization can be calculated using a one-dimensional integral over \(\gamma\) instead [3],

\[\iiint e^{-N_b\,\vartheta(\tilde{\boldsymbol{\gamma}})} \,d^3\tilde{\boldsymbol{\gamma}} = 4\pi\int e^{-N_b\,\vartheta(\tilde{\gamma})} \,\tilde{\gamma}^2 d\tilde{\gamma} .\]

The normalization for the distribution is computed for a given approach once the function is called using that approach.

Parameters:
  • gamma (array_like) – The nondimensional end-to-end length(s).

  • **kwargs – Arbitrary keyword arguments. Passed to vartheta.

Returns:

The nondimensional equilibrium probability density(s).

Return type:

numpy.ndarray

Example

Evaluate the probability density at an end-to-end length of 0.23 using two different approximation methods:

>>> from ufjc import uFJC
>>> model = uFJC(potential='morse', N_b=8)
>>> model.P_eq(23e-2)
array([4.19686303])
>>> model.P_eq(23e-2, gaussian=True)
array([3.86142625])

Example

Verify that the probability density is normalized:

>>> from ufjc import uFJC
>>> model = uFJC()
>>> from scipy.integrate import quad
>>> integrand = lambda gamma: \
...     4*np.pi*gamma**2*model.P_eq(gamma)
>>> P_tot_eq = quad(integrand, 0, np.inf)[0]
>>> np.isclose(P_tot_eq, 1)
True
eta(gamma, **kwargs)[source]

The nondimensional force.

This function computes the scalar nondimensional force as a function of the scalar nondimensional end-to-end length, \(\eta(\gamma)\). In the isometric ensemble, it is given by [7]

\[\gamma(\eta) = -\frac{1}{N_b} \frac{\partial}{\partial\gamma} \,\ln\mathfrak{Q}(\boldsymbol{\gamma}) ,\]

where the nondimensional isometric partition function is given by [2]

\[\mathfrak{Q}(\boldsymbol{\gamma}) = \iiint d^3\boldsymbol{\lambda}_1 \,e^{-\varepsilon\phi(\lambda_1)} \cdots \iiint d^3\boldsymbol{\lambda}_{N_b} \,e^{-\varepsilon\phi(\lambda_{N_b})} ~\delta^3\left(\sum_{j=1}^{N_b}\boldsymbol{\lambda}_j - \boldsymbol{\gamma}\right) .\]

This is analytically intractable for nearly all single-chain models, including the uFJC model. Fortunately, there are both exact numerical methods, such as Monte Carlo calculations, and accurate approximation methods, such as the Legendre transformation method which becomes valid relatively quickly as \(N_b\) increases.

Parameters:
  • gamma (array_like) – The nondimensional end-to-end length(s).

  • **kwargs – Arbitrary keyword arguments.

Returns:

The nondimensional force(s).

Return type:

numpy.ndarray

Example

Compute the nondimensional force at an end-to-end length in the isometric ensemble, which (by default) is done using the Legendre transformation method and the asymptotic approach to get \(\gamma(\eta)\) in the isotensional ensemble. Using the Legendre transformation method with a given isotensional approach for the isometric ensemble is equivalent to using the same approach in the isotensional ensemble:

>>> from ufjc import uFJC
>>> model = uFJC(potential='harmonic', varepsilon=23, N_b=8)
>>> model.eta(1.25508427, ensemble='isometric')
array([8.00000001])
>>> model.eta(1.25508427, ensemble='isotensional')
array([8.00000001])

Example

Verify that the ideal approximation is valid for small nondimensional end-to-end lengths:

>>> import numpy as np
>>> from ufjc import uFJC
>>> model = uFJC()
>>> eta_exact = model.eta(1e-2, approach='exact')
>>> eta_ideal = model.eta(1e-2, ideal=True)
>>> np.abs((eta_exact - eta_ideal)/eta_exact) < 1e-3
array([ True])
g_eq(gamma, **kwargs)[source]

The nondimensional equilibrium radial distribution function.

This function computes the nondimensional radial probability density distribution of nondimensional end-to-end lengths at equilibrium,

\[\mathscr{g}_\mathrm{eq}(\gamma) = 4\pi\gamma^2 \mathscr{P}_\mathrm{eq}(\gamma) .\]
Parameters:
  • gamma (array_like) – The nondimensional end-to-end length(s).

  • **kwargs – Arbitrary keyword arguments. Passed to P_eq.

Returns:

The nondimensional equilibrium radial probability density(s).

Return type:

numpy.ndarray

Example

Verify that the function is normalized by integrating over all permissible end-to-end lengths, which in effect ensures that the total probability is unity:

>>> import numpy as np
>>> from ufjc import uFJC
>>> model = uFJC(potential='morse')
>>> from scipy.integrate import quad
>>> P_tot_eq = quad(model.g_eq, 0, model.lambda_max)[0]
>>> np.isclose(P_tot_eq, 1)
True
gamma(eta, **kwargs)[source]

The nondimensional end-to-end length.

This function computes the scalar nondimensional end-to-end length as a function of the scalar nondimensional force, \(\gamma(\eta)\). In the isotensional ensemble, it is given by [6]

\[\gamma(\eta) = \frac{\partial}{\partial\eta}\,\ln\mathfrak{z}(\eta) ,\]

where the single-link (since \(\gamma(\eta)\) is independent of \(N_b\) in the isotensional ensemble for the uFJC model) nondimensional isotensional partition function is given by [5]

\[\mathfrak{z}(\eta) = 4\pi\int \frac{\sinh(\lambda\eta)}{\lambda\eta} \, e^{-\varepsilon\phi(\lambda)} \, \lambda^2 d\lambda .\]

This function is rarely exactly analytically available, is sometimes accurately approximated using analytically tractable asymptotic relations, and is otherwise numerically calculated using Monte Carlo or quadrature approaches. For the uFJC model, accurate asymptotic relations are available, and for the EFJC model (harmonic link potentials), it is actually possible to exactly find \(\gamma(\eta)\) closed-form [1].

Parameters:
  • eta (array_like) – The nondimensional force(s).

  • **kwargs – Arbitrary keyword arguments.

Returns:

The nondimensional end-to-end length(s).

Return type:

numpy.ndarray

Examples

Create an EFJC model with a nondimensional link energy \(\varepsilon=23\), and calculate the nondimensional end-to-end length under an nondimensional force of 8 using many of the available approaches:

>>> from ufjc import uFJC
>>> model = uFJC(potential='harmonic', varepsilon=23)
>>> model.gamma(8, approach='exact')
array([1.25508427])
>>> model.gamma(8, approach='asymptotic')
array([1.25508427])
>>> model.gamma(8, approach='reduced')
array([1.22282631])
>>> model.gamma(8, approach='quadrature')
array([1.25508427])

Example

Verify that the ideal approximation is valid for small nondimensional forces:

>>> import numpy as np
>>> from ufjc import uFJC
>>> model = uFJC()
>>> gamma_exact = model.gamma(1e-2, approach='exact')
>>> gamma_ideal = model.gamma(1e-2, ideal=True)
>>> np.isclose(gamma_exact, gamma_ideal)
array([ True])

Example

Verify that the single-chain mechanical response is the same in the isotensional ensemble as it is in the isometric ensemble when utilizing the Legendre transformation method:

>>> from ufjc import uFJC
>>> model = uFJC(potential='harmonic', varepsilon=23)
>>> eta = np.random.rand(88)
>>> (np.abs(model.gamma(eta, ensemble='isotensional') -
...         model.gamma(eta, ensemble='isometric',
...                     method='legendre')) < 1e-6).all()
True
k(gamma)[source]

The nondimensional net forward reaction rate coefficient function.

This function computes the net forward reaction rate coefficient as a function of the nondimensional end-to-end length, i.e. \(k(\gamma)=N_bk'(\gamma)\), where \(k'(\gamma)\) is the forward reaction rate coefficient function for a single link. This function is obtained using the axioms of transition state theory [8], and is currently implemented to use the Legendre transformation method and the reduced asymptotic approach [4].

Parameters:

gamma (array_like) – The nondimensional end-to-end length(s).

Returns:

The nondimensional net forward reaction rate coefficient(s).

Return type:

numpy.ndarray

Example

Compute the nondimensional net forward reaction rate coefficients at several nondimensional end-to-end lengths:

>>> from ufjc import uFJC
>>> model = uFJC(potential='morse')
>>> model.k([0.8, 1.01, 1.6])
array([6.23755528e-01, 9.80472713e-01, 6.59105919e+07])
vartheta(gamma, **kwargs)[source]

The nondimensional Helmholtz free energy per link.

This function compute the nondimensional Helmholtz free energy per link in the isometric ensemble, which is generally given by

\[\vartheta(\boldsymbol{\gamma}) = -\ln\left[\mathfrak{Q}(\boldsymbol{\gamma})\right]^{1/N_b} .\]
Parameters:
  • gamma (array_like) – The nondimensional end-to-end length(s).

  • **kwargs – Arbitrary keyword arguments. Passed to vartheta_isometric.

Returns:

The nondimensional Helmholtz free energy(s) per link.

Return type:

numpy.ndarray

Example

Create a uFJC single-chain model with log-squared link potentials and calculate the Helmholtz free energy per link at a nondimensional end-to-end length of a half:

>>> from ufjc import uFJC
>>> model = uFJC(potential='log-squared')
>>> model.vartheta(0.5, method='legendre', approach='reduced')
array([0.39097337])

Example

Verify that the ideal approximation is valid for small nondimensional end-to-end lengths:

>>> import numpy as np
>>> from ufjc import uFJC
>>> model = uFJC()
>>> vartheta = model.vartheta(1e-2, approach='reduced')
>>> vartheta_ideal = model.vartheta(1e-2, ideal=True)
>>> np.abs((vartheta - vartheta_ideal)/vartheta) < 1e-1
array([ True])

References

[1]

NK Balabaev and TN Khazanovich. Extension of chains composed of freely joined elastic segments. Russian Journal of Physical Chemistry B, 3(2):242–246, 2009. doi:10.1134/S1990793109020109.

[2]

Michael R. Buche. Fundamental Theories for the Mechanics of Polymer Chains and Networks. PhD thesis, Cornell University, 2021. URL: https://www.proquest.com/openview/620bd18d1bf93950b88a41aa62ebdb3c.

[3]

Michael R. Buche and Meredith N. Silberstein. Statistical mechanical constitutive theory of polymer networks: the inextricable links between distribution, behavior, and ensemble. Physical Review E, 102(012501):012501, 2020. doi:10.1103/PhysRevE.102.012501.

[4]

Michael R. Buche and Meredith N. Silberstein. Chain breaking in the statistical mechanical constitutive theory of polymer networks. Journal of the Mechanics and Physics of Solids, 156:104593, 2021. doi:10.1016/j.jmps.2021.104593.

[5]

Michael R. Buche, Meredith N. Silberstein, and Scott J. Grutzik. Freely jointed chain models with extensible links. Physical Review E, 106(024502):024502, 2022. doi:10.1103/PhysRevE.106.024502.

[6]

Alessandro Fiasconaro and Fernando Falo. Analytical results of the extensible freely jointed chain model. Physica A: Statistical Mechanics and its Applications, 532:121929, 2019. doi:10.1016/j.physa.2019.121929.

[7]

Fabio Manca, Stefano Giordano, Pier Luca Palla, Rinaldo Zucca, Fabrizio Cleri, and Luciano Colombo. Elasticity of flexible and semiflexible polymers with extensible bonds in the Gibbs and Helmholtz ensembles. The Journal of Chemical Physics, 136(15):154906, 2012. doi:10.1063/1.4704607.

[8]

Robert W. Zwanzig. Nonequilibrium Statistical Mechanics. Oxford university press, 2001. URL: https://global.oup.com/academic/product/nonequilibrium-statistical-mechanics-9780195140187.