Multi-species models & predation

This vignette walks through a script that will generate a gadget3 model, explaining concepts along the way.

Code blocks that make up the gadget3 script are marked in blue, like this:

### Modelling maturity & sex with multiple stocks

When combined they will form a full model, see the appendix for the entire script.

Introduction

A multi-species model can be assembled mostly using concepts from previous vignettes. Within an individual species setup is identical, See previous vignettes, such as vignette('introduction-single-stock'), for how this works.

Initial setup & time-keeping is identical:

library(gadget3)
library(dplyr)

actions <- list()
area_names <- g3_areas(c("a"))

# Create time definitions ####################

actions_time <- list(
  g3a_time(
    1991L, 2020L,
    step_lengths = c(3L, 3L, 3L, 3L)),
  NULL)

actions <- c(actions, actions_time)

Prey stocks

Before we can define the predator, we need something for the predator to eat.

Prey can be fully simulated species as in previous vignettes, or as static “otherfood”.

For this example, we will define vendace with mature and immature substocks, as we did in vignette('multiple-substocks'):

# Create stock definition for vendace ####################
ven_imm <- g3_stock(c(species = "ven", 'imm'), seq(3, 18, 3)) |>
  g3s_age(1L, 5L) |>
  g3s_livesonareas(area_names["a"])

ven_mat <- g3_stock(c(species = "ven", 'mat'), seq(3, 21, 3)) |>
  g3s_age(1L, 10L) |>
  g3s_livesonareas(area_names["a"])
stocks_ven = list(imm = ven_imm, mat = ven_mat)

actions_ven_imm <- list(
  g3a_growmature(ven_imm,
    g3a_grow_impl_bbinom(
      maxlengthgroupgrowth = 4L ),
    # Add maturation
    maturity_f = g3a_mature_continuous(),
    output_stocks = list(ven_mat),
    transition_f = ~TRUE ),
  g3a_naturalmortality(ven_imm),
  g3a_initialconditions_normalcv(ven_imm),
  g3a_renewal_normalparam(ven_imm),
  g3a_age(ven_imm, output_stocks = list(ven_mat)),
  NULL)

actions_ven_mat <- list(
  g3a_growmature(ven_mat,
    g3a_grow_impl_bbinom(
      maxlengthgroupgrowth = 4L )),
  g3a_naturalmortality(ven_mat),
  g3a_initialconditions_normalcv(ven_mat),
  g3a_age(ven_mat),
  NULL)

actions_likelihood_ven <- list(
  g3l_understocking(stocks_ven, nll_breakdown = TRUE),
  NULL)

actions <- c(actions, actions_ven_imm, actions_ven_mat, actions_likelihood_ven)

Other food stocks

Any stomach content surveys will contain data on more preys than you will want to dynamically model.

To account for these we can include “otherfood” stocks, which are always available regardless of predation, as they are renewed to their default level at every step.

There are 2 sources of otherfood in our model, “herring” and “otherfood”:

# Create stock definition for otherfood ##################
her <- g3_stock(c(species = "her"), seq(0, 30, 10)) |> g3s_age(10, 20) |>
  g3s_livesonareas(area_names["a"])

otherfood <- g3_stock("other", 0) |>
  g3s_livesonareas(area_names["a"])
stocks_oth = list(her = her, other = otherfood)

actions_otherfood <- list(
    g3a_otherfood_normalcv(her),
    g3a_otherfood(otherfood) )

actions <- c(actions, actions_otherfood)

Herring is defined with length & age structure, and we use g3a_otherfood_normalcv() to define it’s abundance using a vonB curve, as we do with initialconditions. This means we will have her.t0 & her.Linf parameters to populate, as we do a dynamically modelled stock.

“otherfood” models an always-present level of biomass in the model. It doesn’t have any stock structure, so we don’t break it down further. Here we use g3a_otherfood(), which set a static abundance / mean-weight, via. the other.of_abund.199x, other.of_abund.step.x and other.of_meanwgt parameters. other.of_abund.199x allows abundance to vary year-on-year, and other.of_abund.step.x allows for seasonal variance.

Predator stocks

The predator in our model, the ringed seal, is divided into 2 substocks; males and females.

Again, setup within the ringed seals is identical to previous:

# Create stock definition for ringed seal ################
rin_male <- g3_stock(c(species = "rin", sex = "m"), seq(3, 18, 3)) |>
  g3s_age(0L, 5L) |>
  g3s_livesonareas(area_names["a"])

rin_female <- g3_stock(c(species = "rin", sex = "f"), seq(3, 18, 3)) |>
  g3s_age(0L, 5L) |>
  g3s_livesonareas(area_names["a"])
stocks_rin = list(male = rin_male, female = rin_female)

actions_rin_male <- list(
  g3a_growmature(rin_male,
    g3a_grow_impl_bbinom(
      maxlengthgroupgrowth = 4L )),
  g3a_naturalmortality(rin_male),
  g3a_initialconditions_normalcv(rin_male),
  g3a_renewal_normalparam(rin_male),
  g3a_age(rin_male),
  NULL)

actions_rin_female <- list(
  g3a_growmature(rin_female,
    g3a_grow_impl_bbinom(
      maxlengthgroupgrowth = 4L )),
  g3a_naturalmortality(rin_female),
  g3a_initialconditions_normalcv(rin_female),
  g3a_renewal_normalparam(rin_female),
  g3a_age(rin_female),
  NULL)

actions <- c(actions, actions_rin_male, actions_rin_female)

We could also set up ringed seals with maturity in an identical fashion to vendace, but this is sufficient for the current example.

Predation

Finally, we can add an action for ringed seals to eat the other stocks in the model. The process is nearly identical to adding a fleet to a model; first we generate some random observation data, then we add predation actions, finally we add likelihood components for our prey preference data:

# Predation by ringed seal ##############################
stocks_rin_prey <- c(stocks_ven, stocks_oth)  # Ringed seals will predate on all preys previously set up

stomach_prey_preference <- expand.grid(
    year = 1991:2020,
    step = 1:4,
    stock = c("ven", "her", "other"),
    number = 0 )

# Generate random preference data for each stock
for (s in levels(stomach_prey_preference$stock)) {
    stomach_prey_preference[stomach_prey_preference$stock == s, "number"] <- pmax(0, rnorm(
        sum(stomach_prey_preference$stock == s),
        mean = c(ven = 0.419, her = 0.216, other = 0.492)[[s]],
        sd = 0.1 ))
}

stomach_length_preference <- expand.grid(
    year = 1991:2020,
    step = 1:4,
    predator_length = seq(3, 18, 3),
    length = seq(3, 21, 3), 
    number = 0 )

# Generate random length preference data for each stock
lengths <- unique(stomach_length_preference$length)
lenpref <- table(cut(rnorm(900, mean=14, sd=1.5), c(lengths, Inf)))
for (i in seq_along(lengths)) {
    stomach_length_preference[stomach_length_preference$length == lengths[[i]], "number"] <- pmax(0, rnorm(
        sum(stomach_length_preference$length == lengths[[i]]),
        mean = lenpref[[i]],
        sd = 10 ))
}

actions_rin_predate <- list(
    g3a_predate(
        stocks_rin,
        stocks_rin_prey,
        suitabilities = list(
            # Otherfood uses a fixed suitability
            other = g3_suitability_constant(by_predator='species'),
            # Everything else will use g3_suitability_exponentiall50:
            g3_suitability_exponentiall50(
                alpha=g3_parameterized('alpha', by_stock='species', by_predator='species'),
                l50=g3_parameterized('l50', by_stock='species', by_predator='species') )),
# TODO: Temperature
        catchability_f = g3a_predate_catchability_predator() ),
    NULL )
actions_likelihood_rin_predate <- list(
    g3l_catchdistribution(
        "ppref.rin",
        obs_data = stomach_prey_preference,
        stocks = stocks_rin_prey,
        predators = stocks_rin,
        function_f = g3l_distribution_sumofsquares(),
        area_group = area_names,
        report = TRUE,
        nll_breakdown = TRUE ),
    g3l_catchdistribution(
        "lpref.rin.ven",
        obs_data = stomach_length_preference,
        stocks = stocks_ven,  # NB: Just vendace, not all prey
        predators = stocks_rin,
        function_f = g3l_distribution_sumofsquares(),
        area_group = area_names,
        report = TRUE,
        nll_breakdown = TRUE),
    NULL )

actions <- c(actions, actions_rin_predate, actions_likelihood_rin_predate)

The main differences between defining a fleet and predator are:

  • Instead of a single g3_fleet() we supply all the ringed seal g3_stocks() as a list. This will set up predation identically on males & females.
  • Instead of g3a_predate_catchability_totalfleet(), and providing landings data, we use g3a_predate_catchability_predator() to define the predator/prey relationship, equivalent to consumption in gadget2.
  • We slice our stomach data into 2 likelihood components, giving prey preference and length preference. Using g3l_distribution_sumofsquares() is equivalent to gadget2 scsimple comparison function.

g3a_predate_catchability_predator() will by default create parameters for the configuration settings in gadget2, including half-feeding & energy content. Look at the function documentation for the definitions of each:

names(attr(suppressWarnings(g3_to_r(actions_rin_predate)), 'parameter_template'))
##  [1] "her.rin.alpha"         "her.rin.l50"           "other.rin.suit"        "ven.rin.alpha"         "ven.rin.l50"           "her.energycontent"    
##  [7] "other.energycontent"   "ven_imm.energycontent" "ven_mat.energycontent" "rin_f.halffeeding"     "rin_f.consumption.m0"  "rin_f.consumption.m1" 
## [13] "rin_f.consumption.m2"  "rin_f.consumption.m3"  "rin_m.halffeeding"     "rin_m.consumption.m0"  "rin_m.consumption.m1"  "rin_m.consumption.m2" 
## [19] "rin_m.consumption.m3"

Fleet actions

We add a fleet for vendace, as we have in previous models by generating historical data, adding fleet predation and likelihood components to compare the model to generated observed output.

Note that we can use ven in the stock column to refer to both ven_imm and ven_mat:

# Fleet data for f_surv #################################

# Landings data: For each year/step/area
expand.grid(year = 1991:2020, step = 2, area = 'a', stock = 'ven') |>
    # Generate a random total landings by weight
    mutate(weight = rnorm(n(), mean = 1000, sd = 100)) |>
    # Assign result to landings_f_surv
    identity() -> landings_f_surv

# Length distribution data: Generate 100 random samples in each year/step/area
expand.grid(year = 1991:2020, step = 2, area = 'a', stock = 'ven', length = rep(NA, 100)) |>
  # Generate random lengths for these samples
  mutate(length = rnorm(n(), mean = 50, sd = 20)) |>
  # Save unagggregated data into ldist_f_surv.raw
  identity() -> ldist_f_surv.raw

# Aggregate .raw data
ldist_f_surv.raw |>
  # Group into length bins
  group_by(
      year = year,
      step = step,
      stock = stock,
      length = cut(length, breaks = c(seq(0, 80, 20), Inf), right = FALSE) ) |>
  # Report count in each length bin
  summarise(number = n(), .groups = 'keep') |>
  # Save into ldist_f_surv
  identity() -> ldist_f_surv

# Assume 5 * 5 samples in each year/step/area
expand.grid(year = 1991:2020, step = 2, area = 'a', stock = 'ven', age = rep(NA, 5), length = rep(NA, 5)) |>
  # Generate random lengths/ages for these samples
  mutate(length = rnorm(n(), mean = 50, sd = 20)) |>
  # Generate random whole numbers for age
  mutate(age = floor(runif(n(), min = 1, max = 5))) |>
  # Group into length/age bins
  group_by(
      year = year,
      step = step,
      age = age,
      stock = stock,
      length = cut(length, breaks = c(seq(0, 80, 20), Inf), right = FALSE) ) |>
  # Report count in each length bin
  summarise(number = n(), .groups = 'keep') ->
  aldist_f_surv

For more information on how this works, see vignette("incorporating-observation-data").

Our fleet is defined with the same set of actions as the single-species model:

# Create fleet definition for f_surv ####################
f_surv <- g3_fleet("f_surv") |> g3s_livesonareas(area_names["a"])

actions_f_surv <- list(
  g3a_predate_fleet(
    f_surv,
    stocks_ven,
    suitabilities = g3_suitability_exponentiall50(by_stock = 'species'),
    catchability_f = g3a_predate_catchability_totalfleet(
      g3_timeareadata("landings_f_surv", landings_f_surv, "weight", areas = area_names))),
  NULL)
actions_likelihood_f_surv <- list(
  g3l_catchdistribution(
    "ldist_f_surv",
    obs_data = ldist_f_surv,
    fleets = list(f_surv),
    stocks = stocks_ven,
    function_f = g3l_distribution_sumofsquares(),
    area_group = area_names,
    report = TRUE,
    nll_breakdown = TRUE),
  g3l_catchdistribution(
    "aldist_f_surv",
    obs_data = aldist_f_surv,
    fleets = list(f_surv),
    stocks = stocks_ven,
    function_f = g3l_distribution_sumofsquares(),
    area_group = area_names,
    report = TRUE,
    nll_breakdown = TRUE),
  NULL)

actions <- c(actions, actions_f_surv, actions_likelihood_f_surv)

Survey indices

As in other vignettes, we can add survey indices for each stock. Doing so for a predator is no different:

# Create abundance index for vendace ########################

# Generate random data
expand.grid(year = 1991:2020, step = 3, area = names(area_names)) |>
    # Fill in a weight column with total biomass for the year/step/area combination
    mutate(weight = runif(n(), min = 1e4, max = 1e5)) ->
    dist_si_ven

actions_likelihood_dist_ven <- list(

  g3l_abundancedistribution(
    "dist_si_ven",
    dist_si_ven,

    stocks = stocks_ven,
    function_f = g3l_distribution_surveyindices_log(alpha = NULL, beta = 1),
    area_group = area_names,
    report = TRUE,
    nll_breakdown = TRUE),
  NULL)

actions <- c(actions, actions_likelihood_dist_ven)

# Create abundance index for ringed seal ####################

# Generate random data
expand.grid(year = 1991:2020, step = 3, area = names(area_names)) |>
    # Fill in a weight column with total biomass for the year/step/area combination
    mutate(weight = runif(n(), min = 1e3, max = 1e4)) ->
    dist_si_rin

actions_likelihood_dist_rin <- list(

  g3l_abundancedistribution(
    "dist_si_rin",
    dist_si_rin,

    stocks = stocks_rin,
    function_f = g3l_distribution_surveyindices_log(alpha = NULL, beta = 1),
    area_group = area_names,
    report = TRUE,
    nll_breakdown = TRUE),
  NULL)

actions <- c(actions, actions_likelihood_dist_rin)

Creating model functions and Parameterization

At this point, we are ready to convert our model into code:

# Create model objective function ####################

model_code <- g3_to_tmb(c(actions, list(
    g3a_report_detail(actions),
    g3l_bounds_penalty(actions) )))

Now we should be configuring parameters based on the template:

attr(model_code, "parameter_template") |>

  g3_init_val("her*.Linf", max(g3_stock_def(her, "midlen")), spread = 0.2) |>
  g3_init_val("rin*.Linf", max(sapply(stocks_rin, function (s) max(g3_stock_def(s, "midlen")))), spread = 0.2) |>
  g3_init_val("ven*.Linf", max(sapply(stocks_ven, function (s) max(g3_stock_def(s, "midlen")))), spread = 0.2) |>
  g3_init_val("her.t0", g3_stock_def(her, "minage"), spread = 0.2) |>
  g3_init_val("rin*.t0", min(unlist(g3_stock_def(stocks_rin, "minage"))) - 0.8, spread = 0.2) |>
  g3_init_val("ven*.t0", min(unlist(g3_stock_def(stocks_ven, "minage"))) - 0.8, spread = 0.2) |>
  g3_init_val("*.K", 0.3, lower = 0.04, upper = 1.2) |>

  g3_init_val('other.of_abund.#', 1e5, lower = 0, upper = 1e8, optimise = FALSE) |>
  g3_init_val('other.of_abund.step.#', 1, lower = 0, upper = 1e8, optimise = FALSE) |>
  g3_init_val('other.of_meanwgt', 0.1, lower = 0.01, upper = 3) |>

  g3_init_val('her*.energycontent', 1, lower = 0.01, upper = 3) |>
  g3_init_val('ven*.energycontent', 1, lower = 0.01, upper = 3) |>
  g3_init_val('other.energycontent', 1, lower = 0.01, upper = 3) |>
  g3_init_val('rin*.halffeeding', 0.01, lower = 0.01, upper = 1e5) |>

  g3_init_val("*.rin.alpha", value = 0.5, lower = 0.01, upper = 3) |>
  g3_init_val("her.rin.l50", mean(g3_stock_def(her, "midlen")), lower = 0, upper = max(g3_stock_def(her, "midlen")) ) |>
  g3_init_val("ven.rin.l50", mean(unlist(g3_stock_def(stocks_ven, "midlen"))), lower = 0,
      upper = max(unlist(g3_stock_def(stocks_ven, "midlen"))) ) |>
  g3_init_val("ven.f_surv.l50", mean(unlist(g3_stock_def(stocks_ven, "midlen"))), lower = 0,
      upper = max(unlist(g3_stock_def(stocks_ven, "midlen"))) ) |>
  g3_init_val("other.rin.suit", value = 0.4, lower = 0.1, upper = 1) |>
  g3_init_val('rin_*.consumption.m0', 0.0004, lower = 0.00001, upper = 1e5) |>
  g3_init_val('rin_*.consumption.m1', 0) |>
  g3_init_val('rin_*.consumption.m2', 0) |>
  g3_init_val('rin_*.consumption.m3', 2.57, lower = 0.1, upper = 10) |>

  g3_init_val("ven*.mat.l50", mean(unlist(g3_stock_def(stocks_ven, "midlen"))), lower = 0,
      upper = max(unlist(g3_stock_def(stocks_ven, "midlen"))) ) |>
  g3_init_val("*.rec|init.scalar", 10, lower = 0.001, upper = 200) |>
  g3_init_val("*.init.#", 10, lower = 0.001, upper = 200) |>
  g3_init_val("*.rec.#", 100, lower = 1e-6, upper = 1000) |>
  g3_init_val("*.rec.sd", 5, lower = 4, upper = 20) |>
  g3_init_val("*.M.#", 0.15, lower = 0.001, upper = 1) |>
  g3_init_val("init.F", 0.5, lower = 0.1, upper = 1) |>
  g3_init_val("*.walpha", 0.01, optimise = FALSE) |>
  g3_init_val("*.wbeta", 3, optimise = FALSE) |>
  g3_init_val("*.*.alpha", 0.07, lower = 0.01, upper = 0.2) |>
  g3_init_val("*.bbin", 100, lower = 1e-05, upper = 1000) |>

  identity() -> params.in

Finally we are ready to use gadgetutils::g3_iterative() to optimise based on iterative reweighting:

# Optimise model ################################
obj.fn <- g3_tmb_adfun(model_code, params.in)

params.out <- gadgetutils::g3_iterative(getwd(),
    wgts = "WGTS",
    model = model_code,
    params.in = params.in,
    grouping = list(
        fleet = c("ldist_f_surv", "aldist_f_surv"),
        abund = c("dist_si_ven", "dist_si_rin") ),
    method = "BFGS",
    control = list(maxit = 1000, reltol = 1e-10),
    cv_floor = 0.05)

Once this has finished, we can view the output using gadgetplots::gadget_plots().

# Generate detailed report ######################
fit <- gadgetutils::g3_fit(model_code, params.out)
## Joining with `by = join_by(stock, fleet, length, lower, upper)`
## Warning: There were 480 warnings in `dplyr::mutate()`.
## The first warning was:
## ℹ In argument: `mortality = -log(1 - .data$number_consumed/.data$number)/step_size`.
## ℹ In group 1: `year = 1991`, `step = 1`, `area = "a"`, `stock = "her"`, `fleet = "rin_f"`.
## Caused by warning in `log()`:
## ! NaNs produced
## ℹ Run `dplyr::last_dplyr_warnings()` to see the 479 remaining warnings.
# TODO: Activate plots
#gadgetplots::gadget_plots(fit, "figs", file_type = "html")

Once finished, you can view the output in your web browser:

utils::browseURL("figs/model_output_figures.html")

Appendix: Full model script

For convenience, here is all the sections of the model script above joined together: