Thermodynamics

Burnman has a number of functions and classes which deal with the thermodynamics of single phases and aggregates.

Lattice Vibrations

Debye model

burnman.eos.debye.jit(fn)[source]
burnman.eos.debye.debye_fn(x)[source]

Evaluate the Debye function. Takes the parameter xi = Debye_T/T

burnman.eos.debye.debye_fn_cheb(x)[source]

Evaluate the Debye function using a Chebyshev series expansion coupled with asymptotic solutions of the function. Shamelessly adapted from the GSL implementation of the same function (Itself adapted from Collected Algorithms from ACM). Should give the same result as debye_fn(x) to near machine-precision.

burnman.eos.debye.thermal_energy(T, debye_T, n)[source]

calculate the thermal energy of a substance. Takes the temperature, the Debye temperature, and n, the number of atoms per molecule. Returns thermal energy in J/mol

burnman.eos.debye.molar_heat_capacity_v(T, debye_T, n)[source]

Heat capacity at constant volume. In J/K/mol

burnman.eos.debye.helmholtz_free_energy(T, debye_T, n)[source]

Helmholtz free energy of lattice vibrations in the Debye model. It is important to note that this does NOT include the zero point energy of vibration for the lattice. As long as you are calculating relative differences in F, this should cancel anyways. In Joules.

burnman.eos.debye.entropy(T, debye_T, n)[source]

Entropy due to lattice vibrations in the Debye model [J/K]

Einstein model

burnman.eos.einstein.thermal_energy(T, einstein_T, n)[source]

calculate the thermal energy of a substance. Takes the temperature, the Einstein temperature, and n, the number of atoms per molecule. Returns thermal energy in J/mol

burnman.eos.einstein.molar_heat_capacity_v(T, einstein_T, n)[source]

Heat capacity at constant volume. In J/K/mol

Chemistry parsing and thermodynamics

burnman.tools.chemistry.read_masses()[source]

A simple function to read a file with a two column list of elements and their masses into a dictionary

burnman.tools.chemistry.atomic_masses = {'Ag': 0.107868, 'Al': 0.0269815, 'Ar': 0.039948, 'As': 0.0749216, 'Au': 0.196967, 'B': 0.010811, 'Ba': 0.137327, 'Be': 0.00901218, 'Bi': 0.20898, 'Br': 0.079904, 'C': 0.0120107, 'Ca': 0.040078, 'Cd': 0.112411, 'Ce': 0.140116, 'Cl': 0.035453, 'Co': 0.0589332, 'Cr': 0.0519961, 'Cs': 0.132905, 'Cu': 0.063546, 'Dy': 0.1625, 'Er': 0.167259, 'Eu': 0.151964, 'F': 0.0189984, 'Fe': 0.055845, 'Ga': 0.069723, 'Gd': 0.15725, 'Ge': 0.07264, 'H': 0.00100794, 'He': 0.0040026, 'Hf': 0.17849, 'Hg': 0.20059, 'Ho': 0.16493, 'I': 0.126904, 'In': 0.114818, 'Ir': 0.192217, 'K': 0.0390983, 'Kr': 0.083798, 'La': 0.138905, 'Li': 0.006941, 'Lu': 0.174967, 'Mg': 0.024305, 'Mn': 0.054938, 'Mo': 0.09596, 'N': 0.0140067, 'Na': 0.0229898, 'Nb': 0.0929064, 'Nd': 0.144242, 'Ne': 0.0201797, 'Ni': 0.0586934, 'O': 0.0159994, 'Os': 0.19023, 'P': 0.0309738, 'Pa': 0.231036, 'Pb': 0.2072, 'Pd': 0.10642, 'Pr': 0.140908, 'Pt': 0.195084, 'Rb': 0.0854678, 'Re': 0.186207, 'Rh': 0.102905, 'Ru': 0.10107, 'S': 0.032065, 'Sb': 0.12176, 'Sc': 0.0449559, 'Se': 0.07896, 'Si': 0.0280855, 'Sm': 0.15036, 'Sn': 0.11871, 'Sr': 0.08762, 'Ta': 0.180948, 'Tb': 0.158925, 'Te': 0.1276, 'Th': 0.232038, 'Ti': 0.047867, 'Tl': 0.204383, 'Tm': 0.168934, 'U': 0.238029, 'V': 0.0509415, 'Vc': 0.0, 'W': 0.18384, 'Xe': 0.131293, 'Y': 0.0889058, 'Yb': 0.173054, 'Zn': 0.06538, 'Zr': 0.091224}

IUPAC_element_order provides a list of all the elements. Element order is based loosely on electronegativity, following the scheme suggested by IUPAC, except that H comes after the Group 16 elements, not before them.

burnman.tools.chemistry.dictionarize_formula(formula)[source]

A function to read a chemical formula string and convert it into a dictionary

Parameters
formulastring object

Chemical formula, written in the XnYm format, where the formula has n atoms of element X and m atoms of element Y

Returns
fdictionary object

The same chemical formula, but expressed as a dictionary.

burnman.tools.chemistry.sum_formulae(formulae, amounts=None)[source]

Adds together a set of formulae.

Parameters
formulaelist of dictionary or counter objects

List of chemical formulae

amountslist of floats

List of amounts of each formula

Returns
summed_formulaCounter object

The sum of the user-provided formulae

burnman.tools.chemistry.formula_mass(formula)[source]

A function to take a chemical formula and compute the formula mass.

Parameters
formuladictionary or counter object

A chemical formula

Returns
massfloat

The mass per mole of formula

burnman.tools.chemistry.convert_formula(formula, to_type='mass', normalize=False)[source]

Converts a chemical formula from one type (mass or molar) into the other. Renormalises amounts if normalize=True

Parameters
formuladictionary or counter object

A chemical formula

to_typestring, one of ‘mass’ or ‘molar’

Conversion type

normalizeboolean

Whether or not to normalize the converted formula to 1

Returns
fdictionary

The converted formula

burnman.tools.chemistry.process_solution_chemistry(solution_model)[source]

This function parses a class instance with a “formulas” attribute containing site information, e.g.

[ ‘[Mg]3[Al]2Si3O12’, ‘[Mg]3[Mg1/2Si1/2]2Si3O12’ ]

It outputs the bulk composition of each endmember (removing the site information), and also a set of variables and arrays which contain the site information. These are output in a format that can easily be used to calculate activities and gibbs free energies, given molar fractions of the phases and pressure and temperature where necessary.

Parameters
solution_modelinstance of class

Class must have a “formulas” attribute, containing a list of chemical formulae with site information

Returns
none

Nothing is returned from this function, but the solution_model object gains the following attributes.

solution_formulaelist of dictionaries

List of endmember formulae is output from site formula strings

n_sitesinteger

Number of sites in the solid solution. Should be the same for all endmembers.

siteslist of lists of strings

A list of elements for each site in the solid solution

site_nameslist of strings

A list of elements_site pairs in the solid solution, where each distinct site is given by a unique uppercase letter e.g. [‘Mg_A’, ‘Fe_A’, ‘Al_A’, ‘Al_B’, ‘Si_B’]

n_occupanciesinteger

Sum of the number of possible elements on each of the sites in the solid solution. Example: A binary solution [[A][B],[B][C1/2D1/2]] would have n_occupancies = 5, with two possible elements on Site 1 and three on Site 2

site_multiplicitiesarray of floats

The number of each site per formula unit To simplify computations later, the multiplicities are repeated for each element on each site

endmember_occupancies2d array of floats

A 1D array for each endmember in the solid solution, containing the fraction of atoms of each element on each site.

endmember_noccupancies2d array of floats

A 1D array for each endmember in the solid solution, containing the number of atoms of each element on each site per mole of endmember.

burnman.tools.chemistry.site_occupancies_to_strings(site_species_names, site_multiplicities, endmember_occupancies)[source]

Converts a list of endmember site occupancies into a list of string representations of those occupancies.

Parameters
site_species_names2D list of strings

A list of list of strings, giving the names of the species which reside on each site. List of sites, each of which contains a list of the species occupying each site.

site_multiplicitiesnumpy array of floats

List of floats giving the multiplicity of each site Must be either the same length as the number of sites, or the same length as site_species_names (with an implied repetition of the same number for each species on a given site).

endmember_occupancies2D numpy array of floats

A list of site-species occupancies for each endmember. The first dimension loops over the endmembers, and the second dimension loops over the site-species occupancies for that endmember. The total number and order of occupancies must be the same as the strings in site_species_names.

Returns
site_formulaelist of strings

A list of strings in standard burnman format. For example, [Mg]3[Al]2 would correspond to the classic two-site pyrope garnet.

burnman.tools.chemistry.compositional_array(formulae)[source]
Parameters
formulaelist of dictionaries

List of chemical formulae

Returns
formula_array2D array of floats

Array of endmember formulae

elementsList of strings

List of elements

burnman.tools.chemistry.ordered_compositional_array(formulae, elements)[source]
Parameters
formulaelist of dictionaries

List of chemical formulae

elementsList of strings

List of elements

Returns
formula_array2D array of floats

Array of endmember formulae

burnman.tools.chemistry.formula_to_string(formula)[source]
Parameters
formuladictionary or counter

Chemical formula

Returns
formula_stringstring

A formula string, with element order as given in the list IUPAC_element_order. If one or more keys in the dictionary are not one of the elements in the periodic table, then they are added at the end of the string.

burnman.tools.chemistry.sort_element_list_to_IUPAC_order(element_list)[source]
Parameters
element_listlist

List of elements

Returns
sorted_listlist

List of elements sorted into IUPAC order

burnman.tools.chemistry.convert_fractions(composite, phase_fractions, input_type, output_type)[source]

Takes a composite with a set of user defined molar, volume or mass fractions (which do not have to be the fractions currently associated with the composite) and converts the fractions to molar, mass or volume.

Conversions to and from mass require a molar mass to be defined for all phases. Conversions to and from volume require set_state to have been called for the composite.

Parameters
compositeComposite

Composite for which fractions are to be defined.

phase_fractionslist of floats

List of input phase fractions (of type input_type)

input_typestring

Input fraction type: ‘molar’, ‘mass’ or ‘volume’

output_typestring

Output fraction type: ‘molar’, ‘mass’ or ‘volume’

Returns
output_fractionslist of floats

List of output phase fractions (of type output_type)

burnman.tools.chemistry.chemical_potentials(assemblage, component_formulae)[source]

The compositional space of the components does not have to be a superset of the compositional space of the assemblage. Nor do they have to compose an orthogonal basis.

The components must each be described by a linear mineral combination

The mineral compositions must be linearly independent

Parameters
assemblagelist of classes

List of material classes set_method and set_state should already have been used the composition of the solid solutions should also have been set

component_formulaelist of dictionaries

List of chemical component formula dictionaries No restriction on length

Returns
component_potentialsarray of floats

Array of chemical potentials of components

burnman.tools.chemistry.fugacity(standard_material, assemblage)[source]
Parameters
standard_material: class

Material class set_method and set_state should already have been used material must have a formula as a dictionary parameter

assemblage: list of classes

List of material classes set_method and set_state should already have been used

Returns
fugacityfloat

Value of the fugacity of the component with respect to the standard material

burnman.tools.chemistry.relative_fugacity(standard_material, assemblage, reference_assemblage)[source]
Parameters
standard_material: class

Material class set_method and set_state should already have been used material must have a formula as a dictionary parameter

assemblage: list of classes

List of material classes set_method and set_state should already have been used

reference_assemblage: list of classes

List of material classes set_method and set_state should already have been used

Returns
relative_fugacityfloat

Value of the fugacity of the component in the assemblage with respect to the reference_assemblage

burnman.tools.chemistry.equilibrium_pressure(minerals, stoichiometry, temperature, pressure_initial_guess=100000.0)[source]

Given a list of minerals, their reaction stoichiometries and a temperature of interest, compute the equilibrium pressure of the reaction.

Parameters
mineralslist of minerals

List of minerals involved in the reaction.

stoichiometrylist of floats

Reaction stoichiometry for the minerals provided. Reactants and products should have the opposite signs [mol]

temperaturefloat

Temperature of interest [K]

pressure_initial_guessoptional float

Initial pressure guess [Pa]

Returns
pressurefloat

The equilibrium pressure of the reaction [Pa]

burnman.tools.chemistry.equilibrium_temperature(minerals, stoichiometry, pressure, temperature_initial_guess=1000.0)[source]

Given a list of minerals, their reaction stoichiometries and a pressure of interest, compute the equilibrium temperature of the reaction.

Parameters
mineralslist of minerals

List of minerals involved in the reaction.

stoichiometrylist of floats

Reaction stoichiometry for the minerals provided. Reactants and products should have the opposite signs [mol]

pressurefloat

Pressure of interest [Pa]

temperature_initial_guessoptional float

Initial temperature guess [K]

Returns
temperaturefloat

The equilibrium temperature of the reaction [K]

burnman.tools.chemistry.invariant_point(minerals_r1, stoichiometry_r1, minerals_r2, stoichiometry_r2, pressure_temperature_initial_guess=[1000000000.0, 1000.0])[source]

Given a list of minerals, their reaction stoichiometries and a pressure of interest, compute the equilibrium temperature of the reaction.

Parameters
mineralslist of minerals

List of minerals involved in the reaction.

stoichiometrylist of floats

Reaction stoichiometry for the minerals provided. Reactants and products should have the opposite signs [mol]

pressurefloat

Pressure of interest [Pa]

temperature_initial_guessoptional float

Initial temperature guess [K]

Returns
temperaturefloat

The equilibrium temperature of the reaction [K]

burnman.tools.chemistry.hugoniot(mineral, P_ref, T_ref, pressures, reference_mineral=None)[source]

Calculates the temperatures (and volumes) along a Hugoniot as a function of pressure according to the Hugoniot equation U2-U1 = 0.5*(p2 - p1)(V1 - V2) where U and V are the internal energies and volumes (mass or molar) and U = F + TS

Parameters
mineralmineral

Mineral for which the Hugoniot is to be calculated.

P_reffloat

Reference pressure [Pa]

T_reffloat

Reference temperature [K]

pressuresnumpy array of floats

Set of pressures [Pa] for which the Hugoniot temperature and volume should be calculated

reference_mineralmineral

Mineral which is stable at the reference conditions Provides an alternative U_0 and V_0 when the reference mineral transforms to the mineral of interest at some (unspecified) pressure.

Returns
temperaturesnumpy array of floats

The Hugoniot temperatures at pressure

volumesnumpy array of floats

The Hugoniot volumes at pressure

Equilibrium Thermodynamics

burnman.equilibrate(composition, assemblage, equality_constraints, tol=0.001, store_iterates=False, store_assemblage=True, max_iterations=100.0, verbose=False)[source]

A function that equilibrates an assemblage subject to given bulk composition and equality constraints by solving the equilibrium relations (chemical affinities for feasible reactions in the system should be equal to zero).

Parameters
compositiondictionary

The bulk composition that the assemblage must satisfy

assemblageburnman.Composite object

The assemblage to be equilibrated

equality_constraintslist

A list of equality constraints. Each constraint should have the form: [<constraint type>, <constraint>], where <constraint type> is one of P, T, S, V, X, PT_ellipse, phase_fraction, or phase_composition. The <constraint> object should either be a float or an array of floats for P, T, S, V (representing the desired pressure, temperature, entropy or volume of the material). If the constraint type is X (a generic constraint on the solution vector) then the constraint c is represented by the following equality: np.dot(c[0], x) - c[1]. If the constraint type is PT_ellipse, the equality is given by norm(([P, T] - c[0])/c[1]) - 1. The constraint_type phase_fraction assumes a tuple of the phase object (which must be one of the phases in the burnman.Composite) and a float or vector corresponding to the phase fractions. Finally, a phase_composition constraint has the format (site_names, n, d, v), where n*x/d*x = v and n and d are fixed vectors of site coefficients. So, one could for example choose a constraint ([Mg_A, Fe_A], [1., 0.], [1., 1.], [0.5]) which would correspond to equal amounts Mg and Fe on the A site.

tolfloat

The tolerance for the nonlinear solver.

store_iteratesboolean

Whether to store the parameter values for each iteration in each solution object.

store_assemblageboolean

Whether to store a copy of the assemblage object in each solution object.

max_iterationsinteger

The maximum number of iterations for the nonlinear solver.

verboseboolean

Whether to print output updating the user on the status of equilibration.

Returns
sol_listsingle, list, or 2D list of solver solution objects
prmnamedtuple object

A tuple with attributes n_parameters (the number of parameters for the current equilibrium problem) and phase_amount_indices (the indices of the parameters that correspond to phase amounts).