Log-squared-FJC model thermodynamics (isotensional/asymptotic/reduced)

class LOGSQUAREDFJC(number_of_links, link_length, hinge_mass, link_stiffness)

The log-squared link potential freely-jointed chain (log-squared-FJC) model thermodynamics in the isotensional ensemble approximated using an asymptotic approach.

The number of links in the chain.

The length of each link in the chain in units of nm.

hinge_mass

The mass of each hinge in the chain in units of kg/mol.

The stiffness of each link in the chain in units of J/(mol⋅nm^2).

legendre

The thermodynamic functions of the model in the isotensional ensemble approximated using a reduced asymptotic approach and a Legendre transformation.

end_to_end_length(force, temperature)

The expected end-to-end length as a function of the applied force and temperature,

\[\xi(f, T) = -\frac{\partial\varphi}{\partial f}.\]
Parameters:
  • force (numpy.ndarray) – The force \(f\).

  • temperature (float) – The temperature \(T\).

Returns:

The end-to-end length \(\xi\).

Return type:

numpy.ndarray

The expected end-to-end length per link as a function of the applied force and temperature.

Parameters:
  • force (numpy.ndarray) – The force \(f\).

  • temperature (float) – The temperature \(T\).

Returns:

The end-to-end length per link \(\xi/N_b=\ell_b\gamma\).

Return type:

numpy.ndarray

nondimensional_end_to_end_length(nondimensional_force)

The expected nondimensional end-to-end length as a function of the applied nondimensional force.

Parameters:
  • nondimensional_force (numpy.ndarray) – The nondimensional force \(\eta\equiv\beta f\ell_b\).

  • temperature (float) – The temperature \(T\).

Returns:

The nondimensional end-to-end length \(N_b\gamma=\xi/\ell_b\).

Return type:

numpy.ndarray

The expected nondimensional end-to-end length per link as a function of the applied nondimensional force, given by Buche et al.[1] as

\[\gamma(\eta) \sim \mathcal{L}(\eta) + \Delta\lambda(\eta) \quad \text{for } \varepsilon,\kappa\gg 1,\]

where \(\mathcal{L}(x)=\coth(x)-1/x\) is the Langevin function, and \(\Delta\lambda(\eta)\) is the incremental link stretch,

\[\Delta\lambda(\eta) = e^{-W_0(-\eta/\kappa)} - 1,\]

where \(W_0(\eta)\) is the Lambert \(W\) function.

Parameters:
  • nondimensional_force (numpy.ndarray) – The nondimensional force \(\eta\equiv\beta f\ell_b\).

  • temperature (float) – The temperature \(T\).

Returns:

The nondimensional end-to-end length per link \(\gamma\equiv \xi/N_b\ell_b\).

Return type:

numpy.ndarray

gibbs_free_energy(force, temperature)

The Gibbs free energy as a function of the applied force and temperature,

\[\varphi(f, T) = -kT\ln Z(f, T).\]
Parameters:
  • force (numpy.ndarray) – The force \(f\).

  • temperature (float) – The temperature \(T\).

Returns:

The Gibbs free energy \(\varphi\).

Return type:

numpy.ndarray

The Gibbs free energy per link as a function of the applied force and temperature.

Parameters:
  • force (numpy.ndarray) – The force \(f\).

  • temperature (float) – The temperature \(T\).

Returns:

The Gibbs free energy per link \(\varphi/N_b\).

Return type:

numpy.ndarray

relative_gibbs_free_energy(force, temperature)

The relative Gibbs free energy as a function of the applied force and temperature.

Parameters:
  • force (numpy.ndarray) – The force \(f\).

  • temperature (float) – The temperature \(T\).

Returns:

The relative Gibbs free energy \(\Delta\varphi\equiv\varphi(f,T)-\varphi(0,T)\).

Return type:

numpy.ndarray

The relative Gibbs free energy per link as a function of the applied force and temperature.

Parameters:
  • force (numpy.ndarray) – The force \(f\).

  • temperature (float) – The temperature \(T\).

Returns:

The relative Gibbs free energy per link \(\Delta\varphi/N_b\).

Return type:

numpy.ndarray

nondimensional_gibbs_free_energy(nondimensional_force, temperature)

The nondimensional Gibbs free energy as a function of the applied nondimensional force and temperature.

Parameters:
  • nondimensional_force (numpy.ndarray) – The nondimensional force \(\eta\equiv\beta f\ell_b\).

  • temperature (float) – The temperature \(T\).

Returns:

The nondimensional Gibbs free energy \(\beta\varphi=N_b\varrho\).

Return type:

numpy.ndarray

The nondimensional Gibbs free energy per link as a function of the applied nondimensional force and temperature.

Parameters:
  • nondimensional_force (numpy.ndarray) – The nondimensional force \(\eta\equiv\beta f\ell_b\).

  • temperature (float) – The temperature \(T\).

Returns:

The nondimensional Gibbs free energy per link \(\varrho\equiv\beta\varphi/N_b\).

Return type:

numpy.ndarray

nondimensional_relative_gibbs_free_energy(nondimensional_force)

The nondimensional relative Gibbs free energy as a function of the applied nondimensional force.

Parameters:
  • nondimensional_force (numpy.ndarray) – The nondimensional force \(\eta\equiv\beta f\ell_b\).

  • temperature (float) – The temperature \(T\).

Returns:

The nondimensional relative Gibbs free energy \(\beta\Delta\varphi=N_b\Delta\varrho\).

Return type:

numpy.ndarray

The nondimensional relative Gibbs free energy per link as a function of the applied nondimensional force, given by Buche et al.[1] as

\[\Delta\varrho(\eta) \sim \ln\left[\frac{\eta}{\sinh(\eta)}\right] + \beta u[\lambda(\eta)] - \eta\Delta\lambda(\eta) \quad \text{for } \varepsilon,\kappa\gg 1,\]

where the nondimensional link potential \(\beta u\) is given by

\[\beta u(\lambda) = \frac{\varepsilon}{2}\left[\ln(\lambda)\right]^2,\]

where \(\varepsilon\equiv\beta u_b=\kappa\) is the nondimensional potential energy scale, \(\kappa\equiv\beta k_b\ell_b^2\) is the nondimensional link stiffness, and \(\lambda\equiv\ell/\ell_b\) is the nondimensional link stretch.

Parameters:
  • nondimensional_force (numpy.ndarray) – The nondimensional force \(\eta\equiv\beta f\ell_b\).

  • temperature (float) – The temperature \(T\).

Returns:

The nondimensional relative Gibbs free energy per link \(\Delta\varrho\equiv\beta\Delta\varphi/N_b\).

Return type:

numpy.ndarray


References