estarfm_job: Execute a single self-contained time-series imagefusion job...

Description Usage Arguments Details Value Author(s) References Examples

View source: R/estarfm_job.R

Description

A wrapper function for execute_estarfm_job_cpp. Intended to execute a single job, that is a number of predictions based on the same input pairs. It ensures that all of the arguments passed are of the correct type and creates sensible defaults.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
estarfm_job(
  input_filenames,
  input_resolutions,
  input_dates,
  pred_dates,
  pred_filenames,
  pred_area,
  winsize,
  date1,
  date3,
  n_cores,
  data_range_min,
  data_range_max,
  uncertainty_factor,
  number_classes,
  hightag,
  lowtag,
  MASKIMG_options,
  MASKRANGE_options,
  use_local_tol,
  use_quality_weighted_regression,
  output_masks,
  use_nodata_value,
  verbose = TRUE
)

Arguments

input_filenames

A string vector containing the filenames of the input images

input_resolutions

A string vector containing the resolution-tags (corresponding to the arguments hightag and lowtag, which are by default "high" and "low") of the input images.

input_dates

An integer vector containing the dates of the input images.

pred_dates

An integer vector containing the dates for which images should be predicted.

pred_filenames

A string vector containing the filenames for the predicted images. Must match pred_dates in length and order. Must include an extension relating to one of the drivers supported by GDAL, such as ".tif".

pred_area

(Optional) An integer vector containing parameters in image coordinates for a bounding box which specifies the prediction area. The prediction will only be done in this area. (x_min, y_min, width, height). By default will use the entire area of the first input image.

winsize

(Optional) Window size of the rectangle around the current pixel. Default is 51.

date1

(Optional) Set the date of the first input image pair. By default, will use the pair with the lowest date value.

date3

(Optional) Set the date of the second input image pair. By default, will use the pair with the highest date value.

n_cores

(Optional) Set the number of cores to use when using parallelization. Default is 1.

data_range_min

(Optional) When predicting pixel values ESTARFM can exceed the values that appear in the image. To prevent from writing invalid values (out of a known data range) you can set bounds. By default, the value range will be limited by the natural data range (e. g. -32767 for INT2S).

data_range_max

(Optional) When predicting pixel values ESTARFM can exceed the values that appear in the image. To prevent from writing invalid values (out of a known data range) you can set bounds. By default, the value range will be limited by the natural data range (e. g. 32767 for INT2S).

uncertainty_factor

(Optional) Sets the uncertainty factor. This is multiplied to the upper limit of the high resolution range. The range can be given by a mask range. Default: 0.002 (0.2 per cent)

number_classes

(Optional) The number of classes used for similarity. Note all channels of a pixel are considered for similarity. So this value holds for each channel, e. g. with 3 channels there are n^3 classes. Default: c-th root of 64, where c is the number of channels.

hightag

(Optional) A string which is used in input_resolutions to describe the high-resolution images. Default is "high".

lowtag

(Optional) A string which is used in input_resolutions to describe the low-resolution images. Default is "low".

MASKIMG_options

(Optional) A string containing information for a mask image (8-bit, boolean, i. e. consists of 0 and 255). "For all input images the pixel values at the locations where the mask is 0 is replaced by the mean value." Example: --mask-img=some_image.png

MASKRANGE_options

(Optional) Specify one or more intervals for valid values. Locations with invalid values will be masked out. Ranges should be given in the format '[<float>,<float>]', '(<float>,<float>)', '[<float>,<float>' or '<float>,<float>]'. There are a couple of options:'

  • "–mask-valid-ranges" Intervals which are marked as valid. Valid ranges can excluded from invalid ranges or vice versa, depending on the order of options.

  • "–mask-invalid-ranges" Intervals which are marked as invalid. Invalid intervals can be excluded from valid ranges or vice versa, depending on the order of options.

  • "–mask-high-res-valid-ranges" This is the same as –mask-valid-ranges, but is applied only for the high resolution images.

  • "–mask-high-res-invalid-ranges" This is the same as –mask-invalid-ranges, but is applied only for the high resolution images.

  • "–mask-low-res-valid-ranges" This is the same as –mask-valid-ranges, but is applied only for the low resolution images.

  • "–mask-low-res-invalid-ranges" This is the same as –mask-invalid-ranges, but is applied only for the low resolution images.

use_local_tol

(Optional) This enables the usage of local tolerances to find similar pixels instead of using the global tolerance. When searching similar pixels, a tolerance of 2σ/m is used. This options sets whether is calculated only from the local window region around the central pixel or from the whole image. Default is "false".

use_quality_weighted_regression

(Optional) This enables the smooth weighting of the regression coefficient by its quality. The regression coefficient is not limited strictly by the quality, but linearly blended to 1 in case of bad quality. Default is "false".

output_masks

(Optional) Write mask images to disk? Default is "false".

use_nodata_value

(Optional) Use the nodata value as invalid range for masking? Default is "true".

verbose

(Optional) Print progress updates to console? Default is "true".

Details

Executes the ESTARFM algorithm to create a number of synthetic high-resolution images from two pairs of matching high- and low-resolution images. Assumes that the input images already have matching size. See the original paper for details (Note: There is a difference to the algorithm as described in the paper though. The regression for R is now done with all candidates of one window. This complies to the reference implementation, but not to the paper, since there the regression is done only for the candidates that belong to one single coarse pixel. However, the coarse grid is not known at prediction and not necessarily trivial to find out (e. g. in case of higher order interpolation).

Value

Nothing. Output files are written to disk. The Geoinformation for the output images is adopted from the first input pair images.

Author(s)

Christof Kaufmann (C++)

Johannes Mast (R)

References

Zhu, X., Chen, J., Gao, F., Chen, X., & Masek, J. G. (2010). An enhanced spatial and temporal adaptive reflectance fusion model for complex heterogeneous regions. Remote Sensing of Environment, 114(11), 2610-2623.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
# Load required libraries
library(ImageFusion)
library(raster)
# Get filesnames of high resolution images
landsat <- list.files(
  system.file("landsat/filled",
              package = "ImageFusion"),
  ".tif",
  recursive = TRUE,
  full.names = TRUE
)

# Get filesnames of low resolution images
modis <- list.files(
  system.file("modis",
              package = "ImageFusion"),
  ".tif",
  recursive = TRUE,
  full.names = TRUE
)

#Select the first two landsat images 
landsat_sel <- landsat[1:2]
#Select the corresponding modis images
modis_sel <- modis[1:10]
# Create output directory in temporary folder
out_dir <- file.path(tempdir(),"Outputs")
if(!dir.exists(out_dir)) dir.create(out_dir, recursive = TRUE)

#Run the job, fusing two images
estarfm_job(input_filenames = c(landsat_sel,modis_sel),
            input_resolutions = c("high","high",
                                  "low","low","low",
                                  "low","low","low",
                                  "low","low","low","low"),
            input_dates = c(68,77,68,69,70,71,72,73,74,75,76,77),
            pred_dates = c(72,74),
            pred_filenames = c(file.path(out_dir,"estarfm_72.tif"),
                               file.path(out_dir,"estarfm_74.tif"))
)
# remove the output directory
unlink(out_dir,recursive = TRUE)

ImageFusion documentation built on March 4, 2021, 5:06 p.m.