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

TTQuantiSH. More...

#include <ttquanti.h>

+ Inheritance diagram for TTQuantiSH:
+ Collaboration diagram for TTQuantiSH:

Public Member Functions

 TTQuantiSH (TProtoQuanti *TP)
 
virtual ~TTQuantiSH ()
 
void resetPtrs ()
 
virtual void init ()
 
virtual bool setStatRecorders (std::string &token)
 
void addQuanti (age_t AGE)
 
void addQuantiPerPatch (age_t AGE)
 
void addAvgPerPatch (age_t AGE)
 
void addGenotPerPatch (age_t AGE)
 
void addVarPerPatch (age_t AGE)
 
void addCovarPerPatch (age_t AGE)
 
void addEigenPerPatch (age_t AGE)
 
void addEigenValuesPerPatch (age_t AGE)
 
void addEigenVect1PerPatch (age_t AGE)
 
void addSkewPerPatch (age_t AGE)
 
void setDataTables (age_t AGE)
 
void setAdultStats ()
 
void setOffsprgStats ()
 
void setStats (age_t AGE)
 
double getMeanGenot (unsigned int i)
 
double getMeanPhenot (unsigned int i)
 
double getVa (unsigned int i)
 
double getVg (unsigned int i)
 
double getVb (unsigned int i)
 
double getVp (unsigned int i)
 
double getQst (unsigned int i)
 
double getCovar (unsigned int i)
 
double getMeanGenotPerPatch (unsigned int i, unsigned int p)
 
double getMeanPhenotPerPatch (unsigned int i, unsigned int p)
 
double getVaPerPatch (unsigned int i, unsigned int p)
 
double getVpPerPatch (unsigned int i, unsigned int p)
 
double getEigenValuePerPatch (unsigned int i, unsigned int p)
 
double getCovarPerPatch (unsigned int p, unsigned int i)
 
double getEigenVectorEltPerPatch (unsigned int p, unsigned int v)
 
double getSkewPerPatch (unsigned int i, unsigned int p)
 
vector< double > getSNPalleleFreqInPatch (Patch *patch, const age_idx AGE)
 
vector< double > getVaWithDominance (Patch *curPop, const age_idx AGE)
 computation of the additive genetic variance from the average excess of each allele exact under random mating only More...
 
vector< double > getVaNoDominance (Patch *curPop, const age_idx AGE)
 
- Public Member Functions inherited from TraitStatHandler< TProtoQuanti, TTQuantiSH >
 TraitStatHandler (TProtoQuanti *trait_proto)
 
virtual ~TraitStatHandler ()
 
- Public Member Functions inherited from StatHandler< SH >
 StatHandler ()
 
virtual ~StatHandler ()
 
virtual void clear ()
 Empties the _recorders list, they are destroyed in StatHandlerBase::reset(). More...
 
virtual StatRecorder< SH > * add (std::string Title, std::string Name, age_t AGE, unsigned int ARG1, unsigned int ARG2, double(SH::*getStatNoArg)(void), double(SH::*getStatOneArg)(unsigned int), double(SH::*getStatTwoArg)(unsigned int, unsigned int), void(SH::*setStat)(void))
 Adds a StatRecorder to the list, it is also added to the StatHandlerBase::_stats list. More...
 
- Public Member Functions inherited from StatHandlerBase
 StatHandlerBase ()
 
virtual ~StatHandlerBase ()
 
virtual void reset ()
 Empties the _stats list and calls clear() (defined in the derived class). More...
 
Metapopget_pop_ptr ()
 
void set_service (StatServices *srv)
 
StatServicesget_service ()
 
unsigned int getOccurrence ()
 
unsigned int getNumOccurrences ()
 
unsigned int getCurrentOccurrence ()
 
unsigned int getNbRecorders ()
 
std::list< StatRecBase * > & getStats ()
 
virtual void add (StatRecBase *rec)
 
virtual void update ()
 This function is left empty as the StatServices calls StatRecorder::setVal directly. More...
 
- Public Member Functions inherited from Handler
virtual ~Handler ()
 

Private Attributes

double * _meanP
 
double * _meanG
 
double * _Va
 
double * _Vg
 
double * _Vb
 
double * _Vp
 
double * _covar
 
double ** _pmeanP
 
double ** _pmeanG
 
double ** _pVa
 
double ** _pVp
 
double ** _pcovar
 
double ** _peigval
 
double ** _peigvect
 
bool _eVar
 
unsigned int _num_locus
 
unsigned int _num_trait
 
unsigned int _patchNbr
 
gsl_matrix * _G
 
gsl_matrix * _evec
 
gsl_vector * _eval
 
gsl_eigen_symmv_workspace * _ws
 
DataTable< double > _phenoTable
 
DataTable< double > _genoTable
 
unsigned int _table_set_gen
 
unsigned int _table_set_age
 
unsigned int _table_set_repl
 

Additional Inherited Members

- Protected Types inherited from StatHandler< SH >
typedef std::list< StatRecorder< SH > * >::iterator REC_IT
 
- Protected Attributes inherited from TraitStatHandler< TProtoQuanti, TTQuantiSH >
TProtoQuanti_SHLinkedTrait
 Pointer to a TraitProtoype object. More...
 
int _SHLinkedTraitIndex
 Index of the trait in the Individual::Traits table. More...
 
- Protected Attributes inherited from StatHandler< SH >
std::list< StatRecorder< SH > * > _recorders
 The list of stat recorders. More...
 
- Protected Attributes inherited from StatHandlerBase
Metapop_pop
 Link to the current population, set through the link to the StatService. More...
 

Detailed Description

Constructor & Destructor Documentation

◆ TTQuantiSH()

TTQuantiSH::TTQuantiSH ( TProtoQuanti TP)
inline
753  _meanP(0), _meanG(0), _Va(0), _Vg(0), _Vb(0), _Vp(0), _covar(0),
754  _pmeanP(0), _pmeanG(0), _pVa(0), _pVp(0), _pcovar(0), _peigval(0), _peigvect(0),
755  _eVar(0), _num_locus(0), _num_trait(0),_patchNbr(0), _G(0), _evec(0),_eval(0),_ws(0),
756  _table_set_gen(999999), _table_set_age(999999), _table_set_repl(999999)
757  {}
double * _Vp
Definition: ttquanti.h:736
double ** _pVa
Definition: ttquanti.h:737
double * _Vb
Definition: ttquanti.h:736
unsigned int _num_locus
Definition: ttquanti.h:740
double * _meanP
Definition: ttquanti.h:736
double ** _pVp
Definition: ttquanti.h:737
bool _eVar
Definition: ttquanti.h:738
gsl_matrix * _evec
Definition: ttquanti.h:742
unsigned int _table_set_age
Definition: ttquanti.h:747
double ** _pmeanP
Definition: ttquanti.h:737
double * _meanG
Definition: ttquanti.h:736
gsl_vector * _eval
Definition: ttquanti.h:743
double * _covar
Definition: ttquanti.h:736
gsl_matrix * _G
Definition: ttquanti.h:742
double ** _peigvect
Definition: ttquanti.h:737
unsigned int _table_set_repl
Definition: ttquanti.h:747
gsl_eigen_symmv_workspace * _ws
Definition: ttquanti.h:744
unsigned int _table_set_gen
Definition: ttquanti.h:747
unsigned int _num_trait
Definition: ttquanti.h:740
double * _Va
Definition: ttquanti.h:736
double * _Vg
Definition: ttquanti.h:736
double ** _peigval
Definition: ttquanti.h:737
unsigned int _patchNbr
Definition: ttquanti.h:740
double ** _pcovar
Definition: ttquanti.h:737
double ** _pmeanG
Definition: ttquanti.h:737

◆ ~TTQuantiSH()

virtual TTQuantiSH::~TTQuantiSH ( )
inlinevirtual
759 {resetPtrs();}
void resetPtrs()
Definition: ttquanti.cc:4545

References resetPtrs().

Member Function Documentation

◆ addAvgPerPatch()

void TTQuantiSH::addAvgPerPatch ( age_t  AGE)
4824 {
4825  if (AGE == ALL) {
4828  return;
4829  }
4830 
4831  unsigned int patchNbr = _pop->getPatchNbr();
4832 
4833  string suffix = (AGE == ADULTS ? "adlt.":"off.");
4834  string name;
4835  string patch;
4836  string t1;
4837 
4838  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
4840 
4841  add("Mean phenotype of trait 1 in patch 1", suffix + "q1.p1", AGE, 0, 0, 0, 0,
4843 
4844  for(unsigned int p = 0; p < patchNbr; p++) {
4845  for(unsigned int i = 0; i < _num_trait; i++) {
4846  if(p == 0 && i == 0) continue;
4847  name = "Mean phenotype of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
4848  t1 = "q" + tstring::int2str(i+1);
4849  patch = ".p" + tstring::int2str(p+1);
4850  add(name, suffix + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getMeanPhenotPerPatch, 0);
4851  }
4852  }
4853 
4854 }
unsigned int getPatchNbr()
Definition: metapop.h:276
Metapop * _pop
Link to the current population, set through the link to the StatService.
Definition: stathandler.h:61
virtual StatRecorder< SH > * add(std::string Title, std::string Name, age_t AGE, unsigned int ARG1, unsigned int ARG2, double(SH::*getStatNoArg)(void), double(SH::*getStatOneArg)(unsigned int), double(SH::*getStatTwoArg)(unsigned int, unsigned int), void(SH::*setStat)(void))
Adds a StatRecorder to the list, it is also added to the StatHandlerBase::_stats list.
Definition: stathandler.h:144
TTQuantiSH.
Definition: ttquanti.h:734
void addAvgPerPatch(age_t AGE)
Definition: ttquanti.cc:4823
double getMeanPhenotPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:792
void setAdultStats()
Definition: ttquanti.h:779
void setOffsprgStats()
Definition: ttquanti.h:780
static string int2str(const int i)
Writes an integer value into a string.
Definition: tstring.h:95
#define ALL
All ages age class flag.
Definition: types.h:56
#define ADULTS
Adults age class flag (breeders).
Definition: types.h:54
#define OFFSPRG
Offspring age class flag.
Definition: types.h:50

References _num_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, ALL, getMeanPhenotPerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by addQuantiPerPatch(), and setStatRecorders().

◆ addCovarPerPatch()

void TTQuantiSH::addCovarPerPatch ( age_t  AGE)
4947 {
4948  if(_num_trait < 2) {
4949  // warning("not recording traits covariance with only one \"quanti\" trait!");
4950  return;
4951  }
4952 
4953  if (AGE == ALL) {
4956  return;
4957  }
4958 
4959  string suffix = (AGE == ADULTS ? "adlt.":"off.");
4960  string patch;
4961  string cov;
4962  unsigned int patchNbr = _pop->getPatchNbr();
4963 
4964  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
4966 
4967  add("Genetic covariance of trait 1 and trait 2 in patch 1", suffix + "cov.q12.p1", AGE, 0, 0, 0, 0,
4968  &TTQuantiSH::getCovarPerPatch, setter);
4969 
4970  unsigned int c;
4971  for(unsigned int p = 0; p < patchNbr; p++) {
4972  patch = ".p" + tstring::int2str(p+1);
4973  c = 0;
4974  for(unsigned int t = 0; t < _num_trait; ++t) {
4975  for(unsigned int v = t + 1; v < _num_trait; ++v){
4976  if(p==0 && t==0 && v==1) {c++; continue;} //this on is already recorder above
4977  cov = tstring::int2str((t+1)*10+v+1);
4978  add("", suffix + "cov.q" + cov + patch, AGE, c++, p, 0, 0, &TTQuantiSH::getCovarPerPatch, 0);
4979  }
4980  }
4981  }
4982 
4983 }
void addCovarPerPatch(age_t AGE)
Definition: ttquanti.cc:4946
double getCovarPerPatch(unsigned int p, unsigned int i)
Definition: ttquanti.h:796

References _num_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, ALL, getCovarPerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by addQuantiPerPatch(), and setStatRecorders().

◆ addEigenPerPatch()

void TTQuantiSH::addEigenPerPatch ( age_t  AGE)
4988 {
4989  if(_num_trait < 2) {
4990  warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
4991  return;
4992  }
4993 
4994  if (AGE == ALL) {
4997  return;
4998  }
4999 
5000  unsigned int patchNbr = _pop->getPatchNbr();
5001  string suffix = (AGE == ADULTS ? "adlt.":"off.");
5002  string patch;
5003  unsigned int pv =0;
5004 
5005  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
5007 
5008 
5009  add("First G-matrix eigenvalue in patch 1", suffix + "qeval1.p1", AGE, 0, 0, 0, 0,
5011 
5012  for(unsigned int p = 0; p < patchNbr; ++p) {
5013  patch = ".p" + tstring::int2str(p+1);
5014  for(unsigned int t = 0; t < _num_trait; ++t) {
5015  if(p==0 && t==0) continue;
5016  add("", suffix + "qeval" + tstring::int2str(t+1) + patch, AGE, t, p, 0, 0, &TTQuantiSH::getEigenValuePerPatch,0);
5017  }
5018  }
5019  for(unsigned int p = 0; p < patchNbr; ++p) {
5020  patch = ".p" + tstring::int2str(p+1);
5021  pv = 0;
5022  for(unsigned int t = 0; t < _num_trait; ++t)
5023  for(unsigned int v = 0; v < _num_trait; ++v)
5024  add("", suffix + "qevect" + tstring::int2str((t+1)*10+v+1) + patch, AGE, p, pv++, 0, 0, &TTQuantiSH::getEigenVectorEltPerPatch,0);
5025  }
5026 
5027 }
double getEigenVectorEltPerPatch(unsigned int p, unsigned int v)
Definition: ttquanti.h:797
double getEigenValuePerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:795
void addEigenPerPatch(age_t AGE)
Definition: ttquanti.cc:4987
void warning(const char *str,...)
Definition: output.cc:58

References _num_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, ALL, getEigenValuePerPatch(), getEigenVectorEltPerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, setAdultStats(), setOffsprgStats(), and warning().

Referenced by setStatRecorders().

◆ addEigenValuesPerPatch()

void TTQuantiSH::addEigenValuesPerPatch ( age_t  AGE)
5032 {
5033  if(_num_trait < 2) {
5034  warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
5035  return;
5036  }
5037 
5038  if (AGE == ALL) {
5041  return;
5042  }
5043 
5044  unsigned int patchNbr = _pop->getPatchNbr();
5045  string suffix = (AGE == ADULTS ? "adlt.":"off.");
5046  string patch;
5047 
5048  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
5050 
5051  add("First G-matrix eigenvalue in patch 1", suffix + "qeval1.p1", AGE, 0, 0, 0, 0,
5053 
5054  for(unsigned int p = 0; p < patchNbr; ++p) {
5055  patch = ".p" + tstring::int2str(p+1);
5056  for(unsigned int t = 0; t < _num_trait; ++t) {
5057  if(p==0 && t==0) continue;
5058  add("", suffix + "qeval" + tstring::int2str(t+1) + patch, AGE, t, p, 0, 0, &TTQuantiSH::getEigenValuePerPatch,0);
5059  }
5060  }
5061 }
void addEigenValuesPerPatch(age_t AGE)
Definition: ttquanti.cc:5031

References _num_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, ALL, getEigenValuePerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, setAdultStats(), setOffsprgStats(), and warning().

Referenced by setStatRecorders().

◆ addEigenVect1PerPatch()

void TTQuantiSH::addEigenVect1PerPatch ( age_t  AGE)
5066 {
5067  if(_num_trait < 2) {
5068  warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
5069  return;
5070  }
5071 
5072  if (AGE == ALL) {
5075  return;
5076  }
5077 
5078  unsigned int patchNbr = _pop->getPatchNbr();
5079  string suffix = (AGE == ADULTS ? "adlt.":"off.");
5080  string patch;
5081  unsigned int pv =0;
5082 
5083  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
5086 
5087 
5088  add("First G-matrix eigenvector in patch 1", suffix + "qevect11.p1", AGE, 0, 0, 0, 0,
5090 
5091  for(unsigned int p = 0; p < patchNbr; ++p) {
5092  patch = ".p" + tstring::int2str(p+1);
5093  pv = 0;
5094  // for(unsigned int t = 0; t < _num_trait; ++t)
5095  for(unsigned int v = 0; v < _num_trait; ++v){
5096  if(p==0 && v==0) {pv++; continue;}
5097  add("", suffix + "qevect1" + tstring::int2str(v+1) + patch, AGE, p, pv++, 0, 0, &TTQuantiSH::getEigenVectorEltPerPatch,0);
5098  }
5099  }
5100 }
void addEigenVect1PerPatch(age_t AGE)
Definition: ttquanti.cc:5065

References _num_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, ALL, getEigenVectorEltPerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, setAdultStats(), setOffsprgStats(), and warning().

Referenced by setStatRecorders().

◆ addGenotPerPatch()

void TTQuantiSH::addGenotPerPatch ( age_t  AGE)
4859 {
4860  if (AGE == ALL) {
4863  return;
4864  }
4865 
4866  unsigned int patchNbr = _pop->getPatchNbr();
4867 
4868  string suffix = (AGE == ADULTS ? "adlt.":"off.");
4869  string name;
4870  string patch;
4871  string t1;
4872 
4873  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
4875 
4876  add("Mean genotype of trait 1 in patch 1", suffix + "g1.p1", AGE, 0, 0, 0, 0,
4878 
4879  for(unsigned int p = 0; p < patchNbr; p++) {
4880  for(unsigned int i = 0; i < _num_trait; i++) {
4881  if(p == 0 && i == 0) continue;
4882  name = "Mean genotype of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
4883  t1 = "g" + tstring::int2str(i+1);
4884  patch = ".p" + tstring::int2str(p+1);
4885  add(name, suffix + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getMeanGenotPerPatch, 0);
4886  }
4887  }
4888 
4889 }
double getMeanGenotPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:791
void addGenotPerPatch(age_t AGE)
Definition: ttquanti.cc:4858

References _num_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, ALL, getMeanGenotPerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by setStatRecorders().

◆ addQuanti()

void TTQuantiSH::addQuanti ( age_t  AGE)
4739 {
4740  if (AGE == ALL) {
4741  addQuanti(ADULTS);
4742  addQuanti(OFFSPRG);
4743  return;
4744  }
4745 
4746  string suffix = (AGE == ADULTS ? "adlt.":"off.");
4747  string name = suffix + "q";
4748  string t1, t2;
4749 
4750  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
4752 
4753  add("",suffix + "q1",AGE,0,0,0,&TTQuantiSH::getMeanPhenot,0,setter);
4754 
4755  for(unsigned int i = 1; i < _num_trait; i++) {
4756  t1 = tstring::int2str(i+1);
4757  add("", name + t1,AGE,i,0,0,&TTQuantiSH::getMeanPhenot,0,0);
4758  }
4759 
4760  //Va, no dominance, no epistasis
4762  for(unsigned int i = 0; i < _num_trait; i++) {
4763  t1 = tstring::int2str(i+1);
4764  add("", name + t1 +".Va",AGE,i,0,0,&TTQuantiSH::getVa,0,0);
4765  }
4766  } else {
4767  //if dominance or epistasis, then what we record is Vg, not Va
4768  for(unsigned int i = 0; i < _num_trait; i++) {
4769  t1 = tstring::int2str(i+1);
4770  add("", name + t1 +".Vg",AGE,i,0,0,&TTQuantiSH::getVa,0,0);
4771  }
4772  }
4773 
4774  //Vp
4775  if(_eVar) {
4776  for(unsigned int i = 0; i < _num_trait; i++) {
4777  t1 = tstring::int2str(i+1);
4778  add("", name + t1 +".Vp",AGE,i,0,0,&TTQuantiSH::getVp,0,0);
4779  }
4780  }
4781  //Vb
4782  if(_pop->getPatchNbr() > 1){
4783  for(unsigned int i = 0; i < _num_trait; i++) {
4784  t1 = tstring::int2str(i+1);
4785  add("", name + t1 +".Vb",AGE,i,0,0,&TTQuantiSH::getVb,0,0);
4786  }
4787  //Qst
4788  for(unsigned int i = 0; i < _num_trait; i++) {
4789  t1 = tstring::int2str(i+1);
4790  add("", name + t1 +".Qst",AGE,i,0,0,&TTQuantiSH::getQst,0,0);
4791  }
4792  }
4793 
4794  if (_num_trait > 1) {
4795  unsigned int c = 0;
4796  for(unsigned int t = 0; t < _num_trait; ++t) {
4797  for(unsigned int v = t + 1; v < _num_trait; ++v) {
4798  t1 = tstring::int2str((t+1)*10+(v+1));
4799  add("", name + t1 +".cov", AGE, c++, 0,0,&TTQuantiSH::getCovar,0,0);
4800  }
4801  }
4802  }
4803 }
bool do_epistasis()
Definition: ttquanti.h:567
unsigned int get_dominance_model()
Definition: ttquanti.h:478
void addQuanti(age_t AGE)
Definition: ttquanti.cc:4738
double getVa(unsigned int i)
Definition: ttquanti.h:784
double getVb(unsigned int i)
Definition: ttquanti.h:786
double getMeanPhenot(unsigned int i)
Definition: ttquanti.h:783
double getQst(unsigned int i)
Definition: ttquanti.h:788
double getVp(unsigned int i)
Definition: ttquanti.h:787
double getCovar(unsigned int i)
Definition: ttquanti.h:789
TProtoQuanti * _SHLinkedTrait
Pointer to a TraitProtoype object.
Definition: stathandler.h:171

References _eVar, _num_trait, StatHandlerBase::_pop, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, StatHandler< SH >::add(), ADULTS, ALL, TProtoQuanti::do_epistasis(), TProtoQuanti::get_dominance_model(), getCovar(), getMeanPhenot(), Metapop::getPatchNbr(), getQst(), getVa(), getVb(), getVp(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by setStatRecorders().

◆ addQuantiPerPatch()

void TTQuantiSH::addQuantiPerPatch ( age_t  AGE)
4808 {
4809 
4810  if (AGE == ALL) {
4813  return;
4814  }
4815 
4816  addAvgPerPatch(AGE);
4817  addVarPerPatch(AGE);
4818  addCovarPerPatch(AGE);
4819 }
void addVarPerPatch(age_t AGE)
Definition: ttquanti.cc:4893
void addQuantiPerPatch(age_t AGE)
Definition: ttquanti.cc:4807

References addAvgPerPatch(), addCovarPerPatch(), addVarPerPatch(), ADULTS, ALL, and OFFSPRG.

Referenced by setStatRecorders().

◆ addSkewPerPatch()

void TTQuantiSH::addSkewPerPatch ( age_t  AGE)
5105 {
5106  if (AGE == ALL) {
5109  return;
5110  }
5111 
5112  unsigned int patchNbr = _pop->getPatchNbr();
5113 
5114  string suffix = (AGE == ADULTS ? "adlt.":"off.");
5115  string name;
5116  string patch;
5117  string t1;
5118 
5119  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
5122 
5123  add("Genetic skew of trait 1 in patch 1", suffix + "Sk.q1.p1", AGE, 0, 0, 0, 0,
5124  &TTQuantiSH::getSkewPerPatch, setter);
5125 
5126  for(unsigned int p = 0; p < patchNbr; p++) {
5127  for(unsigned int i = 0; i < _num_trait; i++) {
5128  if(p == 0 && i == 0) continue;
5129  name = "Genetic skew of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
5130  t1 = "q" + tstring::int2str(i+1);
5131  patch = ".p" + tstring::int2str(p+1);
5132  add(name, suffix + "Sk." + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getSkewPerPatch, 0);
5133  }
5134  }
5135 
5136 }
double getSkewPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.cc:5140
void addSkewPerPatch(age_t AGE)
Definition: ttquanti.cc:5104

References _num_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, ALL, Metapop::getPatchNbr(), getSkewPerPatch(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by setStatRecorders().

◆ addVarPerPatch()

void TTQuantiSH::addVarPerPatch ( age_t  AGE)
4894 {
4895  if (AGE == ALL) {
4898  return;
4899  }
4900 
4901  unsigned int patchNbr = _pop->getPatchNbr();
4902 
4903  string suffix = (AGE == ADULTS ? "adlt.":"off.");
4904  string name;
4905  string patch;
4906  string t1;
4907  string Vtype;
4908 
4909 
4911  Vtype = "Vg.";
4912  else
4913  Vtype = "Va.";
4914 
4915  void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
4918 
4919  add("Genetic variance of trait 1 in patch 1", suffix + Vtype + "q1.p1", AGE, 0, 0, 0, 0,
4920  &TTQuantiSH::getVaPerPatch, setter);
4921 
4922  for(unsigned int p = 0; p < patchNbr; p++) {
4923  for(unsigned int i = 0; i < _num_trait; i++) {
4924  if(p == 0 && i == 0) continue;
4925  name = "Genetic variance of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
4926  t1 = "q" + tstring::int2str(i+1);
4927  patch = ".p" + tstring::int2str(p+1);
4928  add(name, suffix + Vtype + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getVaPerPatch, 0);
4929  }
4930  }
4931 
4932  if(_eVar) {
4933  for(unsigned int p = 0; p < patchNbr; p++) {
4934  for(unsigned int i = 0; i < _num_trait; i++) {
4935  name = "Phenotypic variance of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
4936  t1 = "q" + tstring::int2str(i+1);
4937  patch = ".p" + tstring::int2str(p+1);
4938  add(name, suffix + "Vp." + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getVpPerPatch, 0);
4939  }
4940  }
4941  }
4942 }
double getVpPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:794
double getVaPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:793

References _eVar, _num_trait, StatHandlerBase::_pop, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, StatHandler< SH >::add(), ADULTS, ALL, TProtoQuanti::do_epistasis(), TProtoQuanti::get_dominance_model(), Metapop::getPatchNbr(), getVaPerPatch(), getVpPerPatch(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by addQuantiPerPatch(), and setStatRecorders().

◆ getCovar()

double TTQuantiSH::getCovar ( unsigned int  i)
inline
789 {return _covar[i];}

References _covar.

Referenced by addQuanti().

◆ getCovarPerPatch()

double TTQuantiSH::getCovarPerPatch ( unsigned int  p,
unsigned int  i 
)
inline
796 {return _pcovar[p][i];}

References _pcovar.

Referenced by addCovarPerPatch().

◆ getEigenValuePerPatch()

double TTQuantiSH::getEigenValuePerPatch ( unsigned int  i,
unsigned int  p 
)
inline
795 {return _peigval[i][p];}

References _peigval.

Referenced by addEigenPerPatch(), and addEigenValuesPerPatch().

◆ getEigenVectorEltPerPatch()

double TTQuantiSH::getEigenVectorEltPerPatch ( unsigned int  p,
unsigned int  v 
)
inline
797 {return _peigvect[p][v];}

References _peigvect.

Referenced by addEigenPerPatch(), and addEigenVect1PerPatch().

◆ getMeanGenot()

double TTQuantiSH::getMeanGenot ( unsigned int  i)
inline
782 {return _meanG[i];}

References _meanG.

◆ getMeanGenotPerPatch()

double TTQuantiSH::getMeanGenotPerPatch ( unsigned int  i,
unsigned int  p 
)
inline
791 {return _pmeanG[i][p];}

References _pmeanG.

Referenced by addGenotPerPatch().

◆ getMeanPhenot()

double TTQuantiSH::getMeanPhenot ( unsigned int  i)
inline
783 {return _meanP[i];}

References _meanP.

Referenced by addQuanti().

◆ getMeanPhenotPerPatch()

double TTQuantiSH::getMeanPhenotPerPatch ( unsigned int  i,
unsigned int  p 
)
inline
792 {return _pmeanP[i][p];}

References _pmeanP.

Referenced by addAvgPerPatch().

◆ getQst()

double TTQuantiSH::getQst ( unsigned int  i)
inline
788 {return _Vb[i]/(_Vb[i]+2*_Va[i]);}

References _Va, and _Vb.

Referenced by addQuanti().

◆ getSkewPerPatch()

double TTQuantiSH::getSkewPerPatch ( unsigned int  i,
unsigned int  p 
)
5141 {
5142  double skew = 0;
5143 
5144  double *phenot = _phenoTable.getClassWithinGroup(i, p);
5145  unsigned int patch_size = _phenoTable.size(i, p);
5146 
5147  for(unsigned int k = 0; k < patch_size; ++k)
5148  skew += pow( phenot[k] - _pmeanP[i][p], 3); //the mean has been set by setStats()
5149 
5150  return skew / patch_size;
5151 }
unsigned int size(unsigned int i, unsigned int j)
Definition: datatable.h:260
T * getClassWithinGroup(unsigned int group, unsigned int Class)
Accessor to a class array whithin a group.
Definition: datatable.h:225
DataTable< double > _phenoTable
Definition: ttquanti.h:746

References _phenoTable, _pmeanP, DataTable< T >::getClassWithinGroup(), and DataTable< T >::size().

Referenced by addSkewPerPatch().

◆ getSNPalleleFreqInPatch()

vector< double > TTQuantiSH::getSNPalleleFreqInPatch ( Patch patch,
const age_idx  AGE 
)
5348 {
5349  //ONLY IMPLEMENTED FOR DIALLELIC LOCI!
5350  if(_SHLinkedTrait->get_allele_model() > 1)
5351  fatal("TTQuantiSH::getSNPalleleFreqInPatch::SNP allele frequencies can only be computed for di-allelic loci, sorry!\n");
5352 
5353  Individual* ind;
5354  TTQuanti* Qtrait;
5355 
5356  if(!_SHLinkedTrait) fatal("TTQuantiSH::getSNPalleleFreqInPatch:: Pointer to trait prototype is not set!\n");
5357 
5358  unsigned int seqLength = _SHLinkedTrait->get_seq_length();
5359 
5360  //create a null vector of the total sequence size, including alleles affecting all traits (max 2 for di-allelic loci)
5361  vector<double> snp_tab(seqLength, 0.0);
5362 
5363  for(unsigned int s = 0; s < 2; ++s){
5364 
5365  for(unsigned int j = 0, size = patch->size(sex_t(s), AGE); j < size; ++j) {
5366 
5367  ind = patch->get(sex_t(s), AGE, j);
5368 
5369  Qtrait = dynamic_cast<TTQuanti*>(ind->getTrait(_SHLinkedTraitIndex));
5370 
5371  //compute allele freq at all positions, we count allele 0 (freq 'p' is for trait-increasing allele in diallelic model)
5372  for(unsigned int i = 0; i < seqLength; ++i) {
5373  snp_tab[i] += !Qtrait->get_allele_bit(i, 0) + !Qtrait->get_allele_bit(i, 1);
5374  }
5375 
5376  }
5377 
5378  }
5379 
5380  // pre-compute 1/2N :
5381  double inv_size = 1.0/(2.0*patch->size(AGE));
5382 
5383  std::transform(snp_tab.cbegin(), snp_tab.cend(), snp_tab.begin(), [inv_size](double i){return i*inv_size;});
5384 
5385  return snp_tab;
5386 }
This class contains traits along with other individual information (sex, pedigree,...
Definition: individual.h:49
TTrait * getTrait(IDX T)
Trait accessor.
Definition: individual.h:277
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
unsigned int get_allele_model()
Definition: ttquanti.h:434
unsigned int get_seq_length()
Definition: ttquanti.h:427
TTQuanti.
Definition: ttquanti.h:61
virtual bool get_allele_bit(unsigned int position, unsigned int allele) const =0
int _SHLinkedTraitIndex
Index of the trait in the Individual::Traits table.
Definition: stathandler.h:173
void fatal(const char *str,...)
Definition: output.cc:100
sex_t
Sex types, males are always 0 and females 1!!
Definition: types.h:36

References TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, fatal(), Patch::get(), TTQuanti::get_allele_bit(), TProtoQuanti::get_allele_model(), TProtoQuanti::get_seq_length(), Individual::getTrait(), and Patch::size().

Referenced by getVaNoDominance(), and getVaWithDominance().

◆ getVa()

double TTQuantiSH::getVa ( unsigned int  i)
inline
784 {return _Va[i];}

References _Va.

Referenced by addQuanti().

◆ getVaNoDominance()

vector< double > TTQuantiSH::getVaNoDominance ( Patch curPop,
const age_idx  AGE 
)

CYCLE THROUGH ALL TRAITS

5391 {
5392  // function called by LCE_QuantiModifier when setting fixed h^2
5393  unsigned int sizeF = patch->size(FEM, AGE),
5394  sizeM = patch->size(MAL, AGE);
5395  unsigned int size = sizeF + sizeM;
5396 
5397  if(!size) {return vector<double>(_num_trait, -1.0);}
5398 
5399  TTQuanti* trait;
5400 
5401  double G, meanG;
5402 
5403  vector< double > varA(_num_trait, 0.0);
5404 
5405  double *genot = new double[size];
5406 
5408  for(unsigned int TRAIT = 0; TRAIT < _num_trait; ++TRAIT)
5409  {
5410 
5411  meanG = 0;
5412 
5413  // females
5414 
5415  Individual* curInd;
5416  unsigned int gpos = 0;
5417 
5418  for(unsigned int i = 0; i < sizeF && gpos < size; ++i )
5419  {
5420  curInd = patch->get(FEM, AGE, i);
5421 
5422  trait = dynamic_cast<TTQuanti*>(curInd->getTrait(_SHLinkedTraitIndex));
5423 
5424  G = trait->get_additive_genotype(TRAIT); // genotype value
5425 
5426  genot[gpos++] = G;
5427 
5428  meanG += G;
5429 
5430  } // for each female
5431 
5432  // males
5433  for(unsigned int i = 0; i < sizeM && gpos < size; ++i )
5434  {
5435  curInd = patch->get(MAL, AGE, i);
5436 
5437  trait = dynamic_cast<TTQuanti*>(curInd->getTrait(_SHLinkedTraitIndex));
5438 
5439  G = trait->get_additive_genotype(TRAIT); // genotype value
5440 
5441  genot[gpos++] = G;
5442 
5443  meanG += G;
5444 
5445  } // for each male
5446 
5447  assert(gpos == size);
5448 
5449  meanG /= size;
5450 
5451  varA[TRAIT] = my_variance_with_fixed_mean(genot, size, meanG); // not a sample variance
5452 
5453  // record the genotypic variance in the stater
5454  _Vg[TRAIT] = varA[TRAIT];
5455 
5456 #ifdef _DEBUG_
5457  message("TTQuantiSH::getVaNoDominance:: Va = %d\n",varA[TRAIT]);
5458 
5459  const TMatrix& allele_value = _SHLinkedTrait->get_diallele_values(); // nlocus x 2
5460  double Vacheck = 0, a, d,avG;
5461  //compute the allele frequencies:
5462  vector< double > freqs = getSNPalleleFreqInPatch(patch, AGE); //it has num locus x num traits elements, in a 1D array
5463 
5464  for(unsigned int l = 0, pos; l < _num_locus; ++l)
5465  {
5466  pos = l * _num_trait + TRAIT;
5467  a = abs(allele_value.get(l,0) - allele_value.get(l,1))/2;
5468  d = allele_value.get(l,0)+allele_value.get(l,1);
5469  avG = a + d*(1-2*freqs[pos]);
5470  Vacheck += 2*freqs[pos]*(1 - freqs[pos])*avG*avG;
5471  }
5472 
5473  message("TTQuantiSH::getVaNoDominance:: Va check = %d\n",Vacheck);
5474 #endif
5475 
5476  }
5477 
5478  delete [] genot;
5479 
5480  return varA;
5481 }
A class to handle matrix in params, coerces matrix into a vector of same total size.
Definition: tmatrix.h:50
double get(unsigned int i, unsigned int j) const
Accessor to element at row i and column j.
Definition: tmatrix.h:193
const TMatrix & get_diallele_values()
Definition: ttquanti.h:435
vector< double > getSNPalleleFreqInPatch(Patch *patch, const age_idx AGE)
Definition: ttquanti.cc:5347
virtual double get_additive_genotype(const unsigned int trait) const =0
void message(const char *message,...)
Definition: output.cc:40
@ FEM
Definition: types.h:37
@ MAL
Definition: types.h:37
double my_variance_with_fixed_mean(double *data, unsigned int size, double mean)
Definition: utils.cc:63

References _num_locus, _num_trait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, _Vg, FEM, Patch::get(), TMatrix::get(), TTQuanti::get_additive_genotype(), TProtoQuanti::get_diallele_values(), getSNPalleleFreqInPatch(), Individual::getTrait(), MAL, message(), my_variance_with_fixed_mean(), and Patch::size().

Referenced by LCE_QuantiModifier::setParameters(), and LCE_Breed_Quanti::setParameters().

◆ getVaPerPatch()

double TTQuantiSH::getVaPerPatch ( unsigned int  i,
unsigned int  p 
)
inline
793 {return _pVa[i][p];}

References _pVa.

Referenced by addVarPerPatch().

◆ getVaWithDominance()

vector< double > TTQuantiSH::getVaWithDominance ( Patch curPop,
const age_idx  AGE 
)

computation of the additive genetic variance from the average excess of each allele exact under random mating only

CYCLE THROUGH ALL TRAITS

5490 {
5491  // function called by LCE_QuantiModifier when setting fixed h^2
5492  if(_SHLinkedTrait->get_allele_model() > 2)
5493  fatal("Quanti Trait StatHandler::Va with dominance can only be estimated for di-allelic loci, sorry!\n");
5494 
5495  unsigned int size = patch->size(AGE);
5496 
5497 #ifdef _DEBUG_
5498  message("\nTTQuantiSH::getVaWithDominance in patch %i; size = %i\n", patch->getID(), size);
5499 #endif
5500 
5501  if(!size)
5502  return vector<double>(_num_trait, -1.0);
5503 
5504 
5505  //compute the allele frequencies:
5506  vector< double > freqs = getSNPalleleFreqInPatch(patch, AGE); //it has num locus x num traits elements, in a 1D array
5507 
5508  unsigned int numLocus, pos, a1, a2;
5509  TTQuanti* trait;
5510  const TMatrix& allele_value = _SHLinkedTrait->get_diallele_values(); // nlocus x 2
5511 
5512  double G, meanG;
5513  vector< double > genotypeValues[3];
5514  vector< double > genotypeFreq[3];
5515 
5516  double *genot = new double[size];
5517 
5518  vector< double > varA(_num_trait, 0.0);
5519 
5521  for(unsigned int TRAIT = 0; TRAIT < _num_trait; ++TRAIT)
5522  {
5523 
5524  numLocus = _SHLinkedTrait->get_num_locus(TRAIT);
5525 
5526  meanG = 0;
5527 
5528  genotypeValues[0].assign(numLocus, 0.0); //AA
5529  genotypeValues[1].assign(numLocus, 0.0); //Aa/aA
5530  genotypeValues[2].assign(numLocus, 0.0); //aa
5531 
5532  genotypeFreq[0].assign(numLocus, 0.0); //AA
5533  genotypeFreq[1].assign(numLocus, 0.0); //Aa/aA
5534  genotypeFreq[2].assign(numLocus, 0.0); //aa
5535 
5536  for(unsigned int l = 0, locID; l < numLocus; ++l)
5537  {
5538  locID = _SHLinkedTrait->get_locus_ID(l, TRAIT);
5539  genotypeValues[0][l] = _SHLinkedTrait->get_genotype_with_dominance( allele_value.get(locID,0), allele_value.get(locID,0), locID, TRAIT);
5540  genotypeValues[1][l] = _SHLinkedTrait->get_genotype_with_dominance( allele_value.get(locID,0), allele_value.get(locID,1), locID, TRAIT);
5541  genotypeValues[2][l] = _SHLinkedTrait->get_genotype_with_dominance( allele_value.get(locID,1), allele_value.get(locID,1), locID, TRAIT);
5542  }
5543 
5544  // cycle trhough all individuals:
5545 
5546  Individual* curInd;
5547  unsigned int gpos = 0, fpos = 0;
5548 
5549  for(unsigned int s = 0; s < 2; ++s) {
5550  for(unsigned int i = 0; i < patch->size(sex_t(s), AGE) && gpos < size; ++i )
5551  {
5552  curInd = patch->get(sex_t(s), AGE, i);
5553 
5554  trait = dynamic_cast<TTQuanti*>(curInd->getTrait(_SHLinkedTraitIndex));
5555 
5556  G = trait->get_full_genotype(TRAIT); // genotype value
5557 
5558  meanG += G;
5559 
5560  genot[gpos++] = G;
5561 
5562  for(unsigned int l = 0; l < numLocus; ++l)
5563  {
5564  pos = _SHLinkedTrait->get_locus_seq_pos(l, TRAIT); //l * _num_trait + TRAIT;
5565 
5566  a1 = trait->get_allele_bit(pos, 0);
5567  a2 = trait->get_allele_bit(pos, 1);
5568 
5569  fpos = a1 + a2; //which genotype: 0 = AA, 1 = Aa/aA, 2 = aa
5570 
5571  genotypeFreq[fpos][l]++;
5572 
5573  // if( genotypeValues[fpos][l] != _SHLinkedTrait->get_genotype_with_dominance( genes[0][pos], genes[1][pos], locID, TRAIT))
5574  // {
5575  // cout << "assert failure: genot table value: "<<genotypeValues[fpos][l]<<" at "<<fpos
5576  // <<"; genot value: "<<_SHLinkedTrait->get_genotype_with_dominance( genes[0][pos], genes[1][pos], locID, TRAIT)
5577  // <<" ("<< genes[0][pos]<<", "<< genes[1][pos] <<", "<<_SHLinkedTrait->get_dominance(locID, TRAIT)<<")\n";
5578  // }
5579  }
5580 
5581 
5582  } // for each female
5583  }
5584 
5585  assert(gpos == size);
5586 
5587  meanG /= size;
5588 
5589  // record the genotypic variance in the stater, will be needed to comput Ve
5590  _Vg[TRAIT] = my_variance_with_fixed_mean(genot, size, meanG);
5591 
5592  for(unsigned int l = 0; l < numLocus; ++l)
5593  {
5594  genotypeFreq[0][l] /= size;
5595  genotypeFreq[1][l] /= size;
5596  genotypeFreq[2][l] /= size;
5597  }
5598 
5599 
5600  // compute the average excess (alphaStar) which is identical to the additive effects under random mating
5601 
5602  // !!! THIS DOESN'T WORK AT ALL !!!
5603 
5604 // vector< double > alphaStar[2];
5605 // alphaStar[0].assign(_num_locus, 0.0);
5606 // alphaStar[1].assign(_num_locus, 0.0);
5607 //
5608 // double mean_alpha = 0;
5609 //
5610 // for(unsigned int l = 0; l < _num_locus; ++l)
5611 // {
5612 //
5613 // pos = l * _num_trait + TRAIT;
5614 //
5615 // // average excess of each allele; the frequency table contains the frequency of the first allele only
5616 // // here, we condition on the genotype being present in the population, in case HW equilibrium is not reached
5617 // // eq 4.13b Lynch & Walsh 98
5618 // // allele frequency freq[pos] is that of the first allele at [0]
5619 // alphaStar[0][l] = genotypeValues[1][l] * (1.0 - freqs[pos]) + genotypeValues[0][l] * freqs[pos] - meanG; // alpha*_1
5620 //
5621 // alphaStar[1][l] = genotypeValues[1][l] * freqs[pos] + genotypeValues[2][l] * (1.0-freqs[pos]) - meanG; // alpha*_2
5622 //
5626 //
5627 // mean_alpha += (alphaStar[0][l] + alphaStar[1][l])/2;
5628 //
5629 // varA[TRAIT] += alphaStar[0][l]*alphaStar[0][l]*freqs[pos] + alphaStar[1][l]*alphaStar[1][l]*(1-freqs[pos]); // SUM(p_i*alpha_i^2)
5630 //
5631 // }
5632 // varA[TRAIT] *= 2; // for diploids
5633 //
5634 // cout <<"TTQuantiSH::getVaWithDominance:: Va = "<<varA[TRAIT]<< " (mean alpha = "<<mean_alpha/_num_locus<<")"<<endl;
5635 
5636 
5637  // check form general formula Va = 2pq[a + d(q-p)]^2 !! q here is 1-freqs[pos], freq of smaller allele (reducing value) which has index [1]
5638  double Vacheck = 0, a, d,avG, VDcheck = 0, vd, p, q, pq;
5639  for(unsigned int l = 0; l < numLocus; ++l)
5640  {
5641  // the freq vector has same arrangement of loci as the quanti sequence in the genome
5642  pos = _SHLinkedTrait->get_locus_seq_pos(l, TRAIT);
5643  p = freqs[pos];
5644  q = 1 - freqs[pos];
5645  pq= p*q;
5646  a = abs(genotypeValues[0][l] - genotypeValues[2][l])/2.0;
5647  d = genotypeValues[1][l];
5648  avG = a + d*(1-2*p); //'freqs' stores 'p' not 'q', at same location in the array as loci in the sequences
5649  Vacheck += 2*pq*avG*avG;
5650  vd = 2*pq*d;
5651  VDcheck += vd*vd;
5652  }
5653 
5655  varA [TRAIT] = _Vg[TRAIT];
5656  else
5657  varA[TRAIT] = Vacheck;
5658 
5659 #ifdef _DEBUG_
5660  message("TTQuantiSH::getVaWithDominance:: mean genotype = %d\n", meanG);
5661  message("TTQuantiSH::getVaWithDominance:: var genotype Vg = %d\n",_Vg[TRAIT]);
5662  message("TTQuantiSH::getVaWithDominance:: estimates: Va = %d, Vd = %d, Vg = %d\n", Vacheck, VDcheck, Vacheck+VDcheck);
5663 #endif
5664 
5665  } //for each trait
5666 
5667  delete [] genot;
5668 
5669  return varA;
5670 }
unsigned int get_locus_ID(unsigned int locus, unsigned int trait)
Definition: ttquanti.h:455
unsigned int get_locus_seq_pos(unsigned int loc, unsigned int trait)
Definition: ttquanti.h:453
unsigned int get_num_locus()
Definition: ttquanti.h:424
bool get_h2_isBroad()
Definition: ttquanti.h:433
double get_genotype_with_dominance(const double a1, const double a2, const unsigned int locus, const unsigned int trait)
Definition: ttquanti.cc:2539
virtual double get_full_genotype(unsigned int trait)=0

References _num_trait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, _Vg, fatal(), Patch::get(), TMatrix::get(), TTQuanti::get_allele_bit(), TProtoQuanti::get_allele_model(), TProtoQuanti::get_diallele_values(), TTQuanti::get_full_genotype(), TProtoQuanti::get_genotype_with_dominance(), TProtoQuanti::get_h2_isBroad(), TProtoQuanti::get_locus_ID(), TProtoQuanti::get_locus_seq_pos(), TProtoQuanti::get_num_locus(), Patch::getID(), getSNPalleleFreqInPatch(), Individual::getTrait(), message(), my_variance_with_fixed_mean(), and Patch::size().

Referenced by LCE_QuantiModifier::setParameters(), and LCE_Breed_Quanti::setParameters().

◆ getVb()

double TTQuantiSH::getVb ( unsigned int  i)
inline
786 {return _Vb[i];}

References _Vb.

Referenced by addQuanti().

◆ getVg()

double TTQuantiSH::getVg ( unsigned int  i)
inline
785 {return _Vg[i];}

References _Vg.

Referenced by LCE_QuantiModifier::setVefromVa().

◆ getVp()

double TTQuantiSH::getVp ( unsigned int  i)
inline
787 {return _Vp[i];}

References _Vp.

Referenced by addQuanti().

◆ getVpPerPatch()

double TTQuantiSH::getVpPerPatch ( unsigned int  i,
unsigned int  p 
)
inline
794 {return _pVp[i][p];}

References _pVp.

Referenced by addVarPerPatch().

◆ init()

void TTQuantiSH::init ( )
virtual

Reimplemented from StatHandlerBase.

4604 {
4606 
4608 
4609  _eVar = (!_SHLinkedTrait->get_env_var().empty() || !_SHLinkedTrait->get_heritability().empty());
4610 
4611  resetPtrs(); //deallocate anything that was previously allocated
4612 
4613  if(_patchNbr != _pop->getPatchNbr())
4614  _patchNbr = _pop->getPatchNbr();
4615 
4618 
4621 
4622  _meanP = new double [_num_trait];
4623  _meanG = new double [_num_trait];
4624  _Va = new double [_num_trait];
4625  _Vg = new double [_num_trait];
4626  _Vb = new double [_num_trait];
4627  _Vp = new double [_num_trait];
4628 
4629  _pVa = new double* [_num_trait];
4630  _pVp = new double* [_num_trait];
4631  _pmeanP = new double* [_num_trait];
4632  _pmeanG = new double* [_num_trait];
4633 
4634  for(unsigned int i=0; i < _num_trait; ++i) {
4635 
4636  _pVa[i] = new double [_patchNbr];
4637  _pVp[i] = new double [_patchNbr];
4638  _pmeanP[i] = new double [_patchNbr];
4639  _pmeanG[i] = new double [_patchNbr];
4640 
4641  }
4642 
4643  if(_num_trait > 1) {
4644 
4645  _G = gsl_matrix_alloc(_num_trait,_num_trait);
4646  _eval = gsl_vector_alloc (_num_trait);
4647  _evec = gsl_matrix_alloc (_num_trait, _num_trait);
4648  _ws = gsl_eigen_symmv_alloc (_num_trait);
4649 
4650  _peigval = new double* [_num_trait];
4651 // _eigvect = new double* [_num_trait];
4652 
4653  for(unsigned int i=0; i < _num_trait; ++i) {
4654 // _eigvect[i] = new double [_num_trait];
4655  _peigval[i] = new double [_patchNbr];
4656  }
4657 
4658  _peigvect = new double* [_patchNbr];
4659 
4660  for(unsigned int i = 0; i < _patchNbr; i++)
4661  _peigvect[i] = new double [_num_trait*_num_trait];
4662 
4663  _covar = new double [_num_trait*(_num_trait -1)/2];
4664 
4665  _pcovar = new double* [_num_trait*(_num_trait - 1)/2];
4666 
4667  for(unsigned int i = 0; i < _num_trait*(_num_trait - 1)/2; i++)
4668  _pcovar[i] = new double [_patchNbr];
4669 
4670  }
4671 
4672 }
int getTraitIndex(trait_t type)
Gives the index of trait with type.
Definition: indfactory.cc:128
virtual void init()
Definition: stathandler.cc:39
vector< double > get_heritability()
Definition: ttquanti.h:431
virtual trait_t get_type() const
Definition: ttquanti.h:576
unsigned int get_num_traits()
Definition: ttquanti.h:423
vector< double > get_env_var()
Definition: ttquanti.h:430

References _covar, _eval, _eVar, _evec, _G, _meanG, _meanP, _num_locus, _num_trait, _patchNbr, _pcovar, _peigval, _peigvect, _pmeanG, _pmeanP, StatHandlerBase::_pop, _pVa, _pVp, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, _Va, _Vb, _Vg, _Vp, _ws, TProtoQuanti::get_env_var(), TProtoQuanti::get_heritability(), TProtoQuanti::get_num_locus(), TProtoQuanti::get_num_traits(), TProtoQuanti::get_type(), Metapop::getPatchNbr(), IndFactory::getTraitIndex(), StatHandlerBase::init(), and resetPtrs().

◆ resetPtrs()

void TTQuantiSH::resetPtrs ( )
4546 {
4547 
4548  if(_G) gsl_matrix_free(_G);
4549  if(_eval) gsl_vector_free(_eval);
4550  if(_evec) gsl_matrix_free(_evec);
4551  if(_ws) gsl_eigen_symmv_free (_ws);
4552 
4553  if(_meanP) delete [] _meanP; _meanP = NULL;
4554  if(_meanG) delete [] _meanG; _meanG = NULL;
4555  if(_Va) delete [] _Va; _Va = NULL;
4556  if(_Vg) delete [] _Vg; _Vg = NULL;
4557  if(_Vb) delete [] _Vb; _Vb = NULL;
4558  if(_Vp) delete [] _Vp; _Vp = NULL;
4559 
4560 
4561  if(_pVa) {
4562  for(unsigned int i=0; i < _num_trait; ++i) delete [] _pVa[i];
4563  delete [] _pVa; _pVa = NULL;
4564  }
4565  if(_pVp) {
4566  for(unsigned int i=0; i < _num_trait; ++i) delete [] _pVp[i];
4567  delete [] _pVp; _pVp = NULL;
4568  }
4569  if(_pmeanP) {
4570  for(unsigned int i=0; i < _num_trait; ++i) delete [] _pmeanP[i];
4571  delete [] _pmeanP; _pmeanP = NULL;
4572  }
4573  if(_pmeanG) {
4574  for(unsigned int i=0; i < _num_trait; ++i) delete [] _pmeanG[i];
4575  delete [] _pmeanG; _pmeanG = NULL;
4576  }
4577 
4578  if(_peigval) {
4579  for(unsigned int i=0; i < _num_trait; ++i) delete [] _peigval[i];
4580  delete [] _peigval;
4581  _peigval = NULL;
4582  }
4583 
4584  if(_covar) delete [] _covar;
4585  _covar = NULL;
4586 
4587  if(_pcovar) {
4588  for(unsigned int i = 0; i < (_num_trait*(_num_trait - 1)/2); i++) delete [] _pcovar[i];
4589  delete [] _pcovar;
4590  _pcovar = NULL;
4591  }
4592 
4593  if(_peigvect) {
4594  for(unsigned int i = 0; i < _patchNbr; i++) delete [] _peigvect[i];
4595  delete [] _peigvect;
4596  _peigvect = NULL;
4597  }
4598 
4599 }

References _covar, _eval, _evec, _G, _meanG, _meanP, _num_trait, _patchNbr, _pcovar, _peigval, _peigvect, _pmeanG, _pmeanP, _pVa, _pVp, _Va, _Vb, _Vg, _Vp, and _ws.

Referenced by init(), and ~TTQuantiSH().

◆ setAdultStats()

void TTQuantiSH::setAdultStats ( )
inline

◆ setDataTables()

void TTQuantiSH::setDataTables ( age_t  AGE)
5156 {
5157  unsigned int **sizes;
5158  unsigned int nb_patch = _pop->getPatchNbr();
5159 
5160  sizes = new unsigned int * [_num_trait];
5161 
5162  for(unsigned int i = 0; i < _num_trait; ++i) {
5163  sizes[i] = new unsigned int [nb_patch];
5164  for(unsigned int j = 0; j < nb_patch; ++j)
5165  sizes[i][j] = _pop->size(AGE, j);
5166  }
5167 
5168  _phenoTable.update(_num_trait, nb_patch, sizes);
5169  _genoTable.update(_num_trait, nb_patch, sizes);
5170 
5171  for(unsigned int i = 0; i < _num_trait; ++i)
5172  delete [] sizes[i];
5173  delete [] sizes;
5174 
5175  Patch* patch;
5176  age_idx age = (AGE == ADULTS ? ADLTx : OFFSx);
5177 
5178  for(unsigned int i = 0, n; i < nb_patch; i++) {
5179 
5180  patch = _pop->getPatch(i);
5181 
5182  n = 0;
5183 
5184  if ((patch->size(MAL, age)+patch->size(FEM, age)) != _phenoTable.size(0,i)) {
5185  fatal("problem while recording quanti trait values; table size doesn't match patch size.\n");
5186  }
5187 
5188  store_quanti_trait_values(patch, i, patch->size(MAL, age), &n, MAL, age, &_phenoTable, &_genoTable,
5190 
5191  store_quanti_trait_values(patch, i, patch->size(FEM, age), &n, FEM, age, &_phenoTable, &_genoTable,
5193 
5194  if (n != _phenoTable.size(0,i) || n != _genoTable.size(0,i)) {
5195  fatal("problem while recording quanti trait values; size counter doesn't match table size.\n");
5196  }
5197  } // end for patches
5198 }
void update(unsigned int nbgroups, unsigned int nbclasses, unsigned int **classSizes)
Updates the group and classe sizes and re-allocates the table according to its new length.
Definition: datatable.h:154
unsigned int size()
Get the total number of individuals present in the population, all sex and age classes together.
Definition: metapop.h:312
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
DataTable< double > _genoTable
Definition: ttquanti.h:746
void store_quanti_trait_values(Patch *patch, unsigned int patchID, unsigned int size, unsigned int *cntr, sex_t SEX, age_idx AGE, DataTable< double > *ptable, DataTable< double > *gtable, unsigned int nTrait, unsigned int TraitIndex)
Definition: ttquanti.cc:5202
age_idx
Array index of the age classes in the patch sizes and containers arrays.
Definition: types.h:41
@ OFFSx
Definition: types.h:42
@ ADLTx
Definition: types.h:42

References _genoTable, _num_trait, _phenoTable, StatHandlerBase::_pop, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, ADLTx, ADULTS, fatal(), FEM, Metapop::getPatch(), Metapop::getPatchNbr(), MAL, OFFSx, Metapop::size(), Patch::size(), DataTable< T >::size(), store_quanti_trait_values(), and DataTable< T >::update().

Referenced by setStats().

◆ setOffsprgStats()

◆ setStatRecorders()

bool TTQuantiSH::setStatRecorders ( std::string &  token)
virtual

Implements StatHandlerBase.

4677 {
4678 #ifdef _DEBUG_
4679  message("-TTQuantiSH::setStatRecorders ");
4680 #endif
4681  string age_tag = token.substr(0,token.find_first_of("."));
4682  string sub_token;
4683  age_t AGE = ALL;
4684 
4685  if (age_tag.size() != 0 && age_tag.size() != string::npos) {
4686 
4687  if (age_tag == "adlt") AGE = ADULTS;
4688 
4689  else if (age_tag == "off") AGE = OFFSPRG;
4690 
4691  else age_tag = "";
4692 
4693  } else {
4694  age_tag = "";
4695  }
4696 
4697  if (age_tag.size() != 0)
4698  sub_token = token.substr(token.find_first_of(".") + 1, string::npos);
4699  else
4700  sub_token = token;
4701 
4702  if(sub_token == "quanti") {
4703  addQuanti(AGE);
4704  } else if(sub_token == "quanti.patch") {
4705  addQuantiPerPatch(AGE);
4706  } else if(sub_token == "quanti.mean.patch") {
4707  addAvgPerPatch(AGE);
4708  } else if(sub_token == "quanti.pheno.patch") { // same as quanti.mean.patch
4709  addAvgPerPatch(AGE);
4710  } else if(sub_token == "quanti.geno.patch") {
4711  addGenotPerPatch(AGE);
4712  } else if(sub_token == "quanti.var.patch") {
4713  addVarPerPatch(AGE);
4714  } else if(sub_token == "quanti.covar.patch") {
4715  addCovarPerPatch(AGE);
4716  } else if(sub_token == "quanti.eigen.patch") {
4717  addEigenPerPatch(AGE);
4718  } else if(sub_token == "quanti.eigenvalues.patch") {
4720  } else if(sub_token == "quanti.eigenvect1.patch") {
4721  addEigenVect1PerPatch(AGE);
4722  } else if(sub_token == "quanti.skew.patch") {
4723  addSkewPerPatch(AGE);
4724 // } else if((sub_token == "quanti.epistasis") && (_SHLinkedTrait->do_epistasis())) {
4725 // addQuantiEpistasis(AGE);
4726 // _epistats = true;
4727 // } else if(sub_token == "quanti.effect.stats") {
4728 // addEffectStats(AGE);
4729 // _alleleStats = true;
4730  } else
4731  return false;
4732 
4733  return true;
4734 }
unsigned int age_t
Age class flags.
Definition: types.h:46

References addAvgPerPatch(), addCovarPerPatch(), addEigenPerPatch(), addEigenValuesPerPatch(), addEigenVect1PerPatch(), addGenotPerPatch(), addQuanti(), addQuantiPerPatch(), addSkewPerPatch(), addVarPerPatch(), ADULTS, ALL, message(), and OFFSPRG.

◆ setStats()

void TTQuantiSH::setStats ( age_t  AGE)
5227 {
5228  if(_table_set_age == AGE
5231  return;
5232 
5233  unsigned int pop_size = _pop->size(AGE);
5234  unsigned int patch_size;
5235  double *phenot1, *genot1, *genot2;
5236 
5237  unsigned int nb_patch = _pop->getPatchNbr();
5238 
5239  if(nb_patch < _patchNbr) {
5240  warning("increase in patch number detected (in Quanti Stat Handler),");
5241  warning("stats for quanti trait will not be recorded in new patches, patch identity may have changed.\n");
5242  _patchNbr = nb_patch; // record stat only in remaining patches
5243  }
5244 
5245  if(nb_patch > _patchNbr) nb_patch = _patchNbr; //tables have been allocated for _patchNbr patches
5246 
5247  setDataTables(AGE);
5248 
5249 #ifdef HAS_GSL
5250  unsigned int c = 0; //covariance position counter
5251  unsigned int pv = 0; //eigenvector position counter
5252 
5253  //within deme stats:
5254  for(unsigned int j = 0; j < nb_patch; j++) {
5255 
5256  patch_size = _pop->size(AGE, j);
5257 
5258  for(unsigned int t=0; t < _num_trait; t++) {
5259 
5260  phenot1 = _phenoTable.getClassWithinGroup(t,j);
5261  genot1 = _genoTable.getClassWithinGroup(t,j);
5262 
5263  _pmeanP[t][j] = my_mean (phenot1, patch_size );
5264  _pmeanG[t][j] = my_mean (genot1, patch_size );
5265  _pVp[t][j] = my_variance_with_fixed_mean (phenot1, patch_size, _pmeanP[t][j]);
5266  _pVa[t][j] = my_variance_with_fixed_mean (genot1, patch_size, _pmeanG[t][j]);
5267 
5268  }
5269 
5270  if(_num_trait > 1) {
5271  c = 0;
5272  // calculate the covariances and G, need to adjust dimensions in class declaration
5273  for(unsigned int t1 = 0; t1 < _num_trait; t1++) {
5274  // set the diagonal elements of G here
5275 
5276  gsl_matrix_set(_G, t1, t1, _pVa[t1][j]);
5277 
5278  for(unsigned int t2 = t1 + 1; t2 < _num_trait; t2++) {
5279 
5280  genot1 = _genoTable.getClassWithinGroup(t1,j);
5281  genot2 = _genoTable.getClassWithinGroup(t2,j);
5282 
5283  _pcovar[c][j] = gsl_stats_covariance_m (genot1, 1, genot2, 1, patch_size,
5284  _pmeanG[t1][j], _pmeanG[t2][j]);
5285 
5286  gsl_matrix_set(_G, t1, t2, _pcovar[c][j]);
5287  gsl_matrix_set(_G, t2, t1, _pcovar[c++][j]);
5288  }
5289  }
5290 
5291  gsl_eigen_symmv (_G, _eval, _evec, _ws);
5292  gsl_eigen_symmv_sort (_eval, _evec, GSL_EIGEN_SORT_VAL_DESC);
5293 
5294  pv = 0;
5295 
5296  for(unsigned int t = 0; t < _num_trait; t++) {
5297  _peigval[t][j] = gsl_vector_get (_eval, t);
5298  for(unsigned int v = 0; v < _num_trait; v++)
5299  _peigvect[j][pv++] = gsl_matrix_get (_evec, v, t); //read eigenvectors column-wise
5300  }
5301  }
5302  } //_end_ for patch
5303 
5304  double meanGamong1;
5305  c = 0; //reset covariance positioner
5306 
5307  //among demes stats:
5308  for(unsigned int t1 = 0; t1 < _num_trait; t1++) {
5309 
5310  phenot1 = _phenoTable.getGroup(t1);
5311  genot1 = _genoTable.getGroup(t1);
5312 
5313  _meanP[t1] = my_mean (phenot1, pop_size ); //grand mean, pooling all individuals
5314  _meanG[t1] = my_mean (genot1, pop_size ); //grand mean, pooling all individuals
5315 
5316  _Vp[t1] = my_mean_no_nan (_pVp[t1], nb_patch); //mean within patch variance
5317  _Va[t1] = my_mean_no_nan (_pVa[t1], nb_patch); //mean within patch variance
5318 
5319  meanGamong1 = my_mean_no_nan (_pmeanG[t1], nb_patch); //mean of within patch mean genotypic values
5320 
5321  _Vb[t1] = my_variance_with_fixed_mean_no_nan (_pmeanG[t1], nb_patch, meanGamong1); //variance of patch means
5322 
5323  for(unsigned int t2 = t1 + 1; t2 < _num_trait; t2++) {
5324 
5325  _covar[c] = my_covariance_no_nan(_pmeanG[t1], _pmeanG[t2], nb_patch, nb_patch);
5326 
5327  //set covariance to the mean covariance within patch for recording
5328  _covar[c] = my_mean_no_nan (_pcovar[c], nb_patch);
5329 
5330  c++;
5331  }
5332 
5333  } //_end_ for trait
5334 
5335 #else
5336  fatal("install the GSL library to get the quanti stats!\n");
5337 #endif
5338 
5339  _table_set_age = AGE;
5342 
5343 }
T * getGroup(unsigned int group)
Accessor to a group array.
Definition: datatable.h:223
unsigned int getCurrentReplicate()
Definition: metapop.h:295
unsigned int getCurrentGeneration()
Definition: metapop.h:296
void setDataTables(age_t AGE)
Definition: ttquanti.cc:5155
double my_variance_with_fixed_mean_no_nan(double *data, unsigned int size, double mean)
Definition: utils.cc:74
double my_mean(double *data, unsigned int size)
Definition: utils.cc:38
double my_mean_no_nan(double *data, unsigned int size)
Definition: utils.cc:49
double my_covariance_no_nan(double *data1, double *data2, unsigned int size1, unsigned int size2)
Definition: utils.cc:88

References _covar, _eval, _evec, _G, _genoTable, _meanG, _meanP, _num_trait, _patchNbr, _pcovar, _peigval, _peigvect, _phenoTable, _pmeanG, _pmeanP, StatHandlerBase::_pop, _pVa, _pVp, _table_set_age, _table_set_gen, _table_set_repl, _Va, _Vb, _Vp, _ws, fatal(), DataTable< T >::getClassWithinGroup(), Metapop::getCurrentGeneration(), Metapop::getCurrentReplicate(), DataTable< T >::getGroup(), Metapop::getPatchNbr(), my_covariance_no_nan(), my_mean(), my_mean_no_nan(), my_variance_with_fixed_mean(), my_variance_with_fixed_mean_no_nan(), setDataTables(), Metapop::size(), and warning().

Referenced by setAdultStats(), and setOffsprgStats().

Member Data Documentation

◆ _covar

double * TTQuantiSH::_covar
private

Referenced by getCovar(), init(), resetPtrs(), and setStats().

◆ _eval

gsl_vector* TTQuantiSH::_eval
private

Referenced by init(), resetPtrs(), and setStats().

◆ _eVar

bool TTQuantiSH::_eVar
private

Referenced by addQuanti(), addVarPerPatch(), and init().

◆ _evec

gsl_matrix * TTQuantiSH::_evec
private

Referenced by init(), resetPtrs(), and setStats().

◆ _G

gsl_matrix* TTQuantiSH::_G
private

Referenced by init(), resetPtrs(), and setStats().

◆ _genoTable

DataTable< double > TTQuantiSH::_genoTable
private

Referenced by setDataTables(), and setStats().

◆ _meanG

double * TTQuantiSH::_meanG
private

◆ _meanP

double* TTQuantiSH::_meanP
private

◆ _num_locus

unsigned int TTQuantiSH::_num_locus
private

Referenced by getVaNoDominance(), and init().

◆ _num_trait

◆ _patchNbr

unsigned int TTQuantiSH::_patchNbr
private

Referenced by init(), resetPtrs(), and setStats().

◆ _pcovar

double ** TTQuantiSH::_pcovar
private

◆ _peigval

double ** TTQuantiSH::_peigval
private

◆ _peigvect

double ** TTQuantiSH::_peigvect
private

◆ _phenoTable

DataTable< double > TTQuantiSH::_phenoTable
private

◆ _pmeanG

double ** TTQuantiSH::_pmeanG
private

◆ _pmeanP

double** TTQuantiSH::_pmeanP
private

◆ _pVa

double ** TTQuantiSH::_pVa
private

◆ _pVp

double ** TTQuantiSH::_pVp
private

◆ _table_set_age

unsigned int TTQuantiSH::_table_set_age
private

Referenced by setStats().

◆ _table_set_gen

unsigned int TTQuantiSH::_table_set_gen
private

Referenced by setStats().

◆ _table_set_repl

unsigned int TTQuantiSH::_table_set_repl
private

Referenced by setStats().

◆ _Va

double * TTQuantiSH::_Va
private

Referenced by getQst(), getVa(), init(), resetPtrs(), and setStats().

◆ _Vb

double * TTQuantiSH::_Vb
private

Referenced by getQst(), getVb(), init(), resetPtrs(), and setStats().

◆ _Vg

double * TTQuantiSH::_Vg
private

◆ _Vp

double * TTQuantiSH::_Vp
private

Referenced by getVp(), init(), resetPtrs(), and setStats().

◆ _ws

gsl_eigen_symmv_workspace* TTQuantiSH::_ws
private

Referenced by init(), resetPtrs(), and setStats().


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