Methods in the SAFT class
SurfPack - Latest (beta)
The SAFT
class, found in surfpack/saft.py
, is an abstract class, that is inherited
by the SAFT_VR_Mie
and PC_SAFT
classes. It contains some generic utility methods to
compute quantities of interest when investigating SAFT-type functionals.
Table of contents
- Utility methods
- Internal methods
- Profile property interfaces
- Helmholtz energy contributions
- Weighted densities
- Weight functions
- Deprecated methods
Utility methods
Methods for computing specific parameters and contributions to the residual Helmholtz energy for SAFT-type equations of state
Table of contents
get_caching_id(self)
See Functional for docs.
get_eps_div_k(self, ic)
Get the epsilon parameter
Args:
ic (int) :
Component index (zero-indexed)
Returns:
float :
Model epsilon-parameter, divided by Boltzmanns constant [K]
get_sigma(self, ic)
Get the sigma parameter
Args:
ic (int) :
Component index (zero-indexed)
Returns:
float :
Model sigma-parameter [m]
pair_potential(self, i, j, r)
Evaluate the pair potential between component i
and j
at distance r
i (int) :
Component index
j (int) :
Component index
r (float) :
Distance [m]
Returns:
float :
Interaction potential energy [J]
reduce_temperature(self, T, c=0)
Compute the reduced temperature (LJ units)
Args:
T (float) :
Temperature [K]
c (int) :
Component idx to use
Return:
float :
Temperature in LJ units
set_association_active(self, active)
Toggle association contribution on/off
Args:
active (bool) :
Whether association is active
set_chain_active(self, active)
Toggle chain contribution on/off
Args:
active (bool) :
Whether chain is active
set_dispersion_active(self, active)
Toggle dispersion contribution on/off
Args:
active (bool) :
Whether dispersion is active
set_eps_div_k(self, ic, eps_div_k)
Set the model epsilon-parameter
Args:
ic (int) :
Component index (zero-indexed)
eps_div_k (float) :
epsilon-parameter divided by Boltzmanns constant [K]
set_multipole_active(self, active)
Toggle multipole contribution on/off
Args:
active (bool) :
Whether multipole is active
set_pure_assoc_param(self, ic, eps, beta)
Set pure-conponent association parameters
Args:
ic (int) :
Component index (zero indexed)
eps (float) :
Association energy [J / mol]
beta (float) :
Associaiton volume [-]
set_pure_fluid_param(self, ic, m, sigma, eps_div_k, *assoc_param)
Set all pure component parameters
Args:
ic (int) :
Component index (zero indexed)
m (float) :
Segment number [-]
sigma (float) :
sigma-parameter [m]
eps_div_k (float) :
epsilon-parameter [K]
set_segment_number(self, ic, m)
Set the segment number
ic (int) :
Component index (zero indexed)
m (float) :
Segment number
:
return:
set_sigma(self, ic, sigma)
Set the model sigma-parameter
Args:
ic (int) :
Component index (zero-indexed)
sigma (float) :
sigma-parameter [m]
Internal methods
Internal use, of little interest to outside users.
Table of contents
__init__(self, comps, eos, hs_model=<class 'surfpack.saft_hardsphere.SAFT_WhiteBear'>, parameter_ref='default')
This class is inherited by SAFT_VR_Mie and PC_SAFT. The only thing they do is pass the correct eos to the eos argument of this initialiser.
Args:
comps (str) :
Thermopack component string
eos (thermo) :
Thermopack EoS class (not initialised object)
hs_model (SAFT_HardSphere) :
A class of the same type as SAFT_WhiteBear (inheriting from SAFT_HardSphere)
__repr__(self)
Called from inheriting classes.
Returns:
str :
id. string with parameters, active contributions and hard sphere model identifier.
refresh_hs_model(self)
Update hard-sphere model such that parameters are in sync.
Profile property interfaces
Compute properties using a given density profile, without iterating the density profile to ensure equilibrium.
Table of contents
reduced_helmholtz_energy_density(self, rho, T, dphidn=False, bulk=False, asarray=False, dphidT=False)
Compute the reduced helmholtz energy density [1 / Å^3] (see the Functional class for explanation, this is simply an overlaod that puts toghether all the contributions.)
Args:
rho (list[Profile] or list[float]) :
Density profiles [particles / Å^3], indexed as rho[
Takes list[Profile] if bulk is False
, and list[float] if bulk is True
.
T (float) :
Temperature [K]
dphidn (bool) :
Whether to compute differentials wrt. weighted densities
bulk (bool) :
passed on to get_weighted_density, because convolutions are not neccesary for bulk.
If True, rho should be list[float].
asarray (bool) :
Do not set to True (I dont think that works yet)
Intent:
return helmholtz density as a 1d numpy array rather than a Profile
dphidT (bool) :
Compute temperature derivative
Returns:
Profile :
Reduced Helmholtz energy density [1 / Å^3]
Optionally list[Profile] :
derivatives wrt. weighted densities.
Helmholtz energy contributions
Methods for computing various contributions to the Helmholtz energy that are present in all SAFT-based functionals.
Table of contents
dispersion_helmholtz_energy_density(self, rho, T, bulk=False, dphidn=False, dphidT=False, dphidrho=False, n_disp=None)
Args:
rho (list[Profile] or list[float]) :
Particle density for each species
T (float) :
Temperature [K]
bulk (bool) :
Default False. Set to True if rho
is list[float]
dphidn (bool) :
Compute derivative
dphidT (bool) :
Compute derivative
dphidrho (bool) :
Compute derivative
n_disp (list[Profile], optional) :
Pre-computed weighted densities.
Returns:
Profile or float :
The (reduced) dispersion helmholtz energy density [-]
Weighted densities
Methods to compute various weighted densities required by the different Helmholtz energy contributions.
Table of contents
get_dispersion_weighted_density(self, rho_arr, T, bulk=False, dndrho=False, dndT=False)
Get the weighted density for the dispersion term
Args:
rho_arr (list[Profile]) :
Density profiles [particles / Å^3], indexed as rho_arr[
T (float) :
Temperature [K]
bulk (bool) :
For bulk, no convolution is required (this weighted density equals the bulk density)
dndrho (bool) :
Whether to compute the derivatives (all are equal to unity)
dndT (bool) :
Compute the derivatives
Returns:
list[Profile] :
Weighted densities, indexed as n[
optionally 2d array of ones (dndrho)
get_fmt_weighted_densities(self, rho, T, bulk=False, dndrho=False, dndT=False)
Compute weighted densities from FMT model.
Args:
rho (list[Profile] or list[float]) :
Particle densities of each species
T (float) :
Temperature [K]
bulk (bool) :
Default False, set to True if rho
is list[float]
dndrho (bool) :
Also compute derivatives (only for bulk)
dndT (bool) :
Also compute derivatives
Returns:
list[Profile] or list[float] :
The weighted densities and optionally one differential.
get_weighted_densities(self, rho, T, bulk=False, dndrho=False, dndT=False)
Compute the weighted densities.
Remember:
dphidn is a list[Profile], indexed as dphidn[
Where
Remember also :
Each weight gives rise to ONE weighted density, becuase you sum over the components.
The exception is the dispersion weight, which gives rise to ncomps weighted densities,
one for each component.
Weight functions
Get-methods for various weight functions.
Table of contents
get_dispersion_weights(self, T, dwdT=False, d=None, dd_dT=None)
Get the weights for the dispersion term
Args:
T (float) :
Temperature [K]
dwdT (bool, optimal) :
Compute derivative wrt. T? Defaults to False.
d (1d array) :
Pre-computed Barker-Henderson diameters
dd_dT (1d array) :
Pre-computed temperature derivatives of Barker-Henderson diameters.
Returns:
2d array of WeightFunction :
The weights for the dispersion weighted densities, indexed as wts[
get_fmt_weights(self, T, dwdT=False)
Get the FMT weight functions
Args:
T (float) :
Temperature [K]
dwdT (bool) :
Compute derivative instead
Returns:
list[list[WeightFunction]] :
FMT weight functions, indexed as wts[
get_weights(self, T, dwdT=False)
Get all the weights used for weighted densities in a 2D array, indexed as weight[
Args:
T (float) :
Temperature [K], used to get hard sphere diameters
dwdT (bool) :
Compute temperature differentials
Return:
2D array [Analytical] :
Weight functions, indexed as weight[
Deprecated methods
Deprecated methods are not maintained, and may be removed in the future.
Table of contents
pressure_tv(self, rho, T)