View source: R/makeMultiGrid.R
makeMultiGrid | R Documentation |
Constructs a (possibly multimember) multigrid from different (multimember) grids. A multigrid can be considered as a “stack” of grids with similar spatiotemporal extents, useful to handle sets of predictors as a single block.
makeMultiGrid(..., spatial.tolerance = 0.001, skip.temporal.check = FALSE)
... |
Input grids to form the multigrid. These must be compatible in time and space (see details). For flexibility, they can be introduced as a list or directly as consecutive arguments. |
spatial.tolerance |
numeric. Coordinate differences smaller than |
skip.temporal.check |
Logical. Should temporal matching check of input fields be skipped?
Default to |
The function makes a number of checks in order to test the spatiotemporal compatibility of the input multi-member grids.
Regarding the temporal concordance, it can be skipped by setting the argument skip.temporal.check
to TRUE
.
It is implicitly assumed that all temporal data from the different
multimember grids correspond to the same time zone ("GMT"). The time zone itself is not important, as long as it is the
same across datasets, because temporal consistency is checked on a daily basis (not hourly), allowing the inclusion of
predictors with different verification times and temporal aggregations. For instance, instantaneous geopotential at 12:00
is compatible with mean daily surface temperature, always that both variables correspond to the same days. Different time
resolutions are not compatible and will return an error (for instance, 6-hourly data is incompatible with daily values,
because their respective time series for a given season have different lengths).
The spatial consistency of the input grids is also checked. In order to avoid possible errors from the user, the spatial
consistency (i.e., equal XY coordinates) of the input grids must be ensured before attempting the creation of the multigrid,
otherwise giving an error. This can be achieved either through the specification of the same 'lonLim' and 'latLim' argument
values when loading the grids, or using the interpGrid
interpolator in conjuntion with the getGrid
method.
Variable names
When a vertical level of the variable is defined (i.e., is not NA
),
the short name of the variable is modified to the standard format "var@level"
(e.g., "ta@850"
for air temperature at 850mb surface pressure level). This way, the same
variable at different vertical levels can be differentiated at a glance using getVarNames
,
for instance.
A (multimember) multigrid object encompassing the different input (multimember) grids
A multigrid can not be passed to the interpolator interpGrid
directly. Instead, the
multimember grids should be interpolated individually prior to multigrid construction.
J. Bedia
interpGrid
for spatial consistency of input grids.
Other downscaling.helpers:
convert2bin()
,
filterNA()
require(climate4R.datasets)
# Creation of a multigrid from three different grids:
data(NCEP_Iberia_ta850)
data(NCEP_Iberia_hus850)
data(NCEP_Iberia_psl)
# An example of different temporal aggregations, temporally compatible:
# sea-level pressure is a daily mean, while specific humidity and air temperature
# (850 mb surface isobaric pressure level) are instantaneous data verifying at 12:00 UTC:
# air temperature
mf <- makeMultiGrid(NCEP_Iberia_hus850, NCEP_Iberia_psl, NCEP_Iberia_ta850)
# The new object inherits the global attributes from the first grid, as it is assumed
# that all input grids come from the same data source:
attributes(mf)
# The data structure has now one additional dimension ("var"), along which the data arrays
# have been binded:
str(mf$Data)
# Example of multimember multigrid creation from several multimember grids:
# Load three different multimember grids with the same spatiotemporal ranges:
data("CFS_Iberia_tas")
data("CFS_Iberia_hus850")
data("CFS_Iberia_pr")
mm.mf <- makeMultiGrid(CFS_Iberia_tas, CFS_Iberia_hus850, CFS_Iberia_pr)
# Different fields should not be plotted together in the same plot directly, unless
# their units are compatible.
# subsetGrid and visualizeR::spatialPlot can be used to this aim, if needed.
# For instance:
tas <- subsetGrid(mm.mf, var = "tas")
require(visualizeR)
spatialPlot(climatology(tas), backdrop.theme = "coastline", rev.colors = TRUE)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.