spwb_land {medfateland}R Documentation

Watershed simulations

Description

Functions to perform simulations on a watershed described by a set of connected grid cells.

Usage

spwb_land(
  r,
  sf,
  SpParams,
  meteo = NULL,
  dates = NULL,
  CO2ByYear = numeric(0),
  summary_frequency = "years",
  local_control = defaultControl(soilDomains = "single"),
  watershed_control = default_watershed_control(),
  parallelize = FALSE,
  num_cores = detectCores() - 1,
  chunk_size = NULL,
  progress = TRUE
)

growth_land(
  r,
  sf,
  SpParams,
  meteo = NULL,
  dates = NULL,
  CO2ByYear = numeric(0),
  summary_frequency = "years",
  local_control = medfate::defaultControl(soilDomains = "single"),
  watershed_control = default_watershed_control(),
  parallelize = FALSE,
  num_cores = detectCores() - 1,
  chunk_size = NULL,
  progress = TRUE
)

fordyn_land(
  r,
  sf,
  SpParams,
  meteo = NULL,
  dates = NULL,
  CO2ByYear = numeric(0),
  local_control = medfate::defaultControl(soilDomains = "single"),
  watershed_control = default_watershed_control(),
  dispersal_control = default_dispersal_control(),
  management_function = NULL,
  parallelize = FALSE,
  num_cores = detectCores() - 1,
  chunk_size = NULL,
  progress = TRUE
)

cell_neighbors(sf, r)

## S3 method for class 'spwb_land'
summary(object, ...)

## S3 method for class 'growth_land'
summary(object, ...)

Arguments

r

An object of class SpatRaster, defining the raster topology.

sf

An object of class sf with the following columns:

  • geometry: Spatial point geometry corresponding to cell centers.

  • elevation: Elevation above sea level (in m).

  • slope: Slope (in degrees).

  • aspect: Aspect (in degrees).

  • land_cover_type: Land cover type of each grid cell (values should be 'wildland', 'agriculture', 'rock', 'artificial' or 'water').

  • forest: Objects of class forest.

  • soil: Objects of class soil or data frames of physical properties.

  • state: Objects of class spwbInput or growthInput (optional).

  • meteo: Data frames with weather data (required if parameter meteo = NULL).

  • crop_factor: Crop evapo-transpiration factor. Only required for 'agriculture' land cover type.

  • local_control: A list of control parameters (optional). Used to override function parameter local_control for specific cells (values can be NULL for the remaining ones).

  • snowpack: An optional numeric vector with the snow water equivalent content of the snowpack in each cell (in mm). If missing it will be initialized to zero.

  • management_arguments: Lists with management arguments (optional, relevant for fordyn_land only).

  • result_cell: A logical indicating that local model results are desired (optional, relevant for spwb_land and growth_land only). Model results are only produced for wildland and agriculture cells.

When using TETIS watershed model, the following columns are also REQUIRED:

  • depth_to_bedrock: Depth to bedrock (mm).

  • bedrock_conductivity: Bedrock (saturated) conductivity (in m·day-1).

  • bedrock_porosity: Bedrock porosity (the proportion of pore space in the rock).

When using TETIS watershed model, the following columns are OPTIONAL:

  • aquifer: A numeric vector with the water content of the aquifer in each cell (in mm). If missing, it will be initialized to zero.

  • deep_aquifer_loss: A numeric vector with the maximum daily loss to a deeper aquifer (in mm·day-1). If missing all cells take their value from deep_aquifer_loss in default_watershed_control

SpParams

A data frame with species parameters (see SpParamsMED).

meteo

Input meteorological data (see spwb_spatial and details).

dates

A Date object describing the days of the period to be modeled.

CO2ByYear

A named numeric vector with years as names and atmospheric CO2 concentration (in ppm) as values. Used to specify annual changes in CO2 concentration along the simulation (as an alternative to specifying daily values in meteo).

summary_frequency

Frequency in which cell summary will be produced (e.g. "years", "months", ...) (see cut.Date). In fordyn_land summaries are always produced at monthly resolution.

local_control

A list of control parameters (see defaultControl) for function spwb_day or growth_day. By default, parameter soilDomains is set to "single", meaning a single-domain Richards model.

watershed_control

A list of watershed control parameters (see default_watershed_control). Importantly, the sub-model used for lateral water flows - either Francés et al. (2007) or Caviedes-Voullième et al. (2023) - is specified there.

parallelize

Boolean flag to try parallelization (see details).

num_cores

Integer with the number of cores to be used for parallel computation (by default it will use all clusters minus one).

chunk_size

Integer indicating the size of chunks to be sent to different processes (by default, the number of spatial elements divided by the number of cores).

progress

Boolean flag to display progress information for simulations.

dispersal_control

A list of dispersal control parameters (see default_dispersal_control). If NULL, then dispersal is not simulated.

management_function

A function that implements forest management actions (see fordyn). of such lists, one per spatial unit.

object

An object of class spwb_land or groth_land

...

Additional parameters for summary functions

Details

The default soilDomains = "single" means that vertical bulk soil flows are simulated using a single permeability domain with Richards equation.

Two sub-models are available for lateral water transfer processes (overland flow, sub-surface flow, etc.), either "TETIS" (similar to Francés et al. 2007) or "SERGHEI" (Caviedes-Voullième et al. 2023).

IMPORTANT: medfateland needs to be compiled along with SERGHEI model in order to launch simulations with using this distributed hydrological model.

When running fordyn_land, the input 'sf' object has to be in a Universal Transverse Mercator (UTM) coordinate system (or any other projection using meters as length unit) for appropriate behavior of dispersal sub-model.

Parallel computation is only recommended for watersheds with large number of grid cells (e.g. > 10,000 when using transpirationMode = "granier"). In watershed with a small number of cells, parallel computation can result in larger processing times than sequential computation, due to the communication overload.

When dealing with large data sets, weather data included in the 'sf' object will likely be very data hungry. In those cases, it is recommended to resort on weather interpolation (see spwb_spatial). Weather interpolation can be done using a coarser resolution than that of raster 'r', by changing the watershed control parameter called 'weather_aggregation_factor' (see default_watershed_control).

Value

Functions spwb_land, growth_land and fordyn_land return a list of class of the same name as the function with the following elements:

Author(s)

Miquel De Cáceres Ainsa, CREAF.

Maria González-Sanchís, Universitat Politecnica de Valencia.

Daniel Caviedes-Voullième, Forschungszentrum Julich.

Mario Morales-Hernández, Universidad de Zaragoza.

References

Francés, F., Vélez, J.I. & Vélez, J.J. (2007). Split-parameter structure for the automatic calibration of distributed hydrological models. Journal of Hydrology, 332, 226–240.

Caviedes-Voullième, D., Morales-Hernández, M., Norman, M.R. & Ogzen-Xian, I. (2023). SERGHEI (SERGHEI-SWE) v1.0: a performance-portable high-performance parallel-computing shallow-water solver for hydrology and environmental hydraulics. Geoscientific Model Development, 16, 977-1008.

See Also

default_watershed_control, initialize_landscape, spwb_land_day, spwb_day, growth_day, spwb_spatial, fordyn_spatial, dispersal

Examples


# Load example watershed data
data("example_watershed")

# Set crop factor 
example_watershed$crop_factor <- NA
example_watershed$crop_factor[example_watershed$land_cover_type=="agriculture"] <- 0.75

# Set request for daily model results in cells number 3, 6 (outlet) and 9
example_watershed$result_cell <- FALSE
example_watershed$result_cell[c(3,6,9)] <- TRUE

# Get bounding box to determine limits
b <- sf::st_bbox(example_watershed)
b

# Define a raster topology, using terra package, 
# with the same CRS as the watershed. In this example cells have 100 m side.
# Coordinates in the 'sf' object are assumed to be cell centers
r <-terra::rast(xmin = 401380, ymin = 4671820, xmax = 402880, ymax = 4672620, 
                nrow = 8, ncol = 15, crs = "epsg:32631")

# Load example meteo data frame from package meteoland
data("examplemeteo")
  
# Load default medfate parameters
data("SpParamsMED")
  
# Set simulation period
dates <- seq(as.Date("2001-01-01"), as.Date("2001-03-31"), by="day")

# Watershed control parameters (TETIS model; Frances et al. 2007)
ws_control <- default_watershed_control("tetis")

# Launch simulations 
res <- spwb_land(r, example_watershed, SpParamsMED, examplemeteo, 
                 dates = dates, summary_frequency = "month",
                 watershed_control = ws_control)
                 
# Print a summary of water balance components
summary(res)

# Option 'simplify = TRUE' in initialization, may be useful to speed up calculations
example_simplified <- initialize_landscape(example_watershed, SpParams = SpParamsMED,
                                           local_control = defaultControl(soilDomains = "single"), 
                                           simplify = TRUE)

# Launch simulations over simplified landscape (should be considerably faster)
res_simplified <- spwb_land(r, example_simplified, SpParamsMED, examplemeteo, 
                            dates = dates, summary_frequency = "month",
                            watershed_control = ws_control)



[Package medfateland version 2.4.5 Index]