Skip to content

Commit 8e23d8c

Browse files
committed
improve docs
1 parent 0fc7b7c commit 8e23d8c

File tree

2 files changed

+53
-30
lines changed

2 files changed

+53
-30
lines changed

burnman/tools/thermobarometry.py

Lines changed: 52 additions & 29 deletions
Original file line numberDiff line numberDiff line change
@@ -170,13 +170,19 @@ def assemblage_set_state_from_params(assemblage, params):
170170
Set the state of the assemblage (P, T, compositions) from the given
171171
list of parameters.
172172
173-
:param assemblage: The mineral assemblage to set the state for.
173+
:param assemblage: The mineral assemblage for which to set the state.
174+
If there are free compositional vectors for any solution phases
175+
(e.g., to account for an unknown state of order at fixed composition),
176+
the corresponding phases must have the attributes
177+
`baseline_composition` and `free_compositional_vectors` set,
178+
and the length of params must match the number of
179+
free compositional vectors plus two (for P and T).
174180
:type assemblage: Assemblage
175181
176182
:param params: List of parameters, where the first two are pressure (Pa)
177183
and temperature (K), and any additional parameters correspond to
178184
compositional degrees of freedom. Each compositional degree of freedom
179-
corresponds to a free compositional vector defined on one of the
185+
corresponds to a free compositional vector in one of the
180186
solution phases in the assemblage.
181187
:type params: list or np.array
182188
@@ -227,7 +233,7 @@ def assemblage_affinity_covariance_matrix(
227233
228234
:param include_state_uncertainties: If True, includes the contribution
229235
from uncertainties in pressure and temperature to the covariance matrix.
230-
If True, the assemblage must have its attribute `state_covariances` set.
236+
If True, the assemblage must have the attribute `state_covariances`.
231237
:type include_state_uncertainties: bool
232238
233239
:return: Covariance matrix of the affinities of the independent reactions.
@@ -269,6 +275,11 @@ def assemblage_affinity_misfit(
269275
:type dataset_covariances: dict, with keys 'endmember_names' and 'covariance_matrix'.
270276
Default is None, in which case only compositional uncertainties are considered.
271277
278+
:param include_state_uncertainties: If True, includes the contribution
279+
from uncertainties in pressure and temperature to the covariance matrix.
280+
If True, the assemblage must have the attribute `state_covariances`.
281+
:type include_state_uncertainties: bool
282+
272283
:return: Chi-squared misfit value.
273284
:rtype: float
274285
"""
@@ -284,7 +295,8 @@ def assemblage_state_misfit(assemblage):
284295
"""
285296
Compute the objective misfit function (chi-squared) for given P and T
286297
based on prior expectations for P and T and their covariance.
287-
The assemblage must have attributes `state_priors` and `state_inverse_covariances`.
298+
The assemblage must have attributes `state_priors` and
299+
`state_inverse_covariances`.
288300
289301
:param assemblage: The mineral assemblage for which to compute the misfit.
290302
:type assemblage: Assemblage
@@ -314,29 +326,35 @@ def estimate_conditions(
314326
max_it=100,
315327
):
316328
"""
317-
Perform a least-squares inversion to find the optimal pressure and temperature
318-
for a given mineral assemblage. Algorithm modified from Powell and Holland (1994).
319-
320-
:param assemblage: The mineral assemblage for which to perform the inversion.
321-
Each solution phase in the assemblage must have its composition set along with
322-
its compositional covariance matrix (called `compositional_covariances`).
323-
If there are compositional degrees of freedom, they can be added by setting
324-
the attribute `free_compositional_vectors` on the relevant solution phases, where each
325-
row of the array corresponds to a free compositional vector, and the columns correspond
326-
to the amounts of endmembers of that phase in each vector.
329+
Perform a least-squares inversion to find the optimal pressure and
330+
temperature for a given mineral assemblage.
331+
Algorithm modified from Powell and Holland (1994).
332+
333+
:param assemblage: The mineral assemblage for which to perform the
334+
inversion. Each solution phase in the assemblage must have its
335+
composition set along with its compositional covariance matrix
336+
(called `compositional_covariances`).
337+
If there are compositional degrees of freedom, they can be added by
338+
setting the attribute `free_compositional_vectors` on the relevant
339+
solution phases, where each row of the array corresponds to a
340+
free compositional vector, and the columns correspond to the
341+
amounts of endmembers of that phase in each vector.
327342
:type assemblage: Assemblage
328343
329-
:param dataset_covariances: The covariance data from the thermodynamic dataset.
330-
:type dataset_covariances: dict, with keys 'endmember_names' and 'covariance_matrix'.
331-
Default is None, in which case only compositional uncertainties are considered.
344+
:param dataset_covariances: The covariance data from the thermodynamic
345+
dataset.
346+
:type dataset_covariances: dict, with keys 'endmember_names' and
347+
'covariance_matrix'. Default is None, in which case only
348+
compositional uncertainties are considered.
332349
333-
:param include_state_misfit: If True, includes the misfit from prior expectations on P and T.
334-
The assemblage must also have attributes `state_priors` and `state_inverse_covariances`.
350+
:param include_state_misfit: If True, includes the misfit from
351+
prior expectations on P and T. The assemblage must also have
352+
attributes `state_priors` and `state_inverse_covariances`.
335353
:type include_state_misfit: bool
336354
337-
:param guessed_conditions: Initial guess for pressure (Pa) and temperature (K).
338-
If not provided, the initial guess will be taken from the current
339-
state of the assemblage.
355+
:param guessed_conditions: Initial guess for pressure (Pa) and
356+
temperature (K). If not provided, the initial guess will be taken
357+
from the current state of the assemblage.
340358
:type guessed_conditions: np.array of shape (2,), optional, default None
341359
342360
:param pressure_bounds: Bounds for pressure (Pa) during optimization.
@@ -348,17 +366,21 @@ def estimate_conditions(
348366
:param P_scaling: Scaling factor for pressure to improve numerical stability.
349367
:type P_scaling: float
350368
351-
:param small_fraction_tol: If > 0.0, reduces the number of endmembers in solution phases by
352-
transforming to a smaller set of independent endmembers using a greedy algorithm
353-
and excluding those with molar fractions smaller than this value during the inversion.
369+
:param small_fraction_tol: If > 0.0, reduces the number of endmembers in
370+
solution phases by transforming to a smaller set of independent
371+
endmembers using a greedy algorithm and excluding those with
372+
molar fractions smaller than this value during the inversion.
354373
:type small_fraction_tol: float, optional, default 0.0
355374
356375
:param max_it: Maximum number of iterations for the optimization algorithm.
357376
:type max_it: int, optional, default 100
358377
359-
:return: Result object from the optimization containing optimal P and T (x) and other properties
360-
of the solution. These include the covariance matrix of the estimated parameters (xcov),
361-
the correlation coefficient between P and T (xcorr),
378+
:return: Result object from the optimization containing the optimal
379+
conditions (x, which includes P, T, and any free compositional vectors)
380+
and other properties of the solution. These include
381+
the covariance matrix of the estimated parameters (xcov),
382+
the standard deviations of the estimated parameters (var),
383+
the correlation matrix (xcorr),
362384
the affinities of the independent reactions at the optimal conditions (affinities),
363385
the covariance matrix of the affinities (acov),
364386
the affinities weighted by the inverse square root of their covariance matrix (weighted_affinities),
@@ -479,7 +501,8 @@ def chisqr(args):
479501
res.dadx.T.dot(np.linalg.pinv(res.acov)).dot(res.dadx)
480502
+ (assemblage.state_inverse_covariances if include_state_misfit else 0)
481503
)
482-
res.xcorr = res.xcov[0, 1] / np.sqrt(res.xcov[0, 0] * res.xcov[1, 1])
504+
res.var = np.sqrt(np.diag(res.xcov))
505+
res.xcorr = res.xcov / np.outer(res.var, res.var)
483506

484507
res.affinities = R.dot(assemblage.endmember_partial_gibbs)
485508
res.weighted_affinities = np.linalg.pinv(sqrtm(res.acov)).dot(res.affinities)

examples/example_optimal_thermobarometry.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -310,7 +310,7 @@
310310
f" Temperature: {assemblage.temperature:.0f} "
311311
f"+/- {np.sqrt(res.xcov[1, 1]):.0f} K"
312312
)
313-
print(f" Correlation between P and T: {res.xcorr:.2f}")
313+
print(f" Correlation between P and T: {res.xcorr[0][1]:.2f}")
314314
print(f" Number of Reactions: {res.n_reactions}")
315315
print(f" Number of Parameters: {res.n_params}")
316316
print(f" Degrees of Freedom: {res.degrees_of_freedom}")

0 commit comments

Comments
 (0)