CHNOSZ FAQ

library(CHNOSZ)
options(width = 80)
# Use pngquant to optimize PNG images
library(knitr)
knit_hooks$set(pngquant = hook_pngquant)
pngquant <- "--speed=1 --quality=0-25"
if (!nzchar(Sys.which("pngquant"))) pngquant <- NULL
# To make warnings appear in text box 20230619
# https://selbydavid.com/2017/06/18/rmarkdown-alerts/
knitr::knit_hooks$set(
   error = function(x, options) {
     paste('\n\n<div class="alert alert-danger">',
           gsub('##', '\n', gsub('^##\ Error', '**Error:**', x)),
           '</div>', sep = '\n')
   },
   warning = function(x, options) {
     paste('\n\n<div class="alert alert-warning">',
           gsub('##', '\n', gsub('^##\ Warning:', '**Warning:**', x)),
           '</div>', sep = '\n')
   },
   message = function(x, options) {
     paste('\n\n<div class="alert alert-info">',
           gsub('##', '\n', x),
           '</div>', sep = '\n')
   }
)

# Set dpi 20231129
knitr::opts_chunk$set(
  dpi = if(nzchar(Sys.getenv("CHNOSZ_BUILD_LARGE_VIGNETTES"))) 100 else 72
)
# Merge consecutive messages into a single div 20231114
knitr::knit_hooks$set(document = function(x){

  # Not sure why this is needed, but simply computing on 'x' doesn't work
  file <- tempfile()
  writeLines(x, file)
  x <- readLines(file)

  # Line numbers of the document with </div>
  enddiv <- which(x == "</div>")
  # Line numbers with <div class="alert alert-info">
  beginalert <- which(x == '<div class="alert alert-info">')
  # Find </div> followed <div class="alert alert-info"> (skip empty lines)
  removediv <- (enddiv + 2) %in% beginalert
  if(any(removediv)) {
    # Lines to remove
    rmlines1 <- enddiv[removediv]
    rmlines2 <- enddiv[removediv] + 1
    rmlines3 <- enddiv[removediv] + 2
    # Do the removing
    x <- x[-c(rmlines1, rmlines2, rmlines3)]
  }

  x

})
NOTE <- '<span style="background-color: yellow;">NOTE</span>'
# CHNOSZ functions
equilibrate_ <- '<code>equilibrate()</code>'
info_ <- '<code>info()</code>'
thermo.refs_ <- '<code>thermo.refs()</code>'
# Math stuff
logK <- "log&thinsp;<i>K</i>"
Hplus <- "H<sup>+</sup>"
HCO2_ <- "HCO<sub>2</sub><sup>−</sup>"
HCO3_ <- "HCO<sub>3</sub><sup>−</sup>"
O2 <- "O<sub>2</sub>"
S2 <- "S<sub>2</sub>"
log <- "log&thinsp;"
aHCO2_ <- "<i>a</i><sub>HCO<sub>2</sub><sup>−</sup></sub>"
aHCO3_ <- "<i>a</i><sub>HCO<sub>3</sub><sup>−</sup></sub>"
logfO2 <- "log&thinsp;<i>f</i><sub>O<sub>2</sub></sub>"
Ctot <- "C<sub>tot</sub>"
C3H5O2_ <- "C<sub>3</sub>H<sub>5</sub>O<sub>2</sub><sup>−</sup>"
a3HCO3_ <- "<i>a</i><sup>3</sup><sub>HCO<sub>3</sub><sup>−</sup></sub>"
aC3H5O2_ <- "<i>a</i><sub>C<sub>3</sub>H<sub>5</sub>O<sub>2</sub><sup>−</sup></sub>"
a2HCO3_ <- "<i>a</i><sup>2</sup><sub>HCO<sub>3</sub><sup>−</sup></sub>"
logCtot <- "log&thinsp;C<sub>tot</sub>"
CO2 <- "CO<sub>2</sub>"
H2O <- "H<sub>2</sub>O"
S3minus <- "S<sub>3</sub><sup>-</sup>"
H2S <- "H<sub>2</sub>S"
SO4__ <- "SO<sub>4</sub><sup>-2</sup>"
Kplus <- "K<sup>+</sup>"
Naplus <- "Na<sup>+</sup>"
Clminus <- "Cl<sup>-</sup>"
H2 <- "H<sub>2</sub>"

This vignette was compiled on r Sys.Date() with CHNOSZ version r sessionInfo()$otherPkgs$CHNOSZ$Version.

How is 'CHNOSZ' pronounced?

As one syllable that starts with an sh sound and rhymes with Oz. CHNOSZ and schnoz are homophones.

Added on 2023-05-22.

How should CHNOSZ be cited?

info(info("alanine"))[c("ref1", "ref2")]
thermo.refs(info("alanine"))
basis(c("pyrite", "pyrrhotite", "oxygen"))
sres <- subcrt("magnetite", 1)
info(sres$reaction$ispecies)[, 1:6]
thermo.refs(sres)
reset()

Added on 2023-05-27; PPM example added on 2023-11-15.

What thermodynamic models are used in CHNOSZ?

Added to https://chnosz.net/ website on 2018-11-13; moved to FAQ on 2023-05-27; added references for revised HKF on 2023-11-17.

When and why do equal-activity boundaries depend on total activity?

Short answer: When the species have the same number of the conserved element (let's take C for example), their activities are raised to the same exponent in reaction quotient, so the activity ratio in the law of mass action becomes unity. But when the species have different numbers of the conserved element (for example, propanoate with 3 C and bicarbonate with 1 C), their activities are raised to different exponents, and the activity ratio does not become unity even when the activities are equal (except for the specific case where the activities themselves are equal to 1). Therefore, in general, the condition of "equal activity" is not sufficient to define boundaries on a relative stability diagram; instead, we need to say "activity of each species equal to x" or alternatively "total activity equal to y".

Long answer: First, consider a reaction between formate and bicarbonate: r HCO2_ + 0.5 r O2 $\rightleftharpoons$ r HCO3_. The law of mass action (LMA) is r logK $=$ r log(r aHCO3_ / r aHCO2_) $-$ 0.5 r logfO2. The condition of equal activity is r aHCO2_ $=$ r aHCO3_. Then, the LMA simplifies to r logK $=$ $-$ 0.5 r logfO2. The total activity of C is given by r Ctot $=$ r aHCO2_ $+$ r aHCO3_. According to the LMA, r logfO2 is a function only of r logK, so dr logfO2/dr logCtot $=$ 0. In other words, the position of the equal-activity boundary is independent of the value of r Ctot.

Next, consider a reaction between propanoate and bicarbonate: r C3H5O2_ + 72 r O2 $\rightleftharpoons$ 3 r HCO3_ + 2 r Hplus. The LMA is r logK $=$ r log(r a3HCO3_ / r aC3H5O2_) $-$ pH $-$ 72 r logfO2. The condition of equal activity is r aC3H5O2_ $=$ r aHCO3_. Then, the LMA simplifies to r logK $=$ r log``r a2HCO3_ $-$ pH $-$ 72 r logfO2. The total activity of C is given by r Ctot $=$ 3 r aC3H5O2_ $+$ r aHCO3_; combined with the condition of equal activity, this gives r Ctot $=$ 4 r aHCO3_. Substituting this into the LMA gives r logK $=$ r log(r Ctot / 4)2 $-$ pH $-$ 72 r logfO2, which can be rearranged to write r logfO2 $=$ 27 (2 r log``r Ctot $-$ r logK $-$ r log16 $-$ pH). It follows that dr logfO2/dr logCtot $=$ 47, and the position of the equal-activity boundary depends on r Ctot.

According to this analysis, increasing r Ctot from 0.03 to 3 molal (a 2 log-unit increase) would have no effect on the location of the formate-bicarbonate equal-activity boundary, but would raise the propanoate-bicarbonate equal-activity boundary by 87 units on the r logfO2 scale. Because the reaction between bicarbonate and r CO2 does not involve r O2 (but rather r H2O and r Hplus), the same effect should occur on the propanoate-r CO2 equal-activity boundary. The plots below, which are made using r equilibrate_ for species in the Deep Earth Water (DEW) model, illustrate this effect.


Added on 2023-05-17.

How can minerals with polymorphic transitions be added to the database?

The different crystal forms of a mineral are called polymorphs. Many minerals undergo polymorphic transitions upon heating. Each polymorph for a given mineral should have its own entry in OBIGT, including values of the standard thermodynamic properties (ΔG°~f~, ΔH°~f~, and S°) at 25 °C. The 25 °C (or 298.15 K) values for high-temperature polymorphs are often not listed in thermodynamic tables, but they can be calculated. This thermodynamic cycle shows how: we calculate the changes of a thermodyamic property (pictured here as DS1 and DS2) between 298.15 K and the transition temperature (Ttr) for two polymorphs, then combine those with the property of the polymorphic transition (DStr) to obtain the difference of the property between the polymorphs at 298.15 K (DS298).

               DStr              DStr: entropy of transition between polymorphs 1 and 2
      Ttr  o---------->o          Ttr: temperature of transition
           ^           |         
           |           |         
       DS1 |           | DS2      DS1: entropy change of polymorph 1 from 298.15 K to Ttr
           |           |          DS2: entropy change of polymorph 2 from 298.15 K to Ttr
           |           v
 298.15 K  o==========>o    DS298: entropy difference between polymorphs 1 and 2 at 298.15 K
               DS298        DS298 = DS1 + DStr - DS2
Polymorph  1           2

As an example, let's add pyrrhotite (Fe0.877S) from @PMW87. The formula and thermodynamic properties of this pyrrhotite differ from those of FeS from @HDNB78, which is already in OBIGT. We begin by defining all the input values in the next code block. In addition to G, H, S, and the heat capacity coefficients, non-NA values of volume (V) must be provided for the polymorph transitions to be calculated correctly by subcrt().

# The formula of the new mineral and literature reference
formula <- "Fe0.877S"
ref1 <- "PMW87"
# Because the properties from Pankratz et al. (1987) are listed in calories,
# we need to change the output of subcrt() to calories (the default is Joules)
E.units("cal")
# Use temperature in Kelvin for the calculations below
T.units("K")
# Thermodynamic properties of polymorph 1 at 25 °C (298.15 K)
G1 <- -25543
H1 <- -25200
S1 <- 14.531
Cp1 <- 11.922
# Heat capacity coefficients for polymorph 1
a1 <- 7.510
b1 <- 0.014798
# For volume, use the value from Helgeson et al. (1978)
V1 <- V2 <- 18.2
# Transition temperature
Ttr <- 598
# Transition enthalpy (cal/mol)
DHtr <- 95
# Heat capacity coefficients for polymorph 2
a2 <- -1.709
b2 <- 0.011746
c2 <- 3073400
# Maximum temperature of polymorph 2
T2 <- 1800

Use the temperature (Ttr) and enthalpy of transition (DHtr) to calculate the entropy of transition (DStr). Note that the Gibbs energy of transition (DGtr) is zero at Ttr.

DGtr <- 0  # DON'T CHANGE THIS
TDStr <- DHtr - DGtr  # TΔS° = ΔH° - ΔG°
DStr <- TDStr / Ttr

Start new database entries that include basic information, volume, and heat capacity coefficients for each polymorph. @PMW87 don't list C~p~° of the high-temperature polymorph extrapolated to 298.15 K, so leave it out. If the properties were in Joules, we would omit E_units = "cal" or change it to E_units = "J".

mod.OBIGT("pyrrhotite_new", formula = formula, state = "cr", ref1 = ref1,
  E_units = "cal", G = 0, H = 0, S = 0, V = V1, Cp = Cp1,
  a = a1, b = b1, c = 0, d = 0, e = 0, f = 0, lambda = 0, T = Ttr)
mod.OBIGT("pyrrhotite_new", formula = formula, state = "cr2", ref1 = ref1,
  E_units = "cal", G = 0, H = 0, S = 0, V = V2,
  a = a2, b = b2, c = c2, d = 0, e = 0, f = 0, lambda = 0, T = T2)

For the time being, we set G, H, and S (i.e., the properties at 25 °C) to zero in order to easily calculate the temperature integrals of the properties from 298.15 K to Ttr. Values of zero are placeholders that don't satisfy ΔG°~f~ = ΔH°~f~ − TΔS°~f~ for either polymorph (the subscript f represents formation from the elements), as indicated by the following messages. We will check again for consistency of the thermodynamic parameters at the end of the example.

info(info("pyrrhotite_new", "cr"))
info(info("pyrrhotite_new", "cr2"))

In order to calculate the temperature integral of ΔG°~f~, we need not only the heat capacity coefficients but also actual values of S°. Therefore, we start by calculating the entropy changes of each polymorph from 298.15 to r Ttr K (DS1 and DS2) and combining those with the entropy of transition to obtain the difference of entropy between the polymorphs at 298.15 K. For polymorph 1 (with state = "cr") it's advisable to include use.polymorphs = FALSE to prevent subcrt() from trying to identify the most stable polymorph at the indicated temperature.

DS1 <- subcrt("pyrrhotite_new", "cr", P = 1, T = Ttr, use.polymorphs = FALSE)$out[[1]]$S
DS2 <- subcrt("pyrrhotite_new", "cr2", P = 1, T = Ttr)$out[[1]]$S
DS298 <- DS1 + DStr - DS2

Put the values of S° at 298.15 into OBIGT, then calculate the changes of all thermodynamic properties of each polymorph between 298.15 K and Ttr.

mod.OBIGT("pyrrhotite_new", state = "cr", S = S1)
mod.OBIGT("pyrrhotite_new", state = "cr2", S = S1 + DS298)
D1 <- subcrt("pyrrhotite_new", "cr", P = 1, T = Ttr, use.polymorphs = FALSE)$out[[1]]
D2 <- subcrt("pyrrhotite_new", "cr2", P = 1, T = Ttr)$out[[1]]

It's a good idea to check that the entropy of transition is calculated correctly.

stopifnot(all.equal(D2$S - D1$S, DStr))

Now we're ready to add up the contributions to ΔG°~f~ and ΔH°~f~ from the left, top, and right sides of the cycle. This gives us the differences between the polymorphs at 298.15 K, which we use to make the final changes to the database.

DG298 <- D1$G + DGtr - D2$G
DH298 <- D1$H + DHtr - D2$H
mod.OBIGT("pyrrhotite_new", state = "cr", G = G1, H = H1)
mod.OBIGT("pyrrhotite_new", state = "cr2", G = G1 + DG298, H = H1 + DH298)

It's a good idea to check that the values of G, H, and S at 25 °C for a given polymorph are consistent with each other. Here we use check.GHS() to calculate the difference between the value given for G and the value calculated from H and S. The difference of less than 1 r E.units()/mol can probably be attributed to small differences in the entropies of the elements used by @PMW87 and in CHNOSZ. We still get a message that the database value of C~p~° at 25 °C for the high-temperature polymorph is NA; this is OK because the (extrapolated) value can be calculated from the heat capacity coefficients.

cr_parameters <- info(info("pyrrhotite_new", "cr"))
stopifnot(abs(check.GHS(cr_parameters)) < 1)
cr2_parameters <- info(info("pyrrhotite_new", "cr2"))
stopifnot(abs(check.GHS(cr2_parameters)) < 1)

For the curious, here are the parameter values:

cr_parameters
cr2_parameters

Finally, let's look at the thermodynamic properties of the newly added pyrrhotite as a function of temperature around Ttr. Here, we use the feature of subcrt() that identifies the stable polymorph at each temperature. Note that ΔG°~f~ is a continuous function -- visual confirmation that the parameters yield zero for DGtr -- but ΔH°~f~, S°, and C~p~° are discontinuous at the transition temperature.


For additional polymorphs, we could repeat the above procedure using polymorph 2 as the starting point to calculate G, H, and S of polymorph 3, and so on.

Added on 2023-06-23.

How can I make a diagram with the trisulfur radical ion (r S3minus)?

A r logfO2--pH plot for aqueous sulfur species including r S3minus was first presented by @PD11. Later, @PD15 reported parameters in the revised HKF equations of state for r S3minus, which are available in OBIGT.

par(mfrow = c(1, 3))

## BLOCK 1
T <- 350
P <- 5000
res <- 200

## BLOCK 2
wt_percent_S <- 10
wt_permil_S <- wt_percent_S * 10
molar_mass_S <- mass("S") # 32.06
moles_S <- wt_permil_S / molar_mass_S
grams_H2O <- 1000 - wt_permil_S
molality_S <- 1000 * moles_S / grams_H2O
logm_S <- log10(molality_S)

## BLOCK 3
basis(c("Ni", "SiO2", "Fe2O3", "H2S", "H2O", "oxygen", "H+"))
species(c("H2S", "HS-", "SO2", "HSO4-", "SO4-2", "S3-"))
a <- affinity(pH = c(2, 10, res), O2 = c(-34, -22, res), T = T, P = P)
e <- equilibrate(a, loga.balance = logm_S)
diagram(e)

## BLOCK 4
mod.buffer("NNO", c("nickel", "bunsenite"), state = "cr", logact = 0)
for(buffer in c("HM", "QFM", "NNO")) {
  basis("O2", buffer)
  logfO2_ <- affinity(return.buffer = TRUE, T = T, P = P)$O2
  abline(h = logfO2_, lty = 3, col = 2)
  text(8.5, logfO2_, buffer, adj = c(0, 0), col = 2, cex = 0.9)
}

## BLOCK 5
pH <- subcrt(c("H2O", "H+", "OH-"), c(-1, 1, 1), T = T, P = P)$out$logK / -2
abline(v = pH, lty = 2, col = 4)

## BLOCK 6
lTP <- describe.property(c("T", "P"), c(T, P))
lS <- paste(wt_percent_S, "wt% S(aq)")
ltext <- c(lTP, lS)
legend("topright", ltext, bty = "n", bg = "transparent")
title(quote("Parameters for"~S[3]^"-"~"from Pokrovski and Dubessy (2015)"), xpd = NA)

######## Plot 2: Modify Gibbs energy

oldG <- info(info("S3-"))$G
newG <- oldG - 12548
mod.OBIGT("S3-", G = newG)
a <- affinity(pH = c(2, 10, res), O2 = c(-34, -22, res), T = T, P = P)
e <- equilibrate(a, loga.balance = logm_S)
diagram(e)
legend("topright", ltext, bty = "n", bg = "transparent")
title(quote("Modified"~log*italic(K)~"after Pokrovski and Dubrovinsky (2011)"), xpd = NA)
OBIGT()

######## Plot 3: Do it with DEW

T <- 700
P <- 10000 # 10000 bar = 1 GPa
oldwat <- water("DEW")
add.OBIGT("DEW")
info(species()$ispecies)
a <- affinity(pH = c(2, 10, res), O2 = c(-18, -10, res), T = T, P = P)
e <- equilibrate(a, loga.balance = logm_S)
diagram(e)
lTP <- describe.property(c("T", "P"), c(T, P))
ltext <- c(lTP, lS)
legend("topright", ltext, bty = "n", bg = "transparent")
title(quote("Deep Earth Water (DEW)"~"model"))
water(oldwat)
OBIGT()

The blocks of code are commented here:

  1. Set temperature, pressure, and resolution.
  2. Calculate molality of S from given weight percent [this is rather tedious and could be condensed to fewer lines of code].
  3. Define the given weight percent (10 wt% S).
  4. Calculate weight permil S.
  5. Divide by molar mass to calculate moles of S in 1 kg of solution.
  6. Calculate grams of r H2O in 1 kg of solution.
  7. Calculate molality (moles of S per kg of r H2O, not kg of solution).
  8. Calculate decimal logarithm of molality.
  9. Define the basis species and formed species, calculate affinity, equilibrate activities, and make the diagram.
  10. If we didn't want to plot the buffer lines, we could just use basis(c("H2S", "H2O", "oxygen", "H+")).
  11. Basis species with Fe, Si, and Ni are needed for the HM, QFM, and NNO buffers.
  12. Note that "oxygen" matches r O2(gas), not r O2(aq), so the variable on the diagram is r logfO2.
  13. Define Ni-NiO (NNO) buffer and plot buffer lines for HM, QFM, and NNO.
  14. QFM (quartz-fayalite-magnetite) is also known as FMQ.
  15. Calculate and plot pH of neutrality for water.
  16. Add a legend and title.

Why does the published diagram have a much larger stability field for r S3minus?

Let's calculate r logK for the reaction 2 r H2S(aq) + r SO4__ + r Hplus = r S3minus + 0.75 r O2(gas) + 2.5 r H2O.

species <- c("H2S", "SO4-2", "H+", "S3-", "oxygen", "H2O")
coeffs <- c(-2, -1, -1, 1, 0.75, 2.5)
(calclogK <- subcrt(species, coeffs, T = seq(300, 450, 50), P = 5000)$out$logK)
fcalclogK <- formatC(calclogK, format = "f", digits = 1)
reflogK <- -9.6
dlogK <- calclogK[2] - reflogK
# Put in a test so that we don't get surprised by
# future database updates or changes to this vignette
stopifnot(round(dlogK, 4) == -4.4132)

By using the thermodynamic parameters for r S3minus in OBIGT that are taken from @PD15, r logK is calculated to be r paste(paste(fcalclogK[1:3], collapse = ", "), fcalclogK[4], sep = ", and ") at 300, 350, 400, and 450 °C and 5000 bar. In contrast, ref. 22 of @PD11 lists r reflogK for r logK at 350 °C; this is r round(-dlogK, 1) log units higher than the calculated value of r fcalclogK[2]. This corresponds to a difference of Gibbs energy of -2.303 * 1.9872 * (350 + 273.15) * r round(-dlogK, 1) = r formatC(-2.303 * 1.9872 * (350 + 273.15) * -dlogK, format = "f", digits = 0) cal/mol.

In the code below, we use the difference of Gibbs energy to temporarily update the OBIGT entry for r S3minus. Then, we make a new diagram that is more similar to that from @PD11. Finally, we reset the OBIGT database so that the temporary parameters don't interfere with later calculations.

Can I make the diagram using the Deep Earth Water (DEW) model?

Yes! Just set a new temperature and pressure and activate the DEW water model and load the DEW aqueous species. You can also use r info_ to see which species are affected by loading the DEW parameters; it turns out that r SO4__ isn't. Then, use similar commands as above to make the diagram. At the end, reset the water model and OBIGT database.

Here are the three plots that we made:


Added on 2023-09-08.

In OBIGT, what is the meaning of T for solids, liquids, and gases?

The value in this column can be one of the following:

  1. The temperature of transition to the next polymorph of a mineral;
  2. For the highest-temperature (or only) polymorph, if T is positive, it is the phase stability limit (i.e., the temperature of melting or decomposition of a solid or vaporization of a liquid);
  3. For the highest-temperature polymorph, if T is negative, the opposite (positive) value is the T limit for validity of the C~p~ equation. (New feature in development version of CHNOSZ)

These cases are handled by subcrt() as follows. The units of T in OBIGT are Kelvin, but subcrt() by default uses °C:

1. For polymorphic transitions, the properties of specific polymorphs are returned:

subcrt("pyrrhotite", T = c(25, 150, 350), property = "G")$out

Note: In both SUPCRT92 and OBIGT, tin, sulfur, and selenium are listed as minerals with one or more polymorphic transitions, but the highest-temperature polymorph actually represents the liquid state. Furthermore, quicksilver is listed as a mineral whose polymorphs actually correspond to the liquid and gaseous states.

2. For a phase stability limit, ΔG° is set to NA above the temperature limit:

subcrt("pyrite", T = seq(200, 1000, 200), P = 1)

This feature is intended to make it harder to obtain potentially unreliable results at temperatures where a mineral (or an organic solid or liquid) is not stable. If you want the extrapolated ΔG° above the listed phase stability limit, then add exceed.Ttr = TRUE to the function call to subcrt().

OBIGT has a non-exhaustive list of temperatures of melting, decomposition, or other phase change, some of which were taken from SUPCRT92 while others were taken from @RH95. These minerals are listed below:

file <- system.file("extdata/OBIGT/inorganic_cr.csv", package = "CHNOSZ")
dat <- read.csv(file)
# Reverse rows so highest-T polymorph for each mineral is listed first
dat <- dat[nrow(dat):1, ]
# Remove low-T polymorphs
dat <- dat[!duplicated(dat$name), ]
# Remove minerals with no T limit for phase stability (+ve) or Cp equation (-ve)
dat <- dat[!is.na(dat$z.T), ]
# Keep minerals with phase stability limit
dat <- dat[dat$z.T > 0, ]
# Get names of minerals and put into original order
rev(dat$name)

OBIGT now uses the decomposition temperature of covellite [780.5 K from @RH95] in contrast to the previous Tmax from SUPCRT92 [1273 K, which is referenced to a \Cp equation described as "estimated" on p. 62 of @Kel60]. Selected organic solids and liquids have melting or vaporization temperatures listed as well. However, no melting temperatures are listed for minerals that use the Berman model.

3. For a C~p~ equation limit, extrapolated values of ΔG° are shown and a warning is produced:

add.OBIGT("SUPCRT92")
subcrt("muscovite", T = 850, P = 4500)
reset()

The warning is similar to that produced by SUPCRT92 ("CAUTION: BEYOND T LIMIT OF CP COEFFS FOR A MINERAL OR GAS") at temperatures above maximum temperature of validity of the Maier-Kelley equation (Tmax). Notably, SUPCRT92 outputs ΔG° and other standard thermodynamic properties at temperatures higher than Tmax despite the warning.

This is a new feature in CHNOSZ version 2.1.0. In previous versions of CHNOSZ, values of ΔG° above the C~p~ equation limit were set to NA without a warning, as with the phase stability limit described above.

4. Finally, if T is NA or 0, then no upper temerature limit is imposed by subcrt().

Added on 2023-11-15.

How can mineral pH buffers be plotted?

Unlike mineral redox buffers, the K-feldspar–muscovite–quartz (KMQ) and muscovite–kaolinite (MC) pH buffers are known as "sliding scale" buffers because they do not determine pH but rather the activity ratio of r Kplus to r Hplus [@HA05]. To add these buffers to a r logfO2–pH diagram in CHNOSZ, choose basis species that include Al+3 (the least mobile element, which the reactions are balanced on), quartz (this is needed for the KMQ buffer), the mobile ions r Kplus and r Hplus, and the remaining elements in r H2O and r O2; oxygen denotes the gas in OBIGT. The formation reactions for these minerals don't involve r O2, but it must be present so that the number of basis species equals the number of elements +1 (i.e. elements plus charge).

basis(c("Al+3", "quartz", "K+", "H+", "H2O", "oxygen"))
species(c("kaolinite", "muscovite", "K-feldspar"))

We could go right ahead and make a r logfO2–pH diagram, but the implied assumption would be that the r Kplus activity is unity, which may not be valid. Instead, we can obtain an independent estimate for r Kplus activity based on 1) the activity ratio of r Naplus to r Kplus for the reaction between albite and K-feldspar and 2) charge balance among r Naplus, r Kplus, and r Clminus for a given activity of the latter [@HC14]. Using the variables defined below, those conditions are expressed as K_AK = m_Na / m_K and m_Na + m_K = m_Cl, which combine to give m_K = m_Cl / (K_AK + 1). The reason for writing the equations with molality instead of activity is that ionic strength (IS) is provided in the arguments to subcrt(), so the function returns a value of r logK adjusted for ionic strength. Furthermore, it is assumed that this is a chloride-dominated solution, so ionic strength is taken to be equal to the molality of r Clminus.

# Define temperature, pressure, and molality of Cl- (==IS)
T <- 150
P <- 500
IS <- m_Cl <- 1
# Calculate equilibrium constant for Ab-Kfs reaction, corrected for ionic strength
logK_AK <- subcrt(c("albite", "K+", "K-feldspar", "Na+"), c(-1, -1, 1, 1),
  T = T, P = P, IS = IS)$out$logK
K_AK <- 10 ^ logK_AK
# Calculate molality of K+
(m_K <- m_Cl / (K_AK + 1))

This calculation gives a molality of r Kplus that is lower than unity and accordingly makes the buffers less acidic [@HC14]. Now we can apply the calculated molality of r Kplus to the basis species and add the buffer lines to the diagram. The IS argument is also used for affinity() so that activities are replaced by molalities (that is, affinity is calculated with standard Gibbs energies adjusted for ionic strength; this has the same effect as calculating activity coefficients).

par(mfrow = c(1, 2))
basis("K+", log10(m_K))
a <- affinity(pH = c(2, 10), O2 = c(-55, -38), T = T, P = P, IS = IS)
diagram(a, srt = 90)
lTP <- as.expression(c(lT(T), lP(P)))
legend("topleft", legend = lTP, bty = "n", inset = c(-0.05, 0), cex = 0.9)
ltxt <- c(quote("Unit molality of Cl"^"-"), "Quartz saturation")
legend("topright", legend = ltxt, bty = "n", cex = 0.9)
title("Mineral data from Berman (1988)\nand Sverjensky et al. (1991) (OBIGT default)",
  font.main = 1, cex.main = 0.9)
add.OBIGT("SUPCRT92")
a <- affinity(pH = c(2, 10), O2 = c(-55, -38), T = T, P = P, IS = IS)
diagram(a, srt = 90)
title("Mineral data from Helgeson et al. (1978)\n(as used in SUPCRT92)",
  font.main = 1, cex.main = 0.9)
OBIGT()

The gray area, which is automatically drawn by diagram(), is below the reducing stability limit of water; that is, this area is where the equilibrium fugacity of r H2 exceeds unity. NOTE: Although the muscovite–kaolinite (MC) buffer was mentioned by @HC14 in the context of "clay-rich but feldspar-free sediments", this example uses the feldspathic Ab–Kfs reaction for calculating r Kplus molality for both the KMQ and MC buffers. A more appropriate reaction to constrain the Na/K ratio with the MC buffer may be that between paragonite and muscovite [e.g., @Yar05].

The diagram in Fig. 4 of @HC14 shows the buffer lines at somewhat higher pH values of ca. 5 and 6. Removing IS from the code moves the lines to lower rather than higher pH (not shown -- try it yourself!), so the calculation of activity coefficients does not explain the differences. One possible reason for these differences is the use of different thermodynamic data for the minerals. The parameters for these minerals in the default OBIGT database come from @Ber88 and @SHD91.

thermo.refs(species()$ispecies)

CHNOSZ doesn't implement the thermodynamic model for minerals from @HP98, which is one of the sources cited by @HC14. If we use the thermodynamic parameters for minerals from @HDNB78 [these do not include the revisions for aluminosilicates described by @SHD91], we get the lines shown in the second plot above, representing a larger stability field for muscovite. This moves the KMQ buffer closer to the value shown by @HC14, but the MC buffer further away, so this still doesn't explain why we get a different result.


Added on 2023-11-28.

References



Try the CHNOSZ package in your browser

Any scripts or data that you put into this service are public.

CHNOSZ documentation built on Feb. 12, 2024, 3 p.m.