SingleWell< InfectionState > Class Template Reference

CPP API: SingleWell< InfectionState > Class Template Reference

Implementation of diffusive ABM, see dabm::Model. More...

#include <single_well.h>

Classes

struct  Agent
 

Public Types

using Position = Eigen::Vector2d
 
using Status = InfectionState
 

Public Member Functions

ScalarType adoption_rate (const Agent &agent, const Status &new_status) const
 Calculate adoption rate for an Agent. More...
 
std::map< std::tuple< mio::regions::Region, Status, Status >, mio::AdoptionRate< ScalarType, Status > > & get_adoption_rates ()
 
mio::RandomNumberGenerator & get_rng ()
 Get the RandomNumberGenerator used by this Model for random events. More...
 
void move (const ScalarType, const ScalarType dt, Agent &agent)
 Perform one integration step of the diffusion process for a given Agent using the Euler-Maruyama method. More...
 
 SingleWell (const std::vector< Agent > &agents, const std::vector< mio::AdoptionRate< ScalarType, Status >> &rates, ScalarType contact_radius=0.4, ScalarType sigma=0.4, std::vector< Status > non_moving_states={})
 Set up a diffusive ABM using the quadwell potential F(x,y) = (x^4 + y^4)/2. More...
 
Eigen::VectorX< ScalarTypetime_point () const
 Calculate the current system state i.e. More...
 

Static Public Member Functions

static constexpr void adopt (Agent &agent, const Status &new_status)
 Perform infection state adoption for an Agent. More...
 

Public Attributes

std::vector< Agentpopulations
 Vector containing all Agents in the model. More...
 

Private Member Functions

bool is_contact (const Agent &agent, const Agent &contact) const
 Evaluate whether two agents are within their contact radius. More...
 
bool is_in_domain (const Position &p, const ScalarType lower_domain_border=-2, const ScalarType upper_domain_border=2) const
 Restrict domain to [-2, 2]^2 where "escaping" is impossible. More...
 

Static Private Member Functions

static Position grad_U (const Position p)
 Calculate the gradient of the potential at a given position. More...
 

Private Attributes

std::map< std::tuple< mio::regions::Region, Status, Status >, mio::AdoptionRate< ScalarType, Status > > m_adoption_rates
 Map of AdoptionRates according to their region index and their from -> to infection states. More...
 
ScalarType m_contact_radius
 Agents' interaction radius. Within this radius agents are considered as contacts. More...
 
std::vector< Statusm_non_moving_states
 Infection states within which agents do not change their location. More...
 
mio::RandomNumberGenerator m_rng
 Model's random number generator. More...
 
ScalarType m_sigma
 Noise term of the diffusion process. More...
 

Detailed Description

template<class InfectionState>
class SingleWell< InfectionState >

Implementation of diffusive ABM, see dabm::Model.

This implementation defines a diffusion process for the potential F(x,y) = (x^4 + y^4)/2.

Template Parameters
InfectionStateAn infection state enum.

Member Typedef Documentation

◆ Position

template<class InfectionState >
using SingleWell< InfectionState >::Position = Eigen::Vector2d

◆ Status

template<class InfectionState >
using SingleWell< InfectionState >::Status = InfectionState

Constructor & Destructor Documentation

◆ SingleWell()

template<class InfectionState >
SingleWell< InfectionState >::SingleWell ( const std::vector< Agent > &  agents,
const std::vector< mio::AdoptionRate< ScalarType, Status >> &  rates,
ScalarType  contact_radius = 0.4,
ScalarType  sigma = 0.4,
std::vector< Status non_moving_states = {} 
)
inline

Set up a diffusive ABM using the quadwell potential F(x,y) = (x^4 + y^4)/2.

Parameters
[in]agentsA vector of Agents representing the population.
[in]ratesAdoptionRates defining InfectionState adoptions, see mio::AdoptionRate.
[in]contact_radiusContact radius for second-order adoptions.
[in]sigmaNoise term for the diffusion process.
[in]non_moving_stateInfectionStates that are excluded from movement e.g. Dead.

Member Function Documentation

◆ adopt()

template<class InfectionState >
static constexpr void SingleWell< InfectionState >::adopt ( Agent agent,
const Status new_status 
)
inlinestaticconstexpr

Perform infection state adoption for an Agent.

Parameters
[in,out]agentAgent whose infection state is changed.
[in]new_statusAgent's new infection state.

◆ adoption_rate()

template<class InfectionState >
ScalarType SingleWell< InfectionState >::adoption_rate ( const Agent agent,
const Status new_status 
) const
inline

Calculate adoption rate for an Agent.

Parameters
[in]agentAgent for whom the adoption rate is calculated.
[in]new_statusTarget infection state of the adoption rate, see mio::AdoptionRate.to.
Returns
Value of agent-dependent AdoptionRate.

◆ get_adoption_rates()

template<class InfectionState >
std::map<std::tuple<mio::regions::Region, Status, Status>, mio::AdoptionRate<ScalarType, Status> >& SingleWell< InfectionState >::get_adoption_rates ( )
inline

◆ get_rng()

template<class InfectionState >
mio::RandomNumberGenerator& SingleWell< InfectionState >::get_rng ( )
inline

Get the RandomNumberGenerator used by this Model for random events.

Returns
The random number generator.

◆ grad_U()

template<class InfectionState >
static Position SingleWell< InfectionState >::grad_U ( const Position  p)
inlinestaticprivate

Calculate the gradient of the potential at a given position.

Parameters
[in]pPosition at which the gradient is evaluated.
Returns
Value of potential gradient.

◆ is_contact()

template<class InfectionState >
bool SingleWell< InfectionState >::is_contact ( const Agent agent,
const Agent contact 
) const
inlineprivate

Evaluate whether two agents are within their contact radius.

Parameters
[in]agentAgent whose position specifies the contact area.
[in]contactPotential contact of agent.
Returns
Boolean specifying whether are within each others contact radius.

◆ is_in_domain()

template<class InfectionState >
bool SingleWell< InfectionState >::is_in_domain ( const Position p,
const ScalarType  lower_domain_border = -2,
const ScalarType  upper_domain_border = 2 
) const
inlineprivate

Restrict domain to [-2, 2]^2 where "escaping" is impossible.

Parameters
[in]pPosition to check.
Returns
Boolean specifying whether p is in [-2, 2]^2.

◆ move()

template<class InfectionState >
void SingleWell< InfectionState >::move ( const  ScalarType,
const ScalarType  dt,
Agent agent 
)
inline

Perform one integration step of the diffusion process for a given Agent using the Euler-Maruyama method.

Parameters
[in]dtStep size.
[in]agentAgent to be moved.

◆ time_point()

template<class InfectionState >
Eigen::VectorX<ScalarType> SingleWell< InfectionState >::time_point ( ) const
inline

Calculate the current system state i.e.

the populations for each infection state.

Returns
Vector containing the number of agents per infection state.

Member Data Documentation

◆ m_adoption_rates

template<class InfectionState >
std::map<std::tuple<mio::regions::Region, Status, Status>, mio::AdoptionRate<ScalarType, Status> > SingleWell< InfectionState >::m_adoption_rates
private

Map of AdoptionRates according to their region index and their from -> to infection states.

◆ m_contact_radius

template<class InfectionState >
ScalarType SingleWell< InfectionState >::m_contact_radius
private

Agents' interaction radius. Within this radius agents are considered as contacts.

◆ m_non_moving_states

template<class InfectionState >
std::vector<Status> SingleWell< InfectionState >::m_non_moving_states
private

Infection states within which agents do not change their location.

◆ m_rng

template<class InfectionState >
mio::RandomNumberGenerator SingleWell< InfectionState >::m_rng
private

Model's random number generator.

◆ m_sigma

template<class InfectionState >
ScalarType SingleWell< InfectionState >::m_sigma
private

Noise term of the diffusion process.

◆ populations

template<class InfectionState >
std::vector<Agent> SingleWell< InfectionState >::populations

Vector containing all Agents in the model.