ufjc.isotensional

A module for the uFJC single-chain model in the isotensional ensemble.

This module consist of the class uFJCIsotensional which contains methods for computing single-chain quantities in the isotensional (constant end-to-end force) thermodynamic ensemble.

Example

Import and instantiate the class:

>>> from ufjc.isotensional import uFJCIsotensional
>>> class_instance = uFJCIsotensional()
class uFJCIsotensional[source]

Bases: BasicUtility

The uFJC single-chain model class for the isotensional ensemble.

This class contains methods for computing single-chain quantities in the isotensional (constant end-to-end force) thermodynamic ensemble. It inherits all attributes and methods from the BasicUtility class.

beta_Pi_config(config, eta_vector)[source]

The nondimensional total potential energy of a configuration.

This function provides the nondimensional total potential energy \(\beta \Pi = \beta U - N_b\eta\gamma\) given the configuration of the chain, i.e. the vector position of each atom/hinge relative to the first one.

Parameters:

config (numpy.ndarray) – The configuration of the chain, a \((N_b+1)\)-by-3 numpy array.

Returns:

The nondimensional total potential energy \(\beta\Pi\).

Return type:

float

efjc_gamma_isotensional_exact(eta)[source]

The exact approach for the EFJC isotensional \(\gamma(\eta)\).

This function provides the exact, analytical approach for obtaining \(\gamma(\eta)\) in the isotensional ensemble for harmonic link potentials, i.e. for the EFJC model,

\[\gamma(\eta) = \mathcal{L}(\eta) + \frac{\eta}{\kappa}\left[1 + \frac{1 - \mathcal{L}(\eta)\coth(\eta)} {1 + (\eta/\kappa)\coth(\eta)}\right] + \frac{\partial}{\partial\eta}\,\ln\left[1+g(\eta)\right] ,\]

where the function \(g(\eta)\) is given by [2]

\[g(\eta) = \frac{ e^{\eta}\left(\frac{\eta}{\kappa} + 1\right) \,\mathrm{erf}\left(\frac{\eta+\kappa}{\sqrt{2\kappa}}\right) - e^{-\eta}\left(\frac{\eta}{\kappa} - 1\right) \,\mathrm{erf}\left(\frac{\eta-\kappa}{\sqrt{2\kappa}}\right)} {4\sinh(\eta)\left[1 + (\eta/\kappa)\coth(\eta)\right]} - \frac{1}{2} .\]
Parameters:

eta (array_like) – The nondimensional force(s).

Returns:

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

Return type:

numpy.ndarray

Example

Compute the exact nondimensional end-to-end length of the EFJC in the isotensional ensemble under a nondimensional force of 23:

>>> from ufjc import uFJC
>>> model = uFJC(potential='harmonic')
>>> model.gamma(23, approach='exact')
array([1.22689438])
eta_isotensional(gamma, **kwargs)[source]

Invert the isotensional \(\gamma(\eta)\) for the isotensional \(\eta(\gamma)\).

This function obtains the isotensional nondimensional single-chain mechanical response \(\eta(\gamma)\) by inverting the isotensional \(\gamma(\eta)\).

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

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

Returns:

The nondimensional force(s).

Return type:

numpy.ndarray

Example

Check that \(\eta[\gamma(\eta)] = \eta\,\):

>>> import numpy as np
>>> from ufjc import uFJC
>>> model = uFJC(potential='morse')
>>> def check_gamma(eta):
...     gamma_fun = lambda eta: model.gamma_isotensional(eta)
...     eta_fun = lambda gamma: model.eta_isotensional(gamma)
...     return np.isclose(eta_fun(gamma_fun(eta))[0], eta)
>>> check_gamma(np.random.rand())
True
gamma_isotensional(eta, **kwargs)[source]

Main function for the isotensional \(\gamma(\eta)\).

This is the main function utilized to compute the isotensional nondimensional single-chain mechanical response. Keyword arguments specify and are passed onto the approach.

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

  • **kwargs – Arbitrary keyword arguments. Passed to gamma_isotensional_MHMCMC when the Monte Carlo approach is chosen.

Returns:

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

Return type:

numpy.ndarray

Example

Compute the nondimensional end-to-end length of the EFJC in the isotensional ensemble under a nondimensional force of 23 using numerical quadrature:

>>> from ufjc import uFJC
>>> model = uFJC(potential='harmonic')
>>> model.gamma_isotensional(23, approach='quadrature')
array([1.22689438])

Note

An improperly-specified approach will result in nan:

>>> from ufjc import uFJC
>>> uFJC().gamma_isotensional(0.8, approach='blah blah')
array([nan])
gamma_isotensional_MHMCMC(eta, **kwargs)[source]

The Monte Carlo approach for the isotensional \(\gamma(\eta)\).

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

  • **kwargs – Arbitrary keyword arguments. Passed to MHMCMC.gamma_isotensional_MHMCMC.

Returns:

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

Return type:

numpy.ndarray

Example

Compute the nondimensional end-to-end length of the log-squared-FJC in the isotensional ensemble under a nondimensional force of 5 using the Monte Carlo approach and ensure that it is within ten percent of the numerical quadrature result:

>>> from ufjc import uFJC
>>> model = uFJC(potential='log-squared')
>>> gamma_MHMCMC = model.gamma(5, approach='monte carlo')
>>> gamma_quad = model.gamma(5, approach='quadrature')
>>> np.abs(gamma_quad - gamma_MHMCMC)/gamma_quad < 1e-1
array([ True])
gamma_isotensional_asymptotic(eta)[source]

The full asymptotic approach for the isotensional \(\gamma(\eta)\).

This function provides the asymptotic approach for obtaining \(\gamma(\eta)\) in the isotensional ensemble [2], a composite approximation for all \(\eta\) that is valid for \(\varepsilon\gg 1\),

\[\gamma(\eta) \sim \mathcal{L}(\eta) + \frac{\eta}{\kappa}\left[ \frac{1 - \mathcal{L}(\eta)\coth(\eta)} {c + (\eta/\kappa)\coth(\eta)}\right] + \Delta\lambda(\eta) .\]
Parameters:

eta (array_like) – The nondimensional force(s).

Returns:

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

Return type:

numpy.ndarray

Example

Compute the asymptotically-correct nondimensional end-to-end length of the Mie-FJC in the isotensional ensemble under a nondimensional force of 23:

>>> from ufjc import uFJC
>>> model = uFJC(potential='mie')
>>> model.gamma(23, approach='asymptotic')
array([0.96204056])
gamma_isotensional_asymptotic_reduced(eta)[source]

The reduced asymptotic approach for the isotensional \(\gamma(\eta)\).

This function provides the reduced asymptotic approach for obtaining \(\gamma(\eta)\) in the isotensional ensemble, a composite approximation for all \(\eta\) that is valid for \(\varepsilon\gg 1\),

\[\gamma(\eta) \sim \mathcal{L}(\eta) + \Delta\lambda(\eta).\]

The simplification involves neglecting terms that contribute negligibly apart when \(\eta=\mathrm{ord}(\varepsilon)\), a region that tends to vanish for \(\varepsilon\gg 1\). This simplification is then asymptotically valid in the same limit that made the original asymptotic approximation possible. Further, this asymptotic approximation is equivalent to asymptotically matching the low-force entropic mechanical response \(\mathcal{L}(\eta)\) to the high-force enthalpic mechanical response \(\lambda(\eta)\) [1].

Parameters:

eta (array_like) – The nondimensional force(s).

Returns:

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

Return type:

numpy.ndarray

Example

Verify that either asymptotic approximation obtains close-to-exact end-to-end lengths over many nondimensional forces when the nondimensional energy scale is high:

>>> import numpy as np
>>> from ufjc import uFJC
>>> model = uFJC(potential='morse', varepsilon=888)
>>> eta_trial = np.logspace(-2, 1, 88)
>>> gamma_quad = model.gamma(eta_trial, approach='quadrature')
>>> gamma_asymp = model.gamma(eta_trial, approach='asymptotic')
>>> (np.abs(gamma_quad-gamma_asymp)/gamma_quad < 1e-2).all()
True
>>> gamma_simpl = model.gamma(eta_trial, approach='asymptotic')
>>> (np.abs(gamma_quad-gamma_simpl)/gamma_quad < 1e-2).all()
True
gamma_isotensional_exact(eta)[source]

Exact approach for the isotensional \(\gamma(\eta)\).

This function provides the exact, analytical approach for obtaining \(\gamma(\eta)\) in the isotensional ensemble.

Parameters:

eta (array_like) – The nondimensional force(s).

Returns:

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

Return type:

numpy.ndarray

Warning

This is currently only available for harmonic link potentials, and will therefore return nan for non-harmonic potentials:

>>> from ufjc import uFJC
>>> uFJC().gamma([88e-2, 8], approach='exact')
array([0.29518908, 0.97632595])
>>> uFJC(potential='morse').gamma([88e-2, 8], approach='exact')
array([nan, nan])
gamma_isotensional_quadrature(eta)[source]

The numerical quadrature approach for the isotensional \(\gamma(\eta)\).

Parameters:

eta (array_like) – The nondimensional force(s).

Returns:

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

Return type:

numpy.ndarray

Example

Compute the nondimensional end-to-end length of the EFJC model at many nondimensional forces using the numerical quadrature approach, ensuring the exact result is always obtained:

>>> import numpy as np
>>> from ufjc import uFJC
>>> model = uFJC()
>>> eta_trial = np.logspace(-2, 2, 88)
>>> gamma_quad = model.gamma(eta_trial, approach='quadrature')
>>> gamma_exact = model.gamma(eta_trial, approach='exact')
>>> np.isclose(gamma_quad, gamma_exact).all()
True

References

[1]

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.

[2] (1,2)

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.