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 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 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< 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
 
TMatrix _linear_selection_offset
 
vector< double > _mean_trait_value
 
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),
56 _letheq(0),
57 _base_fitness(1),
58 _mean_fitness(0),
59 _max_fitness(0),
61 _is_local(0),
62 _is_absolute(1),
63 _eVariance(0),
65 _getFitness(0),
68 _stater(0),
69 _writer(0)
70 {
71  // mandatory parameter (no default):
72  add_parameter("selection_trait", STR, true, false, 0, 0); //no updaters here, to keep it safe...
73 
76 
77  add_parameter("selection_fitness_model",STR,false,false,0,0, updater);
78 
80 
81  add_parameter("selection_model",STR,true,false,0,0, updater); //mandatory, depends on trait
82  add_parameter("selection_matrix",MAT,false,false,0,0, updater);
83  add_parameter("selection_variance",DBL,false,false,0,0, updater);
84  add_parameter("selection_correlation",DBL,false,false,0,0, updater);
85  add_parameter("selection_trait_dimension",INT,false,false,0,0, updater);
86  add_parameter("selection_base_fitness",DBL,false,true,0,1, updater);
87  add_parameter("selection_lethal_equivalents",DBL,false,false,0,0, updater);
88  add_parameter("selection_pedigree_F",MAT,false,false,0,0, updater);
89  add_parameter("selection_randomize",BOOL,false,false,0,0, updater);
90  add_parameter("selection_environmental_variance", DBL, false, false, 0, 0, updater);
91 
92  // linear selection:
93  add_parameter("selection_linear_strength", DBL, false, false, 0, 0, updater);
94  add_parameter("selection_linear_offset", DBL, false, false, 0, 0, updater);
95 
97  add_parameter("selection_local_optima",DBL,false,false,0,0, updater);
98 
100  add_parameter("selection_rate_environmental_change", DBL, false, false, 0, 0, updater);
101  add_parameter("selection_std_rate_environmental_change", DBL, false, false, 0, 0, updater);
102  add_parameter("selection_std_rate_set_at_generation", INT, false, false, 0, 0, updater);
103  add_parameter("selection_std_rate_reference_patch", INT, false, false, 0, 0, updater);
104 
105 
106  add_parameter("selection_output", BOOL, false, false, 0, 0, 0);
107  add_parameter("selection_output_logtime", INT, false, false, 0, 0, 0);
108  add_parameter("selection_output_dir", STR, false, false, 0, 0, 0);
109 }
int _selectTraitDimension
Number of quantitative traits under selection.
Definition: LCEselection.h:112
double _base_fitness
Definition: LCEselection.h:100
double(LCE_Selection_base::* _getFitness)(Individual *, unsigned int)
Pointer to the function returning the individual fitness.
Definition: LCEselection.h:137
double _max_fitness
Definition: LCEselection.h:101
vector< gsl_matrix * > _gsl_selection_matrix
Definition: LCEselection.h:56
LCE_SelectionFH * _writer
Definition: LCEselection.h:151
gsl_vector * _diffs
Definition: LCEselection.h:57
bool set_param_rate_of_change()
Definition: LCEselection.cc:740
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:142
bool _is_absolute
Definition: LCEselection.h:103
double _eVariance
Evironmental variance.
Definition: LCEselection.h:118
double _scaling_factor
Definition: LCEselection.h:101
bool _do_change_local_opt
Definition: LCEselection.h:61
bool _is_local
Definition: LCEselection.h:102
gsl_vector * _res1
Definition: LCEselection.h:58
bool set_sel_model()
Definition: LCEselection.cc:304
unsigned int _set_std_rate_at_generation
Definition: LCEselection.h:63
double _mean_fitness
Definition: LCEselection.h:101
bool set_fit_model()
Definition: LCEselection.cc:246
LCE_SelectionSH * _stater
Definition: LCEselection.h:150
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:146
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:133
bool set_local_optima()
Definition: LCEselection.cc:712
double _letheq
Definition: LCEselection.h:81
vector< TMatrix * > _selection_matrix
Definition: LCEselection.h:55
bool _rate_of_change_is_std
Definition: LCEselection.h:62
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

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
114 {
115  vector< TMatrix* >::iterator selIT = _selection_matrix.begin();
116  for(; selIT != _selection_matrix.end(); ++selIT)
117  if((*selIT)) delete (*selIT);
118  _selection_matrix.clear();
119 
120  for(unsigned int i = 0; i < _gsl_selection_matrix.size(); ++i)
121  if(_gsl_selection_matrix[i]) gsl_matrix_free(_gsl_selection_matrix[i]);
122  _gsl_selection_matrix.clear();
123 
124  if(_diffs) gsl_vector_free(_diffs);
125  if(_res1) gsl_vector_free(_res1);
126  if(_stater) delete _stater;
127 }

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.

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

References NONE.

◆ addPhenotypicSD()

void LCE_Selection_base::addPhenotypicSD ( unsigned int  deme,
double *  stDev 
)
1313 {
1314  Patch *patch = _popPtr->getPatch(deme);
1315  Individual* ind;
1316  unsigned int table_size = patch->size(OFFSx);
1317  double* _phe;
1318 
1319  double **phenot = new double* [_selectTraitDimension];
1320 
1321  for (int i = 0; i < _selectTraitDimension; ++i) {
1322 
1323  phenot[i] = new double [table_size];
1324 
1325  }
1326 
1327  unsigned int pos = 0;
1328 
1329  for (unsigned int j = 0; j < patch->size(FEM, OFFSx) && pos < table_size; ++j) {
1330 
1331  ind = patch->get(FEM, OFFSx, j);
1332  _phe = (double*)ind->getTraitValue(_LCELinkedTraitIndex);
1333 
1334  for (int i = 0; i < _selectTraitDimension; ++i)
1335  phenot[i][pos] = _phe[i];
1336 
1337  pos++;
1338 
1339  }
1340 
1341  for (unsigned int j = 0; j < patch->size(MAL, OFFSx) && pos < table_size; ++j) {
1342 
1343  ind = patch->get(MAL, OFFSx, j);
1344  _phe = (double*)ind->getTraitValue(_LCELinkedTraitIndex);
1345 
1346  for (int i = 0; i < _selectTraitDimension; ++i)
1347  phenot[i][pos] = _phe[i];
1348 
1349  pos++;
1350 
1351  }
1352 
1353  assert(pos == table_size);
1354 
1355  for (int i = 0; i < _selectTraitDimension; ++i) {
1356  stDev[i] += sqrt( my_variance_with_fixed_mean( phenot[i], table_size, my_mean(phenot[i], table_size) ) );
1357  delete [] phenot[i];
1358  }
1359 
1360  delete [] phenot;
1361 }
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
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
@ OFFSx
Definition: types.h:42
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, 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 ( )
1096 {
1097  for (int i = 0; i < _selectTraitDimension; ++i)
1098  for (unsigned int j = 0; j < _local_optima.nrows(); ++j) {
1100  }
1101 }
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.

1106 {
1107  // check if change in local optima for quantitative traits
1108  if(_do_change_local_opt) {
1109 
1110  if(_popPtr->getCurrentGeneration() == 1) {
1111 
1112  set_local_optima(); //reset local optima to initial values
1113 
1114  if(_rate_of_change_is_std) //reset rate of change relative to SD of that replicate
1116  }
1117 
1118  (this->*_setNewLocalOptima)();
1119  }
1120 }
void set_std_rate_of_change()
Definition: LCEselection.cc:1228
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.

320 {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
1159 {
1160  Individual* ind;
1161  double fitness;
1162  bool not_survived;
1163 
1164  for(unsigned int i = 0; i < patch->size(SEX, AGE); i++) {
1165 
1166  ind = patch->get(SEX, AGE, i);
1167 
1168  fitness = getFitness( ind, p);
1169 
1170  not_survived = ( RAND::Uniform() > fitness );
1171 
1172  if( not_survived ) {
1173 
1174  //this individual dies
1175  patch->remove(SEX, AGE, i);
1176 
1177  _popPtr->recycle(ind);
1178 
1179  i--;
1180 
1181  } //else; this individual stays in the patch
1182 
1183  updateFitnessCounters(fitness, ind->getPedigreeClass(), !not_survived);
1184  }
1185 }
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:213
void updateFitnessCounters(double fitness, unsigned int ped_class, bool survived)
Updates the fitness and survival mean counters.
Definition: LCEselection.cc:1216
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:124

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.

1125 {
1126  Patch * patch;
1127  unsigned int popSize = _popPtr->size(OFFSPRG);
1128 
1129 #ifdef _DEBUG_
1130  message("LCE_Selection_base::execute (tot offsprg nb: %i", popSize);
1131 #endif
1132 
1133  resetCounters();
1134 
1135  // check if change in local optima for quantitative traits
1137 
1138  for(unsigned int p = 0; p < _popPtr->getPatchNbr(); p++) {
1139 
1140  (this->*_setScalingFactor)(OFFSx, p);
1141 
1142  patch = _popPtr->getPatch(p);
1143 
1144  doViabilitySelection(FEM, OFFSx, patch, p);
1145 
1146  doViabilitySelection(MAL, OFFSx, patch, p);
1147  }
1148 
1149  setMeans(popSize);
1150 
1151 #ifdef _DEBUG_
1152  message(", after selection: %i, mean fitness offspring=%f)\n",_popPtr->size(OFFSPRG), _mean_fitness);
1153 #endif
1154 }
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:1158
void checkChangeLocalOptima()
Check is rate of change in local optima has been set and must be applied.
Definition: LCEselection.cc:1105
void setMeans(unsigned int tot_ind)
Computes the average fitness of each pedigree class.
Definition: LCEselection.cc:1200
void resetCounters()
Resets the fitness counters.
Definition: LCEselection.cc:1189
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
214  {
215  return (this->*_getFitness)(ind, patch);
216  }

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
934 {
935  double fitness = 1;
936  //at this point, _getRawFitness.size() == _TraitIndices.size()
937  //and we assume that the functions are "aligned" with the traits
938  for(unsigned int i = 0; i < _getRawFitness.size(); i++)
939  fitness *= (this->*_getRawFitness[i])(ind, patch, _TraitIndices[i]);
940  return fitness;
941 }
vector< unsigned int > _TraitIndices
The indices of the traits under selection.
Definition: LCEselection.h:125

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.

226  {
227  return *(double*)offsprg->getTraitValue(trait);
228  }

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.

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

References _FitnessFixModel, and Individual::getPedigreeClass().

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.

875 {
876  double res2;
877 
878 // if(!ind) fatal("passing NULL ind ptr to LCE_Selection_base::getFitnessMultivariateGaussian!!!\n");
879 
880  double *_phe = (double*)ind->getTraitValue(trait);
881 
882  for( int i = 0; i < _selectTraitDimension; i++)
883  gsl_vector_set(_diffs, i, _phe[i] - _local_optima.get(patch, i));
884 
885  //(diff)T * W * diff:
886  //right partial product:
887  gsl_blas_dsymv(CblasUpper, 1.0, _gsl_selection_matrix[patch], _diffs, 0.0, _res1);
888  //left product:
889  gsl_blas_ddot(_diffs, _res1, &res2);
890 
891  return exp( -0.5 * res2 );
892 }

References _diffs, _gsl_selection_matrix, _local_optima, _res1, _selectTraitDimension, 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.

913 {
914  double res2;
915 
916  double *_phe = (double*)ind->getTraitValue(trait);
917 
918  //add the environmental variance here:
919  for( int i = 0; i < _selectTraitDimension; i++)
920  gsl_vector_set(_diffs, i, _phe[i] + RAND::Gaussian(_eVariance) - _local_optima.get(patch, i));
921 
922  //(diff)T * W * diff:
923  //right partial product:
924  gsl_blas_dsymv(CblasUpper, 1.0, _gsl_selection_matrix[patch], _diffs, 0.0, _res1);
925  //left product:
926  gsl_blas_ddot(_diffs, _res1, &res2);
927 
928  return exp( -0.5 * res2 );
929 }
static double Gaussian(double sigma)
Definition: Uniform.h:261

References _diffs, _eVariance, _gsl_selection_matrix, _local_optima, _res1, _selectTraitDimension, 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
262  {
263  return getFitnessAbsolute(ind, patch) * _scaling_factor;
264  }
double getFitnessAbsolute(Individual *ind, unsigned int patch)
Returns the raw fitness of the individual, without adjustment (absolute fitness).
Definition: LCEselection.cc:933

References _scaling_factor, and getFitnessAbsolute().

Referenced by set_fit_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.

860 {
861  double res2, diff;
862 
863  double *_phe = (double*)ind->getTraitValue(trait);
864 
865  diff = _phe[0] - _local_optima.get(patch, 0);
866 
867  res2 = diff*diff / _selection_variance[patch];
868 
869  return exp( -0.5 * res2 );
870 }
vector< double > _selection_variance
Patch-specific selection variance.
Definition: LCEselection.h:115

References _local_optima, _selection_variance, 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.

897 {
898  double res2, diff;
899 
900  double *_phe = (double*)ind->getTraitValue(trait);
901 
902  //add the environmental variance here:
903  diff = _phe[0] + RAND::Gaussian(_eVariance) - _local_optima.get(patch, 0);
904 
905  res2 = diff*diff / _selection_variance[patch];
906 
907  return exp( -0.5 * res2 );
908 }

References _eVariance, _local_optima, _selection_variance, 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.

826 {
827  double *_phe = (double*)ind->getTraitValue(trait);
828 
829  double w = 1.0;
830 
831 // for(unsigned int t = 0; t < _selectTraitDimension; ++t)
832 // w *= _linear_selection_offset[patch] + _linear_selection_strength[patch] * (_phe[t] - _mean_absolute_fitness[patch]);
833 
834  return w;
835 }

References Individual::getTraitValue().

◆ 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.

840 {
841  double res2, diff, w = 1.0;
842 
843  double *_phe = (double*)ind->getTraitValue(trait);
844 
845  for(unsigned int t = 0; t < _selectTraitDimension; ++t){
846 
847  diff = _phe[t] - _local_optima.get(patch, t);
848 
849  res2 = 1 - diff*diff / _selection_variance[patch];
850 
851  w *= (res2 < 0 ? 0 : res2);
852  }
853 
854  return w;
855 }

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

Referenced by set_sel_model().

◆ getLocalOptima()

const TMatrix& LCE_Selection_base::getLocalOptima ( ) const
inline
161 {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
980 {
981  double max = 0, local_max;
982 
983  for(unsigned int i = 0, npatch = _popPtr->getPatchNbr(); i < npatch; i++) {
984 
985  local_max = getMaxPatchFitness(age, i);
986 
987  max = (local_max > max ? local_max : max);
988  }
989 
990  return max;
991 }
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:995

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
996 {
997  double max = 0, fit;
998  Patch *patch = _popPtr->getPatch(p);
999 
1000  for(unsigned int j = 0, size = patch->size(FEM, age); j < size; j++) {
1001  fit = getFitness( patch->get(FEM, age, j), p);
1002  max = (fit > max ? fit : max);
1003  }
1004  for(unsigned int j = 0, size = patch->size(MAL, age); j < size; j++){
1005  fit = getFitness( patch->get(MAL, age, j), p);
1006  max = (fit > max ? fit : max);
1007  }
1008 
1009  return max;
1010 }

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
946 {
947  double mean = 0;
948  Patch *patch;
949 // age_idx age = (AGE == ADULTS ? ADLTx : OFFSx);
950 
951  for(unsigned int i = 0, npatch = _popPtr->getPatchNbr(); i < npatch; i++) {
952  patch = _popPtr->getPatch(i);
953  for(unsigned int j = 0, size = patch->size(FEM, age); j < size; j++)
954  mean += getFitness( patch->get(FEM, age, j), i);
955  for(unsigned int j = 0, size = patch->size(MAL, age); j < size; j++)
956  mean += getFitness( patch->get(MAL, age, j), i);
957  }
958  return mean/_popPtr->size(age);
959 }

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
964 {
965  double mean = 0;
966  Patch *patch = _popPtr->getPatch(p);
967 
968  for(unsigned int j = 0, size = patch->size(FEM, age); j < size; j++)
969  mean += getFitness( patch->get(FEM, age, j), p);
970 
971  for(unsigned int j = 0, size = patch->size(MAL, age); j < size; j++)
972  mean += getFitness( patch->get(MAL, age, j), p);
973 
974  return mean/patch->size(age);
975 }

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.

140 {
141  if(_paramSet->isSet("selection_output")) {
142 
143  if(_writer == NULL) _writer = new LCE_SelectionFH(this);
144 
145  Param* param = get_parameter("selection_output_logtime");
146 
147  if(!param->isSet()) fatal("parameter \"selection_output_logtime\" is missing\n");
148 
149  if(param->isMatrix()) {
150 
151  TMatrix temp;
152 
153  param->getMatrix(&temp);
154 
155  _writer->set_multi(true, true, 1, &temp, get_parameter("selection_output_dir")->getArg());
156 
157  } else
158  // rpl_per, gen_per, rpl_occ, gen_occ, rank, path, self-ref
159  _writer->set(true, param->isSet(), 1, (param->isSet() ? (int)param->getValue() : 0), 0,
160  get_parameter("selection_output_dir")->getArg(), this);
161 
162 
163  loader->attach(_writer);
164 
165  } else {
166 
167  if(_writer) delete _writer;
168 
169  return;
170  }
171 
172 
173 }
virtual void set(bool rpl_per, bool gen_per, int rpl_occ, int gen_occ, int rank, string path, LCE *event)
Definition: filehandler.h:272
virtual void set_multi(bool rpl_per, bool gen_per, int rpl_occ, TMatrix *Occ, string path)
Definition: filehandler.h:197
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:95
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:96

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.

132 {
133  if(_stater == NULL) _stater = new LCE_SelectionSH(this);
134  loader->attach(_stater);
135 }
friend class LCE_SelectionSH
The StatHandler associated class is a friend.
Definition: LCEselection.h:94
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.

321 {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.

323 {return OFFSPRG;}

References OFFSPRG.

◆ resetCounters()

void LCE_Selection_base::resetCounters ( )

Resets the fitness counters.

1190 {
1191  for(unsigned int i = 0; i < 5; i++) {
1192  _fitness[i] = 0;
1193  _survival[i] = 0;
1194  _ind_cntr[i] = 0;
1195  }
1196 }
double _fitness[5]
Fitness counters, one for each pedigree class.
Definition: LCEselection.h:106
double _ind_cntr[5]
Definition: LCEselection.h:106
double _survival[5]
Definition: LCEselection.h:106

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.

319 {return false;}

◆ set_fit_model()

bool LCE_Selection_base::set_fit_model ( )
247 {
248 
249  if(_paramSet->isSet("selection_fitness_model")) {
250 
251  string fit_model = _paramSet->getArg("selection_fitness_model");
252 
253  if( fit_model == "absolute") {
254 
255  _is_local = false;
256  _is_absolute = true;
259 
260  } else if( fit_model == "relative_local" ) {
261 
262  _is_local = true;
263  _is_absolute = false;
266 
267  } else if(fit_model == "relative_global") {
268 
269  _is_local = false;
270  _is_absolute = false;
273 
274  } else if(fit_model == "relative_max_local") {
275 
276  _is_local = true;
277  _is_absolute = false;
280 
281  } else if(fit_model == "relative_max_global") {
282 
283  _is_local = false;
284  _is_absolute = false;
287 
288  } else {
289  return error("Unknown fitness model \"%s\"", fit_model.c_str());
290  }
291  } //default case:
292  else {
293  _is_local = false;
294  _is_absolute = true;
297  }
298 
299  return true;
300 }
double getFitnessRelative(Individual *ind, unsigned int patch)
Returns the relative fitness of the individual, adjusted by a scaling factor.
Definition: LCEselection.h:261
void setScalingFactorAbsolute(age_idx age, unsigned int p)
Resets the fitness scaling factor equal to one.
Definition: LCEselection.h:299
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:1041
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:1033
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:1023
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:1015
string getArg(string name)
Accessor to the parameters argument string.
Definition: param.h:300
int error(const char *str,...)
Definition: output.cc:77

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 ( )
713 {//set the traits' local optima, _selectTraitDimension must be set before that (= nbr of traits to select on)
714  string model = _paramSet->getArg("selection_model");
715 
716  if( model == "fix" || model == "direct" || model == "linear") return true;
717 
718  if( (model == "gaussian" || model == "quadratic")
719  && !get_parameter("selection_local_optima")->isSet()) {
720 
721  return error("parameter \"selection_local_optima\" must be set to have Gaussian or quadratic selection.\n");
722 
723  } else {
724 
725  TMatrix tmp_mat;
726 
728 
729  _paramSet->getMatrix("selection_local_optima", &tmp_mat);
730 
731  return setSpatialMatrix("selection_local_optima", "\"selection_trait_dimension\"", &tmp_mat, &_local_optima,
732  _selectTraitDimension, _popPtr->getPatchNbr(), _paramSet->isSet("selection_randomize"));
733  }
734 
735  return true;
736 }
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 ( )
741 {
742  _rate_of_change_is_std = false;
743  _do_change_local_opt = false;
744  _setNewLocalOptima = 0;
746 
747  if (!get_parameter("selection_rate_environmental_change")->isSet()
748  && !get_parameter("selection_std_rate_environmental_change")->isSet() ) {
749  return true;
750  }
751 
752  if (get_parameter("selection_rate_environmental_change")->isSet()
753  && get_parameter("selection_std_rate_environmental_change")->isSet() ) {
754  return error("both \"selection_rate_environmental_change\" and \"selection_std_rate_environmental_change\" are set, need only one.\n");
755  }
756 
757  TMatrix tmpMat;
758 
759  if (get_parameter("selection_rate_environmental_change")->isSet()) {
760 
761  if(!get_parameter("selection_rate_environmental_change")->isMatrix()) {
762 
763  double val = get_parameter_value("selection_rate_environmental_change");
764 
765  tmpMat.reset(1, _selectTraitDimension);
766  tmpMat.assign(val);
767 
768  } else {
769  get_parameter("selection_rate_environmental_change")->getMatrix(&tmpMat);
770  }
771 
773 
774  } else if (get_parameter("selection_std_rate_environmental_change")->isSet()){
775 
776  if(!get_parameter("selection_std_rate_environmental_change")->isMatrix()) {
777 
778  double val = get_parameter_value("selection_std_rate_environmental_change");
779 
780  tmpMat.reset(1, _selectTraitDimension);
781  tmpMat.assign(val);
782 
783  } else {
784  get_parameter("selection_std_rate_environmental_change")->getMatrix(&tmpMat);
785  }
786 
787  _rate_of_change_is_std = true;
788 
789  if(get_parameter("selection_std_rate_set_at_generation")->isSet())
790  _set_std_rate_at_generation = (unsigned int)get_parameter_value("selection_std_rate_set_at_generation");
791  else
793 
794  //check if phenotypic SD is to be computed in a single reference patch
795  //is -1 if parameter not set, which corresponds to the whole population then
796  _std_rate_reference_patch = get_parameter_value("selection_std_rate_reference_patch");
797 
799  }
800 
801  if(tmpMat.nrows() > _popPtr->getPatchNbr())
802  return error("The matrix of rate of change in local optima has more rows than the number of patches\n");
803 
804  if((int)tmpMat.ncols() > _selectTraitDimension)
805  return error("The matrix of rate of change in local optima has more columns than the number of quantitative traits\n");
806 
807  _do_change_local_opt = true;
808 
810 
811  unsigned int nCol = tmpMat.ncols(), nRow = tmpMat.nrows();
812 
813  //copy values, with pattern propagation
814  for (int p = 0; p < _popPtr->getPatchNbr(); ++p) {
815  for (int i = 0; i < _selectTraitDimension; ++i) {
816  _rate_of_change_local_optima.set(p, i, tmpMat.get(p % nRow, i % nCol));
817  }
818  }
819 
820  return true;
821 }
void changeLocalOptima()
Definition: LCEselection.cc:1095
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 ( )
305 {
306 
307  _selectTraitDimension = (int)get_parameter_value("selection_trait_dimension");
308 
309  // -------------------------------------------------------------------------------------
310  // ENVIRONMENTAL VARIANCE
311 
312  if(get_parameter("selection_environmental_variance")->isSet())
313  //the variable actually holds the standard dev...
314  _eVariance = sqrt(get_parameter_value("selection_environmental_variance"));
315  else
316  _eVariance = 0;
317 
318  // -------------------------------------------------------------------------------------
319  // SELECTION MODELS
320 
321  // empty containers:
322  _SelectionModels.clear();
323  _getRawFitness.clear();
324 
325  // default case: direct fitness model
326  if(!_paramSet->isSet("selection_model")) {
327 
329 
330  _SelectionModels.push_back("direct");
331 
332  } else {
333 
334  _SelectionModels = get_parameter("selection_model")->getMultiArgs();
335  } //end_if isSet()
336 
337  // check consistency with selection_trait info
338  if (_SelectionModels.size() != _TraitIndices.size())
339  return error("\"selection_trait\" and \"selection_model\" must have the same number of arguments.");
340 
341  // -------------------------------------------------------------------------------------
342  // SET FUNCTION POINTERS TO FITNESS FUNCTIONS
343 
344  string sel_model;
345 
346  for(unsigned int t = 0; t < _SelectionModels.size(); t++) {
347 
348  sel_model = _SelectionModels[t];
349 
350 
351  // FIX -----------------------------------------------------------------------
352  if(sel_model == "fix") {
353 
354  //check if the trait type in the trait table is correct
355  if(_Traits[t] != "fix")
356  return error("the selection model \"fix\" does not match with the trait type \"%s\"\n", _Traits[t].c_str());
357 
358  if(!get_parameter("selection_lethal_equivalents")->isSet()) {
359 
360  return error("\"selection_lethal_equivalents\" parameter is missing with \"fix\" selection!\n");
361 
362  } else
363  _letheq = get_parameter_value("selection_lethal_equivalents");
364 
365  if(!get_parameter("selection_base_fitness")->isSet()) {
366 
367  warning("\"selection_base_fitness\" parameter is missing under fix selection model, setting it to 1.\n");
368 
369  _base_fitness = 1.0;
370 
371  } else {
372 
373  _base_fitness = get_parameter_value("selection_base_fitness");
374 
375  }
376 
377  if(!get_parameter("selection_pedigree_F")->isSet()) {
378 
379  return error("\"selection_pedigree_F\" parameter is missing with \"fix\" selection!\n");
380 
381  } else {
382 
383  TMatrix tmp_mat;
384 
385  get_parameter("selection_pedigree_F")->getMatrix(&tmp_mat);
386 
387  if(tmp_mat.getNbCols() != 5) {
388  return error("\"selection_pedigree_F\" must be an array of size 5.\n");
389  }
390 
391  for(unsigned int i = 0; i < 5; i++) _Fpedigree[i] = tmp_mat.get(0,i);
392  }
393 
394  for(unsigned int i = 0; i < 5; i++)
395  _FitnessFixModel[i] = _base_fitness * exp( -_letheq * _Fpedigree[i] );
396 
397  //everything has been set correctly, add function pointer to container
399 
400  // DIRECT -----------------------------------------------------------------------
401  } else if(sel_model == "direct") {
402 
403  //check if the trait type in the trait table is correct
404  if(_Traits[t] != "delet" && _Traits[t] != "dmi")
405  return error("the selection model \"direct\" does not match with the trait type \"%s\"\n", _Traits[t].c_str());
406 
408 
409  // QUADRATIC --------------------------------------------------------------------
410  } else if(sel_model == "quadratic") {
411 
412  //check if the trait type in the trait table is correct
413  if(_Traits[t] != "quant")
414  return error("the selection model \"quadratic\" works only with the \"quant\" trait, your trait %i has type \"%s\"\n",
415  t, _Traits[t].c_str());
416 
417  if( get_parameter("selection_trait_dimension")->isSet() )
418  _selectTraitDimension = (unsigned int)get_parameter_value("selection_trait_dimension");
419  else
420  _selectTraitDimension = 1; // this can be changed in setSelectionMatrix() from the matrix dimensions
421 
422  //check if the number of quanti traits modeled match with the dimensionality:
423  if(_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") != _selectTraitDimension)
424  return error("The number of quantitative traits set with \"quanti_traits\" does not match with\
425  the selection trait dimensionality set with \"selection_trait_dimension\"\n");
426 
427  //check if the number of quanti traits modeled match with the dimensionality:
428  if((int)_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") > 1)
429  return error("The number of quantitative traits set with \"quanti_traits\" must be one with the \"quadratic\" selection model\n");
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( get_parameter("selection_trait_dimension")->isSet() )
445  _selectTraitDimension = (unsigned int)get_parameter_value("selection_trait_dimension");
446  else
447  _selectTraitDimension = 1; // this can be changed in setSelectionMatrix() from the matrix dimensions
448 
449  //check if the number of quanti traits modeled match with the dimensionality:
450  if(_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") < _selectTraitDimension)
451  return error("The number of quantitative traits set with \"quanti_traits\" does not match with\
452  the selection trait dimensionality set with \"selection_trait_dimension\"\n");
453 
454  // set the selection matrices for multi-stage or spatially variable cases
455  if( !setSelectionMatrix() ) return error("LCE selection::failed to set the selection matrix\n");
456 
457  if(_selectTraitDimension > 1) {
458 
459  if(_eVariance > 0)
461  else
463 
464  } else {
465 
466  if(_eVariance > 0)
468  else
470 
471  }
472 
473 
474  // LINEAR --------------------------------------------------------------------
475 
476  } else if(sel_model == "linear") {
477 
478  return error("The Linear selection model is currently off.\n");
479 //
480 // unsigned int patchNbr = _popPtr->getPatchNbr();
481 // TMatrix tmp_matx;
482 //
483 //
484 // _calc_mean_absolute_fitness = true;
485 //
486 //
487 // if( get_parameter("selection_trait_dimension")->isSet() )
488 // _selectTraitDimension = (unsigned int)get_parameter_value("selection_trait_dimension");
489 // else
490 // _selectTraitDimension = 1; // this can be changed in setSelectionMatrix() from the matrix dimensions
491 //
492 // //check if the number of quanti traits modeled match with the dimensionality:
493 // if(_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") != _selectTraitDimension)
494 // return error("The number of quantitative traits set with \"quanti_traits\" does not match with\
495 //the selection trait dimensionality set with \"selection_trait_dimension\"\n");
496 //
497 //
498 // _linear_selection_strength.reset(patchNbr, _selectTraitDimension, 1.0);
499 //
500 // if(get_parameter("selection_linear_strength")->isSet()) {
501 //
502 // get_parameter("selection_linear_strength")->getMatrix(&tmp_matx);
503 //
504 // if(tmp_matx.getNbRows() > patchNbr)
505 // return error("\"selection_linear_strength\" must have at max as many rows as the number of patches, with selection strengths for each trait as elements.\n");
506 //
507 // if(tmp_matx.ncols() > _selectTraitDimension)
508 // return error("\"selection_linear_strength\" must have at max as many columns as the number of quanti traits, with selection strengths for each patch as elements.\n");
509 //
510 // _linear_selection_strength.copy_recycle(tmp_matx);
511 //
512 // }
513 //
514 //
515 //
516 // _linear_selection_offset.reset(patchNbr, _selectTraitDimension, 1.0);
517 //
518 // if(get_parameter("selection_linear_offset")->isSet()) {
519 //
520 // get_parameter("selection_linear_offset")->getMatrix(&tmp_matx);
521 //
522 // if(tmp_matx.getNbRows() > patchNbr)
523 // return error("\"selection_linear_offset\" must have at max as many rows as the number of patches, with selection strengths for each trait as elements.\n");
524 //
525 // if(tmp_matx.ncols() > _selectTraitDimension)
526 // return error("\"selection_linear_offset\" must have at max as many columns as the number of quanti traits, with selection strengths for each patch as elements.\n");
527 //
528 // _linear_selection_offset.copy_recycle(tmp_matx);
529 // }
530 //
531 // _getRawFitness.push_back(&LCE_Selection_base::getFitnessUnivariateLinear);
532 //
533  } else return error("wrong selection model, must be either \"fix\", \"direct\", \"quadratic\", \"gaussian\", or \"linear\".\n");
534 
535  }
536 
537  if(sel_model != "fix" && !_paramSet->isSet("selection_trait"))
538  return error("trait under selection is not set, please add parameter \"selection_trait\"\n");
539 
540  return true;
541 }
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:896
bool setSelectionMatrix()
Definition: LCEselection.cc:545
double getFitnessFixedEffect(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual in the fixed selection model.
Definition: LCEselection.h:219
double getFitnessDirect(Individual *offsprg, unsigned int patch, unsigned int trait)
Returns the fitness of an individual following the direct selection model.
Definition: LCEselection.h:225
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:839
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:912
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:874
double _Fpedigree[5]
Array of pedigree values used in the fixed-selection model.
Definition: LCEselection.h:84
vector< string > _SelectionModels
The selection models associated with each trait under selection.
Definition: LCEselection.h:128
vector< string > _Traits
The list of trait types under selection.
Definition: LCEselection.h:122
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:859
vector< string > getMultiArgs()
Definition: param.cc:120
unsigned int getNbCols() const
Gives the number of columns.
Definition: tmatrix.h:215
void warning(const char *str,...)
Definition: output.cc:58

References _base_fitness, _eVariance, _FitnessFixModel, _Fpedigree, _getRawFitness, _letheq, SimComponent::_paramSet, LifeCycleEvent::_popPtr, _SelectionModels, _selectTraitDimension, _TraitIndices, _Traits, error(), TMatrix::get(), SimComponent::get_parameter(), SimComponent::get_parameter_value(), getFitnessDirect(), getFitnessFixedEffect(), getFitnessMultivariateGaussian(), getFitnessMultivariateGaussian_VE(), getFitnessUnivariateGaussian(), getFitnessUnivariateGaussian_VE(), getFitnessUnivariateQuadratic(), Param::getMatrix(), Param::getMultiArgs(), TMatrix::getNbCols(), IndFactory::getTraitPrototype(), Param::isSet(), ParamSet::isSet(), setSelectionMatrix(), and warning().

Referenced by LCE_Selection_base(), and setParameters().

◆ set_std_rate_of_change()

void LCE_Selection_base::set_std_rate_of_change ( )
1229 {
1231 
1232  //reset rate_of_change matrix in each new replicate:
1233  TMatrix tmpMat;
1234 
1235  if(!get_parameter("selection_std_rate_environmental_change")->isMatrix()) {
1236 
1237  double val = get_parameter_value("selection_std_rate_environmental_change");
1238 
1240  tmpMat.assign(val);
1241 
1242  } else {
1243  get_parameter("selection_std_rate_environmental_change")->getMatrix(&tmpMat);
1244  }
1245 
1246  unsigned int nCol = tmpMat.ncols(), nRow = tmpMat.nrows();
1247 
1249 
1250  //copy values, with pattern propagation
1251  for (int p = 0; p < _popPtr->getPatchNbr(); ++p) {
1252  for (int i = 0; i < _selectTraitDimension; ++i) {
1253  _rate_of_change_local_optima.set(p, i, tmpMat.get(p % nRow, i % nCol));
1254  }
1255  }
1256 
1257  double* SD = new double [_selectTraitDimension];
1258 
1259  for (int i = 0; i < _selectTraitDimension; ++i) {
1260  SD[i] = 0;
1261  }
1262 
1263  // check if SD is taken in a reference patch instead of the whole pop
1264  if (_std_rate_reference_patch > -1) {
1265 
1268 
1269  } else {
1270  // compute SD as the mean within-patch SD
1271 
1272  unsigned int cnt = 0;
1273 
1274  //get SD only in extant demes, to avoid nans
1275  for(unsigned int i = 0; i < _popPtr->getPatchNbr(); ++i) {
1276  if (_popPtr->getPatch(i)->size(OFFSx) != 0 ) {
1277  cnt++;
1278  addPhenotypicSD(i, SD);
1279  }
1280  }
1281 
1282  //compute mean within-patch phenotypic standard deviation:
1283  for (int i = 0; i < _selectTraitDimension; ++i) {
1284  SD[i]/= cnt;
1285 
1286  }
1287  } //end if
1288 
1289  //multiply the per-trait rates of change by SD:
1290  for (int p = 0; p < _popPtr->getPatchNbr(); ++p) {
1291  for (int i = 0; i < _selectTraitDimension; ++i){
1292  _rate_of_change_local_optima.multi(p, i, SD[i]);
1293  }
1294 
1295  }
1296  //log the rates in the simulation log file:
1297  SIMenv::MainSim->_FileServices.log("#selection_rate_environmental_change " +
1299 
1300  //compute the change of local optima for current generation:
1302 
1303  //now reset the function pointer to changeLocalOptima() for next generation:
1305 
1306  delete [] SD;
1307  }
1308 }
void log(string message)
Write to the parameter logfile.
Definition: fileservices.cc:407
void addPhenotypicSD(unsigned int deme, double *stDev)
Definition: LCEselection.cc:1312
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
207 { 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:945

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
1201 {
1202  _mean_fitness = 0;
1203 
1204  for(unsigned int i = 0; i < 5; i++) {
1205  _mean_fitness += _fitness[i];
1206  _fitness[i] /= _ind_cntr[i];
1207  _survival[i] /= _ind_cntr[i];
1208  _ind_cntr[i] /= tot_ind;
1209  }
1210 
1211  _mean_fitness /= tot_ind;
1212 }

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.

179 {
180  //selection may act on more than one trait at a time
181  _Traits.clear();
182  _Traits = get_parameter("selection_trait")->getMultiArgs();
183 
184  _TraitIndices.clear();
185 
186  for(unsigned int i = 0; i < _Traits.size(); i++) {
187 
188  if(_popPtr->getTraitIndex(_Traits[i].c_str()) == -1) {
189 
190  return error("cannot attach trait \"%s\" to life cycle event \"%s\", trait has not been initiated.\n",
191  _Traits[i].c_str(), _event_name.c_str());
192 
193  } else {
194  _TraitIndices.push_back(_popPtr->getTraitIndex(_Traits[i].c_str()));
195 
196 
197 #ifdef _DEBUG_
198  cout<<"#### attaching trait \""<<_Traits[i].c_str()<<"\" to selection LCE, index = "<<_TraitIndices[i]<<endl;
199 #endif
200 
201  }
202 
203  }
204 
205 #ifdef _DEBUG_
206  cout<<"#### setting fitness model "<<endl;
207 #endif
208 
209  if(!set_fit_model()) return false;
210 
211 #ifdef _DEBUG_
212  cout<<"#### setting selection model "<<endl;
213 #endif
214 
215  if(!set_sel_model()) return false;
216 
217 #ifdef _DEBUG_
218  cout<<"#### setting local optima "<<endl;
219 #endif
220 
221  if(!set_local_optima()) return false;
222 
223 #ifdef _DEBUG_
224  cout<<"#### setting rate of environmental change "<<endl;
225 #endif
226 
227  if(!set_param_rate_of_change()) return false;
228 
230  _scaling_factor = 1;
231  resetCounters();
232 
233  assert(_TraitIndices.size() == _getRawFitness.size());
234 
235 #ifdef _DEBUG_
236  cout<<"#### selection is acting on "<<_TraitIndices.size()<<" traits\n";
237  cout<<"#### selection is using "<< _getRawFitness.size() <<" fitness functions\n";
238 #endif
239 
240 
241  return true;
242 }
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
299 { _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
1051 {
1052 // Patch *patch;
1053 // Individual* ind;
1054 // unsigned int npatch = _popPtr->getPatchNbr();
1055 // double tot_ind;
1056 // double* _phe;
1057 //
1058 // // we have to compute the average phenotype of each quantitative trait, not of any other trait
1059 // _mean_trait_value.assign(_selectTraitDimension, 0.0);
1060 //
1061 // patch = _popPtr->getPatch(p);
1062 //
1063 // unsigned int patchAdultSize = 0, sizeFEM = patch->size(FEM, ADLTx),
1064 // sizeMAL = patch->size(MAL, ADLTx);
1065 //
1066 // for(unsigned int j = 0; j < sizeFEM; j++) {
1067 // tot_ind_fitness = 1;
1068 // _LCELinkedTraitIndex
1069 // for(unsigned int trait = 0; trait < _TraitIndices.size(); trait++) {
1070 // ind = patch->get(FEM, ADLTx, j);
1071 // _phe = (double*)ind->getTraitValue(_TraitIndices[trait]);
1072 // tot_ind_fitness *= _phe[0];
1073 // }
1074 // _mean_absolute_fitness[i] += tot_ind_fitness;
1075 // }
1076 //
1077 // for(unsigned int j = 0; j < sizeMAL; j++) {
1078 // tot_ind_fitness = 1;
1079 // for(unsigned int trait = 0; trait < _TraitIndices.size(); trait++) {
1080 // ind = patch->get(MAL, ADLTx, j);
1081 // _phe = (double*)ind->getTraitValue(_TraitIndices[trait]);
1082 // tot_ind_fitness *= _phe[0];
1083 // }
1084 // _mean_absolute_fitness[i] += tot_ind_fitness;
1085 // }
1086 //
1087 // patchAdultSize += sizeFEM + sizeMAL;
1088 // _mean_absolute_fitness[i] /= patchAdultSize;
1089 
1090 }

◆ 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
1024 {
1025  if(p != 0) return; //we compupte the mean fitness only once
1026  _scaling_factor = 1;
1027  _scaling_factor = 1.0/getMeanFitness(age);
1028 }

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
1016 {
1017  _scaling_factor = 1; //this to have the raw mean fitness below
1018  _scaling_factor = 1.0/getMeanPatchFitness(age, p);
1019 }
double getMeanPatchFitness(age_idx age, unsigned int p)
Computes the mean fitness in a given patch for a given age class.
Definition: LCEselection.cc:963

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
1042 {
1043  if(p != 0) return; //we compute the max fitness only once
1044  _scaling_factor = 1;
1045  _scaling_factor = 1.0/getMaxFitness(age);
1046 }
double getMaxFitness(age_idx age)
Computes the maximum fitness value of the whole population for a given age class.
Definition: LCEselection.cc:979

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
1034 {
1035  _scaling_factor = 1; //this to have the raw fitness when computing fitnesses below
1036  _scaling_factor = 1.0/getMaxPatchFitness(age, p);
1037 }

References _scaling_factor, and getMaxPatchFitness().

Referenced by set_fit_model().

◆ setSelectionMatrix()

bool LCE_Selection_base::setSelectionMatrix ( )
546 {
547  TMatrix tmp_mat;
548  unsigned int patchNbr = _popPtr->getPatchNbr();
549 
550  if(!get_parameter("selection_matrix")->isSet() && !get_parameter("selection_variance")->isSet())
551  return error("\"selection_matrix\" or \"selection_variance\" must be set with selection model = \"gaussian\".\n");
552 
553 
554  if(get_parameter("selection_variance")->isSet() && !get_parameter("selection_trait_dimension")->isSet())
555  return error("parameter \"selection_trait_dimension\" is missing!\n");
556 
557  //clear the selection matrix container
558  vector< TMatrix* >::iterator selIT = _selection_matrix.begin();
559 
560  for(; selIT != _selection_matrix.end(); ++selIT) if((*selIT)) delete (*selIT);
561 
562  _selection_matrix.clear();
563 
564 
565  if(get_parameter("selection_matrix")->isSet()) {
566 
567  //selection matrix provided, same selection surface in each patch
568  _paramSet->getMatrix("selection_matrix", &tmp_mat);
569 
570  if(tmp_mat.getNbCols() != tmp_mat.getNbRows())
571  return error("\"selection_matrix\" must be a square matrix!\n");
572 
573  if( get_parameter("selection_trait_dimension")->isSet() ) {
574 
575  if( (unsigned int)get_parameter_value("selection_trait_dimension") != tmp_mat.ncols() ) {
576  return error("dimensionality of \"selection_matrix\" does not match with \"selection_trait_dimension\"\n");
577  }
578 
579  } else {
580 
581  //check if the number of quanti traits modeled match with the dimensionality of the selection matrix:
582  if(_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") < tmp_mat.ncols())
583  return error("The number of quantitative traits set with \"quanti_traits\" does not match with\
584  the dimensionality of the selection matrix set with \"selection_matrix\"\n");
585  else
586  _selectTraitDimension = tmp_mat.ncols();
587 
588  }
589 
590  //we have one selection matrix per patch, copy it in the container for each patch
591  for(unsigned int i = 0; i < patchNbr; ++i)
592  _selection_matrix.push_back( new TMatrix(tmp_mat) );
593 
594  } else {
595  // selection_matrix not set in input, only selection_variance is set
596 
597  // some double checking in case selParameters wasn't call before
598  if( get_parameter("selection_trait_dimension")->isSet() ) {
599 
600  if(_selectTraitDimension != (unsigned int)get_parameter_value("selection_trait_dimension"))
601  _selectTraitDimension = (unsigned int)get_parameter_value("selection_trait_dimension");
602 
603  //check if the number of quanti traits modeled match with the dimensionality:
604  if(_popPtr->getTraitPrototype("quant")->get_parameter_value("quanti_traits") < (unsigned int)get_parameter_value("selection_trait_dimension"))
605  return error("The number of quantitative traits set with \"quanti_traits\" does not match with\
606 the selection trait dimensionality set with \"selection_trait_dimension\"\n");
607 
608  } // alternative case already checked at the onset
609 
610  TMatrix var_spatmat, corr_spatmat;
611 
612  //setting variance spatial matrix:
613  var_spatmat.reset(patchNbr, (unsigned)_selectTraitDimension);
614 
615  if(get_parameter("selection_variance")->isMatrix()) {
616 
617  _paramSet->getMatrix("selection_variance", &tmp_mat);
618 
619  if( !setSpatialMatrix("selection_variance","\"selection_trait_dimension\"", &tmp_mat, &var_spatmat,
620  (unsigned)_selectTraitDimension, patchNbr, _paramSet->isSet("selection_randomize") ) )
621  return error("LCE selection::failed to set patch-specific selection matrices from \"selection_variance\"\n");
622 
623  } else {
624 
625  var_spatmat.assign(get_parameter_value("selection_variance"));
626  }
627 
628  //setting correlation spatial matrix:
629  if(get_parameter("selection_correlation")->isSet() && _selectTraitDimension > 1)
630  {
631  corr_spatmat.reset(patchNbr, (unsigned)_selectTraitDimension*(_selectTraitDimension-1)/2);
632 
633  if(get_parameter("selection_correlation")->isMatrix()) {
634 
635  _paramSet->getMatrix("selection_correlation", &tmp_mat);
636 
637  if( !setSpatialMatrix("selection_correlation","the num of correlation coefficients", &tmp_mat, &corr_spatmat,
638  (unsigned)_selectTraitDimension*(_selectTraitDimension-1)/2, patchNbr, _paramSet->isSet("selection_randomize") ) )
639  return error("LCE selection::failed to set patch-specific selection matrices from \"selection_correlation\"\n");;
640 
641  } else {
642 
643  corr_spatmat.assign(get_parameter_value("selection_correlation"));
644  }
645  } else {
646  //a single trait or not set in input; matrix has a single row.
647  corr_spatmat.reset(patchNbr, 1, 0.0);
648  }
649 
650 
651  //set the selection matrix:
653  double covar;
654  unsigned int col;
655  for( unsigned int p = 0; p < patchNbr; p++) {
656  col = 0;
657  for( int i = 0; i < _selectTraitDimension; i++) {
658  tmp_mat.set(i, i, var_spatmat.get(p, i));
659  for( int j = i+1; j < _selectTraitDimension; j++) {
660  covar = corr_spatmat.get(p, col) * sqrt( var_spatmat.get(p, i) * var_spatmat.get(p, j) );
661  tmp_mat.set(i, j, covar);
662  tmp_mat.set(j, i, covar);
663  if(corr_spatmat.ncols() > 1) col++;
664  }
665  }
666  _selection_matrix.push_back(new TMatrix(tmp_mat));
667  }
668  }
669 
670  if(_selectTraitDimension > 1) {
671 
672  //selection on more than one trait
673 
674  //inverting the selection matrices:
675  if(_gsl_selection_matrix.size() != 0)
676  for(unsigned int i = 0; i < _gsl_selection_matrix.size(); ++i)
677  if(_gsl_selection_matrix[i] != NULL) gsl_matrix_free( _gsl_selection_matrix[i] );
678 
679  _gsl_selection_matrix.clear();
680 
681  for( unsigned int p = 0; p < patchNbr; p++) {
682  _selection_matrix[p]->inverse();
683 
684  _gsl_selection_matrix.push_back( gsl_matrix_alloc(_selectTraitDimension, _selectTraitDimension) );
685 
686  _selection_matrix[p]->get_gsl_matrix(_gsl_selection_matrix[p]);
687  }
688 
689  //allocate the vectors used by the fitness function:
690  if(_diffs != NULL) gsl_vector_free(_diffs);
691  _diffs = gsl_vector_alloc( _selectTraitDimension );
692 
693  if(_res1 != NULL) gsl_vector_free(_res1);
694  _res1 = gsl_vector_alloc( _selectTraitDimension );
695 
696  } else {
697 
698  //selection on one trait only
699  _selection_variance.clear();
700 
701  for(unsigned int p = 0; p < patchNbr; p++)
702  _selection_variance.push_back( _selection_matrix[p]->get(0, 0) );
703  }
704 
705 
706  return true;
707 }
unsigned int getNbRows() const
Gives the number of rows.
Definition: tmatrix.h:212

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(), Param::isMatrix(), ParamSet::isSet(), TMatrix::ncols(), TMatrix::reset(), TMatrix::set(), and setSpatialMatrix().

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)
1217 {
1218 
1219  _fitness[ped_class] += fitness;
1220  _ind_cntr[ped_class]++;
1221 
1222  if(survived) _survival[ped_class]++;
1223 
1224 }

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

gsl_vector* LCE_Selection_base::_diffs
private

◆ _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_offset

TMatrix LCE_Selection_base::_linear_selection_offset
private

◆ _linear_selection_strength

TMatrix LCE_Selection_base::_linear_selection_strength
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

gsl_vector * LCE_Selection_base::_res1
private

◆ _scaling_factor

◆ _selection_matrix

vector< TMatrix* > LCE_Selection_base::_selection_matrix
private

◆ _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

◆ _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

void(LCE_Selection_base::* LCE_Selection_base::_setScalingFactor) (age_idx, unsigned int)
protected

Pointer to the function used to set the fitness scaling factor when fitness is relative.

May point to LCE_Selection_base::setScalingFactorGlobal, LCE_Selection_base::setScalingFactorLocal, or LCE_Selection_base::setScalingFactorAbsolute (returns 1).

Referenced by LCE_Breed_Selection::do_breed_selection_FecFitness(), execute(), and set_fit_model().

◆ _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().

◆ _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