ensembler.potentials package¶
Submodules¶
ensembler.potentials.ND module¶
- Module: Potential
- This module shall be used to implement subclasses of Potentials that formulate a potential as an Function with N-Dimensions.
This module contains all available potentials.
-
class
ensembler.potentials.ND.
envelopedPotential
(V_is: List[ensembler.potentials._basicPotentials._potentialNDCls] = (<ensembler.potentials.ND.harmonicOscillatorPotential object>, <ensembler.potentials.ND.harmonicOscillatorPotential object>), s: float = 1.0, eoff: List[float] = None, T: float = 1, kb: float = 1)[source]¶ Bases:
ensembler.potentials._basicPotentials._potentialNDCls
This implementation of exponential Coupling for EDS is a more numeric robust and variable implementation, it allows N states. Therefore the computation of energies and the deviation is not symbolic.
Here N-states are coupled by the log-sum-exp resulting in a new reference state $V_R$,
$V_R = -1/{eta} * ln(sum_i^Ne^(-eta*s*(V_i-E^R_i)))$
This potential coupling is for example used in EDS.
-
property
Eoff
¶ The Energy offsets are used to bias the single states in the reference potential by a constant offset. Therefore each state of the enveloping potential has its own energy offset.
- Returns
Eoff
- Return type
t.List[Number]
-
property
Eoff_i
¶ The Energy offsets are used to bias the single states in the reference potential by a constant offset. Therefore each state of the enveloping potential has its own energy offset.
- Returns
Eoff
- Return type
t.List[Number]
-
Eoffis
= Matrix([[Eoff_i]])¶
-
T
= T¶
-
V_functional
: sympy.core.function.Function = -T*kb*log(Sum(exp(-(-Matrix([[Eoff_i]])[i, 0] + Matrix([[V_i]])[i, 0])*Matrix([[s_i]])[i, 0]/(T*kb)), (i, 0, N)))/s_i¶
-
property
V_is
¶ V_is are the state potential classes enveloped by the reference state.
- Returns
V_is
- Return type
t.List[_potential1DCls]
-
Vis
= Matrix([[V_i]])¶
-
__init__
(V_is: List[ensembler.potentials._basicPotentials._potentialNDCls] = (<ensembler.potentials.ND.harmonicOscillatorPotential object>, <ensembler.potentials.ND.harmonicOscillatorPotential object>), s: float = 1.0, eoff: List[float] = None, T: float = 1, kb: float = 1)[source]¶ - __init__
This function constructs a enveloped potential, enveloping all given states.
- Parameters
V_is (List[_potential1DCls], optional) – The states(potential classes) to be enveloped (default: [harmonicOscillatorPotential(), harmonicOscillatorPotential(x_shift=3)])
s (float, optional) – the smoothing parameter, lowering the barriers between the states
eoff (List[float], optional) – the energy offsets of the individual states in the reference potential. These can be used to allow a more uniform sampling. (default: seta ll to 0)
T (float, optional) – the temperature of the reference state (default: 1 = T)
kb (float, optional) – the boltzman constant (default: 1 = kb)
-
_calculate_dvdpos_singlePos_overwrite
(positions: Iterable[float]) → numpy.array[source]¶ - Parameters
positions
-
beta
= 1/(T*kb)¶
-
i
= i¶
-
kb
= kb¶
-
nStates
= N¶
-
name
: str = 'Enveloping Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
property
s
¶
-
property
s_i
¶
-
set_Eoff
(Eoff: Union[numbers.Number, Iterable[numbers.Number]])[source]¶ This function is setting the Energy offsets of the states enveloped by the reference state. :Parameters: Eoff (Union[Number, Iterable[Number]])
-
set_s
(s: Union[numbers.Number, Iterable[numbers.Number]])[source]¶ set_s is a function used to set an smoothing parameter.
- Parameters
s (Union[Number, Iterable[Number]])
-
sis
= Matrix([[s_i]])¶
-
property
-
class
ensembler.potentials.ND.
harmonicOscillatorPotential
(k: numpy.array = array([1.0, 1.0, 1.0]), r_shift: numpy.array = array([0.0, 0.0, 0.0]), Voff: numpy.array = array([0.0, 0.0, 0.0]), nDimensions: int = 3)[source]¶ Bases:
ensembler.potentials._basicPotentials._potentialNDCls
ND harmonic oscillator potential
-
V_dim
= Matrix([[V_off + 0.5*k*(r - r_shift)**2]])¶
-
V_functional
: sympy.core.function.Function = Sum(Matrix([[V_off + 0.5*k*(r - r_shift)**2]])[i, 0], (i, 0, nDimensions))¶
-
Voff
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[V_off]])¶
-
__init__
(k: numpy.array = array([1.0, 1.0, 1.0]), r_shift: numpy.array = array([0.0, 0.0, 0.0]), Voff: numpy.array = array([0.0, 0.0, 0.0]), nDimensions: int = 3)[source]¶ - __init__
Constructs an harmonic Oscillator with an on runtime defined dimensionality.
- Parameters
k (List[float], optional) – force constants, as many as nDim, defaults to [1.0, 1.0, 1.0]
x_shift (List[float], optional) – shift of the minimum in the x Axis, as many as nDim, defaults to [0.0, 0.0, 0.0]
y_shift (List[float], optional) – shift on the y Axis, as many as nDim, defaults to [0.0, 0.0, 0.0]
nDim – dimensionality of the harmoic oscillator object. default: 3
-
i
= i¶
-
k
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[k]])¶
-
nDimensions
: int = nDimensions¶
-
name
: str = 'harmonicOscilator'¶
-
position
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[r]])¶
-
r_shift
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[r_shift]])¶
-
-
class
ensembler.potentials.ND.
lambdaEDSPotential
(V_is: List[ensembler.potentials._basicPotentials._potentialNDCls] = (<ensembler.potentials.ND.harmonicOscillatorPotential object>, <ensembler.potentials.ND.harmonicOscillatorPotential object>), lam: numbers.Number = 0.5, s: float = 1.0, eoff: List[float] = None, T: float = 1, kb: float = 1)[source]¶ Bases:
ensembler.potentials.ND.envelopedPotential
This implementation of exponential Coupling combined with linear compling is called $lambda$-EDS the implementation of function is more numerical robust to the hybrid coupling class.
Here two-states are coupled by the log-sum-exp and weighted by lambda resulting in a new reference state $V_R$,
$V_R = -
rac{1}{eta*s} * ln(lambda * e^(-eta*s*(V_A-E^R_A)) + (1-lambda)*e^(-eta*s*(V_B-E^R_B)))$
This potential coupling is for example used in $lambda$-EDS.
-
Eoffis
= Matrix([[Eoff_i]])¶
-
T
= T¶
-
V_functional
: sympy.core.function.Function = -T*kb*log(Sum(exp(-(-Matrix([[Eoff_i]])[i, 0] + Matrix([[V_i]])[i, 0])*Matrix([[s_i]])[i, 0]/(T*kb))*Matrix([[λ]])[i, 0], (i, 0, N)))/s_i¶
-
Vis
= Matrix([[V_i]])¶
-
beta
= 1/(T*kb)¶
-
i
= i¶
-
kb
= kb¶
-
property
lam
¶
-
property
lam_i
¶
-
lamis
= Matrix([[λ]])¶
-
nStates
= N¶
-
name
: str = 'lambda enveloped Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
sis
= Matrix([[s_i]])¶
-
-
class
ensembler.potentials.ND.
sumPotentials
(potentials: List[ensembler.potentials._basicPotentials._potentialNDCls] = (<ensembler.potentials.ND.harmonicOscillatorPotential object>, <ensembler.potentials.ND.harmonicOscillatorPotential object>))[source]¶ Bases:
ensembler.potentials._basicPotentials._potentialNDCls
- Adds n different potentials.
For adding up wavepotentials, we recommend using the addedwavePotential class.
-
V_functional
: sympy.core.function.Function = Sum(Matrix([[V_x]])[i, 0], (i, 0, N))¶
-
__init__
(potentials: List[ensembler.potentials._basicPotentials._potentialNDCls] = (<ensembler.potentials.ND.harmonicOscillatorPotential object>, <ensembler.potentials.ND.harmonicOscillatorPotential object>))[source]¶ This is the Constructor of an summed Potentials
- Parameters
potentials (List[_potential2DCls], optional) – it uses the 2D potential class to generate its potential, default to (wavePotential(), wavePotential(multiplicity=[3, 3]))
-
_initialize_functions
()[source]¶ converts the symbolic mathematics of sympy to a matrix representation that is compatible with multi-dimentionality.
-
i
= i¶
-
nPotentials
= N¶
-
name
: str = 'Summed Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
potentials
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[V_x]])¶
ensembler.potentials.OneD module¶
Module: Potential This module shall be used to implement subclasses of Potential. This module contains all available potentials.
-
class
ensembler.potentials.OneD.
_metadynamicsPotentialSympy
(origPotential, amplitude=0.1, sigma=0.1, n_trigger=100)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
The metadynamics bias potential adds Gaussian potentials on top of the original potential. The added gaussian potential is centered on the current position. Thereby the valleys of the potential “flooded” and barrier crossing is easier
This implementation uses sympy instead of a grid and is therefore super slow
-
__init__
(origPotential, amplitude=0.1, sigma=0.1, n_trigger=100)[source]¶ This is the Constructor of the metadynamicsPotential class. :Parameters: * origPotential (potential type) – The unbiased potential
amplitude (float) – scaling of the gaussian potential added in the metadynamcis step
sigma (float) – standard deviation of the gaussian potential added in the metadynamcis step
n_trigger (int) – Metadynamics potential will be added after every n_trigger’th steps
-
_update_potential
(curr_position)[source]¶ Is triggered by check_for_metastep(). Adds a gaussian centered on the current position to the potential
- Parameters
curr_position (float) – current x position
-
check_for_metastep
(curr_position)[source]¶ Checks if the bias potential should be added at the current step :Parameters: curr_position (flaot) – current x position
-
name
: str = 'Metadynamics Enhanced Sampling System using sympy'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
_timedependendBias
(origPotential, addPotential, n_trigger)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
The timedependend bias potential adds a user defined potential on top of the original potential.
This implementation uses sympy instead of a grid and is therefore super slow
-
__init__
(origPotential, addPotential, n_trigger)[source]¶ This is the Constructor of the addedPotential class. :Parameters: * origPotential (potential type) – The unbiased potential
addPotential (potential type) – The potential added on top of the unbiased potential to bias the system, usually of gaussian type
n_trigger (int) – Added potential will be added after every n_trigger’th steps
-
_update_potential
()[source]¶ Is triggered by check_for_metastep(). Adds the pre-defined potential on the current position to the potential
-
check_for_metastep
(curr_position)[source]¶ Checks if the bias potential should be added at the current step :Parameters: curr_position (flaot) – current x position
-
name
: str = 'Metadynamics Enhanced Sampling System'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
addedPotentials
(origPotential=<ensembler.potentials.OneD.harmonicOscillatorPotential object>, addPotential=<ensembler.potentials.OneD.gaussPotential object>)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Adds two different potentials on top of each other. Can be used to generate harmonic potential umbrella sampling or scaled potentials
-
__init__
(origPotential=<ensembler.potentials.OneD.harmonicOscillatorPotential object>, addPotential=<ensembler.potentials.OneD.gaussPotential object>)[source]¶ This is the Constructor of the addedPotential class. :Parameters: * origPotential (potential type) – The unbiased potential
addPotential (potential type) – The potential added on top of the unbiased potential to bias the system
-
bias_potential
= True¶
-
name
: str = 'Added Potential Enhanced Sampling System'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
coulombPotential
(q1=1, q2=1, epsilon=1)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Coulomb potential representing the pairwise electrostatic interaction of two charged particles
-
V_functional
: sympy.core.function.Function = q1*q2/(4*pi*e*r)¶
-
__init__
(q1=1, q2=1, epsilon=1)[source]¶ This is the Constructor of the Coulomb potential :Parameters: * q1 (int, optional) – Charge of atom 1, defaults to 1
q2 (int, optional) – Charge of atom 2, defaults to 1
epsilon (int, optional) – Electric Permetivitty, defaults to 1
-
charge1
= q1¶
-
charge2
= q2¶
-
electric_permetivity
= e¶
-
name
: str = 'Coulomb Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
doubleWellPotential
(Vmax=5, a=- 1, b=1)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
This is an implementation of a double Well potential
-
V_functional
: sympy.core.function.Function = V_max*(-b**2 + (-a/2 + r)**2)**2/b**4¶
-
Vmax
= V_max¶
-
__init__
(Vmax=5, a=- 1, b=1)[source]¶ This is the Constructor of the double well Potential
- Parameters
Vmax (int, optional) – Maximal barrier between minima, defaults to 5
a (int, optional) – defines x position of the minimum of the first well, defaults to -1
b (int, optional) – defines x position of the minimum of the second well, defaults to 1
-
a
= a¶
-
b
= b¶
-
name
: str = 'Double Well'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
dummyPotential
(y_shift: float = 0)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
This Dummy potential returns a simple constant value for each position
-
__init__
(y_shift: float = 0)[source]¶ This Class is representing the a dummy potential. It returns a constant value equalling the y_shift parameter.
- Parameters
y_shift (float, optional) – This will be the constant return value, defaults to 0
-
name
: str = 'Dummy Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
y_shift
= Voffset¶
-
-
class
ensembler.potentials.OneD.
envelopedPotential
(V_is: List[ensembler.potentials._basicPotentials._potential1DCls] = (<ensembler.potentials.OneD.harmonicOscillatorPotential object>, <ensembler.potentials.OneD.harmonicOscillatorPotential object>), s: float = 1.0, eoff: List[float] = None, T: float = 1, kb: float = 1)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
This implementation of exponential Coupling for EDS is a more numeric robust and variable implementation, it allows N states. Therefore the computation of energies and the deviation is not symbolic.
Here N-states are coupled by the log-sum-exp resulting in a new reference state $V_R$,
$V_R = -1/{eta} * ln(sum_i^Ne^(-eta*s*(V_i-E^R_i)))$
This potential coupling is for example used in EDS.
-
property
Eoff
¶ The Energy offsets are used to bias the single states in the reference potential by a constant offset. Therefore each state of the enveloping potential has its own energy offset.
- Returns
Eoff
- Return type
t.List[Number]
-
property
Eoff_i
¶ The Energy offsets are used to bias the single states in the reference potential by a constant offset. Therefore each state of the enveloping potential has its own energy offset.
- Returns
Eoff
- Return type
t.List[Number]
-
Eoffis
= Matrix([[Eoff_i]])¶
-
T
= T¶
-
V_functional
: sympy.core.function.Function = -T*kb*log(Sum(exp(-(-Matrix([[Eoff_i]])[i, 0] + Matrix([[V_i]])[i, 0])*Matrix([[s_i]])[i, 0]/(T*kb)), (i, 0, N)))/s_i¶
-
property
V_is
¶ V_is are the state potential classes enveloped by the reference state.
- Returns
V_is
- Return type
t.List[_potential1DCls]
-
Vis
= Matrix([[V_i]])¶
-
__init__
(V_is: List[ensembler.potentials._basicPotentials._potential1DCls] = (<ensembler.potentials.OneD.harmonicOscillatorPotential object>, <ensembler.potentials.OneD.harmonicOscillatorPotential object>), s: float = 1.0, eoff: List[float] = None, T: float = 1, kb: float = 1)[source]¶ - __init__
This function constructs a enveloped potential, enveloping all given states.
- Parameters
V_is (List[_potential1DCls], optional) – The states(potential classes) to be enveloped (default: [harmonicOscillatorPotential(), harmonicOscillatorPotential(x_shift=3)])
s (float, optional) – the smoothing parameter, lowering the barriers between the states
eoff (List[float], optional) – the energy offsets of the individual states in the reference potential. These can be used to allow a more uniform sampling. (default: seta ll to 0)
T (float, optional) – the temperature of the reference state (default: 1 = T)
kb (float, optional) – the boltzman constant (default: 1 = kb)
-
_calculate_dvdpos_singlePos_overwrite
(positions: Iterable[float]) → numpy.array[source]¶ - Parameters
positions
-
beta
= 1/(T*kb)¶
-
i
= i¶
-
kb
= kb¶
-
nStates
= N¶
-
name
: str = 'Enveloping Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
property
s
¶
-
property
s_i
¶
-
set_Eoff
(Eoff: Union[numbers.Number, Iterable[numbers.Number]])[source]¶ This function is setting the Energy offsets of the states enveloped by the reference state. :Parameters: Eoff (Union[Number, Iterable[Number]])
-
set_s
(s: Union[numbers.Number, Iterable[numbers.Number]])[source]¶ set_s is a function used to set an smoothing parameter.
- Parameters
s (Union[Number, Iterable[Number]])
-
sis
= Matrix([[s_i]])¶
-
property
-
class
ensembler.potentials.OneD.
exponentialCoupledPotentials
(Va: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, Vb: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, eoffA: float = 0, eoffB: float = 0, s: float = 1.0, temp: float = 298)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
This implementation of exponential Coupling is the symbolic variant of the more robust eds potential implementation. Here N-states are coupled by the log-sum-exp resulting in a new reference state $V_R$,
$V_R = -1/{eta} * ln(sum_i^Ne^(-eta*s*(V_i-E^R_i)))$
This potential coupling is for example used in EDS.
-
Va
= V_a¶
-
Vb
= V_b¶
-
__init__
(Va: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, Vb: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, eoffA: float = 0, eoffB: float = 0, s: float = 1.0, temp: float = 298)[source]¶ - __init__
This constructor is building a exponential coupled Potential out of two given end-states.
- Parameters
Va (_potential1DCls, optional) – potential function of state A (default: harmonic oscillator)
Vb (_potential1DCls, optional) – potential function of state B (default: harmonic oscillator)
eoffA (float, optional) – Energy offset of state A in the reference potential (default: 0)
eoffB (float, optional) – Energy offset of state B in the reference potential (default: 0)
s (float, optional) – smoothing factor of the reference potential (default: 1.0)
temp (float, optional) – Temperature of the reference state. (default: 298)
-
beta
= 0.008314462618*T¶
-
coupling
= -120.272355044943*log(exp(-0.008314462618*T*V_a*s - eoffJ) + exp(-0.008314462618*T*V_b*s - eoffI))/(T*s)¶
-
eoffA
= eoffI¶
-
eoffB
= eoffJ¶
-
name
: str = 'exponential Coupled System'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
s
= s¶
-
set_Eoff
(eoffA: float = 0, eoffB: float = 0)[source]¶ - set_Eoff
set the energy offsets for the states in the reference state.
- Parameters
eoffA (float, optional) – set a new offset for state A (default: None)
eoffB (float, optional) – set a new E offset for state B in the reference state (default: None)
-
set_s
(s: float)[source]¶ - set_s
sets a new s-value. (please only use this function to change s)
- Parameters
s (float) – the new sval.
-
temp
= T¶
-
-
class
ensembler.potentials.OneD.
flatwellPotential
(x_range: list = 0, 1, y_max: float = 4, y_min: float = 0)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
A flatwell potential returns a simple constant value (y_max) for each position except positions in the range x_min-x_max. In the defined phase space range the potential always returns a second defined value (y_min)
-
__init__
(x_range: list = 0, 1, y_max: float = 4, y_min: float = 0)[source]¶ __init__ This potential is a flatwell potential. The flatwell potential is a function similar to an if case. If a position is inside a the x_range, it returns the y_min val. If a position is outside, the y_max val will be returned.
- Parameters
x_range ((list, range), optional)
range inside this the y_min val will be returned, defaults to (0, 1)
y_max (float, optional)
outside of the range this value will be returned, defaults to 1000
y_min (float, optional)
inside the range this value will be returned, defaults to 0
-
name
: str = 'Flat Well'¶
-
x_max
: float = None¶
-
x_min
: float = None¶
-
y_max
: float = None¶
-
y_min
: float = None¶
-
-
class
ensembler.potentials.OneD.
forceField
[source]¶ Bases:
object
Force field potential energy that combines Coulomb, Lennard Jones and Torsion potentials
-
class
ensembler.potentials.OneD.
fourWellPotential
(Vmax=4, a=1.5, b=4.0, c=7.0, d=9.0, ah=2.0, bh=0.0, ch=0.5, dh=1.0)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Unperturbed four well potential
-
V_functional
: sympy.core.function.Function = -V_max*log(exp(-ah - (-a + r)**2) + exp(-bh - (-b + r)**2) + exp(-ch - (-c + r)**2) + exp(-dh - (-d + r)**2))¶
-
Vmax
= V_max¶
-
__init__
(Vmax=4, a=1.5, b=4.0, c=7.0, d=9.0, ah=2.0, bh=0.0, ch=0.5, dh=1.0)[source]¶ This is the Constructor of the four well Potential
- Parameters
Vmax (float, optional) – scaling of the whole potential
a (float, optional) – x position of the minimum of the first well
b (float, optional) – x position of the minimum of the second well
c (float, optional) – x position of the minimum of the third well
d (float, optional) – x position of the minimum of the fourth well
ah (str, optional) – ah*Vmax = y position of the first well
bh (str, optional) – bh*Vmax = y position of the second well
ch (str, optional) – ch*Vmax = y position of the third well
dh (str, optional) – dh*Vmax = y position of the fourth well
-
a
= a¶
-
ah
= ah¶
-
b
= b¶
-
bh
= bh¶
-
c
= c¶
-
ch
= ch¶
-
d
= d¶
-
dh
= dh¶
-
name
: str = 'Four Well Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
gaussPotential
(A=1.0, mu=0.0, sigma=1.0)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Gaussian like potential, usually used for metadynamics
-
A
= A¶
-
V_functional
: sympy.core.function.Function = A*exp(-(-mu + r)**2/(2*sigma**2))¶
-
__init__
(A=1.0, mu=0.0, sigma=1.0)[source]¶ This is the Constructor of a 1D Gauss Potential
- Parameters
A (float, optional) – scaling of the gauss function, defaults to 1.
mu (float, optional) – mean of the gauss function, defautls to 0.
sigma (float, optional) –
standard deviation of the gauss function, defaults to 1.
TODO: improve numerical stablility
-
_update_functions
()[source]¶ This function is needed to simplyfiy the symbolic equation on the fly and to calculate the position derivateive.
-
mu
= mu¶
-
name
: str = 'Gaussian Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
sigma
= sigma¶
-
-
class
ensembler.potentials.OneD.
harmonicOscillatorPotential
(k: float = 1.0, x_shift: float = 0.0, y_shift: float = 0.0)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Implementation of an 1D harmonic oscillator potential following hooke’s law
-
V_functional
: sympy.core.function.Function = Voffset + 0.5*k*(r - r_0)**2¶
-
__init__
(k: float = 1.0, x_shift: float = 0.0, y_shift: float = 0.0)[source]¶ This is the Constructor of the 1D harmonic oscillator
- Parameters
k (float, optional) – force constant, defaults to 1.0
x_shift (float, optional) – shift of the minimum in the x Axis, defaults to 0.0
y_shift (float, optional) – shift on the y Axis, defaults to 0.0
-
k
= k¶
-
name
: str = 'Harmonic Oscillator'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
x_shift
= r_0¶
-
y_shift
= Voffset¶
-
-
class
ensembler.potentials.OneD.
hybridCoupledPotentials
(Va: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, Vb: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, lam: float = 0.5, s: float = 1.0, temp: float = 298)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DClsPerturbed
This implementation of exponential Coupling combined with linear compling is called $lambda$-EDS this function is the purely symbolic class and therefore note very numerical stable (see lambda EDS).
Here two-states are coupled by the log-sum-exp and weighted by lambda resulting in a new reference state $V_R$,
$V_R = -
rac{1}{eta*s} * ln(lambda * e^(-eta*s*(V_A-E^R_A)) + (1-lambda)*e^(-eta*s*(V_B-E^R_B)))$
This potential coupling is for example used in $lambda$-EDS.
-
T
= T¶
-
Va
= V_a¶
-
Vb
= V_b¶
-
__init__
(Va: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, Vb: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, lam: float = 0.5, s: float = 1.0, temp: float = 298)[source]¶ - __init__
This function constructs a $lambda$-enveloped potential, enveloping all given states and weighting them by $lambda$.
- Parameters
V_is (List[_potential1DCls], optional) – The states(potential classes) to be enveloped (default: [harmonicOscillatorPotential(), harmonicOscillatorPotential(x_shift=3)])
s (float, optional) – the smoothing parameter, lowering the barriers between the states
Eoff_i (List[float], optional) – the energy offsets of the individual states in the reference potential. These can be used to allow a more uniform sampling. (default: seta ll to 0)
T (float, optional) – the temperature of the reference state (default: 1 = T)
kb (float, optional) – the boltzman constant (default: 1 = kb)
-
beta
= 1¶
-
coupling
: sympy.core.function.Function = -log(λ*exp(-V_b*s) + (1 - λ)*exp(-V_a*s))/s¶
-
lam
= λ¶
-
name
: str = 'hybrid Coupled Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
s
= s¶
-
-
class
ensembler.potentials.OneD.
lambdaEDSPotential
(V_is: List[ensembler.potentials._basicPotentials._potential1DCls] = (<ensembler.potentials.OneD.harmonicOscillatorPotential object>, <ensembler.potentials.OneD.harmonicOscillatorPotential object>), lam: numbers.Number = 0.5, s: float = 1.0, eoff: List[float] = None, T: float = 1, kb: float = 1)[source]¶ Bases:
ensembler.potentials.OneD.envelopedPotential
This implementation of exponential Coupling combined with linear compling is called $lambda$-EDS the implementation of function is more numerical robust to the hybrid coupling class.
Here two-states are coupled by the log-sum-exp and weighted by lambda resulting in a new reference state $V_R$,
$V_R = -
rac{1}{eta*s} * ln(lambda * e^(-eta*s*(V_A-E^R_A)) + (1-lambda)*e^(-eta*s*(V_B-E^R_B)))$
This potential coupling is for example used in $lambda$-EDS.
-
Eoffis
= Matrix([[Eoff_i]])¶
-
T
= T¶
-
V_functional
: sympy.core.function.Function = -T*kb*log(Sum(exp(-(-Matrix([[Eoff_i]])[i, 0] + Matrix([[V_i]])[i, 0])*Matrix([[s_i]])[i, 0]/(T*kb))*Matrix([[λ]])[i, 0], (i, 0, N)))/s_i¶
-
Vis
= Matrix([[V_i]])¶
-
beta
= 1/(T*kb)¶
-
i
= i¶
-
kb
= kb¶
-
property
lam
¶
-
property
lam_i
¶
-
lamis
= Matrix([[λ]])¶
-
nStates
= N¶
-
name
: str = 'lambda enveloped Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
sis
= Matrix([[s_i]])¶
-
-
class
ensembler.potentials.OneD.
lennardJonesForceFieldPotential
(c6: float = 0.2, c12: float = 0.0001, x_shift: float = 0, y_shift: float = 0)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
This is a forcefield like implementation of a lennard Jones Potential
-
V_functional
: sympy.core.function.Function = V_off + c12/(r - r_0)**12 - c6/(r - r_0**6)¶
-
__init__
(c6: float = 0.2, c12: float = 0.0001, x_shift: float = 0, y_shift: float = 0)[source]¶ This is the Constructor of the Lennard-Jones Field Potential :Parameters: * c6 (float, optional) – prefactor of the interaction term that scales with **6, defaults to 0.2
c12 (float, optional) – prefactor of the interaction term that scales with **12, defaults to 0.0001
x_shift (float, optional) – shift of potential on x Axis, defaults to 0
y_shift (float, optional) – shift of potential on y Axis, defaults to 0
-
c12
= c12¶
-
c6
= c6¶
-
name
: str = 'Lennard Jones Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
x_shift
= r_0¶
-
y_shift
= V_off¶
-
-
class
ensembler.potentials.OneD.
lennardJonesPotential
(sigma: float = 1.5, epsilon: float = 2, x_shift: float = 0, y_shift=0)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Lennard Jones potential representing the pairwise van-der-Waals interaction of two particles
-
V_functional
: sympy.core.function.Function = V_off + 4*e*(s**12/(r - r_0)**12 - s**6/(r - r_0)**6)¶
-
__init__
(sigma: float = 1.5, epsilon: float = 2, x_shift: float = 0, y_shift=0)[source]¶ This is the Constructor of the Lennard-Jones Potential
- Parameters
sigma (float, optional) – x - Position of the minimum, defaults to 1.5
epsilon (float, optional) – y - position of minimum, defaults to 2
x_shift (float, optional) – shift of potential on x Axis, defaults to 0
y_shift (int, optional) – shift of potential on y Axis, defaults to 0
-
epsilon
= e¶
-
name
: str = 'Lennard Jones Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
sigma
= s¶
-
x_shift
= r_0¶
-
y_shift
= V_off¶
-
-
class
ensembler.potentials.OneD.
linearCoupledPotentials
(Va: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, Vb: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, lam: float = 0.5)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DClsPerturbed
Linear Coupled Potential combines two potential as linear combinations, $ V_{lambda} = lambda * V_a + (1-lambda)*V_b $
This variant of coupling states is used for example in FEP, TI or BAR approaches.
-
Va
= V_a¶
-
Vb
= V_b¶
-
__init__
(Va: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, Vb: ensembler.potentials._basicPotentials._potential1DCls = <ensembler.potentials.OneD.harmonicOscillatorPotential object>, lam: float = 0.5)[source]¶ - __init__
This constructor builds a linear combination of Va and Vb potentials, with lam as a cofactor. Linear Coupled Potentials, like in FEP or TI simulations.]
- Parameters
Va (_potential1DCls, optional) – Potential A that is mixed to the new potential.
Vb (_potential1DCls, optional) – Potential B that is mixed to the new potential.
lam (float) – lam is representing the lambda variable
-
coupling
: sympy.core.function.Function = V_a*(1 - λ) + V_b*λ¶
-
lam
= λ¶
-
name
: str = 'Linear Coupled System'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
metadynamicsPotential
(origPotential=<ensembler.potentials.OneD.harmonicOscillatorPotential object>, amplitude=0.1, sigma=1, n_trigger=100, bias_grid_min=0, bias_grid_max=10, numbins=100)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
The metadynamics bias potential adds 1D Gaussian potentials on top of the original 1D potential. The added gaussian potential is centered on the current position. Thereby the valleys of the potential “flooded” and barrier crossing is easier.
This implementation uses a grid to store the biasing. This is much faster than calculating an ever increasing potential with sympy
-
__init__
(origPotential=<ensembler.potentials.OneD.harmonicOscillatorPotential object>, amplitude=0.1, sigma=1, n_trigger=100, bias_grid_min=0, bias_grid_max=10, numbins=100)[source]¶ This is the Constructor of the metadynamicsPotential class. :Parameters: * origPotential (potential type) – The unbiased potential
amplitude (float) – scaling of the gaussian potential added in the metadynamcis step
sigma (float) – standard deviation of the gaussian potential added in the metadynamcis step
n_trigger (int) – Metadynamics potential will be added after every n_trigger’th steps
bias_grid_min (float) – min value of the bias grid
bias_grid_max (float) – max value of the bias grid
numbins (float) – size of the grid bias and forces are saved in
-
_find_nearest
(array, value)[source]¶ Function that finds position of the closest entry to a given value in an array
- Parameters
array (np.array) – 1D array containing the midpoints of the metadynamics grid
value (int or float) – search value
Returns
Index of the entry closest to the given value
——-
-
_update_potential
(curr_position)[source]¶ Is triggered by check_for_metastep(). Adds a gaussian centered on the current position to the potential
- Parameters
curr_position (float) – current x position
-
bias_potential
= True¶
-
check_for_metastep
(curr_position)[source]¶ Checks if the bias potential should be added at the current step :Parameters: curr_position (tuple) – current x,y position
-
ene
(positions)[source]¶ calculates energy of particle also takes bias into account :Parameters: positions (tuple) – position on 1D potential energy surface
- Returns
- Return type
current energy
-
force
(positions)[source]¶ calculates derivative with respect to position also takes bias into account
- Parameters
positions (tuple) – position on 1D potential energy surface
Returns
current derivative dh/dpos
——-
-
name
: str = 'Metadynamics Enhanced Sampling System using grid bias'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.OneD.
sumPotentials
(potentials: List[ensembler.potentials._basicPotentials._potential1DCls] = (<ensembler.potentials.OneD.harmonicOscillatorPotential object>, <ensembler.potentials.OneD.harmonicOscillatorPotential object>))[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
- Adds n different potentials.
For adding up wavepotentials, we recommend using the addedwavePotential class.
-
V_functional
: sympy.core.function.Function = Sum(Matrix([[V_x]])[i, 0], (i, 0, N))¶
-
__init__
(potentials: List[ensembler.potentials._basicPotentials._potential1DCls] = (<ensembler.potentials.OneD.harmonicOscillatorPotential object>, <ensembler.potentials.OneD.harmonicOscillatorPotential object>))[source]¶ This is the Constructor of an summed Potentials
- Parameters
potentials (List[_potential2DCls], optional) – it uses the 2D potential class to generate its potential, default to (wavePotential(), wavePotential(multiplicity=[3, 3]))
-
_initialize_functions
()[source]¶ converts the symbolic mathematics of sympy to a matrix representation that is compatible with multi-dimentionality.
-
i
= i¶
-
nPotentials
= N¶
-
name
: str = 'Summed Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
potentials
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[V_x]])¶
-
class
ensembler.potentials.OneD.
torsionPotential
(wavePotentials=[<ensembler.potentials.OneD.wavePotential object>, <ensembler.potentials.OneD.wavePotential object>], radians=False)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Torsion potential that represents the energy potential of a torsion angle
-
N
= N¶
-
V_functional
: sympy.core.function.Function = Sum([1][i, 0], (i, 0, N))¶
-
__init__
(wavePotentials=[<ensembler.potentials.OneD.wavePotential object>, <ensembler.potentials.OneD.wavePotential object>], radians=False)[source]¶ This is the Constructor of a Torsion Potential
- Parameters
wavePotentials (list of two potentialTypes, optionel) – Torsion potential use the wave potential class to generate its potential, default to [wavePotential(), wavePotential(multiplicity=3)]
radians (bool, optional) – set potential to radians or degrees, defaults to False
-
i
= i¶
-
name
: str = 'Torsion Potential'¶
-
phase
: float = 1.0¶
-
position
: sympy.core.symbol.Symbol = r¶
-
set_degrees
(degrees: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
degrees (bool, optional,) – if True, output will be given in degrees, otherwise in radians, default: True
-
set_radians
(radians: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
radians (bool, optional,) – if True, output will be given in radians, otherwise in degree, default: True
-
wavePotentials
= [1]¶
-
-
class
ensembler.potentials.OneD.
wavePotential
(amplitude: float = 1.0, multiplicity: float = 1.0, phase_shift: float = 0.0, y_shift: float = 0.0, radians: bool = False)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential1DCls
Simple 1D wave potential consisting of a cosine function with given multiplicity, that can be shifted and elongated
-
V_functional
: sympy.core.function.Function = A*cos(m*(r + w)) + Voff¶
-
__init__
(amplitude: float = 1.0, multiplicity: float = 1.0, phase_shift: float = 0.0, y_shift: float = 0.0, radians: bool = False)[source]¶ This is the Constructor of the 1D wave potential function :Parameters: * amplitude (float, optional) – absolute min and max of the potential, defaults to 1.0
multiplicity (float, optional) – amount of minima in one phase, defaults to 1.0
phase_shift (float, optional) – shift of the potential on the x Axis, defaults to 0.0
y_offset (float, optional) – shift on the y Axis, defaults to 0.0
radians (bool, optional) – in radians or degrees, defaults to False
-
amplitude
= A¶
-
multiplicity
= m¶
-
name
: str = 'Wave Potential'¶
-
phase_shift
= w¶
-
position
: sympy.core.symbol.Symbol = r¶
-
set_degrees
(degrees: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
degrees (bool, optional,) – if True, output will be given in degrees, otherwise in radians, default: True
-
set_radians
(radians: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
radians (bool, optional,) – if True, output will be given in radians, otherwise in degree, default: True
-
y_shift
= Voff¶
-
ensembler.potentials.TwoD module¶
- Module: Potential
This module shall be used to implement subclasses of Potential. This module contains all available potentials.
-
class
ensembler.potentials.TwoD.
addedPotentials
(origPotential=<ensembler.potentials.TwoD.harmonicOscillatorPotential object>, addPotential=<ensembler.potentials.TwoD.gaussPotential object>)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential2DCls
Adds two different potentials on top of each other. Can be used to generate harmonic potential umbrella sampling or scaled potentials
-
__init__
(origPotential=<ensembler.potentials.TwoD.harmonicOscillatorPotential object>, addPotential=<ensembler.potentials.TwoD.gaussPotential object>)[source]¶ This is the Constructor of the addedPotential class. :Parameters: * origPotential (2D potential type) – The unbiased potential
addPotential (2D potential type) – The potential added on top of the unbiased potential to bias the system
-
bias_potential
= True¶
-
name
: str = 'Added Potential Enhanced Sampling System for 2D'¶
-
position
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[r]])¶
-
-
class
ensembler.potentials.TwoD.
addedWavePotential
(wave_potentials: List[ensembler.potentials.TwoD.wavePotential] = (<ensembler.potentials.TwoD.wavePotential object>, <ensembler.potentials.TwoD.wavePotential object>), degrees: bool = True)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential2DCls
Adds two wave potentials
-
V_functional
: sympy.core.function.Function = Sum(Matrix([[V_x]])[i, 0], (i, 0, N))¶
-
__init__
(wave_potentials: List[ensembler.potentials.TwoD.wavePotential] = (<ensembler.potentials.TwoD.wavePotential object>, <ensembler.potentials.TwoD.wavePotential object>), degrees: bool = True)[source]¶ This is the Constructor of an added wave Potential
- Parameters
wavePotentials (list of two 2D potentialTypes, optional) – is uses the 2D wave potential class to generate its potential, default to (wavePotential(), wavePotential(multiplicity=[3, 3]))
radians (bool, optional) – set potential to radians or degrees, defaults to False
-
_initialize_functions
()[source]¶ converts the symbolic mathematics of sympy to a matrix representation that is compatible with multi-dimentionality.
-
i
= i¶
-
nWavePotentials
= N¶
-
name
: str = 'Torsion Potential'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
set_degrees
(degrees: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
degrees (bool, optional,) – if True, output will be given in degrees, otherwise in radians, default: True
-
set_radians
(radians: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
radians (bool, optional,) – if True, output will be given in radians, otherwise in degree, default: True
-
wave_potentials
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[V_x]])¶
-
-
class
ensembler.potentials.TwoD.
gaussPotential
(amplitude=1.0, mu=0.0, 0.0, sigma=1.0, 1.0, negative_sign: bool = False)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential2DCls
Gaussian like potential, usually used for metadynamics
-
V_dim
= Matrix([[A_gauss*exp(-0.5*(-mu + r)**2/sigma**2)]])¶
-
V_functional
: sympy.core.function.Function = Sum(Matrix([[A_gauss*exp(-0.5*(-mu + r)**2/sigma**2)]])[i, 0], (i, 0, nDimensions))¶
-
__init__
(amplitude=1.0, mu=0.0, 0.0, sigma=1.0, 1.0, negative_sign: bool = False)[source]¶ - __init__
This is the Constructor of a 2D Gauss Potential
- Parameters
A (float, optional) – scaling of the gauss function, defaults to 1.
mu (tupel, optional) – mean of the gauss function, defaults to (0., 0.)
sigma (tupel, optional) – standard deviation of the gauss function, defaults to (1., 1.)
negative_sign (bool, optional) – this option is switching the sign of the final potential energy landscape. ==> mu defines the minima location, not maxima location
-
_initialize_functions
()[source]¶ converts the symbolic mathematics of sympy to a matrix representation that is compatible with multi-dimentionality.
-
_update_functions
()[source]¶ This function is needed to simplyfiy the symbolic equation on the fly and to calculate the position derivateive.
-
amplitude
= A_gauss¶
-
i
= i¶
-
mean
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[mu]])¶
-
nDimensions
: sympy.core.symbol.Symbol = nDimensions¶
-
name
: str = 'Gaussian Potential 2D'¶
-
position
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[r]])¶
-
sigma
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[sigma]])¶
-
-
class
ensembler.potentials.TwoD.
harmonicOscillatorPotential
(k: numpy.array = array([1.0, 1.0]), r_shift: numpy.array = array([0.0, 0.0]), Voff: numpy.array = array([0.0, 0.0]))[source]¶ Bases:
ensembler.potentials._basicPotentials._potential2DCls
Implementation of an 2D harmonic oscillator potential following hooke’s law
-
V_dim
= Matrix([[V_off + 0.5*k*(r - r_shift)**2]])¶
-
V_functional
: sympy.core.function.Function = Sum(Matrix([[V_off + 0.5*k*(r - r_shift)**2]])[i, 0], (i, 0, nDimensions))¶
-
Voff
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[V_off]])¶
-
__init__
(k: numpy.array = array([1.0, 1.0]), r_shift: numpy.array = array([0.0, 0.0]), Voff: numpy.array = array([0.0, 0.0]))[source]¶ This is the Constructor of the 2D harmonic oscillator
- Parameters
k (array, optional) – force constants in x and y direction, defaults to [1.0, 1.0]
r_shift (array, optional) – shift of the minimum in the x and y direction, defaults to [0.0, 0.0]
Voff (array, optional) – offset of the minimum, defaults to [0.0, 0.0]
-
_initialize_functions
()[source]¶ converts the symbolic mathematics of sympy to a matrix representation that is compatible with multi-dimentionality.
-
i
= i¶
-
k
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[k]])¶
-
nDimensions
: int = nDimensions¶
-
name
: str = 'harmonicOscilator'¶
-
position
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[r]])¶
-
r_shift
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[r_shift]])¶
-
-
class
ensembler.potentials.TwoD.
metadynamicsPotential
(origPotential=<ensembler.potentials.TwoD.harmonicOscillatorPotential object>, amplitude=1.0, sigma=(1.0, 1.0), n_trigger=100, bias_grid_min=(0, 0), bias_grid_max=(10, 10), numbins=(100, 100))[source]¶ Bases:
ensembler.potentials._basicPotentials._potential2DCls
The metadynamics bias potential adds 2D Gaussian potentials on top of the original 2D potential. The added gaussian potential is centered on the current position. Thereby the valleys of the potential “flooded” and barrier crossing is easier.
This implementation uses a grid to store the biasing. This is much faster than calculating an ever increasing potential with sympy
-
__init__
(origPotential=<ensembler.potentials.TwoD.harmonicOscillatorPotential object>, amplitude=1.0, sigma=(1.0, 1.0), n_trigger=100, bias_grid_min=(0, 0), bias_grid_max=(10, 10), numbins=(100, 100))[source]¶ - Parameters
origPotential (potential 2D type) – The unbiased potential
amplitude (float) – scaling of the gaussian potential added in the metadynamcis step
sigma (tuple) – standard deviation of the gaussian potential in x and y added in the metadynamcis step
n_trigger (int) – Metadynamics potential will be added after every n_trigger’th steps
bias_grid_min (tuple) – min value in x and y direction for the grid
bias_grid_max (tuple) – max value in x and y direction for the grid
numbins (tuple) – size of the grid bias and forces are saved in
-
_find_nearest
(array, value)[source]¶ Function that finds position of the closest entry to a given value in an array
- Parameters
array (np.array) – 1D array containing the midpoints of the metadynamics grid
value (int or float) – search value
Returns
Index of the entry closest to the given value
——-
-
_update_potential
(curr_position)[source]¶ Is triggered by check_for_metastep(). Adds a gaussian centered on the current position to the potential
- Parameters
curr_position (tuple) – current x,y position
-
bias_potential
= True¶
-
check_for_metastep
(curr_position)[source]¶ Checks if the bias potential should be added at the current step :Parameters: curr_position (tuple) – current x,y position
-
ene
(positions)[source]¶ calculates energy of particle also takes bias into account :Parameters: positions (tuple) – position on 2D potential energy surface
- Returns
- Return type
current energy
-
force
(positions)[source]¶ calculates derivative with respect to position also takes bias into account
- Parameters
positions (tuple) – position on 2D potential energy surface
Returns
current derivative dh/dpos
——-
-
name
: str = 'Metadynamics Enhanced Sampling System using grid bias in 2D'¶
-
position
: sympy.core.symbol.Symbol = r¶
-
-
class
ensembler.potentials.TwoD.
wavePotential
(amplitude=1, 1, multiplicity=1, 1, phase_shift=0, 0, y_offset=0, 0, radians: bool = False)[source]¶ Bases:
ensembler.potentials._basicPotentials._potential2DCls
Simple 2D wave potential consisting of cosine functions with given multiplicity, that can be shifted and elongated
-
V_dim
= Matrix([[A*cos(m*(omega + r)) + y_off]])¶
-
V_functional
: sympy.core.function.Function = Sum(Matrix([[A*cos(m*(omega + r)) + y_off]])[i, 0], (i, 0, nDimensions))¶
-
__init__
(amplitude=1, 1, multiplicity=1, 1, phase_shift=0, 0, y_offset=0, 0, radians: bool = False)[source]¶ This is the Constructor of the 2D wave potential function :Parameters: * amplitude (tuple, optional) – absolute min and max of the potential for the cosines in x and y direction, defaults to (1, 1)
multiplicity (tuple, optional) – amount of minima in one phase for the cosines in x and y direction, defaults to (1, 1)
phase_shift (tuple, optional) – position shift of the potential for the cosines in x and y direction, defaults to (0, 0)
y_offset (tuple, optional) – potential shift for the cosines in x and y direction, defaults to (0, 0)
radians (bool, optional) – in radians or degrees, defaults to False
-
_initialize_functions
()[source]¶ converts the symbolic mathematics of sympy to a matrix representation that is compatible with multi-dimentionality.
-
amplitude
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[A]])¶
-
i
= i¶
-
multiplicity
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[m]])¶
-
nDimensions
: sympy.core.symbol.Symbol = nDimensions¶
-
name
: str = 'Wave Potential'¶
-
phase_shift
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[omega]])¶
-
position
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[r]])¶
-
set_degrees
(degrees: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
degrees (bool, optional,) – if True, output will be given in degrees, otherwise in radians, default: True
-
set_radians
(radians: bool = True)[source]¶ Sets output to either degrees or radians
- Parameters
radians (bool, optional,) – if True, output will be given in radians, otherwise in degree, default: True
-
yOffset
: sympy.matrices.dense.MutableDenseMatrix = Matrix([[y_off]])¶
-
Module contents¶
Definition of the potential energy functions, that can be explored by the simulations.