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

#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 addEigen (age_t AGE)
 
void addEigenValues (age_t AGE)
 
void addEigenVect1 (age_t AGE)
 
void addQuantiPerPatch (age_t AGE)
 
void addAvgPerPatch (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 addEigenStatsPerPatcg (age_t AGE)
 
void addSkewPerPatch (age_t AGE)
 
void setDataTables (age_t AGE)
 
void setAdultStats ()
 
void setOffsprgStats ()
 
void setStats (age_t AGE)
 
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 getEigenValue (unsigned int i)
 
double getEigenVectorElt (unsigned int t1, unsigned int t2)
 
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)
 
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 void init ()=0
 Inits state. More...
 
virtual void update ()=0
 Updates the handler state. More...
 
virtual ~Handler ()
 

Private Attributes

double * _meanP
 
double * _meanG
 
double * _Va
 
double * _Vg
 
double * _Vb
 
double * _Vp
 
double * _covar
 
double * _eigval
 
double ** _eigvect
 
double ** _pmeanP
 
double ** _pmeanG
 
double ** _pVa
 
double ** _pVp
 
double ** _pcovar
 
double ** _peigval
 
double ** _peigvect
 
unsigned int _nb_locus
 
unsigned int _nb_trait
 
unsigned int _patchNbr
 
bool _eVar
 
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...
 

Constructor & Destructor Documentation

◆ TTQuantiSH()

TTQuantiSH::TTQuantiSH ( TProtoQuanti TP)
inline
276 _meanP(0), _meanG(0), _Va(0), _Vg(0), _Vb(0), _Vp(0), _covar(0), _eigval(0), _eigvect(0),
277 _pmeanP(0), _pmeanG(0), _pVa(0), _pVp(0), _pcovar(0), _peigval(0), _peigvect(0),
278 _nb_locus(0), _nb_trait(0),_patchNbr(0), _eVar(0),_G(0),_evec(0),_eval(0),_ws(0),
279 _table_set_gen(999999), _table_set_age(999999), _table_set_repl(999999)
280 {}
double * _Vp
Definition: ttquanti.h:259
double ** _pVa
Definition: ttquanti.h:260
double * _Vb
Definition: ttquanti.h:259
unsigned int _nb_trait
Definition: ttquanti.h:262
unsigned int _nb_locus
Definition: ttquanti.h:262
double * _meanP
Definition: ttquanti.h:259
double ** _pVp
Definition: ttquanti.h:260
bool _eVar
Definition: ttquanti.h:263
gsl_matrix * _evec
Definition: ttquanti.h:265
unsigned int _table_set_age
Definition: ttquanti.h:270
double ** _pmeanP
Definition: ttquanti.h:260
double * _meanG
Definition: ttquanti.h:259
gsl_vector * _eval
Definition: ttquanti.h:266
double * _covar
Definition: ttquanti.h:259
gsl_matrix * _G
Definition: ttquanti.h:265
double ** _peigvect
Definition: ttquanti.h:260
unsigned int _table_set_repl
Definition: ttquanti.h:270
gsl_eigen_symmv_workspace * _ws
Definition: ttquanti.h:267
unsigned int _table_set_gen
Definition: ttquanti.h:270
double * _Va
Definition: ttquanti.h:259
double * _Vg
Definition: ttquanti.h:259
double ** _peigval
Definition: ttquanti.h:260
unsigned int _patchNbr
Definition: ttquanti.h:262
double ** _pcovar
Definition: ttquanti.h:260
double ** _pmeanG
Definition: ttquanti.h:260
double * _eigval
Definition: ttquanti.h:259
double ** _eigvect
Definition: ttquanti.h:259
Template class for the trait's StatHandler.
Definition: stathandler.h:168

◆ ~TTQuantiSH()

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

References resetPtrs().

Member Function Documentation

◆ addAvgPerPatch()

void TTQuantiSH::addAvgPerPatch ( age_t  AGE)
2097{
2098 if (AGE == ALL) {
2101 return;
2102 }
2103
2104 unsigned int patchNbr = _pop->getPatchNbr();
2105
2106 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2107 string name;
2108 string patch;
2109 string t1;
2110
2111 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2113
2114 add("Mean phenotype of trait 1 in patch 1", suffix + "q1.p1", AGE, 0, 0, 0, 0,
2116
2117 for(unsigned int p = 0; p < patchNbr; p++) {
2118 for(unsigned int i = 0; i < _nb_trait; i++) {
2119 if(p == 0 && i == 0) continue;
2120 name = "Mean phenotype of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
2121 t1 = "q" + tstring::int2str(i+1);
2122 patch = ".p" + tstring::int2str(p+1);
2123 add(name, suffix + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getMeanPhenotPerPatch, 0);
2124 }
2125 }
2126
2127}
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
Definition: ttquanti.h:257
void addAvgPerPatch(age_t AGE)
Definition: ttquanti.cc:2096
double getMeanPhenotPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:317
void setAdultStats()
Definition: ttquanti.h:304
void setOffsprgStats()
Definition: ttquanti.h:305
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 _nb_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), addAvgPerPatch(), ADULTS, ALL, getMeanPhenotPerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by addAvgPerPatch(), addQuantiPerPatch(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ addCovarPerPatch()

void TTQuantiSH::addCovarPerPatch ( age_t  AGE)
2178{
2179 if(_nb_trait < 2) {
2180 // warning("not recording traits covariance with only one \"quanti\" trait!");
2181 return;
2182 }
2183
2184 if (AGE == ALL) {
2187 return;
2188 }
2189
2190 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2191 string patch;
2192 string cov;
2193 unsigned int patchNbr = _pop->getPatchNbr();
2194
2195 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2197
2198 add("Genetic covariance of trait 1 and trait 2 in patch 1", suffix + "cov.q12.p1", AGE, 0, 0, 0, 0,
2200
2201 unsigned int c;
2202 for(unsigned int p = 0; p < patchNbr; p++) {
2203 patch = ".p" + tstring::int2str(p+1);
2204 c = 0;
2205 for(unsigned int t = 0; t < _nb_trait; ++t) {
2206 for(unsigned int v = t + 1; v < _nb_trait; ++v){
2207 if(p==0 && t==0 && v==1) {c++; continue;} //this on is already recorder above
2208 cov = tstring::int2str((t+1)*10+v+1);
2209 add("", suffix + "cov.q" + cov + patch, AGE, c++, p, 0, 0, &TTQuantiSH::getCovarPerPatch, 0);
2210 }
2211 }
2212 }
2213
2214}
void addCovarPerPatch(age_t AGE)
Definition: ttquanti.cc:2177
double getCovarPerPatch(unsigned int p, unsigned int i)
Definition: ttquanti.h:321

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

Referenced by addCovarPerPatch(), addQuantiPerPatch(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ addEigen()

void TTQuantiSH::addEigen ( age_t  AGE)
1993{
1994 if(_nb_trait < 2) {
1995 warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
1996 return;
1997 }
1998
1999 if (AGE == ALL) {
2002 return;
2003 }
2004
2005 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2006
2007 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2009
2010 add("", suffix + "q.eval1",AGE,0,0,0,&TTQuantiSH::getEigenValue, 0, setter); //this one calls the setter
2011
2012 for(unsigned int t = 1; t < _nb_trait; ++t)
2013 add("", suffix + "q.eval" + tstring::int2str(t+1),AGE,t,0,0,&TTQuantiSH::getEigenValue, 0, 0);
2014
2015 for(unsigned int t = 0; t< _nb_trait; ++t)
2016 for(unsigned int v = 0; v < _nb_trait; ++v)
2017 add("", suffix + "q.evect" + tstring::int2str((t+1)*10+(v+1)),AGE,t,v,0,0,&TTQuantiSH::getEigenVectorElt,0);
2018
2019}
double getEigenVectorElt(unsigned int t1, unsigned int t2)
Definition: ttquanti.h:315
void addEigen(age_t AGE)
Definition: ttquanti.cc:1992
double getEigenValue(unsigned int i)
Definition: ttquanti.h:314
void warning(const char *str,...)
Definition: output.cc:58

References _nb_trait, StatHandler< SH >::add(), addEigen(), ADULTS, ALL, getEigenValue(), getEigenVectorElt(), tstring::int2str(), OFFSPRG, setAdultStats(), setOffsprgStats(), and warning().

Referenced by addEigen(), addEigenValues(), addEigenVect1(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ addEigenPerPatch()

void TTQuantiSH::addEigenPerPatch ( age_t  AGE)
2219{
2220 if(_nb_trait < 2) {
2221 warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
2222 return;
2223 }
2224
2225 if (AGE == ALL) {
2228 return;
2229 }
2230
2231 unsigned int patchNbr = _pop->getPatchNbr();
2232 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2233 string patch;
2234 unsigned int pv =0;
2235
2236 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2238
2239
2240 add("First G-matrix eigenvalue in patch 1", suffix + "qeval1.p1", AGE, 0, 0, 0, 0,
2242
2243 for(unsigned int p = 0; p < patchNbr; ++p) {
2244 patch = ".p" + tstring::int2str(p+1);
2245 for(unsigned int t = 0; t < _nb_trait; ++t) {
2246 if(p==0 && t==0) continue;
2247 add("", suffix + "qeval" + tstring::int2str(t+1) + patch, AGE, t, p, 0, 0, &TTQuantiSH::getEigenValuePerPatch,0);
2248 }
2249 }
2250 for(unsigned int p = 0; p < patchNbr; ++p) {
2251 patch = ".p" + tstring::int2str(p+1);
2252 pv = 0;
2253 for(unsigned int t = 0; t < _nb_trait; ++t)
2254 for(unsigned int v = 0; v < _nb_trait; ++v)
2255 add("", suffix + "qevect" + tstring::int2str((t+1)*10+v+1) + patch, AGE, p, pv++, 0, 0, &TTQuantiSH::getEigenVectorEltPerPatch,0);
2256 }
2257
2258}
double getEigenVectorEltPerPatch(unsigned int p, unsigned int v)
Definition: ttquanti.h:322
double getEigenValuePerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:320
void addEigenPerPatch(age_t AGE)
Definition: ttquanti.cc:2218

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

Referenced by addEigenPerPatch(), addEigenValuesPerPatch(), addEigenVect1PerPatch(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ addEigenStatsPerPatcg()

void TTQuantiSH::addEigenStatsPerPatcg ( age_t  AGE)

◆ addEigenValues()

void TTQuantiSH::addEigenValues ( age_t  AGE)
2024{
2025 if(_nb_trait < 2) {
2026 warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
2027 return;
2028 }
2029
2030 if (AGE == ALL) {
2033 return;
2034 }
2035
2036 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2037
2038 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2040
2041 add("", suffix + "q.eval1",AGE,0,0,0,&TTQuantiSH::getEigenValue, 0, setter);
2042
2043 for(unsigned int t = 1; t < _nb_trait; ++t)
2044 add("", suffix + "q.eval" + tstring::int2str(t+1),AGE,t,0,0,&TTQuantiSH::getEigenValue, 0, 0);
2045
2046
2047}

References _nb_trait, StatHandler< SH >::add(), addEigen(), ADULTS, ALL, getEigenValue(), tstring::int2str(), OFFSPRG, setAdultStats(), setOffsprgStats(), and warning().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ addEigenValuesPerPatch()

void TTQuantiSH::addEigenValuesPerPatch ( age_t  AGE)
2263{
2264 if(_nb_trait < 2) {
2265 warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
2266 return;
2267 }
2268
2269 if (AGE == ALL) {
2272 return;
2273 }
2274
2275 unsigned int patchNbr = _pop->getPatchNbr();
2276 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2277 string patch;
2278
2279 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2281
2282 add("First G-matrix eigenvalue in patch 1", suffix + "qeval1.p1", AGE, 0, 0, 0, 0,
2284
2285 for(unsigned int p = 0; p < patchNbr; ++p) {
2286 patch = ".p" + tstring::int2str(p+1);
2287 for(unsigned int t = 0; t < _nb_trait; ++t) {
2288 if(p==0 && t==0) continue;
2289 add("", suffix + "qeval" + tstring::int2str(t+1) + patch, AGE, t, p, 0, 0, &TTQuantiSH::getEigenValuePerPatch,0);
2290 }
2291 }
2292}

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

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ addEigenVect1()

void TTQuantiSH::addEigenVect1 ( age_t  AGE)
2052{
2053 if(_nb_trait < 2) {
2054 warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
2055 return;
2056 }
2057
2058 if (AGE == ALL) {
2061 return;
2062 }
2063
2064 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2065
2066 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2068
2069 add("", suffix + "q.evect11",AGE,0,0,0,0,&TTQuantiSH::getEigenVectorElt,setter);
2070
2071 for(unsigned int v = 1; v < _nb_trait; ++v)
2072 add("", suffix + "q.evect1" + tstring::int2str(v+1),AGE,0,v,0,0,&TTQuantiSH::getEigenVectorElt,0);
2073
2074}

References _nb_trait, StatHandler< SH >::add(), addEigen(), ADULTS, ALL, getEigenVectorElt(), tstring::int2str(), OFFSPRG, setAdultStats(), setOffsprgStats(), and warning().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ addEigenVect1PerPatch()

void TTQuantiSH::addEigenVect1PerPatch ( age_t  AGE)
2297{
2298 if(_nb_trait < 2) {
2299 warning("not recording G-matrix eigen stats with only one \"quanti\" trait!\n");
2300 return;
2301 }
2302
2303 if (AGE == ALL) {
2306 return;
2307 }
2308
2309 unsigned int patchNbr = _pop->getPatchNbr();
2310 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2311 string patch;
2312 unsigned int pv =0;
2313
2314 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2317
2318
2319 add("First G-matrix eigenvector in patch 1", suffix + "qevect11.p1", AGE, 0, 0, 0, 0,
2321
2322 for(unsigned int p = 0; p < patchNbr; ++p) {
2323 patch = ".p" + tstring::int2str(p+1);
2324 pv = 0;
2325 // for(unsigned int t = 0; t < _nb_trait; ++t)
2326 for(unsigned int v = 0; v < _nb_trait; ++v){
2327 if(p==0 && v==0) {pv++; continue;}
2328 add("", suffix + "qevect1" + tstring::int2str(v+1) + patch, AGE, p, pv++, 0, 0, &TTQuantiSH::getEigenVectorEltPerPatch,0);
2329 }
2330 }
2331}

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

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ addQuanti()

void TTQuantiSH::addQuanti ( age_t  AGE)
1924{
1925 if (AGE == ALL) {
1928 return;
1929 }
1930
1931 string suffix = (AGE == ADULTS ? "adlt.":"off.");
1932 string name = suffix + "q";
1933 string t1, t2;
1934
1935 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
1937
1938 add("",suffix + "q1",AGE,0,0,0,&TTQuantiSH::getMeanPhenot,0,setter);
1939
1940 for(unsigned int i = 1; i < _nb_trait; i++) {
1941 t1 = tstring::int2str(i+1);
1942 add("", name + t1,AGE,i,0,0,&TTQuantiSH::getMeanPhenot,0,0);
1943 }
1944
1945 //Va, no dominance
1947 for(unsigned int i = 0; i < _nb_trait; i++) {
1948 t1 = tstring::int2str(i+1);
1949 add("", name + t1 +".Va",AGE,i,0,0,&TTQuantiSH::getVa,0,0);
1950 }
1951 } else {
1952 //if dominance, then what we record is Vg, not Va
1953 for(unsigned int i = 0; i < _nb_trait; i++) {
1954 t1 = tstring::int2str(i+1);
1955 add("", name + t1 +".Vg",AGE,i,0,0,&TTQuantiSH::getVa,0,0);
1956 }
1957 }
1958
1959 //Vp
1960 if(_eVar) {
1961 for(unsigned int i = 0; i < _nb_trait; i++) {
1962 t1 = tstring::int2str(i+1);
1963 add("", name + t1 +".Vp",AGE,i,0,0,&TTQuantiSH::getVp,0,0);
1964 }
1965 }
1966 //Vb
1967 if(_pop->getPatchNbr() > 1){
1968 for(unsigned int i = 0; i < _nb_trait; i++) {
1969 t1 = tstring::int2str(i+1);
1970 add("", name + t1 +".Vb",AGE,i,0,0,&TTQuantiSH::getVb,0,0);
1971 }
1972 //Qst
1973 for(unsigned int i = 0; i < _nb_trait; i++) {
1974 t1 = tstring::int2str(i+1);
1975 add("", name + t1 +".Qst",AGE,i,0,0,&TTQuantiSH::getQst,0,0);
1976 }
1977 }
1978
1979 if (_nb_trait > 1) {
1980 unsigned int c = 0;
1981 for(unsigned int t = 0; t < _nb_trait; ++t) {
1982 for(unsigned int v = t + 1; v < _nb_trait; ++v) {
1983 t1 = tstring::int2str((t+1)*10+(v+1));
1984 add("", name + t1 +".cov", AGE, c++, 0,0,&TTQuantiSH::getCovar,0,0);
1985 }
1986 }
1987 }
1988}
unsigned int get_dominance_model()
Definition: ttquanti.h:146
void addQuanti(age_t AGE)
Definition: ttquanti.cc:1923
double getVa(unsigned int i)
Definition: ttquanti.h:308
double getVb(unsigned int i)
Definition: ttquanti.h:310
double getMeanPhenot(unsigned int i)
Definition: ttquanti.h:307
double getQst(unsigned int i)
Definition: ttquanti.h:312
double getVp(unsigned int i)
Definition: ttquanti.h:311
double getCovar(unsigned int i)
Definition: ttquanti.h:313
TProtoQuanti * _SHLinkedTrait
Pointer to a TraitProtoype object.
Definition: stathandler.h:171

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

Referenced by addQuanti(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ addQuantiPerPatch()

void TTQuantiSH::addQuantiPerPatch ( age_t  AGE)
2079{
2080
2081 if (AGE == ALL) {
2084 return;
2085 }
2086
2087 addAvgPerPatch(AGE);
2088 addVarPerPatch(AGE);
2089 addCovarPerPatch(AGE);
2090 // addEigenPerPatch(AGE);
2091
2092}
void addVarPerPatch(age_t AGE)
Definition: ttquanti.cc:2131
void addQuantiPerPatch(age_t AGE)
Definition: ttquanti.cc:2078

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

Referenced by addQuantiPerPatch(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ addSkewPerPatch()

void TTQuantiSH::addSkewPerPatch ( age_t  AGE)
2336{
2337 if (AGE == ALL) {
2340 return;
2341 }
2342
2343 unsigned int patchNbr = _pop->getPatchNbr();
2344
2345 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2346 string name;
2347 string patch;
2348 string t1;
2349
2350 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2353
2354 add("Genetic skew of trait 1 in patch 1", suffix + "Sk.q1.p1", AGE, 0, 0, 0, 0,
2356
2357 for(unsigned int p = 0; p < patchNbr; p++) {
2358 for(unsigned int i = 0; i < _nb_trait; i++) {
2359 if(p == 0 && i == 0) continue;
2360 name = "Genetic skew of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
2361 t1 = "q" + tstring::int2str(i+1);
2362 patch = ".p" + tstring::int2str(p+1);
2363 add(name, suffix + "Sk." + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getSkewPerPatch, 0);
2364 }
2365 }
2366
2367}
double getSkewPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.cc:2371
void addSkewPerPatch(age_t AGE)
Definition: ttquanti.cc:2335

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

Referenced by addSkewPerPatch(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ addVarPerPatch()

void TTQuantiSH::addVarPerPatch ( age_t  AGE)
2132{
2133 if (AGE == ALL) {
2136 return;
2137 }
2138
2139 unsigned int patchNbr = _pop->getPatchNbr();
2140
2141 string suffix = (AGE == ADULTS ? "adlt.":"off.");
2142 string name;
2143 string patch;
2144 string t1;
2145
2146 void (TTQuantiSH::* setter) (void) = (AGE == ADULTS ?
2149
2150 add("Genetic variance of trait 1 in patch 1", suffix + "Va.q1.p1", AGE, 0, 0, 0, 0,
2151 &TTQuantiSH::getVaPerPatch, setter);
2152
2153 for(unsigned int p = 0; p < patchNbr; p++) {
2154 for(unsigned int i = 0; i < _nb_trait; i++) {
2155 if(p == 0 && i == 0) continue;
2156 name = "Genetic variance of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
2157 t1 = "q" + tstring::int2str(i+1);
2158 patch = ".p" + tstring::int2str(p+1);
2159 add(name, suffix + "Va." + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getVaPerPatch, 0);
2160 }
2161 }
2162
2163 if(_eVar) {
2164 for(unsigned int p = 0; p < patchNbr; p++) {
2165 for(unsigned int i = 0; i < _nb_trait; i++) {
2166 name = "Phenotypic variance of trait " + tstring::int2str(i+1) + " in patch " + tstring::int2str(p+1);
2167 t1 = "q" + tstring::int2str(i+1);
2168 patch = ".p" + tstring::int2str(p+1);
2169 add(name, suffix + "Vp." + t1 + patch, AGE, i, p, 0, 0, &TTQuantiSH::getVpPerPatch, 0);
2170 }
2171 }
2172 }
2173}
double getVpPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:319
double getVaPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:318

References _eVar, _nb_trait, StatHandlerBase::_pop, StatHandler< SH >::add(), addVarPerPatch(), ADULTS, ALL, Metapop::getPatchNbr(), getVaPerPatch(), getVpPerPatch(), tstring::int2str(), OFFSPRG, setAdultStats(), and setOffsprgStats().

Referenced by addQuantiPerPatch(), addVarPerPatch(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ getCovar()

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

References _covar.

Referenced by addQuanti().

+ Here is the caller graph for this function:

◆ getCovarPerPatch()

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

References _pcovar.

Referenced by addCovarPerPatch().

+ Here is the caller graph for this function:

◆ getEigenValue()

double TTQuantiSH::getEigenValue ( unsigned int  i)
inline
314{return _eigval[i];}

References _eigval.

Referenced by addEigen(), and addEigenValues().

+ Here is the caller graph for this function:

◆ getEigenValuePerPatch()

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

References _peigval.

Referenced by addEigenPerPatch(), and addEigenValuesPerPatch().

+ Here is the caller graph for this function:

◆ getEigenVectorElt()

double TTQuantiSH::getEigenVectorElt ( unsigned int  t1,
unsigned int  t2 
)
inline
315{return _eigvect[t2][t1];}//eigenvectors arranged column-wise

References _eigvect.

Referenced by addEigen(), and addEigenVect1().

+ Here is the caller graph for this function:

◆ getEigenVectorEltPerPatch()

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

References _peigvect.

Referenced by addEigenPerPatch(), and addEigenVect1PerPatch().

+ Here is the caller graph for this function:

◆ getMeanPhenot()

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

References _meanP.

Referenced by addQuanti().

+ Here is the caller graph for this function:

◆ getMeanPhenotPerPatch()

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

References _pmeanP.

Referenced by addAvgPerPatch().

+ Here is the caller graph for this function:

◆ getQst()

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

References _Va, and _Vb.

Referenced by addQuanti().

+ Here is the caller graph for this function:

◆ getSkewPerPatch()

double TTQuantiSH::getSkewPerPatch ( unsigned int  i,
unsigned int  p 
)
2372{
2373 double skew = 0;
2374
2375 double *phenot = _phenoTable.getClassWithinGroup(i, p);
2376 unsigned int patch_size = _phenoTable.size(i, p);
2377
2378 for(unsigned int k = 0; k < patch_size; ++k)
2379 skew += pow( phenot[k] - _pmeanP[i][p], 3); //the mean has been set by setStats()
2380
2381 return skew / patch_size;
2382}
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:269

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

Referenced by addSkewPerPatch().

+ Here is the caller graph for this function:

◆ getSNPalleleFreqInPatch()

vector< double > TTQuantiSH::getSNPalleleFreqInPatch ( Patch patch)
2594{
2595 Individual* ind;
2596 double **genes;
2597 int total_size = 0;
2598
2599 if(!_SHLinkedTrait) fatal("TTQuantiSH::getSNPalleleFreqInPatch:: Pointer to trait prototype is not set!\n");
2600
2601 double **snp_allele = _SHLinkedTrait->get_allele_values();
2602
2603
2604 unsigned int loc;
2605 unsigned int nb_trait = _SHLinkedTrait->get_nb_traits();
2606 unsigned int nb_locus = _SHLinkedTrait->get_nb_locus();
2607
2608
2609
2610 vector<double> snp_tab(nb_trait * nb_locus, 0.0);
2611
2612 //ONLY IMPLEMENTED FOR DIALLELIC LOCI!
2613 for(unsigned int k = 0; k < nb_trait; k++) {
2614 for(unsigned int l = 0; l < nb_locus; l++) {
2615 loc = l * nb_trait + k;
2616 snp_tab [loc] = 0.0;
2617 }
2618 }
2619
2620 total_size = patch->size(FEM, ADLTx) + patch->size(MAL, ADLTx);
2621
2622
2623 for(unsigned int j = 0, size = patch->size(FEM, ADLTx); j < size; j++) {
2624
2625 ind = patch->get(FEM, ADLTx, j);
2626
2627 genes = (double**)ind->getTrait(_SHLinkedTraitIndex)->get_sequence();
2628
2629
2630 for(unsigned int k = 0; k < nb_trait; k++) {
2631 for(unsigned int l = 0; l < nb_locus; l++) {
2632
2633 loc = l * nb_trait + k;
2634
2635 if (genes[0][loc] == snp_allele[l][0]){
2636 snp_tab [loc] += 1;
2637 }
2638
2639 if (genes[1][loc] == snp_allele[l][0]){
2640 snp_tab [loc] += 1;
2641 }
2642 }
2643 }
2644 }
2645
2646
2647 for(unsigned int j = 0, size = patch->size(MAL, ADLTx); j < size; j++) {
2648
2649 ind = patch->get(MAL, ADLTx, j);
2650
2651 genes = (double**)ind->getTrait(_SHLinkedTraitIndex)->get_sequence();
2652
2653 for(unsigned int k = 0; k < nb_trait; k++) {
2654 for(unsigned int l = 0; l < nb_locus; l++) {
2655
2656 loc = l * nb_trait + k;
2657
2658 if (genes[0][loc] == snp_allele[l][0]){
2659 snp_tab [loc] += 1;
2660 }
2661 if (genes[1][loc] == snp_allele[l][0]){
2662 snp_tab [loc] += 1;
2663 }
2664 }
2665 }
2666
2667
2668 }
2669
2670 double inv_size = 1.0/(2.0*total_size);
2671
2672 for(unsigned int k = 0; k < nb_trait; k++) {
2673 for(unsigned int l = 0; l < nb_locus; l++) {
2674 loc = l * nb_trait + k;
2675 snp_tab [loc] *= inv_size;
2676 }
2677 }
2678
2679 return snp_tab;
2680}
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:496
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:532
unsigned int get_nb_locus()
Definition: ttquanti.h:122
unsigned int get_nb_traits()
Definition: ttquanti.h:121
double ** get_allele_values() const
Definition: ttquanti.h:131
virtual void ** get_sequence() const =0
sequence accessor.
int _SHLinkedTraitIndex
Index of the trait in the Individual::Traits table.
Definition: stathandler.h:173
void fatal(const char *str,...)
Definition: output.cc:96
@ FEM
Definition: types.h:37
@ MAL
Definition: types.h:37
@ ADLTx
Definition: types.h:42

References TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, ADLTx, fatal(), FEM, Patch::get(), TProtoQuanti::get_allele_values(), TProtoQuanti::get_nb_locus(), TProtoQuanti::get_nb_traits(), TTrait::get_sequence(), Individual::getTrait(), MAL, and Patch::size().

Referenced by getVaNoDominance(), and getVaWithDominance().

+ Here is the caller graph for this function:

◆ getVa()

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

References _Va.

Referenced by addQuanti().

+ Here is the caller graph for this function:

◆ getVaNoDominance()

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

CYCLE THROUGH ALL TRAITS

2685{
2686 unsigned int sizeF = patch->size(FEM, AGE),
2687 sizeM = patch->size(MAL, AGE);
2688 unsigned int size = sizeF + sizeM;
2689
2690 if(!size) {return vector<double>(_nb_trait, -1.0);}
2691
2692 TTQuanti* trait;
2693 double **allele_value = _SHLinkedTrait->get_allele_values(); // nlocus x 2
2694
2695 double G, meanG;
2696
2697 vector< double > varA(_nb_trait, 0.0);
2698
2699 double *genot = new double[size];
2700
2702 for(unsigned int TRAIT = 0; TRAIT < _nb_trait; ++TRAIT)
2703 {
2704
2705 meanG = 0;
2706
2707 // females
2708
2709 Individual* curInd;
2710 unsigned int gpos = 0;
2711
2712 for(unsigned int i = 0; i < sizeF && gpos < size; ++i )
2713 {
2714 curInd = patch->get(FEM, AGE, i);
2715
2716 trait = dynamic_cast<TTQuanti*>(curInd->getTrait(_SHLinkedTraitIndex));
2717
2718 G = trait->get_additive_genotype(TRAIT); // genotype value
2719
2720 genot[gpos++] = G;
2721
2722 meanG += G;
2723
2724 } // for each female
2725
2726 // males
2727 for(unsigned int i = 0; i < sizeM && gpos < size; ++i )
2728 {
2729 curInd = patch->get(MAL, AGE, i);
2730
2731 trait = dynamic_cast<TTQuanti*>(curInd->getTrait(_SHLinkedTraitIndex));
2732
2733 G = trait->get_additive_genotype(TRAIT); // genotype value
2734
2735 genot[gpos++] = G;
2736
2737 meanG += G;
2738
2739 } // for each male
2740
2741 assert(gpos == size);
2742
2743 meanG /= size;
2744
2745 varA[TRAIT] = my_variance_with_fixed_mean(genot, size, meanG); // not a sample variance
2746
2747 // record the genotypic variance in the stater
2748 _Vg[TRAIT] = varA[TRAIT];
2749
2750#ifdef _DEBUG_
2751 message("TTQuantiSH::getVaNoDominance:: Va = %d\n",varA[TRAIT]);
2752
2753 double Vacheck = 0, a, d,avG;
2754 //compute the allele frequencies:
2755 vector< double > freqs = getSNPalleleFreqInPatch(patch); //it has num locus x num traits elements, in a 1D array
2756
2757 for(unsigned int l = 0, pos; l < _nb_locus; ++l)
2758 {
2759 pos = l * _nb_trait + TRAIT;
2760 a = abs(allele_value[l][0] - allele_value[l][1])/2;
2761 d = allele_value[l][0]+allele_value[l][1];
2762 avG = a + d*(1-2*freqs[pos]);
2763 Vacheck += 2*freqs[pos]*(1 - freqs[pos])*avG*avG;
2764 }
2765
2766 message("TTQuantiSH::getVaNoDominance:: Va check = %d\n",Vacheck);
2767#endif
2768
2769 }
2770
2771 delete [] genot;
2772
2773 return varA;
2774}
vector< double > getSNPalleleFreqInPatch(Patch *patch)
Definition: ttquanti.cc:2593
TTQuanti.
Definition: ttquanti.h:56
double get_additive_genotype(unsigned int trait)
Definition: ttquanti.cc:1664
void message(const char *message,...)
Definition: output.cc:40
double my_variance_with_fixed_mean(double *data, unsigned int size, double mean)
Definition: utils.cc:62

References _nb_locus, _nb_trait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, _Vg, FEM, Patch::get(), TTQuanti::get_additive_genotype(), TProtoQuanti::get_allele_values(), getSNPalleleFreqInPatch(), Individual::getTrait(), MAL, message(), my_variance_with_fixed_mean(), and Patch::size().

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

+ Here is the caller graph for this function:

◆ getVaPerPatch()

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

References _pVa.

Referenced by addVarPerPatch().

+ Here is the caller graph for this function:

◆ 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

2783{
2784 unsigned int sizeF = patch->size(FEM, AGE),
2785 sizeM = patch->size(MAL, AGE);
2786 unsigned int size = sizeF + sizeM;
2787
2788#ifdef _DEBUG_
2789 message("\nTTQuantiSH::getVaWithDominance in patch %i; size = %i + %i\n", patch->getID(), sizeF, sizeM);
2790#endif
2791
2792 if(!size)
2793 return vector<double>(_nb_trait, -1.0);
2794
2795
2796 //compute the allele frequencies:
2797 vector< double > freqs = getSNPalleleFreqInPatch(patch); //it has num locus x num traits elements, in a 1D array
2798
2799 unsigned int l, pos, a1, a2;
2800 double** genes;
2801 TTQuanti* trait;
2802 double **allele_value = _SHLinkedTrait->get_allele_values(); // nlocus x 2
2803
2804 double G, meanG, varG;
2805 vector< double > genotypeValues[3];
2806 vector< double > genotypeFreq[3];
2807
2808 double *genot = new double[size];
2809
2810 vector< double > varA(_nb_trait, 0.0);
2811
2813 for(unsigned int TRAIT = 0; TRAIT < _nb_trait; ++TRAIT)
2814 {
2815
2816 meanG = 0;
2817
2818 genotypeValues[0].assign(_nb_locus, 0.0); //AA
2819 genotypeValues[1].assign(_nb_locus, 0.0); //Aa/aA
2820 genotypeValues[2].assign(_nb_locus, 0.0); //aa
2821
2822 genotypeFreq[0].assign(_nb_locus, 0.0); //AA
2823 genotypeFreq[1].assign(_nb_locus, 0.0); //Aa/aA
2824 genotypeFreq[2].assign(_nb_locus, 0.0); //aa
2825
2826 for(unsigned int l = 0; l < _nb_locus; ++l)
2827 {
2828 genotypeValues[0][l] = _SHLinkedTrait->get_genotype( allele_value[l][0], allele_value[l][0], l, TRAIT);
2829 genotypeValues[1][l] = _SHLinkedTrait->get_genotype( allele_value[l][0], allele_value[l][1], l, TRAIT);
2830 genotypeValues[2][l] = _SHLinkedTrait->get_genotype( allele_value[l][1], allele_value[l][1], l, TRAIT);
2831 }
2832
2833 // females
2834
2835 Individual* curInd;
2836 unsigned int gpos = 0, fpos = 0;
2837
2838 for(unsigned int i = 0; i < sizeF && gpos < size; ++i )
2839 {
2840 curInd = patch->get(FEM, AGE, i);
2841
2842 trait = dynamic_cast<TTQuanti*>(curInd->getTrait(_SHLinkedTraitIndex));
2843
2844 G = trait->get_full_genotype(TRAIT); // genotype value
2845
2846 meanG += G;
2847
2848 genot[gpos++] = G;
2849
2850 genes = (double**)trait->get_sequence();
2851
2852 for(unsigned int l = 0; l < _nb_locus; ++l)
2853 {
2854 pos = l * _nb_trait + TRAIT;
2855
2856 a1 = (genes[0][pos] != allele_value[l][0]); //allele index 1.allele
2857 a2 = (genes[1][pos] != allele_value[l][0]);
2858
2859 fpos = a1 + a2; //which genotype: 0 = AA, 1 = Aa/aA, 2 = aa
2860
2861 genotypeFreq[fpos][l]++;
2862
2863 if( genotypeValues[fpos][l] != _SHLinkedTrait->get_genotype( genes[0][pos], genes[1][pos], l, TRAIT))
2864 {
2865 cout << "assert failure: genot table value: "<<genotypeValues[fpos][l]<<" at "<<fpos
2866 <<"; genot value: "<<_SHLinkedTrait->get_genotype( genes[0][pos], genes[1][pos], l, TRAIT)
2867 <<" ("<< genes[0][pos]<<", "<< genes[1][pos] <<", "<<_SHLinkedTrait->get_dominance(l, TRAIT)<<")\n";
2868 }
2869 }
2870
2871
2872 } // for each female
2873
2874 // males
2875 for(unsigned int i = 0; i < sizeM && gpos < size; ++i )
2876 {
2877 curInd = patch->get(MAL, AGE, i);
2878
2879 trait = dynamic_cast<TTQuanti*>(curInd->getTrait(_SHLinkedTraitIndex));
2880
2881 G = trait->get_full_genotype(TRAIT); // genotype value
2882
2883 meanG += G;
2884
2885 genot[gpos++] = G;
2886
2887 genes = (double**)trait->get_sequence();
2888
2889 for(unsigned int l = 0; l < _nb_locus; ++l)
2890 {
2891 pos = l * _nb_trait + TRAIT;
2892
2893 a1 = (genes[0][pos] != allele_value[l][0]); //allele index 1.allele
2894 a2 = (genes[1][pos] != allele_value[l][0]);
2895
2896 fpos = a1 + a2; //which genotype: 0 = AA, 1 = Aa/aA, 2 = aa
2897
2898 genotypeFreq[fpos][l]++;
2899 }
2900
2901 } // for each male
2902
2903 assert(gpos == size);
2904
2905 meanG /= size;
2906
2907 // record the genotypic variance in the stater, will be needed to comput Ve
2908 _Vg[TRAIT] = my_variance_with_fixed_mean(genot, size, meanG);
2909
2910 for(unsigned int l = 0; l < _nb_locus; ++l)
2911 {
2912 genotypeFreq[0][l] /= size;
2913 genotypeFreq[1][l] /= size;
2914 genotypeFreq[2][l] /= size;
2915 }
2916
2917
2918 // compute the average excess (alphaStar) which is identical to the additive effects under random mating
2919 // !!! THIS DOESN'T WORK AT ALL !!!
2920// vector< double > alphaStar[2];
2921// alphaStar[0].assign(_nb_locus, 0.0);
2922// alphaStar[1].assign(_nb_locus, 0.0);
2923//
2924// double mean_alpha = 0;
2925//
2926// for(unsigned int l = 0; l < _nb_locus; ++l)
2927// {
2928//
2929// pos = l * _nb_trait + TRAIT;
2930//
2931// // average excess of each allele; the frequency table contains the frequency of the first allele only
2932// // here, we condition on the genotype being present in the population, in case HW equilibrium is not reached
2933// // eq 4.13b Lynch & Walsh 98
2934// // allele frequency freq[pos] is that of the first allele at [0]
2935// alphaStar[0][l] = genotypeValues[1][l] * (1.0 - freqs[pos]) + genotypeValues[0][l] * freqs[pos] - meanG; // alpha*_1
2936//
2937// alphaStar[1][l] = genotypeValues[1][l] * freqs[pos] + genotypeValues[2][l] * (1.0-freqs[pos]) - meanG; // alpha*_2
2938//
2942//
2943// mean_alpha += (alphaStar[0][l] + alphaStar[1][l])/2;
2944//
2945// 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)
2946//
2947// }
2948// varA[TRAIT] *= 2; // for diploids
2949//
2950// cout <<"TTQuantiSH::getVaWithDominance:: Va = "<<varA[TRAIT]<< " (mean alpha = "<<mean_alpha/_nb_locus<<")"<<endl;
2951
2952 // 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]
2953 double Vacheck = 0, a, d,avG, VDcheck = 0, vd;
2954 for(unsigned int l = 0; l < _nb_locus; ++l)
2955 {
2956 pos = l * _nb_trait + TRAIT;
2957 a = abs(genotypeValues[0][l] - genotypeValues[2][l])/2.0;
2958 d = genotypeValues[1][l];
2959 avG = a + d*(1-2*freqs[pos]);
2960 Vacheck += 2*freqs[pos]*(1 - freqs[pos])*avG*avG;
2961 vd = 2*freqs[pos]*(1 - freqs[pos])*d;
2962 VDcheck += vd*vd;
2963 }
2964
2966 varA [TRAIT] = _Vg[TRAIT];
2967 else
2968 varA[TRAIT] = Vacheck;
2969
2970#ifdef _DEBUG_
2971 message("TTQuantiSH::getVaWithDominance:: mean genotype = %d\n", meanG);
2972 message("TTQuantiSH::getVaWithDominance:: var genotype Vg = %d\n",_Vg[TRAIT]);
2973 message("TTQuantiSH::getVaWithDominance:: estimates: Va = %d, Vd = %d, Vg = %d\n", Vacheck, VDcheck, Vacheck+VDcheck);
2974#endif
2975
2976 } //for each trait
2977
2978 delete [] genot;
2979
2980 return varA;
2981}
double get_dominance(unsigned int locus, unsigned int trait)
Definition: ttquanti.h:147
double get_genotype(const double a1, const double a2, const unsigned int locus, const unsigned int trait)
Definition: ttquanti.cc:1106
bool get_h2_isBroad()
Definition: ttquanti.h:127
virtual void ** get_sequence() const
Definition: ttquanti.h:79
double get_full_genotype(unsigned int trait)
Definition: ttquanti.cc:1671

References _nb_locus, _nb_trait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTrait, TraitStatHandler< TProtoQuanti, TTQuantiSH >::_SHLinkedTraitIndex, _Vg, FEM, Patch::get(), TProtoQuanti::get_allele_values(), TProtoQuanti::get_dominance(), TTQuanti::get_full_genotype(), TProtoQuanti::get_genotype(), TProtoQuanti::get_h2_isBroad(), TTQuanti::get_sequence(), Patch::getID(), getSNPalleleFreqInPatch(), Individual::getTrait(), MAL, message(), my_variance_with_fixed_mean(), and Patch::size().

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

+ Here is the caller graph for this function:

◆ getVb()

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

References _Vb.

Referenced by addQuanti().

+ Here is the caller graph for this function:

◆ getVg()

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

References _Vg.

Referenced by LCE_QuantiModifier::setVefromVa().

+ Here is the caller graph for this function:

◆ getVp()

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

References _Vp.

Referenced by addQuanti().

+ Here is the caller graph for this function:

◆ getVpPerPatch()

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

References _pVp.

Referenced by addVarPerPatch().

+ Here is the caller graph for this function:

◆ init()

void TTQuantiSH::init ( )
virtual

Reimplemented from StatHandlerBase.

1795{
1797
1799
1800 _eVar = (!_SHLinkedTrait->get_env_var().empty() || !_SHLinkedTrait->get_heritability().empty());
1801
1802 resetPtrs(); //deallocate anything that was previously allocated
1803
1804 if(_patchNbr != _pop->getPatchNbr())
1806
1809
1812
1813 _G = gsl_matrix_alloc(_nb_trait,_nb_trait);
1814 _eval = gsl_vector_alloc (_nb_trait);
1815 _evec = gsl_matrix_alloc (_nb_trait, _nb_trait);
1816 _ws = gsl_eigen_symmv_alloc (_nb_trait);
1817
1818 _meanP = new double [_nb_trait];
1819 _meanG = new double [_nb_trait];
1820 _Va = new double [_nb_trait];
1821 _Vg = new double [_nb_trait];
1822 _Vb = new double [_nb_trait];
1823 _Vp = new double [_nb_trait];
1824 _eigval = new double [_nb_trait];
1825
1826 _eigvect = new double* [_nb_trait];
1827 _pVa = new double* [_nb_trait];
1828 _pVp = new double* [_nb_trait];
1829 _pmeanP = new double* [_nb_trait];
1830 _pmeanG = new double* [_nb_trait];
1831 _peigval = new double* [_nb_trait];
1832
1833 for(unsigned int i=0; i < _nb_trait; ++i) {
1834
1835 _eigvect[i] = new double [_nb_trait];
1836
1837 _pVa[i] = new double [_patchNbr];
1838 _pVp[i] = new double [_patchNbr];
1839 _pmeanP[i] = new double [_patchNbr];
1840 _pmeanG[i] = new double [_patchNbr];
1841 _peigval[i] = new double [_patchNbr];
1842
1843 }
1844
1845 _peigvect = new double* [_patchNbr];
1846 for(unsigned int i = 0; i < _patchNbr; i++)
1847 _peigvect[i] = new double [_nb_trait*_nb_trait];
1848
1849
1850 if(_nb_trait > 1) {
1851
1852 _covar = new double [_nb_trait*(_nb_trait -1)/2];
1853
1854 _pcovar = new double* [_nb_trait*(_nb_trait - 1)/2];
1855
1856 for(unsigned int i = 0; i < _nb_trait*(_nb_trait - 1)/2; i++)
1857 _pcovar[i] = new double [_patchNbr];
1858
1859 }
1860
1861}
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_env_var()
Definition: ttquanti.h:124
vector< double > get_heritability()
Definition: ttquanti.h:125
virtual trait_t get_type() const
Definition: ttquanti.h:179

References _covar, _eigval, _eigvect, _eval, _eVar, _evec, _G, _meanG, _meanP, _nb_locus, _nb_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_nb_locus(), TProtoQuanti::get_nb_traits(), TProtoQuanti::get_type(), Metapop::getPatchNbr(), IndFactory::getTraitIndex(), StatHandlerBase::init(), and resetPtrs().

◆ resetPtrs()

void TTQuantiSH::resetPtrs ( )
1735{
1736
1737 if(_G != NULL) gsl_matrix_free(_G);
1738 if(_eval != NULL) gsl_vector_free(_eval);
1739 if(_evec != NULL) gsl_matrix_free(_evec);
1740 if(_ws != NULL) gsl_eigen_symmv_free (_ws);
1741
1742 if(_meanP != NULL) delete [] _meanP;
1743 if(_meanG != NULL) delete [] _meanG;
1744 if(_Va != NULL) delete [] _Va;
1745 if(_Vg != NULL) delete [] _Vg;
1746 if(_Vb != NULL) delete [] _Vb;
1747 if(_Vp != NULL) delete [] _Vp;
1748 if(_covar != NULL) delete [] _covar;
1749 _covar = 0;
1750 if(_eigval != NULL) delete [] _eigval;
1751
1752
1753 if(_eigvect) {
1754 for(unsigned int i=0; i < _nb_trait; ++i) delete [] _eigvect[i];
1755 delete [] _eigvect;
1756 }
1757 if(_pVa) {
1758 for(unsigned int i=0; i < _nb_trait; ++i) delete [] _pVa[i];
1759 delete [] _pVa;
1760 }
1761 if(_pVp) {
1762 for(unsigned int i=0; i < _nb_trait; ++i) delete [] _pVp[i];
1763 delete [] _pVp;
1764 }
1765 if(_pmeanP) {
1766 for(unsigned int i=0; i < _nb_trait; ++i) delete [] _pmeanP[i];
1767 delete [] _pmeanP;
1768 }
1769 if(_pmeanG) {
1770 for(unsigned int i=0; i < _nb_trait; ++i) delete [] _pmeanG[i];
1771 delete [] _pmeanG;
1772 }
1773 if(_peigval) {
1774 for(unsigned int i=0; i < _nb_trait; ++i) delete [] _peigval[i];
1775 delete [] _peigval;
1776 }
1777
1778 if(_pcovar != NULL) {
1779 for(unsigned int i = 0; i < (_nb_trait*(_nb_trait - 1)/2); i++) delete [] _pcovar[i];
1780 delete [] _pcovar;
1781 _pcovar = 0;
1782 }
1783
1784 if(_peigvect != NULL) {
1785 for(unsigned int i = 0; i < _patchNbr; i++) delete [] _peigvect[i];
1786 delete [] _peigvect;
1787 _peigvect = 0;
1788 }
1789
1790}

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

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

+ Here is the caller graph for this function:

◆ setAdultStats()

void TTQuantiSH::setAdultStats ( )
inline
void setStats(age_t AGE)
Definition: ttquanti.cc:2456

References ADULTS, and setStats().

Referenced by addAvgPerPatch(), addCovarPerPatch(), addEigen(), addEigenPerPatch(), addEigenValues(), addEigenValuesPerPatch(), addEigenVect1(), addEigenVect1PerPatch(), addQuanti(), addSkewPerPatch(), and addVarPerPatch().

+ Here is the caller graph for this function:

◆ setDataTables()

void TTQuantiSH::setDataTables ( age_t  AGE)
2387{
2388 unsigned int **sizes;
2389 unsigned int nb_patch = _pop->getPatchNbr();
2390
2391 sizes = new unsigned int * [_nb_trait];
2392
2393 for(unsigned int i = 0; i < _nb_trait; ++i) {
2394 sizes[i] = new unsigned int [nb_patch];
2395 for(unsigned int j = 0; j < nb_patch; ++j)
2396 sizes[i][j] = _pop->size(AGE, j);
2397 }
2398
2399 _phenoTable.update(_nb_trait, nb_patch, sizes);
2400 _genoTable.update(_nb_trait, nb_patch, sizes);
2401
2402 for(unsigned int i = 0; i < _nb_trait; ++i)
2403 delete [] sizes[i];
2404 delete [] sizes;
2405
2406 Patch* patch;
2407 age_idx age = (AGE == ADULTS ? ADLTx : OFFSx);
2408
2409 for(unsigned int i = 0, n; i < nb_patch; i++) {
2410
2411 patch = _pop->getPatch(i);
2412
2413 n=0;
2414
2415 if ((patch->size(MAL, age)+patch->size(FEM, age)) != _phenoTable.size(0,i)) {
2416 fatal("problem while recording quanti trait values; table size doesn't match patch size.\n");
2417 }
2418 store_quanti_trait_values(patch, i, patch->size(MAL, age), &n, MAL, age, &_phenoTable, &_genoTable,
2420
2421 store_quanti_trait_values(patch, i, patch->size(FEM, age), &n, FEM, age, &_phenoTable, &_genoTable,
2423
2424 if (n != _phenoTable.size(0,i) || n != _genoTable.size(0,i)) {
2425 fatal("problem while recording quanti trait values; size counter doesn't match table size.\n");
2426 }
2427 }
2428}
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:310
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:430
DataTable< double > _genoTable
Definition: ttquanti.h:269
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:2432
age_idx
Array index of the age classes in the patch sizes and containers arrays.
Definition: types.h:41
@ OFFSx
Definition: types.h:42

References _genoTable, _nb_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().

+ Here is the caller graph for this function:

◆ setOffsprgStats()

void TTQuantiSH::setOffsprgStats ( )
inline

◆ setStatRecorders()

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

Implements StatHandlerBase.

1866{
1867#ifdef _DEBUG_
1868 message("-TTQuantiSH::setStatRecorders ");
1869#endif
1870 string age_tag = token.substr(0,token.find_first_of("."));
1871 string sub_token;
1872 age_t AGE = ALL;
1873
1874 if (age_tag.size() != 0 && age_tag.size() != string::npos) {
1875
1876 if (age_tag == "adlt") AGE = ADULTS;
1877
1878 else if (age_tag == "off") AGE = OFFSPRG;
1879
1880 else age_tag = "";
1881
1882 } else {
1883 age_tag = "";
1884 }
1885
1886 if (age_tag.size() != 0)
1887 sub_token = token.substr(token.find_first_of(".") + 1, string::npos);
1888 else
1889 sub_token = token;
1890
1891 if(sub_token == "quanti") {
1892 addQuanti(AGE);
1893 } else if(sub_token == "quanti.eigen") { //based on Vb; among-patch (D) matrix
1894 addEigen(AGE);
1895 } else if(sub_token == "quanti.eigenvalues") { //based on Vb; among-patch (D) matrix
1896 addEigenValues(AGE);
1897 } else if(sub_token == "quanti.eigenvect1") { //based on Vb; among-patch (D) matrix
1898 addEigenVect1(AGE);
1899 } else if(sub_token == "quanti.patch") {
1900 addQuantiPerPatch(AGE);
1901 } else if(sub_token == "quanti.mean.patch") {
1902 addAvgPerPatch(AGE);
1903 } else if(sub_token == "quanti.var.patch") {
1904 addVarPerPatch(AGE);
1905 } else if(sub_token == "quanti.covar.patch") {
1906 addCovarPerPatch(AGE);
1907 } else if(sub_token == "quanti.eigen.patch") {
1908 addEigenPerPatch(AGE);
1909 } else if(sub_token == "quanti.eigenvalues.patch") {
1911 } else if(sub_token == "quanti.eigenvect1.patch") {
1913 } else if(sub_token == "quanti.skew.patch") {
1914 addSkewPerPatch(AGE);
1915 } else
1916 return false;
1917
1918 return true;
1919}
void addEigenVect1PerPatch(age_t AGE)
Definition: ttquanti.cc:2296
void addEigenValuesPerPatch(age_t AGE)
Definition: ttquanti.cc:2262
void addEigenVect1(age_t AGE)
Definition: ttquanti.cc:2051
void addEigenValues(age_t AGE)
Definition: ttquanti.cc:2023
unsigned int age_t
Age class flags.
Definition: types.h:46

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

◆ setStats()

void TTQuantiSH::setStats ( age_t  AGE)
2457{
2458 if(_table_set_age == AGE
2461 return;
2462
2463 unsigned int pop_size = _pop->size(AGE);
2464 unsigned int patch_size;
2465 double *phenot1, *genot1, *genot2;
2466
2467 unsigned int nb_patch = _pop->getPatchNbr();
2468
2469 if(nb_patch < _patchNbr) {
2470 warning("increase in patch number detected (in Quanti Stat Handler),");
2471 warning("stats for quanti trait will not be recorded in new patches, patch identity may have changed.\n");
2472 _patchNbr = nb_patch; // record stat only in remaining patches
2473 }
2474
2475 if(nb_patch > _patchNbr) nb_patch = _patchNbr; //tables have been allocated for _patchNbr patches
2476
2477 setDataTables(AGE);
2478
2479#ifdef HAS_GSL
2480 unsigned int c = 0; //covariance position counter
2481 unsigned int pv = 0; //eigenvector position counter
2482 //within deme stats:
2483 for(unsigned int j = 0; j < nb_patch; j++) {
2484
2485 patch_size = _pop->size(AGE, j);
2486
2487 for(unsigned int t=0; t < _nb_trait; t++) {
2488
2489 phenot1 = _phenoTable.getClassWithinGroup(t,j);
2490 genot1 = _genoTable.getClassWithinGroup(t,j);
2491
2492 _pmeanP[t][j] = my_mean (phenot1, patch_size );
2493 _pmeanG[t][j] = my_mean (genot1, patch_size );
2494 _pVp[t][j] = my_variance_with_fixed_mean (phenot1, patch_size, _pmeanP[t][j]);
2495 _pVa[t][j] = my_variance_with_fixed_mean (genot1, patch_size, _pmeanG[t][j]);
2496
2497 }
2498
2499 if(_nb_trait > 1) {
2500 c = 0;
2501 // calculate the covariances and G, need to adjust dimensions in class declaration
2502 for(unsigned int t1 = 0; t1 < _nb_trait; t1++) {
2503 // set the diagonal elements of G here
2504
2505 gsl_matrix_set(_G, t1, t1, _pVa[t1][j]);
2506
2507 for(unsigned int t2 = t1 + 1; t2 < _nb_trait; t2++) {
2508
2509 genot1 = _genoTable.getClassWithinGroup(t1,j);
2510 genot2 = _genoTable.getClassWithinGroup(t2,j);
2511
2512 _pcovar[c][j] = gsl_stats_covariance_m (genot1, 1, genot2, 1, patch_size,
2513 _pmeanG[t1][j], _pmeanG[t2][j]);
2514
2515 gsl_matrix_set(_G, t1, t2, _pcovar[c][j]);
2516 gsl_matrix_set(_G, t2, t1, _pcovar[c++][j]);
2517 }
2518 }
2519
2520 gsl_eigen_symmv (_G, _eval, _evec, _ws);
2521 gsl_eigen_symmv_sort (_eval, _evec, GSL_EIGEN_SORT_VAL_DESC);
2522
2523 pv = 0;
2524
2525 for(unsigned int t = 0; t < _nb_trait; t++) {
2526 _peigval[t][j] = gsl_vector_get (_eval, t);
2527 for(unsigned int v = 0; v < _nb_trait; v++)
2528 _peigvect[j][pv++] = gsl_matrix_get (_evec, v, t); //read eigenvectors column-wise
2529 }
2530 }
2531 }
2532
2533 double meanGamong1, meanGamong2;
2534 c = 0; //reset covariance positioner
2535 //among demes stats:
2536 for(unsigned int t1 = 0; t1 < _nb_trait; t1++) {
2537
2538 phenot1 = _phenoTable.getGroup(t1);
2539 genot1 = _genoTable.getGroup(t1);
2540
2541 _meanP[t1] = my_mean (phenot1, pop_size ); //grand mean, pooling all individuals
2542 _meanG[t1] = my_mean (genot1, pop_size ); //grand mean, pooling all individuals
2543
2544 _Vp[t1] = my_mean_no_nan (_pVp[t1], nb_patch); //mean within patch variance
2545 _Va[t1] = my_mean_no_nan (_pVa[t1], nb_patch); //mean within patch variance
2546
2547 meanGamong1 = my_mean_no_nan (_pmeanG[t1], nb_patch); //mean of within patch mean genotypic values
2548
2549 _Vb[t1] = my_variance_with_fixed_mean_no_nan (_pmeanG[t1], nb_patch, meanGamong1); //variance of patch means
2550
2551 gsl_matrix_set(_G, t1, t1, _Vb[t1]); //_G here becomes the D-matrix, the among-deme (Difference) covariance matrix
2552
2553 for(unsigned int t2 = t1 + 1; t2 < _nb_trait; t2++) {
2554
2555 meanGamong2 = my_mean (_pmeanG[t2], nb_patch);
2556
2557 _covar[c] = my_covariance_no_nan(_pmeanG[t1], _pmeanG[t2], nb_patch, nb_patch);
2558
2559 gsl_matrix_set(_G, t1, t2, _covar[c]);
2560 gsl_matrix_set(_G, t2, t1, _covar[c]);
2561
2562 //set covariance to the mean covariance within patch for recording
2563 _covar[c] = my_mean_no_nan (_pcovar[c], nb_patch);
2564
2565 c++;
2566 }
2567 }
2568
2569 if(_nb_trait > 1) {
2570 gsl_eigen_symmv (_G, _eval, _evec, _ws);
2571 gsl_eigen_symmv_sort (_eval, _evec, GSL_EIGEN_SORT_VAL_DESC);
2572
2573 for(unsigned int t1 = 0; t1 < _nb_trait; t1++) {
2574 _eigval[t1] = gsl_vector_get (_eval, t1);
2575 for(unsigned int t2 = 0; t2 < _nb_trait; t2++) {
2576 _eigvect[t2][t1] = gsl_matrix_get (_evec, t2, t1);
2577 }
2578 }
2579 }
2580
2581#else
2582 fatal("install the GSL library to get the quanti stats!\n");
2583#endif
2584
2585 _table_set_age = AGE;
2588
2589}
T * getGroup(unsigned int group)
Accessor to a group array.
Definition: datatable.h:223
unsigned int getCurrentReplicate()
Definition: metapop.h:293
unsigned int getCurrentGeneration()
Definition: metapop.h:294
void setDataTables(age_t AGE)
Definition: ttquanti.cc:2386
double my_variance_with_fixed_mean_no_nan(double *data, unsigned int size, double mean)
Definition: utils.cc:73
double my_mean(double *data, unsigned int size)
Definition: utils.cc:37
double my_mean_no_nan(double *data, unsigned int size)
Definition: utils.cc:48
double my_covariance_no_nan(double *data1, double *data2, unsigned int size1, unsigned int size2)
Definition: utils.cc:87

References _covar, _eigval, _eigvect, _eval, _evec, _G, _genoTable, _meanG, _meanP, _nb_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().

+ Here is the caller graph for this function:

Member Data Documentation

◆ _covar

double * TTQuantiSH::_covar
private

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

◆ _eigval

double * TTQuantiSH::_eigval
private

◆ _eigvect

double ** TTQuantiSH::_eigvect
private

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

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

◆ _meanP

double* TTQuantiSH::_meanP
private

◆ _nb_locus

unsigned int TTQuantiSH::_nb_locus
private

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

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

◆ _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.3.56 by  doxygen 1.9.0 -- Nemo is hosted on  Download Nemo

Locations of visitors to this page
Catalogued on GSR