Package 'serofoi'

Title: Estimates the Force-of-Infection of a Given Pathogen from Population Based Seroprevalence Studies
Description: Estimate time or age varying Force-of-Infection from population based seroprevalence studies using a bayesian framework.
Authors: Zulma M. Cucunubá [aut, cre] , Nicolás T. Domínguez [aut] , Ben Lambert [aut], Pierre Nouvellet [aut], Miguel Gámez [ctb], Geraldine Gómez [ctb] , Jaime A. Pavlich-Mariscal [ctb]
Maintainer: Zulma M. Cucunubá <[email protected]>
License: MIT + file LICENSE
Version: 1.0.2
Built: 2024-12-11 10:33:01 UTC
Source: https://github.com/epiverse-trace/serofoi

Help Index


The 'serofoi' package.

Description

A DESCRIPTION OF THE PACKAGE

Author(s)

Maintainer: Zulma M. Cucunubá [email protected] (ORCID)

Authors:

Other contributors:

  • Miguel Gámez [contributor]

  • Geraldine Gómez (ORCID) [contributor]

  • Jaime A. Pavlich-Mariscal (ORCID) [contributor]

References

Stan Development Team (NA). RStan: the R interface to Stan. R package version 2.26.22. https://mc-stan.org

#' @keywords internal


Add bins based on age intervals.

Description

It generates a new column 'group' in the survey_features dataframe, representing the group interval for each row based on the age_min and age_max columns.

Usage

add_age_bins(survey_features)

Arguments

survey_features

A dataframe containing age_min and age_max columns representing the minimum and maximum age boundaries for each group.

Value

A dataframe with an additional 'group' column representing the group interval for each row based on the age_min and age_max columns.


Adds age group marker to serosurvey

Description

Adds age group marker to serosurvey

Usage

add_age_group_to_serosurvey(serosurvey)

Arguments

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group


Builds stan data for sampling depending on the selected model

Description

Builds stan data for sampling depending on the selected model

Usage

build_stan_data(
  serosurvey,
  model_type = "constant",
  foi_prior = sf_uniform(),
  foi_index = NULL,
  is_log_foi = FALSE,
  foi_sigma_rw = sf_none(),
  is_seroreversion = FALSE,
  seroreversion_prior = sf_none()
)

Arguments

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

model_type

Type of the model. Either "constant", "age" or "time"

foi_prior

Force-of-infection distribution specified by means of the helper functions. Currently available options are:

sf_normal

Function to set normal distribution priors

sf_uniform

Function to set uniform distribution priors

foi_index

Integer vector specifying the age-groups for which force-of-infection values will be estimated. It can be specified by means of get_foi_index

is_log_foi

Boolean to set logarithmic scale in the FOI

foi_sigma_rw

Prior distribution for the standard deviation of the force-of-infection. Currently available options are:

sf_normal

Function to set normal distribution prior. Available for time models in the log-scale

sf_cauchy

Function to set Cauchy distribution prior. Available for time models in regular scale.

is_seroreversion

Boolean specifying whether to include seroreversion rate estimation in the model

seroreversion_prior

seroreversion distribution specified by means of the helper functions. Currently available options are:

sf_normal

Function to set normal distribution priors

sf_uniform

Function to set uniform distribution priors

sf_none

Function to set no prior distribution

Value

List with necessary data for sampling the specified model


Chagas seroprevalence data in serofoi

Description

Datasets that measure the seroprevalence of IgG antibodies against Trypanosoma cruzi infection in rural areas of Colombia corresponding to a serosurvey conducted in 2012 for a rural indigenous community known to have long-term endemic transmission, where some control interventions have taken place over the years.

Usage

data(chagas2012)

Format

chagas2012

A ⁠<data.frame>⁠ with 4 rows and 5 columns:

survey_year

Year in which the serosurvey was conducted

n_sample

Number of collected samples per age group

n_seropositive

Number of positive samples per age group

age_min

Age group minimal age

age_max

Age group maximal age

Examples

data(chagas2012)

Chikungunya seroprevalence data in serofoi

Description

Datasets that measure the seroprevalence of IgG antibodies against the Chikungunya virus conducted in Bahia, Brazil in October-December 2015 by Dias et al. (2018). The survey was conducted immediately after a large Chikungunya epidemic in the area.

Usage

data(chik2015)

Format

chik2015

A ⁠<data.frame>⁠ with 4 rows and 5 columns:

survey_year

Year in which the serosurvey was conducted

n_sample

Number of collected samples per age group

n_seropositive

Number of positive samples per age group

age_min

Age group minimal age

age_max

Age group maximal age

Examples

data(chik2015)

Extracts central estimates from stan_fit object for specified parameter

Description

Extracts central estimates from stan_fit object for specified parameter

Usage

extract_central_estimates(
  seromodel,
  serosurvey,
  alpha = 0.05,
  par_name = "foi_vector"
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

alpha

1 - alpha indicates the credibility level to be used

par_name

String specifying the parameter to be extracted from seromodel

Value

A dataframe with the following columns

median

Median of the samples computed as the 0.5 quantile

lower

Lower quantile alpha

upper

Upper quantile 1 - alpha

Examples

data(veev2012)
seromodel <- fit_seromodel(veev2012, iter = 100)
central_estimates <- extract_central_estimates(
  seromodel,
  veev2012,
  par_name = "foi"
)

Runs specified stan model for the force-of-infection

Description

Runs specified stan model for the force-of-infection

Usage

fit_seromodel(
  serosurvey,
  model_type = "constant",
  is_log_foi = FALSE,
  foi_prior = sf_normal(),
  foi_sigma_rw = sf_none(),
  foi_index = NULL,
  foi_init = NULL,
  is_seroreversion = FALSE,
  seroreversion_prior = sf_normal(),
  ...
)

Arguments

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

model_type

Type of the model. Either "constant", "age" or "time"

is_log_foi

Boolean to set logarithmic scale in the FOI

foi_prior

Force-of-infection distribution specified by means of the helper functions. Currently available options are:

sf_normal

Function to set normal distribution priors

sf_uniform

Function to set uniform distribution priors

foi_sigma_rw

Prior distribution for the standard deviation of the force-of-infection. Currently available options are:

sf_normal

Function to set normal distribution prior. Available for time models in the log-scale

sf_cauchy

Function to set Cauchy distribution prior. Available for time models in regular scale.

foi_index

Integer vector specifying the age-groups for which force-of-infection values will be estimated. It can be specified by means of get_foi_index

foi_init

Initialization function for sampling. If null, default is chosen depending on the foi-scale of the model

is_seroreversion

Boolean specifying whether to include seroreversion rate estimation in the model

seroreversion_prior

seroreversion distribution specified by means of the helper functions. Currently available options are:

sf_normal

Function to set normal distribution priors

sf_uniform

Function to set uniform distribution priors

sf_none

Function to set no prior distribution

...

Additional parameters for rstan

Value

stan_fit object with force-of-infection and seroreversion (when applicable) samples

Examples

data(chagas2012)
seromodel <- fit_seromodel(
  serosurvey = chagas2012,
  model_type = "time",
  foi_index = data.frame(
    year = 1935:2011,
    foi_index = c(rep(1, 46), rep(2, 31))
  ),
  iter = 100
)

Generate random sample sizes for each age group.

Description

This function generates random sample sizes for each age group based on the overall sample size and the distribution of individuals across age groups. It uses multinomial sampling to allocate the total sample size to each age group proportionally.

Usage

generate_random_sample_sizes(survey_df_long)

Arguments

survey_df_long

A dataframe with columns 'age', 'group' and 'overall_sample_size'.

Value

A dataframe with random sample sizes generated for each age based on the overall sample size.


Construct age-group variable from age column

Description

Generates age intervals of length step in the interval spanned by age_min and age_max in a serosurvey. In cases where max(age_max)%%(step+1)!=0, the last age interval is truncated and will have a different length than the others.

Usage

get_age_intervals(serosurvey, step)

Arguments

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

step

step used to split the age interval

Value

Serosurvey with addition factor variable grouping age_intervals. The interval is taken as closed to the right and to the left.


Generates force-of-infection indexes for heterogeneous age groups

Description

Generates a list of integers indexing together the time/age intervals for which FOI values will be estimated in fit_seromodel. The max value in foi_index corresponds to the number of FOI values to be estimated when sampling. The serofoi approach to fitting serological data currently supposes that FOI is piecewise-constant across either groups of years or ages, and this function creates a Data Frame that communicates this grouping to the Stan model

Usage

get_foi_index(serosurvey, group_size, model_type)

Arguments

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

group_size

Age groups size

model_type

Type of the model. Either "age" or "time"

Value

A Data Frame which describes the grouping of years or ages (dependent on model) into pieces within which the FOI is assumed constant when performing model fitting. A single FOI value will be estimated for ages/years assigned with the same index

Examples

data(chagas2012)
foi_index <- get_foi_index(chagas2012, group_size = 25, model_type = "time")

Generate random sample sizes using multinomial sampling.

Description

This function generates random sample sizes for each age group using multinomial sampling. It takes the total sample size and the number of age groups as input and returns a vector of sample sizes for each age group.

Usage

multinomial_sampling_group(n_sample, n_ages)

Arguments

n_sample

The total sample size to be distributed among age groups.

n_ages

The number of age groups.

Value

A vector containing random sample sizes for each age group.


Plots force-of-infection central estimates

Description

Plots force-of-infection central estimates

Usage

plot_foi_estimates(
  seromodel,
  serosurvey,
  alpha = 0.05,
  foi_df = NULL,
  foi_max = NULL,
  size_text = 11
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

alpha

1 - alpha indicates the credibility level to be used

foi_df

Dataframe with columns

year/age

Year/Age (depending on the model)

foi

Force-of-infection values by year/age

foi_max

Max force-of-infection value for plotting

size_text

Size of text for plotting (base_size in ggplot2)

Value

ggplot object with estimated force-of-infection

Examples

data(chagas2012)
seromodel <- fit_seromodel(
  serosurvey = chagas2012,
  model_type = "time",
  foi_index = data.frame(
    year = 1935:2011,
    foi_index = c(rep(1, 46), rep(2, 31))
  ),
  iter = 100,
  chains = 2
)
plot_foi_estimates(seromodel, chagas2012)

Plot r-hats convergence criteria for the specified model

Description

Plot r-hats convergence criteria for the specified model

Usage

plot_rhats(seromodel, serosurvey, par_name = "foi_expanded", size_text = 11)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

par_name

String specifying the parameter to be extracted from seromodel

size_text

Size of text for plotting (base_size in ggplot2)

Value

ggplot object showing the r-hats of the model to be compared with the convergence criteria (horizontal dashed line)

Examples

data(chagas2012)
seromodel <- fit_seromodel(
  serosurvey = chagas2012,
  model_type = "time",
  foi_index = data.frame(
    year = 1935:2011,
    foi_index = c(rep(1, 46), rep(2, 31))
  ),
  iter = 100,
  chains = 2
)
plot_rhats(seromodel, chagas2012)

Visualise results of the provided model

Description

Visualise results of the provided model

Usage

plot_seromodel(
  seromodel,
  serosurvey,
  alpha = 0.05,
  bin_serosurvey = FALSE,
  bin_step = 5,
  foi_df = NULL,
  foi_max = NULL,
  loo_estimate_digits = 1,
  central_estimate_digits = 2,
  seroreversion_digits = 2,
  rhat_digits = 2,
  size_text = 11
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

alpha

1 - alpha indicates the credibility level to be used

bin_serosurvey

If TRUE, serodata is binned by means of prepare_bin_serosurvey. Otherwise, age groups are kept as originally input.

bin_step

Integer specifying the age groups bin size to be used when bin_serosurvey is set to TRUE.

foi_df

Dataframe with columns

year/age

Year/Age (depending on the model)

foi

Force-of-infection values by year/age

foi_max

Max force-of-infection value for plotting

loo_estimate_digits

Number of loo estimate digits

central_estimate_digits

Number of central estimate digits

seroreversion_digits

Number of seroreversion rate digits

rhat_digits

Number of rhat estimate digits

size_text

Size of text for plotting (base_size in ggplot2)

Value

seromodel summary plot

Examples

data(veev2012)
seromodel <- fit_seromodel(veev2012, iter = 100)
plot_seromodel(seromodel, veev2012)

Plot seroprevalence estimates on top of the serosurvey

Description

Plot seroprevalence estimates on top of the serosurvey

Usage

plot_seroprevalence_estimates(
  seromodel,
  serosurvey,
  alpha = 0.05,
  size_text = 11,
  bin_serosurvey = FALSE,
  bin_step = 5
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

alpha

1 - alpha indicates the credibility level to be used

size_text

Size of text for plotting (base_size in ggplot2)

bin_serosurvey

If TRUE, serodata is binned by means of prepare_bin_serosurvey. Otherwise, age groups are kept as originally input.

bin_step

Integer specifying the age groups bin size to be used when bin_serosurvey is set to TRUE.

Value

ggplot object with seroprevalence estimates and serosurveys plots

Examples

data(veev2012)
seromodel <- fit_seromodel(veev2012, iter = 100)
plot_seroprevalence_estimates(seromodel, veev2012)

Plots seroprevalence from the given serosurvey

Description

Plots seroprevalence from the given serosurvey

Usage

plot_serosurvey(
  serosurvey,
  size_text = 11,
  bin_serosurvey = FALSE,
  bin_step = 5
)

Arguments

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

size_text

Size of text for plotting (base_size in ggplot2)

bin_serosurvey

If TRUE, serodata is binned by means of prepare_bin_serosurvey. Otherwise, age groups are kept as originally input.

bin_step

Integer specifying the age groups bin size to be used when bin_serosurvey is set to TRUE.

Value

ggplot object with seroprevalence plot

Examples

# Chikungunya example serosurvey
data(chik2015)
plot_serosurvey(chik2015)

# VEEV example serosurvey
data(veev2012)
plot_serosurvey(veev2012)

# Chagas disease example serosurvey
data(chagas2012)
plot_serosurvey(chagas2012, bin_serosurvey = TRUE)

Plots model summary

Description

Plots model summary

Usage

plot_summary(
  seromodel,
  serosurvey,
  loo_estimate_digits = 1,
  central_estimate_digits = 2,
  rhat_digits = 2,
  size_text = 11
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

loo_estimate_digits

Number of loo estimate digits

central_estimate_digits

Number of central estimate digits

rhat_digits

Number of rhat estimate digits

size_text

Size of text for plotting (base_size in ggplot2)

Value

ggplot object with a summary of the specified model

Examples

data(veev2012)
seromodel <- fit_seromodel(veev2012, iter = 100)
plot_summary(seromodel, veev2012)

Prepares serosurvey for plotting

Description

Adds seroprevalence values with corresponding binomial confidence interval

Usage

prepare_serosurvey_for_plotting(serosurvey, alpha = 0.05)

Arguments

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

alpha

1 - alpha indicates the confidence level to be used

Value

serosurvey with additional columns:

seroprev

Seroprevalence computed as the proportion of positive cases n_seropositive in the number of samples n_sample for each age group

seroprev_lower

Lower limit of the binomial confidence interval of seroprev

seroprev_upper

Upper limit of the binomial confidence interval of seroprev


Computes the probability of being seropositive when FOIs vary by age

Description

Computes the probability of being seropositive when FOIs vary by age

Usage

probability_exact_age_varying(ages, fois, seroreversion_rate = 0)

Arguments

ages

Integer indicating the ages of the exposed cohorts

fois

Numeric atomic vector corresponding to the age-varying force-of-infection to simulate from

seroreversion_rate

Non-negative seroreversion rate. Default is 0.

Value

vector of probabilities of being seropositive for age-varying FoI including seroreversion (ordered from youngest to oldest individuals)


Computes the probability of being seropositive when FOIs vary by time

Description

Computes the probability of being seropositive when FOIs vary by time

Usage

probability_exact_time_varying(years, fois, seroreversion_rate = 0)

Arguments

years

Integer indicating the years covering the birth ages of the sample

fois

Numeric atomic vector corresponding to the age-varying force-of-infection to simulate from

seroreversion_rate

Non-negative seroreversion rate. Default is 0.

Value

vector of probabilities of being seropositive for age-varying FoI including seroreversion (ordered from youngest to oldest individuals)


Generate probabilities of seropositivity by age based on an age-and-time varying FOI model.

Description

This function calculates the probabilities of seropositivity by age based on an age-and-time-varying FOI model. It takes into account the FOI and the rate of seroreversion.

Usage

probability_seropositive_age_and_time_model_by_age(foi, seroreversion_rate)

Arguments

foi

A dataframe containing the force of infection (FOI) values for different ages. It should have three columns: 'year', 'age' and 'foi'.

seroreversion_rate

A non-negative numeric value representing the rate of seroreversion.

Value

A dataframe with columns 'age' and 'seropositivity'.


Generate probabilities of seropositivity by age based on an age-varying FOI model.

Description

This function calculates the probabilities of seropositivity by age based on an age-varying FOI model. It takes into account the FOI and the rate of seroreversion.

Usage

probability_seropositive_age_model_by_age(foi, seroreversion_rate)

Arguments

foi

A dataframe containing the force of infection (FOI) values for different ages. It should have two columns: 'age' and 'foi'.

seroreversion_rate

A non-negative numeric value representing the rate of seroreversion.

Value

A dataframe with columns 'age' and 'seropositivity'.


Generate probabilities of seropositivity by age based user's choice of model.

Description

This function generates seropositivity probabilities based on either a time-varying FOI model, an age-varying FOI model, or an age-and-time-varying FOI model. In all cases, it is possible to optionally include seroreversion.

Usage

probability_seropositive_by_age(model, foi, seroreversion_rate = 0)

Arguments

model

A string specifying the model type which can be one of 'age', 'time', 'age-time'.

foi

A dataframe containing the force of infection (FOI) values. For time-varying models the columns should be 'year', 'foi'. For age-varying models the columns should be 'age', 'foi'. For age-and-time-varying models the columns should be 'age', 'time', 'foi'.

seroreversion_rate

A non-negative value determining the rate of seroreversion (per year). Default is 0.

Value

A dataframe with columns 'age' and 'seropositivity'.

Examples

probability_seropositive_by_age(
  model = "age",
  foi = data.frame(
    age = 1:80,
    foi = rep(0.01, 80)
  )
)

Generate probabilities of seropositivity by age based on a general FOI model.

Description

This function calculates the probabilities of seropositivity by age based on an abstract model of the serocatalytic system.

Usage

probability_seropositive_general_model_by_age(
  construct_A_fn,
  calculate_seropositivity_function,
  initial_conditions,
  max_age,
  ...
)

Arguments

construct_A_fn

A function that constructs a matrix that defines the multiplier term in the linear ODE system.

calculate_seropositivity_function

A function which takes the state vector and returns the seropositive fraction.

initial_conditions

The initial state vector proportions for each birth cohort.

max_age

The maximum age to simulate seropositivity for.

...

Additional parameters for sum_of_A

Value

A dataframe with columns 'age' and 'seropositivity'.

Examples

# define age- and time-specific multipliers
foi_df_time <- data.frame(
  year = seq(1946, 2025, 1),
  foi = c(rep(0, 40), rep(1, 40))
)

foi_df_age <- data.frame(
  age = 1:80,
  foi = 2 * dlnorm(1:80, meanlog = 3.5, sdlog = 0.5)
)

u <- foi_df_age$foi
v <- foi_df_time$foi

# function to construct A matrix for one piece
construct_A <- function(t, tau, u, v) {
  u_bar <- u[t - tau]
  v_bar <- v[t]

  A <- diag(-1, ncol = 12, nrow = 12)
  A[row(A) == (col(A) + 1)] <- 1
  A[1, 1] <- -u_bar * v_bar
  A[2, 1] <- u_bar * v_bar
  A[12, 12] <- 0

  A
}

# determines the sum of seropositive compartments of those still alive
calculate_seropositivity_fn <- function(Y) {
  sum(Y[2:11]) / (1 - Y[12])
}

# initial conditions in 12D state vector
initial_conditions <- rep(0, 12)
initial_conditions[1] <- 1

# calculate probability
seropositive_hiv <- probability_seropositive_general_model_by_age(
  construct_A,
  calculate_seropositivity_fn,
  initial_conditions,
  max_age = 80,
  u,
  v
)

Generate probabilities of seropositivity by age based on a time-varying FOI model.

Description

This function calculates the probabilities of seropositivity by age based on a time-varying FOI model. It takes into account the FOI and the rate of seroreversion.

Usage

probability_seropositive_time_model_by_age(foi, seroreversion_rate)

Arguments

foi

A dataframe containing the force of infection (FOI) values for different years. It should have two columns: 'year' and 'foi'.

seroreversion_rate

A non-negative numeric value representing the rate of seroreversion.

Value

A dataframe with columns 'age' and 'seropositivity'.


Generate random sample sizes for each individual age based on survey features.

Description

This function generates random sample sizes for each individual age based on the provided survey features. It first creates age bins, assigns each individual in the survey features to an age bin, calculates the overall sample size by group, and then generates random sample sizes for each age group. Finally, it returns a dataframe with the random sample sizes for each individual age.

Usage

sample_size_by_individual_age_random(survey_features)

Arguments

survey_features

A dataframe containing information about individuals' age ranges and sample sizes.

Value

A dataframe with random sample sizes generated for each individual age based on the provided survey features.


Sets initialization function for sampling

Description

Sets initialization function for sampling

Usage

set_foi_init(foi_init, is_log_foi, foi_index)

Arguments

foi_init

Initialization function for sampling. If null, default is chosen depending on the foi-scale of the model

is_log_foi

Boolean to set logarithmic scale in the FOI

foi_index

Integer vector specifying the age-groups for which force-of-infection values will be estimated. It can be specified by means of get_foi_index

Examples

data(chagas2012)
foi_index <- get_foi_index(chagas2012, group_size = 5, model_type = "age")
foi_init <- set_foi_init(
  foi_init = NULL,
  is_log_foi = FALSE,
  foi_index = foi_index
)

Set stan data defaults for sampling

Description

Set stan data defaults for sampling

Usage

set_stan_data_defaults(stan_data, is_log_foi = FALSE, is_seroreversion = FALSE)

Arguments

stan_data

List to be passed to rstan

is_log_foi

Boolean to set logarithmic scale in the FOI

is_seroreversion

Boolean specifying whether to include seroreversion rate estimation in the model

Value

List with default values of stan data for sampling


Sets Cauchy distribution parameters for sampling

Description

Sets Cauchy distribution parameters for sampling

Usage

sf_cauchy(location = 0, scale = 1)

Arguments

location

Location of the Cauchy distribution

scale

Scale of the Cauchy distribution

Value

List with specified statistics and name of the model

Examples

my_prior <- sf_cauchy()

Sets empty distribution

Description

Sets empty distribution

Usage

sf_none()

Sets normal distribution parameters for sampling

Description

Sets normal distribution parameters for sampling

Usage

sf_normal(mean = 0, sd = 1)

Arguments

mean

Mean of the normal distribution

sd

Standard deviation of the normal distribution

Value

List with specified statistics and name of the model

Examples

my_prior <- sf_normal()

Sets uniform distribution parameters for sampling

Description

Sets uniform distribution parameters for sampling

Usage

sf_uniform(min = 0, max = 10)

Arguments

min

Minimum value of the random variable of the uniform distribution

max

Maximum value of the random variable of the uniform distribution

Value

List with specified statistics and name of the model

Examples

my_prior <- sf_uniform()

Simulate serosurvey data based on various FOI models.

Description

This function generates binned serosurvey data based on either a time-varying FOI model, an age-varying FOI model, or an age-and-time-varying FOI model. In all cases, it is possible to optionally include seroreversion. This function allows construction of serosurveys with binned age groups, and it generates uncertainty in the distribution of a sample size within an age bin through multinomial sampling.

Usage

simulate_serosurvey(model, foi, survey_features, seroreversion_rate = 0)

Arguments

model

A string specifying the model type which can be one of 'age', 'time', 'age-time'.

foi

A dataframe containing the force of infection (FOI) values. For time-varying models the columns should be 'year', 'foi'. For age-varying models the columns should be 'age', 'foi'. For age-and-time-varying models the columns should be 'age', 'time', 'foi'.

survey_features

A dataframe containing information about the binned age groups and sample sizes for each. It should contain columns: 'age_min', 'age_max', 'n_sample'.

seroreversion_rate

A non-negative value determining the rate of seroreversion (per year). Default is 0.

Value

A dataframe with simulated serosurvey data, including age group information, overall sample sizes, the number of seropositive individuals, and other survey features.

Examples

# time-varying model
foi_df <- data.frame(
  year = seq(1990, 2009, 1),
  foi = rnorm(20, 0.1, 0.01)
)
survey_features <- data.frame(
  age_min = c(1, 3, 15),
  age_max = c(2, 14, 20),
  n_sample = c(1000, 2000, 1500))
serosurvey <- simulate_serosurvey(
model = "time",
foi = foi_df,
survey_features = survey_features)

# age-varying model
foi_df <- data.frame(
  age = seq(1, 20, 1),
  foi = rnorm(20, 0.1, 0.01)
)
survey_features <- data.frame(
  age_min = c(1, 3, 15),
  age_max = c(2, 14, 20),
  n_sample = c(1000, 2000, 1500))
serosurvey <- simulate_serosurvey(
model = "age",
foi = foi_df,
survey_features = survey_features)

# age-and-time varying model
foi_df <- expand.grid(
  year = seq(1990, 2009, 1),
  age = seq(1, 20, 1)
)
foi_df$foi <- rnorm(20 * 20, 0.1, 0.01)
survey_features <- data.frame(
  age_min = c(1, 3, 15),
  age_max = c(2, 14, 20),
  n_sample = c(1000, 2000, 1500))
serosurvey <- simulate_serosurvey(
model = "age-time",
foi = foi_df,
survey_features = survey_features)

Simulate serosurvey data based on an age-and-time-varying FOI model.

Description

This function generates binned serosurvey data based on an age-and-time-varying FOI model, optionally including seroreversion. This function allows construction of serosurveys with binned age groups, and it generates uncertainty in the distribution of a sample size within an age bin through multinomial sampling.

Usage

simulate_serosurvey_age_and_time_model(
  foi,
  survey_features,
  seroreversion_rate = 0
)

Arguments

foi

A dataframe containing the force of infection (FOI) values for different ages. It should have two columns: 'year', 'age' and 'foi'.

survey_features

A dataframe containing information about the binned age groups and sample sizes for each. It should contain columns: 'age_min', 'age_max', 'n_sample'.

seroreversion_rate

A non-negative value determining the rate of seroreversion (per year). Default is 0.

Value

A dataframe with simulated serosurvey data, including age group information, overall sample sizes, the number of seropositive individuals, and other survey features.

Examples

# specify FOIs for each year
foi_df <- expand.grid(
  year = seq(1990, 2009, 1),
  age = seq(1, 20, 1)
)
foi_df$foi <- rnorm(20 * 20, 0.1, 0.01)
survey_features <- data.frame(
  age_min = c(1, 3, 15),
  age_max = c(2, 14, 20),
  n_sample = c(1000, 2000, 1500))
serosurvey <- simulate_serosurvey_age_and_time_model(
foi_df, survey_features)

Simulate serosurvey data based on an age-varying FOI model.

Description

This function generates binned serosurvey data based on an age-varying FOI model, optionally including seroreversion. This function allows construction of serosurveys with binned age groups, and it generates uncertainty in the distribution of a sample size within an age bin through multinomial sampling.

Usage

simulate_serosurvey_age_model(foi, survey_features, seroreversion_rate = 0)

Arguments

foi

A dataframe containing the force of infection (FOI) values for different ages. It should have two columns: 'age' and 'foi'.

survey_features

A dataframe containing information about the binned age groups and sample sizes for each. It should contain columns: 'age_min', 'age_max', 'n_sample'.

seroreversion_rate

A non-negative value determining the rate of seroreversion (per year). Default is 0.

Value

A dataframe with simulated serosurvey data, including age group information, overall sample sizes, the number of seropositive individuals, and other survey features.

Examples

# specify FOIs for each year
foi_df <- data.frame(
  age = seq(1, 20, 1),
  foi = rnorm(20, 0.1, 0.01)
)
survey_features <- data.frame(
  age_min = c(1, 3, 15),
  age_max = c(2, 14, 20),
  n_sample = c(1000, 2000, 1500))
serosurvey <- simulate_serosurvey_age_model(
foi_df, survey_features)

Simulate serosurvey data based on general serocatalytic model.

Description

This simulation method assumes only that the model system can be written as a piecewise-linear ordinary differential equation system.

Usage

simulate_serosurvey_general_model(
  construct_A_fn,
  calculate_seropositivity_function,
  initial_conditions,
  survey_features,
  ...
)

Arguments

construct_A_fn

A function that constructs a matrix that defines the multiplier term in the linear ODE system.

calculate_seropositivity_function

A function which takes the state vector and returns the seropositive fraction.

initial_conditions

The initial state vector proportions for each birth cohort.

survey_features

A dataframe containing information about the binned age groups and sample sizes for each. It should contain columns: 'age_min', 'age_max', 'n_sample'.

...

Additional parameters for sum_of_A

Value

A dataframe with simulated serosurvey data, including age group information, overall sample sizes, the number of seropositive individuals, and other survey features.

Examples

foi_df_time <- data.frame(
  year = seq(1946, 2025, 1),
  foi = c(rep(0, 40), rep(1, 40))
)

foi_df_age <- data.frame(
  age = 1:80,
  foi = 2 * dlnorm(1:80, meanlog = 3.5, sdlog = 0.5)
)

# generate age and time dependent FOI from multipliers
foi_age_time <- expand.grid(
  year = foi_df_time$year,
  age = foi_df_age$age
) |>
  dplyr::left_join(foi_df_age, by = "age") |>
  dplyr::rename(foi_age = foi) |>
  dplyr::left_join(foi_df_time, by = "year") |>
  dplyr::rename(foi_time = foi) |>
  dplyr::mutate(foi = foi_age * foi_time) |>
  dplyr::select(-c("foi_age", "foi_time"))

# create survey features for simulating
max_age <- 80
n_sample <- 50
survey_features <- data.frame(
  age_min = seq(1, max_age, 5),
  age_max = seq(5, max_age, 5)) |>
  dplyr::mutate(n_sample = rep(n_sample, length(age_min))
  )

# simulate survey from age and time FOI
serosurvey <- simulate_serosurvey(
  model = "age-time",
  foi = foi_age_time,
  survey_features = survey_features
)

Simulate serosurvey data based on a time-varying FOI model.

Description

This function generates binned serosurvey data based on a time-varying FOI model, optionally including seroreversion. This function allows construction of serosurveys with binned age groups, and it generates uncertainty in the distribution of a sample size within an age bin through multinomial sampling.

Usage

simulate_serosurvey_time_model(foi, survey_features, seroreversion_rate = 0)

Arguments

foi

A dataframe containing the force of infection (FOI) values for different years. It should have two columns: 'year' and 'foi'.

survey_features

A dataframe containing information about the binned age groups and sample sizes for each. It should contain columns: 'age_min', 'age_max', 'n_sample'.

seroreversion_rate

A non-negative value determining the rate of seroreversion (per year). Default is 0.

Value

A dataframe with simulated serosurvey data, including age group information, overall sample sizes, the number of seropositive individuals, and other survey features.

Examples

# specify FOIs for each year
foi_df <- data.frame(
  year = seq(1990, 2009, 1),
  foi = rnorm(20, 0.1, 0.01)
)
survey_features <- data.frame(
  age_min = c(1, 3, 15),
  age_max = c(2, 14, 20),
  n_sample = c(1000, 2000, 1500))
serosurvey <- simulate_serosurvey_time_model(
foi_df, survey_features)

Summarise central estimate

Description

Summarise central estimate

Usage

summarise_central_estimate(
  seromodel,
  serosurvey,
  alpha,
  par_name = "seroreversion_rate",
  central_estimate_digits = 2
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

alpha

1 - alpha indicates the credibility level to be used

par_name

String specifying the parameter to be extracted from seromodel

central_estimate_digits

Number of central estimate digits

Value

Text summarising specified central estimate

Examples

data(veev2012)
seromodel <- fit_seromodel(veev2012, iter = 100)
summarise_central_estimate(
  seromodel,
  veev2012,
  alpha = 0.05,
  par_name = "foi"
  )

Extract specified loo estimate

Description

Extract specified loo estimate

Usage

summarise_loo_estimate(
  seromodel,
  par_loo_estimate = "elpd_loo",
  loo_estimate_digits = 2
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

par_loo_estimate

Name of the loo estimate to be extracted. Available options are:

"elpd_loo"

Expected log pointwise predictive density

"p_loo"

Effective number of parameters

"looic"

Leave-one-out cross-validation information criteria

For additional information refer to loo.

loo_estimate_digits

Number of loo estimate digits

Value

Text summarising specified loo estimate

Examples

data(veev2012)
seromodel <- fit_seromodel(veev2012, iter = 100)
summarise_loo_estimate(seromodel)

Summarise specified model

Description

Summarise specified model

Usage

summarise_seromodel(
  seromodel,
  serosurvey,
  alpha = 0.05,
  par_loo_estimate = "elpd_loo",
  loo_estimate_digits = 1,
  central_estimate_digits = 2,
  rhat_digits = 2
)

Arguments

seromodel

stan_fit object obtained from sampling a model with fit_seromode

serosurvey
survey_year

Year in which the survey took place (only needed to plot time models)

age_min

Floor value of the average between age_min and age_max

age_max

The size of the sample

n_sample

Number of samples for each age group

n_seropositive

Number of positive samples for each age group

alpha

1 - alpha indicates the credibility level to be used

par_loo_estimate

Name of the loo estimate to be extracted. Available options are:

"elpd_loo"

Expected log pointwise predictive density

"p_loo"

Effective number of parameters

"looic"

Leave-one-out cross-validation information criteria

For additional information refer to loo.

loo_estimate_digits

Number of loo estimate digits

central_estimate_digits

Number of central estimate digits

rhat_digits

Number of rhat estimate digits

Value

A list summarising the specified model

model_name

Name of the model

elpd

elpd and its standard deviation

foi

Estimated foi with credible interval (for 'constant' model)

foi_rhat

foi rhat value (for 'constant' model)

seroreversion_rate

Estimated seroreversion rate

seroreversion_rate_rhat

Seroreversion rate rhat value

Examples

data(veev2012)
seromodel <- fit_seromodel(veev2012, iter = 100)
summarise_seromodel(seromodel, veev2012)

Create a survey dataframe with per individual age information.

Description

Create a survey dataframe with per individual age information.

Usage

survey_by_individual_age(survey_features, age_df)

Arguments

survey_features

A dataframe containing information about age groups and sample sizes.

age_df

A dataframe containing 'age' and 'group'.

Value

A dataframe with overall sample sizes calculated by joining survey_features and age_df. This dataframe has columns including 'age' and 'overall_sample_size'.


Venezuelan Equine Encephalitis Virus (VEEV) seroprevalence data in serofoi

Description

Datasets that measure the seroprevalence of IgG antibodies against VEEV in a rural village in Panamá in 2012. carrera2020

Usage

data(veev2012)

Format

veev2012

A ⁠<data.frame>⁠ with 4 rows and 5 columns:

survey_year

Year in which the serosurvey was conducted

n_sample

Number of collected samples per age group

n_seropositive

Number of positive samples per age group

age_min

Age group minimal age

age_max

Age group maximal age

Examples

data(veev2012)