polyStatsLevel2BVPM: Compute descriptive statistics of GEDI Canopy Cover and...

View source: R/polyStatsLevel2BVPM.R

polyStatsLevel2BVPMR Documentation

Compute descriptive statistics of GEDI Canopy Cover and Vertical Profile Metrics

Description

Computes a Series of Statistics of GEDI derived Canopy Cover and Vertical Profile metrics within a given area defined or not by a polygon

Usage

polyStatsLevel2BVPM(level2BVPM, func, id=NULL)

Arguments

level2BVPM

A GEDI Level2BVPM object (output of getLevel2BVPM() function). An S4 object of class "data.table".

func

The function to be applied for computing the defined statistics

id

A vector containing the polygon id for each GEDI observation. Default is NULL

Value

Returns an S4 object of class data.table::data.table Containing Statistics of GEDI level2BVPM defined metrics

See Also

https://lpdaac.usgs.gov/products/gedi02_bv002/

Examples

# Specifying the path to GEDI level2B data (zip file)
outdir <- tempdir()
level2B_fp_zip <- system.file("extdata",
  "GEDI02_B_2019108080338_O01964_T05337_02_001_01_sub.zip",
  package = "rGEDI"
)

# Unzipping GEDI level2A data
level2Bpath <- unzip(level2B_fp_zip, exdir = outdir)

# Reading GEDI level2B data (h5 file)
level2b <- readLevel2B(level2Bpath = level2Bpath)

# Specifying the path to shapefile
polygon_filepath <- system.file("extdata", "stands_cerrado.shp", package = "rGEDI")

# Reading shapefile as sf object
library(sf)
polygon <- sf::st_read(polygon_filepath)

# Extracting GEDI Canopy Cover and Vertical Profile Metrics
level2BVPM <- getLevel2BVPM(level2b)
head(level2BVPM)

# Clipping GEDI data by geometry
level2BVPM_clip <- clipLevel2BVPMGeometry(level2BVPM, polygon, split_by = "id")

# Define your own function
mySetOfMetrics <- function(x) {
  metrics <- list(
    min = min(x), # Min of x
    max = max(x), # Max of x
    mean = mean(x), # Mean of x
    sd = sd(x) # Sd of x
  )
  return(metrics)
}

# Computing the max of the Total Plant Area Index
pai_max <- polyStatsLevel2BVPM(level2BVPM_clip, func = max(pai), id = NULL)
pai_max

# Computing the max of the Total Plant Area Index stratified by polygon
pai_max_poly <- polyStatsLevel2BVPM(level2BVPM_clip, func = max(pai), id = "poly_id")
head(pai_max_poly)

# Computing the serie of statistics of canopy cover stratified by polygon
cover_metrics <- polyStatsLevel2BVPM(level2BVPM_clip,
  func = mySetOfMetrics(cover),
  id = level2BVPM_clip$id
)
head(cover_metrics)
close(level2b)

carlos-alberto-silva/rGEDI documentation built on Oct. 18, 2024, 4:46 a.m.