[CHNOSZ-commits] r236 - in pkg/CHNOSZ: . R demo inst man

noreply at r-forge.r-project.org noreply at r-forge.r-project.org
Mon Oct 2 13:35:04 CEST 2017


Author: jedick
Date: 2017-10-02 13:35:03 +0200 (Mon, 02 Oct 2017)
New Revision: 236

Added:
   pkg/CHNOSZ/demo/lambda.R
Modified:
   pkg/CHNOSZ/DESCRIPTION
   pkg/CHNOSZ/R/berman.R
   pkg/CHNOSZ/R/examples.R
   pkg/CHNOSZ/R/util.expression.R
   pkg/CHNOSZ/demo/00Index
   pkg/CHNOSZ/demo/DEW.R
   pkg/CHNOSZ/inst/NEWS
   pkg/CHNOSZ/man/berman.Rd
   pkg/CHNOSZ/man/examples.Rd
Log:
add transition calculations to berman() and add demo labmda.R


Modified: pkg/CHNOSZ/DESCRIPTION
===================================================================
--- pkg/CHNOSZ/DESCRIPTION	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/DESCRIPTION	2017-10-02 11:35:03 UTC (rev 236)
@@ -1,6 +1,6 @@
 Date: 2017-10-02
 Package: CHNOSZ
-Version: 1.1.0-34
+Version: 1.1.0-35
 Title: Thermodynamic Calculations for Geobiochemistry
 Author: Jeffrey Dick
 Maintainer: Jeffrey Dick <j3ffdick at gmail.com>

Modified: pkg/CHNOSZ/R/berman.R
===================================================================
--- pkg/CHNOSZ/R/berman.R	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/R/berman.R	2017-10-02 11:35:03 UTC (rev 236)
@@ -1,10 +1,16 @@
 # CHNOSZ/berman.R 20170930
-# calculate thermodynamic properties of minerals using Berman formulation
+# calculate thermodynamic properties of minerals using equations from:
+#   Berman, R. G. (1988) Internally-consistent thermodynamic data for minerals
+#      in the system Na2O-K2O-CaO-MgO-FeO-Fe2O3-Al2O3-SiO2-TiO2-H2O-CO2.
+#      J. Petrol. 29, 445-522. https://doi.org/10.1093/petrology/29.2.445
 
-berman <- function(name, T = 298.15, P = 1, thisinfo=NULL, check.G=FALSE) {
+berman <- function(name, T = 298.15, P = 1, thisinfo=NULL, check.G=FALSE, calc.transition=TRUE, calc.disorder=TRUE, units="cal") {
   # reference temperature and pressure
   Pr <- 1
   Tr <- 298.15
+  # the number of conditions we have
+  ncond <- max(length(T), length(P))
+  # get thermodynamic parameters
   file <- system.file("extdata/Berman/Ber88.csv", package="CHNOSZ")
   dat <- read.csv(file, as.is=TRUE)
   # remove the multipliers
@@ -17,10 +23,10 @@
   dat[, 2:27] <- t(t(dat[, 2:27]) / 10^multexp)
   # which row has data for this mineral?
   irow <- which(dat$name == name)
-  # only the immediately following assign() call is needed for the function to work,
+  # the function works fine with just the following assign() call,
   # but an explicit dummy assignment here is used to avoid "Undefined global functions or variables" in R CMD check
-  GfPrTr <- HfPrTr <- SPrTr <- Tmax <- Tmin <- VPrTr <-
-    d0 <- d1 <- d2 <- d3 <- d4 <- d5 <- k0 <- k1 <- k2 <- k3 <- v1 <- v2 <- v3 <- v4 <- NA
+  GfPrTr <- HfPrTr <- SPrTr <- Tlambda <- Tmax <- Tmin <- Tref <- VPrTr <-
+    d0 <- d1 <- d2 <- d3 <- d4 <- d5 <- dTdP <- k0 <- k1 <- k2 <- k3 <- l1 <- l2 <- v1 <- v2 <- v3 <- v4 <- NA
   # assign values to the variables used below
   for(i in 1:ncol(dat)) assign(colnames(dat)[i], dat[irow, i])
   # get the entropy of the elements using the chemical formula in thermo$obigt
@@ -58,8 +64,46 @@
   intdVdT <- -VPrTr*(v3 + v4*(-2*Tr + 2*T)) + P*VPrTr*(v3 + v4*(-2*Tr + 2*T))
   S <- SPrTr + intCpoverT - intdVdT
 
+  ### polymorphic transition properties ***
+  if(!is.na(Tlambda) & !is.na(Tref) & any(T > Tref) & calc.transition) {
+    # starting transition contributions are 0
+    Cptr <- Htr <- Str <- Gtr <- numeric(ncond)
+    ## Ber88 Eq. 8: Cp at 1 bar
+    #Cplambda_1bar <- T * (l1 + l2 * T)^2
+    # Eq. 9: Tlambda at P
+    Tlambda_P <- Tlambda + dTdP * (P - 1)
+    # Eq. 8a: Cp at P
+    Td <- Tlambda - Tlambda_P
+    Tprime <- T + Td
+    # with the condition that Tref < Tprime < Tlambda(1bar)
+    iTprime <- Tref < Tprime & Tprime < Tlambda
+    Tprime <- Tprime[iTprime]
+    Cptr[iTprime] <- Tprime * (l1 + l2 * Tprime)^2
+    # we got Cp, now calculate the integrations for H and S
+    # the lower integration limit is Tref
+    iTtr <- T > Tref
+    Ttr <- T[iTtr]
+    # the upper integration limit is Tlambda_P
+    Ttr[Ttr > Tlambda_P] <- Tlambda_P
+    # derived variables
+    tref <- Tref - Td
+    x1 <- l1^2 * Td + 2 * l1 * l2 * Td^2 + l2^2 * Td^3
+    x2 <- l1^2 + 4 * l1 * l2 * Td + 3 * l2^2 * Td^2
+    x3 <- 2 * l1 * l2 + 3 * l2^2 * Td
+    x4 <- l2 ^ 2
+    # Eqs. 10, 11, 12
+    Htr[iTtr] <- x1 * (Ttr - tref) + x2 / 2 * (Ttr^2 - tref^2) + x3 / 3 * (Ttr^3 - tref^3) + x4 / 4 * (Ttr^4 - tref^4)
+    Str[iTtr] <- x1 * (log(Ttr) - log(tref)) + x2 * (Ttr - tref) + x3 / 2 * (Ttr^2 - tref^2) + x4 / 3 * (Ttr^3 - tref^3)
+    Gtr <- Htr - T * Str
+    # apply the transition contributions
+    Ga <- Ga + Gtr
+    Ha <- Ha + Htr
+    S <- S + Str
+    Cp <- Cp + Cptr
+  }
+
   ### disorder thermodynamic properties ###
-  if(!is.na(Tmin) & !is.na(Tmax) & any(T > Tmin)) {
+  if(!is.na(Tmin) & !is.na(Tmax) & any(T > Tmin) & calc.disorder) {
     # starting disorder contributions are 0
     Cpds <- Hds <- Sds <- Vds <- Gds <- 0
     # the lower integration limit is Tmin
@@ -95,11 +139,16 @@
   ### thermodynamic and unit conventions used in SUPCRT ###
   # use entropy of the elements in calculation of G --> Benson-Helgeson convention (DG = DH - T*DS)
   Gf <- Ga + Tr * SPrTr_elements
+  # the output will just have "G" and "H"
+  G <- Gf
+  H <- Ha
   # convert J to cal
-  G <- convert(Gf, "cal")
-  H <- convert(Ha, "cal")
-  S <- convert(S, "cal")
-  Cp <- convert(Cp, "cal")
+  if(grepl("cal", units)) {
+    G <- convert(Gf, "cal")
+    H <- convert(Ha, "cal")
+    S <- convert(S, "cal")
+    Cp <- convert(Cp, "cal")
+  }
   # convert J/bar to cm^3/mol
   V <- V * 10
 

Modified: pkg/CHNOSZ/R/examples.R
===================================================================
--- pkg/CHNOSZ/R/examples.R	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/R/examples.R	2017-10-02 11:35:03 UTC (rev 236)
@@ -29,7 +29,7 @@
 demos <- function(which=c("sources", "protein.equil", "affinity", "NaCl", "density", 
   "ORP", "revisit", "findit", "ionize", "buffer", "protbuff", "yeastgfp", "mosaic",
   "copper", "solubility", "wjd", "dehydration", "bugstab", "Shh", "activity_ratios",
-  "adenine", "DEW"), to.file=FALSE) {
+  "adenine", "DEW", "lambda"), to.file=FALSE) {
   # run one or more demos from CHNOSZ with ask=FALSE, and return the value of the last one
   for(i in 1:length(which)) {
     # say something so the user sees where we are

Modified: pkg/CHNOSZ/R/util.expression.R
===================================================================
--- pkg/CHNOSZ/R/util.expression.R	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/R/util.expression.R	2017-10-02 11:35:03 UTC (rev 236)
@@ -96,9 +96,11 @@
     # D for greek Delta
     # p for subscript italic P (in Cp)
     # 0 for degree sign (but not immediately following a number e.g. 2.303)
+    # l for subscript small lambda
     if(thischar=='D') thisexpr <- substitute(Delta)
     if(thischar=='p') thisexpr <- substitute(a[italic(P)], list(a=""))
     if(thischar=='0' & !can.be.numeric(prevchar)) thisexpr <- substitute(degree)
+    if(thischar=='l') thisexpr <- substitute(a[lambda], list(a=""))
     # put it together
     expr <- substitute(a*b, list(a=expr, b=thisexpr))
   }

Modified: pkg/CHNOSZ/demo/00Index
===================================================================
--- pkg/CHNOSZ/demo/00Index	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/demo/00Index	2017-10-02 11:35:03 UTC (rev 236)
@@ -20,3 +20,4 @@
 activity_ratios mineral stability plots with activity ratios on the axes
 adenine         HKF parameters regressed from heat capacity and volume of aqueous adenine
 DEW             Deep Earth Water (DEW) model for high pressures
+lambda          Effects of lambda transtion on thermodynamic properties of quartz

Modified: pkg/CHNOSZ/demo/DEW.R
===================================================================
--- pkg/CHNOSZ/demo/DEW.R	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/demo/DEW.R	2017-10-02 11:35:03 UTC (rev 236)
@@ -24,9 +24,9 @@
 PT20.0 <- data.frame(P=20000, T=seq(200, 800, 10))
 PT <- rbind(PT0.5, PT1.0, PT2.0, PT5.0, PT10.0, PT20.0)
 # reaction 1: quartz = SiO2(aq) [equivalent to quartz + 3 H2O = Si(OH)4]
-SiO2_logK <- subcrt(c("quartz", "SiO2"), c(-1, 1), P=PT$P, T=PT$T)$out$logK
+SiO2_logK <- subcrt(c("quartz", "SiO2"), c("cr_Berman", "aq"), c(-1, 1), P=PT$P, T=PT$T)$out$logK
 # reaction 2: 2 quartz = Si2O4(aq) [equivalent to 2 quartz + 3 H2O = Si2O(OH)6]
-Si2O4_logK <- subcrt(c("quartz", "Si2O4"), c(-2, 1), P=PT$P, T=PT$T)$out$logK
+Si2O4_logK <- subcrt(c("quartz", "Si2O4"), c("cr_Berman", "aq"), c(-2, 1), P=PT$P, T=PT$T)$out$logK
 # plot the sum of molalities (== activities) for each pressure
 plot(c(200, 1000), c(-2.5, 0.5), type="n", xlab=axis.label("T"), ylab="log molality")
 for(P in unique(PT$P)) {

Added: pkg/CHNOSZ/demo/lambda.R
===================================================================
--- pkg/CHNOSZ/demo/lambda.R	                        (rev 0)
+++ pkg/CHNOSZ/demo/lambda.R	2017-10-02 11:35:03 UTC (rev 236)
@@ -0,0 +1,61 @@
+# plot effects of lambda transition in quartz
+# after Berman 1988 Figs. 1 and 2
+layout(matrix(c(1, 4:2, 1, 7:5), nrow=4), heights=c(0.7, 3, 3, 3))
+# plot title first
+par(mar=c(0, 0, 0, 0))
+plot.new()
+text(0.5, 0.5, "Effects of lambda transition in quartz, after Berman (1988) Figs. 1 and 2", cex=1.8)
+opar <- par(mar=c(4, 4.5, 1, 0.5), cex=0.8)
+
+T <- convert(seq(0, 1400, 1), "K")
+
+labplot <- function(x) label.plot(x, xfrac=0.9, yfrac=0.1, paren=TRUE)
+# this sets the units used for making the axis labels
+E.units("J")
+Cplab <- axis.label("Cp")
+Vlab <- axis.label("V")
+Tlab <- axis.label("T")
+
+# calculate properties at 1 kbar with and without transition
+Qz_1bar <- berman("quartz", T=T, units="J")
+Qz_1bar_notrans <- berman("quartz", T=T, calc.transition=FALSE, units="J")
+# Fig. 1a: volume
+plot(T, Qz_1bar$V, type="l", xlab=Tlab, ylab=Vlab, ylim=c(22.5, 24))
+legend("topleft", legend="1 bar", bty="n")
+labplot("a")
+# TODO: why don't we get the curvature his plot for V shows?
+# Should it be in the v4 parameter (but it's zero)??
+
+# Fig. 1b: heat capacity
+plot(T, Qz_1bar$Cp, type="l", xlab=Tlab, ylab=Cplab)
+lines(T, Qz_1bar_notrans$Cp, lty=3)
+legend("topleft", legend="1 bar", bty="n")
+labplot("b")
+
+# calculate properties at 10 kbar with and without transition
+Qz_10bar <- berman("quartz", T=T, P=10000, units="J")
+Qz_10bar_notrans <- berman("quartz", T=T, P=10000, calc.transition=FALSE, units="J")
+# Fig. 1c: heat capacity
+plot(T, Qz_10bar$Cp, type="l", xlab=Tlab, ylab=Cplab)
+lines(T, Qz_10bar_notrans$Cp, lty=3)
+legend("topleft", legend="10 kb", bty="n")
+labplot("c")
+
+# like Ber88 Fig. 2
+Tlambda <- 848 # Kelvin
+dTdP <- 0.0237
+Pkb <- seq(1, 50, 1)
+P <- 1000 * Pkb
+T <- Tlambda + dTdP * (P - 1)
+Qz_withtrans <- berman("quartz", T=T, P=P, units="J")
+Qz_notrans <- berman("quartz", T=T, P=P, calc.transition=FALSE, units="J")
+Qz_lambda <- Qz_withtrans - Qz_notrans
+Plab <- expression(list(italic(P), "kb"))
+plot(Pkb, Qz_lambda$G, type="l", ylim=c(-300, -50), ylab=axis.label("DlG"), xlab=Plab)
+labplot("d")
+plot(Pkb, Qz_lambda$H, type="l", ylim=c(1200, 1800), ylab=axis.label("DlH"), xlab=Plab)
+labplot("e")
+plot(Pkb, Qz_lambda$S, type="l", ylim=c(0, 3), ylab=axis.label("DlS"), xlab=Plab)
+labplot("f")
+
+par(opar)

Modified: pkg/CHNOSZ/inst/NEWS
===================================================================
--- pkg/CHNOSZ/inst/NEWS	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/inst/NEWS	2017-10-02 11:35:03 UTC (rev 236)
@@ -1,4 +1,4 @@
-CHANGES IN CHNOSZ 1.1.0-33 (2017-10-02)
+CHANGES IN CHNOSZ 1.1.0-35 (2017-10-02)
 ---------------------------------------
 
 MAJOR CHANGES:
@@ -14,8 +14,11 @@
 - Usage of the DEW model is demonstrated in the new demo DEW.R.
 
 - Add berman() function and extdata/Berman/*.csv files for calculating
-  the thermodynamic properties using equations of Berman, 1988.
+  thermodynamic properties of mineral using equations of Berman, 1988.
 
+- Calculations related to Berman's (1988) Figs. 1 and 2 for the lambda
+  transition of quartz are presented in the new demo lambda.R.
+
 - Implement SUPCRT92's handling of variable volume for quartz and
   coesite. Calculations for other minerals still assume constant
   volume of each phase.

Modified: pkg/CHNOSZ/man/berman.Rd
===================================================================
--- pkg/CHNOSZ/man/berman.Rd	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/man/berman.Rd	2017-10-02 11:35:03 UTC (rev 236)
@@ -7,7 +7,8 @@
 }
 
 \usage{
-  berman(name, T = 298.15, P = 1, thisinfo = NULL, check.G = FALSE)
+  berman(name, T = 298.15, P = 1, thisinfo = NULL, check.G = FALSE,
+         calc.transition = TRUE, calc.disorder = TRUE, units = "cal")
 }
 
 \arguments{
@@ -16,13 +17,16 @@
   \item{P}{numeric, pressure(s) at which to calculate properties (bar)}
   \item{thisinfo}{dataframe, row for mineral from \code{thermo$obigt}}
   \item{check.G}{logical, check consistency of G, H, and S?}
+  \item{calc.transition}{logical, include calculation of polymorphic transition properties?}
+  \item{calc.disorder}{logical, include calculation of disordering properties?}
+  \item{units}{character, energy units, \samp{cal} or \samp{J}}
 }
 
 \details{
 This function calculates the thermodynamic properties of minerals at high \P and \T using equations given by Berman (1988).
 The \code{name} refers to a mineral that must be listed in \code{thermo$obigt} with the state \samp{cr_Berman}.
 This file also holds the chemical formula, which is required for calculating the entropies of the elements in the mineral.
-These entropies are used to convert the apparent Gibbs energies from the Berman-Brown convention to the the Benson-Helgeson convention.
+These entropies are used to convert the apparent Gibbs energies from the Berman-Brown convention to the the Benson-Helgeson convention (cf. Anderson, 2005).
 
 Becuase they use a different set of parameters than Helgeson et al., 1978 (see \code{\link{cgl}}), the standard state thermodynamic properties and parameters for the calculations are stored in files under \code{extdata/Berman}.
 
@@ -63,6 +67,8 @@
 }
 
 \references{
+Anderson, G. M. (2005) \emph{Thermodynamics of Natural Systems}, 2nd ed., Cambridge University Press, 648 p. \url{http://www.worldcat.org/oclc/474880901}
+
 Berman, R. G. (1988) Internally-consistent thermodynamic data for minerals in the system Na{\s2}O–K{\s2}O–CaO–MgO–FeO–Fe{\s2}O{\s3}–Al{\s2}O{\s3}–SiO{\s2}–TiO{\s2}–H{\s2}O–CO{\s2}. \emph{J. Petrol.} \bold{29}, 445-522. \url{https://doi.org/10.1093/petrology/29.2.445}
 
 Helgeson, H. C., Delany, J. M., Nesbitt, H. W. and Bird, D. K. (1978) Summary and critique of the thermodynamic properties of rock-forming minerals. \emph{Am. J. Sci.} \bold{278-A}, 1--229. \url{http://www.worldcat.org/oclc/13594862}

Modified: pkg/CHNOSZ/man/examples.Rd
===================================================================
--- pkg/CHNOSZ/man/examples.Rd	2017-10-02 06:22:50 UTC (rev 235)
+++ pkg/CHNOSZ/man/examples.Rd	2017-10-02 11:35:03 UTC (rev 236)
@@ -18,7 +18,7 @@
     "density", "ORP", "revisit", "findit", "ionize", "buffer",
     "protbuff", "yeastgfp", "mosaic", "copper", "solubility",
     "wjd", "dehydration", "bugstab", "Shh", "activity_ratios",
-    "adenine", "DEW"),
+    "adenine", "DEW", "lambda"),
     to.file=FALSE)
 }
 
@@ -56,34 +56,37 @@
     \code{activity_ratios} \tab mineral stability plots with activity ratios on the axes \cr
     \code{adenine} \tab HKF parameters regressed from heat capacity and volume of aqueous adenine (Lowe et al., 2017) \cr
     \code{adenine} \tab Deep Earth Water (DEW) model for high pressures (Sverjensky et al., 2014a and 2014b) \cr
+    \code{lambda} \tab Effects of lambda transtion on thermodynamic properties of quartz (Berman, 1988) \cr
   }
 
 }
 
 \references{
-  Aksu, S. and Doyle, F. M. (2001) Electrochemistry of copper in aqueous glycine solutions. \emph{J. Electrochem. Soc.} \bold{148}, B51--B57. \url{https://doi.org/10.1149/1.1344532}
+Aksu, S. and Doyle, F. M. (2001) Electrochemistry of copper in aqueous glycine solutions. \emph{J. Electrochem. Soc.} \bold{148}, B51--B57. \url{https://doi.org/10.1149/1.1344532}
 
-  Dick, J. M. (2009) Calculation of the relative metastabilities of proteins in subcellular compartments of \emph{Saccharomyces cerevisiae}. \emph{BMC Syst. Biol.} \bold{3}:75. \url{https://doi.org/10.1186/1752-0509-3-75}
+Berman, R. G. (1988) Internally-consistent thermodynamic data for minerals in the system Na{\s2}O–K{\s2}O–CaO–MgO–FeO–Fe{\s2}O{\s3}–Al{\s2}O{\s3}–SiO{\s2}–TiO{\s2}–H{\s2}O–CO{\s2}. \emph{J. Petrol.} \bold{29}, 445-522. \url{https://doi.org/10.1093/petrology/29.2.445}
 
-  Dick, J. M. and Shock, E. L. (2011) Calculation of the relative chemical stabilities of proteins as a function of temperature and redox chemistry in a hot spring. \emph{PLoS ONE} \bold{6}, e22782. \url{https://doi.org/10.1371/journal.pone.0022782}
+Dick, J. M. (2009) Calculation of the relative metastabilities of proteins in subcellular compartments of \emph{Saccharomyces cerevisiae}. \emph{BMC Syst. Biol.} \bold{3}:75. \url{https://doi.org/10.1186/1752-0509-3-75}
 
-  Dick, J. M. (2016) Proteomic indicators of oxidation and hydration state in colorectal cancer. \emph{PeerJ} \bold{4}:e2238. \url{https://doi.org/10.7717/peerj.2238}
+Dick, J. M. and Shock, E. L. (2011) Calculation of the relative chemical stabilities of proteins as a function of temperature and redox chemistry in a hot spring. \emph{PLoS ONE} \bold{6}, e22782. \url{https://doi.org/10.1371/journal.pone.0022782}
 
-  Garrels, R. M. and Christ, C. L. (1965) \emph{Solutions, Minerals, and Equilibria}, Harper & Row, New York, 450 p. \url{http://www.worldcat.org/oclc/517586}
+Dick, J. M. (2016) Proteomic indicators of oxidation and hydration state in colorectal cancer. \emph{PeerJ} \bold{4}:e2238. \url{https://doi.org/10.7717/peerj.2238}
 
-  Lowe, A. R., Cox, J. S. and Tremaine, P. R. (2017) Thermodynamics of aqueous adenine: Standard partial molar volumes and heat capacities of adenine, adeninium chloride, and sodium adeninate from \emph{T} = 278.15 K to 393.15 K. \emph{J. Chem. Thermodyn.} \bold{112}, 129--145. \url{https://doi.org/10.1016/j.jct.2017.04.005}
+Garrels, R. M. and Christ, C. L. (1965) \emph{Solutions, Minerals, and Equilibria}, Harper & Row, New York, 450 p. \url{http://www.worldcat.org/oclc/517586}
+
+Lowe, A. R., Cox, J. S. and Tremaine, P. R. (2017) Thermodynamics of aqueous adenine: Standard partial molar volumes and heat capacities of adenine, adeninium chloride, and sodium adeninate from \emph{T} = 278.15 K to 393.15 K. \emph{J. Chem. Thermodyn.} \bold{112}, 129--145. \url{https://doi.org/10.1016/j.jct.2017.04.005}
   
-  Manning, C. E., Shock, E. L. and Sverjensky, D. A. (2013) The chemistry of carbon in aqueous fluids at crustal and upper-mantle conditions: Experimental and theoretical constraints. \emph{Rev. Mineral. Geochem.} \bold{75}, 109--148. \url{https://doi.org/10.2138/rmg.2013.75.5}
+Manning, C. E., Shock, E. L. and Sverjensky, D. A. (2013) The chemistry of carbon in aqueous fluids at crustal and upper-mantle conditions: Experimental and theoretical constraints. \emph{Rev. Mineral. Geochem.} \bold{75}, 109--148. \url{https://doi.org/10.2138/rmg.2013.75.5}
 
-  Schulte, M. D. and Shock, E. L. (1995) Thermodynamics of Strecker synthesis in hydrothermal systems. \emph{Orig. Life Evol. Biosph.} \bold{25}, 161--173. \url{https://doi.org/10.1007/BF01581580}
+Schulte, M. D. and Shock, E. L. (1995) Thermodynamics of Strecker synthesis in hydrothermal systems. \emph{Orig. Life Evol. Biosph.} \bold{25}, 161--173. \url{https://doi.org/10.1007/BF01581580}
 
-  Shock, E. L., Oelkers, E. H., Johnson, J. W., Sverjensky, D. A. and Helgeson, H. C. (1992) Calculation of the thermodynamic properties of aqueous species at high pressures and temperatures: Effective electrostatic radii, dissociation constants and standard partial molal properties to 1000 \eqn{^{\circ}}{°}C and 5 kbar. \emph{J. Chem. Soc. Faraday Trans.} \bold{88}, 803--826. \url{https://doi.org/10.1039/FT9928800803}
+Shock, E. L., Oelkers, E. H., Johnson, J. W., Sverjensky, D. A. and Helgeson, H. C. (1992) Calculation of the thermodynamic properties of aqueous species at high pressures and temperatures: Effective electrostatic radii, dissociation constants and standard partial molal properties to 1000 \eqn{^{\circ}}{°}C and 5 kbar. \emph{J. Chem. Soc. Faraday Trans.} \bold{88}, 803--826. \url{https://doi.org/10.1039/FT9928800803}
 
-  Stumm, W. and Morgan, J. J. (1996) \emph{Aquatic Chemistry: Chemical Equilibria and Rates in Natural Waters}, John Wiley & Sons, New York, 1040 p. \url{http://www.worldcat.org/oclc/31754493}
+Stumm, W. and Morgan, J. J. (1996) \emph{Aquatic Chemistry: Chemical Equilibria and Rates in Natural Waters}, John Wiley & Sons, New York, 1040 p. \url{http://www.worldcat.org/oclc/31754493}
 
-  Sverjensky, D. A., Harrison, B. and Azzolini, D. (2014a) Water in the deep Earth: The dielectric constant and the solubilities of quartz and corundum to 60 kb and 1,200 \degC. \emph{Geochim. Cosmochim. Acta} \bold{129}, 125--145. \url{https://doi.org/10.1016/j.gca.2013.12.019}
+Sverjensky, D. A., Harrison, B. and Azzolini, D. (2014a) Water in the deep Earth: The dielectric constant and the solubilities of quartz and corundum to 60 kb and 1,200 \degC. \emph{Geochim. Cosmochim. Acta} \bold{129}, 125--145. \url{https://doi.org/10.1016/j.gca.2013.12.019}
 
-  Sverjensky, D. A., Stagno, V. and Huang, F. (2014b) Important role for organic carbon in subduction-zone fluids in the deep carbon cycle. \emph{Nat. Geosci.} \bold{7}, 909--913. \url{https://doi.org/10.1038/ngeo2291}
+Sverjensky, D. A., Stagno, V. and Huang, F. (2014b) Important role for organic carbon in subduction-zone fluids in the deep carbon cycle. \emph{Nat. Geosci.} \bold{7}, 909--913. \url{https://doi.org/10.1038/ngeo2291}
 }
 
 



More information about the CHNOSZ-commits mailing list