View source: R/cliBrtSunDurFrc.R
cliBrtSunDurFrcPoints | R Documentation |
Estimates monthly averages for daily fraction of bright sunshine duration, for a given geographical location (latitude, longitude, and elevation) and year, by using the monthly time series of temperature and precipitation.
cliBrtSunDurFrcPoints(
temp,
prec,
lat,
lon,
elv,
year = 2000,
aprchSIM = c("Solar123", "SPLASH")
)
temp |
'numeric' R object with one-year time series of monthly mean air temperature (in °C) |
prec |
'numeric' R object with one-year time series of monthly precipitation sum (in mm) |
lat |
'numeric' vector with the latitude coordinates (in decimal degrees) |
lon |
'numeric' vector with the longitude coordinates (in decimal degrees) |
elv |
'numeric' vector with the elevation values (in meters above sea level) |
year |
'numeric' vector with values of the year (using astronomical year numbering) |
aprchSIM |
'character' vector of length 1 that indicates the formula used to estimate the value of solar
irradiance/irradiation for a specific day. Valid values are as follows: |
To estimate the monthly averages of relative sunlight duration, the approach presented by Yin (1999) is
implemented here. Many variables in this estimation scheme can be easily and unambiguously determined, but the
approach uses two important quantities, the calculation method of which can be chosen here depending on the
purpose of the investigations. One of them is the estimated value of the mean hourly solar irradiance under
cloudless-sky conditions. This quantity can be estimated in this implementation of the approach with the
original method (aprchSIM = 'Solar123'
) or with the solar radiation model used in the SPLASH algorithm,
considering the variability of orbital parameters of the Earth over time (aprchSIM = 'SPLASH'
). The
latter is recommended for paleo-climatological and paleo-environmental studies. These solar radiation models
is also applied to calculate the daylength, whose monthly averages are used to estimate monthly averages of
daily potential evapotranspiration (Eqs. A10 and A11 in Yin (1998)).
The procedure proposed by Yin (1999) requires the calculation of several regional factors (see Eq 3.3 in Yin
(1999)). Each regional factor is activated as a function of latitude and longitude. However, it is important
to note that in this implementation, these factors are activated with the current configuration of continents
and islands. Continents and regions are classified using the medium-resolution world map of the
rnaturalearthdata
, if the high-resolution world map of the
rnaturalearthhires is not available. In checking whether
or not a given geographic location can be defined as an island, the high-resolution world maps (version 5.1.1)
of the Natural Earth are applied.
A 12-column matrix with monthly averages of relative sunshine duration.
As with any function with a point mode, a set of basic input data is defined here. In this case, they are as
follows: 'temp'
(one-year time series of monthly mean air temperature), and 'prec'
(one-year time
series of monthly precipitation sum). The objects 'temp'
and 'prec'
must be either 12-length
vectors or 12-column matrices. The first dimensions of these matrices have to be the same length. The function
automatically converts vectors into single-row matrices during the error handling, and then uses these
matrices. The first dimensions of these matrices determines the number of rows in the result matrix. In the
case of arguments that do not affect the course of the calculation procedure or the structure of the return
object, scalar values (i.e., 'numeric' vector of length 1) may also be allowed. In this case, they are as
follows: 'lat'
(latitude coordinates in decimal degrees), 'lon'
(longitude coordinates in decimal
degrees), 'elv'
(elevation in meters above sea level), and 'year'
(year using astronomical year
numbering). These scalars are converted to vectors by the function during the error handling, and these vectors
are applied in the further calculations. If these data are stored in vectors of length at least 2, their length
must be the same size of first dimension of the matrices containing the basic data.
Allen RG (1996) Assessing integrity of weather data for reference evapotranspiration estimation. J Irrig Drain Eng 122(2):97–106. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1061/(ASCE)0733-9437(1996)122:2(97)")}
Berger A, Loutre MF (1991) Insolation values for the climate of the last 10 million years. Quat Sci Rev 10(4):297-317. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1016/0277-3791(91)90033-Q")}
Brock TD (1981) Calculating solar radiation for ecological studies. Ecol Model 14(1–2):1-19. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1016/0304-3800(81)90011-9")}
Davis TW, Prentice IC, Stocker BD, Thomas RT, Whitley RJ, Wang H, Evans BJ, Gallego-Sala AV, Sykes MT, Cramer W (2017) Simple process-led algorithms for simulating habitats (SPLASH v.1.0): robust indices of radiation, evapotranspiration and plant-available moisture. Geosci Model Dev 10(2):689–708. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.5194/gmd-10-689-2017")}
Duffie JA, Beckman WA (1991) Solar Engineering of Thermal Processes. Second Edition. Wiley-Interscience, New York, NY
Yin X (1997a) Calculating daytime mean relative air mass. Agric For Meteorol 87(2-3):85-90. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1016/S0168-1923(97)00029-4")}
Yin X (1997b) Optical Air Mass: Daily Integration and its Applications. Meteorol Atmos Phys 63(3-4):227-233. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1007/BF01027387")}
Yin X (1998) The Albedo of Vegetated Land Surfaces: Systems Analysis and Mathematical Modeling. Theor Appl Climatol 60(1–4):121–140. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1007/s007040050038")}
Yin X (1999) Bright Sunshine Duration in Relation to Precipitation, Air Temperature and Geographic Location. Theor Appl Climatol 64(1–2):61–68. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1007/s007040050111")}
library (graphics)
# Loading mandatory data for the Example 'Points'
data(inp_exPoints)
# Measured and estimated one-year time series of the monthly mean relative sunshine duration,
# at a grid cell near Szeged, Hungary (46.3N, 20.2E), in the year 2010
with(inp_exPoints, {
bsdf01 <- matrix(nrow = 0, ncol = 12, dimnames = list(NULL, month.abb))
bsdf01 <- rbind(bsdf01, "Measured" = bsdf["2010", ])
bsdf01 <- rbind(bsdf01, "Solar123" = cliBrtSunDurFrcPoints(temp["2010", ], prec["2010", ],
lat, lon, elv, year = 2010))
bsdf01 <- rbind(bsdf01, "SPLASH" = cliBrtSunDurFrcPoints(temp["2010", ], prec["2010", ],
lat, lon, elv, year = 2010, aprchSIM = "SPLASH"))
cols <- c("black", "green", "blue")
matplot(t(bsdf01), type = "l", lwd = 2, col = cols, xaxt = "n", xlab = "Month",
ylab = "Average relative sunshine duration (unitless)")
axis(1, at = seq(1, ncol(bsdf01)), labels = colnames(bsdf01))
legend(1, 0.7, legend = rownames(bsdf01), col = cols, lty = 1 : 2, lwd = 2, xpd = TRUE)
})
# Relative root mean square error between measured and estimated values for the 'bsdf',
# at a grid cell near Szeged, Hungary (46.3N, 20.2E), in the period 1981-2010
with(inp_exPoints, {
years <- seq(1981, 2010)
bsdf02 <- cliBrtSunDurFrcPoints(temp, prec, lat, lon, elv, year = years)
rrmse <- function(pre, obs) { (sqrt(mean((pre - obs) ^ 2.)) / mean(obs)) * 100. }
rrmse_bsdf <- sapply(1 : 12, function(i) { rrmse(bsdf02[, i], bsdf[, i]) })
cols <- c("black", "green")
plot(rrmse_bsdf, type = "l", lwd = 2, col = cols, xaxt = "n", xlab = "Month",
ylab = "Relative root mean square error (%)")
axis(1, at = 1 : 12, labels = month.abb)
})
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.