GeneDrive.ConstantTemperatureType
mutable struct ConstantTemperature <: Temperature end

Data for simulation that features a single constant temperature in °C. Generated internally.

Arguments

  • value::Float64: Temperature value in °C; re-used by model at every time step.
source
GeneDrive.DensityType
mutable struct Density{D <: DensityDependence}
    model::Type{D}
    param::Float64
end

Data for density dependence model.

Arguments

  • model::Type{D}: User-selected density dependence formulation. Includes LogisticDensity, LinearDensity, and NoDensity.
  • param::Float64: Calculated internally by init function.
source
GeneDrive.DriveType
mutable struct Drive{G <: Genotype}
    genotype::Type{G}
    likelihood_slice::Array{Float64,2}
    s::Float64
    τ::Array{Float64,2}
    ϕ::Float64
    ξ_m::Float64
    ξ_f::Float64
    ω::Float64
    β::Float64
    η::Float64
    wildtype::Int64
    modified::Int64
end

Data for individual genotypes.

Fields

  • genotype::Type{G}: Single genotype.
  • likelihood_slice::Array{Float64,2}: Offspring likelihoods for this genotype.
  • s::Float64: Multiplicative fertility modifier.
  • τ::Array{Float64,2}: Offspring viability.
  • ϕ::Float64: Male to female emergence ratio.
  • ξ_m::Float64: Male pupatory success.
  • ξ_f::Float64: Female pupatory success.
  • ω::Float64: Multiplicative adult mortality modifier.
  • β::Float64: Female fecundity.
  • η::Float64: Male mating fitness.
  • wildtype::Int64: Boolean demarcates homozygous recessive.
  • modified::Int64: Boolean demarcates homozygous modified.
source
GeneDrive.EggType
struct Egg <: LifeStage end

Egg life stage. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.EggDurationRossiType
mutable struct EggDurationRossi <: TemperatureResponse
    a::Float64
    b::Float64
    c::Float64
    d::Float64
    e::Float64
    f::Float64
end

Data for q_temperature_response. Species: AedesAegypti, LifeStage:, egg stage. Source: Rossi et al (2014).

source
GeneDrive.EggMortalityRossiType
mutable struct EggMortalityRossi <: TemperatureResponse
    a::Float64
    b::Float64
end

Data for μ_temperature_response. Species: AedesAegypti, LifeStage: Egg. Source: Rossi et al (2014).

Arguments

  • a::Float64: Death rate at low temperatures, validation range: (0, nothing)
  • b::Float64: Influence factor of temperature, validation range: (0, nothing)
source
GeneDrive.ExogenousInputsType
struct ExogenousInputs
    intervention::Dict{Symbol, Dict}
    temperature::Dict{Symbol, Float64}
    function ExogenousInputs(;
        intervention = Dict{Symbol, Dict}(),
        temperature = Dict{Symbol, Float64}())
        new(intervention,
        temperature)
    end
end

Data for perturbations to the simulated system, including (1) biological control interventions and (2) externally defined temperature series/shocks.

source
GeneDrive.FemaleType
struct Female <: LifeStage end

Female life stage. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.GeneticsType
mutable struct Genetics

    all_genotypes::Array{Drive{<:Genotype}}
    likelihood::Array{Float64, 3}
    S::Vector{Float64}
    Τ::Array{Float64,3}
    Φ::Vector{Float64}
    Ξ_m::Vector{Float64}
    Ξ_f::Vector{Float64}
    Ω::Vector{Float64}
    Β::Vector{Float64}
    Η::Vector{Float64}
    all_wildtypes::Vector{Int64}
    all_modified::Vector{Int64}

        function Genetics(all_genotypes::Array{Drive{<:Genotype}})

            gN = length(all_genotypes)
            likelihood = Array{Float64, 3}(undef, gN, gN, gN)
            S = Vector{Float64}(undef, gN)
            Τ = Array{Float64,3}(undef, gN, gN, gN)
            Φ = Vector{Float64}(undef, gN)
            Ξ_m = Vector{Float64}(undef, gN)
            Ξ_f = Vector{Float64}(undef, gN)
            Ω = Vector{Float64}(undef, gN)
            Β = Vector{Float64}(undef, gN)
            Η = Vector{Float64}(undef, gN)
            all_wildtypes = Vector{Int64}(undef, gN)
            all_modified = Vector{Int64}(undef, gN)

            for (index, g) in enumerate(all_genotypes)
                likelihood[:,:,index] = g.likelihood_slice
                S[index] = g.s
                Τ[:,:,index] = g.τ
                Φ[index] = g.ϕ
                Ξ_m[index] = g.ξ_m
                Ξ_f[index] = g.ξ_f
                Ω[index] = g.ω
                Β[index] = g.β
                Η[index] = g.η
                all_wildtypes[index] = g.wildtype
                all_modified[index] = g.modified
            end

            new(all_genotypes, likelihood, S, Τ, Φ, Ξ_m, Ξ_f, Ω, Β, Η,
            all_wildtypes, all_modified)

        end

end

Data for all genotypes in a population.

Fields

  • all_genotypes::Array{Drive{<:Genotype}}: All genotypes in a population.
  • likelihood::Array{Float64, 3}: Offspring likelihoods per genotype.
  • S::Vector{Float64}: Multiplicative fertility modifier genedata_splitdriveper genotype, applied to oviposition.
  • Τ::Array{Float64,3}: Offspring viability per genotype, applied to oviposition.
  • Φ::Vector{Float64}: Male to female emergence ratio per genotype.
  • Ξ_m::Vector{Float64}: Male pupatory success per genotype.
  • Ξ_f::Vector{Float64}: Female pupatory success per genotype.
  • Ω::Vector{Float64}: Multiplicative adult mortality modifier per genotype.
  • Β::Vector{Float64}: Female fecundity per genotype (count of eggs laid).
  • Η::Vector{Float64}: Male mating fitness per genotype.
  • all_wildtypes::Vector{Int64}: Collect homozygous wildtype booleans.
  • all_modified::Vector{Int64}: Collect homozygous modified booleans.
source
GeneDrive.LarvaType
struct Larva <: LifeStage end

Larva life stage. Also referred to as "nymph" life stage. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.MaleType
struct Male <: LifeStage end

Male life stage. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.NetworkType
struct Network
    name::Symbol
    nodes::DataStructures.OrderedDict{Symbol,Node}
    migration::DataStructures.OrderedDict{DataType},Array{Matrix{Float64},2}}
    locations_key_map
end

Data for multiple interconnected spatial nodes.

Arguments

  • name::Symbol: Name of network, usually location-relevant.
  • nodes::DataStructures.OrderedDict{Symbol,Node}: Dictionary of data for nodes included in network (metapopulation).
  • migration::DataStructures.OrderedDict{DataType}, Array{Matrix{Float64},2}}: Data defining species, genotype, and stage-specific transition rates. Entries default to zero.
  • locations_key_map: Mapping of node location information. Used to assign transition rates.
source
GeneDrive.NetworkMethod
Network(name::Symbol, nodes::DataStructures.OrderedDict{Symbol,Node}...)

Return Network instance containing specified metapopulation information.

source
GeneDrive.NoResponseType
mutable struct NoResponse <: TemperatureResponse
    baseline_value::Float64
end

Data for model without temperature response.

Arguments:

  • baseline_value::Float64: Literature-sourced (rather than dynamically calculated) vital rate parameter.
source
GeneDrive.NodeType
mutable struct Node
    name::Symbol
    organisms::DataStructures.OrderedDict{Type{<:Species}, Organism}
    temperature::Temperature
    location::Tuple{Float64,Float64}
end

Data for a single spatial node.

Arguments

  • name::Symbol: Name of node, usually location-relevant.
  • organisms::DataStructures.OrderedDict{Type{<:Species}, Organism}: Dictionary containing data for all species inhabiting node.
  • temperature::Temperature: Climatic specification for temperature.
  • location::Tuple{Float64,Float64}: Geographic location denoted by coordinates.
source
GeneDrive.OrganismType
mutable struct Organism{S <: Species}
    gene_data::Genetics
    all_stages::DataStructures.OrderedDict{Type{<:LifeStage}, Stage}
end

Generic data container for species-specific genetic and life stage information.

Arguments

  • gene_data::Genetics: Species and modification-specific genetic data.
  • all_stages::DataStructures.OrderedDict{Type{<:LifeStage}, Stage}: Dictionary of wildtype species-specific life stage data.
source
GeneDrive.ProportionalReleaseType
mutable struct ProportionalRelease <: Intervention
    node::Symbol
    organism::Type{<:Species}
    stage::Type{<:LifeStage}
    gene_index::Int64
    times::Vector{Float64}
    proportion::Float64
    callbacks::Vector
    adults_counting::String
end

Data for biological control interventions predicated on releasing modified organisms.

Arguments

  • node::Symbol: Node where releases occur.
  • organism::Type{<:Species}: Species of organism being released.
  • stage::Type{<:LifeStage}: LifeStage of organism being released.
  • gene_index::Int64: Genotype being released to implement the intervention.
  • times::Vector{Float64}: Release time points.
  • proportion::Float64: Proportion of modified organisms released during each time period.
  • callbacks::Vector
  • adults_to_count::String: The LifeStage of the standing population against which to measure proportional release size. Choices include Male, Female, or All.
source
GeneDrive.ProportionalReleaseMethod
ProportionalRelease(node::Node, organism, stage::Type{T}, gene_index, times::Vector, proportion::Float64, adult_counting::String) where T <: LifeStage

Return ProportionalRelease object specifying details of biological control intervention where release size is predicated on real-time wild population level.

source
GeneDrive.PupaType
struct Pupa <: LifeStage end

Pupa life stage. Applicable to holometabolous (complete) metamorphosing insect species.

source
GeneDrive.PupaDurationRossiType

Data for q_temperature_response. Species: AedesAegypti, LifeStage: Pupa. Source: Rossi et al (2014) and Poletti et al (2011) Table 1.

source
GeneDrive.ReleaseType
mutable struct Release <: Intervention
    node::Symbol
    organism::Type{<:Species}
    stage::Type{<:LifeStage}
    gene_index::Int64
    times::Vector{Float64}
    values::Vector{Float64}
    callbacks::Vector
end

Data for biological control interventions predicated on releasing modified organisms. Applicable to both suppression and replacement techniques.

Arguments

  • node::Symbol: Node where releases occur.
  • organism::Type{<:Species}: Species of organism being released.
  • stage::Type{<:LifeStage}: Life stage of organism being released.
  • gene_index::Int64: Genotype being released to implement the intervention.
  • times::Vector{Float64}: Release time points.
  • values::Union{Float64, Vector{Float64}}: Number of modified organisms released during each time period. Variably sized releases permitted.
  • callbacks::Vector
source
GeneDrive.ReleaseMethod
Release(node::Node, organism, stage, gene, times::Vector, fixed_release::Float64)

Return Release object specifying details of biological control intervention where release size is fixed over time.

source
GeneDrive.ReleaseMethod
Release(node::Node, organism, stage::Type{T}, gene, times::Vector, variable_release::Vector{Float64}) where T <: LifeStage

Return Release object specifying details of biological control intervention where release size is variable over time.

source
GeneDrive.ReleaseStrategyType
mutable struct ReleaseStrategy
    release_this_gene_index::Union{Nothing, Int64}=nothing
    release_this_life_stage=nothing
    release_location_force::Union{Nothing, Bool}=nothing
    release_start_time::Union{Nothing,Int64}=nothing
    release_end_time::Union{Nothing,Int64}=nothing
    release_time_interval::Int64=1
    release_size_min_per_timestep::Union{Int64, Float64}=0.0
    release_size_max_per_timestep::Union{Int64,Float64}=9e9
    release_max_over_timehorizon::Union{Int64,Float64}=9e9

end

Data defining the operational constraints for each node.

Arguments

  • release_this_gene_index: Gene index to be released. Generally defined by get_homozygous_modified.
  • release_this_life_stage: Lifestage to be released. Varies according to genetic technology.
  • release_location_force: Specify locations where releases are obligatory; only applicable when decision model is being run as an MINLP.
  • release_start_time: Timestep (day) that releases are permitted to start.
  • release_end_time: Timestep (day) that releases are required to end.
  • release_time_interval: Minimum timestep interval (in days) permitted for releases.
  • release_size_min_per_timestep: Minimum number of organisms that may be released on a single day.
  • release_size_max_per_timestep: Maximum number of organisms that may be released on a single day.
  • release_max_over_timehorizon: Maximum number of organisms that may be released over the problem horizon.
source
GeneDrive.ScenarioTemperatureType
mutable struct ScenarioTemperature <: Temperature end

Data for simulation that uses ensembles of temperature time series in °C.

Arguments

  • data::Matrix{Float64}: Array of ensemble members, each column of which is a time series of temperature values in °C.
  • probability::Vector{Float64}: Vector of probabilities with which given scenarios are expected to occur (one probability applicable per ensemble member).
  • selected_scenario::Int: Index to run dynamic model over selected scenarios; Int argument required. For decision model, use nothing.
source
GeneDrive.SinusoidalTemperatureType
mutable struct SinusoidalTemperature <: Temperature
    a::Float64
    b::Int64
    c::Int64
    d::Float64
end

Data for simulation that features sinusoidal temperature fluctuation in °C. Uses cosine implementation applicable to the Southern Hemisphere. Generated internally.

Arguments

  • a::Float64: Amplitude.
  • b::Float64: Periodicity coefficient.
  • c::Float64: Time period (days).
  • d::Float64: Mean.
source
GeneDrive.StageType
mutable struct Stage{L <: LifeStage}
    μ_temperature_response::TemperatureResponse
    q_temperature_response::TemperatureResponse
    n::Union{Nothing, Int64}
    density::Density{<: DensityDependence}
    dependency::Union{Nothing, Type{<:LifeStage}}
    N0::Int64
end

Data for life stages. Applies to any organism represented by stage-structured population equations.

Arguments

  • μ_temperature_response::TemperatureResponse: Mortality rate. Responsive to temperature.
  • q_temperature_response::TemperatureResponse: Developmental rate. 1/total time (days or portion thereof) spent in stage. Responsive to temperature.
  • n::Union{Nothing, Int64}: Number of bins allocated to stage (parameter, Erlang distribution).
  • density::Density: Specify density dependence model.
  • dependency::Union{Nothing, Type{<:LifeStage}}: Organism internal reference, do not modify.
  • N0::Int64: Initial stage-specific population count per node. Specify "0" for all stages except Female life stage.
source
GeneDrive.StageMethod
Stage{Female}(μ::Float64, n::Int64, density, dependency, N0::Int64)

Return Female life stage populated with input data. Not dynamically responsive to temperature. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.StageMethod
Stage{Female}(μ::TemperatureResponse, n::Int64, density, dependency, N0::Int64)

Return Female life stage populated with input data. Dynamically responsive to temperature. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.StageMethod
Stage{L}(μ::Float64, n::Int64, density, dependency, N0::Int64) where {L <: LifeStage}

Return juvenile life stage populated with input data. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.StageMethod
Stage{Male}(μ::Float64, n::Int64, density, dependency, N0::Int64)

Return Male life stage populated with input data. Not dynamically responsive to temperature. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.StageMethod
Stage{Male}(μ::TemperatureResponse, n::Int64, density, dependency, N0::Int64)

Return Male life stage populated with input data. Dynamically responsive to temperature. Applicable to holometabolous (complete) or hemimetabolous (partial) metamorphosing insect species.

source
GeneDrive.TemperatureSeriesDataType
mutable struct TemperatureSeriesData <: ExogenousChange
    node::Symbol
    times::Vector{Float64}
    values::Vector{Float64}
    set::diffeq.CallbackSet
end

Data for time series of temperature.

Arguments

  • node::Symbol: Node where temperature is realized.
  • times::Vector{Float64}: Time step (day) of realized temperature value.
  • values::Vector{Float64}: Temperature in °C.
  • set::diffeq.CallbackSet
source
GeneDrive.TemperatureShockDataType
mutable struct TemperatureShockData <: ExogenousChange
    node::Symbol
    times::Vector{Tuple{Float64, Float64}}
    values::Union{Float64, Vector{Float64}}
    set::diffeq.CallbackSet
end

Data for temperature shocks.

Arguments

  • node::Symbol: Node where temperature shock occurs.
  • times::Vector{Tuple{Float64,Float64}}: Vector of tuples defining shock start/stop time points. Multiple discrete shock periods may be included.
  • values::Union{Float64, Vector{Float64}}: Single fixed temperature value in °C or vector of variable temperature values in °C added to the baseline temperature during the specified time period. Multiple discrete temperature changes may be included; values may be positive or negative.
  • set::diffeq.CallbackSet
source
GeneDrive.TemperatureShockDataMethod
function TemperatureShockData(node::Node, times::Vector{Tuple{Float64,Float64}}, fixed_shock::Float64)

Return TemperatureShockData instance specifying data for fixed (all the same) shock values.

source
GeneDrive.TemperatureShockDataMethod
TemperatureShockData(node::Node, times::Vector{Tuple{Float64,Float64}}, variable_shock::Vector{Float64})

Return TemperatureShockData instance specifying data for variable shock values.

source
GeneDrive.TimeSeriesTemperatureType
mutable struct TimeSeriesTemperature <: Temperature end

Data for simulation that uses temperature time series in °C.

Arguments

  • value::Float64: Time series of temperature values in °C.
  • probability::Float64: Probability with which timeseries expected to be observed.
  • selected_scenario::Int
source
GeneDrive.assign_migration!Method
assign_migration!(network::Network, migration_data::Dict, species::Type{<:Species})

Return Network instance with migration field populated by user-specified transition rates.

Note to user:

  • Input data must be formatted as a nested dictionary. First level denotes relevant life stage and gene, second level includes to/from nodes and transition rate.
  • Stage and gene combinations not specified by input data retain a default transition rate of zero.
source
GeneDrive.compute_densityMethod
compute_density(data::Density{LinearDensity}, stage)

Return the effect of LinearDensity on life stage instance. Calculated each time step.

source
GeneDrive.compute_densityMethod
compute_density(data::Density{LogisticDensity}, stage)

Return the effect of LogisticDensity on life stage instance. Calculated each time step.

source
GeneDrive.compute_densityMethod
compute_density(data::Density{NoDensity}, stage)

Return the effect of NoDensity on life stage instance. Calculated each time step.

source
GeneDrive.count_genotypesMethod
count_genotypes(network::Network, node::Symbol, species::Type{<:Species})

Return the total count of genotypes from Genetics for Species in Node of Network.

source
GeneDrive.count_genotypesMethod
count_genotypes(node::Node, species::Type{<:Species})

Return the total count of genotypes from Genetics for Species in Node.

source
GeneDrive.count_substagesMethod
count_substages(network::Network, node::Symbol, species::Type{<:Species})

Return count of the total substages for Species in the specified Node of Network.

source
GeneDrive.count_substagesMethod
count_substages(node::Node, species::Type{<:Species}, stage::Type{<:LifeStage})

Return count of the substages in LifeStage for Species in Node.

source
GeneDrive.count_substagesMethod
count_substages(node::Node, species::Type{<:Species}, stage::Type{Female})

Return count of the substages for Species in Node. Specific to the Female lifestage.

source
GeneDrive.count_temperature_scenariosMethod
count_temperature_scenarios(temperature_model::ScenarioTemperature)

Return the total count of temperature scenarios in ScenarioTemperature or TimeSeriesTemperature object.

source
GeneDrive.create_decision_modelMethod
create_decision_model(network::Network, tspan; node_strategy::Dict, species::Type{<:Species}=AedesAegypti,do_binary::Bool=false, optimizer=nothing)

Build mathematical program. Problem created as an NLP (dobinary=false) or MINLP (dobinary=true).

source
GeneDrive.create_decision_modelMethod
create_decision_model(node::Node, tspan; node_strategy::Dict, species::Type{<:Species}=AedesAegypti,do_binary::Bool=false, optimizer=nothing)

Build mathematical program. Problem created as an NLP (dobinary=false) or MINLP (dobinary=true). NB: Node is recreated as a Network object internally; this does not change the problem but is relevant for data exploration as it adds one index layer to the formatted results.

source
GeneDrive.get_densityMethod
get_density(node::Node, species::Type{<:Species}, life_stage::Type{<:LifeStage})

Return the density dependence model and parameterization for LifeStage of Species in Node.

source
GeneDrive.get_durationMethod
get_duration(node::Node, species::Type{<:Species}, life_stage::Type{<:LifeStage})

Return the temperature-sensitive duration (q_temperature_response) for the LifeStage of Species in Node.

source
GeneDrive.get_exogenous_interventionMethod
    function get_exogenous_intervention(inputs::ExogenousInputs, node::Node, organism, stage, gene)

Returns biological control intervention relevant to specified `Node`, `Organism`, `LifeStage`, and `Genotype`.
source
GeneDrive.get_exogenous_interventionMethod
get_exogenous_intervention(inputs::ExogenousInputs, node::Node, organism, stage::Type{Female}, gene, female_index)

Return biological control intervention relevant to Female.

source
GeneDrive.get_geneticsMethod
get_genetics(network::Network, node::Symbol, species::Type{<:Species})

Return Genetics data for Species in Node of Network.

source
GeneDrive.get_genotypesMethod
get_genotypes(node::Node, species::Type{<:Species})

Return genotype data from Genetics for Species in Node.

source
GeneDrive.get_homozygous_modifiedMethod
get_homozygous_modified(network::Network, node::Node, species::Type{<:Species})

Return the index of the homozygous modified genotype in Genetics for Species in Node of Network.

source
GeneDrive.get_homozygous_modifiedMethod
get_homozygous_modified(node::Node, species::Type{<:Species})

Return the index of the homozygous modified genotype in Genetics for Species in Node.

source
GeneDrive.get_lifestageMethod
get_lifestage(node::Node, species::Type{<:Species}, life_stage::Type{<:LifeStage})

Return data for specified LifeStage of Species in Node.

source
GeneDrive.get_migrationMethod
get_migration(network::Network, species::Type{<:Species})

Return the migration characterizing each genotype and Lifestage for Species in Network.

source
GeneDrive.get_mortalityMethod
get_mortality(node::Node, species::Type{<:Species}, life_stage::Type{<:LifeStage})

Return the temperature-sensitive mortality (μ_temperature_response) for the LifeStage of Species in Node.

source
GeneDrive.get_namesMethod
get_names(network::Network)

Return the names (symbols) of the nodes within the Network object.

source
GeneDrive.get_organismsMethod
get_organisms(network::Network, node::Symbol)

Return vector of species populating the specified Node of Network.

source
GeneDrive.get_previous_lifestageMethod
get_previous_lifestage(node::Node, species::Type{<:Species}, ::Stage{T}) where T <: LifeStage

Show LifeStage dependency: return data for the LifeStage previous to the specified LifeStage of Species in Node.

source
GeneDrive.get_probabilityMethod
get_probability(temperature_model::ScenarioTemperature)

Return the probability with which each temperature scenario occurs from ScenarioTemperature or TimeSeriesTemperature object.

source
GeneDrive.get_release_dataMethod
get_release_data(optimized_schedule::Vector{Float64})

Return re-formatted schedule of optimal release times and values produced by decision model, for simplified application in dynamic model.

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::AdultMortalityAbiodun, ::Float64)

Return q_temperature_response. Species: AnophelesGambiae, LifeStage: Male, Female. Source: Abiodun et al (2016).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::AdultMortalityMoustaid, ::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage: Male, Female. Source: El Moustaid et al (2019).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::AdultMortalityRossi, ::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage: Male, Female. Source: Rossi et al (2014).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::EggDurationAbiodun, ::Float64)

Return q_temperature_response. Species: AnophelesGambiae, LifeStage: Egg. Source: Abiodun et al (2016).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::EggDurationMoustaid, ::Float64)

Return q_temperature_response. Species: AedesAegypti, LifeStage: Egg. Source: El Moustaid et al (2019).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::EggDurationRossi, ::Float64)

Return q_temperature_response. Species: AedesAegypti, LifeStage: Egg. Source: Rossi et al (2014).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::EggMortalityAbiodun, ::Float64)

Return q_temperature_response. Species: AnophelesGambiae, LifeStage: Egg. Source: Abiodun et al (2016).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::EggMortalityMoustaid, duration::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage: Egg. Source: El Moustaid et al (2019).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::EggMortalityRossi, ::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage: Egg. Source: Rossi et al (2014).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::LarvaDurationAbiodun, ::Float64)

Return q_temperature_response. Species: AnophelesGambiae, LifeStage: Larva. Source: Abiodun et al (2016).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::LarvaDurationMoustaid, ::Float64)

Return q_temperature_response. Species: AedesAegypti, LifeStage: Larva. Source: El Moustaid et al (2019).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::LarvaDurationRossi, ::Float64)

Return q_temperature_response. Species: AedesAegypti, LifeStage: Larva. Source: Rossi et al (2014).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::LarvaMortalityAbiodun, ::Float64)

Return q_temperature_response. Species: AnophelesGambiae, LifeStage: Larva. Source: Abiodun et al (2016).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::LarvaMortalityMoustaid, duration::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage: Larva. Source: El Moustaid et al (2019).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::LarvaMortalityRossi, ::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage:, larval stage. Source: Rossi et al (2014).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(::Float64, response::NoResponse, ::Float64) = response.baseline_value

Return NoResponse data (model without temperature response).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::PupaDurationAbiodun, ::Float64)

Return q_temperature_response. Species: AnophelesGambiae, LifeStage: Pupa. Source: Abiodun et al (2016).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::PupaDurationMoustaid, ::Float64)

Return q_temperature_response. Species: AedesAegypti, LifeStage: Pupa. Source: El Moustaid et al (2019).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::PupaDurationRossi, ::Float64)

Return q_temperature_response. Species: AedesAegypti, LifeStage: Pupa. Source: Rossi et al (2014) and Poletti et al (2011) Table 1.

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::PupaMortalityAbiodun, ::Float64)

Return q_temperature_response. Species: AnophelesGambiae, LifeStage: Pupa. Source: Abiodun et al (2016).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::PupaMortalityMoustaid, duration::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage: Pupa. Source: El Moustaid et al (2019).

source
GeneDrive.get_temperature_responseMethod
get_temperature_response(ctemp::Float64, response::PupaMortalityRossi, ::Float64)

Return μ_temperature_response. Species: AedesAegypti, LifeStage: Pupa. Source: Rossi et al (2014).

source
GeneDrive.get_temperature_scenariosMethod
get_temperature_scenarios(temperature_model::ScenarioTemperature)

Return the values of all temperature scenarios in ScenarioTemperature or TimeSeriesTemperature object.

source
GeneDrive.get_temperature_valueMethod
get_temperature_value(temperature_model::ConstantTemperature, temp_value_from_inputs::Float64, t)

Return value of temperature in °C for specifed time step. Accounts for perturbations to baseline temperature model where applicable.

source
GeneDrive.get_temperature_valueMethod
get_temperature_value(temperature_model::ScenarioTemperature, temp_value_from_inputs::Float64, t)

Return value of temperature in °C for selected scenario and specifed time step. Perturbations to baseline temperature, where applicable, should be directly added to time series prior to running model.

source
GeneDrive.get_temperature_valueMethod
get_temperature_value(temperature_model::SinusoidalTemperature, temp_value_from_inputs::Float64, t)

Return value of temperature in °C for specifed time step. Accounts for perturbations to baseline temperature model where applicable.

source
GeneDrive.get_temperature_valueMethod
get_temperature_value(temperature_model::TimeSeriesTemperature, temp_value_from_inputs::Float64, t)

Return value of temperature in °C for specifed time step. Perturbations to baseline temperature, where applicable, should be directly added to time series prior to running model.

source
GeneDrive.get_wildtypeMethod
get_wildtype(network::Network, node::Node, species::Type{<:Species})

Return the index of the wildtype genotype in Genetics for Species in Node of Network.

source
GeneDrive.get_wildtypeMethod
get_wildtype(node::Node, species::Type{<:Species})

Return the index of the wildtype genotype in Genetics for Species in Node.

source
GeneDrive.make_organismsMethod
make_organisms(species::Type{<:Species},genetics::Genetics,stages::DataStructures.OrderedDict)

Return Organism object.

source
GeneDrive.perturb_temperature_timeseries!Method
perturb_temperature_timeseries(current_temperature::TimeSeriesTemperature, perturbation)

Alter daily values across entire temperature timeseries by the size of perturbation input. Perturbation may be positive or negative.

source
GeneDrive.set_scenario!Method
set_scenario!(data::ScenarioTemperature, selected_scenario::Int)

Update the selected_scenario of Temperature for ScenarioTemperature in Node.

source
GeneDrive.set_scenario!Method
set_scenario!(data::ScenarioTemperature, selected_scenario::Int)

Update the selected_scenario of Temperature for ScenarioTemperature.

source
GeneDrive.solve_decision_modelFunction
solve_decision_model(model::JuMP.Model, objective_function::Nothing=nothing; kwargs...)

Solve mathematical program using default objective function @objective(model, Min, 0). This permits comparison to dynamic population model output.

source
GeneDrive.solve_dynamic_modelMethod
solve_dynamic_model(network::Network, releases::Union{Vector{Release}, Vector{ProportionalRelease}}, algorithm, tspan)

Return ODE model solution for network problem with releases.

source
GeneDrive.solve_dynamic_modelMethod
solve_dynamic_model(network::Network, releases::Union{Vector{Release},Vector{ProportionalRelease}}, shocks::Vector{TemperatureShockData}, algorithm, tspan)

Return ODE model solution for network problem with releases and temperature shocks.

source
GeneDrive.solve_dynamic_modelMethod
solve_dynamic_model(network::Network, shocks::Vector{TemperatureShockData}, algorithm, tspan)

Return ODE model solution for network problem with temperature shocks.

source
GeneDrive.solve_dynamic_modelMethod
solve_dynamic_model(node::Node, shocks::TemperatureShockData, algorithm, tspan)

Return ODE model solution for single node problem with temperature shocks.

source
GeneDrive.solve_dynamic_modelMethod
solve_dynamic_model(node::Node, releases::Union{Vector{Release},Vector{ProportionalRelease}}, algorithm, tspan)

Return ODE model solution for single node problem with releases.

source
GeneDrive.solve_dynamic_modelMethod
solve_dynamic_model(node::Node, releases::Union{Vector{Release},Vector{ProportionalRelease}}, shocks::TemperatureShockData, algorithm, tspan)

Return ODE model solution for single node problem with releases and temperature shocks.

source
GeneDrive.solve_scenarios_dynamic_modelMethod
solve_scenarios_dynamic_model(network::Network, algorithm, tspan, scenarios_of_interest::Vector{Int})

Return ODE model solution for network problem with scenarions and no releases.

source
GeneDrive.solve_scenarios_dynamic_modelMethod
solve_scenarios_dynamic_model(network::Network, releases, algorithm, tspan, scenarios_of_interest::Vector{Int})

Return ODE model solution for node problem with scenarios and releases.

source
GeneDrive.solve_scenarios_dynamic_modelMethod
solve_scenarios_dynamic_model(node::Node, algorithm, tspan, scenarios_of_interest::Vector{Int})

Return ODE model solution for node problem with scenarios and no releases.

source
GeneDrive.solve_scenarios_dynamic_modelMethod
solve_scenarios_dynamic_model(node::Node, releases, algorithm, tspan, scenarios_of_interest::Vector{Int})

Return ODE model solution for node problem with scenarios and releases.

source
GeneDrive.temperature_effectMethod
temperature_effect(ctemp::Float64, stage)

Return effect of temperature on organism vital rates (mortality, duration = μ_temperature_response, q_temperature_response).

source
GeneDrive.update_density!Method
update_density!(node::Node, species::Type{<:Species}, ::Type{T}; new_density::Density) where T <: LifeStage

Update both the functional form and parameterization of the density dependence model for LifeStage of Species in Node.

source
GeneDrive.update_density!Method
update_density!(stages::DataStructures.OrderedDict,lifestage::Type{T},new_density::Density) where {T <: LifeStage}

Update both the functional form and parameterization of the density dependence model for LifeStage.

source
GeneDrive.update_density_modelMethod
update_density_model(node::Node, species::Type{<:Species}, ::Type{T}; new_density_model::Density) where T <: LifeStage

Update the functional form of the density dependence model for LifeStage of Species in Node.

source
GeneDrive.update_density_parameterMethod
update_density_parameter(node::Node, species::Type{<:Species}, ::Type{T}; new_param_value::Float64) where T <: LifeStage

Update the parameterization of the density dependence model for LifeStage of Species in Node.

source
GeneDrive.update_durationMethod
update_duration(node::Node, species::Type{<:Species}, life_stage::Type{<:LifeStage}, new_q)

Update the temperature-sensitive duration (q_temperature_response) for the LifeStage of Species in Node.

source
GeneDrive.update_egg_durationMethod
update_egg_duration(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated duration for Egg stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_egg_mortalityMethod
update_egg_mortality(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated mortality for Egg stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_female_mortalityMethod
update_female_mortality(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated mortality for Female stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_genetics_SMethod
update_genetics_S(node::Node, species::Type{<:Species}, new_sigma::Array{Float64,1})

Update S data in Genetics for Species in Node.

source
GeneDrive.update_genetics_ΒMethod
update_genetics_Β(node::Node, species::Type{<:Species}, new_beta::Array{Float64,1})

Update Β data in Genetics for Species in Node.

source
GeneDrive.update_genetics_ΗMethod
update_genetics_Η(node::Node, species::Type{<:Species}, new_eta::Array{Float64,1})

Update Η data in Genetics for Species in Node.

source
GeneDrive.update_genetics_ΩMethod
update_genetics_Ω(node::Node, species::Type{<:Species}, new_omega::Array{Float64,1})

Update Ω data in Genetics for Species in Node.

source
GeneDrive.update_larva_durationMethod
update_larva_duration(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated duration for Larva stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_larva_mortalityMethod
update_larva_mortality(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated mortality for Larva stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_male_mortalityMethod
update_male_mortality(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated mortality for Male stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_migrationMethod
update_migration(network::Network, species::Type{<:Species}, new_migration)

Update the migration characterizing each genotype and Lifestage for Species in Network.

source
GeneDrive.update_mortalityMethod
update_mortality(node::Node, species::Type{<:Species}, life_stage::Type{<:LifeStage}, new_μ)

Update the temperature-sensitive mortality (μ_temperature_response) for the LifeStage of Species in Node.

source
GeneDrive.update_organismMethod
update_organism(network::Network, node::Symbol, new_species)

Update the species populating the specified Node of Network.

source
GeneDrive.update_population_sizeMethod
update_population_size(stages::DataStructures, Stage}, new_popsize::Int64)

Return updated population size. Note: new_popsize argument refers specifically to the Female population; if e.g. new_popsize = 500, the full adult population (Females and Males) will be 500*2.

source
GeneDrive.update_pupa_durationMethod
update_pupa_duration(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated duration for Pupa stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_pupa_mortalityMethod
update_pupa_mortality(stages::DataStructures.OrderedDict, vital_rate::Float64)

Return updated mortality for Pupa stage. Note: Exclusively applicable to NoResponse temperature type.

source
GeneDrive.update_temperature!Method
update_temperature!(node::Node, temp_type::Type{<:ConstantTemperature}, new_temperature::Float64)

Update the values of ConstantTemperature for Node.

source
GeneDrive.update_temperature!Method
update_temperature!(node::Node, temp_type::Type{<:SinusoidalTemperature}, new_temperature::Vector{Float64})

Update the values of Temperature for SinusoidalTemperature in Node.

source
GeneDrive.update_temperature!Method
update_temperature!(node::Node, temp_type::Type{<:TimeSeriesTemperature}, new_temperature::Vector{Float64})

Update the values of Temperature for TimeSeriesTemperature in Node.

source

API Reference