This document describes how to run a water and energy balance model that uses a rather complex approach for hydraulics and stomatal regulation. This document is meant to teach users to run the simulation model within R. All the details of the model design and formulation are given in a separate document that can be found on the website of the R package (http://vegmod.ctfc.cat/medfateweb).
Any water and energy balance model needs information on climate, vegetation and soils of the forest stand to be simulated. Moreover, since the soil water balance in medfate
differentiates between species, information on species-specific model parameters is also needed. In this subsection we explain the different steps (including optional pathways) to prepare all the data needed to run function spwb()
.
Simulation models in medfate
require information on the physical attributes of soil, namely soil depth, texture, bulk density and rock fragment content. Soil information needs to be entered as a data frame with soil layers in rows and physical attributes in columns. The model accepts different layer definitions (from one to five layers). Soil physical attributes can be initialized to default values, for a given number of layers, using function defaultSoilParams()
:
spar = defaultSoilParams(2)
print(spar)
## widths clay sand om bd rfc
## 1 300 25 25 NA 1.5 20
## 2 700 25 25 NA 1.5 40
where widths
are soil layer widths in mm; clay
and sand
are the percentage of clay and sand, in percent of dry weight, om
stands for organic matter, bd
is bulk density (in g·cm\(^{-3}\)) and rfc
the percentage of rock fragments. Because soil properties vary strongly at fine spatial scales, ideally soil physical attributes should be measured on samples taken at the forest stand to be simulated. For those users lacking such data, soil properties modelled at larger scales are available via soilgrids.org (see function soilgridsParams()
).
The soil input for function spwb()
is actually an object of class soil
that is created using a function with the same name:
examplesoil = soil(spar, VG_PTF = "Toth")
names(examplesoil)
## [1] "SoilDepth" "W" "SWE" "Temp"
## [5] "Ksoil" "Gsoil" "dVec" "sand"
## [9] "clay" "om" "usda_Type" "VG_alpha"
## [13] "VG_n" "VG_theta_res" "VG_theta_sat" "macro"
## [17] "rfc"
In addition to the physical soil description, this object contains soil parameters needed for soil water balance simulations. For example, macro
specifies the macroporosity of each layer; Gsoil
and Ksoil
are parameters needed to model infiltration of water into the soil. With VG_PTF = "Toth"
, we specify that van Genuchten parameters are estimated from texture using the pedotransfer functions of Toth et al (2015). The details of all elements in the soil object can be found in the help page for function soil()
.
The soil object is also used to store the moisture degree of each soil layer. In particular, W
contains the state variable that represents moisture (i.e. the proportion of moisture relative to field capacity), which is normally initialized to 1 for each layer:
examplesoil$W
## [1] 1 1
Analogously, Temp
contains the temperature (in degrees) of soil layers, initialized to missing values:
examplesoil$Temp
## [1] NA NA
It is important to remember that, unlike normal objects in R, any water balance simulation will modify the values of state variables in soil object. That is, the state of the soil at the end of the simulated process (i.e. W
) will be stored (the same for Temp
). Hence, one can use the same object to simulate water balance sequentially and the final state of one simulation is the initial state of the next.
At any time, one can print the status of the soil object using its print
function:
print(examplesoil, model = "SX")
## Soil depth (mm): 1000
##
## Layer 1 [ 0 to 300 mm ]
## clay (%): 25 silt (%): 50 sand (%): 25 organic matter (%): NA [ Silt loam ]
## Rock fragment content (%): 20 Macroporosity (%): 5
## Theta WP (%): 14 Theta FC (%): 30 Theta SAT (%): 49 Theta current (%) 30
## Vol. WP (mm): 34 Vol. FC (mm): 73 Vol. SAT (mm): 118 Vol. current (mm): 73
## Temperature (Celsius): NA
##
## Layer 2 [ 300 to 1000 mm ]
## clay (%): 25 silt (%): 50 sand (%): 25 organic matter (%): NA [ Silt loam ]
## Rock fragment content (%): 40 Macroporosity (%): 5
## Theta WP (%): 14 Theta FC (%): 30 Theta SAT (%): 49 Theta current (%) 30
## Vol. WP (mm): 60 Vol. FC (mm): 127 Vol. SAT (mm): 207 Vol. current (mm): 127
## Temperature (Celsius): NA
##
## Total soil saturated capacity (mm): 325
## Total soil water holding capacity (mm): 200
## Total soil extractable water (mm): 106
## Total soil current Volume (mm): 200
## Water table depth (mm): 1000
##
## Snow pack water equivalent (mm): 0
The modelled moisture content of the soil depends on the water retention curve used to represent the relationship between soil volumetric water content (\(\theta\); %) and soil water potential (\(\Psi\); MPa). By default the Saxton equations are used as water retention curve, but the user may choose to follow Van Genuchten - Mualem equations, which will give different values for the same texture:
print(examplesoil, model="VG")
## Soil depth (mm): 1000
##
## Layer 1 [ 0 to 300 mm ]
## clay (%): 25 silt (%): 50 sand (%): 25 organic matter (%): NA [ Silt loam ]
## Rock fragment content (%): 20 Macroporosity (%): 5
## Theta WP (%): 13 Theta FC (%): 29 Theta SAT (%): 42 Theta current (%) 29
## Vol. WP (mm): 32 Vol. FC (mm): 68 Vol. SAT (mm): 102 Vol. current (mm): 68
## Temperature (Celsius): NA
##
## Layer 2 [ 300 to 1000 mm ]
## clay (%): 25 silt (%): 50 sand (%): 25 organic matter (%): NA [ Silt loam ]
## Rock fragment content (%): 40 Macroporosity (%): 5
## Theta WP (%): 13 Theta FC (%): 30 Theta SAT (%): 42 Theta current (%) 30
## Vol. WP (mm): 54 Vol. FC (mm): 127 Vol. SAT (mm): 178 Vol. current (mm): 127
## Temperature (Celsius): NA
##
## Total soil saturated capacity (mm): 280
## Total soil water holding capacity (mm): 196
## Total soil extractable water (mm): 111
## Total soil current Volume (mm): 196
## Water table depth (mm): 1000
##
## Snow pack water equivalent (mm): 0
While Saxton equations use texture and organic matter as inputs, the Van Genuchten-Mualem equations need other parameters, which are estimated using pedotransfer functions and their names start with VG_
. The following figure illustrates the difference between the two water retention models:
par(mar=c(4,4,1,1))
# Plot Saxton's retention curve
psi = seq(-0.01, -2.0, by=-0.001)
plot(-psi, lapply(as.list(psi), FUN=soil.psi2thetaSX, clay=25, sand=25),
type="l", ylim=c(0,0.5),ylab="Water content (prop. volume)", xlab = "Soil water potential (-MPa)")
#Add Van Genuchten retention curve
lines(-psi, lapply(as.list(psi), FUN=soil.psi2thetaVG,
alpha=examplesoil$VG_alpha[1], n = examplesoil$VG_n[1],
theta_res = examplesoil$VG_theta_res[1],
theta_sat = examplesoil$VG_theta_sat[1]), lty=2)
legend("topright", legend=c("Saxton", "Van Genuchten"), lty=c(1,2), bty="n")
Functions soil.psi2thetaSX()
and soil.psi2thetaVG()
(and their counterparts) can be used to calculate volumetric soil moisture from the water potential using the two models. When simulating soil water balance, the user can choose among the two models (see control parameters below). The soil water balance model described in this vignette uses the van Genuchten-Mualem equations for water retention curves and rhizosphere conductance.
Simulation models in medfate
require a data frame with species parameter values. The package provides a default data set of parameter values for 89 Mediterranean species (rows), resulting from bibliographic search, fit to empirical data or expert-based guesses:
data("SpParamsMED")
These species commonly occur in the Spanish forest inventory of Catalonia, but may not be sufficient for other areas. A large number of parameters (columns) can be found in SpParamsMED
:
names(SpParamsMED)
## [1] "Name" "IFNcodes" "SpIndex"
## [4] "Group" "Order" "Family"
## [7] "GrowthForm" "TreeType" "Hmed"
## [10] "Hmax" "Z50" "Z95"
## [13] "Zmax" "a_ash" "a_bsh"
## [16] "b_bsh" "cr" "a_fbt"
## [19] "b_fbt" "c_fbt" "d_fbt"
## [22] "a_cr" "b_1cr" "b_2cr"
## [25] "b_3cr" "c_1cr" "c_2cr"
## [28] "a_cw" "b_cw" "fHDmin"
## [31] "fHDmax" "SLA" "LeafDensity"
## [34] "r635" "pDead" "maxFMC"
## [37] "minFMC" "LeafPI0" "LeafEPS"
## [40] "LeafAF" "StemPI0" "StemEPS"
## [43] "StemAF" "LeafDuration" "LigninPercent"
## [46] "ParticleDensity" "LeafLitterFuelType" "Flammability"
## [49] "SAV" "HeatContent" "albedo"
## [52] "k" "g" "Sgdd"
## [55] "Psi_Extract" "WUE" "pRootDisc"
## [58] "Gwmin" "Gwmax" "VCleaf_kmax"
## [61] "VCleaf_c" "VCleaf_d" "Kmax_stemxylem"
## [64] "VCstem_c" "VCstem_d" "Kmax_rootxylem"
## [67] "VCroot_c" "VCroot_d" "LeafWidth"
## [70] "Vmax298" "Jmax298" "Al2As"
## [73] "WoodDensity" "WoodC" "RGRmax"
## [76] "Cstoragepmax"
Not all parameters are needed for all models. The user can find parameter definitions in the help page of this data set. However, to fully understand the role of parameters in the model, the user should read the details of model design and formulation (http://vegmod.ctfc.cat/medfateweb).
Models included in medfate
were primarily designed to be ran on forest inventory plots. In this kind of data, the vegetation of a sampled area is described in terms of woody plants (trees and shrubs) along with their size and species identity. Forest plots in medfate
are assumed to be in a format that follows closely the Spanish forest inventory. Each forest plot is represented in an object of class forest
, a list that contains several elements. Among them, the most important items are two data frames, treeData
(for trees) and shrubData
for shrubs:
data(exampleforest)
exampleforest
## $ID
## [1] "1"
##
## $patchsize
## [1] 10000
##
## $treeData
## Species N DBH Height Z50 Z95
## 1 54 168 37.55 800 750 3000
## 2 68 384 14.60 660 750 3000
##
## $shrubData
## Species Cover Height Z50 Z95
## 1 65 3.75 30 50 1000
##
## $herbCover
## [1] 10
##
## $herbHeight
## [1] 20
##
## $seedBank
## Species Abundance
## 5 54 100
## 8 65 100
## 9 68 100
##
## attr(,"class")
## [1] "forest" "list"
Trees are expected to be primarily described in terms of species, diameter (DBH) and height, whereas shrubs are described in terms of species, percent cover and mean height.
Because the forest plot format is rather specific, simulation functions in medfate
allow starting in a more general way using two data frames, one with aboveground information (i.e. the leave area and size of plants) and the other with belowground information (i.e. root distribution). The aboveground data frame does not distinguish between trees and shrubs. It includes, for each plant cohort to be considered in rows, its species identity, height and leaf area index (LAI). While users can build their input data themselves, we use function forest2aboveground()
on the object exampleforest
to show how should the data look like:
above = forest2aboveground(exampleforest, SpParamsMED)
above
## SP N DBH Cover H CR LAI_live LAI_expanded
## T1_54 54 168.000 37.55 NA 800 0.7150421 0.81630007 0.81630007
## T2_68 68 384.000 14.60 NA 660 0.6055507 0.79744714 0.79744714
## S1_65 65 5504.183 NA 3.75 30 0.9740000 0.08911235 0.08911235
## LAI_dead
## T1_54 0
## T2_68 0
## S1_65 0
Note that the call to forest2aboveground()
included species parameters, because species-specific values are needed to calculate leaf area from tree diameters or shrub cover. Columns N
, DBH
and Cover
are required for simulating growth, but not for soil water balance, which only requires columns SP
, H
(in cm), CR
(i.e. the crown ratio), LAI_live
, LAI_expanded
and LAI_dead
. Here plant cohorts are given unique codes that tell us whether they correspond to trees or shrubs, but the user can use other row identifiers as long as they are unique. In practice, the user only needs to worry to calculate the values for LAI_live
. LAI_live
and LAI_expanded
can contain the same LAI values, and LAI_dead
is normally zero. This is so because models update LAI_expanded
and LAI_dead
according to the leaf phenology of species.
Regarding belowground information, a matrix describing for each plant cohort, the proportion of fine roots in each soil layer. As before, we use function forest2belowground()
on the object exampleforest
to show how should the data look like:
below = forest2belowground(exampleforest, examplesoil, SpParamsMED)
below
## 1 2
## T1_54 0.1933638 0.8066362
## T2_68 0.1933638 0.8066362
## S1_65 0.8981075 0.1018925
Function forest2belowground()
internally takes values of Z50
and Z95
and calls function root.ldrDistribution()
to estimate the distribution of fine roots according to the linear dose response model. For example the first row is:
root.ldrDistribution(exampleforest$treeData$Z50[1],
exampleforest$treeData$Z95[1],
examplesoil$dVec)
## [,1] [,2]
## [1,] 0.1933638 0.8066362
An analogous function root.conicDistribution()
can be used to estimate fine root distribution according to a cone. The user is free to build any numeric matrix for root distribution, as long as values in rows sum always one (i.e. we have proportions). Otherwise, functions root.conicDistribution()
and root.ldrDistribution()
can be used to calculate root distribution under specific assumptions.
Soil water simulations require daily weather inputs. The weather variables that are required depend on the complexity of the soil water balance model we are using. The complex simulation model requires precipitation, radiation, wind speed, min/max temparature and relative humitidy. Here we show an example of meteorological forcing data.
data(examplemeteo)
head(examplemeteo)
## MeanTemperature MinTemperature MaxTemperature Precipitation
## 2001-01-01 3.57668969 -0.5934215 6.287950 4.869109
## 2001-01-02 1.83695972 -2.3662458 4.569737 2.498292
## 2001-01-03 0.09462563 -3.8541036 2.661951 0.000000
## 2001-01-04 1.13866156 -1.8744860 3.097705 5.796973
## 2001-01-05 4.70578690 0.3288287 7.551532 1.884401
## 2001-01-06 4.57036721 0.5461322 7.186784 13.359801
## MeanRelativeHumidity MinRelativeHumidity MaxRelativeHumidity
## 2001-01-01 78.73709 65.15411 100.00000
## 2001-01-02 69.70800 57.43761 94.71780
## 2001-01-03 70.69610 58.77432 94.66823
## 2001-01-04 76.89156 66.84256 95.80950
## 2001-01-05 76.67424 62.97656 100.00000
## 2001-01-06 89.01940 74.25754 100.00000
## Radiation WindSpeed WindDirection PET
## 2001-01-01 12.89251 2.000000 172 1.3212770
## 2001-01-02 13.03079 7.662544 278 2.2185985
## 2001-01-03 16.90722 2.000000 141 1.8045176
## 2001-01-04 11.07275 2.000000 172 0.9200627
## 2001-01-05 13.45205 7.581347 321 2.2914449
## 2001-01-06 12.84841 6.570501 141 1.7255058
Simulation models in medfate
have been designed to work along with data generated from package meteoland
. The user is strongly recommended to resort to this package to obtain suitable weather input for soil water balance simulations.
Apart from data inputs, the behaviour of simulation models can be controlled using a set of global parameters. The default parameterization is obtained using function defaultControl()
:
control = defaultControl()
The names of all control parameters are:
names(control)
## [1] "verbose" "subdailyResults"
## [3] "soilFunctions" "snowpack"
## [5] "drainage" "transpirationMode"
## [7] "hydraulicCostFunction" "verticalLayerSize"
## [9] "nStemSegments" "capacitance"
## [11] "cavitationRefill" "klat"
## [13] "taper" "numericParams"
## [15] "averageFracRhizosphereResistance" "Catm"
## [17] "ndailysteps" "thermalCapacityLAI"
## [19] "defaultWindSpeed" "storagePool"
Most of these parameters should normally be left to their default value. However, there are some that deserve explanation here:
verbose = FALSE
.soilFunctions
.transpirationMode
.ndailysteps
.To use the complex soil water balance model we must change the values of transpirationMode
(to switch from “Simple” to “Complex”) and soilFunctions
(to switch from Saxton’s retention curve, “SX”, to Van Genuchten’s retention curve, “VG”):
control$transpirationMode = "Complex"
control$soilFunctions = "VG"
A last step is needed before calling simulation functions. It consists in the compilation of all aboveground and belowground parameters and the specification of additional parameter values for each plant cohort, such as their light extinction coefficient or their response to drought. This is done by calling function spwbInput()
and taking species parameter values from species parameter data:
x = spwbInput(above, below, examplesoil, SpParamsMED, control)
If one has a forest
object, the spwbInput
object can be generated in directly from it, avoiding the need to explicitly build aboveground and belowground data frames:
x = forest2spwbInput(exampleforest, examplesoil, SpParamsMED, control)
All the input information for forest data and species parameter values can be inspected by printing different elements of the input object, which are:
names(x)
## [1] "control" "canopy" "cohorts"
## [4] "above" "below" "paramsBase"
## [7] "paramsAnatomy" "paramsTransp" "paramsWaterStorage"
## [10] "Transpiration" "Photosynthesis" "PLCstem"
## [13] "RWCsympstem" "RWCsympleaf" "Einst"
## [16] "psiRhizo" "psiRoot" "psiStem"
## [19] "psiLeaf"
First, information about the cohort species is found in element cohorts
(i.e. code, species and name):
x$cohorts
## SP Name
## T1_54 54 Pinus halepensis
## T2_68 68 Quercus ilex
## S1_65 65 Quercus coccifera
Element canopy
contains state variables of the whole canopy, which include growth degree days and canopy temperature:
x$canopy
## $gdd
## [1] 0
##
## $Temp
## [1] NA
Element above
contains the aboveground structure data that we already know:
x$above
## H CR LAI_live LAI_expanded LAI_dead
## T1_54 800 0.7150421 0.81630007 0.81630007 0
## T2_68 660 0.6055507 0.79744714 0.79744714 0
## S1_65 30 0.9740000 0.08911235 0.08911235 0
As with the soil input object, the spwbInput
object may be modified during simulations. In the case of soil water balance, these modifications concern LAI_expanded
and LAI_dead
in element above
, as well as canopy variables.
Aboveground parameters related to plant transpiration can be seen in paramsTransp
:
x$paramsTransp
## Gwmin Gwmax Vmax298 Jmax298 Kmax_stemxylem Kmax_rootxylem
## T1_54 0.00203 0.3766111 20.0 50.0000 0.1500000 0.890000
## T2_68 0.00450 0.4518345 66.2 130.0000 0.7735834 3.094334
## S1_65 0.01045 0.4518345 100.0 163.6253 0.2900000 1.160000
## VCleaf_kmax VCleaf_c VCleaf_d VCstem_kmax VCstem_c VCstem_d
## T1_54 6 2.204804 -2.466667 1.257285 2.204804 -3.7
## T2_68 8 2.188300 -4.138200 4.299811 1.158507 -4.4
## S1_65 8 2.931921 -5.266667 12.080826 2.931921 -7.9
## VCroot_kmax VCroot_c VCroot_d pRootDisc
## T1_54 35.08225 1.056900 -1.240000 0
## T2_68 63.95962 1.494100 -2.130000 0
## S1_65 26.37917 2.500429 -5.423197 0
Belowground parameters can be seen in below
and include root distribution as well as maximum root and rhizosphere conductances by soil layer:
x$below
## $V
## 1 2
## T1_54 0.1933638 0.8066362
## T2_68 0.1933638 0.8066362
## S1_65 0.8981075 0.1018925
##
## $VGrhizo_kmax
## 1 2
## T1_54 28167894 42557532
## T2_68 55150059 85968578
## S1_65 2046525891 86760581
##
## $VCroot_kmax
## 1 2
## T1_54 22.71909 12.36316
## T2_68 41.41993 22.53970
## S1_65 11.36508 15.01409
Finally, note that one can play with plant-specific parameters for soil water balance (instead of using species-level values) by modifying manually the parameter values in this object.
Before using the soil water balance model, is important to understand its sub-model components. Package medfate
provides many functions corresponding to sub-models (light extinction, hydraulics, transpiration, photosynthesis…). In addition, there are high-level plotting functions that allow examining several aspects of these processes.
Given a spwbInput
object, we can use function hydraulics.vulnerabilityCurvePlot()
to inspect vulnerability curves, i.e. how hydraulic conductance changes with the water potential, for each plant cohort and each of the different elements of the soil-plant hydraulic network: rhizosphere, roots, stems and leaves:
par(mar=c(5,5,1,1), mfrow=c(2,2))
hydraulics.vulnerabilityCurvePlot(x, type="leaf")
hydraulics.vulnerabilityCurvePlot(x, type="stem")
hydraulics.vulnerabilityCurvePlot(x, type="root")
hydraulics.vulnerabilityCurvePlot(x, examplesoil, type="rhizo")
Note that the last call to hydraulics.vulnerabilityCurvePlot()
includes a soil
object. This is because the van Genuchten parameters that define the vulnerability curve for the rhizosphere are stored in this object.
The vulnerability curves conformng the hydraulic network are used in the model to build the supply function, which relates water flow (i.e. transpiration) with the drop of water potential along the continuum. The supply function contains not only these two variables, but also the water potential of intermediate nodes in the the hydraulic network. Function hydraulics.supplyFunctionPlot()
can be used to inspect any of this variables:
par(mar=c(5,5,1,1), mfrow=c(2,2))
hydraulics.supplyFunctionPlot(x, examplesoil, type="E")
hydraulics.supplyFunctionPlot(x, examplesoil, type="ERhizo")
hydraulics.supplyFunctionPlot(x, examplesoil, type="dEdP")
hydraulics.supplyFunctionPlot(x, examplesoil, type="psiStem")
Calls to hydraulics.supplyFunctionPlot()
always need both a spwbInput
object and a soil
object. The soil moisture state (i.e. its water potential) is the starting point for the calculation of the supply function, so different curves will be obtained for different values of soil moisture.
The soil water balance model determines stomatal conductance and transpiration separately for sunlit and shade leaves. Stomatal conductance is determined after building a photosynthesis function corresponding to the supply function and finding the value of stomatal conductance that maximizes carbon revenue while avoiding hydraulic damage. Given a meteorological input and a chosen moment in time (day and timestep), the R function transp.stomatalRegulationPlot()
allows displaying the supply and photosynthesis curves for sunlit and shade leaves, along with an indication of the values corresponding to the chosen stomatal conductance:
d = 100
transp.stomatalRegulationPlot(x, examplesoil, examplemeteo, day = d, timestep=12,
latitude = 41.82592, elevation = 100)
Soil water balance simulations will normally span periods of several months or years, but since the model operates at a daily and subdaily temporal scales, it is possible to perform soil water balance for one day only. This is done using function spwb.day()
. In the following code we select the same day as before from the meteorological input data and perform soil water balance for that day only:
sd1<-spwb.day(x, examplesoil, rownames(examplemeteo)[d],
examplemeteo$MinTemperature[d], examplemeteo$MaxTemperature[d],
examplemeteo$MinRelativeHumidity[d], examplemeteo$MaxRelativeHumidity[d],
examplemeteo$Radiation[d], examplemeteo$WindSpeed[d],
latitude = 41.82592, elevation = 100,
slope= 0, aspect = 0, prec = examplemeteo$Precipitation[d])
The output of spwb.day()
is a list with five elements:
names(sd1)
## [1] "cohorts" "WaterBalance" "EnergyBalance" "Soil"
## [5] "SoilInst" "RhizoPsi" "Plants" "PlantsInst"
Element WaterBalance
contains the soil water balance flows of the day (precipitation, infiltration, transpiration, …)
sd1$WaterBalance
## Rain Snow NetRain
## 0.0000000 0.0000000 0.0000000
## Runon Infiltration Runoff
## 0.0000000 0.0000000 0.0000000
## DeepDrainage SoilEvaporation PlantExtraction
## 0.0000000 0.1795397 0.9797028
## Transpiration HydraulicRedistribution LAIcell
## 0.9797028 0.0000000 1.7028596
## LAIcelldead Cm Lground
## 0.0000000 1.2373017 0.4082981
And Soil
contains water evaporated from each soil layer, water transpired from each soil layer and the final soil water potential:
sd1$Soil
## SoilEvaporation HydraulicInput HydraulicOutput PlantExtraction
## 1 1.795396e-01 0 0.6064550 0.6064550
## 2 5.492161e-08 0 0.3732478 0.3732478
## psi
## 1 -0.03493828
## 2 -0.03346174
Element EnergyBalance
contains subdaily variation in atmosphere, canopy and soil temperatures, as well as canopy and soil energy balance components.
names(sd1$EnergyBalance)
## [1] "Temperature" "CanopyEnergyBalance" "SoilEnergyBalance"
Package medfate
provides a plot
function for objects of class spwb.day
that can be used to inspect the results of the simulation. We use this function to display subdaily dynamics in plant, soil and canopy variables. For example, we can use it to display temperature variations (only the temperature of the topmost soil layer is drawn):
par(mar=c(4,4,1,1))
plot(sd1, type = "Temperature")
Element Plants
contains output values by plant cohort and subdaily time step. Many output variables can be inspected in this element.
names(sd1$Plants)
## [1] "LAI" "Extraction" "Transpiration" "RootPsi"
## [5] "StemPsi" "StemPLC" "LeafPsiMin" "LeafPsiMax"
## [9] "dEdP" "DDS" "StemRWC" "LeafRWC"
With the following code, we use function plot.spwb.day()
to draw plots of daily variation per species of plant transpiration (L·m\(^{-2}\)), transpiration per leaf area (in L·m\(^{-2}\)·LAI\(^{-1}\)), net photosynthesis per leaf area (in g C·m\(^{-2}\)·LAI\(^{-1}\)), and plant water potential (in MPa):
par(mar=c(5,5,1,1), mfrow=c(2,2))
plot(sd1, type = "PlantTranspiration", bySpecies = T)
plot(sd1, type = "LeafTranspiration", bySpecies = T)
plot(sd1, type = "LeafPhotosynthesis", bySpecies = T)
plot(sd1, type = "LeafPsi", bySpecies = T)
We can also use the same function to plot variations in temperature for sunlit and shade leaves (note that sunlit leaves of some species reach temperatures higher than the canopy):
plot(sd1, type = "LeafTemperature", bySpecies=TRUE)
Or we can plot variations in stomatal conductance for sunlit and shade leaves:
plot(sd1, type = "LeafStomatalConductance", bySpecies=TRUE)
Users will often use function spwb()
to run the soil water balance model for several days. This function requires the spwbInput
object, the soil
object and the meteorological data frame. However, running spwb.day()
modified the input objects. In particular, the soil moisture at the end of the simulation was:
examplesoil$W
## [1] 0.9885153 0.9970670
And the temperature of soil layers:
examplesoil$Temp
## [1] 8.305204 3.157197
We can also see the current state of canopy variables:
x$canopy
## $gdd
## [1] 1.232373
##
## $Temp
## [1] 5.738397
We simply use function spwb.resetInputs()
to reset state variables to their default values, so that the new simulation is not affected by the end state of the previous simulation:
spwb.resetInputs(x, examplesoil)
examplesoil$W
## [1] 1 1
examplesoil$Temp
## [1] NA NA
x$canopy
## $gdd
## [1] 0
##
## $Temp
## [1] NA
Now we are ready to call function spwb()
. In this example, we only simulate 61 days to save computational time:
S = spwb(x, examplesoil, examplemeteo[110:170,], latitude = 41.82592, elevation = 100)
## W1i:1 W2i:1
## Daily balance:.............................................................done
## Precipitation (mm) 36
## Rain (mm) 36 Snow (mm) 0
## Interception (mm) 11 Net rainfall (mm) 25
## Infiltration (mm) 25 Runoff (mm) 0 Deep drainage (mm) 11
## Soil evaporation (mm) 3 Transpiration (mm) 78
## Plant extraction from soil (mm) 78 Hydraulic redistribution (mm) 6
## W1f:0.67 W2f:0.66
## Final soil water content (mm): 131
##
## Building SPWB output ...done.
Function spwb()
returns an object of class spwb, actually a list:
class(S)
## [1] "spwb" "list"
If we inspect its elements, we realize that the output is arranged differently than in spwb.day()
:
names(S)
## [1] "spwbInput" "soilInput" "WaterBalance"
## [4] "Soil" "EnergyBalance" "Temperature"
## [7] "PlantLAI" "PlantAbsorbedSWR" "PlantAbsorbedLWR"
## [10] "PlantTranspiration" "PlantPhotosynthesis" "dEdP"
## [13] "LeafPsiMin" "LeafPsiMax" "LeafRWC"
## [16] "StemPsi" "StemPLC" "StemRWC"
## [19] "RootPsi" "RhizoPsi" "PlantStress"
## [22] "subdaily"
In particular, element spwbInput
contains a copy of the input parameters that were used to run the model:
names(S$spwbInput)
## [1] "control" "canopy" "cohorts"
## [4] "above" "below" "paramsBase"
## [7] "paramsAnatomy" "paramsTransp" "paramsWaterStorage"
## [10] "Transpiration" "Photosynthesis" "PLCstem"
## [13] "RWCsympstem" "RWCsympleaf" "Einst"
## [16] "psiRhizo" "psiRoot" "psiStem"
## [19] "psiLeaf"
As before, WaterBalance
contains water balance components, but in this case in form of a data frame with days in rows:
head(S$WaterBalance)
## GDD LAIcell LAIcelldead Cm Lground PET
## 2001-04-20 0.000000 1.70286 0 1.237302 0.4082981 0
## 2001-04-21 0.000000 1.70286 0 1.237302 0.4082981 0
## 2001-04-22 0.000000 1.70286 0 1.237302 0.4082981 0
## 2001-04-23 0.463942 1.70286 0 1.237302 0.4082981 0
## 2001-04-24 5.058209 1.70286 0 1.237302 0.4082981 0
## 2001-04-25 7.532502 1.70286 0 1.237302 0.4082981 0
## Precipitation Rain Snow NetRain Infiltration Runoff
## 2001-04-20 2.1625135 2.1625135 0 0.88295013 0.88295013 0
## 2001-04-21 3.7992356 3.7992356 0 2.24524902 2.24524902 0
## 2001-04-22 4.1962782 4.1962782 0 2.59530545 2.59530545 0
## 2001-04-23 1.4698434 1.4698434 0 0.60013427 0.60013427 0
## 2001-04-24 0.1538991 0.1538991 0 0.06283671 0.06283671 0
## 2001-04-25 0.1317847 0.1317847 0 0.05380744 0.05380744 0
## DeepDrainage Evapotranspiration Interception SoilEvaporation
## 2001-04-20 0.8829501 2.123090 1.27956335 0.25423061
## 2001-04-21 1.4017222 2.433345 1.55398659 0.24460077
## 2001-04-22 1.7159467 2.555043 1.60097277 0.22075970
## 2001-04-23 0.0000000 1.909326 0.86970916 0.16697936
## 2001-04-24 0.0000000 1.199514 0.09106239 0.08817851
## 2001-04-25 0.0000000 1.321388 0.07797725 0.05928082
## PlantExtraction Transpiration HydraulicRedistribution
## 2001-04-20 0.5892962 0.5892962 0
## 2001-04-21 0.6347580 0.6347580 0
## 2001-04-22 0.7333101 0.7333101 0
## 2001-04-23 0.8726371 0.8726371 0
## 2001-04-24 1.0202731 1.0202731 0
## 2001-04-25 1.1841298 1.1841298 0
Elements PlantLAI
, PlantStress
, LeafPsi
, PlantTranspiration
and PlantPhotosynthesis
contain daily values by plant cohorts, for example plant water potentials are:
head(S$LeafPsi)
## NULL
Package medfate
also provides a plot
function for objects of class spwb
. It can be used to show the meteorological input. Additionally, it can also be used to draw soil and plant variables. In the code below we draw water fluxes, soil water potentials, plant transpiration and plant (mid-day) water potential:
par(mar=c(5,5,1,1), mfrow=c(2,2))
plot(S, type="Evapotranspiration", bySpecies = TRUE)
plot(S, type="SoilPsi", bySpecies = TRUE)
plot(S, type="PlantTranspiration", bySpecies = TRUE)
plot(S, type="LeafPsiMin", bySpecies = TRUE)
While the simulation model uses daily steps, users may be interested in outputs at larger time scales. The package provides a summary
for objects of class spwb
. This function can be used to summarize the model’s output at different temporal steps (i.e. weekly, annual, …). For example, to obtain the average soil moisture and water potentials by months one can use:
summary(S, freq="months",FUN=mean, output="Soil")
## W.1 W.2 ML.1 ML.2 MLTot WTD SWE
## 2001-04-01 0.9716282 0.9887216 66.49646 125.82252 192.3190 1000 0
## 2001-05-01 0.9041201 0.9353835 61.87634 119.03483 180.9112 1000 0
## 2001-06-01 0.7344578 0.7373630 50.26495 93.83518 144.1001 1000 0
## PlantExt.1 PlantExt.2 HydraulicInput.1 HydraulicInput.2
## 2001-04-01 0.4976419 0.4475738 0.01694259 0
## 2001-05-01 0.4801647 0.7666154 0.09464791 0
## 2001-06-01 0.4161221 1.0859666 0.17146613 0
## psi.1 psi.2
## 2001-04-01 -0.03954751 -0.03538743
## 2001-05-01 -0.05892868 -0.04785538
## 2001-06-01 -0.15437834 -0.13649334
Parameter output
is used to indicate the element of the spwb
object for which we desire summaries. Similarly, it is possible to calculate the average stress of plant cohorts by months:
summary(S, freq="months",FUN=mean, output="PlantStress")
## T1_54 T2_68 S1_65
## 2001-04-01 0.007192883 0.01599768 0.0004713716
## 2001-05-01 0.013106652 0.02475831 0.0013300643
## 2001-06-01 0.024550669 0.04282354 0.0027279681
The summary
function can be also used to aggregate the output by species. In this case, the values of plant cohorts belonging to the same species will be averaged using LAI values as weights. For example, we may average the daily drought stress across cohorts of the same species (here there is only one cohort by species, so this does not modify the output):
head(summary(S, freq="day", output="PlantStress", bySpecies = TRUE))
## Pinus halepensis Quercus coccifera Quercus ilex
## 2001-04-20 0.003332744 9.900714e-05 0.008806105
## 2001-04-21 0.003995170 1.824316e-04 0.008874972
## 2001-04-22 0.004233576 1.388487e-04 0.011669000
## 2001-04-23 0.006375923 2.595975e-04 0.014082426
## 2001-04-24 0.008223007 4.496005e-04 0.017051551
## 2001-04-25 0.010196826 6.701092e-04 0.020149900
Or we can combine the aggregation by species with a temporal aggregation (here monthly averages):
summary(S, freq="month", FUN = mean, output="PlantStress", bySpecies = TRUE)
## Pinus halepensis Quercus coccifera Quercus ilex
## 2001-04-01 0.007192883 0.0004713716 0.01599768
## 2001-05-01 0.013106652 0.0013300643 0.02475831
## 2001-06-01 0.024550669 0.0027279681 0.04282354
Function spwb()
does not store subdaily values in the output, because otherwise a lot of memory would be consumed. This has the drawback that, while subdaily steps are simulated, there is no trace of them. If we detect something unusual in the results at the daily scale, we may decide to inspect this in detail by re-running a particular day. Say we want to inspect day 5 of the previous simulation, which correspond to day 144 of the whole year. We start by setting the state variables to the end of the day before:
spwb.resetInputs(x, examplesoil, from = S, day = 4)
Now we call function spwb.day()
for day 144 of the year:
d = 144
sd1<-spwb.day(x, examplesoil, rownames(examplemeteo)[d],
examplemeteo$MinTemperature[d], examplemeteo$MaxTemperature[d],
examplemeteo$MinRelativeHumidity[d],
examplemeteo$MaxRelativeHumidity[d],
examplemeteo$Radiation[d], examplemeteo$WindSpeed[d],
latitude = 41.82592, elevation = 100,
slope= 0, aspect = 0, prec = examplemeteo$Precipitation[d])
We use function plot.spwb.day()
to draw plots of daily variation per species of plant transpiration (L·m\(^{-2}\)), transpiration per leaf area (in L·m\(^{-2}\)·LAI\(^{-1}\)), net photosynthesis per leaf area (in g C·m\(^{-2}\)·LAI\(^{-1}\)), and plant water potential (in MPa):
par(mar=c(5,5,1,1), mfrow=c(2,2))
plot(sd1, type = "PlantTranspiration", bySpecies = T)
plot(sd1, type = "LeafTranspiration", bySpecies = T)
plot(sd1, type = "LeafPhotosynthesis", bySpecies = T)
plot(sd1, type = "LeafPsi", bySpecies = T)