Introduction to radiocarbon models
Radiocarbon is a useful tracer of the global carbon cycle. In soils, radiocarbon measurements help to estimate how fast new carbon is stored and respired again by microorganisms. It also helps to determine how stable different soil carbon fractions are. Radiocarbon measurements in respired CO2 from soils are closely related to the mean transit time of carbon. Measurements of radiocarbon in bulk soil are closely related to the mean age of carbon. Therefore, radiocarbon is a useful indicator of the timescales of carbon cycling in soils and it’s often used as a surrogate of the concept of soil organic matter persistence.
Radiocarbon dynamics in a one-pool model
Radiocarbon atoms, 14C, in the Earth system follow the same dynamics as the more abundant 12C atoms. The main difference is that radiocarbon atoms go through radioactive decay with a mean lifetime of 8267 years. Although there are isotopic fractionation effects related to the higher atomic weight of radiocarbon, these are corrected in their reporting.
A one pool model for carbon atoms \(C\) in a soil is generally described as $$ \frac{dC}{dt} = I - k \cdot C, $$ where \(I\) are carbon inputs to soil, and \(k\) is a decomposition rate constant. An equivalent model for radiocarbon can be written as $$ \frac{d^{14}C}{dt} = I_{atm} + k \cdot ^{14}C - \lambda \cdot ^{14}C $$ where \(I_{atm}\) is the inputs of radiocarbon from the atmosphere, and \(\lambda\) is the radioactive decay constant, \(\lambda = 1/8267\).
Radiocarbon is commonly reported as a fraction of 14C to 12C atoms with respect to the same ratio of a standard material, corrected for isotopic fractionation and radioactive decay. This fraction is expressed as \(F\), therefore we can modify the differential equation for radiocarbon as $$ \frac{d(F\cdot C)}{dt} = I \cdot F_{atm} + (k-\lambda)(F \cdot C), $$ where \(F_{atm}\) is the fraction of radiocarbon atoms in the atmosphere.
Implementation in SoilR
The one-pool radiocarbon model is implemented in SoilR in the function
OnepModel14
, which requires as main arguments: a vector of times where we want
to compute the solution of the differential equation, a decomposition rate \(k\),
and intial amount of carbon in the system at the initial time together with its
radiocarbon fraction, a vector of carbon inputs to the soil, and the fraction of
radiocarbon in the atmosphere during the time period of the simulation. Additional
details can be found in the help of the function typping ?OnepModel14
.
In the following example, we will run a one pool model for a system with about 200 gC m^-2^ yr^-1^ of litter inputs and a decomposition rate of \(k = 0.25\) yr^-1^. We will model the dynamics of radiocarbon in this system during the period from 1950 until 2019 using the standard atmospheric radiocarbon curves of Hua et al. (2021), assuming that the soil is located inside the Norther Hemisphere zone 2. In SoilR, radiocarbon is reported in \(\Delta^{14}\)C notation in permil units.
years<-seq(1950, 2019, by=1/12)
LitterInput=200
k1=0.25
initialC<-1000
initialF<- -27 # in Delta14C
Fatm<-Hua2021$NHZone2[,1:2] # data.frame with time and atmospheric Delta14C
With these elements, we can now call the OnepModel14
function, which creates a
model that can be used later to compute the radiocarbon dynamics in the soil.
Model1<-OnepModel14(t=years, k=k1, C0=initialC, F0_Delta14C = initialF, In=LitterInput,
inputFc = Fatm)
We can compute now the radiocarbon content of the soil using the function
getF14
. In a one pool system, the dynamics of the bulk soil and the respired
CO2 are the same, so this function give us both.
C14_1<-getF14(Model1)
We can now plot the results and compare the dynamics of radiocarbon in the soil versus the atmosphere.
plot(Fatm, type="l", col="blue", xlab="Calendar year", ylab="Delta14C in permil",
ylim=c(-30,1000), bty="n")
lines(years, C14_1)
legend("topright", c("Atmosphere", "Soil"), lty=1, col=c("blue",1), bty="n")
We see from this example that radiocarbon in this soil responds slowly to the inputs of radiocarbon from the atmosphere. It takes a few years for this soil to assimilate the peak of the bomb spike, and then it incorporates radiocarbon with a time lag.
Speed of decomposition and radiocarbon dynamics
We will run now two additional simulations in which we vary the decomposition rate to observe how differences in the speed of decay affects the incorporation of radiocarbon in a soil.
k2=4 # very fast rate
k3=1/64 # very slow rate
Model2<-OnepModel14(t=years, k=k2, C0=initialC, F0_Delta14C = initialF, In=LitterInput,
inputFc = Fatm)
Model3<-OnepModel14(t=years, k=k3, C0=initialC, F0_Delta14C = initialF, In=LitterInput,
inputFc = Fatm)
C14_2<-getF14(Model2)
C14_3<-getF14(Model3)
We can now plot the results and compare the three systems
pal=rainbow(4)
plot(Fatm, type="l", col="blue", xlab="Calendar year", ylab="Delta14C in permil",
ylim=c(-30,1000), bty="n")
lines(years, C14_2, col=pal[1])
lines(years, C14_1, col=pal[2])
lines(years, C14_3, col=pal[3])
legend("topright", c("Atmosphere", "Fast rate", "Intermediate rate", "Slow rate"),
lty=1, col=c("blue",pal), bty="n")
We can clearly see that a fast system follows very closely the radiocarbon signature of the atmosphere. In this case we would say that the soil is in isotopic equilibrium with the atmosphere. The slow decomposition rate shows a trend that differs considerably from the atmosphere because the soil slowly incorporates atmospheric radiocarbon with a long time lag. In this case we would say that the soil has a large istopic disequilibrium with the atmosphere.