Solution models

SolidSolution objects in Burnman (type SolidSolution) take one of several methods which define the properties of the solution.

Inheritance diagram of burnman.classes.solutionmodel.MechanicalSolution, burnman.classes.solutionmodel.IdealSolution, burnman.classes.solutionmodel.AsymmetricRegularSolution, burnman.classes.solutionmodel.SymmetricRegularSolution, burnman.classes.solutionmodel.SubregularSolution

Base class

class burnman.SolidSolution(name=None, solution_type=None, endmembers=None, energy_interaction=None, volume_interaction=None, entropy_interaction=None, energy_ternary_terms=None, volume_ternary_terms=None, entropy_ternary_terms=None, alphas=None, molar_fractions=None)[source]

Bases: burnman.classes.mineral.Mineral

This is the base class for all solid solutions. Site occupancies, endmember activities and the constant and pressure and temperature dependencies of the excess properties can be queried after using set_composition() States of the solid solution can only be queried after setting the pressure, temperature and composition using set_state().

This class is available as burnman.SolidSolution. It uses an instance of burnman.SolutionModel to calculate interaction terms between endmembers.

All the solid solution parameters are expected to be in SI units. This means that the interaction parameters should be in J/mol, with the T and P derivatives in J/K/mol and m^3/mol.

The parameters are relevant to all solution models. Please see the documentation for individual models for details about other parameters.

Parameters
namestring

Name of the solid solution

solution_typestring

String determining which SolutionModel to use. One of ‘mechanical’, ‘ideal’, ‘symmetric’, ‘asymmetric’ or ‘subregular’.

endmemberslist of lists

List of endmembers in this solid solution. The first item of each list should be a burnman.Mineral object. The second item should be a string with the site formula of the endmember.

molar_fractionsnumpy array (optional)

The molar fractions of each endmember in the solid solution. Can be reset using the set_composition() method.

property name

Human-readable name of this material.

By default this will return the name of the class, but it can be set to an arbitrary string. Overriden in Mineral.

get_endmembers()[source]
set_composition(molar_fractions)[source]

Set the composition for this solid solution. Resets cached properties.

Parameters
molar_fractions: list of float

molar abundance for each endmember, needs to sum to one.

set_method(method)[source]

Set the equation of state to be used for this mineral. Takes a string corresponding to any of the predefined equations of state: ‘bm2’, ‘bm3’, ‘mgd2’, ‘mgd3’, ‘slb2’, ‘slb3’, ‘mt’, ‘hp_tmt’, or ‘cork’. Alternatively, you can pass a user defined class which derives from the equation_of_state base class. After calling set_method(), any existing derived properties (e.g., elastic parameters or thermodynamic potentials) will be out of date, so set_state() will need to be called again.

set_state(pressure, temperature)[source]

(copied from set_state):

Set the material to the given pressure and temperature.

Parameters
pressurefloat

The desired pressure in [Pa].

temperaturefloat

The desired temperature in [K].

property formula

Returns molar chemical formula of the solid solution.

property activities

Returns a list of endmember activities [unitless].

property activity_coefficients

Returns a list of endmember activity coefficients (gamma = activity / ideal activity) [unitless].

property molar_internal_energy

Returns molar internal energy of the mineral [J/mol]. Aliased with self.energy

property excess_partial_gibbs

Returns excess partial molar gibbs free energy [J/mol]. Property specific to solid solutions.

property excess_partial_volumes

Returns excess partial volumes [m^3]. Property specific to solid solutions.

property excess_partial_entropies

Returns excess partial entropies [J/K]. Property specific to solid solutions.

property partial_gibbs

Returns excess partial molar gibbs free energy [J/mol]. Property specific to solid solutions.

property partial_volumes

Returns excess partial volumes [m^3]. Property specific to solid solutions.

property partial_entropies

Returns excess partial entropies [J/K]. Property specific to solid solutions.

property excess_gibbs

Returns molar excess gibbs free energy [J/mol]. Property specific to solid solutions.

property gibbs_hessian

Returns an array containing the second compositional derivative of the Gibbs free energy [J]. Property specific to solid solutions.

property entropy_hessian

Returns an array containing the second compositional derivative of the entropy [J/K]. Property specific to solid solutions.

property volume_hessian

Returns an array containing the second compositional derivative of the volume [m^3]. Property specific to solid solutions.

property molar_gibbs

Returns molar Gibbs free energy of the solid solution [J/mol]. Aliased with self.gibbs.

property molar_helmholtz

Returns molar Helmholtz free energy of the solid solution [J/mol]. Aliased with self.helmholtz.

property molar_mass

Returns molar mass of the solid solution [kg/mol].

property excess_volume

Returns excess molar volume of the solid solution [m^3/mol]. Specific property for solid solutions.

property molar_volume

Returns molar volume of the solid solution [m^3/mol]. Aliased with self.V.

property density

Returns density of the solid solution [kg/m^3]. Aliased with self.rho.

property excess_entropy

Returns excess molar entropy [J/K/mol]. Property specific to solid solutions.

property molar_entropy

Returns molar entropy of the solid solution [J/K/mol]. Aliased with self.S.

property excess_enthalpy

Returns excess molar enthalpy [J/mol]. Property specific to solid solutions.

property molar_enthalpy

Returns molar enthalpy of the solid solution [J/mol]. Aliased with self.H.

property isothermal_bulk_modulus

Returns isothermal bulk modulus of the solid solution [Pa]. Aliased with self.K_T.

property adiabatic_bulk_modulus

Returns adiabatic bulk modulus of the solid solution [Pa]. Aliased with self.K_S.

property isothermal_compressibility

Returns isothermal compressibility of the solid solution. (or inverse isothermal bulk modulus) [1/Pa]. Aliased with self.K_T.

property adiabatic_compressibility

Returns adiabatic compressibility of the solid solution. (or inverse adiabatic bulk modulus) [1/Pa]. Aliased with self.K_S.

property shear_modulus

Returns shear modulus of the solid solution [Pa]. Aliased with self.G.

property p_wave_velocity

Returns P wave speed of the solid solution [m/s]. Aliased with self.v_p.

property bulk_sound_velocity

Returns bulk sound speed of the solid solution [m/s]. Aliased with self.v_phi.

property shear_wave_velocity

Returns shear wave speed of the solid solution [m/s]. Aliased with self.v_s.

property grueneisen_parameter

Returns grueneisen parameter of the solid solution [unitless]. Aliased with self.gr.

property thermal_expansivity

Returns thermal expansion coefficient (alpha) of the solid solution [1/K]. Aliased with self.alpha.

property molar_heat_capacity_v

Returns molar heat capacity at constant volume of the solid solution [J/K/mol]. Aliased with self.C_v.

property C_p

Alias for molar_heat_capacity_p()

property C_v

Alias for molar_heat_capacity_v()

property G

Alias for shear_modulus()

property H

Alias for molar_enthalpy()

property K_S

Alias for adiabatic_bulk_modulus()

property K_T

Alias for isothermal_bulk_modulus()

property P

Alias for pressure()

property S

Alias for molar_entropy()

property T

Alias for temperature()

property V

Alias for molar_volume()

property alpha

Alias for thermal_expansivity()

property beta_S

Alias for adiabatic_compressibility()

property beta_T

Alias for isothermal_compressibility()

copy()
debug_print(indent='')

Print a human-readable representation of this Material.

property energy

Alias for molar_internal_energy()

evaluate(vars_list, pressures, temperatures)

Returns an array of material properties requested through a list of strings at given pressure and temperature conditions. At the end it resets the set_state to the original values. The user needs to call set_method() before.

Parameters
vars_listlist of strings

Variables to be returned for given conditions

pressuresndlist or ndarray of float

n-dimensional array of pressures in [Pa].

temperaturesndlist or ndarray of float

n-dimensional array of temperatures in [K].

Returns
outputarray of array of float

Array returning all variables at given pressure/temperature values. output[i][j] is property vars_list[j] and temperatures[i] and pressures[i].

property gibbs

Alias for molar_gibbs()

property gr

Alias for grueneisen_parameter()

property helmholtz

Alias for molar_helmholtz()

property molar_heat_capacity_p

Returns molar heat capacity at constant pressure of the solid solution [J/K/mol]. Aliased with self.C_p.

property pressure

Returns current pressure that was set with set_state().

Returns
pressurefloat

Pressure in [Pa].

Notes

  • Aliased with P().

print_minerals_of_current_state()

Print a human-readable representation of this Material at the current P, T as a list of minerals. This requires set_state() has been called before.

reset()

Resets all cached material properties.

It is typically not required for the user to call this function.

property rho

Alias for density()

property temperature

Returns current temperature that was set with set_state().

Returns
temperaturefloat

Temperature in [K].

Notes

  • Aliased with T().

to_string()

Returns the name of the mineral class

unroll()

Unroll this material into a list of burnman.Mineral and their molar fractions. All averaging schemes then operate on this list of minerals. Note that the return value of this function may depend on the current state (temperature, pressure).

Returns
fractionslist of float

List of molar fractions, should sum to 1.0.

mineralslist of burnman.Mineral

List of minerals.

Notes

Needs to be implemented in derived classes.

property v_p

Alias for p_wave_velocity()

property v_phi

Alias for bulk_sound_velocity()

property v_s

Alias for shear_wave_velocity()

property stoichiometric_matrix[source]

A sympy Matrix where each element M[i,j] corresponds to the number of atoms of element[j] in endmember[i].

property stoichiometric_array[source]

An array where each element arr[i,j] corresponds to the number of atoms of element[j] in endmember[i].

property reaction_basis[source]

An array where each element arr[i,j] corresponds to the number of moles of endmember[j] involved in reaction[i].

property n_reactions[source]

The number of reactions in reaction_basis.

property independent_element_indices[source]

A list of an independent set of element indices. If the amounts of these elements are known (element_amounts), the amounts of the other elements can be inferred by -compositional_null_basis[independent_element_indices].dot(element_amounts).

property dependent_element_indices[source]

The element indices not included in the independent list.

property compositional_null_basis[source]

An array N such that N.b = 0 for all bulk compositions that can be produced with a linear sum of the endmembers in the solid solution.

property endmember_formulae[source]

A list of formulae for all the endmember in the solid solution.

property endmember_names[source]

A list of names for all the endmember in the solid solution.

property n_endmembers[source]

The number of endmembers in the solid solution.

property elements[source]

A list of the elements which could be contained in the solid solution, returned in the IUPAC element order.

class burnman.SolutionModel[source]

Bases: object

This is the base class for a solution model, intended for use in defining solid solutions and performing thermodynamic calculations on them. All minerals of type burnman.SolidSolution use a solution model for defining how the endmembers in the solid solution interact.

A user wanting a new solution model should define the functions included in the base class. All of the functions in the base class return zero, so if the user-defined solution model does not implement them, they essentially have no effect, and the Gibbs free energy and molar volume of a solid solution will be equal to the weighted arithmetic averages of the different endmember values.

excess_gibbs_free_energy(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess Gibbs free energy of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
G_excessfloat

The excess Gibbs free energy

excess_volume(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess volume of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
V_excessfloat

The excess volume of the solution

excess_entropy(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess entropy of the solution. The base class implementation assumes that the excess entropy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
S_excessfloat

The excess entropy of the solution

excess_enthalpy(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess enthalpy of the solution. The base class implementation assumes that the excess enthalpy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
H_excessfloat

The excess enthalpy of the solution

excess_partial_gibbs_free_energies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess Gibbs free energy for each endmember of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_G_excessnumpy array

The excess Gibbs free energy of each endmember

excess_partial_entropies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess entropy for each endmember of the solution. The base class implementation assumes that the excess entropy is zero (true for mechanical solutions).

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_S_excessnumpy array

The excess entropy of each endmember

excess_partial_volumes(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess volume for each endmember of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_V_excessnumpy array

The excess volume of each endmember

Mechanical solution

class burnman.classes.solutionmodel.MechanicalSolution(endmembers)[source]

Bases: burnman.classes.solutionmodel.SolutionModel

An extremely simple class representing a mechanical solution model. A mechanical solution experiences no interaction between endmembers. Therefore, unlike ideal solutions there is no entropy of mixing; the total gibbs free energy of the solution is equal to the dot product of the molar gibbs free energies and molar fractions of the constituent materials.

excess_gibbs_free_energy(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess Gibbs free energy of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
G_excessfloat

The excess Gibbs free energy

excess_volume(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess volume of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
V_excessfloat

The excess volume of the solution

excess_entropy(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess entropy of the solution. The base class implementation assumes that the excess entropy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
S_excessfloat

The excess entropy of the solution

excess_enthalpy(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess enthalpy of the solution. The base class implementation assumes that the excess enthalpy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
H_excessfloat

The excess enthalpy of the solution

excess_partial_gibbs_free_energies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess Gibbs free energy for each endmember of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_G_excessnumpy array

The excess Gibbs free energy of each endmember

excess_partial_volumes(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess volume for each endmember of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_V_excessnumpy array

The excess volume of each endmember

excess_partial_entropies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess entropy for each endmember of the solution. The base class implementation assumes that the excess entropy is zero (true for mechanical solutions).

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_S_excessnumpy array

The excess entropy of each endmember

activity_coefficients(pressure, temperature, molar_fractions)[source]
activities(pressure, temperature, molar_fractions)[source]

Ideal solution

class burnman.classes.solutionmodel.IdealSolution(endmembers)[source]

Bases: burnman.classes.solutionmodel.SolutionModel

A very simple class representing an ideal solution model. Calculate the excess gibbs free energy and entropy due to configurational entropy, excess volume is equal to zero.

excess_partial_gibbs_free_energies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess Gibbs free energy for each endmember of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_G_excessnumpy array

The excess Gibbs free energy of each endmember

excess_partial_entropies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess entropy for each endmember of the solution. The base class implementation assumes that the excess entropy is zero (true for mechanical solutions).

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_S_excessnumpy array

The excess entropy of each endmember

excess_partial_volumes(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess volume for each endmember of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_V_excessnumpy array

The excess volume of each endmember

gibbs_hessian(pressure, temperature, molar_fractions)[source]
entropy_hessian(pressure, temperature, molar_fractions)[source]
volume_hessian(pressure, temperature, molar_fractions)[source]
activity_coefficients(pressure, temperature, molar_fractions)[source]
activities(pressure, temperature, molar_fractions)[source]
excess_enthalpy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess enthalpy of the solution. The base class implementation assumes that the excess enthalpy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
H_excessfloat

The excess enthalpy of the solution

excess_entropy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess entropy of the solution. The base class implementation assumes that the excess entropy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
S_excessfloat

The excess entropy of the solution

excess_gibbs_free_energy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess Gibbs free energy of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
G_excessfloat

The excess Gibbs free energy

excess_volume(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess volume of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
V_excessfloat

The excess volume of the solution

Asymmetric regular solution

class burnman.classes.solutionmodel.AsymmetricRegularSolution(endmembers, alphas, energy_interaction, volume_interaction=None, entropy_interaction=None)[source]

Bases: burnman.classes.solutionmodel.IdealSolution

Solution model implementing the asymmetric regular solution model formulation as described in [HollandPowell03].

The excess nonconfigurational Gibbs energy is given by the expression:

\[\mathcal{G}_{\textrm{excess}} = \alpha^T p (\phi^T W \phi)\]

\(\alpha\) is a vector of van Laar parameters governing asymmetry in the excess properties.

\[\phi_i = \frac{\alpha_i p_i}{\sum_{k=1}^{n} \alpha_k p_k}, W_{ij} = \frac{2 w_{ij}}{\alpha_i + \alpha_j} \textrm{for i<j}\]
excess_partial_gibbs_free_energies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess Gibbs free energy for each endmember of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_G_excessnumpy array

The excess Gibbs free energy of each endmember

excess_partial_entropies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess entropy for each endmember of the solution. The base class implementation assumes that the excess entropy is zero (true for mechanical solutions).

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_S_excessnumpy array

The excess entropy of each endmember

excess_partial_volumes(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess volume for each endmember of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_V_excessnumpy array

The excess volume of each endmember

gibbs_hessian(pressure, temperature, molar_fractions)[source]
entropy_hessian(pressure, temperature, molar_fractions)[source]
volume_hessian(pressure, temperature, molar_fractions)[source]
activity_coefficients(pressure, temperature, molar_fractions)[source]
activities(pressure, temperature, molar_fractions)[source]
excess_enthalpy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess enthalpy of the solution. The base class implementation assumes that the excess enthalpy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
H_excessfloat

The excess enthalpy of the solution

excess_entropy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess entropy of the solution. The base class implementation assumes that the excess entropy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
S_excessfloat

The excess entropy of the solution

excess_gibbs_free_energy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess Gibbs free energy of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
G_excessfloat

The excess Gibbs free energy

excess_volume(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess volume of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
V_excessfloat

The excess volume of the solution

Symmetric regular solution

class burnman.classes.solutionmodel.SymmetricRegularSolution(endmembers, energy_interaction, volume_interaction=None, entropy_interaction=None)[source]

Bases: burnman.classes.solutionmodel.AsymmetricRegularSolution

Solution model implementing the symmetric regular solution model. This is a special case of the burnman.solutionmodel.AsymmetricRegularSolution class.

activities(pressure, temperature, molar_fractions)
activity_coefficients(pressure, temperature, molar_fractions)
entropy_hessian(pressure, temperature, molar_fractions)
excess_enthalpy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess enthalpy of the solution. The base class implementation assumes that the excess enthalpy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
H_excessfloat

The excess enthalpy of the solution

excess_entropy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess entropy of the solution. The base class implementation assumes that the excess entropy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
S_excessfloat

The excess entropy of the solution

excess_gibbs_free_energy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess Gibbs free energy of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
G_excessfloat

The excess Gibbs free energy

excess_partial_entropies(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess entropy for each endmember of the solution. The base class implementation assumes that the excess entropy is zero (true for mechanical solutions).

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_S_excessnumpy array

The excess entropy of each endmember

excess_partial_gibbs_free_energies(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess Gibbs free energy for each endmember of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_G_excessnumpy array

The excess Gibbs free energy of each endmember

excess_partial_volumes(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess volume for each endmember of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_V_excessnumpy array

The excess volume of each endmember

excess_volume(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess volume of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
V_excessfloat

The excess volume of the solution

gibbs_hessian(pressure, temperature, molar_fractions)
volume_hessian(pressure, temperature, molar_fractions)

Subregular solution

class burnman.classes.solutionmodel.SubregularSolution(endmembers, energy_interaction, volume_interaction=None, entropy_interaction=None, energy_ternary_terms=None, volume_ternary_terms=None, entropy_ternary_terms=None)[source]

Bases: burnman.classes.solutionmodel.IdealSolution

Solution model implementing the subregular solution model formulation as described in [HW89]. The excess conconfigurational Gibbs energy is given by the expression:

\[\mathcal{G}_{\textrm{excess}} = \sum_i \sum_{j > i} (p_i p_j^2 W_{ij} + p_j p_i^2 W_{ji} + \sum_{k > j > i} p_i p_j p_k W_{ijk})\]

Interaction parameters are inserted into a 3D interaction matrix during initialization to make use of numpy vector algebra.

Parameters
endmemberslist of lists

A list of all the independent endmembers in the solution. The first item of each list gives the Mineral object corresponding to the endmember. The second item gives the site-species formula.

energy_interactionlist of list of lists

The binary endmember interaction energies. Each interaction[i, j-i-1, 0] corresponds to W(i,j), while interaction[i, j-i-1, 1] corresponds to W(j,i).

volume_interactionlist of list of lists

The binary endmember interaction volumes. Each interaction[i, j-i-1, 0] corresponds to W(i,j), while interaction[i, j-i-1, 1] corresponds to W(j,i).

entropy_interactionlist of list of lists

The binary endmember interaction entropies. Each interaction[i, j-i-1, 0] corresponds to W(i,j), while interaction[i, j-i-1, 1] corresponds to W(j,i).

energy_ternary_termslist of lists

The ternary interaction energies. Each list should contain four entries: the indices i, j, k and the value of the interaction.

volume_ternary_termslist of lists

The ternary interaction volumes. Each list should contain four entries: the indices i, j, k and the value of the interaction.

entropy_ternary_termslist of lists

The ternary interaction entropies. Each list should contain four entries: the indices i, j, k and the value of the interaction.

excess_partial_gibbs_free_energies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess Gibbs free energy for each endmember of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_G_excessnumpy array

The excess Gibbs free energy of each endmember

excess_partial_entropies(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess entropy for each endmember of the solution. The base class implementation assumes that the excess entropy is zero (true for mechanical solutions).

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_S_excessnumpy array

The excess entropy of each endmember

excess_partial_volumes(pressure, temperature, molar_fractions)[source]

Given a list of molar fractions of different phases, compute the excess volume for each endmember of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
partial_V_excessnumpy array

The excess volume of each endmember

gibbs_hessian(pressure, temperature, molar_fractions)[source]
entropy_hessian(pressure, temperature, molar_fractions)[source]
volume_hessian(pressure, temperature, molar_fractions)[source]
activity_coefficients(pressure, temperature, molar_fractions)[source]
excess_enthalpy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess enthalpy of the solution. The base class implementation assumes that the excess enthalpy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
H_excessfloat

The excess enthalpy of the solution

excess_entropy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess entropy of the solution. The base class implementation assumes that the excess entropy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
S_excessfloat

The excess entropy of the solution

excess_gibbs_free_energy(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess Gibbs free energy of the solution. The base class implementation assumes that the excess gibbs free energy is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
G_excessfloat

The excess Gibbs free energy

excess_volume(pressure, temperature, molar_fractions)

Given a list of molar fractions of different phases, compute the excess volume of the solution. The base class implementation assumes that the excess volume is zero.

Parameters
pressurefloat

Pressure at which to evaluate the solution model. [Pa]

temperaturefloat

Temperature at which to evaluate the solution. [K]

molar_fractionslist of floats

List of molar fractions of the different endmembers in solution

Returns
V_excessfloat

The excess volume of the solution

activities(pressure, temperature, molar_fractions)[source]

Solution tools

burnman.tools.solution.transform_solution_to_new_basis(solution, new_basis, n_mbrs=None, solution_name=None, endmember_names=None, molar_fractions=None)[source]

Transforms a solution model from one endmember basis to another. Returns a new SolidSolution object.

Parameters
solutionburnman.SolidSolution object

The original solution object.

new_basis2D numpy array

The new endmember basis, given as amounts of the old endmembers.

n_mbrsfloat (optional)

The number of endmembers in the new solution (defaults to the length of new_basis)

solution_namestring (optional)

A name corresponding to the new solution

endmember_nameslist of strings (optional)

A list corresponding to the names of the new endmembers.

molar_fractionsnumpy array (optional)

Fractions of the new endmembers in the new solution.

Returns
solutionburnman.SolidSolution object

The transformed solid solution