Nemo  2.4.0b
Simulate forward-in-time genetic evolution in a spatially explicit, individual-based stochastic simulator
LCE_Selection_base Class Reference

Base class performing (viability) selection on an arbitrary trait. More...

#include <LCEselection.h>

+ Inheritance diagram for LCE_Selection_base:
+ Collaboration diagram for LCE_Selection_base:

Public Member Functions

 LCE_Selection_base ()
 
virtual ~LCE_Selection_base ()
 
bool setSelectionMatrix ()
 
bool setSelectTraitMapping (unsigned int num_quanti_traits)
 
bool setSelectionOffset (double default_val, double min_val)
 
bool setLocalOptima ()
 
const TMatrixgetLocalOptima () const
 
bool set_sel_model ()
 
bool set_fit_model ()
 
bool set_local_optima ()
 
bool set_param_rate_of_change ()
 
void set_std_rate_of_change ()
 
void addPhenotypicSD (unsigned int deme, double *stDev)
 
void changeLocalOptima ()
 
void resetCounters ()
 Resets the fitness counters. More...
 
void setMeans (unsigned int tot_ind)
 Computes the average fitness of each pedigree class. More...
 
void updateFitnessCounters (double fitness, unsigned int ped_class, bool survived)
 Updates the fitness and survival mean counters. More...
 
double getMeanFitness (age_idx age)
 Computes the mean fitness of the whole population for a given age class. More...
 
double getMeanPatchFitness (age_idx age, unsigned int p)
 Computes the mean fitness in a given patch for a given age class. More...
 
double getMaxFitness (age_idx age)
 Computes the maximum fitness value of the whole population for a given age class. More...
 
double getMaxPatchFitness (age_idx age, unsigned int p)
 Computes the maximum fitness value in a given patch for a given age class. More...
 
double setMeanFitness (age_idx age)
 Sets the _mean_fitness variable to the value of the mean population fitness. More...
 
double getFitness (Individual *ind, unsigned int patch)
 Calls the fitness function according to the fitness model. More...
 
double getFitnessFixedEffect (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual in the fixed selection model. More...
 
double getFitnessDirect (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following the direct selection model. More...
 
double getFitnessUnivariateQuadratic (Individual *ind, unsigned int patch, unsigned int trait)
 Quadratic fitness surface, approximates the Gaussian model for weak selection and/or small deviation from the optimum. More...
 
double getFitnessUnivariateLinear (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following a Linear selection model with a single trait under selection. More...
 
double getFitnessUnivariateDisruptive (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following a Disruptive (inverted Gaussian) selection model with one trait under selection. More...
 
double getFitnessMultivariateDisruptive (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following a Disruptive (inverted Gaussian) selection model with several traits under selection. More...
 
double getFitnessTruncation (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following a Truncation selection model. More...
 
double getFitnessMultivariateGaussian (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following the Gaussian selection model with one trait under selection. More...
 
double getFitnessUnivariateGaussian (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following the Gaussian selection model with several traits under selection. More...
 
double getFitnessMultivariateGaussian_VE (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following the Gaussian selection model with one trait under selection and environmental variance. More...
 
double getFitnessUnivariateGaussian_VE (Individual *offsprg, unsigned int patch, unsigned int trait)
 Returns the fitness of an individual following the Gaussian selection model with several traits under selection and environmental variance. More...
 
double getFitnessAbsolute (Individual *ind, unsigned int patch)
 Returns the raw fitness of the individual, without adjustment (absolute fitness). More...
 
double getFitnessRelative (Individual *ind, unsigned int patch)
 Returns the relative fitness of the individual, adjusted by a scaling factor. More...
 
void setScalingFactorLocal (age_idx age, unsigned int p)
 Sets the fitness scaling factor equal to the inverse of the mean local patch fitness. More...
 
void setScalingFactorGlobal (age_idx age, unsigned int p)
 Sets the fitness scaling factor equal to the inverse of the mean population fitness. More...
 
void setScalingFactorMaxLocal (age_idx age, unsigned int p)
 Sets the fitness scaling factor equal to the inverse of the maximum local patch fitness value. More...
 
void setScalingFactorMaxGlobal (age_idx age, unsigned int p)
 Sets the fitness scaling factor equal to the inverse of the maximum population fitness value. More...
 
void setScalingFactorForLinearSelection (age_idx age, unsigned int p)
 Sets the fitness scaling factor as the mean trait value in the patch to compute the fitness value in the linear selection model. More...
 
void setScalingFactorAbsolute (age_idx age, unsigned int p)
 Resets the fitness scaling factor equal to one. More...
 
void doViabilitySelection (sex_t SEX, age_idx AGE, Patch *patch, unsigned int p)
 Selectively removes individuals in the population depending on their fitness. More...
 
void checkChangeLocalOptima ()
 Check is rate of change in local optima has been set and must be applied. More...
 
Implementations
virtual bool setParameters ()
 
virtual void execute ()
 
virtual void loadStatServices (StatServices *loader)
 
virtual void loadFileServices (FileServices *loader)
 
virtual bool resetParameterFromSource (std::string param, SimComponent *cmpt)
 
virtual LifeCycleEventclone ()
 
virtual age_t removeAgeClass ()
 
virtual age_t addAgeClass ()
 
virtual age_t requiredAgeClass ()
 
- Public Member Functions inherited from LifeCycleEvent
 LifeCycleEvent (const char *name, const char *trait_link)
 Cstor. More...
 
virtual ~LifeCycleEvent ()
 
virtual void init (Metapop *popPtr)
 Sets the pointer to the current Metapop and the trait link if applicable. More...
 
virtual bool attach_trait (string trait)
 
virtual void set_paramset (std::string name, bool required, SimComponent *owner)
 
virtual void set_event_name (std::string &name)
 Set the name of the event (name of the ParamSet) and add the corresponding parameter to the set. More...
 
virtual void set_event_name (const char *name)
 
virtual string & get_event_name ()
 Accessor to the LCE's name. More...
 
virtual int get_rank ()
 Accessor to the LCE rank in the life cycle. More...
 
virtual void set_pop_ptr (Metapop *popPtr)
 Accessors for the population pointer. More...
 
virtual Metapopget_pop_ptr ()
 
- Public Member Functions inherited from SimComponent
 SimComponent ()
 
virtual ~SimComponent ()
 
virtual void loadUpdaters (UpdaterServices *loader)
 Loads the parameters and component updater onto the updater manager. More...
 
virtual void set_paramset (ParamSet *paramset)
 Sets the ParamSet member. More...
 
virtual void set_paramsetFromCopy (const ParamSet &PSet)
 Reset the set of parameters from a another set. More...
 
virtual ParamSetget_paramset ()
 ParamSet accessor. More...
 
virtual void add_parameter (Param *param)
 Interface to add a parameter to the set. More...
 
virtual void add_parameter (std::string Name, param_t Type, bool isRequired, bool isBounded, double low_bnd, double up_bnd)
 Interface to add a parameter to the set. More...
 
virtual void add_parameter (std::string Name, param_t Type, bool isRequired, bool isBounded, double low_bnd, double up_bnd, ParamUpdaterBase *updater)
 Interface to add a parameter and its updater to the set. More...
 
virtual Paramget_parameter (std::string name)
 Param getter. More...
 
virtual double get_parameter_value (std::string name)
 Param value getter. More...
 
virtual string get_name ()
 Returnd the name of the ParamSet, i.e. More...
 
virtual bool has_parameter (std::string name)
 Param getter. More...
 

Protected Attributes

double _base_fitness
 
double _mean_fitness
 
double _max_fitness
 
double _scaling_factor
 
bool _is_local
 
bool _is_absolute
 
double _fitness [5]
 Fitness counters, one for each pedigree class. More...
 
double _survival [5]
 
double _ind_cntr [5]
 
vector< string > _Traits
 The list of trait types under selection. More...
 
vector< unsigned int > _TraitIndices
 The indices of the traits under selection. More...
 
vector< string > _SelectionModels
 The selection models associated with each trait under selection. More...
 
LCE_SelectionSH_stater
 
LCE_SelectionFH_writer
 
Quantitative traits parameters
int _selectTraitDimension
 Number of quantitative traits under selection. More...
 
vector< unsigned int > _selectTraitMapping
 Mapping from selection dimension index to quanti trait index (0-based). More...
 
vector< double > _selection_variance
 Patch-specific selection variance. More...
 
double _eVariance
 Evironmental variance. More...
 
Pointers to functions, set in LCE_Selection_base::setParameters
vector< double(LCE_Selection_base::*)(Individual *, unsigned int, unsigned int) > _getRawFitness
 A vector containing pointers to fitness function related to each trait under selection. More...
 
double(LCE_Selection_base::* _getFitness )(Individual *, unsigned int)
 Pointer to the function returning the individual fitness. More...
 
void(LCE_Selection_base::* _setScalingFactor )(age_idx, unsigned int)
 Pointer to the function used to set the fitness scaling factor when fitness is relative. More...
 
void(LCE_Selection_base::* _setNewLocalOptima )(void)
 Pointer to the function used to change the local phenotypic optima or its rate of change. More...
 
- Protected Attributes inherited from LifeCycleEvent
std::string _event_name
 The param name to be read in the init file. More...
 
Metapop_popPtr
 The ptr to the current Metapop. More...
 
std::string _LCELinkedTraitType
 The name of the linked trait. More...
 
int _LCELinkedTraitIndex
 The index in the individual's trait table of the linked trait. More...
 
- Protected Attributes inherited from SimComponent
ParamSet_paramSet
 The parameters container. More...
 

Private Attributes

Gaussian selection variables
vector< TMatrix * > _selection_matrix
 
vector< gsl_matrix * > _gsl_selection_matrix
 
gsl_vector * _diffs
 
gsl_vector * _res1
 
TMatrix _local_optima
 
TMatrix _rate_of_change_local_optima
 
bool _do_change_local_opt
 
bool _rate_of_change_is_std
 
unsigned int _set_std_rate_at_generation
 
int _std_rate_reference_patch
 
Linear selection parameters

Array to temporarily hold the phenotypic values of an individual, is never allocated.

TMatrix _linear_selection_strength
 
vector< double > _mean_trait_value
 
unsigned int _linearSelTraitIndex
 
age_idx _linearSelAgeClass
 
Shared selection offset (used by linear, disruptive, and truncation models)
TMatrix _selection_offset
 
Truncation selection parameters
TMatrix _truncation_threshold
 
bool _truncation_upper
 
Fixed selection (lethal-equivalents-based) parameters
double _letheq
 
double _Fpedigree [5]
 Array of pedigree values used in the fixed-selection model. More...
 
double _FitnessFixModel [5]
 Absolute fitness values of the five pedigree class for the fixed selection model (lethal equivalents model). More...
 

Friends

class LCE_SelectionSH
 The StatHandler associated class is a friend. More...
 
class LCE_SelectionFH
 
class LCE_Breed_Selection
 

Detailed Description

Base class performing (viability) selection on an arbitrary trait.

Implements several types of selection models, from fixed inbreeding depression-type model to Gaussian stabilizing selection on an adaptive landscape. This class is the base-class for the composite LCEs breed_selection and breed_selection_disperse.

Constructor & Destructor Documentation

◆ LCE_Selection_base()

LCE_Selection_base::LCE_Selection_base ( )
46  : LifeCycleEvent("viability_selection", ""),
49 _diffs(0),
50 _res1(0),
58 _truncation_upper(true),
59 _letheq(0),
60 _base_fitness(1),
61 _mean_fitness(0),
62 _max_fitness(0),
64 _is_local(0),
65 _is_absolute(1),
66 _eVariance(0),
68 _getFitness(0),
71 _stater(0),
72 _writer(0)
73 {
74  // mandatory parameter (no default):
75  add_parameter("selection_trait", STR, true, false, 0, 0); //no updaters here, to keep it safe...
76 
79 
80  add_parameter("selection_fitness_model",STR,false,false,0,0, updater);
81 
83 
84  add_parameter("selection_model",STR,true,false,0,0, updater); //mandatory, depends on trait
85  add_parameter("selection_matrix",MAT,false,false,0,0, updater);
86  add_parameter("selection_variance",DBL,false,false,0,0, updater);
87  add_parameter("selection_correlation",DBL,false,false,0,0, updater);
88  add_parameter("selection_trait_dimension",INT,false,false,0,0, updater);
89  add_parameter("selection_base_fitness",DBL,false,true,0,1, updater);
90  add_parameter("selection_lethal_equivalents",DBL,false,false,0,0, updater);
91  add_parameter("selection_pedigree_F",MAT,false,false,0,0, updater);
92  add_parameter("selection_randomize",BOOL,false,false,0,0, updater);
93  add_parameter("selection_environmental_variance", DBL, false, false, 0, 0, updater);
94 
95  // linear selection:
96  add_parameter("selection_linear_strength", DBL, false, false, 0, 0, updater);
97  add_parameter("selection_offset", DBL, false, false, 0, 0, updater);
98  add_parameter("selection_at_stage", STR, false, false, 0, 0, updater);
99 
100  // truncation selection:
101  add_parameter("selection_truncation_threshold", DBL, false, false, 0, 0, updater);
102  add_parameter("selection_truncation_direction", STR, false, false, 0, 0, updater);
103 
105  add_parameter("selection_local_optima",DBL,false,false,0,0, updater);
106 
108  add_parameter("selection_rate_environmental_change", DBL, false, false, 0, 0, updater);
109  add_parameter("selection_std_rate_environmental_change", DBL, false, false, 0, 0, updater);
110  add_parameter("selection_std_rate_set_at_generation", INT, false, false, 0, 0, updater);
111  add_parameter("selection_std_rate_reference_patch", INT, false, false, 0, 0, updater);
112 
113 
114  add_parameter("selection_output", BOOL, false, false, 0, 0, 0);
115  add_parameter("selection_output_logtime", INT, false, false, 0, 0, 0);
116  add_parameter("selection_output_dir", STR, false, false, 0, 0, 0);
117 }
int _selectTraitDimension
Number of quantitative traits under selection.
Definition: LCEselection.h:124
double _base_fitness
Definition: LCEselection.h:112
double(LCE_Selection_base::* _getFitness)(Individual *, unsigned int)
Pointer to the function returning the individual fitness.
Definition: LCEselection.h:152
double _max_fitness
Definition: LCEselection.h:113
vector< gsl_matrix * > _gsl_selection_matrix
Definition: LCEselection.h:56
LCE_SelectionFH * _writer
Definition: LCEselection.h:166
unsigned int _linearSelTraitIndex
Definition: LCEselection.h:75
gsl_vector * _diffs
Definition: LCEselection.h:57
bool set_param_rate_of_change()
Definition: LCEselection.cc:888
void(LCE_Selection_base::* _setScalingFactor)(age_idx, unsigned int)
Pointer to the function used to set the fitness scaling factor when fitness is relative.
Definition: LCEselection.h:157
bool _is_absolute
Definition: LCEselection.h:115
double _eVariance
Evironmental variance.
Definition: LCEselection.h:133
double _scaling_factor
Definition: LCEselection.h:113
bool _do_change_local_opt
Definition: LCEselection.h:61
bool _is_local
Definition: LCEselection.h:114
gsl_vector * _res1
Definition: LCEselection.h:58
bool set_sel_model()
Definition: LCEselection.cc:312
unsigned int _set_std_rate_at_generation
Definition: LCEselection.h:63
double _mean_fitness
Definition: LCEselection.h:113
bool set_fit_model()
Definition: LCEselection.cc:254
LCE_SelectionSH * _stater
Definition: LCEselection.h:165
int _std_rate_reference_patch
Definition: LCEselection.h:64
void(LCE_Selection_base::* _setNewLocalOptima)(void)
Pointer to the function used to change the local phenotypic optima or its rate of change.
Definition: LCEselection.h:161
vector< double(LCE_Selection_base::*)(Individual *, unsigned int, unsigned int) > _getRawFitness
A vector containing pointers to fitness function related to each trait under selection.
Definition: LCEselection.h:148
bool _truncation_upper
Definition: LCEselection.h:87
bool set_local_optima()
Definition: LCEselection.cc:860
double _letheq
Definition: LCEselection.h:93
vector< TMatrix * > _selection_matrix
Definition: LCEselection.h:55
bool _rate_of_change_is_std
Definition: LCEselection.h:62
age_idx _linearSelAgeClass
Definition: LCEselection.h:76
LifeCycleEvent(const char *name, const char *trait_link)
Cstor.
Definition: lifecycleevent.h:98
Implementation of the ParamUpdaterBase interface.
Definition: param.h:363
virtual void add_parameter(Param *param)
Interface to add a parameter to the set.
Definition: simcomponent.h:112
@ DBL
Definition: types.h:78
@ MAT
Definition: types.h:78
@ BOOL
Definition: types.h:78
@ STR
Definition: types.h:78
@ INT
Definition: types.h:78
@ OFFSx
Definition: types.h:42

References SimComponent::add_parameter(), BOOL, DBL, INT, MAT, set_fit_model(), set_local_optima(), set_param_rate_of_change(), set_sel_model(), and STR.

Referenced by clone().

◆ ~LCE_Selection_base()

LCE_Selection_base::~LCE_Selection_base ( )
virtual
122 {
123  vector< TMatrix* >::iterator selIT = _selection_matrix.begin();
124  for(; selIT != _selection_matrix.end(); ++selIT)
125  if((*selIT)) delete (*selIT);
126  _selection_matrix.clear();
127 
128  for(unsigned int i = 0; i < _gsl_selection_matrix.size(); ++i)
129  if(_gsl_selection_matrix[i]) gsl_matrix_free(_gsl_selection_matrix[i]);
130  _gsl_selection_matrix.clear();
131 
132  if(_diffs) gsl_vector_free(_diffs);
133  if(_res1) gsl_vector_free(_res1);
134  if(_stater) delete _stater;
135 }

References _diffs, _gsl_selection_matrix, _res1, _selection_matrix, and _stater.

Member Function Documentation

◆ addAgeClass()

virtual age_t LCE_Selection_base::addAgeClass ( )
inlinevirtual

Implements LifeCycleEvent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

349 {return NONE;}
#define NONE
No age flag.
Definition: types.h:48

References NONE.

◆ addPhenotypicSD()

void LCE_Selection_base::addPhenotypicSD ( unsigned int  deme,
double *  stDev 
)
1513 {
1514  Patch *patch = _popPtr->getPatch(deme);
1515  Individual* ind;
1516  unsigned int table_size = patch->size(OFFSx);
1517  double* _phe;
1518 
1519  double **phenot = new double* [_selectTraitDimension];
1520 
1521  for (int i = 0; i < _selectTraitDimension; ++i) {
1522 
1523  phenot[i] = new double [table_size];
1524 
1525  }
1526 
1527  unsigned int pos = 0;
1528 
1529  for (unsigned int j = 0; j < patch->size(FEM, OFFSx) && pos < table_size; ++j) {
1530 
1531  ind = patch->get(FEM, OFFSx, j);
1532  _phe = (double*)ind->getTraitValue(_LCELinkedTraitIndex);
1533 
1534  for (int i = 0; i < _selectTraitDimension; ++i)
1535  phenot[i][pos] = _phe[_selectTraitMapping[i]];
1536 
1537  pos++;
1538 
1539  }
1540 
1541  for (unsigned int j = 0; j < patch->size(MAL, OFFSx) && pos < table_size; ++j) {
1542 
1543  ind = patch->get(MAL, OFFSx, j);
1544  _phe = (double*)ind->getTraitValue(_LCELinkedTraitIndex);
1545 
1546  for (int i = 0; i < _selectTraitDimension; ++i)
1547  phenot[i][pos] = _phe[_selectTraitMapping[i]];
1548 
1549  pos++;
1550 
1551  }
1552 
1553  assert(pos == table_size);
1554 
1555  for (int i = 0; i < _selectTraitDimension; ++i) {
1556  stDev[i] += sqrt( my_variance_with_fixed_mean( phenot[i], table_size, my_mean(phenot[i], table_size) ) );
1557  delete [] phenot[i];
1558  }
1559 
1560  delete [] phenot;
1561 }
This class contains traits along with other individual information (sex, pedigree,...
Definition: individual.h:49
void * getTraitValue(IDX T)
Accessor to the value (phenotype) of a particular trait.
Definition: individual.h:271
vector< unsigned int > _selectTraitMapping
Mapping from selection dimension index to quanti trait index (0-based).
Definition: LCEselection.h:127
int _LCELinkedTraitIndex
The index in the individual's trait table of the linked trait.
Definition: lifecycleevent.h:89
Metapop * _popPtr
The ptr to the current Metapop.
Definition: lifecycleevent.h:81
Patch * getPatch(unsigned int i)
Patch accessor, return the ith+1 patch in the metapop.
Definition: metapop.h:257
Second class in the metapopulation design structure, between the Metapop and Individual classes.
Definition: metapop.h:432
unsigned int size(age_t AGE)
Returns the size of the container of the appropriate age class(es) for both sexes.
Definition: metapop.h:498
Individual * get(sex_t SEX, age_idx AGE, unsigned int at)
Returns a pointer to the individual sitting at the index passed.
Definition: metapop.h:534
@ FEM
Definition: types.h:37
@ MAL
Definition: types.h:37
double my_mean(double *data, unsigned int size)
Definition: utils.cc:38
double my_variance_with_fixed_mean(double *data, unsigned int size, double mean)
Definition: utils.cc:63

References LifeCycleEvent::_LCELinkedTraitIndex, LifeCycleEvent::_popPtr, _selectTraitDimension, _selectTraitMapping, FEM, Patch::get(), Metapop::getPatch(), Individual::getTraitValue(), MAL, my_mean(), my_variance_with_fixed_mean(), OFFSx, and Patch::size().

Referenced by set_std_rate_of_change().

◆ changeLocalOptima()

void LCE_Selection_base::changeLocalOptima ( )
1296 {
1297  for (int i = 0; i < _selectTraitDimension; ++i)
1298  for (unsigned int j = 0; j < _local_optima.nrows(); ++j) {
1300  }
1301 }
TMatrix _rate_of_change_local_optima
Definition: LCEselection.h:60
TMatrix _local_optima
Definition: LCEselection.h:59
double get(unsigned int i, unsigned int j) const
Accessor to element at row i and column j.
Definition: tmatrix.h:193
void plus(unsigned int i, unsigned int j, double value)
Adds a value to an element of the matrix.
Definition: tmatrix.h:256
unsigned int nrows() const
Definition: tmatrix.h:213

References _local_optima, _rate_of_change_local_optima, _selectTraitDimension, TMatrix::get(), TMatrix::nrows(), and TMatrix::plus().

Referenced by set_param_rate_of_change(), and set_std_rate_of_change().

◆ checkChangeLocalOptima()

void LCE_Selection_base::checkChangeLocalOptima ( )

Check is rate of change in local optima has been set and must be applied.

1306 {
1307  // check if change in local optima for quantitative traits
1308  if(_do_change_local_opt) {
1309 
1310  if(_popPtr->getCurrentGeneration() == 1) {
1311 
1312  set_local_optima(); //reset local optima to initial values
1313 
1314  if(_rate_of_change_is_std) //reset rate of change relative to SD of that replicate
1316  }
1317 
1318  (this->*_setNewLocalOptima)();
1319  }
1320 }
void set_std_rate_of_change()
Definition: LCEselection.cc:1428
unsigned int getCurrentGeneration()
Definition: metapop.h:296

References _do_change_local_opt, LifeCycleEvent::_popPtr, _rate_of_change_is_std, _setNewLocalOptima, Metapop::getCurrentGeneration(), set_local_optima(), and set_std_rate_of_change().

Referenced by LCE_Breed_Selection::execute(), LCE_Breed_Selection_Disperse::execute(), and execute().

◆ clone()

virtual LifeCycleEvent* LCE_Selection_base::clone ( )
inlinevirtual

Implements LifeCycleEvent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

347 {return new LCE_Selection_base();}
LCE_Selection_base()
Definition: LCEselection.cc:46

References LCE_Selection_base().

◆ doViabilitySelection()

void LCE_Selection_base::doViabilitySelection ( sex_t  SEX,
age_idx  AGE,
Patch patch,
unsigned int  p 
)

Selectively removes individuals in the population depending on their fitness.

Calls the fitness function for each individual. Updates the fitness counters. The fitness scaling factor is set outside this function, in the execute() procedure.

Parameters
SEXthe sex class
AGEthe age-class index
patchthe focal patch
pthe focal patch index
1359 {
1360  Individual* ind;
1361  double fitness;
1362  bool not_survived;
1363 
1364  for(unsigned int i = 0; i < patch->size(SEX, AGE); i++) {
1365 
1366  ind = patch->get(SEX, AGE, i);
1367 
1368  fitness = getFitness( ind, p);
1369 
1370  not_survived = ( RAND::Uniform() > fitness );
1371 
1372  if( not_survived ) {
1373 
1374  //this individual dies
1375  patch->remove(SEX, AGE, i);
1376 
1377  _popPtr->recycle(ind);
1378 
1379  i--;
1380 
1381  } //else; this individual stays in the patch
1382 
1383  updateFitnessCounters(fitness, ind->getPedigreeClass(), !not_survived);
1384  }
1385 }
void recycle(Individual *ind)
Put an individual in the recycling pool.
Definition: indfactory.h:62
unsigned int getPedigreeClass()
Returns the pedigree class of the individual, as set during offspring creation.
Definition: individual.h:179
double getFitness(Individual *ind, unsigned int patch)
Calls the fitness function according to the fitness model.
Definition: LCEselection.h:230
void updateFitnessCounters(double fitness, unsigned int ped_class, bool survived)
Updates the fitness and survival mean counters.
Definition: LCEselection.cc:1416
Individual * remove(sex_t SEX, age_idx AGE, unsigned int at)
Removes the individual sitting at the given index in the appropriate container.
Definition: metapop.h:588
static double Uniform()
Generates a random number from [0.0, 1.0[ uniformly distributed.
Definition: Uniform.h:127

References LifeCycleEvent::_popPtr, Patch::get(), getFitness(), Individual::getPedigreeClass(), IndFactory::recycle(), Patch::remove(), Patch::size(), RAND::Uniform(), and updateFitnessCounters().

Referenced by execute().

◆ execute()

void LCE_Selection_base::execute ( )
virtual

Implements LifeCycleEvent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

1325 {
1326  Patch * patch;
1327  unsigned int popSize = _popPtr->size(OFFSPRG);
1328 
1329 #ifdef _DEBUG_
1330  message("LCE_Selection_base::execute (tot offsprg nb: %i", popSize);
1331 #endif
1332 
1333  resetCounters();
1334 
1335  // check if change in local optima for quantitative traits
1337 
1338  for(unsigned int p = 0; p < _popPtr->getPatchNbr(); p++) {
1339 
1340  (this->*_setScalingFactor)(OFFSx, p);
1341 
1342  patch = _popPtr->getPatch(p);
1343 
1344  doViabilitySelection(FEM, OFFSx, patch, p);
1345 
1346  doViabilitySelection(MAL, OFFSx, patch, p);
1347  }
1348 
1349  setMeans(popSize);
1350 
1351 #ifdef _DEBUG_
1352  message(", after selection: %i, mean fitness offspring=%f)\n",_popPtr->size(OFFSPRG), _mean_fitness);
1353 #endif
1354 }
void doViabilitySelection(sex_t SEX, age_idx AGE, Patch *patch, unsigned int p)
Selectively removes individuals in the population depending on their fitness.
Definition: LCEselection.cc:1358
void checkChangeLocalOptima()
Check is rate of change in local optima has been set and must be applied.
Definition: LCEselection.cc:1305
void setMeans(unsigned int tot_ind)
Computes the average fitness of each pedigree class.
Definition: LCEselection.cc:1400
void resetCounters()
Resets the fitness counters.
Definition: LCEselection.cc:1389
unsigned int size()
Get the total number of individuals present in the population, all sex and age classes together.
Definition: metapop.h:312
unsigned int getPatchNbr()
Definition: metapop.h:276
void message(const char *message,...)
Definition: output.cc:40
#define OFFSPRG
Offspring age class flag.
Definition: types.h:50

References _mean_fitness, LifeCycleEvent::_popPtr, _setScalingFactor, checkChangeLocalOptima(), doViabilitySelection(), FEM, Metapop::getPatch(), Metapop::getPatchNbr(), MAL, message(), OFFSPRG, OFFSx, resetCounters(), setMeans(), and Metapop::size().

◆ getFitness()

double LCE_Selection_base::getFitness ( Individual ind,
unsigned int  patch 
)
inline

Calls the fitness function according to the fitness model.

The fitness model can be "absolute", "relative_local" or "relative_global".

Parameters
indthe focal indvidual, we want to know its fitness
patchthe index of the patch of the focal individual
231  {
232  return (this->*_getFitness)(ind, patch);
233  }

References _getFitness.

Referenced by LCE_Breed_Selection::do_breed_selection_FecFitness(), doViabilitySelection(), getMaxPatchFitness(), getMeanFitness(), getMeanPatchFitness(), LCE_Breed_Selection::makeOffspringWithSelection(), LCE_SelectionSH::setDataTable(), LCE_Breed_Selection::setReproScaledFitness_gsl(), and LCE_Breed_Selection::setReproScaledFitness_sum().

◆ getFitnessAbsolute()

double LCE_Selection_base::getFitnessAbsolute ( Individual ind,
unsigned int  patch 
)

Returns the raw fitness of the individual, without adjustment (absolute fitness).

Calls the fitness function according to the right selection model.

Parameters
indthe focal indvidual, we want to know its fitness
patchthe index of the patch of the focal individual
1135 {
1136  double fitness = 1;
1137  //at this point, _getRawFitness.size() == _TraitIndices.size()
1138  //and we assume that the functions are "aligned" with the traits
1139  for(unsigned int i = 0; i < _getRawFitness.size(); i++)
1140  fitness *= (this->*_getRawFitness[i])(ind, patch, _TraitIndices[i]);
1141  return fitness;
1142 }
vector< unsigned int > _TraitIndices
The indices of the traits under selection.
Definition: LCEselection.h:140

References _getRawFitness, and _TraitIndices.

Referenced by getFitnessRelative(), and set_fit_model().

◆ getFitnessDirect()

double LCE_Selection_base::getFitnessDirect ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)
inline

Returns the fitness of an individual following the direct selection model.

243  {
244  return *(double*)offsprg->getTraitValue(trait);
245  }

References Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessFixedEffect()

double LCE_Selection_base::getFitnessFixedEffect ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)
inline

Returns the fitness of an individual in the fixed selection model.

237  {
238  return _FitnessFixModel[ offsprg->getPedigreeClass() ];
239  }
double _FitnessFixModel[5]
Absolute fitness values of the five pedigree class for the fixed selection model (lethal equivalents ...
Definition: LCEselection.h:100

References _FitnessFixModel, and Individual::getPedigreeClass().

Referenced by set_sel_model().

◆ getFitnessMultivariateDisruptive()

double LCE_Selection_base::getFitnessMultivariateDisruptive ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following a Disruptive (inverted Gaussian) selection model with several traits under selection.

1003 {
1004  double res2;
1005 
1006  double *_phe = (double*)ind->getTraitValue(trait);
1007 
1008  // W = offset * (1 - exp(-0.5 * dT * Sigma^-1 * d)) inverted multivariate Gaussian
1009  for( int i = 0; i < _selectTraitDimension; i++)
1010  gsl_vector_set(_diffs, i, _phe[_selectTraitMapping[i]] - _local_optima.get(patch, i));
1011 
1012  // (diff)T * W * diff:
1013  gsl_blas_dsymv(CblasUpper, 1.0, _gsl_selection_matrix[patch], _diffs, 0.0, _res1);
1014  gsl_blas_ddot(_diffs, _res1, &res2);
1015 
1016  return _selection_offset.get(patch, 0) * (1.0 - exp( -0.5 * res2 ));
1017 }
TMatrix _selection_offset
Definition: LCEselection.h:81

References _diffs, _gsl_selection_matrix, _local_optima, _res1, _selection_offset, _selectTraitDimension, _selectTraitMapping, TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessMultivariateGaussian()

double LCE_Selection_base::getFitnessMultivariateGaussian ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following the Gaussian selection model with one trait under selection.

1078 {
1079  double res2;
1080 
1081  double *_phe = (double*)ind->getTraitValue(trait);
1082 
1083  for( int i = 0; i < _selectTraitDimension; i++)
1084  gsl_vector_set(_diffs, i, _phe[_selectTraitMapping[i]] - _local_optima.get(patch, i));
1085 
1086  //(diff)T * W * diff:
1087  //right partial product:
1088  gsl_blas_dsymv(CblasUpper, 1.0, _gsl_selection_matrix[patch], _diffs, 0.0, _res1);
1089  //left product:
1090  gsl_blas_ddot(_diffs, _res1, &res2);
1091 
1092  return exp( -0.5 * res2 );
1093 }

References _diffs, _gsl_selection_matrix, _local_optima, _res1, _selectTraitDimension, _selectTraitMapping, TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessMultivariateGaussian_VE()

double LCE_Selection_base::getFitnessMultivariateGaussian_VE ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following the Gaussian selection model with one trait under selection and environmental variance.

1114 {
1115  double res2;
1116 
1117  double *_phe = (double*)ind->getTraitValue(trait);
1118 
1119  //add the environmental variance here:
1120  for( int i = 0; i < _selectTraitDimension; i++)
1121  gsl_vector_set(_diffs, i, _phe[_selectTraitMapping[i]] + RAND::Gaussian(_eVariance) - _local_optima.get(patch, i));
1122 
1123  //(diff)T * W * diff:
1124  //right partial product:
1125  gsl_blas_dsymv(CblasUpper, 1.0, _gsl_selection_matrix[patch], _diffs, 0.0, _res1);
1126  //left product:
1127  gsl_blas_ddot(_diffs, _res1, &res2);
1128 
1129  return exp( -0.5 * res2 );
1130 }
static double Gaussian(double sigma)
Definition: Uniform.h:264

References _diffs, _eVariance, _gsl_selection_matrix, _local_optima, _res1, _selectTraitDimension, _selectTraitMapping, RAND::Gaussian(), TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessRelative()

double LCE_Selection_base::getFitnessRelative ( Individual ind,
unsigned int  patch 
)
inline

Returns the relative fitness of the individual, adjusted by a scaling factor.

The scaling factor is set according to the relative fitness model (local or global) outside of this function. If the scaling factor is one, this is the absolute fitness of the individual. Calls the fitness function according to the right selection model.

Parameters
indthe focal indvidual, we want to know its fitness
patchthe index of the patch of the focal individual
289  {
290  return getFitnessAbsolute(ind, patch) * _scaling_factor;
291  }
double getFitnessAbsolute(Individual *ind, unsigned int patch)
Returns the raw fitness of the individual, without adjustment (absolute fitness).
Definition: LCEselection.cc:1134

References _scaling_factor, and getFitnessAbsolute().

Referenced by set_fit_model().

◆ getFitnessTruncation()

double LCE_Selection_base::getFitnessTruncation ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following a Truncation selection model.

Fitness = 1 + offset if phenotype passes the threshold test, offset otherwise.

1022 {
1023  double *_phe = (double*)ind->getTraitValue(trait);
1024 
1025  double w = 1.0;
1026 
1027  // W = prod_t [ (1 + offset) if z passes threshold, offset otherwise ]
1028  for(int t = 0; t < _selectTraitDimension; ++t) {
1029 
1030  bool passes = _truncation_upper
1031  ? (_phe[_selectTraitMapping[t]] >= _truncation_threshold.get(patch, t))
1032  : (_phe[_selectTraitMapping[t]] <= _truncation_threshold.get(patch, t));
1033 
1034  w *= passes ? (1.0 + _selection_offset.get(patch, t)) : _selection_offset.get(patch, t);
1035  }
1036 
1037  return w;
1038 }
TMatrix _truncation_threshold
Definition: LCEselection.h:86

References _selection_offset, _selectTraitDimension, _selectTraitMapping, _truncation_threshold, _truncation_upper, TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessUnivariateDisruptive()

double LCE_Selection_base::getFitnessUnivariateDisruptive ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following a Disruptive (inverted Gaussian) selection model with one trait under selection.

991 {
992  double *_phe = (double*)ind->getTraitValue(trait);
993 
994  // W = offset * (1 - exp(-d^2 / (2*Vs))) inverted Gaussian, single trait
995  double diff = _phe[_selectTraitMapping[0]] - _local_optima.get(patch, 0);
996 
997  return _selection_offset.get(patch, 0) * (1.0 - exp( -0.5 * diff * diff / _selection_variance[patch] ));
998 }
vector< double > _selection_variance
Patch-specific selection variance.
Definition: LCEselection.h:130

References _local_optima, _selection_offset, _selection_variance, _selectTraitMapping, TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessUnivariateGaussian()

double LCE_Selection_base::getFitnessUnivariateGaussian ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following the Gaussian selection model with several traits under selection.

1063 {
1064  double res2, diff;
1065 
1066  double *_phe = (double*)ind->getTraitValue(trait);
1067 
1068  diff = _phe[_selectTraitMapping[0]] - _local_optima.get(patch, 0);
1069 
1070  res2 = diff*diff / _selection_variance[patch];
1071 
1072  return exp( -0.5 * res2 );
1073 }

References _local_optima, _selection_variance, _selectTraitMapping, TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessUnivariateGaussian_VE()

double LCE_Selection_base::getFitnessUnivariateGaussian_VE ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following the Gaussian selection model with several traits under selection and environmental variance.

1098 {
1099  double res2, diff;
1100 
1101  double *_phe = (double*)ind->getTraitValue(trait);
1102 
1103  //add the environmental variance here:
1104  diff = _phe[_selectTraitMapping[0]] + RAND::Gaussian(_eVariance) - _local_optima.get(patch, 0);
1105 
1106  res2 = diff*diff / _selection_variance[patch];
1107 
1108  return exp( -0.5 * res2 );
1109 }

References _eVariance, _local_optima, _selection_variance, _selectTraitMapping, RAND::Gaussian(), TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessUnivariateLinear()

double LCE_Selection_base::getFitnessUnivariateLinear ( Individual offsprg,
unsigned int  patch,
unsigned int  trait 
)

Returns the fitness of an individual following a Linear selection model with a single trait under selection.

974 {
975  double *_phe = (double*)ind->getTraitValue(trait);
976 
977  double w = 1.0;
978 
979  // W = prod_t [ offset(patch,t) + strength(patch,t) * (z_t - z_bar_t) ]
980  // With defaults: W = prod_t [ 1 + beta * (z_t - z_bar_t) ] (Carter et al. 2005)
981  for(int t = 0; t < _selectTraitDimension; ++t)
982  w *= _selection_offset.get(patch, t)
984 
985  return (w < 0 ? 0 : w);
986 }
vector< double > _mean_trait_value
Definition: LCEselection.h:74
TMatrix _linear_selection_strength
Definition: LCEselection.h:73

References _linear_selection_strength, _mean_trait_value, _selection_offset, _selectTraitDimension, _selectTraitMapping, TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getFitnessUnivariateQuadratic()

double LCE_Selection_base::getFitnessUnivariateQuadratic ( Individual ind,
unsigned int  patch,
unsigned int  trait 
)

Quadratic fitness surface, approximates the Gaussian model for weak selection and/or small deviation from the optimum.

1043 {
1044  double res2, diff, w = 1.0;
1045 
1046  double *_phe = (double*)ind->getTraitValue(trait);
1047 
1048  for(unsigned int t = 0; t < _selectTraitDimension; ++t){
1049 
1050  diff = _phe[_selectTraitMapping[t]] - _local_optima.get(patch, t);
1051 
1052  res2 = 1 - diff*diff / _selection_variance[patch];
1053 
1054  w *= (res2 < 0 ? 0 : res2);
1055  }
1056 
1057  return w;
1058 }

References _local_optima, _selection_variance, _selectTraitDimension, _selectTraitMapping, TMatrix::get(), and Individual::getTraitValue().

Referenced by set_sel_model().

◆ getLocalOptima()

const TMatrix& LCE_Selection_base::getLocalOptima ( ) const
inline
178 {return _local_optima;}

References _local_optima.

◆ getMaxFitness()

double LCE_Selection_base::getMaxFitness ( age_idx  age)

Computes the maximum fitness value of the whole population for a given age class.

Parameters
agethe age-class index
1181 {
1182  double max = 0, local_max;
1183 
1184  for(unsigned int i = 0, npatch = _popPtr->getPatchNbr(); i < npatch; i++) {
1185 
1186  local_max = getMaxPatchFitness(age, i);
1187 
1188  max = (local_max > max ? local_max : max);
1189  }
1190 
1191  return max;
1192 }
double getMaxPatchFitness(age_idx age, unsigned int p)
Computes the maximum fitness value in a given patch for a given age class.
Definition: LCEselection.cc:1196

References LifeCycleEvent::_popPtr, getMaxPatchFitness(), and Metapop::getPatchNbr().

Referenced by setScalingFactorMaxGlobal().

◆ getMaxPatchFitness()

double LCE_Selection_base::getMaxPatchFitness ( age_idx  age,
unsigned int  p 
)

Computes the maximum fitness value in a given patch for a given age class.

Parameters
agethe age-class index
pthe patch index
1197 {
1198  double max = 0, fit;
1199  Patch *patch = _popPtr->getPatch(p);
1200 
1201  for(unsigned int j = 0, size = patch->size(FEM, age); j < size; j++) {
1202  fit = getFitness( patch->get(FEM, age, j), p);
1203  max = (fit > max ? fit : max);
1204  }
1205  for(unsigned int j = 0, size = patch->size(MAL, age); j < size; j++){
1206  fit = getFitness( patch->get(MAL, age, j), p);
1207  max = (fit > max ? fit : max);
1208  }
1209 
1210  return max;
1211 }

References LifeCycleEvent::_popPtr, FEM, Patch::get(), getFitness(), Metapop::getPatch(), MAL, and Patch::size().

Referenced by getMaxFitness(), and setScalingFactorMaxLocal().

◆ getMeanFitness()

double LCE_Selection_base::getMeanFitness ( age_idx  age)

Computes the mean fitness of the whole population for a given age class.

Parameters
agethe age-class index
1147 {
1148  double mean = 0;
1149  Patch *patch;
1150 // age_idx age = (AGE == ADULTS ? ADLTx : OFFSx);
1151 
1152  for(unsigned int i = 0, npatch = _popPtr->getPatchNbr(); i < npatch; i++) {
1153  patch = _popPtr->getPatch(i);
1154  for(unsigned int j = 0, size = patch->size(FEM, age); j < size; j++)
1155  mean += getFitness( patch->get(FEM, age, j), i);
1156  for(unsigned int j = 0, size = patch->size(MAL, age); j < size; j++)
1157  mean += getFitness( patch->get(MAL, age, j), i);
1158  }
1159  return mean/_popPtr->size(age);
1160 }

References LifeCycleEvent::_popPtr, FEM, Patch::get(), getFitness(), Metapop::getPatch(), Metapop::getPatchNbr(), MAL, Metapop::size(), and Patch::size().

Referenced by LCE_Breed_Selection_Disperse::breed_selection_disperse(), LCE_Breed_Selection_Disperse::execute(), setMeanFitness(), and setScalingFactorGlobal().

◆ getMeanPatchFitness()

double LCE_Selection_base::getMeanPatchFitness ( age_idx  age,
unsigned int  p 
)

Computes the mean fitness in a given patch for a given age class.

Parameters
agethe age-class index
pthe patch index
1165 {
1166  double mean = 0;
1167  Patch *patch = _popPtr->getPatch(p);
1168 
1169  for(unsigned int j = 0, size = patch->size(FEM, age); j < size; j++)
1170  mean += getFitness( patch->get(FEM, age, j), p);
1171 
1172  for(unsigned int j = 0, size = patch->size(MAL, age); j < size; j++)
1173  mean += getFitness( patch->get(MAL, age, j), p);
1174 
1175  return mean/patch->size(age);
1176 }

References LifeCycleEvent::_popPtr, FEM, Patch::get(), getFitness(), Metapop::getPatch(), MAL, and Patch::size().

Referenced by LCE_Breed_Selection_Disperse::breed_selection_disperse(), and setScalingFactorLocal().

◆ loadFileServices()

void LCE_Selection_base::loadFileServices ( FileServices loader)
virtual

Implements SimComponent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

148 {
149  if(_paramSet->isSet("selection_output")) {
150 
151  if(_writer == NULL) _writer = new LCE_SelectionFH(this);
152 
153  Param* param = get_parameter("selection_output_logtime");
154 
155  if(!param->isSet()) fatal("parameter \"selection_output_logtime\" is missing\n");
156 
157  if(param->isMatrix()) {
158 
159  TMatrix temp;
160 
161  param->getMatrix(&temp);
162 
163  _writer->set_multi(true, true, 1, &temp, get_parameter("selection_output_dir")->getArg());
164 
165  } else
166  // rpl_per, gen_per, rpl_occ, gen_occ, rank, path, self-ref
167  _writer->set(true, param->isSet(), 1, (param->isSet() ? (int)param->getValue() : 0), 0,
168  get_parameter("selection_output_dir")->getArg(), this);
169 
170 
171  loader->attach(_writer);
172 
173  } else {
174 
175  if(_writer) delete _writer;
176 
177  return;
178  }
179 
180 
181 }
virtual void set(bool rpl_per, bool gen_per, int rpl_occ, int gen_occ, int rank, string path, LCE *event)
Definition: filehandler.h:276
virtual void set_multi(bool rpl_per, bool gen_per, int rpl_occ, TMatrix *Occ, string path)
Definition: filehandler.h:201
virtual void attach(Handler *FH)
Attaches the FileHandler to the current list (_writers) of the FileServices.
Definition: fileservices.cc:61
friend class LCE_SelectionFH
Definition: LCEselection.h:107
bool isSet()
Accessor to the status flag.
Definition: param.h:288
This structure stores one parameter, its definition and its string argument.
Definition: param.h:54
double getValue()
Returns the argument value according to its type.
Definition: param.cc:368
bool isMatrix()
Checks if the argument is of matrix type.
Definition: param.h:172
void getMatrix(TMatrix *mat)
Sets the matrix from the argument string if the parameter is set and of matrix type.
Definition: param.cc:378
bool isSet()
Definition: param.h:140
virtual Param * get_parameter(std::string name)
Param getter.
Definition: simcomponent.h:139
ParamSet * _paramSet
The parameters container.
Definition: simcomponent.h:48
A class to handle matrix in params, coerces matrix into a vector of same total size.
Definition: tmatrix.h:50
void fatal(const char *str,...)
Definition: output.cc:100

References SimComponent::_paramSet, _writer, FileServices::attach(), fatal(), SimComponent::get_parameter(), Param::getMatrix(), Param::getValue(), Param::isMatrix(), Param::isSet(), ParamSet::isSet(), LCE_SelectionFH, EventFileHandler< LCE >::set(), and FileHandler::set_multi().

◆ loadStatServices()

void LCE_Selection_base::loadStatServices ( StatServices loader)
virtual

Implements SimComponent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

140 {
141  if(_stater == NULL) _stater = new LCE_SelectionSH(this);
142  loader->attach(_stater);
143 }
friend class LCE_SelectionSH
The StatHandler associated class is a friend.
Definition: LCEselection.h:106
virtual void attach(Handler *H)
attach the StatHandler to the current list (_statHandlers) of the StatServices
Definition: statservices.cc:177

References _stater, StatServices::attach(), and LCE_SelectionSH.

Referenced by LCE_Breed_Selection::loadStatServices(), and LCE_Breed_Selection_Disperse::loadStatServices().

◆ removeAgeClass()

virtual age_t LCE_Selection_base::removeAgeClass ( )
inlinevirtual

Implements LifeCycleEvent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

348 {return NONE;}

References NONE.

◆ requiredAgeClass()

virtual age_t LCE_Selection_base::requiredAgeClass ( )
inlinevirtual

Implements LifeCycleEvent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

350 {return OFFSPRG;}

References OFFSPRG.

◆ resetCounters()

void LCE_Selection_base::resetCounters ( )

Resets the fitness counters.

1390 {
1391  for(unsigned int i = 0; i < 5; i++) {
1392  _fitness[i] = 0;
1393  _survival[i] = 0;
1394  _ind_cntr[i] = 0;
1395  }
1396 }
double _fitness[5]
Fitness counters, one for each pedigree class.
Definition: LCEselection.h:118
double _ind_cntr[5]
Definition: LCEselection.h:118
double _survival[5]
Definition: LCEselection.h:118

References _fitness, _ind_cntr, and _survival.

Referenced by LCE_Breed_Selection::execute(), LCE_Breed_Selection_Disperse::execute(), execute(), and setParameters().

◆ resetParameterFromSource()

virtual bool LCE_Selection_base::resetParameterFromSource ( std::string  param,
SimComponent cmpt 
)
inlinevirtual

Implements SimComponent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

346 {return false;}

◆ set_fit_model()

bool LCE_Selection_base::set_fit_model ( )
255 {
256 
257  if(_paramSet->isSet("selection_fitness_model")) {
258 
259  string fit_model = _paramSet->getArg("selection_fitness_model");
260 
261  if( fit_model == "absolute") {
262 
263  _is_local = false;
264  _is_absolute = true;
267 
268  } else if( fit_model == "relative_local" ) {
269 
270  _is_local = true;
271  _is_absolute = false;
274 
275  } else if(fit_model == "relative_global") {
276 
277  _is_local = false;
278  _is_absolute = false;
281 
282  } else if(fit_model == "relative_max_local") {
283 
284  _is_local = true;
285  _is_absolute = false;
288 
289  } else if(fit_model == "relative_max_global") {
290 
291  _is_local = false;
292  _is_absolute = false;
295 
296  } else {
297  return error("Unknown fitness model \"%s\"", fit_model.c_str());
298  }
299  } //default case:
300  else {
301  _is_local = false;
302  _is_absolute = true;
305  }
306 
307  return true;
308 }
double getFitnessRelative(Individual *ind, unsigned int patch)
Returns the relative fitness of the individual, adjusted by a scaling factor.
Definition: LCEselection.h:288
void setScalingFactorAbsolute(age_idx age, unsigned int p)
Resets the fitness scaling factor equal to one.
Definition: LCEselection.h:326
void setScalingFactorMaxGlobal(age_idx age, unsigned int p)
Sets the fitness scaling factor equal to the inverse of the maximum population fitness value.
Definition: LCEselection.cc:1242
void setScalingFactorMaxLocal(age_idx age, unsigned int p)
Sets the fitness scaling factor equal to the inverse of the maximum local patch fitness value.
Definition: LCEselection.cc:1234
void setScalingFactorGlobal(age_idx age, unsigned int p)
Sets the fitness scaling factor equal to the inverse of the mean population fitness.
Definition: LCEselection.cc:1224
void setScalingFactorLocal(age_idx age, unsigned int p)
Sets the fitness scaling factor equal to the inverse of the mean local patch fitness.
Definition: LCEselection.cc:1216
string getArg(string name)
Accessor to the parameters argument string.
Definition: param.h:300
int error(const char *str,...)
Definition: output.cc:79

References _getFitness, _is_absolute, _is_local, SimComponent::_paramSet, _setScalingFactor, error(), ParamSet::getArg(), getFitnessAbsolute(), getFitnessRelative(), ParamSet::isSet(), setScalingFactorAbsolute(), setScalingFactorGlobal(), setScalingFactorLocal(), setScalingFactorMaxGlobal(), and setScalingFactorMaxLocal().

Referenced by LCE_Selection_base(), and setParameters().

◆ set_local_optima()

bool LCE_Selection_base::set_local_optima ( )
861 {//set the traits' local optima, _selectTraitDimension must be set before that (= nbr of traits to select on)
862  string model = _paramSet->getArg("selection_model");
863 
864  if( model == "fix" || model == "direct" || model == "linear" || model == "truncation") return true;
865 
866  if( (model == "gaussian" || model == "quadratic" || model == "disruptive")
867  && !get_parameter("selection_local_optima")->isSet()) {
868 
869  return error("parameter \"selection_local_optima\" must be set to have Gaussian, quadratic, or disruptive selection.\n");
870 
871  } else {
872 
873  TMatrix tmp_mat;
874 
876 
877  _paramSet->getMatrix("selection_local_optima", &tmp_mat);
878 
879  return setSpatialMatrix("selection_local_optima", "\"selection_trait_dimension\"", &tmp_mat, &_local_optima,
880  _selectTraitDimension, _popPtr->getPatchNbr(), _paramSet->isSet("selection_randomize"));
881  }
882 
883  return true;
884 }
void getMatrix(string name, TMatrix *mat)
Accessor to the parameters matrix.
Definition: param.h:304
void reset(unsigned int rows, unsigned int cols)
Re-allocate the existing matrix with assigned rows and cols dimensions and all elements to 0.
Definition: tmatrix.h:161
bool setSpatialMatrix(string param, string numColCondition, TMatrix *inMat, TMatrix *outMat, unsigned int nVal, unsigned int patchNbr, bool doRandomize)
Definition: utils.cc:115

References _local_optima, SimComponent::_paramSet, LifeCycleEvent::_popPtr, _selectTraitDimension, error(), SimComponent::get_parameter(), ParamSet::getArg(), ParamSet::getMatrix(), Metapop::getPatchNbr(), ParamSet::isSet(), TMatrix::reset(), and setSpatialMatrix().

Referenced by checkChangeLocalOptima(), LCE_Selection_base(), and setParameters().

◆ set_param_rate_of_change()

bool LCE_Selection_base::set_param_rate_of_change ( )
889 {
890  _rate_of_change_is_std = false;
891  _do_change_local_opt = false;
892  _setNewLocalOptima = 0;
894 
895  if (!get_parameter("selection_rate_environmental_change")->isSet()
896  && !get_parameter("selection_std_rate_environmental_change")->isSet() ) {
897  return true;
898  }
899 
900  if (get_parameter("selection_rate_environmental_change")->isSet()
901  && get_parameter("selection_std_rate_environmental_change")->isSet() ) {
902  return error("both \"selection_rate_environmental_change\" and \"selection_std_rate_environmental_change\" are set, need only one.\n");
903  }
904 
905  TMatrix tmpMat;
906 
907  if (get_parameter("selection_rate_environmental_change")->isSet()) {
908 
909  if(!get_parameter("selection_rate_environmental_change")->isMatrix()) {
910 
911  double val = get_parameter_value("selection_rate_environmental_change");
912 
913  tmpMat.reset(1, _selectTraitDimension);
914  tmpMat.assign(val);
915 
916  } else {
917  get_parameter("selection_rate_environmental_change")->getMatrix(&tmpMat);
918  }
919 
921 
922  } else if (get_parameter("selection_std_rate_environmental_change")->isSet()){
923 
924  if(!get_parameter("selection_std_rate_environmental_change")->isMatrix()) {
925 
926  double val = get_parameter_value("selection_std_rate_environmental_change");
927 
928  tmpMat.reset(1, _selectTraitDimension);
929  tmpMat.assign(val);
930 
931  } else {
932  get_parameter("selection_std_rate_environmental_change")->getMatrix(&tmpMat);
933  }
934 
935  _rate_of_change_is_std = true;
936 
937  if(get_parameter("selection_std_rate_set_at_generation")->isSet())
938  _set_std_rate_at_generation = (unsigned int)get_parameter_value("selection_std_rate_set_at_generation");
939  else
941 
942  //check if phenotypic SD is to be computed in a single reference patch
943  //is -1 if parameter not set, which corresponds to the whole population then
944  _std_rate_reference_patch = get_parameter_value("selection_std_rate_reference_patch");
945 
947  }
948 
949  if(tmpMat.nrows() > _popPtr->getPatchNbr())
950  return error("The matrix of rate of change in local optima has more rows than the number of patches\n");
951 
952  if((int)tmpMat.ncols() > _selectTraitDimension)
953  return error("The matrix of rate of change in local optima has more columns than the number of quantitative traits\n");
954 
955  _do_change_local_opt = true;
956 
958 
959  unsigned int nCol = tmpMat.ncols(), nRow = tmpMat.nrows();
960 
961  //copy values, with pattern propagation
962  for (int p = 0; p < _popPtr->getPatchNbr(); ++p) {
963  for (int i = 0; i < _selectTraitDimension; ++i) {
964  _rate_of_change_local_optima.set(p, i, tmpMat.get(p % nRow, i % nCol));
965  }
966  }
967 
968  return true;
969 }
void changeLocalOptima()
Definition: LCEselection.cc:1295
virtual double get_parameter_value(std::string name)
Param value getter.
Definition: simcomponent.h:143
unsigned int ncols() const
Definition: tmatrix.h:216
void set(unsigned int i, unsigned int j, double val)
Sets element at row i and column j to value val.
Definition: tmatrix.h:103
void assign(double val)
Assigns a value to all element of the matrix.
Definition: tmatrix.h:155

References _do_change_local_opt, LifeCycleEvent::_popPtr, _rate_of_change_is_std, _rate_of_change_local_optima, _selectTraitDimension, _set_std_rate_at_generation, _setNewLocalOptima, _std_rate_reference_patch, TMatrix::assign(), changeLocalOptima(), error(), TMatrix::get(), SimComponent::get_parameter(), SimComponent::get_parameter_value(), Param::getMatrix(), Metapop::getPatchNbr(), TMatrix::ncols(), TMatrix::nrows(), TMatrix::reset(), TMatrix::set(), and set_std_rate_of_change().

Referenced by LCE_Selection_base(), and setParameters().

◆ set_sel_model()

bool LCE_Selection_base::set_sel_model ( )
313 {
314 
315  // _selectTraitDimension is set per-model in gaussian/quadratic via setSelectTraitMapping()
317  _selectTraitMapping.clear();
318  _selectTraitMapping.push_back(0);
319 
320  // -------------------------------------------------------------------------------------
321  // ENVIRONMENTAL VARIANCE
322 
323  if(get_parameter("selection_environmental_variance")->isSet())
324  //the variable actually holds the standard dev...
325  _eVariance = sqrt(get_parameter_value("selection_environmental_variance"));
326  else
327  _eVariance = 0;
328 
329  // -------------------------------------------------------------------------------------
330  // SELECTION MODELS
331 
332  // empty containers:
333  _SelectionModels.clear();
334  _getRawFitness.clear();
335 
336  // default case: direct fitness model
337  if(!_paramSet->isSet("selection_model")) {
338 
340 
341  _SelectionModels.push_back("direct");
342 
343  } else {
344 
345  _SelectionModels = get_parameter("selection_model")->getMultiArgs();
346  } //end_if isSet()
347 
348  // check consistency with selection_trait info
349  if (_SelectionModels.size() != _TraitIndices.size())
350  return error("\"selection_trait\" and \"selection_model\" must have the same number of arguments.");
351 
352  // -------------------------------------------------------------------------------------
353  // SET FUNCTION POINTERS TO FITNESS FUNCTIONS
354 
355  string sel_model;
356 
357  for(unsigned int t = 0; t < _SelectionModels.size(); t++) {
358 
359  sel_model = _SelectionModels[t];
360 
361 
362  // FIX -----------------------------------------------------------------------
363  if(sel_model == "fix") {
364 
365  //check if the trait type in the trait table is correct
366  if(_Traits[t] != "fix")
367  return error("the selection model \"fix\" does not match with the trait type \"%s\"\n", _Traits[t].c_str());
368 
369  if(!get_parameter("selection_lethal_equivalents")->isSet()) {
370 
371  return error("\"selection_lethal_equivalents\" parameter is missing with \"fix\" selection!\n");
372 
373  } else
374  _letheq = get_parameter_value("selection_lethal_equivalents");
375 
376  if(!get_parameter("selection_base_fitness")->isSet()) {
377 
378  warning("\"selection_base_fitness\" parameter is missing under fix selection model, setting it to 1.\n");
379 
380  _base_fitness = 1.0;
381 
382  } else {
383 
384  _base_fitness = get_parameter_value("selection_base_fitness");
385 
386  }
387 
388  if(!get_parameter("selection_pedigree_F")->isSet()) {
389 
390  return error("\"selection_pedigree_F\" parameter is missing with \"fix\" selection!\n");
391 
392  } else {
393 
394  TMatrix tmp_mat;
395 
396  get_parameter("selection_pedigree_F")->getMatrix(&tmp_mat);
397 
398  if(tmp_mat.getNbCols() != 5) {
399  return error("\"selection_pedigree_F\" must be an array of size 5.\n");
400  }
401 
402  for(unsigned int i = 0; i < 5; i++) _Fpedigree[i] = tmp_mat.get(0,i);
403  }
404 
405  for(unsigned int i = 0; i < 5; i++)
406  _FitnessFixModel[i] = _base_fitness * exp( -_letheq * _Fpedigree[i] );
407 
408  //everything has been set correctly, add function pointer to container
410 
411  // DIRECT -----------------------------------------------------------------------
412  } else if(sel_model == "direct") {
413 
414  //check if the trait type in the trait table is correct
415  if(_Traits[t] != "delet" && _Traits[t] != "dmi")
416  return error("the selection model \"direct\" does not match with the trait type \"%s\"\n", _Traits[t].c_str());
417 
419 
420  // QUADRATIC --------------------------------------------------------------------
421  } else if(sel_model == "quadratic") {
422 
423  //check if the trait type in the trait table is correct
424  if(_Traits[t] != "quant")
425  return error("the selection model \"quadratic\" works only with the \"quant\" trait, your trait %i has type \"%s\"\n",
426  t, _Traits[t].c_str());
427 
428  if(!setSelectTraitMapping((unsigned int)_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits")))
429  return false;
430 
431  //set the selection variance params
432  if(!setSelectionMatrix()) return error("LCE selection::failed to set the selection matrix\n");
433 
434  //now add the pointer to fitness function:
436 
437  // GAUSSIAN --------------------------------------------------------------------
438  } else if(sel_model == "gaussian") {
439 
440  //check if the trait type in the trait table is correct
441  if(_Traits[t] != "quant")
442  return error("the selection model \"gaussian\" works only with the \"quant\" trait, your trait type is \"%s\"\n", _Traits[t].c_str());
443 
444  if(!setSelectTraitMapping((unsigned int)_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits")))
445  return false;
446 
447  // set the selection matrices for multi-stage or spatially variable cases
448  if( !setSelectionMatrix() ) return error("LCE selection::failed to set the selection matrix\n");
449 
450  if(_selectTraitDimension > 1) {
451 
452  if(_eVariance > 0)
454  else
456 
457  } else {
458 
459  if(_eVariance > 0)
461  else
463 
464  }
465 
466 
467  // LINEAR --------------------------------------------------------------------
468 
469  } else if(sel_model == "linear") {
470 
471  //check if the trait type in the trait table is correct
472  if(_Traits[t] != "quant")
473  return error("the selection model \"linear\" works only with the \"quant\" trait, your trait type is \"%s\"\n", _Traits[t].c_str());
474 
475  if(!setSelectTraitMapping((unsigned int)_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits")))
476  return false;
477 
478  unsigned int patchNbr = _popPtr->getPatchNbr();
479  TMatrix tmp_matx;
480 
481  // store the quanti trait index for mean phenotype computation
483 
484  // selection gradient beta, default 1.0 (Carter et al. 2005)
486 
487  if(get_parameter("selection_linear_strength")->isSet()) {
488 
489  if(get_parameter("selection_linear_strength")->isMatrix()) {
490  get_parameter("selection_linear_strength")->getMatrix(&tmp_matx);
491 
492  if(tmp_matx.getNbRows() > patchNbr)
493  return error("\"selection_linear_strength\" must have at max as many rows as the number of patches.\n");
494 
495  if(tmp_matx.ncols() > (unsigned)_selectTraitDimension)
496  return error("\"selection_linear_strength\" must have at max as many columns as the number of selected traits.\n");
497 
499  } else {
500  // scalar value: uniform strength across patches and traits
501  _linear_selection_strength.assign(get_parameter_value("selection_linear_strength"));
502  }
503  }
504 
505  // offset defaults to 1.0, giving W = 1 + beta*(z - z_bar) (Carter et al. 2005)
506  if(!setSelectionOffset(1.0, 0.0)) return false;
507 
509 
510  // age class for computing mean phenotype: offspring (default) or adults
512  if(get_parameter("selection_at_stage")->isSet()) {
513  string age_str = get_parameter("selection_at_stage")->getArg();
514  if(age_str == "adults")
516  else if(age_str == "offspring")
518  else
519  return error("\"selection_at_stage\" must be \"offspring\" or \"adults\".\n");
520  }
521 
522  // override scaling factor to compute per-patch mean phenotypes before fitness evaluation
524 
526 
527  // DISRUPTIVE ----------------------------------------------------------------
528 
529  } else if(sel_model == "disruptive") {
530 
531  //check if the trait type in the trait table is correct
532  if(_Traits[t] != "quant")
533  return error("the selection model \"disruptive\" works only with the \"quant\" trait, your trait type is \"%s\"\n", _Traits[t].c_str());
534 
535  if(!setSelectTraitMapping((unsigned int)_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits")))
536  return false;
537 
538  // reuse the selection matrix/variance infrastructure (same as Gaussian)
539  if(!setSelectionMatrix()) return error("LCE selection::failed to set the selection matrix\n");
540 
541  unsigned int patchNbr = _popPtr->getPatchNbr();
542  TMatrix tmp_matx;
543 
544  // offset = asymptotic fitness far from optimum, default 1.0
545  // W = offset * (1 - exp(-d^2 / (2*Vs)))
546  if(!setSelectionOffset(1.0, 0.0)) return false;
547 
548  if(_selectTraitDimension > 1)
550  else
552 
553  // TRUNCATION ----------------------------------------------------------------
554 
555  } else if(sel_model == "truncation") {
556 
557  if(_Traits[t] != "quant")
558  return error("the selection model \"truncation\" works only with the \"quant\" trait, your trait type is \"%s\"\n", _Traits[t].c_str());
559 
560  if(!setSelectTraitMapping((unsigned int)_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits")))
561  return false;
562 
563  unsigned int patchNbr = _popPtr->getPatchNbr();
564  TMatrix tmp_matx;
565 
566  // threshold is mandatory
567  if(!get_parameter("selection_truncation_threshold")->isSet())
568  return error("\"selection_truncation_threshold\" must be set with selection model \"truncation\".\n");
569 
571 
572  if(get_parameter("selection_truncation_threshold")->isMatrix()) {
573  get_parameter("selection_truncation_threshold")->getMatrix(&tmp_matx);
574 
575  if(tmp_matx.getNbRows() > patchNbr)
576  return error("\"selection_truncation_threshold\" must have at max as many rows as the number of patches.\n");
577 
578  if(tmp_matx.ncols() > (unsigned)_selectTraitDimension)
579  return error("\"selection_truncation_threshold\" must have at max as many columns as the number of selected traits.\n");
580 
582  } else {
583  _truncation_threshold.assign(get_parameter_value("selection_truncation_threshold"));
584  }
585 
586  // direction: "upper" (default) selects for z >= threshold; "lower" selects for z <= threshold
587  _truncation_upper = true;
588  if(get_parameter("selection_truncation_direction")->isSet()) {
589  string dir_str = get_parameter("selection_truncation_direction")->getArg();
590  if(dir_str == "lower")
591  _truncation_upper = false;
592  else if(dir_str == "upper")
593  _truncation_upper = true;
594  else
595  return error("\"selection_truncation_direction\" must be \"upper\" or \"lower\".\n");
596  }
597 
598  // offset: default 0 (hard truncation: W=0 below, W=1 above)
599  if(!setSelectionOffset(0.0, 0.0)) return false;
600 
602 
603  } else return error("wrong selection model, must be either \"fix\", \"direct\", \"quadratic\", \"gaussian\", \"linear\", \"disruptive\", or \"truncation\".\n");
604 
605  }
606 
607  if(sel_model != "fix" && !_paramSet->isSet("selection_trait"))
608  return error("trait under selection is not set, please add parameter \"selection_trait\"\n");
609 
610  return true;
611 }
TraitPrototype * getTraitPrototype(trait_t type)
Accessor to a TraitPrototype.
Definition: indfactory.cc:140
double getFitnessUnivariateGaussian_VE(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following the Gaussian selection model with several traits under...
Definition: LCEselection.cc:1097
double getFitnessMultivariateDisruptive(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following a Disruptive (inverted Gaussian) selection model with ...
Definition: LCEselection.cc:1002
bool setSelectionMatrix()
Definition: LCEselection.cc:661
double getFitnessFixedEffect(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual in the fixed selection model.
Definition: LCEselection.h:236
double getFitnessDirect(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following the direct selection model.
Definition: LCEselection.h:242
double getFitnessUnivariateQuadratic(Individual *ind, unsigned int patch, unsigned int trait)
Quadratic fitness surface, approximates the Gaussian model for weak selection and/or small deviation ...
Definition: LCEselection.cc:1042
double getFitnessTruncation(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following a Truncation selection model.
Definition: LCEselection.cc:1021
bool setSelectionOffset(double default_val, double min_val)
Definition: LCEselection.cc:824
double getFitnessMultivariateGaussian_VE(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following the Gaussian selection model with one trait under sele...
Definition: LCEselection.cc:1113
double getFitnessMultivariateGaussian(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following the Gaussian selection model with one trait under sele...
Definition: LCEselection.cc:1077
void setScalingFactorForLinearSelection(age_idx age, unsigned int p)
Sets the fitness scaling factor as the mean trait value in the patch to compute the fitness value in ...
Definition: LCEselection.cc:1251
double _Fpedigree[5]
Array of pedigree values used in the fixed-selection model.
Definition: LCEselection.h:96
bool setSelectTraitMapping(unsigned int num_quanti_traits)
Definition: LCEselection.cc:615
vector< string > _SelectionModels
The selection models associated with each trait under selection.
Definition: LCEselection.h:143
double getFitnessUnivariateLinear(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following a Linear selection model with a single trait under sel...
Definition: LCEselection.cc:973
vector< string > _Traits
The list of trait types under selection.
Definition: LCEselection.h:137
double getFitnessUnivariateGaussian(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following the Gaussian selection model with several traits under...
Definition: LCEselection.cc:1062
double getFitnessUnivariateDisruptive(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following a Disruptive (inverted Gaussian) selection model with ...
Definition: LCEselection.cc:990
string getArg()
Definition: param.h:138
vector< string > getMultiArgs()
Definition: param.cc:120
unsigned int getNbRows() const
Gives the number of rows.
Definition: tmatrix.h:212
unsigned int getNbCols() const
Gives the number of columns.
Definition: tmatrix.h:215
void copy_recycle(const TMatrix &mat)
Copy elements of 'mat', recycling elements of 'mat' if its size is smaller than current matrix.
Definition: tmatrix.h:90
void warning(const char *str,...)
Definition: output.cc:58
@ ADLTx
Definition: types.h:42

References _base_fitness, _eVariance, _FitnessFixModel, _Fpedigree, _getRawFitness, _letheq, _linear_selection_strength, _linearSelAgeClass, _linearSelTraitIndex, _mean_trait_value, SimComponent::_paramSet, LifeCycleEvent::_popPtr, _SelectionModels, _selectTraitDimension, _selectTraitMapping, _setScalingFactor, _TraitIndices, _Traits, _truncation_threshold, _truncation_upper, ADLTx, TMatrix::assign(), TMatrix::copy_recycle(), error(), TMatrix::get(), SimComponent::get_parameter(), SimComponent::get_parameter_value(), Param::getArg(), getFitnessDirect(), getFitnessFixedEffect(), getFitnessMultivariateDisruptive(), getFitnessMultivariateGaussian(), getFitnessMultivariateGaussian_VE(), getFitnessTruncation(), getFitnessUnivariateDisruptive(), getFitnessUnivariateGaussian(), getFitnessUnivariateGaussian_VE(), getFitnessUnivariateLinear(), getFitnessUnivariateQuadratic(), Param::getMatrix(), Param::getMultiArgs(), TMatrix::getNbCols(), TMatrix::getNbRows(), Metapop::getPatchNbr(), IndFactory::getTraitPrototype(), ParamSet::isSet(), TMatrix::ncols(), OFFSx, TMatrix::reset(), setScalingFactorForLinearSelection(), setSelectionMatrix(), setSelectionOffset(), setSelectTraitMapping(), and warning().

Referenced by LCE_Selection_base(), and setParameters().

◆ set_std_rate_of_change()

void LCE_Selection_base::set_std_rate_of_change ( )
1429 {
1431 
1432  //reset rate_of_change matrix in each new replicate:
1433  TMatrix tmpMat;
1434 
1435  if(!get_parameter("selection_std_rate_environmental_change")->isMatrix()) {
1436 
1437  double val = get_parameter_value("selection_std_rate_environmental_change");
1438 
1440  tmpMat.assign(val);
1441 
1442  } else {
1443  get_parameter("selection_std_rate_environmental_change")->getMatrix(&tmpMat);
1444  }
1445 
1446  unsigned int nCol = tmpMat.ncols(), nRow = tmpMat.nrows();
1447 
1449 
1450  //copy values, with pattern propagation
1451  for (int p = 0; p < _popPtr->getPatchNbr(); ++p) {
1452  for (int i = 0; i < _selectTraitDimension; ++i) {
1453  _rate_of_change_local_optima.set(p, i, tmpMat.get(p % nRow, i % nCol));
1454  }
1455  }
1456 
1457  double* SD = new double [_selectTraitDimension];
1458 
1459  for (int i = 0; i < _selectTraitDimension; ++i) {
1460  SD[i] = 0;
1461  }
1462 
1463  // check if SD is taken in a reference patch instead of the whole pop
1464  if (_std_rate_reference_patch > -1) {
1465 
1468 
1469  } else {
1470  // compute SD as the mean within-patch SD
1471 
1472  unsigned int cnt = 0;
1473 
1474  //get SD only in extant demes, to avoid nans
1475  for(unsigned int i = 0; i < _popPtr->getPatchNbr(); ++i) {
1476  if (_popPtr->getPatch(i)->size(OFFSx) != 0 ) {
1477  cnt++;
1478  addPhenotypicSD(i, SD);
1479  }
1480  }
1481 
1482  //compute mean within-patch phenotypic standard deviation:
1483  for (int i = 0; i < _selectTraitDimension; ++i) {
1484  SD[i]/= cnt;
1485 
1486  }
1487  } //end if
1488 
1489  //multiply the per-trait rates of change by SD:
1490  for (int p = 0; p < _popPtr->getPatchNbr(); ++p) {
1491  for (int i = 0; i < _selectTraitDimension; ++i){
1492  _rate_of_change_local_optima.multi(p, i, SD[i]);
1493  }
1494 
1495  }
1496  //log the rates in the simulation log file:
1497  SIMenv::MainSim->_FileServices.log("#selection_rate_environmental_change " +
1499 
1500  //compute the change of local optima for current generation:
1502 
1503  //now reset the function pointer to changeLocalOptima() for next generation:
1505 
1506  delete [] SD;
1507  }
1508 }
void log(string message)
Write to the parameter logfile.
Definition: fileservices.cc:436
void addPhenotypicSD(unsigned int deme, double *stDev)
Definition: LCEselection.cc:1512
static SimRunner * MainSim
Definition: simenv.h:42
FileServices _FileServices
Definition: simulation.h:100
void multi(unsigned int i, unsigned int j, double value)
Multiply an element of the matrix by a value.
Definition: tmatrix.h:281
string to_string()
Writes the matrix into a string in Nemo's matrix input format.
Definition: tmatrix.h:338

References SimRunner::_FileServices, LifeCycleEvent::_popPtr, _rate_of_change_local_optima, _selectTraitDimension, _set_std_rate_at_generation, _setNewLocalOptima, _std_rate_reference_patch, addPhenotypicSD(), TMatrix::assign(), changeLocalOptima(), TMatrix::get(), SimComponent::get_parameter(), SimComponent::get_parameter_value(), Metapop::getCurrentGeneration(), Param::getMatrix(), Metapop::getPatch(), Metapop::getPatchNbr(), FileServices::log(), SIMenv::MainSim, TMatrix::multi(), TMatrix::ncols(), TMatrix::nrows(), OFFSx, TMatrix::reset(), TMatrix::set(), Patch::size(), and TMatrix::to_string().

Referenced by checkChangeLocalOptima(), and set_param_rate_of_change().

◆ setLocalOptima()

bool LCE_Selection_base::setLocalOptima ( )

◆ setMeanFitness()

double LCE_Selection_base::setMeanFitness ( age_idx  age)
inline

Sets the _mean_fitness variable to the value of the mean population fitness.

Parameters
agethe age-class index
224 { return (_mean_fitness = getMeanFitness(age));}
double getMeanFitness(age_idx age)
Computes the mean fitness of the whole population for a given age class.
Definition: LCEselection.cc:1146

References _mean_fitness, and getMeanFitness().

◆ setMeans()

void LCE_Selection_base::setMeans ( unsigned int  tot_ind)

Computes the average fitness of each pedigree class.

Called after selection has been done in the whole population.

Parameters
tot_indcount of the total number of individuals produced
1401 {
1402  _mean_fitness = 0;
1403 
1404  for(unsigned int i = 0; i < 5; i++) {
1405  _mean_fitness += _fitness[i];
1406  _fitness[i] /= _ind_cntr[i];
1407  _survival[i] /= _ind_cntr[i];
1408  _ind_cntr[i] /= tot_ind;
1409  }
1410 
1411  _mean_fitness /= tot_ind;
1412 }

References _fitness, _ind_cntr, _mean_fitness, and _survival.

Referenced by LCE_Breed_Selection::execute(), LCE_Breed_Selection_Disperse::execute(), and execute().

◆ setParameters()

bool LCE_Selection_base::setParameters ( )
virtual

Implements SimComponent.

Reimplemented in LCE_Breed_Selection_Disperse, and LCE_Breed_Selection.

187 {
188  //selection may act on more than one trait at a time
189  _Traits.clear();
190  _Traits = get_parameter("selection_trait")->getMultiArgs();
191 
192  _TraitIndices.clear();
193 
194  for(unsigned int i = 0; i < _Traits.size(); i++) {
195 
196  if(_popPtr->getTraitIndex(_Traits[i].c_str()) == -1) {
197 
198  return error("cannot attach trait \"%s\" to life cycle event \"%s\", trait has not been initiated.\n",
199  _Traits[i].c_str(), _event_name.c_str());
200 
201  } else {
202  _TraitIndices.push_back(_popPtr->getTraitIndex(_Traits[i].c_str()));
203 
204 
205 #ifdef _DEBUG_
206  cout<<"#### attaching trait \""<<_Traits[i].c_str()<<"\" to selection LCE, index = "<<_TraitIndices[i]<<endl;
207 #endif
208 
209  }
210 
211  }
212 
213 #ifdef _DEBUG_
214  cout<<"#### setting fitness model "<<endl;
215 #endif
216 
217  if(!set_fit_model()) return false;
218 
219 #ifdef _DEBUG_
220  cout<<"#### setting selection model "<<endl;
221 #endif
222 
223  if(!set_sel_model()) return false;
224 
225 #ifdef _DEBUG_
226  cout<<"#### setting local optima "<<endl;
227 #endif
228 
229  if(!set_local_optima()) return false;
230 
231 #ifdef _DEBUG_
232  cout<<"#### setting rate of environmental change "<<endl;
233 #endif
234 
235  if(!set_param_rate_of_change()) return false;
236 
238  _scaling_factor = 1;
239  resetCounters();
240 
241  assert(_TraitIndices.size() == _getRawFitness.size());
242 
243 #ifdef _DEBUG_
244  cout<<"#### selection is acting on "<<_TraitIndices.size()<<" traits\n";
245  cout<<"#### selection is using "<< _getRawFitness.size() <<" fitness functions\n";
246 #endif
247 
248 
249  return true;
250 }
int getTraitIndex(trait_t type)
Gives the index of trait with type.
Definition: indfactory.cc:128
std::string _event_name
The param name to be read in the init file.
Definition: lifecycleevent.h:77

References LifeCycleEvent::_event_name, _getRawFitness, _max_fitness, _mean_fitness, LifeCycleEvent::_popPtr, _scaling_factor, _TraitIndices, _Traits, error(), SimComponent::get_parameter(), Param::getMultiArgs(), IndFactory::getTraitIndex(), resetCounters(), set_fit_model(), set_local_optima(), set_param_rate_of_change(), and set_sel_model().

Referenced by LCE_Breed_Selection::setParameters(), and LCE_Breed_Selection_Disperse::setParameters().

◆ setScalingFactorAbsolute()

void LCE_Selection_base::setScalingFactorAbsolute ( age_idx  age,
unsigned int  p 
)
inline

Resets the fitness scaling factor equal to one.

Parameters
agethe age-class index
pthe focal patch index
326 { _scaling_factor = 1; }

References _scaling_factor.

Referenced by set_fit_model().

◆ setScalingFactorForLinearSelection()

void LCE_Selection_base::setScalingFactorForLinearSelection ( age_idx  age,
unsigned int  p 
)

Sets the fitness scaling factor as the mean trait value in the patch to compute the fitness value in the linear selection model.

Parameters
agethe age-class index
pthe focal patch index
1252 {
1253  // use the configured age class (overrides the caller's age argument)
1254  age = _linearSelAgeClass;
1255 
1256  // compute the mean phenotypic value per selected trait dimension from the specified age class
1258 
1259  Patch *patch = _popPtr->getPatch(p);
1260  Individual* ind;
1261  double* phe;
1262 
1263  unsigned int sizeFEM = patch->size(FEM, age),
1264  sizeMAL = patch->size(MAL, age);
1265  unsigned int patchSize = sizeFEM + sizeMAL;
1266 
1267  if(patchSize == 0) {
1268  _scaling_factor = 1;
1269  return;
1270  }
1271 
1272  for(unsigned int j = 0; j < sizeFEM; j++) {
1273  ind = patch->get(FEM, age, j);
1274  phe = (double*)ind->getTraitValue(_linearSelTraitIndex);
1275  for(int i = 0; i < _selectTraitDimension; i++)
1276  _mean_trait_value[i] += phe[_selectTraitMapping[i]];
1277  }
1278 
1279  for(unsigned int j = 0; j < sizeMAL; j++) {
1280  ind = patch->get(MAL, age, j);
1281  phe = (double*)ind->getTraitValue(_linearSelTraitIndex);
1282  for(int i = 0; i < _selectTraitDimension; i++)
1283  _mean_trait_value[i] += phe[_selectTraitMapping[i]];
1284  }
1285 
1286  for(int i = 0; i < _selectTraitDimension; i++)
1287  _mean_trait_value[i] /= patchSize;
1288 
1289  _scaling_factor = 1;
1290 }

References _linearSelAgeClass, _linearSelTraitIndex, _mean_trait_value, LifeCycleEvent::_popPtr, _scaling_factor, _selectTraitDimension, _selectTraitMapping, FEM, Patch::get(), Metapop::getPatch(), Individual::getTraitValue(), MAL, and Patch::size().

Referenced by set_sel_model().

◆ setScalingFactorGlobal()

void LCE_Selection_base::setScalingFactorGlobal ( age_idx  age,
unsigned int  p 
)

Sets the fitness scaling factor equal to the inverse of the mean population fitness.

Function exits on p != 0; the mean population fitness is computed only once in the execute() procedure, at the beginning of the patch loop.

Parameters
agethe age-class index
pthe focal patch index
1225 {
1226  if(p != 0) return; //we compupte the mean fitness only once
1227  _scaling_factor = 1;
1228  _scaling_factor = 1.0/getMeanFitness(age);
1229 }

References _scaling_factor, and getMeanFitness().

Referenced by set_fit_model(), and LCE_SelectionSH::setDataTable().

◆ setScalingFactorLocal()

void LCE_Selection_base::setScalingFactorLocal ( age_idx  age,
unsigned int  p 
)

Sets the fitness scaling factor equal to the inverse of the mean local patch fitness.

Parameters
agethe age-class index
pthe focal patch index
1217 {
1218  _scaling_factor = 1; //this to have the raw mean fitness below
1219  _scaling_factor = 1.0/getMeanPatchFitness(age, p);
1220 }
double getMeanPatchFitness(age_idx age, unsigned int p)
Computes the mean fitness in a given patch for a given age class.
Definition: LCEselection.cc:1164

References _scaling_factor, and getMeanPatchFitness().

Referenced by set_fit_model(), and LCE_SelectionSH::setDataTable().

◆ setScalingFactorMaxGlobal()

void LCE_Selection_base::setScalingFactorMaxGlobal ( age_idx  age,
unsigned int  p 
)

Sets the fitness scaling factor equal to the inverse of the maximum population fitness value.

Function exits on p != 0; the mean population fitness is computed only once in the execute() procedure, at the beginning of the patch loop.

Parameters
agethe age-class index
pthe focal patch index
1243 {
1244  if(p != 0) return; //we compute the max fitness only once
1245  _scaling_factor = 1;
1246  _scaling_factor = 1.0/getMaxFitness(age);
1247 }
double getMaxFitness(age_idx age)
Computes the maximum fitness value of the whole population for a given age class.
Definition: LCEselection.cc:1180

References _scaling_factor, and getMaxFitness().

Referenced by set_fit_model().

◆ setScalingFactorMaxLocal()

void LCE_Selection_base::setScalingFactorMaxLocal ( age_idx  age,
unsigned int  p 
)

Sets the fitness scaling factor equal to the inverse of the maximum local patch fitness value.

Parameters
agethe age-class index
pthe focal patch index
1235 {
1236  _scaling_factor = 1; //this to have the raw fitness when computing fitnesses below
1237  _scaling_factor = 1.0/getMaxPatchFitness(age, p);
1238 }

References _scaling_factor, and getMaxPatchFitness().

Referenced by set_fit_model().

◆ setSelectionMatrix()

bool LCE_Selection_base::setSelectionMatrix ( )
662 {
663  TMatrix tmp_mat;
664  unsigned int patchNbr = _popPtr->getPatchNbr();
665 
666  if(!get_parameter("selection_matrix")->isSet() && !get_parameter("selection_variance")->isSet())
667  return error("\"selection_matrix\" or \"selection_variance\" must be set with selection model = \"gaussian\".\n");
668 
669 
670  if(get_parameter("selection_variance")->isSet() && !get_parameter("selection_trait_dimension")->isSet())
671  return error("parameter \"selection_trait_dimension\" is missing!\n");
672 
673  //clear the selection matrix container
674  vector< TMatrix* >::iterator selIT = _selection_matrix.begin();
675 
676  for(; selIT != _selection_matrix.end(); ++selIT) if((*selIT)) delete (*selIT);
677 
678  _selection_matrix.clear();
679 
680 
681  if(get_parameter("selection_matrix")->isSet()) {
682 
683  //selection matrix provided, same selection surface in each patch
684  _paramSet->getMatrix("selection_matrix", &tmp_mat);
685 
686  if(tmp_mat.getNbCols() != tmp_mat.getNbRows())
687  return error("\"selection_matrix\" must be a square matrix!\n");
688 
689  if( get_parameter("selection_trait_dimension")->isSet() ) {
690 
691  if( (unsigned int)get_parameter_value("selection_trait_dimension") != tmp_mat.ncols() ) {
692  return error("dimensionality of \"selection_matrix\" does not match with \"selection_trait_dimension\"\n");
693  }
694 
695  } else {
696 
697  //check if the number of quanti traits modeled match with the dimensionality of the selection matrix:
698  if(_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") < tmp_mat.ncols())
699  return error("The number of quantitative traits set with \"quanti_traits\" does not match with\
700  the dimensionality of the selection matrix set with \"selection_matrix\"\n");
701  else
702  _selectTraitDimension = tmp_mat.ncols();
703 
704  }
705 
706  //we have one selection matrix per patch, copy it in the container for each patch
707  for(unsigned int i = 0; i < patchNbr; ++i)
708  _selection_matrix.push_back( new TMatrix(tmp_mat) );
709 
710  } else {
711  // selection_matrix not set in input, only selection_variance is set
712 
713  // _selectTraitDimension and _selectTraitMapping already set by setSelectTraitMapping()
714  if( get_parameter("selection_trait_dimension")->isSet() ) {
715 
716  //check if the number of quanti traits modeled match with the dimensionality:
717  if((int)_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") < _selectTraitDimension)
718  return error("The number of quantitative traits set with \"quanti_traits\" does not match with\
719  the selection trait dimensionality set with \"selection_trait_dimension\"\n");
720 
721  } // alternative case already checked at the onset
722 
723  TMatrix var_spatmat, corr_spatmat;
724 
725  //setting variance spatial matrix:
726  var_spatmat.reset(patchNbr, (unsigned)_selectTraitDimension);
727 
728  if(get_parameter("selection_variance")->isMatrix()) {
729 
730  _paramSet->getMatrix("selection_variance", &tmp_mat);
731 
732  if( !setSpatialMatrix("selection_variance","\"selection_trait_dimension\"", &tmp_mat, &var_spatmat,
733  (unsigned)_selectTraitDimension, patchNbr, _paramSet->isSet("selection_randomize") ) )
734  return error("LCE selection::failed to set patch-specific selection matrices from \"selection_variance\"\n");
735 
736  } else {
737 
738  var_spatmat.assign(get_parameter_value("selection_variance"));
739  }
740 
741  //setting correlation spatial matrix:
742  if(get_parameter("selection_correlation")->isSet() && _selectTraitDimension > 1)
743  {
744  corr_spatmat.reset(patchNbr, (unsigned)_selectTraitDimension*(_selectTraitDimension-1)/2);
745 
746  if(get_parameter("selection_correlation")->isMatrix()) {
747 
748  _paramSet->getMatrix("selection_correlation", &tmp_mat);
749 
750  if( !setSpatialMatrix("selection_correlation","the num of correlation coefficients", &tmp_mat, &corr_spatmat,
751  (unsigned)_selectTraitDimension*(_selectTraitDimension-1)/2, patchNbr, _paramSet->isSet("selection_randomize") ) )
752  return error("LCE selection::failed to set patch-specific selection matrices from \"selection_correlation\"\n");;
753 
754  } else {
755 
756  corr_spatmat.assign(get_parameter_value("selection_correlation"));
757  }
758  } else {
759  //a single trait or not set in input; matrix has a single row.
760  corr_spatmat.reset(patchNbr, 1, 0.0);
761  }
762 
763 
764  //set the selection matrix:
766  double covar;
767  unsigned int col;
768  for( unsigned int p = 0; p < patchNbr; p++) {
769  col = 0;
770  for( int i = 0; i < _selectTraitDimension; i++) {
771  tmp_mat.set(i, i, var_spatmat.get(p, i));
772  for( int j = i+1; j < _selectTraitDimension; j++) {
773  covar = corr_spatmat.get(p, col) * sqrt( var_spatmat.get(p, i) * var_spatmat.get(p, j) );
774  tmp_mat.set(i, j, covar);
775  tmp_mat.set(j, i, covar);
776  if(corr_spatmat.ncols() > 1) col++;
777  }
778  }
779  _selection_matrix.push_back(new TMatrix(tmp_mat));
780  }
781  }
782 
783  if(_selectTraitDimension > 1) {
784 
785  //selection on more than one trait
786 
787  //inverting the selection matrices:
788  if(_gsl_selection_matrix.size() != 0)
789  for(unsigned int i = 0; i < _gsl_selection_matrix.size(); ++i)
790  if(_gsl_selection_matrix[i] != NULL) gsl_matrix_free( _gsl_selection_matrix[i] );
791 
792  _gsl_selection_matrix.clear();
793 
794  for( unsigned int p = 0; p < patchNbr; p++) {
795  _selection_matrix[p]->inverse();
796 
797  _gsl_selection_matrix.push_back( gsl_matrix_alloc(_selectTraitDimension, _selectTraitDimension) );
798 
799  _selection_matrix[p]->get_gsl_matrix(_gsl_selection_matrix[p]);
800  }
801 
802  //allocate the vectors used by the fitness function:
803  if(_diffs != NULL) gsl_vector_free(_diffs);
804  _diffs = gsl_vector_alloc( _selectTraitDimension );
805 
806  if(_res1 != NULL) gsl_vector_free(_res1);
807  _res1 = gsl_vector_alloc( _selectTraitDimension );
808 
809  } else {
810 
811  //selection on one trait only
812  _selection_variance.clear();
813 
814  for(unsigned int p = 0; p < patchNbr; p++)
815  _selection_variance.push_back( _selection_matrix[p]->get(0, 0) );
816  }
817 
818 
819  return true;
820 }

References _diffs, _gsl_selection_matrix, SimComponent::_paramSet, LifeCycleEvent::_popPtr, _res1, _selection_matrix, _selection_variance, _selectTraitDimension, TMatrix::assign(), error(), TMatrix::get(), SimComponent::get_parameter(), SimComponent::get_parameter_value(), ParamSet::getMatrix(), TMatrix::getNbCols(), TMatrix::getNbRows(), Metapop::getPatchNbr(), IndFactory::getTraitPrototype(), ParamSet::isSet(), TMatrix::ncols(), TMatrix::reset(), TMatrix::set(), and setSpatialMatrix().

Referenced by set_sel_model().

◆ setSelectionOffset()

bool LCE_Selection_base::setSelectionOffset ( double  default_val,
double  min_val 
)
825 {
826  unsigned int patchNbr = _popPtr->getPatchNbr();
827 
828  _selection_offset.reset(patchNbr, _selectTraitDimension, default_val);
829 
830  if(get_parameter("selection_offset")->isSet()) {
831 
832  if(get_parameter("selection_offset")->isMatrix()) {
833  TMatrix tmp_matx;
834  get_parameter("selection_offset")->getMatrix(&tmp_matx);
835 
836  if(tmp_matx.getNbRows() > patchNbr)
837  return error("\"selection_offset\" must have at max as many rows as the number of patches.\n");
838 
839  if(tmp_matx.ncols() > (unsigned)_selectTraitDimension)
840  return error("\"selection_offset\" must have at max as many columns as the number of selected traits.\n");
841 
843  } else {
844  _selection_offset.assign(get_parameter_value("selection_offset"));
845  }
846  }
847 
848  for(unsigned int i = 0; i < patchNbr; ++i)
849  for(int j = 0; j < _selectTraitDimension; ++j)
850  if(_selection_offset.get(i, j) < min_val)
851  return error("\"selection_offset\" value %g at patch %i, trait %i is below the minimum allowed value (%g).\n",
852  _selection_offset.get(i, j), i+1, j+1, min_val);
853 
854  return true;
855 }

References LifeCycleEvent::_popPtr, _selection_offset, _selectTraitDimension, TMatrix::assign(), TMatrix::copy_recycle(), error(), TMatrix::get(), SimComponent::get_parameter(), SimComponent::get_parameter_value(), Param::getMatrix(), TMatrix::getNbRows(), Metapop::getPatchNbr(), TMatrix::ncols(), and TMatrix::reset().

Referenced by set_sel_model().

◆ setSelectTraitMapping()

bool LCE_Selection_base::setSelectTraitMapping ( unsigned int  num_quanti_traits)
616 {
617  _selectTraitMapping.clear();
618 
619  if(get_parameter("selection_trait_dimension")->isSet()) {
620 
621  if(get_parameter("selection_trait_dimension")->isMatrix()) {
622 
623  // matrix mode: array of 1-based trait indices
624  TMatrix tmp;
625  get_parameter("selection_trait_dimension")->getMatrix(&tmp);
626 
627  if(tmp.nrows() > 1)
628  return error("\"selection_trait_dimension\" must be a single row array when specifying trait indices.\n");
629 
630  for(unsigned int i = 0; i < tmp.ncols(); ++i) {
631  int idx = (int)tmp.get(0, i);
632  if(idx < 1 || idx > (int)num_quanti_traits)
633  return error("\"selection_trait_dimension\" trait index %i is out of range [1, %i].\n", idx, num_quanti_traits);
634  _selectTraitMapping.push_back((unsigned int)idx - 1); // 1-based to 0-based
635  }
636 
637  } else {
638 
639  // integer mode: number of traits, identity mapping
640  unsigned int N = (unsigned int)get_parameter_value("selection_trait_dimension");
641 
642  if(N > num_quanti_traits)
643  return error("\"selection_trait_dimension\" (%i) exceeds the number of quantitative traits (%i).\n", N, num_quanti_traits);
644 
645  for(unsigned int i = 0; i < N; ++i)
646  _selectTraitMapping.push_back(i);
647  }
648 
649  } else {
650  // default: single trait
651  _selectTraitMapping.push_back(0);
652  }
653 
655 
656  return true;
657 }

References _selectTraitDimension, _selectTraitMapping, error(), TMatrix::get(), SimComponent::get_parameter(), SimComponent::get_parameter_value(), Param::getMatrix(), TMatrix::ncols(), and TMatrix::nrows().

Referenced by set_sel_model().

◆ updateFitnessCounters()

void LCE_Selection_base::updateFitnessCounters ( double  fitness,
unsigned int  ped_class,
bool  survived 
)

Updates the fitness and survival mean counters.

Parameters
fitnessthe fitness of an individual
ped_classthe pedigree class of the individual
survivedboolean telling whether the individual survived (true) or not (false)
1417 {
1418 
1419  _fitness[ped_class] += fitness;
1420  _ind_cntr[ped_class]++;
1421 
1422  if(survived) _survival[ped_class]++;
1423 
1424 }

References _fitness, _ind_cntr, and _survival.

Referenced by doViabilitySelection(), and LCE_Breed_Selection::makeOffspringWithSelection().

Friends And Related Function Documentation

◆ LCE_Breed_Selection

friend class LCE_Breed_Selection
friend

◆ LCE_SelectionFH

friend class LCE_SelectionFH
friend

Referenced by loadFileServices().

◆ LCE_SelectionSH

friend class LCE_SelectionSH
friend

The StatHandler associated class is a friend.

Referenced by loadStatServices().

Member Data Documentation

◆ _base_fitness

◆ _diffs

◆ _do_change_local_opt

bool LCE_Selection_base::_do_change_local_opt
private

◆ _eVariance

double LCE_Selection_base::_eVariance
protected

◆ _fitness

double LCE_Selection_base::_fitness[5]
protected

◆ _FitnessFixModel

double LCE_Selection_base::_FitnessFixModel[5]
private

Absolute fitness values of the five pedigree class for the fixed selection model (lethal equivalents model).

Fitness is computed as: _FitnessFixModel[i] = _base_fitness * exp( -_letheq * _Fpedigree[i] );

Referenced by getFitnessFixedEffect(), and set_sel_model().

◆ _Fpedigree

double LCE_Selection_base::_Fpedigree[5]
private

Array of pedigree values used in the fixed-selection model.

Referenced by set_sel_model().

◆ _getFitness

double(LCE_Selection_base::* LCE_Selection_base::_getFitness) (Individual *, unsigned int)
protected

Pointer to the function returning the individual fitness.

May point to LCE_Selection_base::getFitnessAbsolute or LCE_Selection_base::getFitnessRelative.

Referenced by getFitness(), and set_fit_model().

◆ _getRawFitness

vector< double (LCE_Selection_base::* ) (Individual*, unsigned int, unsigned int) > LCE_Selection_base::_getRawFitness
protected

A vector containing pointers to fitness function related to each trait under selection.

Referenced by getFitnessAbsolute(), LCE_SelectionFH::print(), set_sel_model(), and setParameters().

◆ _gsl_selection_matrix

vector< gsl_matrix* > LCE_Selection_base::_gsl_selection_matrix
private

◆ _ind_cntr

◆ _is_absolute

bool LCE_Selection_base::_is_absolute
protected

◆ _is_local

bool LCE_Selection_base::_is_local
protected

◆ _letheq

double LCE_Selection_base::_letheq
private

Referenced by set_sel_model().

◆ _linear_selection_strength

TMatrix LCE_Selection_base::_linear_selection_strength
private

◆ _linearSelAgeClass

age_idx LCE_Selection_base::_linearSelAgeClass
private

◆ _linearSelTraitIndex

unsigned int LCE_Selection_base::_linearSelTraitIndex
private

◆ _local_optima

◆ _max_fitness

double LCE_Selection_base::_max_fitness
protected

Referenced by setParameters().

◆ _mean_fitness

◆ _mean_trait_value

vector<double> LCE_Selection_base::_mean_trait_value
private

◆ _rate_of_change_is_std

bool LCE_Selection_base::_rate_of_change_is_std
private

◆ _rate_of_change_local_optima

TMatrix LCE_Selection_base::_rate_of_change_local_optima
private

◆ _res1

◆ _scaling_factor

◆ _selection_matrix

vector< TMatrix* > LCE_Selection_base::_selection_matrix
private

◆ _selection_offset

◆ _selection_variance

vector< double > LCE_Selection_base::_selection_variance
protected

◆ _SelectionModels

vector< string > LCE_Selection_base::_SelectionModels
protected

The selection models associated with each trait under selection.

Referenced by set_sel_model().

◆ _selectTraitDimension

◆ _selectTraitMapping

◆ _set_std_rate_at_generation

unsigned int LCE_Selection_base::_set_std_rate_at_generation
private

◆ _setNewLocalOptima

void(LCE_Selection_base::* LCE_Selection_base::_setNewLocalOptima) (void)
protected

Pointer to the function used to change the local phenotypic optima or its rate of change.

May point to LCE_Selection_base::changeLocalOptima or LCE_Selection_base::set_std_rate_of_change.

Referenced by checkChangeLocalOptima(), set_param_rate_of_change(), and set_std_rate_of_change().

◆ _setScalingFactor

◆ _stater

LCE_SelectionSH* LCE_Selection_base::_stater
protected

◆ _std_rate_reference_patch

int LCE_Selection_base::_std_rate_reference_patch
private

◆ _survival

◆ _TraitIndices

vector< unsigned int > LCE_Selection_base::_TraitIndices
protected

◆ _Traits

vector< string > LCE_Selection_base::_Traits
protected

The list of trait types under selection.

Referenced by set_sel_model(), and setParameters().

◆ _truncation_threshold

TMatrix LCE_Selection_base::_truncation_threshold
private

◆ _truncation_upper

bool LCE_Selection_base::_truncation_upper
private

◆ _writer

LCE_SelectionFH* LCE_Selection_base::_writer
protected

Referenced by loadFileServices().


The documentation for this class was generated from the following files:

Generated for Nemo v2.4.0b by  doxygen 1.9.1 -- Nemo is hosted on  Download Nemo

Locations of visitors to this page
Catalogued on GSR