Title: | Multi-Site Auto-Regressive Weather GENerator |
---|---|
Description: | S3 and S4 functions are implemented for spatial multi-site stochastic generation of daily time series of temperature and precipitation. These tools make use of Vector AutoRegressive models (VARs). The weather generator model is then saved as an object and is calibrated by daily instrumental "Gaussianized" time series through the 'vars' package tools. Once obtained this model, it can it can be used for weather generations and be adapted to work with several climatic monthly time series. |
Authors: | Emanuele Cordano, Emanuele Eccel |
Maintainer: | Emanuele Cordano <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.3.9.3 |
Built: | 2024-11-10 03:13:48 UTC |
Source: | https://github.com/ecor/rmawgen |
Multi-site autoregressive Models for Daily Weather Generation. The modeling in climate change applications for agricultural or hydrological purposes often requires daily time-series of precipitation and temperature. This is the case of downscaled series from monthly or seasonal predictions of Global Climate Models (GCMs). The R package RMAWGEN (R Multi-Sites Auto regressive Weather GENerator) is built to generate daily temperature and precipitation time series in several sites by using the theory of vectorial autoregressive models (VAR). The VAR model is used because it is able to maintain the temporal and spatial correlations among the several series. In particular, observed time series of daily maximum and minimum temperature and precipitation are used to calibrate the parameters of a VAR model (saved as ”GPCAvarest2” or ”varest2” classes, which inherit the "varest" S3 class defined in the package vars [Pfaff, 2008]). Therefore the VAR model, coupled with monthly mean weather variables downscaled by GCM predictions, allows to generate several stochastic daily scenarios. The structure of the package consists in functions that transform precipitation and temperature time series into Gaussian-distributed random variables through deseasonalization and Principal Component Analysis. Then a VAR model is calibrated on transformed time series. The time series generated by VAR are then inversely re transformed into precipitation and/or temperature series. An application dateset is included in the RMAWGEN package as an example; it is presented by using a dataset with daily weather time series recorded in 59 different sites of Trentino (Italy) and its neighborhoods for the period 1958-2007. The software is distributed as a Free Software with General Public License (GPL) and is available on CRAN website. (https://cran.r-project.org/package=RMAWGEN) . A presentation of the package is available on https://docs.google.com/file/d/0B8xDtMCnW3dJU2JIemVqMnpKTHc/edit. Example script files about package usage are available on https://github.com/ecor/RMAWGENCodeCorner.
Package: | RMAWGEN |
Type: | Package |
Version: | 1.3.6 |
Date: | 2019-11-13 |
License: | GPL (>= 2) |
LazyLoad: | yes |
Depends: R(>=2.12),time,chron,vars | |
First release of RMAWGEN was created in the frame of ACE-SAP (http://www.ace-sap.it/) and ENVIROCHANGE (http://www.envirochange.cc/) projects funded by Provincia Autonoma di Trento (http://www.provincia.tn.it/).
RMAWGEN is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
RMAWGEN is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/.
Emanuele Cordano [email protected], Emanuele Eccel [email protected]
Cordano E. and Eccel E. (2016), Tools for stochastic weather series generation in R environment, Italian Journal of Agrometeorology doi:10.19199/2016.3.2038-5625.031 https://doi.org/10.19199/2016.3.2038-5625.031 Pfaff B. (2008). VAR, SVAR and SVEC Models: Implementation Within R Package vars. Journal of Statistical Software 27(4). https://www.jstatsoft.org/v27/i04/(doi:10.18637/jss.v027.i04)
Useful links:
https://docs.google.com/file/d/0B66otCUk3Bv6V3RPbm1mUG4zVHc/edit
http://presentations.copernicus.org/EGU2012-14026_presentation.pdf
http://presentations.copernicus.org/EGU2012-5404_presentation.pdf
Plots the auto- and cross- covariance functions between measured and simulated data for several stations
acvWGEN(measured, simulated, titles = c("Sim.", "Mes."), station = NULL)
acvWGEN(measured, simulated, titles = c("Sim.", "Mes."), station = NULL)
measured |
matrix containing measured time series |
simulated |
matrix containing simulated time series |
titles |
title suffixes for the simulated and measured data respectively c("Sim.","Mes.") |
station |
string vector containing the IDs of the meteorological stations where the autocovariance is calculated.
If it is |
0 in case of success
It uses acf
function
Inserts three columns (year,month,day) passing dates to a matrix or to a dataframe
adddate(data, origin = "1961-1-1")
adddate(data, origin = "1961-1-1")
data |
matrix of daily data |
origin |
character string containing the date of the first row of |
a data frame with dates and data
values
Adds suffixes for daily maximum and minimum temperature to the names of a column data frame
addsuffixes( names = c("T0001", "T0099", "T0001", "T0099"), suffix = c("_Tx", "_Tn"), sep = "" )
addsuffixes( names = c("T0001", "T0099", "T0001", "T0099"), suffix = c("_Tx", "_Tn"), sep = "" )
names |
a character string vector with column names |
suffix |
suffixes to add to the first and second groups of column names respectively |
sep |
separation element |
This function is used for data frames with duplicated field names
the vector of names with suffixes added
names <- addsuffixes()
names <- addsuffixes()
arch.test
function for varest2
objectarch.test
function for varest2
object
arch_test(object, interval = NULL, overlap = 20, list.output = FALSE, ...)
arch_test(object, interval = NULL, overlap = 20, list.output = FALSE, ...)
object |
a |
interval |
string or subset interval of time (e.g. days) or length of this subset interval to which the ARCH test is applied (see Note). Default is |
overlap |
number of time instants (e.g. days) which are overlapped on two different subsequent intervals. Default is 20. It is used only if |
list.output |
logical value. If |
... |
further arguments for |
This function is a wrapper of arch.test
. It can compute the test also for some subsets (intervals) of the time-series or for all the time-series divided in overlapping intervals. The intervals considered for the ARCH test are defined with the argument interval
. If interval
is an integer number instead of a vector, it indicates the length of the intervals in which the time-series is split. If interval
is set to NULL
, the test is done on the comprehensive residual time-series without splitting.
One object or a list of objects with class attribute varcheck
as reported in arch.test
The comprehensive Precipitation Generator
ComprehensivePrecipitationGenerator( station = c("T0001", "T0010", "T0099"), prec_all, mean_climate_prec = NULL, year_max = 1990, year_min = 1961, leap = TRUE, nmonth = 12, cpf = NULL, verbose = TRUE, p = 1, type = "none", lag.max = NULL, ic = "AIC", activateVARselect = FALSE, exogen = NULL, exogen_sim = NULL, is_exogen_gaussian = FALSE, year_max_sim = year_max, year_min_sim = year_min, mean_climate_prec_sim = NULL, onlygeneration = FALSE, varmodel = NULL, type_quantile = 3, qnull = NULL, valmin = 0.5, step = 0, n_GPCA_iteration = 0, n_GPCA_iteration_residuals = n_GPCA_iteration, sample = NULL, extremes = TRUE, exogen_all = NULL, exogen_all_col = station, no_spline = FALSE, nscenario = 1, seed = NULL, noise = NULL )
ComprehensivePrecipitationGenerator( station = c("T0001", "T0010", "T0099"), prec_all, mean_climate_prec = NULL, year_max = 1990, year_min = 1961, leap = TRUE, nmonth = 12, cpf = NULL, verbose = TRUE, p = 1, type = "none", lag.max = NULL, ic = "AIC", activateVARselect = FALSE, exogen = NULL, exogen_sim = NULL, is_exogen_gaussian = FALSE, year_max_sim = year_max, year_min_sim = year_min, mean_climate_prec_sim = NULL, onlygeneration = FALSE, varmodel = NULL, type_quantile = 3, qnull = NULL, valmin = 0.5, step = 0, n_GPCA_iteration = 0, n_GPCA_iteration_residuals = n_GPCA_iteration, sample = NULL, extremes = TRUE, exogen_all = NULL, exogen_all_col = station, no_spline = FALSE, nscenario = 1, seed = NULL, noise = NULL )
station |
character vector of the IDs of the considered meteorological stations |
prec_all |
data frame containing daily precipitation of all meteorological stations. See |
mean_climate_prec |
a matrix containing monthly mean daily precipitation for the considered station. If it is |
year_max |
start year of the recorded (calibration) period |
year_min |
end year of the recorded (calibration) period |
leap |
logical variables. If it is |
nmonth |
number of months in one year (default is 12) |
cpf |
|
verbose |
logical variable |
p , type , lag.max , ic , activateVARselect
|
see respective input parameter on |
exogen |
data frame or matrix containing the (normalized or not) exogenous variables (predictors) for the recorded (calibration) period. |
exogen_sim |
data frame or matrix containing the (normalized or not) exogenous variables (predictors) for the simulation period. Default is |
is_exogen_gaussian |
logical value. If |
year_max_sim |
last year of the simulation period. Default is equal to |
year_min_sim |
first year of the simulation period. Default is equal to |
mean_climate_prec_sim |
a matrix containing monthly mean daily precipitation for the simulation period. If is |
onlygeneration |
logical value. If |
varmodel |
the comprehensinve VAR model as a |
type_quantile |
see |
step |
see |
n_GPCA_iteration |
number of iterations of Gaussianization process for data. Default is 0 (no Gaussianization) |
n_GPCA_iteration_residuals |
number of iterations of Gaussianization process for VAR residuals. Default is 0 (no Gaussianization) |
sample , extremes , qnull , valmin
|
|
exogen_all |
data frame containing exogenous variable formatted like |
exogen_all_col |
vector of considered columns of |
no_spline |
logical value. See |
nscenario |
number of generated scenarios for daily maximum and minimum temperature |
seed |
seed for stochastic random generation see |
noise |
stochastic noise to add for variabile generation. Default is |
A list of the following variables:
prec_mes
matrix containing measured daily precipitation (the data is copied by the measured data given as input for the period and the station considered for varmodel
estimation)
prec_spline
matrix containing climatic "spline-interpolated" daily preciptation from mean_climate_prec
data_prec
matrix containing normalized measured precipitation variable
prec_gen
matrix containing generated daily precipitation [mm]
prec_spline_sim
matrix containing climatic "spline-interpolated" daily preciptation from mean_climate_prec_sim
data_prec_gen
matrix containing normalized generated precipitation variable
mean_climate_prec
matrix containing monthly means of daily precipitation (historical scenario)
mean_climate_prec_sim
matrix containing monthly means of daily precipitation (predicted/simulated scenario)
var
a varest object containing the used VAR model
It pre-processes and generates a multi-site precipitation fields. It uses getVARmodel
. Detailed examples can be viewed of this function in this presentation.
Unfortunately, using this approach, the spatial correlations are underestimated. This is due to the persinstence of zeros in the precipitation records.
This problem is known in literature and can be solved in the future versions of RMAWGEN.
See the R code for further details
Emanuele Cordano, Emanuele Eccel
splineInterpolateMonthlytoDailyforSeveralYears
data(trentino) set.seed(1222) # set the seed for random generations! year_max <- 1990 year_min <- 1961 year_max_sim <- 1982 year_min_sim <- 1981 n_GPCA_iter <- 2 p <- 1 nscenario=1 station <- c("T0090","T0083") ## Not Run: the call to ComprehensivePrecipitationGenerator may elapse too ## long time (more than 5 eseconds) and is not executed by default CRAN check. ## Please uncomment the following line to run the example on your own PC. generation00 <- ComprehensivePrecipitationGenerator(station=station, prec_all=PRECIPITATION,year_min=year_min,year_max=year_max, year_min_sim=year_min_sim,year_max_sim=year_max_sim,p=p, n_GPCA_iteration=n_GPCA_iter,n_GPCA_iteration_residuals=0, sample="monthly",nscenario=nscenario,no_spline=TRUE)
data(trentino) set.seed(1222) # set the seed for random generations! year_max <- 1990 year_min <- 1961 year_max_sim <- 1982 year_min_sim <- 1981 n_GPCA_iter <- 2 p <- 1 nscenario=1 station <- c("T0090","T0083") ## Not Run: the call to ComprehensivePrecipitationGenerator may elapse too ## long time (more than 5 eseconds) and is not executed by default CRAN check. ## Please uncomment the following line to run the example on your own PC. generation00 <- ComprehensivePrecipitationGenerator(station=station, prec_all=PRECIPITATION,year_min=year_min,year_max=year_max, year_min_sim=year_min_sim,year_max_sim=year_max_sim,p=p, n_GPCA_iteration=n_GPCA_iter,n_GPCA_iteration_residuals=0, sample="monthly",nscenario=nscenario,no_spline=TRUE)
The Comprehensive Temperature Generator
ComprehensiveTemperatureGenerator( station = c("T0001", "T0010", "T0099"), Tx_all, Tn_all, mean_climate_Tn = NULL, mean_climate_Tx = NULL, Tx_spline = NULL, Tn_spline = NULL, year_max = 1990, year_min = 1961, leap = TRUE, nmonth = 12, verbose = TRUE, p = 1, type = "none", lag.max = NULL, ic = "AIC", activateVARselect = FALSE, year_max_sim = year_max, year_min_sim = year_min, mean_climate_Tn_sim = NULL, mean_climate_Tx_sim = NULL, Tn_spline_sim = NULL, Tx_spline_sim = NULL, onlygeneration = FALSE, varmodel = NULL, normalize = TRUE, type_quantile = 3, sample = NULL, extremes = TRUE, option = 2, yearly = FALSE, yearly_sim = yearly, n_GPCA_iteration = 0, n_GPCA_iteration_residuals = n_GPCA_iteration, exogen = NULL, exogen_sim = exogen, is_exogen_gaussian = FALSE, exogen_all = NULL, exogen_all_col = station, nscenario = 1, seed = NULL, noise = NULL )
ComprehensiveTemperatureGenerator( station = c("T0001", "T0010", "T0099"), Tx_all, Tn_all, mean_climate_Tn = NULL, mean_climate_Tx = NULL, Tx_spline = NULL, Tn_spline = NULL, year_max = 1990, year_min = 1961, leap = TRUE, nmonth = 12, verbose = TRUE, p = 1, type = "none", lag.max = NULL, ic = "AIC", activateVARselect = FALSE, year_max_sim = year_max, year_min_sim = year_min, mean_climate_Tn_sim = NULL, mean_climate_Tx_sim = NULL, Tn_spline_sim = NULL, Tx_spline_sim = NULL, onlygeneration = FALSE, varmodel = NULL, normalize = TRUE, type_quantile = 3, sample = NULL, extremes = TRUE, option = 2, yearly = FALSE, yearly_sim = yearly, n_GPCA_iteration = 0, n_GPCA_iteration_residuals = n_GPCA_iteration, exogen = NULL, exogen_sim = exogen, is_exogen_gaussian = FALSE, exogen_all = NULL, exogen_all_col = station, nscenario = 1, seed = NULL, noise = NULL )
station |
see respective input parameter on |
Tx_all , Tn_all , mean_climate_Tn , mean_climate_Tx , Tx_spline , Tn_spline
|
see respective input parameter on |
year_max , year_min , leap , nmonth , verbose
|
see respective input parameter on |
p , type , lag.max , ic , activateVARselect
|
see respective input parameter on |
year_max_sim |
last year of the simulation period. Default is equal to |
year_min_sim |
first year of the simulation period. Default is equal to |
mean_climate_Tn_sim |
monthly averaged daily minimum temperatures for the simulated scenario and used by the random generator . Default is |
mean_climate_Tx_sim |
monthly averaged daily maximum temperatures for the simulated scenario and used by the random generator . Default is |
Tn_spline_sim |
daily timeseries (from the first day of |
Tx_spline_sim |
daily timeseries (from the first day of |
onlygeneration |
logical variable. If |
varmodel |
the comprehensinve VAR model as a |
normalize , sample , extremes
|
see |
type_quantile |
see |
option |
integer value. If 1, the generator works with minimun and maximum temperature, if 2 (default) it works with the average value between maximum and minimum temparature and the respective daily thermal range. |
yearly |
logical value. If |
yearly_sim |
logical value. If |
n_GPCA_iteration |
number of iterations of Gaussianization process for data. Default is 0 (no Gaussianization) |
n_GPCA_iteration_residuals |
number of iterations of Gaussianization process for VAR residuals. Default is 0 (no Gaussianization) |
exogen |
data frame or matrix containing the (normalized or not) exogenous variables (predictors) for the recorded (calibration) period. Default is |
exogen_sim |
data frame or matrix containing the (normalized or not) exogenous variables (predictors) for the simulation period. Default is |
is_exogen_gaussian |
logical value, If |
exogen_all |
data frame containing exogenous variable formatted like |
exogen_all_col |
vector of considered columns of |
nscenario |
number of generated scenarios for daily maximum and minimum temperature |
seed |
seed for stochastic random generation see |
noise |
stochastic noise to add for variabile generation. Default is |
A list of the following variables:
input
list of variables returned by setComprehensiveTemperatureGeneratorParameters
var
varest object containing the used VAR model (if useVAR is true), NULL
(otherwise)
output
list variables returned by generateTemperatureTimeseries
(i.e. generated timeseries)
It pre-processes series and generates multi-site temperature fields by using setComprehensiveTemperatureGeneratorParameters
,getVARmodel
and generateTemperatureTimeseries
. Detailed examples can be viewed of this function in this presentation.
Emanuele Cordano, Emanuele Eccel
setComprehensiveTemperatureGeneratorParameters
, generateTemperatureTimeseries
,generateTemperatureTimeseries
,splineInterpolateMonthlytoDailyforSeveralYears
.
data(trentino) set.seed(1222) # set the seed for random generations! year_min <- 1961 year_max <- 1990 year_min_sim <- 1982 year_max_sim <- 1983 n_GPCA_iter <- 5 n_GPCA_iteration_residuals <- 5 p <- 1 vstation <- c("B2440","B6130","B8570","B9100","LAVIO","POLSA","SMICH","T0001", "T0010","T0014","T0018","T0032","T0064","T0083","T0090","T0092", "T0094","T0099","T0102","T0110","T0129","T0139","T0147","T0149", "T0152","T0157","T0168","T0179","T0189","T0193","T0204","T0210", "T0211","T0327","T0367","T0373") ## Not Run: the call to ComprehensiveTemperatureGenerator may elapse ## too long time (more than 5 eseconds) and is not executed by CRAN check. ## Please uncomment the following line to run the example on your own PC. # generation00 <-ComprehensiveTemperatureGenerator(station=vstation[16], # Tx_all=TEMPERATURE_MAX,Tn_all=TEMPERATURE_MIN,year_min=year_min,year_max=year_max, # p=p,n_GPCA_iteration=n_GPCA_iter,n_GPCA_iteration_residuals=n_GPCA_iteration_residuals, # sample="monthly",year_min_sim=year_min_sim,year_max_sim=year_max_sim)
data(trentino) set.seed(1222) # set the seed for random generations! year_min <- 1961 year_max <- 1990 year_min_sim <- 1982 year_max_sim <- 1983 n_GPCA_iter <- 5 n_GPCA_iteration_residuals <- 5 p <- 1 vstation <- c("B2440","B6130","B8570","B9100","LAVIO","POLSA","SMICH","T0001", "T0010","T0014","T0018","T0032","T0064","T0083","T0090","T0092", "T0094","T0099","T0102","T0110","T0129","T0139","T0147","T0149", "T0152","T0157","T0168","T0179","T0189","T0193","T0204","T0210", "T0211","T0327","T0367","T0373") ## Not Run: the call to ComprehensiveTemperatureGenerator may elapse ## too long time (more than 5 eseconds) and is not executed by CRAN check. ## Please uncomment the following line to run the example on your own PC. # generation00 <-ComprehensiveTemperatureGenerator(station=vstation[16], # Tx_all=TEMPERATURE_MAX,Tn_all=TEMPERATURE_MIN,year_min=year_min,year_max=year_max, # p=p,n_GPCA_iteration=n_GPCA_iter,n_GPCA_iteration_residuals=n_GPCA_iteration_residuals, # sample="monthly",year_min_sim=year_min_sim,year_max_sim=year_max_sim)
Calculates the continuity ratio of a set of precipitation measured or generated data in several sites as defined by Wilks, 1998 (see reference link)
continuity_ratio(data, lag = 0, valmin = 0.5)
continuity_ratio(data, lag = 0, valmin = 0.5)
data |
containing daily precipitation time series for several gauges (one gauge time series per column) |
lag |
numeric lag (expressed as number of days) used for computation for "cross" continuity ratio and joint probability of prercipitation (no)occurrence. |
valmin |
threshold precipitation value [mm] for wet/dry day indicator.
If precipitation is lower than |
A list containing the following matrices:
continuity_ratio
: lag
-day lagged continuity ratio ,
occurrence
: joint probability of lag
-day lagged precipitation occurrence
nooccurrence
: joint probability of lag
-day lagged no precipitation occurrence.
nooccurrence_occurrence
: joint probability of lag
-day lagged no precipitation and precipitation occurrence respectively.
occurrence_nooccurrence
: joint probability of lag
-day lagged precipitation and no precipitation occurrence respectively.
probability_continuity_ratio
: lag
-day lagged ratio about precipitation probability contitioned to no precipitation/preciitation occurrence in the other site
If lag==0
the function returns the continuity ratio and joint probability as described by Wilks, 1998. Otherwise the precipitation values for each couple of rain gauges are taken with lag
-day lag.
see the following URL references: http://onlinelibrary.wiley.com/doi/10.1002/joc.2305/abstract and http://www.sciencedirect.com/science/article/pii/S0022169498001863
data(trentino) year_min <- 1961 year_max <- 1990 origin <- paste(year_min,1,1,sep="-") period <- PRECIPITATION$year>=year_min & PRECIPITATION$year<=year_max station <- names(PRECIPITATION)[!(names(PRECIPITATION) %in% c("day","month","year"))] prec_mes <- PRECIPITATION[period,station] ## removing nonworking stations (e.g. time series with NA) accepted <- array(TRUE,length(names(prec_mes))) names(accepted) <- names(prec_mes) for (it in names(prec_mes)) { accepted[it] <- (length(which(!is.na(prec_mes[,it])))==length(prec_mes[,it])) } prec_mes <- prec_mes[,accepted] ## the dateset is reduced!!! prec_mes <- prec_mes[,1:2] continuity_ratio <-continuity_ratio(data=prec_mes,lag=0,valmin=0.5) continuity_ratio1 <-continuity_ratio(data=prec_mes,lag=-1,valmin=0.5)
data(trentino) year_min <- 1961 year_max <- 1990 origin <- paste(year_min,1,1,sep="-") period <- PRECIPITATION$year>=year_min & PRECIPITATION$year<=year_max station <- names(PRECIPITATION)[!(names(PRECIPITATION) %in% c("day","month","year"))] prec_mes <- PRECIPITATION[period,station] ## removing nonworking stations (e.g. time series with NA) accepted <- array(TRUE,length(names(prec_mes))) names(accepted) <- names(prec_mes) for (it in names(prec_mes)) { accepted[it] <- (length(which(!is.na(prec_mes[,it])))==length(prec_mes[,it])) } prec_mes <- prec_mes[,accepted] ## the dateset is reduced!!! prec_mes <- prec_mes[,1:2] continuity_ratio <-continuity_ratio(data=prec_mes,lag=0,valmin=0.5) continuity_ratio1 <-continuity_ratio(data=prec_mes,lag=-1,valmin=0.5)
data
counts NAs in each row of data
countNAs(data)
countNAs(data)
data |
a data input matrix @export |
the vector with numbers of NA values for each data
column
x
Calculates the covariance matrix of the normally standardized variables obtained from the columns of x
covariance( x, data = x, cpf = NULL, mean = 0, sd = 1, step = NULL, prec = 10^-4, use = "pairwise.complete.obs", type = 3, extremes = TRUE, sample = NULL, origin_x = NULL, origin_data = origin_x )
covariance( x, data = x, cpf = NULL, mean = 0, sd = 1, step = NULL, prec = 10^-4, use = "pairwise.complete.obs", type = 3, extremes = TRUE, sample = NULL, origin_x = NULL, origin_data = origin_x )
x |
variable |
data |
a sample of data on which a non-parametric pghjjrobability distribution is estimated |
cpf |
cumulative probability distribution. If |
mean |
mean (expected value) of the normalized random variable. Default is 0. |
sd |
standard deviation of the normalized random variable. Default is 1. |
step |
vector of values in which step discontinuities of the cumulative probability function occur. Default is |
prec |
amplitude of the neighbourhood of the step discontinuities where cumulative probability function is treated as non continuous. |
use |
see |
type |
see |
extremes |
logical variable.
If
where |
sample |
information about sample or probability distribution. Default is |
origin_x |
date corresponding to the first row of |
origin_data |
date corresponding to the first row of |
a matrix with the normalized variable or its inverse
Emanuele Cordano, Emanuele Eccel
normalizeGaussian_severalstations
,normalizeGaussian
@note It applies normalizeGaussian_severalstations
to x
and data
and then calculates the covariances among the column.
See the R code for further details
Extracts the elevation of a meteorological station expressed in meters above a reference (sea level)
ElevationOf(name, station_names, elevation)
ElevationOf(name, station_names, elevation)
name |
character ID of the station |
station_names |
vector of the IDs (characters) of the considered meteorological stations. An example is |
elevation |
vector of the elevation of the considered meteorological stations. An example is |
the elevation given the vectors of station IDs and the respective elevations
data(trentino) ElevationOf("T0099",station_names=STATION_NAMES,elevation=ELEVATION)
data(trentino) ElevationOf("T0099",station_names=STATION_NAMES,elevation=ELEVATION)
Extracts the rows of a matrix corresponding to the requested days (expressed as dates YYYY-MM-DD) given the date (origin) of the first row
extractdays( data = array(1:ndim_max, dim = c(ndim_max, 1)), ndim_max = 1e+05, when = "1990-1-1", origin = "1961-1-1", nday = 1 )
extractdays( data = array(1:ndim_max, dim = c(ndim_max, 1)), ndim_max = 1e+05, when = "1990-1-1", origin = "1961-1-1", nday = 1 )
data |
an input data matrix where each row corresponds to a daily record |
ndim_max |
maximum (integer) number of rows in |
when |
desired dates for which the data are requested |
origin |
date corresponding to the first row of |
nday |
(optional) number of days since |
a matrix containing the requested rows
It uses julian
extractdays()
extractdays()
Extracts the rows of a matrix corresponding to requested months of a year given the date (origin) of the first row
extractmonths( data = array(1:ndim_max, dim = c(ndim_max, 1)), ndim_max = 1e+05, when = c("Dec", "Jan", "Feb"), year = NULL, origin = "1961-1-1" )
extractmonths( data = array(1:ndim_max, dim = c(ndim_max, 1)), ndim_max = 1e+05, when = c("Dec", "Jan", "Feb"), year = NULL, origin = "1961-1-1" )
data |
an input data matrix where each row corresponds to a daily record |
ndim_max |
maximum (integer) number of rows in |
when |
character vactor of months for which the data are required.
It must be a subset of |
year |
year(s) when data must be extracted |
origin |
date corresponding to the first row of |
a matrix containing the requested rows
Emanuele Cordano, Emanuele Eccel
extractmonths() data(trentino) dates <- sprintf("%02d-%02d-%02d",TEMPERATURE_MAX$year,TEMPERATURE_MAX$month,TEMPERATURE_MAX$day) origin <- dates[1] out <- extractmonths(data=TEMPERATURE_MAX,origin=origin)
extractmonths() data(trentino) dates <- sprintf("%02d-%02d-%02d",TEMPERATURE_MAX$year,TEMPERATURE_MAX$month,TEMPERATURE_MAX$day) origin <- dates[1] out <- extractmonths(data=TEMPERATURE_MAX,origin=origin)
generateTemperatureTimeseries
functionExtracts generated time series of Daily Minimum Temperature from a random multi-realization obtained by generateTemperatureTimeseries
function
extractTnFromAnomalies(res_multigen, std, SplineAdv)
extractTnFromAnomalies(res_multigen, std, SplineAdv)
res_multigen |
matrix containing standardized values of daily temperature as returned by |
std |
vector containing standard deviation for each minimun temperature anomalies |
SplineAdv |
matrix containing the averaged daily values of minimum temperature obtained by a spline interpolation of the monthly climate |
a matrix with generated minimum temperature
Emanuele Cordano, Emanuele Eccel
generateTemperatureTimeseries
functionExtracts generated time series of Daily Maximum Temperature from a random multi-realization obtained by generateTemperatureTimeseries
function
extractTxFromAnomalies(res_multigen, std, SplineAdv)
extractTxFromAnomalies(res_multigen, std, SplineAdv)
res_multigen |
matrix containing standardized values of daily temperature as returned by |
std |
vector containing standard deviation for each maximum temperature anomalies |
SplineAdv |
matrix containing the averaged values of maximum temperature obtained by a spline interpolation of monthly climate |
a matrix with generated maximum temperature
Emanuele Cordano, Emanuele Eccel
year_min
and year_max
for the stations listed in station
Extracts the elements of a data frame corresponding to a period between year_min
and year_max
for the stations listed in station
extractyears( data, year_min = 1961, year_max = 1990, station = c("T0001", "T0014", "T0129") )
extractyears( data, year_min = 1961, year_max = 1990, station = c("T0001", "T0014", "T0129") )
data |
a dataframe containing daily data. |
year_min |
start year |
year_max |
end year |
station |
character vector of the IDs of the station where the data are required |
a matrix containing the requested daily data where each day corresponds to a row and each station corresponds to a column
The input data frame data
must have the following fields: year,month,day,variables_ID1,variables_ID2,...
where the fields ,variables_ID1,variables_ID2,...
contain the daily variables referred to the respective stations and the field names are replaced with the respective station ID.
Finds the date corresponding a row index of a matrix given the date (origin) of the first row
findDate( k, origin = "1961-1-1", data.frame = TRUE, decimal = FALSE, character = FALSE )
findDate( k, origin = "1961-1-1", data.frame = TRUE, decimal = FALSE, character = FALSE )
k |
integer or decimal value corresponding to number of days since |
origin |
origin date. See also |
data.frame |
logical variable. If |
decimal |
logical variable. If |
character |
logical variable. It is used if |
the date(s) corresponding to k
under different formats
It uses functions of time
package. It works like an inverse functions of extractdays
.
If k
is a vector, the function returns several dates for each element of k
findDate <- findDate(100,origin="1961-1-1",data.frame=FALSE,character=TRUE)
findDate <- findDate(100,origin="1961-1-1",data.frame=FALSE,character=TRUE)
Forecasts the expected value of a VAR realization given the prievious one
forecastEV(var, xprev = NULL, exogen = NULL)
forecastEV(var, xprev = NULL, exogen = NULL)
var |
A VAR model represented by a |
xprev |
previous status of the random variable |
exogen |
vector containing the values of the "exogen" variables (predictor) for the generation |
a vector of values
@export
Forecasts the residual value of a VAR realization given the white noise covariance matrix
forecastResidual(var, xprev = NULL, B = NULL)
forecastResidual(var, xprev = NULL, B = NULL)
var |
A VAR model represented by a |
xprev |
previous status of the random variable, in this case the "current instant"white-noise". Default is |
B |
matrix of coefficients for the vectorial white-noise component |
a vector of values
Emanuele Cordano, Emanuele Eccel
forecastEV
,NewVAReventRealization
newVARmultieventRealization
. This function is called by ComprehensiveTemperatureGenerator
.Returns time series of Daily Maximum and Minimum with a random multi-realization obtained by using newVARmultieventRealization
. This function is called by ComprehensiveTemperatureGenerator
.
generateTemperatureTimeseries( std_tn, std_tx, SplineTx, SplineTn, SplineTm, SplineDeltaT, std_tm, var = NULL, exogen = NULL, normalize = TRUE, type = 3, extremes = TRUE, sample = NULL, option = 1, original_data, origin_x = NULL, origin_data = NULL, noise = NULL )
generateTemperatureTimeseries( std_tn, std_tx, SplineTx, SplineTn, SplineTm, SplineDeltaT, std_tm, var = NULL, exogen = NULL, normalize = TRUE, type = 3, extremes = TRUE, sample = NULL, option = 1, original_data, origin_x = NULL, origin_data = NULL, noise = NULL )
std_tn |
vector containing standard deviation of daily minimum temperature anomalies. |
std_tx |
vector containing standard deviation of daily maximum temperature anomalies. |
SplineTx |
matrix containing the averaged daily maximum temperature obtained by a spline interpolation of monthly means . |
SplineTn |
matrix containing the averaged daily minimum temperature obtained by a spline interpolation of monthly means . |
SplineTm |
matrix containing the averaged daily "mean" temperature obtained by a spline interpolation of monthly means . |
SplineDeltaT |
matrix containing the rescaled averaged daily temperature range obtained by a spline interpolation of monthly means. |
std_tm |
vector containing standard deviation of daily "mean" temperature anomalies. |
var |
A VAR model represented by a |
exogen |
see |
normalize |
logical variable If |
type |
see |
sample , origin_x , origin_data , extremes
|
|
option |
integer value. If 1, the generator works with minimum and maximum temperature, if 2 (Default) it works with th average value between maximum and minimum temparature and the respective daily Thermal Range. |
original_data |
matrix containing the measured standardized temperature anomalies |
noise |
stochastic noise to add for variabile generation. Default is |
This function returns a list of the following variables:
res_multigen
matrix containing standardized values of daily maximum and minimum temperature anomalies
Tx_spline
matrix containing climatic "spline-interpolated" daily maximum temperature
Tn_spine
matrix containing climatic "spline-interpolated" daily minimum temperature
Tx_gen
matrix containing generated daily maximum daily temperature ()
Tn_gen
matrix containing generated daily minimum daily temperature ()
Tm_gen
matrix containing generated "mean" daily temperature defined as
DeltaT_gen
matrix containing generated daily thermal range defined as
See the R code for further details
Emanuele Cordano, Emanuele Eccel
newVARmultieventRealization
,normalizeGaussian_severalstations
year_min
and year_max
for stations listed in station
Calculates the daily means of a range of days around each date of a data frame corresponding to a period between year_min
and year_max
for stations listed in station
getDailyMean( data, year_min = 1961, year_max = 1990, station = c("T0001", "T0010"), origin = "1961-1-1", lag = 5 )
getDailyMean( data, year_min = 1961, year_max = 1990, station = c("T0001", "T0010"), origin = "1961-1-1", lag = 5 )
data |
a data frame containing daily data. |
year_min |
start year |
year_max |
end year |
station |
character vector of the IDs of the station where the data are requested |
origin |
origin date of time-series |
lag |
lag (number of days) on which daily mean is calculated. The mean is calculated considereing |
a matrix containing the requested daily mean data where each day corresponds to a row and each station corresponds to a column
The input data frame data
must have the following fields: year,month,day,variables_ID1,variables_ID2,...
where the fields ,variables_ID1,variables_ID2,...
contain the daily variables referred to the respective stations and the field names are replaced with the respective station ID.
Emanuele Cordano, Emanuele Eccel
year_min
and year_max
for stations listed in station
@author Emanuele Cordano, Emanuele Eccel
getMonthlyMean( data, year_min = 1961, year_max = 1990, station = names(data), no_date = FALSE, origin = "1961-1-1", yearly = FALSE )
getMonthlyMean( data, year_min = 1961, year_max = 1990, station = names(data), no_date = FALSE, origin = "1961-1-1", yearly = FALSE )
data |
a dataframe containing daily data. |
year_min |
start year |
year_max |
end year |
station |
character vector of the IDs of the station where the data are requested |
no_date |
logical value if |
origin |
date corresponding to the first row |
yearly |
logical value. If |
a matrix containing the requested monthly means where each month corresponds to a row and each station corresponds to a column or a list of such matrices in case the monthly mean values are calculated separately for each year (if yearly
is TRUE
)
The input data frame data
must have the following fields: year,month,day,variables_ID1,variables_ID2,...
where the fields ,variables_ID1,variables_ID2,...
contain the daily variables referred to the respective stations and the field names are replaced with the respective station ID. In case yearly
is TRUE
the returned output is a list of matrices whose names are the corresponding year.
vars
packageEither creates a VAR model or chooses a VAR model by using VAR or VARselect commands of vars
package
getVARmodel( data, suffix = c("_Tx", "_Tn"), sep = "", p = 1, type = "none", season = NULL, exogen = NULL, lag.max = NULL, ic = "AIC", activateVARselect = FALSE, na.rm = TRUE, n_GPCA_iteration = 0, n_GPCA_iteration_residuals = n_GPCA_iteration, extremes = TRUE )
getVARmodel( data, suffix = c("_Tx", "_Tn"), sep = "", p = 1, type = "none", season = NULL, exogen = NULL, lag.max = NULL, ic = "AIC", activateVARselect = FALSE, na.rm = TRUE, n_GPCA_iteration = 0, n_GPCA_iteration_residuals = n_GPCA_iteration, extremes = TRUE )
data |
see |
suffix |
see |
sep |
separator element. See |
p |
lag considered for the auto-regression see |
type |
see |
season |
see |
exogen |
see |
lag.max |
see |
ic |
see |
activateVARselect |
logical variables. If |
na.rm |
logical variables. If |
n_GPCA_iteration |
number of iterations of Gaussianization process for data. Default is 0 (no Gaussianization) |
n_GPCA_iteration_residuals |
number of iterations of Gaussianization process for data. Default is 0 (no Gaussianization) |
extremes |
a varest2
or GPCAvarest2
object representing a VAR model or a GPCA-varest
object which also contains the GPCA transformation parameters
It inherits input parameters of VAR
, VARselect
and addsuffixes
. The variable data
contains the measured data on which the vector auto-regressive models is estimated.
It is a matrix where each row is a realization of the vector random variable.
In some application of this package, the random variables may be the daily maximum and minimum temperature anomalies for different stations.
Often the the columns of data
are called with the IDs of the stations whithout specifying the type of variable (e.g. minimun or maximum temperature anomalies).
This means that two or more columns may have the same name. Therefore the function addsuffixes
, which is called from this function, adds suitable suffixes to the column names.
Emanuele Cordano, Emanuele Eccel
GPCA_iteration
)This function makes a Gaussianization procedure based on PCA iteration ( see GPCA_iteration
)
GPCA(x_prev, n = 30, extremes = TRUE)
GPCA(x_prev, n = 30, extremes = TRUE)
x_prev |
previous set of the random variable |
n |
number of reiterations |
extremes |
A GPCA-class
S3 object returned by GPCA_iteration
at each iteration
and the final results of the G-PCA procedure (matrix final_results
)
This function re-iterates the equation (1) of "PCA Gaussianization for One-Class Remote Sensing Image" by V. Laparra et al., http://www.uv.es/lapeva/papers/SPIE09_one_class.pdf,http://www.uv.es/vista/vistavalencia/papers/SPIE_09_Gaussianization_presentation.pdf
Emanuele Cordano
GPCA
,GPCA_iteration
,inv_GPCA_iteration
,inv_GPCA
,GPCA-class
for 'GPCA' S3 class
library(RMAWGEN) set.seed(1222) nIterations <- 30 N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA(df,n=nIterations,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA(dfn,n=nIterations,extremes=TRUE)
library(RMAWGEN) set.seed(1222) nIterations <- 30 N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA(df,n=nIterations,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA(dfn,n=nIterations,extremes=TRUE)
This function makes an iteration of PCA-Gaussianization process
GPCA_iteration(x_prev, extremes = TRUE)
GPCA_iteration(x_prev, extremes = TRUE)
x_prev |
previous set of random variable |
extremes |
A GPCA_iteration
S3 object which contains the following objects:
x_prev
Previous set of random variable, x_prev
input variable
x_gauss_prev
Marginal Gaussianization of x_prev
obtained through normalizeGaussian_severalstations
B_prev
rotation matrix (i. e. eigenvector matrix of the covariance matrix of x_gauss_prev
x_next
results obtained by multiplying B_prev
by x_gauss_prev
(see equation 1 of the reference)
This function is based on equation (1) of "PCA Gaussianization for One-Class Remote Sensing Image" by V. Laparra et al., https://www.uv.es/lapeva/papers/SPIE09_one_class.pdf and http://ieeexplore.ieee.org/document/5413808/
Emanuele Cordano
GPCA
,GPCA_iteration
,inv_GPCA_iteration
,inv_GPCA
library(RMAWGEN) set.seed(1222) N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA_iteration(df,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA_iteration(dfn,extremes=TRUE)
library(RMAWGEN) set.seed(1222) N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA_iteration(df,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA_iteration(dfn,extremes=TRUE)
GPCA
S3 class returned by GPCA
GPCA_iteration
subsequent GPCA iterations
final_results
data.frame or matrix of the "gaussianized" data
Formal definition with setOldClass
for the S3 class GPCA
Emanuele Cordano
showClass("GPCA")
showClass("GPCA")
GPCAiteration
S3 class returned by GPCA_iteration
x_prev
Previous set of random variable, x_prev
input variable of GPCA_iteration
x_gauss_prev
Marginal Gaussianization of x_prev
obtained through normalizeGaussian_severalstations
B_prev
rotation matrix (i. e. eigenvector matrix of the covariance matrix of x_gauss_prev
)
x_next
results obtained by multiplying B_prev
by x_gauss_prev
(see equation 1 of the reference in GPCA_iteration
)
Formal definition with setOldClass
for the S3 class GPCAiteration
Emanuele Cordano
showClass("GPCAiteration")
showClass("GPCAiteration")
This class inherits varest2
and contains all information about GPCA (GPCA
transformation.
GPCA_data
:A "GPCA"
S3 object containing the parameters of the Multi-variate Gaussianization of the time series, it is the result of GPCA
function applied to the input data of getVARmodel
GPCA_residuals
:A "GPCA"
S3 object containing the parameters of the Multi-variate Gaussianization of the residuals of the VAR model contained in the VAR
slot; it is NULL
if no Gaussiatization of residuals is applied.
Object of class "list"
VAR
:S3 Object of class "varest"
#' @note A GPCAvarest2
object can be created by new("GPCAvarest2", ...)
or returned by the function getVARmodel
Emanuele Cordano
showClass("GPCAvarest2")
showClass("GPCAvarest2")
inv_GPCA_iteration
This function makes an inverse Gaussianization procedure besad on PCA iteration ( see inv_GPCA_iteration
inv_GPCA(x = NULL, GPCA_param, type = 3, extremes = TRUE)
inv_GPCA(x = NULL, GPCA_param, type = 3, extremes = TRUE)
x |
gaussian random variable to transform |
GPCA_param |
|
type |
|
extremes |
the non-Gaussian random variable
This function re-iterates the inverse of equation (1) of "PCA Gaussianization for One-Class Remote Sensing Image" by V. Laparra et al., http://ieeexplore.ieee.org/document/5413808/
Emanuele Cordano
GPCA
,GPCA_iteration
,inv_GPCA_iteration
,inv_GPCA
library(RMAWGEN) set.seed(1222) nIterations <- 30 N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA(df,n=nIterations,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA(dfn,n=nIterations,extremes=TRUE) df_out <- inv_GPCA(GPCA_param=GPCA,extremes=TRUE) dfn_out <- inv_GPCA(GPCA_param=GPCAn,extremes=TRUE)
library(RMAWGEN) set.seed(1222) nIterations <- 30 N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA(df,n=nIterations,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA(dfn,n=nIterations,extremes=TRUE) df_out <- inv_GPCA(GPCA_param=GPCA,extremes=TRUE) dfn_out <- inv_GPCA(GPCA_param=GPCAn,extremes=TRUE)
This function makes an inverse iteration of PCA-Gaussianization process
inv_GPCA_iteration( x = GPCA_iter_param$x_next, GPCA_iter_param, type = 3, extremes = TRUE )
inv_GPCA_iteration( x = GPCA_iter_param$x_next, GPCA_iter_param, type = 3, extremes = TRUE )
x |
matrix of gaussian random variale to transform |
GPCA_iter_param |
|
type |
|
extremes |
the non-Gaussian random variable
This function is based on the inverse of the equation (1) of "PCA Gaussianization for One-Class Remote Sensing Image" by V. Laparra et al., http://ieeexplore.ieee.org/document/5413808/
GPCA
,GPCA_iteration
,inv_GPCA_iteration
,inv_GPCA
,GPCA-class
for 'GPCA' S3 class
library(RMAWGEN) set.seed(1222) N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA_iteration(df,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA_iteration(dfn,extremes=TRUE) df_out <- inv_GPCA_iteration(GPCA_iter_param=GPCA,extremes=TRUE) dfn_out <- inv_GPCA_iteration(GPCA_iter_param=GPCAn,extremes=TRUE)
library(RMAWGEN) set.seed(1222) N <- 20 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) GPCA <- GPCA_iteration(df,extremes=TRUE) x <- rnorm(N) y <- x+rnorm(N) dfn <- data.frame(x=x,y=y) GPCAn <- GPCA_iteration(dfn,extremes=TRUE) df_out <- inv_GPCA_iteration(GPCA_iter_param=GPCA,extremes=TRUE) dfn_out <- inv_GPCA_iteration(GPCA_iter_param=GPCAn,extremes=TRUE)
setComprehensiveTemperatureGeneratorParameters
.Verifies if 'climate' represents the monthly climatology in one year, i.e 'climate' is monthly.climate type matrix whose rows represent months and each column represents a station. It is also used in setComprehensiveTemperatureGeneratorParameters
.
is.monthly.climate(climate, nstation = 3, nmonth = 12, verbose = TRUE)
is.monthly.climate(climate, nstation = 3, nmonth = 12, verbose = TRUE)
climate |
matrix containing the 'monthly climatology' data |
nstation |
number of variable measurement stations (columns of the matrix 'climate') |
nmonth |
number of months in one year (it can be different if climate is represented by seasonal avarages or others), Default is 12 (recommended). (it can be different if climate is represented by seasonal averages, in this case 4) |
verbose |
Prints output and warining messagrs only if is |
A logical variable if the matrix 'climate' is monthly.climate type
Emanuele Cordano, Emanuele Eccel
setComprehensiveTemperatureGeneratorParameters
months REPLACEMANT
months_f(x, ...)
months_f(x, ...)
x |
an object. See |
... |
arguments |
Generates a new realization of a VAR model
NewVAReventRealization(var, xprev, noise, exogen = NULL, B = NULL)
NewVAReventRealization(var, xprev, noise, exogen = NULL, B = NULL)
var |
A VAR model represented by a |
xprev |
previous status of the random variable |
noise |
uncorrelated or white noise (residual). Default is |
exogen |
vector containing the values of the "exogen" variables (predictor) for the generation |
B |
matrix of coefficients for the vectorial white-noise component |
a vector of values
Emanuele Cordano, Emanuele Eccel
Generates several realizations of a VAR model
newVARmultieventRealization( var, xprev = rnorm(var@VAR$K * var@VAR$p), exogen = NULL, nrealization = 10, B = t(chol(cov(residuals(var)))), extremes = TRUE, type = 3, noise = NULL )
newVARmultieventRealization( var, xprev = rnorm(var@VAR$K * var@VAR$p), exogen = NULL, nrealization = 10, B = t(chol(cov(residuals(var)))), extremes = TRUE, type = 3, noise = NULL )
var |
A VAR model represented by a |
xprev |
previous status of the random variable |
exogen |
matrix containing the values of the "exogen" variables (predictor) for the generation |
nrealization |
number of realization (e.g. days to simulate). If |
B |
matrix of coefficients for the vector white-noise component |
extremes , type
|
see |
noise |
stochastic noise to add for variabile generation. Default is |
a matrix of values
Emanuele Cordano, Emanuele Eccel
normality.test
method for varest2
objectnormality.test
method for varest2
object
normality_test(object, ...)
normality_test(object, ...)
object |
a |
... |
passed arguments |
x
extracted by a population represented by the sample data
or sample
to a normally-distributed variable with assigned mean and standard deviation or vice versa in case inverse
is TRUE
Converts a random variable x
extracted by a population represented by the sample data
or sample
to a normally-distributed variable with assigned mean and standard deviation or vice versa in case inverse
is TRUE
normalizeGaussian( x = 0, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, step = NULL, prec = 10^-4, type = 3, extremes = TRUE, sample = NULL )
normalizeGaussian( x = 0, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, step = NULL, prec = 10^-4, type = 3, extremes = TRUE, sample = NULL )
x |
value or vector of values to be converted |
data |
a sample of data on which a non-parametric probability distribution is estimated |
cpf |
cumulative probability distribution. If |
mean |
mean (expected value) of the normalized random variable. Default is 0. |
sd |
standard deviation of the normalized random variable. Default is 1. |
inverse |
logical value. If |
step |
vector of values in which step discontinuities of the cumulative probability function occur. Default is |
prec |
amplitude of the neighbourhood of the step discontinuities where cumulative probability function is treated as non-continuous. |
type |
see |
extremes |
logical variable.
If
where |
sample |
a character string or |
the normalized variable or its inverse
@note This function makes a Marginal Gaussianization. See the R code for further details
Emanuele Cordano, Emanuele Eccel
inverse
is TRUE
Converts precipitation values to "Gaussinized" normally-distributed values taking into account the probability of no precipitation occurrences. values
or vice versa in case inverse
is TRUE
normalizeGaussian_prec( x = 0, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, type = 3, extremes = TRUE, sample = NULL, qnull = 0, valmin = 1 )
normalizeGaussian_prec( x = 0, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, type = 3, extremes = TRUE, sample = NULL, qnull = 0, valmin = 1 )
x |
value or vector of values to be converted |
data |
a sample of data on which a non-parametric probability distribution is estimated |
cpf |
cumulative probability distribution. If |
mean |
mean (expected value) of the normalized random variable. Default is 0. |
sd |
standard deviation of the normalized random variable. Default is 1. |
inverse |
logical value. If |
type |
see |
extremes |
logical variable.
If
where |
sample |
a character string or |
qnull |
probability of no precipitation occurrence |
valmin |
minimum value of precipitation to consider a wet day |
the normalized variable or its inverse
In the version 1.2.5 of RMAWGEN This function is deprecated and not used.
Emanuele Cordano, Emanuele Eccel
library(RMAWGEN) NDATA <- 1000 occurrence <- as.logical(runif(NDATA)>0.5) prec <- rexp(NDATA,rate=1/3) prec[!occurrence] <- 0 valmin <- 0.5 #0.01 x <- normalizeGaussian_prec(x=prec,valmin=valmin) prec2 <- normalizeGaussian_prec(x=x,data=prec,valmin=valmin,inverse=TRUE) qqplot(prec,prec2) occurrence3 <- as.logical(runif(NDATA)>0.5) prec3 <- rexp(NDATA,rate=1/3) prec3[!occurrence3] <- 0 x3 <- normalizeGaussian_prec(x=prec3,valmin=valmin) qqplot(x,x3) abline(0,1)
library(RMAWGEN) NDATA <- 1000 occurrence <- as.logical(runif(NDATA)>0.5) prec <- rexp(NDATA,rate=1/3) prec[!occurrence] <- 0 valmin <- 0.5 #0.01 x <- normalizeGaussian_prec(x=prec,valmin=valmin) prec2 <- normalizeGaussian_prec(x=x,data=prec,valmin=valmin,inverse=TRUE) qqplot(prec,prec2) occurrence3 <- as.logical(runif(NDATA)>0.5) prec3 <- rexp(NDATA,rate=1/3) prec3[!occurrence3] <- 0 x3 <- normalizeGaussian_prec(x=prec3,valmin=valmin) qqplot(x,x3) abline(0,1)
x
random variable extracted by populations represented by the columns of data
respectively or sample
to a normally-distributed samples with assinged mean and standard deviation or vice versa in case inverse
is TRUE
Converts several samples x
random variable extracted by populations represented by the columns of data
respectively or sample
to a normally-distributed samples with assinged mean and standard deviation or vice versa in case inverse
is TRUE
normalizeGaussian_severalstations( x, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, step = NULL, prec = 10^-4, type = 3, extremes = TRUE, sample = NULL, origin_x = NULL, origin_data = origin_x )
normalizeGaussian_severalstations( x, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, step = NULL, prec = 10^-4, type = 3, extremes = TRUE, sample = NULL, origin_x = NULL, origin_data = origin_x )
x |
value to be converted |
data |
a sample of data on which a non-parametric probability distribution is estimated |
cpf |
cumulative probability distribution. If |
mean |
mean (expected value) of the normalized random variable. Default is 0. |
sd |
standard deviation of the normalized random variable. Default is 1. |
inverse |
logical value. If |
step |
vector of values in which step discontinuities of the cumulative probability function occur. Default is |
prec |
amplitude of the neighbourhood of the step discontinuities where cumulative probability function is treated as non-continuous. |
type |
see |
extremes |
logical variable.
If
where |
sample |
information on how to sample |
origin_x |
date corresponding to the first row of |
origin_data |
date corresponding to the first row of |
a matrix with the normalized variable or its inverse
It applies normalizeGaussian
for each column of x
and data
.
See the R code for further details
Emanuele Cordano, Emanuele Eccel
## Not run: library(RMAWGEN) set.seed(1234) N <- 30 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfg <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE,inverse=FALSE) dfi <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE,inverse=TRUE) N <- 365*2 origin <- "1981-01-01" x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfgm <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE, inverse=FALSE,origin_x=origin,origin_data=origin,sample="monthly") dfim <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE, inverse=TRUE,origin_x=origin,origin_data=origin,sample="monthly") ## Compatibility with 'lubridate' package library(lubridate) N <- 30 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfg <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE,inverse=FALSE) dfi <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE,inverse=TRUE) N <- 365*2 origin <- "1981-01-01" x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfgm <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE, inverse=FALSE,origin_x=origin,origin_data=origin,sample="monthly") dfim <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE, inverse=TRUE,origin_x=origin,origin_data=origin,sample="monthly") ## End(Not run)
## Not run: library(RMAWGEN) set.seed(1234) N <- 30 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfg <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE,inverse=FALSE) dfi <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE,inverse=TRUE) N <- 365*2 origin <- "1981-01-01" x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfgm <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE, inverse=FALSE,origin_x=origin,origin_data=origin,sample="monthly") dfim <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE, inverse=TRUE,origin_x=origin,origin_data=origin,sample="monthly") ## Compatibility with 'lubridate' package library(lubridate) N <- 30 x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfg <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE,inverse=FALSE) dfi <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE,inverse=TRUE) N <- 365*2 origin <- "1981-01-01" x <- rexp(N) y <- x+rnorm(N) df <- data.frame(x=x,y=y) dfgm <- normalizeGaussian_severalstations(df,data=df,extremes=TRUE, inverse=FALSE,origin_x=origin,origin_data=origin,sample="monthly") dfim <- normalizeGaussian_severalstations(dfg,data=df,extremes=TRUE, inverse=TRUE,origin_x=origin,origin_data=origin,sample="monthly") ## End(Not run)
x
random variable (daily precipitation values) extracted by populations represented by the columns of data
respectively or sample
to a normally-distributed samples with assinged mean and standard deviation or vice versa in case inverse
is TRUE
using the function normalizeGaussian_prec
DEPRECATED Converts several samples x
random variable (daily precipitation values) extracted by populations represented by the columns of data
respectively or sample
to a normally-distributed samples with assinged mean and standard deviation or vice versa in case inverse
is TRUE
using the function normalizeGaussian_prec
normalizeGaussian_severalstations_prec( x, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, qnull = NULL, valmin = 0.5, type = 3, extremes = TRUE, sample = NULL, origin_x = NULL, origin_data = NULL )
normalizeGaussian_severalstations_prec( x, data = x, cpf = NULL, mean = 0, sd = 1, inverse = FALSE, qnull = NULL, valmin = 0.5, type = 3, extremes = TRUE, sample = NULL, origin_x = NULL, origin_data = NULL )
x |
value to be converted |
data |
a sample of data on which a non-parametric probability distribution is estimated |
cpf |
cumulative probability distribution. If |
mean |
mean (expected value) of the normalized random variable. Default is 0. |
sd |
standard deviation of the normalized random variable. Default is 1. |
inverse |
logical value. If |
qnull |
probability of no precipitation occurrence. (It can be a matrix in case |
valmin |
minimum value of precipitation to consider a wet day |
type |
see |
extremes |
logical variable.
If
where |
sample |
information about sample or probability distribution. Default is |
origin_x |
date corresponding to the first row of |
origin_data |
date corresponding to the first row of |
a matrix or a data.frame with the normalized variable or its inverse
In the version 1.2.5 of RMAWGEN This function is deprecated and not used.
Emanuele Cordano, Emanuele Eccel
x
and y
It makes a plot by sampling (e.g. monthly) the variables x
and y
plot_sample( x, y = normalizeGaussian_severalstations(x = as.data.frame(x), data = as.data.frame(data), origin_x = origin_x, origin_data = origin_data, sample = sample, step = step, prec = prec)[, 1], xlim = range(x, na.rm = TRUE), legend_position = "topleft", ylim = range(y, na.rm = TRUE), pch = 1, col = 1, col_max = 0.9, col_min = 0.1, origin, sample = NULL, xhist = hist(x, breaks = breaks, plot = FALSE), yhist = hist(y, breaks = breaks, plot = FALSE), axes = FALSE, step = NULL, prec = 1e-04, breaks = 50, origin_x = origin, origin_data = origin, data = x, xlab = "", ylab = "", color = FALSE, gray = TRUE, sort = FALSE, valmin_x = valmin, valmin_y = valmin, valmin = -9999, abline = c(0, 1), ... )
plot_sample( x, y = normalizeGaussian_severalstations(x = as.data.frame(x), data = as.data.frame(data), origin_x = origin_x, origin_data = origin_data, sample = sample, step = step, prec = prec)[, 1], xlim = range(x, na.rm = TRUE), legend_position = "topleft", ylim = range(y, na.rm = TRUE), pch = 1, col = 1, col_max = 0.9, col_min = 0.1, origin, sample = NULL, xhist = hist(x, breaks = breaks, plot = FALSE), yhist = hist(y, breaks = breaks, plot = FALSE), axes = FALSE, step = NULL, prec = 1e-04, breaks = 50, origin_x = origin, origin_data = origin, data = x, xlab = "", ylab = "", color = FALSE, gray = TRUE, sort = FALSE, valmin_x = valmin, valmin_y = valmin, valmin = -9999, abline = c(0, 1), ... )
x |
vector of input data |
y |
vector of second input data. Default is |
xlim , ylim , xlab , ylab
|
see |
legend_position |
legend position. Default is |
pch |
integer single or multi values for |
col |
integer single or multi values for |
col_max |
maximum value for color scale to apply to |
col_min |
minimum value for color scale to apply to |
origin |
date of the first row of |
sample |
string character containg informatio how to sample |
xhist |
frequency histogram for |
yhist |
frequency histogram for |
axes |
see |
step , prec
|
|
breaks |
see |
origin_x |
see |
origin_data |
|
data |
|
color |
logical value. If |
gray |
logical value. If |
sort |
logical value. If |
valmin_x |
numerical threshold value over which the variable |
valmin_y |
numerical threshold value over which the variable |
valmin |
numerical threshold value for |
abline |
arguments for |
... |
see graphical parametes on @usage plot_sample(x, y = normalizeGaussian_severalstations(x = as.data.frame(x), data = as.data.frame(data), origin_x = origin_x, origin_data = origin_data, sample = sample, step = step, prec = prec)[, 1], xlim = range(x, na.rm = TRUE), legend_position = "topleft", ylim = range(y, na.rm = TRUE), pch = 1, col = 1, col_max = 0.9, col_min = 0.1, origin, sample = NULL, xhist = hist(x, breaks = breaks, plot = FALSE), yhist = hist(y, breaks = breaks, plot = FALSE), axes = FALSE, step = NULL, prec = 1e-04, breaks = 50, origin_x = origin, origin_data = origin, data = x, xlab = "", ylab = "", color = FALSE, gray = TRUE, sort = FALSE, valmin_x = valmin, valmin_y = valmin, valmin = -9999, abline = c(0, 1), ...) |
0 in case of success
It makes a plot betwee x
and y
and shows thair respective probibilty histograms.
If y
is missing, it is automatically calculated as one-dimensional Gaussianization of x
through the function normalizeGaussian_severalstations
.
plot.default
,extractmonths
, see normalizeGaussian_severalstations
## Not run: library(lubridate) data(trentino) plot_sample(x=TEMPERATURE_MIN$T0090,sample="monthly", origin="1958-1-1",axes=FALSE,xlab="Tn [ degC]", ylab="x") set.seed(123456) z <- rexp(10000,rate=0.5) x <- normalizeGaussian(x=z,data=z) plot_sample(x=z,xlab="z",ylab="x") ## End(Not run)
## Not run: library(lubridate) data(trentino) plot_sample(x=TEMPERATURE_MIN$T0090,sample="monthly", origin="1958-1-1",axes=FALSE,xlab="Tn [ degC]", ylab="x") set.seed(123456) z <- rexp(10000,rate=0.5) x <- normalizeGaussian(x=z,data=z) plot_sample(x=z,xlab="z",ylab="x") ## End(Not run)
Plots daily climatology through one year
plotDailyClimate( data, title = "Daily_Avereged_Temperture_in_one_year", origin = "1961-1-1", when = "1979-1-1", ylab = "Temperature [degC]", xlab = "Time [days]", nday = 365, bicolor = FALSE, col = "black", lwd = 1 )
plotDailyClimate( data, title = "Daily_Avereged_Temperture_in_one_year", origin = "1961-1-1", when = "1979-1-1", ylab = "Temperature [degC]", xlab = "Time [days]", nday = 365, bicolor = FALSE, col = "black", lwd = 1 )
data |
matrix whose columns contain daily-averaged climatic series of variables (e.g. maximum or minum daily averaged temperature obtained by spline interpolation of monthly climatology) |
title , xlab , ylab , col , lwd
|
see |
origin |
origin date corresponding to the first row of |
when |
start day for daily climatology plot |
nday |
number of days in one year. Default is 365. |
bicolor |
logical variable. If |
a matrix containing the plotted variables
Emanuele Cordano, Emanuele Eccel
@author Emanuele Cordano, Emanuele Eccel
PrecipitationEndDay(name, station_names, end_day)
PrecipitationEndDay(name, station_names, end_day)
name |
charcacter ID of the station |
station_names |
vector containing the IDs (characters) of the considered meteorological stations. An example is |
end_day |
vector containing the measurement end day. An example is |
the precipitation measurement end day given the vectors of station IDs and the precipitation measurement end days
data(trentino) PrecipitationEndDay("T0099",station_names=STATION_NAMES,end_day=PRECIPITATION_MEASUREMENT_END_DAY)
data(trentino) PrecipitationEndDay("T0099",station_names=STATION_NAMES,end_day=PRECIPITATION_MEASUREMENT_END_DAY)
@author Emanuele Cordano
PrecipitationStartDay(name, station_names, start_day)
PrecipitationStartDay(name, station_names, start_day)
name |
character ID of the station |
station_names |
vector containing the IDs (characters) of the considered meteorological stations. An example is |
start_day |
vector containing the precipitation measurement start day. An example is |
the precipitation measurement start day given the vectors of station IDs and the respective precipitation measurement start days
data(trentino) PrecipitationStartDay("T0099", station_names=STATION_NAMES, start_day=PRECIPITATION_MEASUREMENT_START_DAY)
data(trentino) PrecipitationStartDay("T0099", station_names=STATION_NAMES, start_day=PRECIPITATION_MEASUREMENT_START_DAY)
print
S3 method for GPCA
or GPCA_iteration
objectprint
S3 method for GPCA
or GPCA_iteration
object
## S3 method for class 'GPCA' print(x, rmin = 1, rmax = 4, cmin = rmin, cmax = rmax, ...) ## S3 method for class 'GPCAiteration' print(x, rmin = 1, rmax = 4, cmin = rmin, cmax = rmax, ...)
## S3 method for class 'GPCA' print(x, rmin = 1, rmax = 4, cmin = rmin, cmax = rmax, ...) ## S3 method for class 'GPCAiteration' print(x, rmin = 1, rmax = 4, cmin = rmin, cmax = rmax, ...)
x |
a |
rmin , rmax , cmin , cmax
|
maximum and minimum rows and columns to be printed |
... |
passed arguments |
It makes the Q-Q plots observed vs generated time series of daily maximum, minimum temperature and daily thermal range for a list of collected stochastic generations
qqplot_RMAWGEN_Tx( Tx_mes, Tx_gen, Tn_gen, Tn_mes, Tx_spline = NULL, Tn_spline = NULL, xlab = "observed", ylab = "simulated", when = 1:nrow(Tx_mes), main = names(Tx_gen), station, pdf = NULL, xlim = range(Tx_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1 ) qqplot_RMAWGEN_Tn( Tx_mes, Tx_gen, Tn_gen, Tn_mes, Tx_spline = NULL, Tn_spline = NULL, xlab = "observed", ylab = "simulated", when = 1:nrow(Tn_mes), main = names(Tn_gen), station, pdf = NULL, xlim = range(Tn_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1 ) qqplot_RMAWGEN_deltaT( Tx_mes, Tx_gen, Tn_gen, Tn_mes, xlab = "observed", ylab = "simulated", when = 1:nrow(Tx_mes), main = names(Tx_gen), station, pdf = NULL, xlim = range(Tx_mes - Tn_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1 ) qqplot_RMAWGEN_prec( prec_mes, prec_gen, xlab = "observed", ylab = "simulated", when = 1:nrow(prec_mes), main = names(prec_gen), station, pdf = NULL, xlim = range(prec_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1, lag = 1 )
qqplot_RMAWGEN_Tx( Tx_mes, Tx_gen, Tn_gen, Tn_mes, Tx_spline = NULL, Tn_spline = NULL, xlab = "observed", ylab = "simulated", when = 1:nrow(Tx_mes), main = names(Tx_gen), station, pdf = NULL, xlim = range(Tx_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1 ) qqplot_RMAWGEN_Tn( Tx_mes, Tx_gen, Tn_gen, Tn_mes, Tx_spline = NULL, Tn_spline = NULL, xlab = "observed", ylab = "simulated", when = 1:nrow(Tn_mes), main = names(Tn_gen), station, pdf = NULL, xlim = range(Tn_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1 ) qqplot_RMAWGEN_deltaT( Tx_mes, Tx_gen, Tn_gen, Tn_mes, xlab = "observed", ylab = "simulated", when = 1:nrow(Tx_mes), main = names(Tx_gen), station, pdf = NULL, xlim = range(Tx_mes - Tn_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1 ) qqplot_RMAWGEN_prec( prec_mes, prec_gen, xlab = "observed", ylab = "simulated", when = 1:nrow(prec_mes), main = names(prec_gen), station, pdf = NULL, xlim = range(prec_mes), ylim = xlim, cex = 0.4, cex.main = 1, cex.lab = 1, cex.axis = 1, lag = 1 )
Tx_mes |
data frame containing measured daily maximum temperature |
Tx_gen |
data frame containing generated daily maximum temperature |
Tn_gen |
data frame containing generated daily minimum temperature |
Tn_mes |
data frame containing measured daily minimum temperature |
Tx_spline |
data frame containing spline-interpolated daily maximum temperature. Default is |
Tn_spline |
data frame containing spline-interpolated daily minimum temperature Default is |
xlab , ylab
|
lables of |
when |
day indices on which the data frame are extracted for Q-Q plot. Default is |
main |
main titles for each plot. Default is |
station |
identification name (ID) of the station used for the Q-Q plot |
pdf |
name of pdf file if output is written in a pdf file |
xlim |
see |
ylim , cex , cex.main , cex.lab , cex.axis
|
|
prec_mes |
data frame containing measured daily precipitation (in millimeters) |
prec_gen |
data frame containing generated daily precipitation (in millimeters) |
lag |
lag (current index included) on whose value the precipitation addition is made. See |
Tx_gen
,Tn_gen
and main
must have an even number of elements.
Emanuele Cordano
lag
-lag moving cumulative addition of the values in the samples x,y,z
This function creates a Q-Q plot of the lag
-lag moving cumulative addition of the values in the samples x,y,z
qqplot.lagged( x = rnorm(1000), y = rnorm(1000), z = NULL, when = 1:length(x), lag = 1, pch = 1, ... )
qqplot.lagged( x = rnorm(1000), y = rnorm(1000), z = NULL, when = 1:length(x), lag = 1, pch = 1, ... )
x , y
|
samples. If |
z |
further samples organized as a list |
when |
(integer) inidices of |
lag |
lag (current index included) on whose value the addition is made. |
pch |
a vector of plotting characters or symbols: see |
... |
further arguments for |
the Q-Q plot
Makes a qqplot of measured and simulated data for several stations.
qqplotprecWGEN( measured, simulated, xlab = "simulated[mm]", ylab = "measured[mm]", title = "daily precipitation", station = NULL, diff = FALSE, quantile = 0 )
qqplotprecWGEN( measured, simulated, xlab = "simulated[mm]", ylab = "measured[mm]", title = "daily precipitation", station = NULL, diff = FALSE, quantile = 0 )
measured |
matrix containing measured data (each station corresponds to a column) |
simulated |
matrix containing respective generated data (each station corresponds to a column) |
xlab , ylab
|
|
title |
title |
station |
character vector containing IDs of analyzed stations. If |
diff , quantile
|
see |
0 in case of success
It uses qqplotWGEN
and makes a figure for each pair of columns from measured
and simulated
. See the R code for further details.
Emanuele Cordano, Emanuele Eccel
Makes four seasonal qqplots (winter, spring, summer and autumn) of measured and simulated data for several stations.
qqplotprecWGEN_seasonal( measured, simulated, origin = "1961-1-1", xlab = "simulated[mm]", ylab = "measured[mm]", title = "daily_precipitation", directorypdf, station = names(simulated) )
qqplotprecWGEN_seasonal( measured, simulated, origin = "1961-1-1", xlab = "simulated[mm]", ylab = "measured[mm]", title = "daily_precipitation", directorypdf, station = names(simulated) )
measured |
matrix containing measured data (each station corresponds to a column) |
simulated |
matrix containing respective generated data (each station corresponds to a column) |
origin |
first day of data, see |
xlab , ylab
|
|
title |
title |
directorypdf |
name of the directory (path included) where to seva the outputs |
station |
character vector containing IDs of analyzed stations. If |
0 in case of success
Uses qqplotprecWGEN
for each season of collected data and saves the output on pdf files. See the R code for further details.
Emanuele Cordano, Emanuele Eccel
Makes a qqplot of measured and simulated data for several stations.
qqplotTnTxWGEN( measured, simulated, xlab = "simulated[degC]", ylab = "measured[degC]", titles = c("Q-Qplot_An._Tx", "Q-Qplot_An._Tn"), station = NULL, diff = FALSE, quantile = 0 )
qqplotTnTxWGEN( measured, simulated, xlab = "simulated[degC]", ylab = "measured[degC]", titles = c("Q-Qplot_An._Tx", "Q-Qplot_An._Tn"), station = NULL, diff = FALSE, quantile = 0 )
measured |
matrix containing measured data (each station corresponds to a column) |
simulated |
matrix containing respective generated data (each station corresponds to a column) |
xlab , ylab
|
|
titles |
titles that will be added to |
station |
character vector containing IDs of analyzed station. If |
diff , quantile
|
see |
0 in case of success
It uses qqplotWGEN
and makes a figure for each pair of columns from measured
and simulated
. See the R code for further details.
Emanuele Cordano, Emanuele Eccel
Makes four seasonal qqplots (winter, spring, summer and autumn) of measured and simulated data for several stations.
qqplotTnTxWGEN_seasonal( measured, simulated, origin = "1961-1-1", xlab = "simulated[degC]", ylab = "measured[degC]", titles = c("Q-Qplot_An._Tx", "Q-Qplot_An._Tn"), directorypdf, station = NULL )
qqplotTnTxWGEN_seasonal( measured, simulated, origin = "1961-1-1", xlab = "simulated[degC]", ylab = "measured[degC]", titles = c("Q-Qplot_An._Tx", "Q-Qplot_An._Tn"), directorypdf, station = NULL )
measured |
matrix containing measured data (each station corresponds to a column) |
simulated |
matrix containing respective generated data (each station corresponds to a column) |
origin |
first day of data, see |
xlab , ylab
|
|
titles |
titles that will be added |
directorypdf |
name of the directory (path included) where to seva the outputs |
station |
character vector containing IDs of analyzed station. If |
0 in case of success
Uses qqplotTnTxWGEN
for each seasons of collected data and saves the output on pdf files. See the R code for further details.
Emanuele Cordano, Emanuele Eccel
val
Makes a qqplot and Wilcoxon test between the two columns of val
qqplotWGEN( val, xlab = "simulated", ylab = "measured", main = "title", ylim = c(min(val), max(val)), xlim = c(min(val), max(val)), diff = FALSE, quantile = 0 )
qqplotWGEN( val, xlab = "simulated", ylab = "measured", main = "title", ylim = c(min(val), max(val)), xlim = c(min(val), max(val)), diff = FALSE, quantile = 0 )
val |
a matrix with two columns containing the two samples to be compared |
xlab , ylab , main
|
see |
xlim , ylim
|
see |
diff |
logical variable, if |
quantile |
quantile value on which data samples in |
Wilcoxon test between the two columns of 'val'
Emanuele Cordano, Emanuele Eccel
Replaces each entry of the rows containing NA values with NA
removeNAs(data)
removeNAs(data)
data |
a matrix @author Emanuele Cordano, Emanuele Eccel |
the matrix data
with the modified rows of NA values
In getVARmodel
,
when using VAR
or VARselect
, all NAs will be removed
This function adjusts the monthly mean to a daily weather dataset (e. g. spline-interpolated temperature)
rescaling_monthly(data, val, origin = "1961-1-1")
rescaling_monthly(data, val, origin = "1961-1-1")
data |
data frame of wheather variables) |
val |
monthly means returned by |
origin |
character string containing the gregorian date of the first day of |
A data frame with data of data
rescaled with val
for each month
Emanuele Cordano
@export
residuals
S3 method for varest2
objectresiduals
S3 method for varest2
object
## S3 method for class 'varest2' residuals(object, squared = FALSE, ...)
## S3 method for class 'varest2' residuals(object, squared = FALSE, ...)
object |
a |
squared |
logical value. Default is |
... |
passed arguments |
residuals of object
as a data frame. In case squared=TRUE
, the squared residauls are returned, otherwise simple residuals are returned. The squared residuals can be useful in case of ARCH analysis.
Emanuele Cordano
serial.test
function for varest2
objectserial.test
function for varest2
object
serial_test(object, ...)
serial_test(object, ...)
object |
a |
... |
passed arguments |
ComprehensiveTemperatureGenerator
.Computes climatic and correlation information useful for creating an auto-regeressive random generation of maximum and minimun daily temparature. This function is called by ComprehensiveTemperatureGenerator
.
setComprehensiveTemperatureGeneratorParameters( station, Tx_all, Tn_all, mean_climate_Tn = NULL, mean_climate_Tx = NULL, Tx_spline = NULL, Tn_spline = NULL, year_max = 1990, year_min = 1961, leap = TRUE, nmonth = 12, verbose = FALSE, cpf = NULL, normalize = TRUE, sample = NULL, option = 2, yearly = FALSE )
setComprehensiveTemperatureGeneratorParameters( station, Tx_all, Tn_all, mean_climate_Tn = NULL, mean_climate_Tx = NULL, Tx_spline = NULL, Tn_spline = NULL, year_max = 1990, year_min = 1961, leap = TRUE, nmonth = 12, verbose = FALSE, cpf = NULL, normalize = TRUE, sample = NULL, option = 2, yearly = FALSE )
station |
character vector of the IDs of the considered meteorological stations |
Tx_all |
data frame containing daily maximum temperature of all meteorological station. See |
Tn_all |
data frame containing daily minimum temperature of all meteorological station. See |
mean_climate_Tn |
a matrix containing monthly mean minimum daily temperature for the considered station or an object as returned by |
mean_climate_Tx |
a matrix containing monthly mean maximum daily temperature for the considered station or an object as returned by |
Tx_spline |
daily timeseries (from the first day of |
Tn_spline |
daily timeseries (from the first day of |
year_max |
start year of the recorded (calibration) period |
year_min |
end year of the recorded (calibration) period |
leap |
logical variables. It is |
nmonth |
number of months in one year. Default is 12. |
verbose |
logical variable |
cpf |
|
normalize |
logical variable If |
sample |
|
option |
integer value. If 1, the generator works with minimum and maximum temperature, if 2 (default) it works with the average value between maximum and minimum temperature and the respective daily thermal range. |
yearly |
logical value. If |
This function creates and returns the following gloabal variables:
data_original
matrix containing normalized and standardized data (i.e. data_original
)
data_for_var
matrix returned from normalizeGaussian_severalstations
by processing data_original
if normalize
is TRUE
), otherwise it is equal to data_original
.
Tn_mes
matrix containing measured minimum daily temperature in the analyzed time period ( )
Tx_mes
matrix containing measured maximum daily temperature in the analyzed time period ( )
Tm_mes
matrix calculated as to
DeltaT_mes
matrix corresponding to
monthly_mean_Tn
matrix containing monthly means of minimum daily temperature for the considered station. It is calculated according to the input format is.monthly.climate
if saveMonthlyClimate
is TRUE
.
monthly_mean_Tx
matrix containing monthly means of maximum daily temperature for the considered station. It is calculated according to the input format is.monthly.climate
if saveMonthlyClimate
is TRUE
.
Tx_spline
matrix containing the averaged daily values of maximimum temperature obtained by a spline interpolation of the monthly climate monthly_mean_Tx
or mean_climate_Tx
using splineInterpolateMonthlytoDailyforSeveralYears
( )
Tn_spline
matrix containing the averaged daily values of minimun temperature obtained by a spline interpolation of the monthly climate monthly_mean_Tn
or mean_climate_Tn
using splineInterpolateMonthlytoDailyforSeveralYears
( )
SplineAdvTm
matrix calculated as
SplineAdvDeltaT
, matrix corresponding to
stdTn
vector containing the standard deviation of minimum temperature anomalies (
)
stdTx
vector containing the standard deviation of maximum temperature anomalies (
)
stdTm
vector containing the standard deviation of "mean" temperature anomalies (
)
Tn_mes_res
standard core (standardization) of obtained
by solving column by column the expression
Tx_mes_res
standard core (standardization) of obtained
by solving column-by-column the expression
Tm_mes_res
standard core (standardization) of obtained
by solving column-by-column the expression
DeltaT_mes_res
equal to DeltaT_mes
data_original
matrix obtained as cbind(Tx_mes_res,Tn_mes_res)
if option
==1, or cbind(Tm_mes_res,DeltaT_mes_res)
if option
==2
See the R code for further details.
Emanuele Cordano, Emanuele Eccel
splineInterpolateMonthlytoDailyforSeveralYears
,ComprehensiveTemperatureGenerator
spline
and preserving monthly mean valuesInterpolates monthly data to daily data using spline
and preserving monthly mean values
splineInterpolateMonthlytoDaily( nday = 365, val = as.matrix(cbind(1 * (0.5:11.5) * nday/12, 2 * (0.5:11.5) * nday/12)), origin = "1961-1-1", first_row = 1, last_row = nday, no_spline = FALSE, no_mean = FALSE )
splineInterpolateMonthlytoDaily( nday = 365, val = as.matrix(cbind(1 * (0.5:11.5) * nday/12, 2 * (0.5:11.5) * nday/12)), origin = "1961-1-1", first_row = 1, last_row = nday, no_spline = FALSE, no_mean = FALSE )
nday |
number of days on which the daily data is requested, e.g. number of days in one year |
val |
matrix containing monthly mean data |
origin |
date corresponding to the first row of the returned matrix |
first_row |
row corresponding the first day of time interval where montlhy mean conservation is applied |
last_row |
corresponding the last day of time interval where montlhy mean conservation is applied |
no_spline |
logical value. If |
no_mean |
logical value. Default is |
a matrix or data frame with interpolated daily data
Emanuele Cordano, Emanuele Eccel
spline
,splineInterpolateMonthlytoDailyforSeveralYears
splineInterpolateMonthlytoDaily
for several yearsInterpolates monthly data to daily data using splineInterpolateMonthlytoDaily
for several years
splineInterpolateMonthlytoDailyforSeveralYears( val, start_year = 2010, nyear = 1, leap = TRUE, offset = 2, no_spline = FALSE, yearly = FALSE )
splineInterpolateMonthlytoDailyforSeveralYears( val, start_year = 2010, nyear = 1, leap = TRUE, offset = 2, no_spline = FALSE, yearly = FALSE )
val |
matrix containing monthly mean data for one year |
start_year |
first year |
nyear |
number of years since |
leap |
logical variable If |
offset |
integer values. Default is 2. Number of years considered beyond the extremes in order to avoid edge errors |
no_spline |
logical value. If |
yearly |
logical value. If @return a matrix or data frame with interpolated daily data |
Emanuele Cordano, Emanuele Eccel
spline
,splineInterpolateMonthlytoDaily
Gets the last day in a temperature time series, expressed as decimal julian days since 1970-1-1 00:00 UTC
TemperatureEndDay(name, station_names, end_day)
TemperatureEndDay(name, station_names, end_day)
name |
character ID of the station |
station_names |
vector containing the IDs (characters) of the considered meteorological stations. An example is |
end_day |
vector containing the measurement end day. An example is |
the temperature measurement end day given the vectors of station IDs and the temperature measurement end days
Emanuele Cordano, Emanuele Eccel
data(trentino) TemperatureEndDay("T0099",station_names=STATION_NAMES,end_day=TEMPERATURE_MEASUREMENT_END_DAY)
data(trentino) TemperatureEndDay("T0099",station_names=STATION_NAMES,end_day=TEMPERATURE_MEASUREMENT_END_DAY)
@author Emanuele Cordano, Emanuele Eccel
TemperatureStartDay(name, station_names, start_day)
TemperatureStartDay(name, station_names, start_day)
name |
character ID of the station |
station_names |
vector containing the IDs (characters) of the considered meteorological stations. An example is |
start_day |
vector containing the temperature measurement start day. Default is @export |
the temperature measurement start day given the vectors of station IDs and the respective temperature measurement start days
@examples data(trentino) TemperatureStartDay("T0099",station_names=STATION_NAMES,start_day=TEMPERATURE_MEASUREMENT_START_DAY)
It contains the following variables:
TEMPERATURE_MIN
Data frame containing year
,month
, day
and daily minimum temperature in 59 stations in Trentino region
TEMPERATURE_MAX
Data frame containing year
,month
, day
and daily maximum temperature in 59 stations in Trentino region
PRECIPITATION
Data frame containing year
,month
, day
and daily precipitation in 59 stations in Trentino region
STATION_NAMES
Vector containing the names of the meteorological stations
ELEVATION
Vector containing the elevations of the meteorological stations respectively
STATION_LATLON
Matrix containing the latitude and longitude coordinates, respectively, of the meteorological stations
LOCATION
Vector containing the names of the location of each meteorological station
TEMPERATURE_MEASUREMENT_START_DAY
Vector containing the first days referred to midday (expressed as decimal julian day since 1970-1-1 00:00 UTC) of temperature measurement of each meteorological station
TEMPERATURE_MEASUREMENT_END_DAY
Vector containing the last days referred to midday (expressed as decimal julian day since 1-1-1970 00:00 UTC) of temperature measurement of each meteorological station
PRECIPITATION_MEASUREMENT_START_DAY
Vector containing the first days referred to midday (expressed as decimal julian day since 1-1-1970 00:00 UTC) of precipitation measurement of each meteorological station
PRECIPITATION_MEASUREMENT_END_DAY
Vector containing the last days referred to midday (expressed as decimal julian day since 1-1-1970) of precipitation measurement of each meteorological station
data(trentino)
data(trentino)
Data frames and vectors
This dataset stores all information about meteorological stations and instrumental timeseries. The user can easily use the package with his/her own data after replacing the values of such variables.
Original data are provided by Provincia Autonoma di Trento (https://www.meteotrentino.it/), Fondazione Edmund Mach (https://www.fmach.it), Provincia Autonama di Bolzano/Autome Provinz Bozen (http://www.provincia.bz.it/meteo), ARPA Lombardia (https://www.arpalombardia.it/), ARPA Veneto (https://www.arpa.veneto.it/previsioni/it/html/index.php).
This dataset is intended for research purposes only, being distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY.
VAR
function allowing to describe white-noise as VAR-(0) model (i. e. varest
objects)Modified version of VAR
function allowing to describe white-noise as VAR-(0) model (i. e. varest
objects)
VAR_mod( y, p = 1, type = c("const", "trend", "both", "none"), season = NULL, exogen = NULL, lag.max = NULL, ic = c("AIC", "HQ", "SC", "FPE") )
VAR_mod( y, p = 1, type = c("const", "trend", "both", "none"), season = NULL, exogen = NULL, lag.max = NULL, ic = c("AIC", "HQ", "SC", "FPE") )
y , p , type , season , exogen , lag.max , ic
|
see |
a Vector Auto-Regeressive model (VAR) as varest
object
varest
S3 class (formal definition) see VAR
The details of the class are reported on VAR
documentation in "vars" package
Formal definition with setOldClass
for the S3 class varest
Bernhard Pfaff
showClass("varest")
showClass("varest")
This class derives from a varest
S3 class which is a list of objects describing a Vectorial AutoRegressive Model (see VAR
)
VAR
:a varest
S3 object created by VAR
A varest2
object can be created by new("varest2", ...)
or returned by the function getVARmodel
Emanuele Cordano
showClass("varest2")
showClass("varest2")
Gets the toponym where a meteorological station is located
WhereIs(name, station_names, location)
WhereIs(name, station_names, location)
name |
character ID of the station |
station_names |
vector containing the IDs (characters) of the considered meteorological stations. An example is |
location |
vector containing the toponyms. An example is |
the location toponym given the vectors of station IDs and the respective location toponyms
Emanuele Cordano, Emanuele Eccel
data(trentino) WhereIs("T0099",station_names=STATION_NAMES,location=LOCATION)
data(trentino) WhereIs("T0099",station_names=STATION_NAMES,location=LOCATION)