Title: | Find, Download and Process MODIS Land Products Data |
---|---|
Description: | Allows automating the creation of time series of rasters derived from MODIS satellite land products data. It performs several typical preprocessing steps such as download, mosaicking, reprojecting and resizing data acquired on a specified time period. All processing parameters can be set using a user-friendly GUI. Users can select which layers of the original MODIS HDF files they want to process, which additional quality indicators should be extracted from aggregated MODIS quality assurance layers and, in the case of surface reflectance products, which spectral indexes should be computed from the original reflectance bands. For each output layer, outputs are saved as single-band raster files corresponding to each available acquisition date. Virtual files allowing access to the entire time series as a single file are also created. Command-line execution exploiting a previously saved processing options file is also possible, allowing users to automatically update time series related to a MODIS product whenever a new image is available. For additional documentation refer to the following article: Busetto and Ranghetti (2016) <doi:10.1016/j.cageo.2016.08.020>. |
Authors: | Lorenzo Busetto [aut] , Luigi Ranghetti [aut, cre] , Leah Wasser [rev] (Leah Wasser reviewed the package for rOpenSci, see https://github.com/ropensci/software-review/issues/184), Jeff Hanson [rev] (Jeff Hanson reviewed the package for rOpenSci, see https://github.com/ropensci/software-review/issues/184), Babak Naimi [ctb] (Babak Naimi wrote the function ModisDownload(), on which some MODIStsp internal functions are based) |
Maintainer: | Luigi Ranghetti <[email protected]> |
License: | GPL-3 |
Version: | 2.1.0 |
Built: | 2024-12-09 06:32:49 UTC |
Source: | https://github.com/ropensci/MODIStsp |
MODIStsp allows automating the creation of time series of rasters derived from MODIS Satellite Land Products data. It performs several typical preprocessing steps such as download, mosaicking, reprojection and resize of data acquired on a specified time period. All processing parameters can be set using a user-friendly GUI. Users can select which layers of the original MODIS HDF files they want to process, which additional Quality Indicators should be extracted from aggregated MODIS Quality Assurance layers and, in the case of Surface Reflectance products , which Spectral Indexes should be computed from the original reflectance bands. For each output layer, outputs are saved as single-band raster files corresponding to each available acquisition date. Virtual files allowing access to the entire time series as a single file are also created. Command-line execution exploiting a previously saved processing options file is also possible, allowing to automatically update time series related to a MODIS product whenever a new image is available.
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015-2017)
https://docs.ropensci.org/MODIStsp/
https://github.com/ropensci/MODIStsp
Helper function used to retrieve the bounding box of a specified spatial file
recognized by sf
or raster
: the function reads the extent using sf::st_bbox()
bbox_from_file(file_path, crs_out)
bbox_from_file(file_path, crs_out)
file_path |
|
crs_out |
( |
License: GPL 3.0
Lorenzo Busetto, phD (2017)
Luigi Ranghetti, phD (2017)
Accessory function used to see if all expected out files for the selected date are already present in the output folder. If all expected out files are already present, check_files is set to TRUE, and the date is skipped in MODIStsp_process.
check_files_existence( out_prod_folder, file_prefix, yy, DOY, bandnames, bandsel_orig_choice, indexes_bandnames, indexes_bandsel, quality_bandnames, quality_bandsel, out_format )
check_files_existence( out_prod_folder, file_prefix, yy, DOY, bandnames, bandsel_orig_choice, indexes_bandnames, indexes_bandsel, quality_bandnames, quality_bandsel, out_format )
out_prod_folder |
|
file_prefix |
|
yy |
|
DOY |
|
bandnames |
|
bandsel_orig_choice |
|
indexes_bandnames |
|
indexes_bandsel |
|
quality_bandnames |
|
quality_bandsel |
|
out_format |
|
check - logical = 1 if all expected output files are already existing
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
helper function used to check that the input projection (passed as UTM zone, EPSG code, WKT string) is a valid projection for MODIStsp.
check_projection(projection, abort = FALSE, verbose = TRUE) ## Default S3 method: check_projection(projection, abort = FALSE, verbose = TRUE) ## S3 method for class 'numeric' check_projection(projection, abort = FALSE, verbose = TRUE) ## S3 method for class 'character' check_projection(projection, abort = FALSE, verbose = TRUE) ## S3 method for class 'crs' check_projection(projection, abort = FALSE, verbose = TRUE)
check_projection(projection, abort = FALSE, verbose = TRUE) ## Default S3 method: check_projection(projection, abort = FALSE, verbose = TRUE) ## S3 method for class 'numeric' check_projection(projection, abort = FALSE, verbose = TRUE) ## S3 method for class 'character' check_projection(projection, abort = FALSE, verbose = TRUE) ## S3 method for class 'crs' check_projection(projection, abort = FALSE, verbose = TRUE)
projection |
|
abort |
|
verbose |
|
character
proj4string of the object or file
This function was forked from package sprawl
, version 0.3.0.
Lorenzo Busetto, phD (2017)
Luigi Ranghetti, phD (2017)
## Not run: check_projection("32632") check_projection("32631") check_projection(32633) check_projection(30, abort = FALSE) check_projection("example of invalid string", abort = FALSE) proj_wkt <- sf::st_as_text(sf::st_crs(32632)) check_projection(proj_wkt) ## End(Not run)
## Not run: check_projection("32632") check_projection("32631") check_projection(32633) check_projection(30, abort = FALSE) check_projection("example of invalid string", abort = FALSE) proj_wkt <- sf::st_as_text(sf::st_crs(32632)) check_projection(proj_wkt) ## End(Not run)
Accessory function to find the folders corresponding to the requested dates period within the full list retrieved by get_moddirs
get_mod_dates(dates, date_dirs)
get_mod_dates(dates, date_dirs)
dates |
2- element string array specifying start/end dates (yyyy.mm.dd) for which the http addresses of folders in lpdaac should be retrieved (e.g., c("2015.1.1", "2015.12.31) |
date_dirs |
data frame full list of folders in lpdaac archive for product of interest |
array of folder names containing data for the MODIS product acquired in the period specified by "dates"
License: GPL 3.0
Luigi Ranghetti, phD (2016)
Lorenzo Busetto, phD (2017)
Accessory function to get the full list of directories on the lpdaac http site containing data included in the time range selected for processing (modified after Barry Rowlingson function):
get_mod_dirs( http, download_server, user, password, yy, n_retries, gui, out_folder_mod )
get_mod_dirs( http, download_server, user, password, yy, n_retries, gui, out_folder_mod )
http |
|
download_server |
|
user |
|
password |
|
yy |
|
n_retries |
|
gui |
'logical“ indicates if processing was called from the GUI environment or not. If not, processing messages are sent to a log file instead than to the console/GTK progress windows. |
out_folder_mod |
|
character array
listing all available folders (a.k.a. dates) for
the requested MODIS product on lpdaac http archive, for the years
included in the time range selected for processing.
License: GPL 3.0
Original code by Babak Naimi (.getModisList
, in
ModisDownload.R)
modified to adapt it to MODIStsp scheme and to http archive (instead than old
FTP) by:
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2016-2017)
Accessory function to find the names of HDF images corresponding to a given date and interval of spatial tiles within the lpdaac archive.
get_mod_filenames( http, used_server, user, password, n_retries, date_dir, v, h, tiled, out_folder_mod, gui )
get_mod_filenames( http, used_server, user, password, n_retries, date_dir, v, h, tiled, out_folder_mod, gui )
http |
|
used_server |
|
user |
|
password |
|
n_retries |
|
date_dir |
|
v |
|
h |
|
tiled |
|
out_folder_mod |
|
gui |
|
character array
containing names of HDF images corresponding to the
requested tiles available for the product in the selected date
License: GPL 3.0
Original code by Babak Naimi (.getModisList
, in
ModisDownload.R)
modified to adapt it to MODIStsp scheme and to http archive (instead than old
FTP) by:
Lorenzo Busetto, phD (2014-2016)
Luigi Ranghetti, phD (2016)
Helper function used in MODIStsp_process to identify which MODIS hdf layers are required for the current process. The required layers include all MODIS original layers selected by the user, plus all those required to compute the Spectral Indexes and Quality Indicators selected by the user
get_reqbands( bands_indexes_matrix, indexes_bandsel, indexes_bandnames, quality_bandsel, quality_bandnames, out_prod_folder, file_prefix, yy, DOY, out_format, reprocess )
get_reqbands( bands_indexes_matrix, indexes_bandsel, indexes_bandnames, quality_bandsel, quality_bandnames, out_prod_folder, file_prefix, yy, DOY, out_format, reprocess )
bands_indexes_matrix |
|
indexes_bandsel |
|
indexes_bandnames |
names of all indexes available for the product being processed |
quality_bandsel |
|
quality_bandnames |
names of all quality indicators available for the product being processed |
out_prod_folder |
|
file_prefix |
File prefix corresponding to the MODIS product being processed. Used to check if a given processed image already exists. |
yy |
Year corresponding to the image being processed. Used to check if a given processed image already exists. |
DOY |
DOY corresponding to the image being processed. Used to check if a given processed image already exists.. Used to check if a given processed image already exists. |
out_format |
|
reprocess |
|
req_bands_indexes
Lorenzo Busetto, phD (2017)
helper function needed to identify the ranges of dates to
be processed for a given year as a function of download_range
selection
and starting/ending dates and years
get_yeardates(download_range, yy, start_year, end_year, start_date, end_date)
get_yeardates(download_range, yy, start_year, end_year, start_date, end_date)
download_range |
|
yy |
|
start_year |
|
end_year |
|
start_date |
|
end_date |
|
OUTPUT_DESCRIPTION
Lorenzo Busetto, phD (2017)
Function which allows to use MODIStsp in batch mode by creating links
install_MODIStsp_launcher( bin_dir = NA, rscript_dir = NA, desktop_dir = NA, desktop_shortcut = TRUE, sudo = FALSE )
install_MODIStsp_launcher( bin_dir = NA, rscript_dir = NA, desktop_dir = NA, desktop_shortcut = TRUE, sudo = FALSE )
bin_dir |
|
rscript_dir |
|
desktop_dir |
|
desktop_shortcut |
|
sudo |
(Linux only) |
MODIStsp can be used also as a stand-alone tool (i.e., without opening RStudio or R-GUI) by launching a bash/batch script, which is stored in the installation folder (/ExtData/Launcher) To allow to easily find it, this function creates a desktop entry and a symbolic link to the bash script (on Linux) or a link in the Start Menu to the batch script and a shortcut on the desktop (on Windows). Note that, if the packages MODIStsp is installed in a version-dependent directory (as the default one is), this function should be re-executed after an R upgrade, otherwise the links would continue to point to the old package version!
The function is called for its side effects.
License: GPL 3.0
Luigi Ranghetti, phD (2015)
# Linux: common installation (script in /usr/bin, # desktop entry in /usr/share/applications) # (requires administrator permissions) ## Not run: # the administrator password is asked interactively install_MODIStsp_launcher(sudo = TRUE) ## End(Not run) # Linux: installation in a directory which does not require administrator # permissions ## Not run: install_MODIStsp_launcher(bin_dir = "~/bin", desktop_dir = "~/Desktop") ## End(Not run) # Windows: common installation # (script in the Start Menu and shortcut on the desktop) ## Not run: install_MODIStsp_launcher() ## End(Not run)
# Linux: common installation (script in /usr/bin, # desktop entry in /usr/share/applications) # (requires administrator permissions) ## Not run: # the administrator password is asked interactively install_MODIStsp_launcher(sudo = TRUE) ## End(Not run) # Linux: installation in a directory which does not require administrator # permissions ## Not run: install_MODIStsp_launcher(bin_dir = "~/bin", desktop_dir = "~/Desktop") ## End(Not run) # Windows: common installation # (script in the Start Menu and shortcut on the desktop) ## Not run: install_MODIStsp_launcher() ## End(Not run)
FUNCTION_DESCRIPTION
load_prodopts()
load_prodopts()
Load characteristics of the different MODIS products from prodopts_file
OUTPUT_DESCRIPTION
Lorenzo Busetto, phD (2017)
Main function for the MODIS Time Series Processing Tool (MODIStsp)
MODIStsp( ..., gui = TRUE, out_folder = NULL, out_folder_mod = NULL, opts_file = NULL, selprod = NULL, prod_version = NULL, bandsel = NULL, quality_bandsel = NULL, indexes_bandsel = NULL, sensor = NULL, download_server = NULL, downloader = NULL, user = NULL, password = NULL, download_range = NULL, start_date = NULL, end_date = NULL, spatmeth = NULL, start_x = NULL, end_x = NULL, start_y = NULL, end_y = NULL, bbox = NULL, spafile = NULL, out_projsel = NULL, output_proj = NULL, out_res_sel = NULL, out_res = NULL, resampling = NULL, reprocess = NULL, delete_hdf = NULL, nodata_change = NULL, scale_val = NULL, ts_format = NULL, out_format = NULL, compress = NULL, test = NULL, n_retries = 5, verbose = TRUE, parallel = TRUE )
MODIStsp( ..., gui = TRUE, out_folder = NULL, out_folder_mod = NULL, opts_file = NULL, selprod = NULL, prod_version = NULL, bandsel = NULL, quality_bandsel = NULL, indexes_bandsel = NULL, sensor = NULL, download_server = NULL, downloader = NULL, user = NULL, password = NULL, download_range = NULL, start_date = NULL, end_date = NULL, spatmeth = NULL, start_x = NULL, end_x = NULL, start_y = NULL, end_y = NULL, bbox = NULL, spafile = NULL, out_projsel = NULL, output_proj = NULL, out_res_sel = NULL, out_res = NULL, resampling = NULL, reprocess = NULL, delete_hdf = NULL, nodata_change = NULL, scale_val = NULL, ts_format = NULL, out_format = NULL, compress = NULL, test = NULL, n_retries = 5, verbose = TRUE, parallel = TRUE )
... |
not used for values, forces later arguments to bind by name |
gui |
|
out_folder |
|
out_folder_mod |
|
opts_file |
|
selprod |
|
prod_version |
Version of the selected MODIS product.
Currently versions |
bandsel |
|
quality_bandsel |
|
indexes_bandsel |
|
sensor |
|
download_server |
|
downloader |
download_server |
user |
|
password |
|
download_range |
|
start_date |
|
end_date |
|
spatmeth |
|
start_x |
|
end_x |
|
start_y |
|
end_y |
|
bbox |
|
spafile |
|
out_projsel |
|
output_proj |
|
out_res_sel |
|
out_res |
|
resampling |
|
reprocess |
|
delete_hdf |
|
nodata_change |
|
scale_val |
|
ts_format |
|
out_format |
|
compress |
|
test |
|
n_retries |
|
verbose |
|
parallel |
|
The function is used to:
initialize the processing (folder names, packages, etc.);
launch the GUI (MODIStsp_GUI()
) on interactive
execution, or load an options file to set processing arguments and/or
retrieve CLI inputs and run processing on non-interactive execution;
launch the routines for downloading and processing the requested datasets.
(MODIStsp_process()
)
launching the function with GUI = FALSE and without specifying a opts_file initializes arguments with default values. This allows making a test run.
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015-2017)
MODIStsp_GUI()
, MODIStsp_process()
#' # - Running the tool using the GUI # Running the tool without any option will start the GUI with the default or # last used settings, in interactive mode (i.e., with gui = TRUE). if (interactive()) { MODIStsp() } #' # - Running the tool specifying processing arguments in the call # **NOTE** Output files of examples are saved to file.path(tempdir(), "MODIStsp"). # Here we process layers __NDVI__ and __EVI__ and quality indicator __usefulness__ # of product __M*D13Q1__, considering both Terra and Aqua platforms, for dates # comprised between 2020-06-01 and 2020-06-15 and saves output to R tempdir # --> See name and available layers for product M*D13Q1. # Note that this example (as well as the following ones) is run in single # core to follow CRAN policies, by setting parallel = FALSE. # Users can exploit multicore functionalities skipping to set this argument. # The following check is performed in order not to provide errors # running the examples if HDF4 is not supported. is_hdf4_supported <- "HDF4" %in% sf::st_drivers("raster")$name MODIStsp_get_prodlayers("M*D13A2") if (is_hdf4_supported) { MODIStsp( gui = FALSE, out_folder = "$tempdir", selprod = "Vegetation_Indexes_16Days_1Km (M*D13A2)", bandsel = c("EVI", "NDVI"), quality_bandsel = "QA_usef", indexes_bandsel = "SR", user = "mstp_test" , password = "MSTP_test_01", start_date = "2020.06.01", end_date = "2020.06.15", verbose = FALSE, parallel = FALSE ) } #' # - Running the tool using the settings previously saved in a specific options file # **NOTE** Output files of examples are saved to file.path(tempdir(), "MODIStsp"). # You can run the examples with `gui = TRUE` to set a different output folder! # Here we use a test json file saved in MODIStsp installation folder which # downloads and processed 3 MOD13A2 images over the Como Lake (Lombardy, Italy) # and retrieves NDVI and EVI data, plus the Usefulness Index Quality Indicator. opts_file <- system.file("testdata/test_MOD13A2.json", package = "MODIStsp") if (is_hdf4_supported) { MODIStsp(gui = FALSE, opts_file = opts_file, verbose = TRUE, parallel = FALSE) } # Running the tool using the settings previously saved in a specific option file # and specifying the extent from a spatial file allows to re-use the same # processing settings to perform download and reprocessing on a different area opts_file <- system.file("testdata/test_MOD13A2.json", package = "MODIStsp") spatial_file <- system.file("testdata/lakeshapes/garda_lake.shp", package = "MODIStsp") if (is_hdf4_supported) { MODIStsp( gui = FALSE, opts_file = opts_file, spatmeth = "file", spafile = spatial_file, verbose = TRUE, parallel = FALSE ) } # Running the tool using the settings previously saved in a # specific options file and specifying each time the extent from a different # spatial file (e.g., to perform the same processing on several extents) # Note that you can also put all your extent files in a specific folder and # create the extent list using for example. extent_list = list.files( system.file("testdata/lakeshapes/", package = "MODIStsp"), "\\.shp$", full.names = TRUE ) extent_list opts_file <- system.file("testdata/test_MOD13A2.json", package = "MODIStsp") if (is_hdf4_supported) { for (single_shape in extent_list) { MODIStsp( gui = FALSE, opts_file = opts_file, spatmeth = "file", spafile = single_shape, verbose = TRUE, parallel = FALSE ) } } # output files are placed in separate folders: outfiles_garda <- list.files( file.path(tempdir(), "MODIStsp/garda_lake/VI_16Days_1Km_v61/NDVI"), full.names = TRUE ) outfiles_garda require(raster) if (length(outfiles_garda) > 0) { plot(raster(outfiles_garda[1] )) } outfiles_iseo <- list.files( file.path(tempdir(), "MODIStsp/iseo_lake/VI_16Days_1Km_v61/NDVI"), full.names = TRUE ) outfiles_iseo if (length(outfiles_garda) > 0) { plot(raster(outfiles_iseo[1])) } # See also https://docs.ropensci.org/MODIStsp/articles/noninteractive_execution.html
#' # - Running the tool using the GUI # Running the tool without any option will start the GUI with the default or # last used settings, in interactive mode (i.e., with gui = TRUE). if (interactive()) { MODIStsp() } #' # - Running the tool specifying processing arguments in the call # **NOTE** Output files of examples are saved to file.path(tempdir(), "MODIStsp"). # Here we process layers __NDVI__ and __EVI__ and quality indicator __usefulness__ # of product __M*D13Q1__, considering both Terra and Aqua platforms, for dates # comprised between 2020-06-01 and 2020-06-15 and saves output to R tempdir # --> See name and available layers for product M*D13Q1. # Note that this example (as well as the following ones) is run in single # core to follow CRAN policies, by setting parallel = FALSE. # Users can exploit multicore functionalities skipping to set this argument. # The following check is performed in order not to provide errors # running the examples if HDF4 is not supported. is_hdf4_supported <- "HDF4" %in% sf::st_drivers("raster")$name MODIStsp_get_prodlayers("M*D13A2") if (is_hdf4_supported) { MODIStsp( gui = FALSE, out_folder = "$tempdir", selprod = "Vegetation_Indexes_16Days_1Km (M*D13A2)", bandsel = c("EVI", "NDVI"), quality_bandsel = "QA_usef", indexes_bandsel = "SR", user = "mstp_test" , password = "MSTP_test_01", start_date = "2020.06.01", end_date = "2020.06.15", verbose = FALSE, parallel = FALSE ) } #' # - Running the tool using the settings previously saved in a specific options file # **NOTE** Output files of examples are saved to file.path(tempdir(), "MODIStsp"). # You can run the examples with `gui = TRUE` to set a different output folder! # Here we use a test json file saved in MODIStsp installation folder which # downloads and processed 3 MOD13A2 images over the Como Lake (Lombardy, Italy) # and retrieves NDVI and EVI data, plus the Usefulness Index Quality Indicator. opts_file <- system.file("testdata/test_MOD13A2.json", package = "MODIStsp") if (is_hdf4_supported) { MODIStsp(gui = FALSE, opts_file = opts_file, verbose = TRUE, parallel = FALSE) } # Running the tool using the settings previously saved in a specific option file # and specifying the extent from a spatial file allows to re-use the same # processing settings to perform download and reprocessing on a different area opts_file <- system.file("testdata/test_MOD13A2.json", package = "MODIStsp") spatial_file <- system.file("testdata/lakeshapes/garda_lake.shp", package = "MODIStsp") if (is_hdf4_supported) { MODIStsp( gui = FALSE, opts_file = opts_file, spatmeth = "file", spafile = spatial_file, verbose = TRUE, parallel = FALSE ) } # Running the tool using the settings previously saved in a # specific options file and specifying each time the extent from a different # spatial file (e.g., to perform the same processing on several extents) # Note that you can also put all your extent files in a specific folder and # create the extent list using for example. extent_list = list.files( system.file("testdata/lakeshapes/", package = "MODIStsp"), "\\.shp$", full.names = TRUE ) extent_list opts_file <- system.file("testdata/test_MOD13A2.json", package = "MODIStsp") if (is_hdf4_supported) { for (single_shape in extent_list) { MODIStsp( gui = FALSE, opts_file = opts_file, spatmeth = "file", spafile = single_shape, verbose = TRUE, parallel = FALSE ) } } # output files are placed in separate folders: outfiles_garda <- list.files( file.path(tempdir(), "MODIStsp/garda_lake/VI_16Days_1Km_v61/NDVI"), full.names = TRUE ) outfiles_garda require(raster) if (length(outfiles_garda) > 0) { plot(raster(outfiles_garda[1] )) } outfiles_iseo <- list.files( file.path(tempdir(), "MODIStsp/iseo_lake/VI_16Days_1Km_v61/NDVI"), full.names = TRUE ) outfiles_iseo if (length(outfiles_garda) > 0) { plot(raster(outfiles_iseo[1])) } # See also https://docs.ropensci.org/MODIStsp/articles/noninteractive_execution.html
Function used to add a user-defined Spectral Index to the default list of computable spectral indexes. Execution without the GUI (i.e., to add a new index from a script) is also possible (see examples).
MODIStsp_addindex( new_indexbandname = "", new_indexfullname = "", new_indexformula = "", new_indexnodata_out = "32767" )
MODIStsp_addindex( new_indexbandname = "", new_indexfullname = "", new_indexformula = "", new_indexnodata_out = "32767" )
new_indexbandname |
|
new_indexfullname |
|
new_indexformula |
|
new_indexnodata_out |
|
The function asks the user to provide the info related to the new desired Spectral Index, checks for correctness of provided information (e.g., correct bandnames, computable formula, etc...). If the index is legit, it modifies the MODIStsp_addindex.json file so to allow computation of the additional index within MODIStsp for all products containing the required reflectance bands.
To remove all custom-added spectral indexes, run MODIStsp_resetindexes()
The function is called for its side effects. On success, the MODIStsp_indexes.json is modified so to allow computation of the additional indexes.
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
# Run the GUI to interactively define a new index ## Not run: MODIStsp_addindex() ## End(Not run) # Define the new index in non-interactive execution ## Not run: MODIStsp_addindex(new_indexbandname = "SSI", new_indexfullname = "Simple Useless Index", new_indexformula = "b2_NIR+b1_Red") ## End(Not run)
# Run the GUI to interactively define a new index ## Not run: MODIStsp_addindex() ## End(Not run) # Define the new index in non-interactive execution ## Not run: MODIStsp_addindex(new_indexbandname = "SSI", new_indexfullname = "Simple Useless Index", new_indexformula = "b2_NIR+b1_Red") ## End(Not run)
Internal function dealing with download of MODIS hdfs from http remote server for a given date.
MODIStsp_download( modislist, out_folder_mod, download_server, http, n_retries, use_aria, date_dir, year, DOY, user, password, sens_sel, date_name, gui, verbose )
MODIStsp_download( modislist, out_folder_mod, download_server, http, n_retries, use_aria, date_dir, year, DOY, user, password, sens_sel, date_name, gui, verbose )
modislist |
|
out_folder_mod |
|
download_server |
|
http |
|
n_retries |
|
use_aria |
|
date_dir |
|
year |
|
DOY |
|
user |
|
password |
|
sens_sel |
|
date_name |
|
gui |
|
verbose |
|
The function is called for its side effects
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
function used to extract time series data from rts files created by MODIStsp on spatial locations provided in the form of "R" spatial objects (SpatialPoints, SpatialPolygons, etc.)
MODIStsp_extract( in_rts, sf_object, start_date = NULL, end_date = NULL, id_field = NULL, FUN = "mean", out_format = "xts", small = TRUE, small_method = "centroids", na.rm = TRUE, verbose = FALSE )
MODIStsp_extract( in_rts, sf_object, start_date = NULL, end_date = NULL, id_field = NULL, FUN = "mean", out_format = "xts", small = TRUE, small_method = "centroids", na.rm = TRUE, verbose = FALSE )
in_rts |
A |
sf_object |
"sf" object OR name of an GDAL-readable vector file specifying the "area" from which data has to be extracted.
|
start_date |
object of class |
end_date |
object of class |
id_field |
|
FUN |
function to summarize the values (e.g. mean) on polygon data frames. The function should take a single numeric vector as argument and return a single value (e.g. mean, min or max), and accept a na.rm argument. Thus, standard R functions not including an na.rm argument must be wrapped as in this example: fun=function(x,...)length(x). Defaults to "mean" |
out_format |
|
small |
|
small_method |
|
na.rm |
|
verbose |
|
The function takes as input a RasterStack object containing time information
in the "z" attribute (set by raster::setZ
), a starting and ending date
and a standard "R" spatial object, and returns the time series for the spatial locations
specified in the spatial object in the form of a "R" xts object OR a plain data.frame
with a "date" column in first position.
If the input spatial object is a "point" or "line" one, the output object
contains one column for each specified point, or for each cell intersecting
the line, and one line for each date. If the input spatial object is a "polygon"
one, the output object contains one column for each polygon, containing values
obtained applying the function specified as the FUN argument over all pixels
belonging to the polygon, and one line for each date.
data.frame or xts object. Each column of data corresponds to one point or one polygon, each row to a date.
License: GPL 3.0
Lorenzo Busetto, phD (2015 - 2017) email: [email protected]
## Not run: # Extract average and standard deviation values from a rts object created by # MODIStsp for each polygon of a shapefile, for each date in the period # between 2001-01-01 and 2014-12-31 # The example uses tif files in testdata/VI_16Days_500m_v61 to build # a MODIStsp rasterStack corresponding to the 2016 time series of the NDVI index # over the Como Lake (Italy). It then extracts data on polygons corresponding # to different land cover classes saved in testdata/extract_polys.shp # First, prepare the test dataset. # __NOTE__ To avoid redownloading, here we copy some test data from MODIStsp # installation folder to tempdir and use it to create a test time series. test_zip <- system.file("testdata/VI_16Days_500m_v6/NDVI.zip", package = "MODIStsp") dir.create(file.path(tempdir(), "MODIStsp/VI_16Days_500m_v61"), showWarnings = FALSE, recursive = TRUE) utils::unzip(test_zip, exdir = file.path(tempdir(), "MODIStsp/VI_16Days_500m_v61")) opts_file <- system.file("testdata/test_extract.json", package = "MODIStsp") MODIStsp(opts_file = opts_file, gui = FALSE, verbose = FALSE) # Now load the MODIStsp stack: This is a MODIS NDVI time series ranging between # 2016-01-01 and 2016-12-18 # __NOTE__: MODIStsp rasterStack files are always saved in the "Time_Series\/RData" # subfolder of your main output folder - see # "https://docs.ropensci.org/MODIStsp/articles/output.html") # Specify the filename of the RData RasterStack of interest stack_file <- file.path(tempdir(), "MODIStsp/VI_16Days_500m_v61/Time_Series/RData/Terra/NDVI", "MOD13A1_NDVI_1_2016_353_2016_RData.RData") basename(stack_file) ts_data <- get(load(stack_file)) ts_data # Now load a shapefile containing polygons from which we want to extract data polygons <- sf::st_read(system.file("testdata/extract_polys.shp", package = "MODIStsp"), quiet = TRUE) polygons # Finally, extract the average values for each polygon and date and plot the # results out_dataavg <- suppressMessages(MODIStsp_extract(ts_data, polygons, id_field = "lc_type", small = FALSE)) head(out_dataavg) plot(out_dataavg, legend.loc = "topleft") # use a different summarization function out_datasd <- MODIStsp_extract(ts_data, polygons, id_field = "lc_type", FUN = "sd", small = FALSE) head(out_datasd) # (See also https://docs.ropensci.org/MODIStsp/articles/Analyze.html for a # worked-out example) ## End(Not run)
## Not run: # Extract average and standard deviation values from a rts object created by # MODIStsp for each polygon of a shapefile, for each date in the period # between 2001-01-01 and 2014-12-31 # The example uses tif files in testdata/VI_16Days_500m_v61 to build # a MODIStsp rasterStack corresponding to the 2016 time series of the NDVI index # over the Como Lake (Italy). It then extracts data on polygons corresponding # to different land cover classes saved in testdata/extract_polys.shp # First, prepare the test dataset. # __NOTE__ To avoid redownloading, here we copy some test data from MODIStsp # installation folder to tempdir and use it to create a test time series. test_zip <- system.file("testdata/VI_16Days_500m_v6/NDVI.zip", package = "MODIStsp") dir.create(file.path(tempdir(), "MODIStsp/VI_16Days_500m_v61"), showWarnings = FALSE, recursive = TRUE) utils::unzip(test_zip, exdir = file.path(tempdir(), "MODIStsp/VI_16Days_500m_v61")) opts_file <- system.file("testdata/test_extract.json", package = "MODIStsp") MODIStsp(opts_file = opts_file, gui = FALSE, verbose = FALSE) # Now load the MODIStsp stack: This is a MODIS NDVI time series ranging between # 2016-01-01 and 2016-12-18 # __NOTE__: MODIStsp rasterStack files are always saved in the "Time_Series\/RData" # subfolder of your main output folder - see # "https://docs.ropensci.org/MODIStsp/articles/output.html") # Specify the filename of the RData RasterStack of interest stack_file <- file.path(tempdir(), "MODIStsp/VI_16Days_500m_v61/Time_Series/RData/Terra/NDVI", "MOD13A1_NDVI_1_2016_353_2016_RData.RData") basename(stack_file) ts_data <- get(load(stack_file)) ts_data # Now load a shapefile containing polygons from which we want to extract data polygons <- sf::st_read(system.file("testdata/extract_polys.shp", package = "MODIStsp"), quiet = TRUE) polygons # Finally, extract the average values for each polygon and date and plot the # results out_dataavg <- suppressMessages(MODIStsp_extract(ts_data, polygons, id_field = "lc_type", small = FALSE)) head(out_dataavg) plot(out_dataavg, legend.loc = "topleft") # use a different summarization function out_datasd <- MODIStsp_extract(ts_data, polygons, id_field = "lc_type", FUN = "sd", small = FALSE) head(out_datasd) # (See also https://docs.ropensci.org/MODIStsp/articles/Analyze.html for a # worked-out example) ## End(Not run)
Function used to retrieve the names of original MODIS layers, quality layers and eventually available spectral indexes given a MODIS product code. It is useful to identify the names of the layers to be processed when launching MODIStsp from the CLI.
MODIStsp_get_prodlayers(prodname, prodver = "061")
MODIStsp_get_prodlayers(prodname, prodver = "061")
prodname |
character containing the code of the desired MODIS product. NOTE: for products available separately for Terra and Aqua (e.g., MOD13Q1/MYD13Q1), use the code MD_code_ (e.g., MD13Q1) |
prodver |
character containing the version of the desired MODIS product.
Currently versions |
list, containing the slots: prodname
, bandnames
, quality_bandnames
and
indexes_bandnames
, band_fullnames
, quality_fullnames
, indexes_fullnames
License: GPL 3.0
Lorenzo Busetto, phD (2014-2020)
Luigi Ranghetti, phD (2021)
# Get layers of product M*13Q1 based on code MODIStsp_get_prodlayers("M*13Q1") # Get layers of product M*13Q1 based on full name MODIStsp_get_prodlayers("Vegetation Indexes_16Days_250m (M*D13Q1)") # Get indexes names of product M*13Q1 based on full name MODIStsp_get_prodlayers("MCD43C4")$indexes_bandnames
# Get layers of product M*13Q1 based on code MODIStsp_get_prodlayers("M*13Q1") # Get layers of product M*13Q1 based on full name MODIStsp_get_prodlayers("Vegetation Indexes_16Days_250m (M*D13Q1)") # Get indexes names of product M*13Q1 based on full name MODIStsp_get_prodlayers("MCD43C4")$indexes_bandnames
Function used to retrieve the names of available MODIS products
MODIStsp_get_prodnames()
MODIStsp_get_prodnames()
character array of product names
License: GPL 3.0
Lorenzo Busetto, phD (2014-2020)
# Get MODIStsp product names MODIStsp_get_prodnames()
# Get MODIStsp product names MODIStsp_get_prodnames()
Function used to generate and handle the GUI used to allow selection of MODIStsp processing parameters.
MODIStsp_GUI()
MODIStsp_GUI()
the function is called for its side effects - opening the GUI and allowing to set, save, load options and eventually launch the processing.
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
Main processing function of MODIStsp. Takes as input processing parameters specified by the user and performs all required processing.
MODIStsp_process(proc_opts, n_retries, verbose = TRUE, parallel = TRUE)
MODIStsp_process(proc_opts, n_retries, verbose = TRUE, parallel = TRUE)
proc_opts |
|
n_retries |
|
verbose |
|
parallel |
|
After retrieving the input processing options, the function
Accesses NASA http archive to determine the list of files to be
downloaded/processed (or, in case of offline processing, get the list
of already available hdf files present in out_mod_folder
);
Performs all required processing steps on each date (download, reprojection, resize, mosaicing, Spectral Indexes and Quality indicators computation);
Creates virtual files of the processed time series.
Reprojection and resize is dealt with by accessing gdal routines through the
gdalUtilities
package.
Extraction of bitfields from Quality layers is done though bitwise computation
Checks are done in order to not re-download already existing HDF images, and not
reprocess already processed dates (if the user did not specify that)
The function is called for its side effects.
Thanks Tomislav Hengl and Babak Naimi, whose scripts made the starting point for development of this function (ModisDownload; Download_and_resampling_of_MODIS_images)
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
Internal function used to perform the required spatial
processing on MODIS original hdf layers (reprojection, resizing, resampling,
mosaicing, computation of scaling factors). The function is based on the
use of gdal
routines.
MODIStsp_process_bands( out_folder_mod, modislist, outproj_str, mod_proj_str, sens_sel, band, bandname, date_name, datatype, nodata_in, nodata_out, full_ext, bbox, scale_val, scale_factor, offset, out_format, outrep_file, compress, out_res_sel, out_res, resampling, nodata_change, gui, verbose, parallel )
MODIStsp_process_bands( out_folder_mod, modislist, outproj_str, mod_proj_str, sens_sel, band, bandname, date_name, datatype, nodata_in, nodata_out, full_ext, bbox, scale_val, scale_factor, offset, out_format, outrep_file, compress, out_res_sel, out_res, resampling, nodata_change, gui, verbose, parallel )
out_folder_mod |
|
modislist |
|
outproj_str |
|
mod_proj_str |
|
sens_sel |
|
band |
|
bandname |
|
date_name |
|
datatype |
|
nodata_in |
|
nodata_out |
|
full_ext |
|
bbox |
|
scale_val |
|
scale_factor |
|
offset |
|
out_format |
|
outrep_file |
|
compress |
|
out_res_sel |
|
out_res |
|
resampling |
|
nodata_change |
|
gui |
|
verbose |
|
parallel |
|
The function is called for its side effects
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
function used to compute spectral indexes, given the index formula
MODIStsp_process_indexes( out_filename, out_prod_folder, formula, bandnames, nodata_out, indexes_nodata_out, file_prefix, compress, yy, out_format, DOY, scale_val )
MODIStsp_process_indexes( out_filename, out_prod_folder, formula, bandnames, nodata_out, indexes_nodata_out, file_prefix, compress, yy, out_format, DOY, scale_val )
out_filename |
|
out_prod_folder |
|
formula |
|
bandnames |
|
nodata_out |
|
indexes_nodata_out |
|
file_prefix |
|
compress |
|
yy |
|
out_format |
|
DOY |
|
scale_val |
|
the function parses the index formula to identify the required bands. On the basis of identified bands, it retrieves the reflectance bands required, gets the data into R raster objects, performs the computation and stores results in a GeoTiff or ENVI raster file
NULL - new raster file saved in out_filename
License: GPL 3.0
Lorenzo Busetto, phD (2017)
Luigi Ranghetti, phD (2017)
function used to extract quality indicator from MODIS aggregated quality layers
MODIStsp_process_QA_bits( out_filename, in_source_file, bitN, out_format, nodata_source, nodata_qa_in, nodata_qa_out, compress )
MODIStsp_process_QA_bits( out_filename, in_source_file, bitN, out_format, nodata_source, nodata_qa_in, nodata_qa_out, compress )
out_filename |
|
in_source_file |
|
bitN |
|
out_format |
output format (ENVI or GTiff) |
nodata_source |
|
nodata_qa_in |
|
nodata_qa_out |
|
compress |
|
On the basis of the name of the image containing the aggregated quality information
(in_source_file``) and of the position of the bit fields corresponding to the QI of interest in the bitfield representation (
bitN“), the function extracts the correct information exploiting
bitwise operators, and save the result in a new raster image
License: GPL 3.0
Based on the modis.qc.R
script by Yann Chemin (2008) (https://goo.gl/7Fhreo)
license GPL 3.0
Lorenzo Busetto, phD (2017)
Luigi Ranghetti, phD (2017)
function used to parse the XML file used to store the characteristics of MODIS Land Products and store them in the "prod_opts" data frame
MODIStsp_read_xml(prodopts_file, xml_file)
MODIStsp_read_xml(prodopts_file, xml_file)
prodopts_file |
string filename of the RData in which to store the data parsed from the XML file |
xml_file |
string filename of the XML file containing the MODIS products characteristics |
The function parses the XML file product by product, stores data in a data frame and saves the data frame within the "MODIStsp_previous" RData file as a list of lists
NULL - retrieved data are stored in the specified RData file
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
Function used to remove all user-defined Spectral Indexes from MODIStsp, thus resetting the list of available indexes to the default ones.
MODIStsp_resetindexes()
MODIStsp_resetindexes()
The function is called for its side effects. On success, the MODIStsp_indexes.json file is modified so to remove all previously custom-specified Spectral Indexes.
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017) [email protected]
## Not run: # Remove all custom-defined spectral indexes from an options file # Add a custom index for testing purposes library(jsonlite) opts_jsfile = system.file("testdata/test_addindex.json", package = "MODIStsp") MODIStsp_addindex( opts_jsfile = opts_jsfile, gui = FALSE, new_indexbandname = paste0("Index_", as.character(sample(10000, 1))), new_indexformula = "b1_Red - b2_NIR", new_indexfullname = paste0("Index_", as.character(sample(10000, 1))) ) opts <- jsonlite::fromJSON(indexes_file) opts$custom_indexes[1] # Now remove all custom indexes MODIStsp_resetindexes() opts <- jsonlite::fromJSON(opts_jsfile) opts$custom_indexes[1] ## End(Not run)
## Not run: # Remove all custom-defined spectral indexes from an options file # Add a custom index for testing purposes library(jsonlite) opts_jsfile = system.file("testdata/test_addindex.json", package = "MODIStsp") MODIStsp_addindex( opts_jsfile = opts_jsfile, gui = FALSE, new_indexbandname = paste0("Index_", as.character(sample(10000, 1))), new_indexformula = "b1_Red - b2_NIR", new_indexfullname = paste0("Index_", as.character(sample(10000, 1))) ) opts <- jsonlite::fromJSON(indexes_file) opts$custom_indexes[1] # Now remove all custom indexes MODIStsp_resetindexes() opts <- jsonlite::fromJSON(opts_jsfile) opts$custom_indexes[1] ## End(Not run)
Function used to create virtual files from time series of single-band files corresponding to different acquisition dates. The function takes as input the folder in which the single-band files are stored, and creates a ENVI Meta file and/or a GDAL vrt file that allows access to the full time series as if it was a single physical file. Created virtual files are stored in the "Time Series" subfolder of 'out_prod_folder“
MODIStsp_vrt_create( sensor, out_prod_folder, bandnames, bandsel, nodata_out, indexes_bandnames, indexes_bandsel, indexes_nodata_out, quality_bandnames, quality_bandsel, quality_nodata_out, file_prefixes, ts_format, out_format, verbose )
MODIStsp_vrt_create( sensor, out_prod_folder, bandnames, bandsel, nodata_out, indexes_bandnames, indexes_bandsel, indexes_nodata_out, quality_bandnames, quality_bandsel, quality_nodata_out, file_prefixes, ts_format, out_format, verbose )
sensor |
|
out_prod_folder |
|
bandnames |
names of all layers available for the product being processed |
bandsel |
|
nodata_out |
|
indexes_bandnames |
names of all indexes available for the product being processed |
indexes_bandsel |
|
indexes_nodata_out |
nodata value for indexes vrts |
quality_bandnames |
names of all quality indicators available for the product being processed |
quality_bandsel |
|
quality_nodata_out |
nodata value for quality vrts |
file_prefixes |
|
ts_format |
|
out_format |
|
verbose |
|
NULL - the function is called for its side effects
License: GPL 3.0
Lorenzo Busetto, phD (2014-2017)
Luigi Ranghetti, phD (2015)
helper function to provide processing messages
process_message(mess_text, verbose = TRUE)
process_message(mess_text, verbose = TRUE)
mess_text |
|
verbose |
|
The function is called for its side effects
Lorenzo Busetto, phD (2017)
Helper function used to reproject bounding boxes; setting the parameter 'enlarge' allows to choose if the new one would be the one which completely includes the original extent in the output projection, or if is simply the one obtained by reprojecting the upper-left and the lower-right corners.
reproj_bbox(bbox, in_proj, out_proj, enlarge = TRUE)
reproj_bbox(bbox, in_proj, out_proj, enlarge = TRUE)
bbox |
The input bounding box (it can be a matrix obtained from |
in_proj |
( |
out_proj |
|
enlarge |
'logical“ if TRUE, the reprojected bounding box is the one which completely include the original one; if FALSE, it is simply the one obtained by reprojecting the upper-left and the lower-right corners. |
License: GPL 3.0
Luigi Ranghetti, phD (2015)
FUNCTION_DESCRIPTION
set_bandind_matrix( bandnames, bandsel, indexes_bandnames, indexes_bandsel, indexes_formula, quality_bandnames, quality_bandsel, quality_source )
set_bandind_matrix( bandnames, bandsel, indexes_bandnames, indexes_bandsel, indexes_formula, quality_bandnames, quality_bandsel, quality_source )
bandnames |
names of all layers available for the product being processed |
bandsel |
|
indexes_bandnames |
names of all indexes available for the product being processed |
indexes_bandsel |
|
indexes_formula |
formulas of all indexes available for the product being processed |
quality_bandnames |
names of all quality indicators available for the product being processed |
quality_bandsel |
|
quality_source |
sources of data (original layers) of all quality indicators available for the product being processed |
matrix
containing info on which bands are needed for computing
each available QI or SI
Lorenzo Busetto, phD (2017)
Internal functions: split_nodata_values splits the ranges of NODATA saved in the xml product file to a readable vector of NoData values; create_nodata_rcl creates the matrix for the reclassification of NODATA values to be used with raster::reclassify function.
split_nodata_values(nodata_in, take_all = TRUE) create_nodata_rcl(nodata_in, nodata_out)
split_nodata_values(nodata_in, take_all = TRUE) create_nodata_rcl(nodata_in, nodata_out)
nodata_in |
Character vector corresponding to input NoData values as saved in the xml product file (one or more values per band). |
take_all |
Logical: if TRUE (default), all the NoData values are considered; if FALSE, only the last one is taken. See "details" for the meaning of this parameter. |
nodata_out |
Character vector corresponding to output NoData values as saved in the xml product file (one single value per band). |
MODIS products can have more than one NoData values (sometimes
with different meanings, e.g. 255 = "fill" and 254 = "detector saturated"
in MOD09A1 product). By setting
"Change NoData values" to "Yes" in the GUI, all the NoData values are
coerced to one single new NoData value; conversely, setting it to "No"
only one value is assumed to be NoData.
The parameter take_all
is assumed to be used in this way, by using this
function with take_all = TRUE
with "Change NoData values" = "Yes" and
take_all = FALSE
with "Change NoData values" = "No".
In the xml product file, NoData ranges are set as:
x
for products with single NoData values;
x,y,z
for products with a vector of NoData values;
x:y
for products with a range of NoData values;
x:y,z
for a combination of NoData ranges and/or values.
In split_nodata_values NoData values are assumed to be integer: this means that intervals are splitted in integer values (e.g. "250:255" becomes "250 251 252 253 254 255"). Conversely, function create_nodata_rcl creates intervals, so it can also manage float values (in practice, this should not make difference within MODIS products, since NoData values are always integer values).
This function interprets these strings and convert them in vectors with single values. Notice that the first NoData value is the only one which is considered if 'Change NoData values' was set to 'No'.
split_nodata_values returns a list with the same length
of nodata_in
vector, in which each element
is a vector with all the NoData values.
create_nodata_rcl returns a list of matrices in the format
specified for parameter rcl
in raster::reclassify.
The parameter right
is intended to be used as right = NA
.
Luigi Ranghetti, phD (2018)
MODIStsp:::create_nodata_rcl(c("255","250,254:255"), c("255","255"))
MODIStsp:::create_nodata_rcl(c("255","250,254:255"), c("255","255"))