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

The stat handler for neutral markers. More...

#include <ttneutralgenes.h>

+ Inheritance diagram for TTNeutralGenesSH:
+ Collaboration diagram for TTNeutralGenesSH:

Public Member Functions

 TTNeutralGenesSH (TProtoNeutralGenes *TP)
 
virtual ~TTNeutralGenesSH ()
 
virtual void init ()
 
virtual bool setStatRecorders (std::string &token)
 
void setFreqRecorders (age_t AGE)
 
void setFreqRecordersPerPatch (age_t AGE)
 
void setFstatRecorders (age_t AGE)
 
void setFstat2Recorders (age_t AGE)
 
void setFstatWCRecorders (age_t AGE)
 
void setCoaMatrixRecorders (age_t AGE, unsigned char dim)
 
void setFstMatrixRecorders (age_t AGE, unsigned char dim)
 
void setNeiGeneticDistanceRecorders (age_t AGE, bool pairwise)
 
void setDxyRecorders (age_t AGE, bool patchwise)
 
Allele and genotype frequencies:
void setAdultAlleleFreq ()
 
void setOffspringAlleleFreq ()
 
void setHeterozygosity (age_t AGE)
 
void setAdultHeterozygosity ()
 
void setOffspringHeterozygosity ()
 
double getGlobalAlleleFreq (unsigned int loc, unsigned int all)
 
double getHeterozygosity (unsigned int loc)
 
F-stats:
void setAlleleTables (age_t AGE)
 
void setHeteroTable (age_t AGE)
 
void allocateTables (unsigned int loci, unsigned int all)
 
DataTable< double > * getAlleleFreqTable ()
 Accessor to the table of allele frequencies, per patch. More...
 
DataTable< unsigned int > * getAlleleCountTable ()
 
DataTable< double > * getHeteroTable ()
 
TMatrixgetGlobalFreqs ()
 Accessor to the table of allele frequencies in the whole population. More...
 
void setFstMatrix (age_t AGE, unsigned char dim)
 Computes the weighted within and between patch Fst's as well as the overall Fst (Theta). More...
 
void setAdultsFstMatrix ()
 
void setAdultsFstWithin ()
 
void setAdultsFstBetween ()
 
void setOffsprgFstMatrix ()
 
void setOffsprgFstWithin ()
 
void setOffsprgFstBetween ()
 
double getWeightedFst ()
 Returns the weighted Fst using Weir & Hill (2002) method. More...
 
double getFst_ij (unsigned int i)
 Accessor to the Fst matrix as set by setFstMatrix(). More...
 
void setFst_li (unsigned int N, unsigned int L, double **array)
 Computes the per-locus per-patch Fst values using Weir&Hill 2002 approach. More...
 
void setFstat (age_t AGE)
 Computes the F-statistics following Nei & Chesser (1983). More...
 
void setOffsprgFstat ()
 
void setAdultsFstat ()
 
double setHo (age_idx age_pos)
 
double setHs (age_idx age_pos)
 
double setHt (age_idx age_pos)
 
double getHsnei ()
 
double getHtnei ()
 
double getHo ()
 
double getHs ()
 
double getHt ()
 
double getFst ()
 
double getFis ()
 
double getFit ()
 
void setFstat2 (age_t AGE)
 New version of Nei & Chesser. More...
 
void setOffsprgFstat2 ()
 
void setAdultsFstat2 ()
 
deque< double > setHo2 (age_idx age_pos)
 
deque< double > setHs2 (age_idx age_pos)
 
deque< double > setHt2 (age_idx age_pos)
 
void setFstatWeirCockerham (age_t AGE)
 Computes the Weir & Cockerham (1984) Fstat values (Theta, F, and f). More...
 
void setFstatWeirCockerham_MS (age_t AGE)
 
void setOffspringFstatWeirCockerham ()
 
void setAdultsFstatWeirCockerham ()
 
double getFstWC ()
 
double getFisWC ()
 
double getFitWC ()
 
void setLociDivCounter (age_t AGE)
 Sets the allelic diversity counters. More...
 
double getNbAllLocal ()
 
double getNbAllGlobal ()
 
double getFixLocLocal ()
 
double getFixLocGlobal ()
 
Coancestries
double Coancestry (void **ind1, void **ind2, unsigned int nb_locus)
 Gives the coancestry (probability of identity by state) of two gene sequences. More...
 
void setCoaMatrix (age_idx age_pos, unsigned char dim)
 Computes the within and between patches coancestry coefficients. More...
 
void setAdultsCoaMatrix ()
 
void setOffsprgCoaMatrix ()
 
void setAdultsCoaWithin ()
 
void setOffsprgCoaWithin ()
 
void setAdultsCoaBetween ()
 
void setOffsprgCoaBetween ()
 
void setAdults_Theta ()
 
double getCoa (unsigned int i)
 Gets the given coancestry coefficient from the coancestry matrix. More...
 
double getMeanTheta ()
 
double getMeanAlpha ()
 
double getTheta_FF ()
 Gives the mean within females coancestry coefficient. More...
 
double getTheta_MM ()
 Gives the mean within males coancestry coefficient. More...
 
double getTheta_FM ()
 Gives the mean between males and females coancestry coefficient. More...
 
void setSibStats ()
 
void setSibCoa (Individual *I1, Individual *I2)
 
double getSibProportions (unsigned int i)
 
double getSibCoaMeans (unsigned int i)
 
Nei's genetic distance:
void setAdltNeiGeneticDistance ()
 
void setOffsprgNeiGeneticDistance ()
 
void setNeiGeneticDistance (age_t AGE)
 
double getNeiGeneticDistance (unsigned int i)
 
double getMeanNeiGeneticDistance ()
 
Sequence divergence: Dxy (Nei & Li 1979; Nei 1987)
double getDxyOffspringPerPatch (unsigned int patch1, unsigned patch2)
 
double getDxyAdultPerPatch (unsigned int patch1, unsigned patch2)
 
double getDxyPerPatch (age_idx age, unsigned int patch1, unsigned patch2)
 
double getDxy (unsigned int age_class)
 
- Public Member Functions inherited from TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >
 TraitStatHandler (TProtoNeutralGenes *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

DataTable< unsigned int_alleleCountTable
 
DataTable< double > _alleleFreqTable
 
DataTable< double > _heteroTable
 
TMatrix _globalAlleleFreq
 
unsigned int _table_set_gen
 
unsigned int _table_set_age
 
unsigned int _table_set_repl
 
double Theta_FF
 
double Theta_MM
 
double Theta_FM
 
double _mean_theta
 
double _mean_alpha
 
TMatrix_coa_matrix
 
double _sib_prop [4]
 Kinship classes proportions. More...
 
double _sib_coa [4]
 
double _ho
 F-statistics. More...
 
double _hs
 
double _ht
 
double _hsnei
 
double _htnei
 
double _nb_all_local
 
double _nb_all_global
 
double _fst
 
double _fis
 
double _fit
 
double _fix_loc_local
 
double _fix_loc_global
 
double _fst_WH
 Weir & Hill (2002) F-stat estimates. More...
 
double _fst_WC
 Weir & Cockerham (1984) F-stat estimates. More...
 
double _fis_WC
 
double _fit_WC
 
double * _fst_WC_loc
 Per-locus F-stats (Weir&Cockerham). More...
 
double * _fis_WC_loc
 
double * _fit_WC_loc
 
double _fst_W1
 
double _fst_W2
 
TMatrix_fst_matrix
 Pairwise Fst matrix. More...
 
TMatrix_D
 
double _meanD
 

Additional Inherited Members

- Protected Types inherited from StatHandler< SH >
typedef std::list< StatRecorder< SH > * >::iterator REC_IT
 
- Protected Attributes inherited from TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >
TProtoNeutralGenes_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

The stat handler for neutral markers.

Constructor & Destructor Documentation

◆ TTNeutralGenesSH()

TTNeutralGenesSH::TTNeutralGenesSH ( TProtoNeutralGenes TP)
inline
316 _D(0)
317 { }
double * _fis_WC_loc
Definition: ttneutralgenes.h:301
TMatrix * _fst_matrix
Pairwise Fst matrix.
Definition: ttneutralgenes.h:305
unsigned int _table_set_gen
Definition: ttneutralgenes.h:283
TMatrix * _D
Definition: ttneutralgenes.h:308
double * _fst_WC_loc
Per-locus F-stats (Weir&Cockerham).
Definition: ttneutralgenes.h:301
double * _fit_WC_loc
Definition: ttneutralgenes.h:301
TMatrix * _coa_matrix
Definition: ttneutralgenes.h:287
unsigned int _table_set_age
Definition: ttneutralgenes.h:283
unsigned int _table_set_repl
Definition: ttneutralgenes.h:283
Template class for the trait's StatHandler.
Definition: stathandler.h:168

◆ ~TTNeutralGenesSH()

virtual TTNeutralGenesSH::~TTNeutralGenesSH ( )
inlinevirtual
320 {
321 if(_coa_matrix != NULL) delete _coa_matrix;
322 if(_fst_matrix != NULL) delete _fst_matrix;
323 if(_fst_WC_loc) delete[]_fst_WC_loc;
324 if(_fis_WC_loc) delete[]_fis_WC_loc;
325 if(_fit_WC_loc) delete[]_fit_WC_loc;
326 if(_D != NULL) delete _D;
327 }

References _coa_matrix, _D, _fis_WC_loc, _fit_WC_loc, _fst_matrix, and _fst_WC_loc.

Member Function Documentation

◆ allocateTables()

void TTNeutralGenesSH::allocateTables ( unsigned int  loci,
unsigned int  all 
)
44{
45 unsigned int nb_patch = _pop->getPatchNbr();
46 unsigned int **sizes;
47
48 sizes = new unsigned int * [nb_patch];
49
50 for(unsigned int i = 0; i < nb_patch; ++i) {
51 sizes[i] = new unsigned int [loci];
52 for(unsigned int j = 0; j < loci; ++j)
53 sizes[i][j] = all;
54 }
55
56 _alleleCountTable.allocate(nb_patch, loci, sizes);
57
58 _alleleFreqTable.allocate(nb_patch, loci, sizes);
59
60 _globalAlleleFreq.reset(loci, all);
61
62 for(unsigned int i = 0; i < nb_patch; ++i)
63 delete [] sizes[i];
64 delete [] sizes;
65
66 //reset the time info, to force recalculation after allocate
70}
void allocate(unsigned int nbgroups, unsigned int nbclasses, unsigned int **classSizes)
Creates a table of size given by the sum of all classes passed by the 'classSizes' matrix.
Definition: datatable.h:107
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
void reset(unsigned int rows, unsigned int cols)
Re-allocate the existing matrix with assigned rows and cols dimensions.
Definition: tmatrix.h:116
DataTable< unsigned int > _alleleCountTable
Definition: ttneutralgenes.h:279
TMatrix _globalAlleleFreq
Definition: ttneutralgenes.h:282
DataTable< double > _alleleFreqTable
Definition: ttneutralgenes.h:280
#define NONE
No age flag.
Definition: types.h:48

References _alleleCountTable, _alleleFreqTable, _globalAlleleFreq, StatHandlerBase::_pop, _table_set_age, _table_set_gen, _table_set_repl, DataTable< T >::allocate(), Metapop::getPatchNbr(), NONE, and TMatrix::reset().

Referenced by init(), setAlleleTables(), setFstat(), setFstat2(), TTNeutralGenesFH::write_Fst_i(), and TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ Coancestry()

double TTNeutralGenesSH::Coancestry ( void **  ind1,
void **  ind2,
unsigned int  nb_locus 
)

Gives the coancestry (probability of identity by state) of two gene sequences.

The probability returned is the average probability of having two identical alleles at a locus between the two sequences.

Parameters
ind1first _sequence, treated as of type (unsigned char**)
ind2second _sequence, treated as of type (unsigned char**)
nb_locusnumber of loci present in each _sequence
44{
45 unsigned int p = 0;
46 unsigned char **seq1, **seq2;
47
48 seq1 = (unsigned char**)ind1;
49 seq2 = (unsigned char**)ind2;
50
51 for (unsigned int k = 0; k < nb_locus; ++k)
52 p += !(seq1[0][k]^seq2[0][k]) + !(seq1[0][k]^seq2[1][k]) + !(seq1[1][k]^seq2[0][k]) + !(seq1[1][k]^seq2[1][k]);
53
54 return (double)p/(4.0*nb_locus);
55}

Referenced by setAdults_Theta(), setCoaMatrix(), and setSibCoa().

+ Here is the caller graph for this function:

◆ getAlleleCountTable()

DataTable< unsigned int > * TTNeutralGenesSH::getAlleleCountTable ( )
inline
372{return &_alleleCountTable;}

References _alleleCountTable.

Referenced by TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ getAlleleFreqTable()

DataTable< double > * TTNeutralGenesSH::getAlleleFreqTable ( )
inline

Accessor to the table of allele frequencies, per patch.

370{return &_alleleFreqTable;}

References _alleleFreqTable.

Referenced by TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ getCoa()

double TTNeutralGenesSH::getCoa ( unsigned int  i)
inline

Gets the given coancestry coefficient from the coancestry matrix.

Parameters
icombination of the row and column indexes (see setCoaMatrixRecorders()).
Note
the upper half and the diagonal of the matrix are filled, other positions are set to 0.
478 {
479 unsigned int scale = (unsigned int)pow( 10.0, (int)log10((float)_coa_matrix->getNbCols()) + 1 );
480 return _coa_matrix->get(i/scale, i%scale);
481 }
double get(unsigned int i, unsigned int j)
Accessor to element at row i and column j.
Definition: tmatrix.h:147
unsigned int getNbCols()
Gives the number of columns.
Definition: tmatrix.h:169

References _coa_matrix, TMatrix::get(), and TMatrix::getNbCols().

Referenced by setCoaMatrixRecorders().

+ Here is the caller graph for this function:

◆ getDxy()

double TTNeutralGenesSH::getDxy ( unsigned int  age_class)
1211{
1212 double D = 0;
1213 unsigned int p = 0;
1214
1215 age_t AGE = (static_cast<age_idx>(age_class) == OFFSx ? OFFSPRG : ADULTS);
1216
1217 for (unsigned int p1 = 0; p1 < _pop->getPatchNbr(); ++p1) {
1218
1219 if( !_pop->size( AGE , p1) ) continue;
1220
1221 for (unsigned int p2 = p1 + 1; p2 < _pop->getPatchNbr(); ++p2) {
1222
1223 if( !_pop->size( AGE , p2) ) continue;
1224
1225 D += getDxyPerPatch(static_cast<age_idx>(age_class), p1, p2);
1226
1227 p++;
1228 }
1229 }
1230
1231 return (p != 0 ? D/p : nanf("NULL"));
1232}
unsigned int size()
Get the total number of individuals present in the population, all sex and age classes together.
Definition: metapop.h:310
double getDxyPerPatch(age_idx age, unsigned int patch1, unsigned patch2)
Definition: stats_fstat.cc:1236
unsigned int age_t
Age class flags.
Definition: types.h:46
#define ADULTS
Adults age class flag (breeders).
Definition: types.h:54
#define OFFSPRG
Offspring age class flag.
Definition: types.h:50
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 StatHandlerBase::_pop, ADULTS, getDxyPerPatch(), Metapop::getPatchNbr(), OFFSPRG, OFFSx, and Metapop::size().

Referenced by setDxyRecorders().

+ Here is the caller graph for this function:

◆ getDxyAdultPerPatch()

double TTNeutralGenesSH::getDxyAdultPerPatch ( unsigned int  patch1,
unsigned  patch2 
)
inline
512{return getDxyPerPatch(ADLTx, patch1, patch2);}
@ ADLTx
Definition: types.h:42

References ADLTx, and getDxyPerPatch().

Referenced by setDxyRecorders().

+ Here is the caller graph for this function:

◆ getDxyOffspringPerPatch()

double TTNeutralGenesSH::getDxyOffspringPerPatch ( unsigned int  patch1,
unsigned  patch2 
)
inline
511{return getDxyPerPatch(OFFSx, patch1, patch2);}

References getDxyPerPatch(), and OFFSx.

Referenced by setDxyRecorders().

+ Here is the caller graph for this function:

◆ getDxyPerPatch()

double TTNeutralGenesSH::getDxyPerPatch ( age_idx  age,
unsigned int  patch1,
unsigned  patch2 
)
1237{
1238 double D = 0;
1239 Patch* patch_1 = _pop->getPatchPtr(patch1);
1240 Patch* patch_2 = _pop->getPatchPtr(patch2);
1241 unsigned int size_1 = patch_1->size(FEM, age);
1242 unsigned int size_2 = patch_2->size(FEM, age);
1243 unsigned int N = 0;//4*size_1 * size_2;// (2N)^2, tot num of haplotype comparisons
1244
1245
1246 unsigned int num_loc = _SHLinkedTrait->get_locus_num();
1247
1248 unsigned char **seq_1, **seq_2;
1249
1250 //females:
1251 for (unsigned int i = 0; i < size_1; ++i) {
1252
1253 seq_1 = (unsigned char**)patch_1->get(FEM, age, i)->getTrait(_SHLinkedTraitIndex)->get_sequence();
1254
1255 for (unsigned int j = 0; j < size_2; ++j) {
1256
1257 seq_2 = (unsigned char**)patch_2->get(FEM, age, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
1258
1259 for (unsigned int l = 0; l < num_loc; ++l) {
1260 D += (seq_1[0][l] != seq_2[0][l]);
1261 D += (seq_1[1][l] != seq_2[1][l]);
1262 D += (seq_1[0][l]!=seq_2[1][l]);
1263 D += (seq_1[1][l]!=seq_2[0][l]);
1264 }
1265 N += 4;
1266 }
1267
1268 for (unsigned int j = 0; j < patch_2->size(MAL, age); ++j) {
1269
1270 seq_2 = (unsigned char**)patch_2->get(MAL, age, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
1271
1272 for (unsigned int l = 0; l < num_loc; ++l) {
1273 D += (seq_1[0][l]!=seq_2[0][l]);
1274 D += (seq_1[1][l]!=seq_2[1][l]);
1275 D += (seq_1[0][l]!=seq_2[1][l]);
1276 D += (seq_1[1][l]!=seq_2[0][l]);
1277 }
1278 N += 4;
1279
1280 }
1281 }
1282
1283 //same for the males:
1284 size_1 = patch_1->size(MAL, age);
1285 size_2 = patch_2->size(MAL, age);
1286
1287// N += 4 * size_1 * size_2;
1288
1289// if( N == 0) return 0; //needless to go further
1290
1291 for (unsigned int i = 0; i < size_1; ++i) {
1292
1293 seq_1 = (unsigned char**)patch_1->get(MAL, age, i)->getTrait(_SHLinkedTraitIndex)->get_sequence();
1294
1295 for (unsigned int j = 0; j < size_2; ++j) {
1296
1297 seq_2 = (unsigned char**)patch_2->get(MAL, age, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
1298
1299 for (unsigned int l = 0; l < num_loc; ++l) {
1300 D += seq_1[0][l]!=seq_2[0][l];
1301 D += seq_1[1][l]!=seq_2[1][l];
1302 D += seq_1[0][l]!=seq_2[1][l];
1303 D += seq_1[1][l]!=seq_2[0][l];
1304 }
1305 N += 4;
1306 }
1307
1308 for (unsigned int j = 0; j < patch_2->size(FEM, age); ++j) {
1309
1310 seq_2 = (unsigned char**)patch_2->get(FEM, age, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
1311
1312 for (unsigned int l = 0; l < num_loc; ++l) {
1313 D += seq_1[0][l]!=seq_2[0][l];
1314 D += seq_1[1][l]!=seq_2[1][l];
1315 D += seq_1[0][l]!=seq_2[1][l];
1316 D += seq_1[1][l]!=seq_2[0][l];
1317 }
1318 N += 4;
1319
1320 }
1321 }
1322
1323 if( N == 0) return 0;
1324
1325 return D/N;
1326}
TTrait * getTrait(IDX T)
Trait accessor.
Definition: individual.h:277
Patch * getPatchPtr(unsigned int patch)
A secure version of the getPatch() method.
Definition: metapop.h:260
Second class in the metapopulation design structure, between the Metapop and Individual classes.
Definition: metapop.h:430
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_locus_num()
Definition: ttneutralgenes.h:190
virtual void ** get_sequence() const =0
sequence accessor.
int _SHLinkedTraitIndex
Index of the trait in the Individual::Traits table.
Definition: stathandler.h:173
TProtoNeutralGenes * _SHLinkedTrait
Pointer to a TraitProtoype object.
Definition: stathandler.h:171
@ FEM
Definition: types.h:37
@ MAL
Definition: types.h:37

References StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, FEM, Patch::get(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Metapop::getPatchPtr(), Individual::getTrait(), MAL, and Patch::size().

Referenced by getDxy(), getDxyAdultPerPatch(), and getDxyOffspringPerPatch().

+ Here is the caller graph for this function:

◆ getFis()

double TTNeutralGenesSH::getFis ( )
inline
422{return _fis;}
double _fis
Definition: ttneutralgenes.h:295

References _fis.

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ getFisWC()

double TTNeutralGenesSH::getFisWC ( )
inline
438{return _fis_WC;}
double _fis_WC
Definition: ttneutralgenes.h:299

References _fis_WC.

Referenced by setFstatWCRecorders().

+ Here is the caller graph for this function:

◆ getFit()

double TTNeutralGenesSH::getFit ( )
inline
423{return _fit;}
double _fit
Definition: ttneutralgenes.h:295

References _fit.

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ getFitWC()

double TTNeutralGenesSH::getFitWC ( )
inline
439{return _fit_WC;}
double _fit_WC
Definition: ttneutralgenes.h:299

References _fit_WC.

Referenced by setFstatWCRecorders().

+ Here is the caller graph for this function:

◆ getFixLocGlobal()

double TTNeutralGenesSH::getFixLocGlobal ( )
inline
445{return _fix_loc_global;}
double _fix_loc_global
Definition: ttneutralgenes.h:295

References _fix_loc_global.

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ getFixLocLocal()

double TTNeutralGenesSH::getFixLocLocal ( )
inline
444{return _fix_loc_local;}
double _fix_loc_local
Definition: ttneutralgenes.h:295

References _fix_loc_local.

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ getFst()

double TTNeutralGenesSH::getFst ( )
inline
421{return _fst;}
double _fst
Definition: ttneutralgenes.h:295

References _fst.

Referenced by setFstat2Recorders(), and setFstatRecorders().

+ Here is the caller graph for this function:

◆ getFst_ij()

double TTNeutralGenesSH::getFst_ij ( unsigned int  i)
inline

Accessor to the Fst matrix as set by setFstMatrix().

400 {
401 unsigned int scale = (unsigned int)pow( 10.0, (int)log10((float)_fst_matrix->getNbCols()) + 1 );
402 return _fst_matrix->get(i/scale, i%scale);
403 }

References _fst_matrix, TMatrix::get(), and TMatrix::getNbCols().

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ getFstWC()

double TTNeutralGenesSH::getFstWC ( )
inline
437{return _fst_WC;}
double _fst_WC
Weir & Cockerham (1984) F-stat estimates.
Definition: ttneutralgenes.h:299

References _fst_WC.

Referenced by setFstatWCRecorders().

+ Here is the caller graph for this function:

◆ getGlobalAlleleFreq()

double TTNeutralGenesSH::getGlobalAlleleFreq ( unsigned int  loc,
unsigned int  all 
)
inline
351 {
352 return _globalAlleleFreq.get(loc, all);
353 }

References _globalAlleleFreq, and TMatrix::get().

Referenced by setFreqRecorders().

+ Here is the caller graph for this function:

◆ getGlobalFreqs()

TMatrix * TTNeutralGenesSH::getGlobalFreqs ( )
inline

Accessor to the table of allele frequencies in the whole population.

377{return &_globalAlleleFreq;}

References _globalAlleleFreq.

Referenced by TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ getHeteroTable()

DataTable< double > * TTNeutralGenesSH::getHeteroTable ( )
inline
374{return &_heteroTable;}
DataTable< double > _heteroTable
Definition: ttneutralgenes.h:281

References _heteroTable.

Referenced by TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ getHeterozygosity()

double TTNeutralGenesSH::getHeterozygosity ( unsigned int  loc)
inline
355 {
356 double het = 0;
357 for(unsigned int i = 0; i < _heteroTable.getNumGroups(); ++i )
358 het += _heteroTable.get(i, loc, 0);
359 return het/_heteroTable.getNumGroups(); //mean per patch heterozygosity
360 }
T get(unsigned int group, unsigned int Class, unsigned int elmnt)
Returns value stored of the element 'elmnt' of the class 'Class' in the group 'group'.
Definition: datatable.h:228
unsigned int getNumGroups()
Definition: datatable.h:258

References _heteroTable, DataTable< T >::get(), and DataTable< T >::getNumGroups().

Referenced by setFreqRecorders().

+ Here is the caller graph for this function:

◆ getHo()

double TTNeutralGenesSH::getHo ( )
inline
418{return _ho;}
double _ho
F-statistics.
Definition: ttneutralgenes.h:294

References _ho.

Referenced by setFstat2Recorders(), and setFstatRecorders().

+ Here is the caller graph for this function:

◆ getHs()

double TTNeutralGenesSH::getHs ( )
inline
419{return _hs;}
double _hs
Definition: ttneutralgenes.h:294

References _hs.

◆ getHsnei()

double TTNeutralGenesSH::getHsnei ( )
inline
416{return _hsnei;}
double _hsnei
Definition: ttneutralgenes.h:294

References _hsnei.

Referenced by setFstat2Recorders(), and setFstatRecorders().

+ Here is the caller graph for this function:

◆ getHt()

double TTNeutralGenesSH::getHt ( )
inline
420{return _ht;}
double _ht
Definition: ttneutralgenes.h:294

References _ht.

◆ getHtnei()

double TTNeutralGenesSH::getHtnei ( )
inline
417{return _htnei;}
double _htnei
Definition: ttneutralgenes.h:294

References _htnei.

Referenced by setFstat2Recorders(), and setFstatRecorders().

+ Here is the caller graph for this function:

◆ getMeanAlpha()

double TTNeutralGenesSH::getMeanAlpha ( )
inline
483{return _mean_alpha;}
double _mean_alpha
Definition: ttneutralgenes.h:286

References _mean_alpha.

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ getMeanNeiGeneticDistance()

double TTNeutralGenesSH::getMeanNeiGeneticDistance ( )
inline
506{return _meanD;}
double _meanD
Definition: ttneutralgenes.h:309

References _meanD.

Referenced by setNeiGeneticDistanceRecorders().

+ Here is the caller graph for this function:

◆ getMeanTheta()

double TTNeutralGenesSH::getMeanTheta ( )
inline
482{return _mean_theta;}
double _mean_theta
Definition: ttneutralgenes.h:286

References _mean_theta.

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ getNbAllGlobal()

double TTNeutralGenesSH::getNbAllGlobal ( )
inline
443{return _nb_all_global;}
double _nb_all_global
Definition: ttneutralgenes.h:294

References _nb_all_global.

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ getNbAllLocal()

double TTNeutralGenesSH::getNbAllLocal ( )
inline
442{return _nb_all_local;}
double _nb_all_local
Definition: ttneutralgenes.h:294

References _nb_all_local.

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ getNeiGeneticDistance()

double TTNeutralGenesSH::getNeiGeneticDistance ( unsigned int  i)
inline
502 {
503 unsigned int scale = (unsigned int)pow( 10.0, (int)log10((float)_D->getNbCols()) + 1 );
504 return _D->get(i/scale,i%scale);
505 }

References _D, TMatrix::get(), and TMatrix::getNbCols().

Referenced by setNeiGeneticDistanceRecorders().

+ Here is the caller graph for this function:

◆ getSibCoaMeans()

double TTNeutralGenesSH::getSibCoaMeans ( unsigned int  i)
inline
493{return _sib_coa[i];}
double _sib_coa[4]
Definition: ttneutralgenes.h:291

References _sib_coa.

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ getSibProportions()

double TTNeutralGenesSH::getSibProportions ( unsigned int  i)
inline
492{return _sib_prop[i];}
double _sib_prop[4]
Kinship classes proportions.
Definition: ttneutralgenes.h:290

References _sib_prop.

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ getTheta_FF()

double TTNeutralGenesSH::getTheta_FF ( )
inline

Gives the mean within females coancestry coefficient.

485{return Theta_FF;}
double Theta_FF
Definition: ttneutralgenes.h:285

References Theta_FF.

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ getTheta_FM()

double TTNeutralGenesSH::getTheta_FM ( )
inline

Gives the mean between males and females coancestry coefficient.

489{return Theta_FM;}
double Theta_FM
Definition: ttneutralgenes.h:285

References Theta_FM.

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ getTheta_MM()

double TTNeutralGenesSH::getTheta_MM ( )
inline

Gives the mean within males coancestry coefficient.

487{return Theta_MM;}
double Theta_MM
Definition: ttneutralgenes.h:285

References Theta_MM.

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ getWeightedFst()

double TTNeutralGenesSH::getWeightedFst ( )
inline

Returns the weighted Fst using Weir & Hill (2002) method.

This Fst is set by a previous call to setFstMatrix().

397{return _fst_WH;}
double _fst_WH
Weir & Hill (2002) F-stat estimates.
Definition: ttneutralgenes.h:297

References _fst_WH.

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ init()

void TTNeutralGenesSH::init ( )
virtual

Reimplemented from StatHandlerBase.

1734{
1736
1738}
virtual void init()=0
Inits state.
unsigned int get_allele_num()
Definition: ttneutralgenes.h:191
void allocateTables(unsigned int loci, unsigned int all)
Definition: stats_fstat.cc:43

References TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, allocateTables(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), and Handler::init().

◆ setAdltNeiGeneticDistance()

void TTNeutralGenesSH::setAdltNeiGeneticDistance ( )
inline
void setNeiGeneticDistance(age_t AGE)
Definition: stats_fstat.cc:1129

References ADULTS, and setNeiGeneticDistance().

Referenced by setNeiGeneticDistanceRecorders().

+ Here is the caller graph for this function:

◆ setAdultAlleleFreq()

void TTNeutralGenesSH::setAdultAlleleFreq ( )
inline
void setAlleleTables(age_t AGE)
Definition: stats_fstat.cc:74

References ADULTS, and setAlleleTables().

Referenced by setFreqRecorders().

+ Here is the caller graph for this function:

◆ setAdultHeterozygosity()

void TTNeutralGenesSH::setAdultHeterozygosity ( )
inline
void setHeterozygosity(age_t AGE)
Definition: stats_fstat.cc:225

References ADULTS, and setHeterozygosity().

Referenced by setFreqRecorders().

+ Here is the caller graph for this function:

◆ setAdults_Theta()

void TTNeutralGenesSH::setAdults_Theta ( )
204{
205 unsigned int Fsize,Msize,i,j,k, FFsize, MMsize, FMsize;
206 unsigned int patchNbr = this->_pop->getPatchNbr();
207 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
208 Patch *P1;
209 double mean = 0, grand_mean = 0;
210
211 Theta_FF = 0;
212 Theta_MM = 0;
213 Theta_FM = 0;
214
215 _mean_theta = 0;
216
217 for(i = 0; i < patchNbr; ++i) {
218
219 P1 = _pop->getPatch(i);
220 Fsize = P1->size(FEM, ADLTx);
221 Msize = P1->size(MAL, ADLTx);
222
223 FFsize = Fsize*(Fsize-1)/2;
224 MMsize = Msize*(Msize-1)/2;
225 FMsize = Fsize*Msize;
226
227 grand_mean = 0;
228
229 if(Fsize != 0 || Msize != 0) {
230
231 if(Fsize != 0) {
232 mean = 0;
233 for(j = 0; j < Fsize-1;++j)
234 for(k = j+1; k < Fsize;++k)
237 nb_locus);
238 Theta_FF += mean/FFsize;
239 grand_mean += mean;
240 }
241
242 if(Msize != 0) {
243 mean = 0;
244 for(j = 0; j < Msize-1;++j)
245 for(k = j+1; k < Msize;++k)
248 nb_locus);
249 Theta_MM += mean/MMsize;
250 grand_mean += mean;
251 }
252
253
254 if(Fsize != 0 && Msize != 0) {
255 mean = 0;
256 for(j = 0; j < Fsize;++j)
257 for(k = 0; k < Msize;++k)
260 nb_locus);
261 Theta_FM += mean/FMsize;
262 grand_mean += mean;
263 }
264 _mean_theta += grand_mean / (FFsize + MMsize + FMsize);
265 }
266 }
267
268 _mean_theta /= patchNbr;
269 Theta_FF /= patchNbr;
270 Theta_MM /= patchNbr;
271 Theta_FM /= patchNbr;
272}
Patch * getPatch(unsigned int i)
Patch accessor, return the ith+1 patch in the metapop.
Definition: metapop.h:257
double Coancestry(void **ind1, void **ind2, unsigned int nb_locus)
Gives the coancestry (probability of identity by state) of two gene sequences.
Definition: stats_coa.cc:43

References _mean_theta, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, ADLTx, Coancestry(), FEM, Patch::get(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), MAL, Patch::size(), Theta_FF, Theta_FM, and Theta_MM.

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setAdultsCoaBetween()

void TTNeutralGenesSH::setAdultsCoaBetween ( )
inline
469{setCoaMatrix(ADLTx, 2);}
void setCoaMatrix(age_idx age_pos, unsigned char dim)
Computes the within and between patches coancestry coefficients.
Definition: stats_coa.cc:59

References ADLTx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ setAdultsCoaMatrix()

void TTNeutralGenesSH::setAdultsCoaMatrix ( )
inline
465{setCoaMatrix(ADLTx, 3);}

References ADLTx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ setAdultsCoaWithin()

void TTNeutralGenesSH::setAdultsCoaWithin ( )
inline
467{setCoaMatrix(ADLTx, 1);}

References ADLTx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ setAdultsFstat()

void TTNeutralGenesSH::setAdultsFstat ( )
inline
void setFstat(age_t AGE)
Computes the F-statistics following Nei & Chesser (1983).
Definition: stats_fstat.cc:254

References ADULTS, and setFstat().

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ setAdultsFstat2()

void TTNeutralGenesSH::setAdultsFstat2 ( )
inline
void setFstat2(age_t AGE)
New version of Nei & Chesser.
Definition: stats_fstat.cc:482

References ADULTS, and setFstat2().

Referenced by setFstat2Recorders().

+ Here is the caller graph for this function:

◆ setAdultsFstatWeirCockerham()

void TTNeutralGenesSH::setAdultsFstatWeirCockerham ( )
inline
void setFstatWeirCockerham(age_t AGE)
Computes the Weir & Cockerham (1984) Fstat values (Theta, F, and f).
Definition: stats_fstat.cc:771

References ADULTS, and setFstatWeirCockerham().

Referenced by setFstatWCRecorders().

+ Here is the caller graph for this function:

◆ setAdultsFstBetween()

void TTNeutralGenesSH::setAdultsFstBetween ( )
inline
391{setFstMatrix(ADULTS, 2);}
void setFstMatrix(age_t AGE, unsigned char dim)
Computes the weighted within and between patch Fst's as well as the overall Fst (Theta).
Definition: stats_fstat.cc:665

References ADULTS, and setFstMatrix().

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ setAdultsFstMatrix()

void TTNeutralGenesSH::setAdultsFstMatrix ( )
inline
389{setFstMatrix(ADULTS, 3);}

References ADULTS, and setFstMatrix().

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ setAdultsFstWithin()

void TTNeutralGenesSH::setAdultsFstWithin ( )
inline
390{setFstMatrix(ADULTS, 1);}

References ADULTS, and setFstMatrix().

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ setAlleleTables()

void TTNeutralGenesSH::setAlleleTables ( age_t  AGE)
75{
76
77 if(_table_set_age == AGE
80 return;
81
82 unsigned char** genes;
83 unsigned int nb_locus = _SHLinkedTrait->get_locus_num(), nb_allele = _SHLinkedTrait->get_allele_num();
84 unsigned int patch_size, patchNbr = _pop->getPatchNbr();
85 age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
86 Patch* crnt_patch;
87
88 //check if the population size has changed:
89 if(_alleleCountTable.getNumGroups() != patchNbr)
91
93
94 //counting the copies of each allele present in each patch:
95 for(unsigned int k = 0; k < patchNbr; ++k) {
96
97 crnt_patch = _pop->getPatch(k);
98
99 for(unsigned int i = 0, size = crnt_patch->size(FEM, age_pos); i < size; ++i) {
100
101 genes = (unsigned char**)crnt_patch->get(FEM, age_pos, i)->getTrait(_SHLinkedTraitIndex)->get_sequence();
102
103 for (unsigned int j = 0; j < nb_locus; ++j)
104 _alleleCountTable.increment(k, j, genes[0][j]);
105
106 for (unsigned int j = 0; j < nb_locus; ++j)
107 _alleleCountTable.increment(k, j, genes[1][j]);
108
109 }
110
111 for(unsigned int i = 0, size = crnt_patch->size(MAL, age_pos); i < size; ++i) {
112
113 genes = (unsigned char**)crnt_patch->get(MAL, age_pos, i)->getTrait(_SHLinkedTraitIndex)->get_sequence();
114
115 for (unsigned int j = 0; j < nb_locus; ++j)
116 _alleleCountTable.increment(k, j, genes[0][j]);
117
118 for (unsigned int j = 0; j < nb_locus; ++j)
119 _alleleCountTable.increment(k, j, genes[1][j]);
120
121 }
122 }
123
124 //allelic frequencies:
125 for (unsigned int i = 0; i < patchNbr; ++i) {
126
127 patch_size = _pop->getPatch(i)->size(age_pos) * 2;
128
129 if (patch_size) {
130
131 for (unsigned int l = 0; l < nb_locus; ++l)
132 for (unsigned int u = 0; u < nb_allele; ++u)
133 _alleleFreqTable.set(i, l, u, (double)_alleleCountTable.get(i, l, u) / (double)patch_size);
134
135 } else {
136 for (unsigned int l = 0; l < nb_locus; ++l)
137 for (unsigned int u = 0; u < nb_allele; ++u)
138 _alleleFreqTable.set(i, l, u, 0 );
139 }
140
141 }
142
143 //population-wide allelic frequencies:
144 unsigned int tot_size = _pop->size(AGE) * 2;
145
147
148 for(unsigned int i = 0; i < patchNbr; i++)
149 for (unsigned int l = 0; l < nb_locus; ++l)
150 for (unsigned int u = 0; u < nb_allele; ++u)
152
153 for (unsigned int l = 0; l < nb_locus; ++l)
154 for (unsigned int u = 0; u < nb_allele; ++u)
155 _globalAlleleFreq.divide(l, u, tot_size);
156
157 _table_set_age = AGE;
158
160
162}
void set(unsigned int group, unsigned int Class, unsigned int elmnt, T val)
Sets the element 'elmnt' of the class 'Class' in the group 'group' to the value 'val'.
Definition: datatable.h:232
void init(T val)
Sets all elements of the table to value 'val'.
Definition: datatable.h:256
void increment(unsigned int group, unsigned int Class, unsigned int elmnt)
Increments 'elmnt' of the class 'Class' in the group 'group' by one.
Definition: datatable.h:236
unsigned int getCurrentReplicate()
Definition: metapop.h:293
unsigned int getCurrentGeneration()
Definition: metapop.h:294
void assign(double val)
Assigns a value to all element of the matrix.
Definition: tmatrix.h:110
void divide(unsigned int i, unsigned int j, double value)
Divide an element of the matrix by a value.
Definition: tmatrix.h:243
void plus(unsigned int i, unsigned int j, double value)
Adds a value to an element of the matrix.
Definition: tmatrix.h:210

References _alleleCountTable, _alleleFreqTable, _globalAlleleFreq, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, _table_set_age, _table_set_gen, _table_set_repl, ADLTx, ADULTS, allocateTables(), TMatrix::assign(), TMatrix::divide(), FEM, Patch::get(), DataTable< T >::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Metapop::getCurrentGeneration(), Metapop::getCurrentReplicate(), DataTable< T >::getNumGroups(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), DataTable< T >::increment(), DataTable< T >::init(), MAL, OFFSx, TMatrix::plus(), DataTable< T >::set(), Metapop::size(), and Patch::size().

Referenced by setAdultAlleleFreq(), setFst_li(), setFstat(), setFstat2(), setFstatWeirCockerham(), setFstatWeirCockerham_MS(), setFstMatrix(), setNeiGeneticDistance(), setOffspringAlleleFreq(), and TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ setCoaMatrix()

void TTNeutralGenesSH::setCoaMatrix ( age_idx  age_pos,
unsigned char  dim 
)

Computes the within and between patches coancestry coefficients.

Parameters
age_posthe age class index
dimthe dimension of the matrix to fill:
  • 1 = the diagonal (i.e. the wihtin patch coancestries or theta's)
  • 2 = the upper half (i.e. the between patch coancestries or alpha's)
  • 3 = both
60{
61 unsigned int Fsize,Msize,tot_size,size_i,size_l, wt;
62 unsigned int nb_coeff = 0;
63 unsigned int i, j, k, l;
64 unsigned int patchNbr = this->_pop->getPatchNbr();
65 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
66 Patch *P1, *P2;
67
68 if(_coa_matrix == NULL)
69
70 _coa_matrix = new TMatrix(patchNbr, patchNbr);
71
72 else if( _coa_matrix->length() != patchNbr * patchNbr)
73
74 _coa_matrix->reset(patchNbr, patchNbr);
75
77
78 if(dim & 1) {
79
80 _mean_theta = 0;
81
82 wt = 0;
83
84 //first fill the diagonale: within deme coancestry (theta)
85 for(i = 0; i < patchNbr; ++i) {
86
87 P1 = _pop->getPatch(i);
88 Fsize = P1->size(FEM, age_pos);
89 Msize = P1->size(MAL, age_pos);
90
91 tot_size = Fsize + Msize;
92
93 if(tot_size != 0) {
94
95 //fem-fem coa
96 for(j = 0; j < Fsize; ++j)
97 for(k = j+1; k < Fsize; ++k)
99 P1->get(FEM, age_pos, k)->getTrait(_SHLinkedTraitIndex)->get_sequence(),
100 nb_locus));
101 //mal-mal coa
102 for(j = 0; j < Msize; ++j)
103 for(k = j+1; k < Msize; ++k)
105 P1->get(MAL, age_pos, k)->getTrait(_SHLinkedTraitIndex)->get_sequence(),
106 nb_locus));
107 //fem-mal coa
108 for(j = 0; j < Fsize; ++j)
109 for(k = 0; k < Msize; ++k)
111 P1->get(MAL, age_pos, k)->getTrait(_SHLinkedTraitIndex)->get_sequence(),
112 nb_locus));
113
114 _coa_matrix->divide(i, i, tot_size*(tot_size -1)/2.0);
115 }//end if
116
117 _mean_theta += tot_size * _coa_matrix->get(i, i);
118
119 wt += tot_size;
120
121 }//end for patchNbr
122
123 //weighted average:
124 _mean_theta /= wt;
125
126 } //end if diag
127
128 if(dim & 2) {
129 //fill the first upper half of the matrix: between deme coancestry (alpha)
130
131 _mean_alpha = 0;
132
133 wt = 0;
134
135 for(i = 0; i < patchNbr-1; ++i) {
136
137 P1 = _pop->getPatch(i);
138
139 for(l = i+1; l < patchNbr; ++l) {
140
141 P2 = _pop->getPatch(l);
142
143 tot_size = P1->size(age_pos) * P2->size(age_pos);
144
145 if(tot_size != 0) {
146 //females i vs. females l
147 size_i = P1->size(FEM, age_pos);
148 size_l = P2->size(FEM, age_pos);
149 nb_coeff = size_i * size_l;
150
151 for(j = 0; j < size_i; ++j)
152 for(k = 0; k < size_l; ++k)
154 P2->get(FEM, age_pos, k)->getTrait(_SHLinkedTraitIndex)->get_sequence(),
155 nb_locus));
156 //females i vs. males l
157 size_l = P2->size(MAL, age_pos);
158 nb_coeff += size_i * size_l;
159
160 for(j = 0; j < size_i; ++j)
161 for(k = 0; k < size_l; ++k)
163 P2->get(MAL, age_pos, k)->getTrait(_SHLinkedTraitIndex)->get_sequence(),
164 nb_locus));
165 //males i vs. males l
166 size_i = P1->size(MAL, age_pos);
167 size_l = P2->size(MAL, age_pos);
168 nb_coeff += size_i * size_l;
169
170 for(j = 0; j < size_i; ++j)
171 for(k = 0; k < size_l; ++k)
173 P2->get(MAL, age_pos, k)->getTrait(_SHLinkedTraitIndex)->get_sequence(),
174 nb_locus));
175 //males i vs. females l
176 size_l = P2->size(FEM, age_pos);
177 nb_coeff += size_i * size_l;
178
179 for(j = 0; j < size_i; ++j)
180 for(k = 0; k < size_l; ++k)
182 P2->get(FEM, age_pos, k)->getTrait(_SHLinkedTraitIndex)->get_sequence(),
183 nb_locus));
184 _coa_matrix->divide(i, l, nb_coeff);
185 }//endif
186
187 _mean_alpha += tot_size * _coa_matrix->get(i, l);
188
189 wt += tot_size;
190
191 }//end for P2
192 }//end for P1
193
194 //weighted average:
195 _mean_alpha /= wt;
196
197 }//end if upper half
198
199}
A class to handle matrix in params, coerces matrix into a vector of same total size.
Definition: tmatrix.h:49
unsigned int length()
Returns the number of elements in the matrix.
Definition: tmatrix.h:172

References _coa_matrix, _mean_alpha, _mean_theta, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, TMatrix::assign(), Coancestry(), TMatrix::divide(), FEM, Patch::get(), TMatrix::get(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), TMatrix::length(), MAL, TMatrix::plus(), TMatrix::reset(), and Patch::size().

Referenced by setAdultsCoaBetween(), setAdultsCoaMatrix(), setAdultsCoaWithin(), setOffsprgCoaBetween(), setOffsprgCoaMatrix(), and setOffsprgCoaWithin().

+ Here is the caller graph for this function:

◆ setCoaMatrixRecorders()

void TTNeutralGenesSH::setCoaMatrixRecorders ( age_t  AGE,
unsigned char  dim 
)
1996{
1997 std::ostringstream name, sub_name;
1998
1999 void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ?
2006
2007 const char *prefix = (AGE == ADULTS ? "adlt." : "off.");
2008
2009 unsigned int nbpatch = _pop->getPatchNbr();
2010 unsigned int scale = (unsigned int)pow(10.0, (int)log10((float)nbpatch) + 1);
2011
2012 if(dim & 1) {
2013 name<<"Wtn Patch Coancestry ("<<prefix<<")";
2014 sub_name<< prefix << "theta";
2015 add(name.str(), sub_name.str(), AGE,0, 0, &TTNeutralGenesSH::getMeanTheta, 0, 0, setter);
2016 name.str("");
2017 sub_name.str("");
2018 }
2019 if(dim & 2){
2020 name<<"Btn Patch Coancestry ("<<prefix<<")";
2021 sub_name<< prefix << "alpha";
2022 add(name.str(), sub_name.str(), AGE, 0, 0, &TTNeutralGenesSH::getMeanAlpha, 0, 0, (!(dim&1) ? setter : 0) );
2023 name.str("");
2024 sub_name.str("");
2025 }
2026 if(dim & 1) {
2027 for(unsigned int i = 0; i < nbpatch; ++i) {
2028 name<<"coancestry "<<i+1<<"."<<i+1;
2029 sub_name<< prefix << "coa" << i+1 << "." << i+1;
2030 add(name.str(), sub_name.str(), AGE, i*scale + i, 0, 0, &TTNeutralGenesSH::getCoa, 0, 0);
2031 name.str("");
2032 sub_name.str("");
2033 }
2034 }
2035 if(dim & 2){
2036 for(unsigned int i = 0; i < nbpatch; ++i) {
2037 for(unsigned int j = i+1; j < nbpatch; ++j) {
2038 name<<"coancestry "<<i+1<< "." <<j+1;
2039 sub_name<< prefix << "coa" << i+1 << "." << j+1;
2040 add(name.str(), sub_name.str(), AGE, i*scale + j, 0, 0, &TTNeutralGenesSH::getCoa, 0, 0);
2041 name.str("");
2042 sub_name.str("");
2043 }
2044 }
2045 }
2046
2047}
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
The stat handler for neutral markers.
Definition: ttneutralgenes.h:277
void setAdultsCoaWithin()
Definition: ttneutralgenes.h:467
void setOffsprgCoaWithin()
Definition: ttneutralgenes.h:468
double getCoa(unsigned int i)
Gets the given coancestry coefficient from the coancestry matrix.
Definition: ttneutralgenes.h:477
void setOffsprgCoaBetween()
Definition: ttneutralgenes.h:470
double getMeanTheta()
Definition: ttneutralgenes.h:482
double getMeanAlpha()
Definition: ttneutralgenes.h:483
void setOffsprgCoaMatrix()
Definition: ttneutralgenes.h:466
void setAdultsCoaBetween()
Definition: ttneutralgenes.h:469
void setAdultsCoaMatrix()
Definition: ttneutralgenes.h:465

References StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, getCoa(), getMeanAlpha(), getMeanTheta(), Metapop::getPatchNbr(), setAdultsCoaBetween(), setAdultsCoaMatrix(), setAdultsCoaWithin(), setOffsprgCoaBetween(), setOffsprgCoaMatrix(), and setOffsprgCoaWithin().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setDxyRecorders()

void TTNeutralGenesSH::setDxyRecorders ( age_t  AGE,
bool  patchwise 
)
2238{
2239 string prefix = (AGE == OFFSPRG ? "off." : "adlt.");
2240 age_idx age = (AGE == OFFSPRG ? OFFSx : ADLTx );
2241 string name = "Sequence divergence - Dxy";
2242 string sub_name = "Dxy";
2243
2244
2245 if (!patchwise) {
2246
2247 add(name, prefix + sub_name, AGE, (unsigned int)age, 0, 0, &TTNeutralGenesSH::getDxy, 0, 0);
2248
2249 } else {
2250
2251 for (unsigned int p1 = 0; p1 < _pop->getPatchNbr(); ++p1) {
2252 for (unsigned int p2 = p1 + 1; p2 < _pop->getPatchNbr(); ++p2) {
2253 add(name, prefix + sub_name + ".p" + tstring::int2str(p1+1) + "p" + tstring::int2str(p2+1),
2254 AGE, p1, p2,
2257 0);
2258 }
2259 }
2260
2261 }
2262}
double getDxyAdultPerPatch(unsigned int patch1, unsigned patch2)
Definition: ttneutralgenes.h:512
double getDxy(unsigned int age_class)
Definition: stats_fstat.cc:1210
double getDxyOffspringPerPatch(unsigned int patch1, unsigned patch2)
Definition: ttneutralgenes.h:511
static string int2str(const int i)
Writes an integer value into a string.
Definition: tstring.h:95

References StatHandlerBase::_pop, StatHandler< SH >::add(), ADLTx, getDxy(), getDxyAdultPerPatch(), getDxyOffspringPerPatch(), Metapop::getPatchNbr(), tstring::int2str(), OFFSPRG, and OFFSx.

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setFreqRecorders()

void TTNeutralGenesSH::setFreqRecorders ( age_t  AGE)
2052{
2053 string prefix = (AGE == OFFSPRG ? "off." : "adlt.");
2054
2055 void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ? &TTNeutralGenesSH::setAdultAlleleFreq :
2057
2058 unsigned int nb_allele = _SHLinkedTrait->get_allele_num();
2059 unsigned int nb_locus = _SHLinkedTrait->get_locus_num();
2060
2061 for (unsigned int l = 0; l < nb_locus; ++l) {
2062 for (unsigned int u = 0; u < nb_allele-1; ++u) {
2063 add("", prefix + "ntrl.l" + tstring::int2str(l+1) + ".a" + tstring::int2str(u+1), AGE, l, u,
2065 }
2066 }
2067
2070
2071 for (unsigned int l = 0; l < nb_locus; ++l) {
2072 add("", prefix + "ntrl.l" + tstring::int2str(l+1) + ".Het", AGE, l, 0, 0,
2074 // add("", prefix + "ntrl.l" + tstring::int2str(l) + ".Hom", AGE, l, 0, 0,
2075 // &TTNeutralGenesSH::getHomozygosity, 0, 0);
2076 }
2077
2078}
void setAdultHeterozygosity()
Definition: ttneutralgenes.h:348
double getGlobalAlleleFreq(unsigned int loc, unsigned int all)
Definition: ttneutralgenes.h:351
double getHeterozygosity(unsigned int loc)
Definition: ttneutralgenes.h:355
void setAdultAlleleFreq()
Definition: ttneutralgenes.h:345
void setOffspringHeterozygosity()
Definition: ttneutralgenes.h:349
void setOffspringAlleleFreq()
Definition: ttneutralgenes.h:346

References TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, StatHandler< SH >::add(), ADULTS, TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), getGlobalAlleleFreq(), getHeterozygosity(), tstring::int2str(), OFFSPRG, setAdultAlleleFreq(), setAdultHeterozygosity(), setOffspringAlleleFreq(), and setOffspringHeterozygosity().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setFreqRecordersPerPatch()

void TTNeutralGenesSH::setFreqRecordersPerPatch ( age_t  AGE)

◆ setFst_li()

void TTNeutralGenesSH::setFst_li ( unsigned int  N,
unsigned int  L,
double **  array 
)

Computes the per-locus per-patch Fst values using Weir&Hill 2002 approach.

1058{
1059 unsigned int patchNbr = this->_pop->getPatchNbr();
1060 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
1061 unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
1062
1063 assert(N == patchNbr);
1064 assert(L == nb_locus);
1065
1067
1068 double *pop_weights = new double[patchNbr];
1069 double *pop_sizes = new double[patchNbr];
1070 double *denominator = new double[nb_locus];
1071 double sum_weights = 0;
1072
1073 double tot_size = _pop->size(ADULTS) * 2; //diploids!
1074
1075 for(unsigned int i = 0; i < patchNbr; ++i) {
1076 pop_sizes[i] = _pop->size(ADULTS, i) * 2;
1077 pop_weights[i] = pop_sizes[i] - (pop_sizes[i] * pop_sizes[i] / tot_size); //n_ic in Weir & Hill 2002
1078 sum_weights += pop_weights[i];
1079 for(unsigned int j = 0; j < nb_locus; ++j)
1080 array[i][j] = nanf("NULL");
1081 }
1082
1083 for(unsigned int l = 0; l < nb_locus; ++l)
1084 denominator[l] = 0;
1085
1086 double p, pq, var;
1087
1088 for (unsigned int i = 0; i < patchNbr; ++i) {
1089
1090 if( !pop_sizes[i] ) continue;
1091
1092 for (unsigned int l = 0; l < nb_locus; ++l) {
1093
1094 array[i][l] = 0;
1095
1096 for (unsigned int u = 0; u < nb_allele; ++u) {
1097
1098 p = _alleleFreqTable.get(i, l, u); //p_liu
1099
1100 pq = p * (1 - p);
1101
1102 var = p - _globalAlleleFreq.get(l, u); //(p_liu - pbar_u)
1103
1104 var *= var;
1105
1106 array[i][l] += pq * pop_sizes[i] / (pop_sizes[i] -1);
1107
1108 denominator[l] += pop_sizes[i] * var + pop_weights[i] * pq;
1109
1110 } // end for allele
1111 }// end for locus
1112 }//end for pop
1113
1114
1115 for (unsigned int i = 0; i < patchNbr; ++i) {
1116 if( !pop_sizes[i] ) continue;
1117 for (unsigned int l = 0; l < nb_locus; ++l)
1118 array[i][l] = 1 - ( array[i][l] * sum_weights / denominator[l]);
1119 }
1120
1121 delete [] pop_weights;
1122 delete [] pop_sizes;
1123 delete [] denominator;
1124}

References _alleleFreqTable, _globalAlleleFreq, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, ADULTS, DataTable< T >::get(), TMatrix::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatchNbr(), setAlleleTables(), and Metapop::size().

Referenced by TTNeutralGenesFH::write_Fst_i().

+ Here is the caller graph for this function:

◆ setFstat()

void TTNeutralGenesSH::setFstat ( age_t  AGE)

Computes the F-statistics following Nei & Chesser (1983).

255{
256 double harmonic = 0, nbpatch = 0, nbind;
257 unsigned int patchNbr = _pop->getPatchNbr();
258 age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
259
260 //check if the population size has changed:
261 if(_alleleCountTable.getNumGroups() != patchNbr)
263
264 setAlleleTables(AGE);
265
267
268 //harmonic mean of patch sizes:
269 for (unsigned int i = 0; i < patchNbr; ++i){
270
271 nbind = _pop->size(AGE, i);
272 if(nbind != 0){
273 nbpatch++;
274 harmonic += 1.0/nbind;
275 }
276 }
277
278 harmonic = nbpatch / harmonic;
279
280 _ho = setHo(age_pos);
281 _hs = setHs(age_pos);
282 _ht = setHt(age_pos);
283 //Nei's corrections:
284 _hsnei = (nbpatch != 0 ? harmonic/(harmonic-1.0)*(_hs-(_ho/(2.0*harmonic))) : nanf("NULL") );
285 _htnei = (nbpatch != 0 ? _ht + (_hsnei/(harmonic*nbpatch))
286 -(_ho/(2.0*harmonic*nbpatch)) : nanf("NULL") );
287
288 _fis = ( _hsnei ? 1.0-(_ho/_hsnei) : nanf("NULL") );
289 _fit = ( _htnei ? 1.0-(_ho/_htnei) : nanf("NULL") );
290 _fst = ( _htnei ? 1.0-(_hsnei/_htnei) : nanf("NULL") );
291}
double setHo(age_idx age_pos)
Definition: stats_fstat.cc:363
double setHt(age_idx age_pos)
Definition: stats_fstat.cc:455
void setLociDivCounter(age_t AGE)
Sets the allelic diversity counters.
Definition: stats_fstat.cc:295
double setHs(age_idx age_pos)
Definition: stats_fstat.cc:410

References _alleleCountTable, _fis, _fit, _fst, _ho, _hs, _hsnei, _ht, _htnei, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, ADLTx, ADULTS, allocateTables(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), DataTable< T >::getNumGroups(), Metapop::getPatchNbr(), OFFSx, setAlleleTables(), setHo(), setHs(), setHt(), setLociDivCounter(), and Metapop::size().

Referenced by setAdultsFstat(), and setOffsprgFstat().

+ Here is the caller graph for this function:

◆ setFstat2()

void TTNeutralGenesSH::setFstat2 ( age_t  AGE)

New version of Nei & Chesser.

483{
484 double harmonic = 0, nbpatch = 0, nbind;
485 unsigned int patchNbr = _pop->getPatchNbr();
486 age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
487 unsigned int nloc = _SHLinkedTrait->get_locus_num();
488
489 //check if the population size has changed:
490 if(_alleleCountTable.getNumGroups() != patchNbr)
492
493 setAlleleTables(AGE);
494
496
497 //harmonic mean of patch sizes:
498 for (unsigned int i = 0; i < patchNbr; ++i){
499
500 nbind = _pop->size(AGE, i);
501 if(nbind != 0){
502 nbpatch++;
503 harmonic += 1.0/nbind;
504 }
505 }
506
507 harmonic = nbpatch / harmonic;
508
509 deque<double> ho = setHo2(age_pos);
510 deque<double> hs = setHs2(age_pos);
511 deque<double> ht = setHt2(age_pos);
512
513 _fis = _fit = _fst = 0;
514
515 double hsnei, htnei;
516
517 for (unsigned int l = 0; l < nloc; ++l) {
518
519 //Nei's corrections:
520 hsnei = (nbpatch != 0 ? harmonic/(harmonic-1.0)*(hs[l]-(ho[l]/(2.0*harmonic))) : nanf("NULL") );
521 htnei = (nbpatch != 0 ? ht[l] + (hsnei/(harmonic*nbpatch))
522 -(ho[l]/(2.0*harmonic*nbpatch)) : 0 );
523 _hsnei += hsnei;
524 _htnei += htnei;
525// _hsnei += hs[l];
526// _htnei += ht[l];
527 _ho += ho[l];
528// _fis = ( _hsnei ? 1.0-(_ho/_hsnei) : nanf("NULL") );
529// _fit = ( _htnei ? 1.0-(_ho/_htnei) : nanf("NULL") );
530 if(ht[l] == 0) nloc--; //would mean locus is fixed
531 _fst += ( ht[l] != 0 ? 1.0-(hsnei/htnei) : 0);
532
533// if(ht[l] == 0) nloc--;
534//
535// _fst += (ht[l] != 0 ? 1 - (hs[l]/ht[l]) : 0);
536 }
537 _hsnei /= nloc;
538 _htnei /= nloc;
539 _ho /= nloc;
540 _fst /= nloc;
541}
deque< double > setHo2(age_idx age_pos)
Definition: stats_fstat.cc:545
deque< double > setHt2(age_idx age_pos)
Definition: stats_fstat.cc:639
deque< double > setHs2(age_idx age_pos)
Definition: stats_fstat.cc:589

References _alleleCountTable, _fis, _fit, _fst, _ho, _hsnei, _htnei, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, ADLTx, ADULTS, allocateTables(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), DataTable< T >::getNumGroups(), Metapop::getPatchNbr(), OFFSx, setAlleleTables(), setHo2(), setHs2(), setHt2(), setLociDivCounter(), and Metapop::size().

Referenced by setAdultsFstat2(), and setOffsprgFstat2().

+ Here is the caller graph for this function:

◆ setFstat2Recorders()

void TTNeutralGenesSH::setFstat2Recorders ( age_t  AGE)
2118{
2119 if(AGE & OFFSPRG) {
2120 add("Ho (offsprg)","off.ho2",OFFSPRG,0,0,&TTNeutralGenesSH::getHo,0,0,
2122 add("Hs-Nei (offsprg)","off.hsnei2",OFFSPRG,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2123 add("Ht-Nei (offsprg)","off.htnei2",OFFSPRG,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2124 add("Fst (offsprg)","off.fst2",OFFSPRG,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2125 }
2126
2127 if(AGE & ADULTS) {
2128 add("Ho (adult)","adlt.ho2",ADULTS,0,0,&TTNeutralGenesSH::getHo,0,0,
2130 add("Hs-Nei (adult)","adlt.hsnei2",ADULTS,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2131 add("Ht-Nei (adult)","adlt.htnei2",ADULTS,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2132 add("Fst (adult)","adlt.fst2",ADULTS,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2133 }
2134}
void setAdultsFstat2()
Definition: ttneutralgenes.h:427
double getHsnei()
Definition: ttneutralgenes.h:416
double getHo()
Definition: ttneutralgenes.h:418
double getHtnei()
Definition: ttneutralgenes.h:417
void setOffsprgFstat2()
Definition: ttneutralgenes.h:426
double getFst()
Definition: ttneutralgenes.h:421

References StatHandler< SH >::add(), ADULTS, getFst(), getHo(), getHsnei(), getHtnei(), OFFSPRG, setAdultsFstat2(), and setOffsprgFstat2().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setFstatRecorders()

void TTNeutralGenesSH::setFstatRecorders ( age_t  AGE)
2084{
2085 if(AGE & OFFSPRG) {
2086 add("Nbr of Alleles per Locus - local (offsprg)","off.allnbp",OFFSPRG,0,0,
2088 add("Nbr of Alleles per Locus - global (offsprg)","off.allnb",OFFSPRG,0,0,&TTNeutralGenesSH::getNbAllGlobal,0,0,0);
2089 add("Nbr of Fixed Loci per Patch (offsprg)","off.fixlocp",OFFSPRG,0,0,&TTNeutralGenesSH::getFixLocLocal,0,0,0);
2090 add("Nbr of Fixed Loci in the Pop (offsprg)","off.fixloc",OFFSPRG,0,0,&TTNeutralGenesSH::getFixLocGlobal,0,0,0);
2091 add("Ho (offsprg)","off.ho",OFFSPRG,0,0,&TTNeutralGenesSH::getHo,0,0,0);
2092 add("Hs-Nei (offsprg)","off.hsnei",OFFSPRG,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2093 add("Ht-Nei (offsprg)","off.htnei",OFFSPRG,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2094 add("Fis (offsprg)","off.fis",OFFSPRG,0,0,&TTNeutralGenesSH::getFis,0,0,0);
2095 add("Fst (offsprg)","off.fst",OFFSPRG,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2096 add("Fit (offsprg)","off.fit",OFFSPRG,0,0,&TTNeutralGenesSH::getFit,0,0,0);
2097 }
2098
2099 if(AGE & ADULTS) {
2100 add("Nbr of Alleles per Locus - local (adult)","adlt.allnbp",ADULTS,0,0,
2102 add("Nbr of Alleles per Locus - global (adult)","adlt.allnb",ADULTS,0,0,&TTNeutralGenesSH::getNbAllGlobal,0,0,0);
2103 add("Nbr of Fixed Loci per Patch (adult)","adlt.fixlocp",ADULTS,0,0,&TTNeutralGenesSH::getFixLocLocal,0,0,0);
2104 add("Nbr of Fixed Loci in the Pop (adult)","adlt.fixloc",ADULTS,0,0,&TTNeutralGenesSH::getFixLocGlobal,0,0,0);
2105 add("Ho (adult)","adlt.ho",ADULTS,0,0,&TTNeutralGenesSH::getHo,0,0,0);
2106 add("Hs-Nei (adult)","adlt.hsnei",ADULTS,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2107 add("Ht-Nei (adult)","adlt.htnei",ADULTS,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2108 add("Fis (adult)","adlt.fis",ADULTS,0,0,&TTNeutralGenesSH::getFis,0,0,0);
2109 add("Fst (adult)","adlt.fst",ADULTS,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2110 add("Fit (adult)","adlt.fit",ADULTS,0,0,&TTNeutralGenesSH::getFit,0,0,0);
2111 }
2112}
double getFixLocGlobal()
Definition: ttneutralgenes.h:445
double getFit()
Definition: ttneutralgenes.h:423
void setAdultsFstat()
Definition: ttneutralgenes.h:412
double getNbAllGlobal()
Definition: ttneutralgenes.h:443
double getFis()
Definition: ttneutralgenes.h:422
double getNbAllLocal()
Definition: ttneutralgenes.h:442
double getFixLocLocal()
Definition: ttneutralgenes.h:444
void setOffsprgFstat()
Definition: ttneutralgenes.h:411

References StatHandler< SH >::add(), ADULTS, getFis(), getFit(), getFixLocGlobal(), getFixLocLocal(), getFst(), getHo(), getHsnei(), getHtnei(), getNbAllGlobal(), getNbAllLocal(), OFFSPRG, setAdultsFstat(), and setOffsprgFstat().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setFstatWCRecorders()

void TTNeutralGenesSH::setFstatWCRecorders ( age_t  AGE)
2140{
2141 if(AGE & OFFSPRG) {
2142
2143 add("Fis Weir & Cockerham","off.fis.WC",OFFSPRG,0,0,&TTNeutralGenesSH::getFisWC,0,0,&TTNeutralGenesSH::setOffspringFstatWeirCockerham);
2144 add("Fst Weir & Cockerham","off.fst.WC",OFFSPRG,0,0,&TTNeutralGenesSH::getFstWC,0,0,0);
2145 add("Fit Weir & Cockerham","off.fit.WC",OFFSPRG,0,0,&TTNeutralGenesSH::getFitWC,0,0,0);
2146 }
2147
2148 if(AGE & ADULTS) {
2149
2150 add("Fis Weir & Cockerham","adlt.fis.WC",ADULTS,0,0,&TTNeutralGenesSH::getFisWC,0,0,&TTNeutralGenesSH::setAdultsFstatWeirCockerham);
2151 add("Fst Weir & Cockerham","adlt.fst.WC",ADULTS,0,0,&TTNeutralGenesSH::getFstWC,0,0,0);
2152 add("Fit Weir & Cockerham","adlt.fit.WC",ADULTS,0,0,&TTNeutralGenesSH::getFitWC,0,0,0);
2153
2154 }
2155}
void setAdultsFstatWeirCockerham()
Definition: ttneutralgenes.h:436
void setOffspringFstatWeirCockerham()
Definition: ttneutralgenes.h:435
double getFitWC()
Definition: ttneutralgenes.h:439
double getFstWC()
Definition: ttneutralgenes.h:437
double getFisWC()
Definition: ttneutralgenes.h:438

References StatHandler< SH >::add(), ADULTS, getFisWC(), getFitWC(), getFstWC(), OFFSPRG, setAdultsFstatWeirCockerham(), and setOffspringFstatWeirCockerham().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setFstatWeirCockerham()

void TTNeutralGenesSH::setFstatWeirCockerham ( age_t  AGE)

Computes the Weir & Cockerham (1984) Fstat values (Theta, F, and f).

772{
773 unsigned int patchNbr = this->_pop->getPatchNbr();
774 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
775 unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
776 //age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
777
778 setAlleleTables(AGE);
779 setHeteroTable(AGE);
780
781 //init
782 double *pop_sizes = new double [patchNbr];
783 double tot_size, inv_ntot;
784 double sum_weights = 0;
785 double nbar, nc, inv_nbar;
786 unsigned int extantPs = 0;
787
788 tot_size = _pop->size(AGE);
789
790 for(unsigned int i = 0; i < patchNbr; i++) {
791 pop_sizes[i] = _pop->size(AGE, i);
792 if(pop_sizes[i]) {
793 extantPs++;
794 sum_weights += (pop_sizes[i] * pop_sizes[i] / tot_size);
795 }
796 }
797 nbar = tot_size/extantPs;
798 nc = (tot_size - sum_weights)/(extantPs-1);
799 inv_nbar = 1/(nbar - 1);
800 inv_ntot = 1/tot_size;
801
802 double var;
803 double s2, pbar, hbar;
804 double s2_denom = 1.0/((extantPs-1)*nbar),
805 r = (double)(extantPs-1)/extantPs,
806 hbar_factor=(2*nbar-1)/(4*nbar);
807 double a = 0, b = 0, c = 0, x;
808
809 for (unsigned int l = 0; l < nb_locus; ++l) {
810
811 for (unsigned int u = 0; u < nb_allele; ++u) {
812
813 s2 = pbar = hbar = 0;
814
815 for (unsigned int i = 0; i < patchNbr; ++i) {
816
817 var = _alleleFreqTable.get(i, l, u) - _globalAlleleFreq.get(l, u); //(p_liu - pbar_u)^2
818
819 var *= var;
820
821 s2 += var * pop_sizes[i];
822
823 hbar += _heteroTable.get(i, l, u);
824
825 }//end for pop
826
827 s2 *= s2_denom;
828 pbar = _globalAlleleFreq.get(l, u);
829 hbar *= inv_ntot;
830
831 x = pbar * (1 - pbar) - r * s2;
832 a += s2 - inv_nbar*( x - 0.25 * hbar);
833 b += x - hbar_factor * hbar;
834 c += hbar;
835
836
837 } // end for allele
838
839 }// end for locus
840
841 a *= nbar/nc;
842 b *= nbar/(nbar - 1);
843 c *= 0.5;
844
845 _fst_WC = a / (a + b + c);
846 _fit_WC = (a + b) / (a + b + c);
847 _fis_WC = b / (b + c);
848
849 delete [] pop_sizes;
850}
void setHeteroTable(age_t AGE)
Definition: stats_fstat.cc:166

References _alleleFreqTable, _fis_WC, _fit_WC, _fst_WC, _globalAlleleFreq, _heteroTable, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, DataTable< T >::get(), TMatrix::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatchNbr(), setAlleleTables(), setHeteroTable(), and Metapop::size().

Referenced by setAdultsFstatWeirCockerham(), and setOffspringFstatWeirCockerham().

+ Here is the caller graph for this function:

◆ setFstatWeirCockerham_MS()

void TTNeutralGenesSH::setFstatWeirCockerham_MS ( age_t  AGE)

Computes W&C F-stats using the Mean Squares approach, similar to the implementation in Hierfstat. This code gives the exact same results as the other version.

855{
858 unsigned int patchNbr = this->_pop->getPatchNbr();
859 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
860 unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
861 //age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
862
863 setAlleleTables(AGE);
864 setHeteroTable(AGE);
865
866 //init
867 double *pop_sizes = new double [patchNbr];
868 double tot_size;
869 double sum_weights = 0;
870 double nc;
871 unsigned int extantPs = 0;
872
873 tot_size = _pop->size(AGE);
874
875 for(unsigned int i = 0; i < patchNbr; i++) {
876 pop_sizes[i] = _pop->size(AGE, i);
877 if(pop_sizes[i]) {
878 extantPs++;
879 sum_weights += (pop_sizes[i] * pop_sizes[i] / tot_size);
880 }
881 }
882
883 nc = (tot_size - sum_weights)/(extantPs-1);
884
885// unsigned int np = extantPs;
886 unsigned int npl = extantPs; //all loci typed in all patches
887
888 //p = _alleleFreqTable
889 //pb = _globalAlleleFreq
890
891 unsigned int *alploc = new unsigned int [nb_locus];
892
893 unsigned int **alploc_table = new unsigned int* [nb_locus];
894
895 for(unsigned int i = 0; i < nb_locus; ++i)
896 alploc_table[i] = new unsigned int[nb_allele];
897
898 unsigned int tot_num_allele = 0;
899
900 for(unsigned int l = 0; l < nb_locus; ++l){
901
902 alploc[l] = 0;
903
904 for(unsigned int cnt, a = 0; a < nb_allele; ++a) {
905
906 cnt=0;
907
908 for(unsigned int i = 0; i < patchNbr; i++) {
909
910 cnt += _alleleCountTable.get(i,l,a);
911
912 }
913 alploc_table[l][a] = (cnt != 0);
914 alploc[l] += (cnt != 0);
915 }
916
917 tot_num_allele += alploc[l];
918 }
919
920 //n, and nal are given by pop_sizes, same num ind typed at all loci in each patch
921 //nc is the same for each locus
922 //nt is given by tot_size, same tot num of ind typed for all loci
923
924 //SSG: het/2 for each allele
925 double *SSG = new double[tot_num_allele];
926 double *SSP = new double[tot_num_allele];
927 double *SSi = new double[tot_num_allele];
928
929 unsigned int all_cntr = 0;
930
931 double het, freq, var;
932
933 for(unsigned int l = 0; l < nb_locus; ++l) {
934
935 for(unsigned int a = 0; a < nb_allele & all_cntr < tot_num_allele; ++a) {
936
937 if(alploc_table[l][a] == 0) continue; //do not consider alleles not present in the pop
938
939 SSG[all_cntr] = 0;
940 SSi[all_cntr] = 0;
941 SSP[all_cntr] = 0;
942
943 for(unsigned int p = 0; p < patchNbr; ++p){
944
945 if(!_pop->size(AGE, p)) continue; //skip empty patches
946
947 het = _heteroTable.get(p, l, a);
948
949 freq = _alleleFreqTable.get(p, l, a);
950
951 var = freq - _globalAlleleFreq.get(l, a); //(p_liu - pbar_u)^2
952
953 var *= var;
954
955 SSG[all_cntr] += het;
956
957 SSi[all_cntr] += 2*pop_sizes[p]*freq*(1-freq) - het/2;
958
959 SSP[all_cntr] += 2*pop_sizes[p]*var;
960 }
961
962 all_cntr++;
963 }
964
965 }
966
967
968 assert(all_cntr == tot_num_allele);
969
970 double *MSG = new double[tot_num_allele];
971 double *MSP = new double[tot_num_allele];
972 double *MSI = new double[tot_num_allele];
973 double *sigw = new double[tot_num_allele];
974 double *siga = new double[tot_num_allele];
975 double *sigb = new double[tot_num_allele];
976
977// double *FST_pal = new double[tot_num_allele];
978// double *FIS_pal = new double[tot_num_allele];
979
980 double SIGA = 0, SIGB = 0, SIGW = 0;
981
982 for(unsigned int i = 0; i < tot_num_allele; ++i){
983
984 MSG[i] = SSG[i] / (2 * tot_size);
985 sigw[i] = MSG[i]; //wasted!
986
987 MSP[i] = SSP[i] / (npl-1);
988
989 MSI[i] = SSi[i]/ (tot_size - npl);
990
991 sigb[i] = 0.5*(MSI[i] - MSG[i]);
992
993 siga[i] = (MSP[i] - MSI[i])/(2*nc);
994
995// FST_pal[i] = siga[i]/(siga[i]+sigb[i]+sigw[i]);
996// FIS_pal[i] = sigb[i]/(sigb[i]+sigw[i]);
997
998 SIGA += siga[i];
999 SIGB += sigb[i];
1000 SIGW += sigw[i];
1001 }
1002
1003 //per locus stats:
1004 if(_fst_WC_loc) delete [] _fst_WC_loc; _fst_WC_loc = new double [nb_locus];
1005 if(_fis_WC_loc) delete [] _fis_WC_loc; _fis_WC_loc = new double [nb_locus];
1006 if(_fit_WC_loc) delete [] _fit_WC_loc; _fit_WC_loc = new double [nb_locus];
1007
1008 double lsiga, lsigb, lsigw;
1009
1010// cout<<" computing sigma per locus\n";
1011
1012 for(unsigned int allcntr = 0, i = 0; i < nb_locus; ++i) {
1013
1014 lsiga = lsigb = lsigw = 0;
1015
1016 for(unsigned int l = 0; l < alploc[i]; ++l) {
1017
1018 lsiga += siga[allcntr];
1019 lsigb += sigb[allcntr];
1020 lsigw += sigw[allcntr];
1021
1022 allcntr++;
1023
1024 }
1025
1026 _fst_WC_loc[i] = lsiga /(lsiga + lsigb + lsigw);
1027 _fis_WC_loc[i] = lsigb /(lsigb + lsigw);
1028 _fit_WC_loc[i] = (lsiga +lsigb) /(lsiga + lsigb + lsigw);
1029
1030 }
1031
1032
1033 // Total F-stats
1034 _fst_WC = SIGA / (SIGA + SIGB + SIGW);
1035 _fit_WC = (SIGA + SIGB) / (SIGA + SIGB + SIGW);
1036 _fis_WC = SIGB / (SIGB + SIGW);
1037
1038 delete[]pop_sizes;
1039 delete[]alploc;
1040 for(unsigned int i = 0; i < nb_locus; ++i)
1041 delete[]alploc_table[i];
1042 delete[]alploc_table;
1043 delete[]SSG;
1044 delete[]SSi;
1045 delete[]SSP;
1046 delete[]MSG;
1047 delete[]MSI;
1048 delete[]MSP;
1049 delete[]sigw;
1050 delete[]siga;
1051 delete[]sigb;
1052
1053}

References _alleleCountTable, _alleleFreqTable, _fis_WC, _fis_WC_loc, _fit_WC, _fit_WC_loc, _fst_WC, _fst_WC_loc, _globalAlleleFreq, _heteroTable, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, DataTable< T >::get(), TMatrix::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatchNbr(), setAlleleTables(), setHeteroTable(), and Metapop::size().

◆ setFstMatrix()

void TTNeutralGenesSH::setFstMatrix ( age_t  AGE,
unsigned char  dim 
)

Computes the weighted within and between patch Fst's as well as the overall Fst (Theta).

The method used here is that of Weir & Hill 2002, Ann. Rev. Genet. 36:721-750. The weighting is done for samples (patches) of unequal sizes.

Parameters
AGEthe age class
dimthe dimension of the matrix to fill:
  • 1 = the diagonal (i.e. the wihtin patch Fst or theta_ii)
  • 2 = the upper half (i.e. the between patch Fst or theta_ii')
  • 3 = both
666{
667 unsigned int patchNbr = this->_pop->getPatchNbr();
668 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
669 unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
670
671 setAlleleTables(AGE);
672
673 if(_fst_matrix == NULL)
674
675 _fst_matrix = new TMatrix(patchNbr, patchNbr);
676
677 else if( _fst_matrix->length() != patchNbr * patchNbr)
678
679 _fst_matrix->reset(patchNbr, patchNbr);
680
681 _fst_matrix->assign(nanf("NULL"));
682
683 //init
684 double *pop_weights = new double[patchNbr];
685 double *pop_sizes = new double[patchNbr];
686 double **numerator = new double*[patchNbr];
687 for(unsigned int i = 0; i < patchNbr; i++) numerator[i] = new double [patchNbr];
688 double tot_size;
689 double numerator_W = 0;
690 double denominator = 0;
691 double sum_weights = 0;
692
693 tot_size = _pop->size(AGE) * 2;
694
695 for(unsigned int i = 0; i < patchNbr; ++i) {
696 pop_sizes[i] = _pop->size(AGE, i) * 2;
697 pop_weights[i] = pop_sizes[i] - (pop_sizes[i] * pop_sizes[i] / tot_size); //n_ic in Weir & Hill 2002
698 sum_weights += pop_weights[i];
699 for(unsigned int j = 0; j < patchNbr; j++)
700 numerator[i][j] = 0;
701 }
702
703 double p, pq, var, num;
704
705 for (unsigned int i = 0; i < patchNbr; ++i) {
706
707 if( !pop_sizes[i] ) continue;
708
709 for (unsigned int l = 0; l < nb_locus; ++l) {
710
711 for (unsigned int u = 0; u < nb_allele; ++u) {
712
713 p = _alleleFreqTable.get(i, l, u); //p_liu
714
715 pq = p * (1 - p);
716
717 var = p - _globalAlleleFreq.get(l, u); //(p_liu - pbar_u)^2
718
719 var *= var;
720
721 num = pq * pop_sizes[i] / (pop_sizes[i] -1);
722
723 numerator[i][i] += num;
724
725 numerator_W += num * pop_sizes[i]; //see equ. 9, Weir & Hill 2002
726
727 denominator += pop_sizes[i] * var + pop_weights[i] * pq; //common denominator
728
729 } // end for allele
730 }// end for locus
731 }//end for pop
732
733 for (unsigned int i = 0; i < patchNbr; ++i) {
734 if( !pop_sizes[i] ) continue;
735 _fst_matrix->set(i, i, 1 - (numerator[i][i] * sum_weights / denominator) );
736 }
737 _fst_WH = 1 - ((numerator_W * sum_weights) / (denominator * tot_size)); //equ. 9 Weir & Hill 2002
738
739 //pairwise Fst:
740 if(dim & 2) {
741 double pi, pj;
742 for (unsigned int l = 0; l < nb_locus; ++l)
743 for (unsigned int u = 0; u < nb_allele; ++u)
744 for (unsigned int i = 0; i < patchNbr - 1; ++i) {
745 if( !pop_sizes[i] ) continue;
746 for (unsigned int j = i + 1; j < patchNbr; ++j) {
747 if( !pop_sizes[j] ) continue;
748 pi = _alleleFreqTable.get(i, l, u);
749 pj = _alleleFreqTable.get(j, l, u);
750 numerator[i][j] += pi * (1 - pj) + pj * (1 - pi); //equ. 7 of Weir & Hill 2002
751 }
752 }
753 for (unsigned int i = 0; i < patchNbr - 1; ++i){
754 if( !pop_sizes[i] ) continue;
755 for (unsigned int j = i + 1; j < patchNbr; ++j){
756 if( !pop_sizes[j] ) continue;
757 _fst_matrix->set(i, j, 1 - ( (numerator[i][j] * sum_weights) / (2 * denominator)) );
758 }
759 }
760 }
761
762
763 delete [] pop_weights;
764 delete [] pop_sizes;
765 for(unsigned int i = 0; i < patchNbr; i++) delete [] numerator[i];
766 delete [] numerator;
767}
void set(unsigned int i, unsigned int j, double val)
Sets element at row i and column j to value val.
Definition: tmatrix.h:102

References _alleleFreqTable, _fst_matrix, _fst_WH, _globalAlleleFreq, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TMatrix::assign(), DataTable< T >::get(), TMatrix::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatchNbr(), TMatrix::length(), TMatrix::reset(), TMatrix::set(), setAlleleTables(), and Metapop::size().

Referenced by setAdultsFstBetween(), setAdultsFstMatrix(), setAdultsFstWithin(), setOffsprgFstBetween(), setOffsprgFstMatrix(), and setOffsprgFstWithin().

+ Here is the caller graph for this function:

◆ setFstMatrixRecorders()

void TTNeutralGenesSH::setFstMatrixRecorders ( age_t  AGE,
unsigned char  dim 
)
2160{
2161 std::ostringstream name, sub_name;
2162
2163 void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ?
2170
2171 const char *prefix = (AGE == ADULTS ? "adlt." : "off.");
2172
2173 unsigned int nbpatch = _pop->getPatchNbr();
2174 unsigned int scale = (unsigned int)pow(10.0, (int)log10((float)nbpatch) + 1);
2175
2176 name<<"Weir&Hill weighted Fst ("<<prefix<<")";
2177 sub_name<< prefix << "fst.WH";
2178 add(name.str(), sub_name.str(), AGE, 0, 0, &TTNeutralGenesSH::getWeightedFst, 0, 0, setter);
2179 name.str("");
2180 sub_name.str("");
2181
2182 if(dim & 1) {
2183 for(unsigned int i = 0; i < nbpatch; ++i) {
2184 name<<"Weighted Fst "<<i+1<<"."<<i+1;
2185 sub_name<< prefix << "fst" << i+1 << "." << i+1;
2186 add(name.str(), sub_name.str(), AGE, i*scale + i, 0, 0, &TTNeutralGenesSH::getFst_ij, 0, 0);
2187 name.str("");
2188 sub_name.str("");
2189 }
2190 }
2191 if(dim & 2){
2192 for(unsigned int i = 0; i < nbpatch; ++i) {
2193 for(unsigned int j = i+1; j < nbpatch; ++j) {
2194 name<<"Weighted Fst "<<i+1<< "." <<j+1;
2195 sub_name<< prefix << "fst" << i+1 << "." << j+1;
2196 add(name.str(), sub_name.str(), AGE, i*scale + j, 0, 0, &TTNeutralGenesSH::getFst_ij, 0, 0);
2197 name.str("");
2198 sub_name.str("");
2199 }
2200 }
2201 }
2202
2203}
double getWeightedFst()
Returns the weighted Fst using Weir & Hill (2002) method.
Definition: ttneutralgenes.h:397
void setAdultsFstBetween()
Definition: ttneutralgenes.h:391
void setOffsprgFstBetween()
Definition: ttneutralgenes.h:394
void setAdultsFstMatrix()
Definition: ttneutralgenes.h:389
void setOffsprgFstWithin()
Definition: ttneutralgenes.h:393
double getFst_ij(unsigned int i)
Accessor to the Fst matrix as set by setFstMatrix().
Definition: ttneutralgenes.h:399
void setAdultsFstWithin()
Definition: ttneutralgenes.h:390
void setOffsprgFstMatrix()
Definition: ttneutralgenes.h:392

References StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, getFst_ij(), Metapop::getPatchNbr(), getWeightedFst(), setAdultsFstBetween(), setAdultsFstMatrix(), setAdultsFstWithin(), setOffsprgFstBetween(), setOffsprgFstMatrix(), and setOffsprgFstWithin().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setHeteroTable()

void TTNeutralGenesSH::setHeteroTable ( age_t  AGE)
167{
168 unsigned int patchNbr = this->_pop->getPatchNbr();
169 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
170 unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
171 Patch* patch;
172 unsigned int isHetero;
173 unsigned char** genes;
174 age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
175
176 unsigned int **sizes;
177
178 sizes = new unsigned int * [patchNbr];
179
180 for(unsigned int i = 0; i < patchNbr; ++i) {
181 sizes[i] = new unsigned int [nb_locus];
182 for(unsigned int j = 0; j < nb_locus; ++j)
183 sizes[i][j] = nb_allele;
184 }
185
186 _heteroTable.update(patchNbr, nb_locus, sizes);
188
189
190 for (unsigned int i = 0; i < patchNbr; ++i) {
191
192 patch = _pop->getPatch(i);
193
194 if(!patch->size(age_pos)) continue;
195
196 for(unsigned int j = 0, size = patch->size(FEM, age_pos); j < size; ++j) {
197
198 genes = (unsigned char**)patch->get(FEM, age_pos, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
199
200 for (unsigned int l = 0; l < nb_locus; ++l) {
201 isHetero = (genes[0][l] != genes[1][l]);
202 _heteroTable.plus(i, l, genes[0][l], isHetero);
203 _heteroTable.plus(i, l, genes[1][l], isHetero);
204 }
205 }
206
207 for(unsigned int j = 0, size = patch->size(MAL, age_pos); j < size; ++j) {
208
209 genes = (unsigned char**)patch->get(MAL, age_pos, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
210
211 for (unsigned int l = 0; l < nb_locus; ++l) {
212 isHetero = (genes[0][l] != genes[1][l]);
213 _heteroTable.plus(i, l, genes[0][l], isHetero);
214 _heteroTable.plus(i, l, genes[1][l], isHetero);
215 }
216 }
217 }
218 for(unsigned int i = 0; i < patchNbr; ++i) delete [] sizes[i];
219 delete [] sizes;
220
221}
void plus(unsigned int group, unsigned int Class, unsigned int elmnt, T val)
Adds 'val' to 'elmnt' in the class 'Class' in the group 'group'.
Definition: datatable.h:240
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

References _heteroTable, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, ADLTx, ADULTS, FEM, Patch::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), DataTable< T >::init(), MAL, OFFSx, DataTable< T >::plus(), Patch::size(), and DataTable< T >::update().

Referenced by setFstatWeirCockerham(), setFstatWeirCockerham_MS(), setHeterozygosity(), and TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ setHeterozygosity()

void TTNeutralGenesSH::setHeterozygosity ( age_t  AGE)
226{
227 unsigned int patchNbr = this->_pop->getPatchNbr();
228 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
229 unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
230 Patch* patch;
231 age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
232
233 setHeteroTable(AGE);
234
235 for (unsigned int i = 0, size; i < patchNbr; ++i) {
236
237 patch = _pop->getPatch(i);
238 size = patch->size(age_pos);
239
240 if(!size) continue;
241
242 for (unsigned int l = 0; l < nb_locus; ++l) {
243 for (unsigned int u = 0; u < nb_allele; ++u) {
244 _heteroTable.divide(i, l , u, size);
245 }
246 }
247 }
248
249
250}
void divide(unsigned int group, unsigned int Class, unsigned int elmnt, T val)
Subdivide 'elmnt' of the class 'Class' in the group 'group' by 'val'.
Definition: datatable.h:252

References _heteroTable, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, ADLTx, ADULTS, DataTable< T >::divide(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatch(), Metapop::getPatchNbr(), OFFSx, setHeteroTable(), and Patch::size().

Referenced by setAdultHeterozygosity(), and setOffspringHeterozygosity().

+ Here is the caller graph for this function:

◆ setHo()

double TTNeutralGenesSH::setHo ( age_idx  age_pos)
364{
365 unsigned int nloc = _SHLinkedTrait->get_locus_num(), nbpatch = _pop->getPatchNbr(), psize;
366 double indloc = 0, hetero = 0;
367 unsigned char** genes;
368 Patch* current_patch;
369
370 for (unsigned int i = 0; i < nbpatch; ++i) {
371
372 current_patch = _pop->getPatch(i);
373
374 psize = current_patch->size(FEM, age_pos);
375
376 indloc += psize;
377
378 for(unsigned int j = 0; j < psize; ++j) {
379
380 genes = (unsigned char**)current_patch->get(FEM, age_pos, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
381
382 for (unsigned int k = 0; k < nloc; ++k) hetero += (genes[0][k] != genes[1][k]);
383 }
384 }
385 //doing some unrolling, loops for the males here:
386 for (unsigned int i = 0; i < nbpatch; ++i) {
387
388 current_patch = _pop->getPatch(i);
389
390 psize = current_patch->size(MAL, age_pos);
391
392 indloc += psize;
393
394 for(unsigned int j = 0; j < psize; ++j) {
395
396 genes = (unsigned char**)current_patch->get(MAL, age_pos, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
397
398 for (unsigned int k = 0; k < nloc; ++k) hetero += (genes[0][k] != genes[1][k]);
399 }
400 }
401
402 indloc *= nloc;
403
404 return (indloc != 0 ? hetero/indloc : 0.0);
405}

References StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, FEM, Patch::get(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), MAL, and Patch::size().

Referenced by setFstat().

+ Here is the caller graph for this function:

◆ setHo2()

deque< double > TTNeutralGenesSH::setHo2 ( age_idx  age_pos)
546{
547 unsigned int nloc = _SHLinkedTrait->get_locus_num(), nbpatch = _pop->getPatchNbr();
548 unsigned int msize, fsize, psize;
549 deque<double> hetero(nloc,0);
550 unsigned char** genes;
551 Patch* current_patch;
552
553 for (unsigned int i = 0; i < nbpatch; ++i) {
554
555 current_patch = _pop->getPatch(i);
556
557 fsize = current_patch->size(FEM, age_pos);
558
559 for(unsigned int j = 0; j < fsize; ++j) {
560
561 genes = (unsigned char**)current_patch->get(FEM, age_pos, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
562
563 for (unsigned int k = 0; k < nloc; ++k) hetero[k] += (genes[0][k] != genes[1][k]);
564 }
565
566 msize = current_patch->size(MAL, age_pos);
567
568 for(unsigned int j = 0; j < msize; ++j) {
569
570 genes = (unsigned char**)current_patch->get(MAL, age_pos, j)->getTrait(_SHLinkedTraitIndex)->get_sequence();
571
572 for (unsigned int k = 0; k < nloc; ++k) hetero[k] += (genes[0][k] != genes[1][k]);
573 }
574
575 }
576
577 psize = _pop->size(age_pos);
578
579 if(psize != 0)
580 for (unsigned int k = 0; k < nloc; ++k)
581 hetero[k] /= psize;
582
583 return hetero;
584}

References StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, FEM, Patch::get(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), MAL, Metapop::size(), and Patch::size().

Referenced by setFstat2(), and TTNeutralGenesFH::write_varcompWC().

+ Here is the caller graph for this function:

◆ setHs()

double TTNeutralGenesSH::setHs ( age_idx  age_pos)
411{
412 unsigned int i, k, x, nb_patch=0;
413 double hs = 0;
414 unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
415 unsigned int patchNbr = _pop->getPatchNbr();
416 deque<double>genehet(nbLoc,1);
417 deque<double>hspop(patchNbr,0);
418 double freq;
419 Patch* current_patch;
420
421 //compute the expected heterozygosity for each patch and return average:
422 for (i = 0; i < patchNbr; ++i) {
423
424 current_patch = _pop->getPatch(i);
425
426 genehet.assign(nbLoc, 1);
427
428 if(current_patch->size(age_pos) != 0) {
429
430 nb_patch++;
431
432 for (k = 0; k < nbLoc; ++k) {
433
434 for (x = 0; x < nbAll; ++x) {
435
436 freq = _alleleFreqTable.get(i, k, x);
437
438 freq *= freq; //squared frequencies (expected _homozygosity)
439
440 genehet[k] -= freq; //1 - sum of p2 = expected heterozygosity
441 }
442 //expected heterozygosity:
443 hspop[i] += genehet[k];
444 }
445 }//end_if size!=0
446
447 hs += hspop[i];
448 }//end patch loop
449
450 return (nb_patch !=0 ? hs/(nbLoc*nb_patch) : 0.0);
451}
void assign(sex_t SEX, age_idx AGE, unsigned int n)
Assigns a new container of given size for the sex and age class passed, sets all values to NULL.
Definition: metapop.h:561

References _alleleFreqTable, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, Patch::assign(), DataTable< T >::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatch(), Metapop::getPatchNbr(), and Patch::size().

Referenced by setFstat().

+ Here is the caller graph for this function:

◆ setHs2()

deque< double > TTNeutralGenesSH::setHs2 ( age_idx  age_pos)
590{
591 unsigned int i, k, x, nb_patch = 0;
592 unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
593 unsigned int patchNbr = _pop->getPatchNbr();
594 deque<double>genehet(nbLoc,1);
595 deque<double>hsloc(nbLoc,0);
596 double freq;
597 Patch* current_patch;
598
599 for (i = 0; i < patchNbr; ++i) {
600 if(_pop->size(age_pos, i) != 0)
601 nb_patch++;
602 }
603
604 //compute the expected heterozygosity for each patch and return average:
605 for (k = 0; k < nbLoc; ++k) {
606
607 for (i = 0; i < patchNbr; ++i) {
608
609 current_patch = _pop->getPatch(i);
610
611 genehet[k] = 1;
612
613 if(current_patch->size(age_pos) != 0) {
614
615 for (x = 0; x < nbAll; ++x) {
616
617 freq = _alleleFreqTable.get(i, k, x);
618
619 freq *= freq; //squared frequencies (expected _homozygosity)
620
621 genehet[k] -= freq; //1 - sum of p2 = expected heterozygosity
622 }
623
624 }//end_if size!=0
625
626 hsloc[k] += genehet[k];
627
628 }//end patch loop
629
630 if(nb_patch !=0) hsloc[k] /= nb_patch;
631
632 }//end loci loop
633
634 return hsloc;
635}

References _alleleFreqTable, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, DataTable< T >::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatch(), Metapop::getPatchNbr(), Metapop::size(), and Patch::size().

Referenced by setFstat2().

+ Here is the caller graph for this function:

◆ setHt()

double TTNeutralGenesSH::setHt ( age_idx  age_pos)
456{
457 double ht = 0;
458 unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
459 deque<double> genehet(nbLoc, 1);
460 double freq;
461
462 //get global allele frequencies per locus
463 for (unsigned int l = 0; l < nbLoc; ++l) {
464 for (unsigned int u = 0; u < nbAll; ++u) {
465
466 freq = _globalAlleleFreq.get(l, u);
467
468 freq *= freq; //squared frequencies
469
470 genehet[l] -= freq; //1 - sum of p2 = expected heterozygosity
471 }
472
473 ht += genehet[l];
474 }
475
476 return (ht/nbLoc);
477}

References _globalAlleleFreq, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TMatrix::get(), TProtoNeutralGenes::get_allele_num(), and TProtoNeutralGenes::get_locus_num().

Referenced by setFstat().

+ Here is the caller graph for this function:

◆ setHt2()

deque< double > TTNeutralGenesSH::setHt2 ( age_idx  age_pos)
640{
641 unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
642 deque<double> genehet(nbLoc, 1);
643 double freq;
644
645 //get global allele frequencies per locus
646 for (unsigned int l = 0; l < nbLoc; ++l) {
647
648 for (unsigned int u = 0; u < nbAll; ++u) {
649
650 freq = _globalAlleleFreq.get(l, u);
651
652 freq *= freq; //squared frequencies
653
654 genehet[l] -= freq; //1 - sum of p2 = expected heterozygosity
655 }
656 }
657
658 return genehet;
659}

References _globalAlleleFreq, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TMatrix::get(), TProtoNeutralGenes::get_allele_num(), and TProtoNeutralGenes::get_locus_num().

Referenced by setFstat2().

+ Here is the caller graph for this function:

◆ setLociDivCounter()

void TTNeutralGenesSH::setLociDivCounter ( age_t  AGE)

Sets the allelic diversity counters.

296{
297 unsigned int i, j, k;
298 unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
299 unsigned int nbpatch = 0, patchNbr = _pop->getPatchNbr();
300 double patch_mean, pop_mean = 0;
301 bool **pop_div;
302
303 //number of alleles per locus, Patch and pop counters:
304 pop_div = new bool * [nbLoc];
305
306 for(i = 0; i < nbLoc; ++i) {
307 pop_div[i] = new bool [nbAll];
308 for(j = 0; j < nbAll;++j)
309 pop_div[i][j] = 0;
310 }
311
312 for(k = 0; k < patchNbr; ++k)
313 nbpatch += (_pop->size(AGE, k) != 0);
314
315 for(k = 0; k < patchNbr; ++k) {
316
317 patch_mean = 0;
318
319 for(i = 0; i < nbLoc; ++i)
320 for(j = 0; j < nbAll;++j) {
321 patch_mean += (_alleleCountTable.get(k,i,j) != 0);
322 pop_div[i][j] |= (_alleleCountTable.get(k,i,j) != 0);
323 }
324 //add mean nb of alleles per locus for Patch k to the pop mean
325 pop_mean += patch_mean/nbLoc;
326 }
327
328 _nb_all_local = (nbpatch ? pop_mean/nbpatch : nanf("NULL"));
329
330 _nb_all_global = 0;
331
332 for(i = 0; i < nbLoc; ++i)
333 for(j = 0; j < nbAll;++j)
334 _nb_all_global += pop_div[i][j];
335
336 _nb_all_global /= nbLoc;
337
338 for(i = 0; i < nbLoc; ++i)
339 delete [] pop_div[i];
340 delete [] pop_div;
341
342 //number of fixed loci, local and global counters:
343 _fix_loc_local = 0;
344
345 for (k = 0; k < patchNbr; ++k)
346 for (i = 0; i < nbLoc; ++i)
347 for (j = 0; j < nbAll; ++j)
348 _fix_loc_local += ( _alleleFreqTable.get(k, i, j) == 1 );
349
350 _fix_loc_local /= nbpatch;
351
352 _fix_loc_global = 0;
353
354 //globally:
355 for (i = 0; i < nbLoc; ++i)
356 for (j = 0; j < nbAll; ++j)
357 _fix_loc_global += ( _globalAlleleFreq.get(i, j) == 1 );
358
359}

References _alleleCountTable, _alleleFreqTable, _fix_loc_global, _fix_loc_local, _globalAlleleFreq, _nb_all_global, _nb_all_local, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, DataTable< T >::get(), TMatrix::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatchNbr(), and Metapop::size().

Referenced by setFstat(), and setFstat2().

+ Here is the caller graph for this function:

◆ setNeiGeneticDistance()

void TTNeutralGenesSH::setNeiGeneticDistance ( age_t  AGE)
1130{
1131 //see: Nei, M. 1978. Genetics 89:583-590
1132 unsigned int patchNbr = this->_pop->getPatchNbr();
1133 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
1134 unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
1135
1136 setAlleleTables(AGE);
1137
1138 if(_D == NULL)
1139
1140 _D = new TMatrix(patchNbr, patchNbr);
1141
1142 else if( _D->length() != patchNbr*patchNbr )
1143
1144 _D->reset(patchNbr, patchNbr);
1145
1146 _D->assign(nanf("NULL"));
1147
1148 double num, denom, p1, p2, sum;
1149 double *sum_square = new double [patchNbr];
1150 double *sample_size = new double [patchNbr];
1151
1152 for(unsigned int i = 0; i < patchNbr; ++i) {
1153
1154 sample_size[i] = 2 * _pop->size(AGE, i);
1155
1156 sum_square[i] = 0;
1157
1158 if( !sample_size[i] ) continue;
1159
1160 for(unsigned int l = 0; l < nb_locus; ++l) {
1161
1162 sum = 0;
1163
1164 for(unsigned int u = 0; u < nb_allele ; ++u) {
1165 p1 = _alleleFreqTable.get(i, l, u);
1166 sum += p1*p1;
1167 }
1168
1169 sum_square[i] += (sample_size[i] * sum -1) / (sample_size[i] -1); //unbiased estimate, equ. 6 in Nei 1978 (Genetics)
1170 }
1171 }
1172
1173 unsigned int pairs = 0;
1174 double Dpair;
1175 _meanD = 0;
1176
1177 for(unsigned int i = 0; i < patchNbr-1; ++i) {
1178
1179 if( !sample_size[i] ) continue;
1180
1181 for(unsigned int j = i + 1; j < patchNbr; ++j) {
1182
1183 if( !sample_size[j] ) continue;
1184
1185 num = 0;
1186 pairs++;
1187
1188 for(unsigned int l = 0; l < nb_locus; ++l) {
1189 for(unsigned int u = 0; u < nb_allele ; ++u) {
1190 p1 = _alleleFreqTable.get(i, l, u);
1191 p2 = _alleleFreqTable.get(j, l, u);
1192 num += p1*p2;
1193 }
1194 }
1195
1196 denom = sqrt(sum_square[i] * sum_square[j]);
1197 Dpair = -log(num/denom);
1198 _meanD += Dpair;
1199 _D->set(i, j, Dpair);
1200 }
1201 }
1202 _meanD /= pairs;
1203
1204 delete [] sum_square;
1205 delete [] sample_size;
1206}

References _alleleFreqTable, _D, _meanD, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TMatrix::assign(), DataTable< T >::get(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatchNbr(), TMatrix::length(), TMatrix::reset(), TMatrix::set(), setAlleleTables(), and Metapop::size().

Referenced by setAdltNeiGeneticDistance(), and setOffsprgNeiGeneticDistance().

+ Here is the caller graph for this function:

◆ setNeiGeneticDistanceRecorders()

void TTNeutralGenesSH::setNeiGeneticDistanceRecorders ( age_t  AGE,
bool  pairwise 
)
2208{
2209 string name, sub_name;
2210 unsigned int nbpatch = _pop->getPatchNbr();
2211 string prefix = (AGE == ADULTS ? "adlt." : "off.");
2212 unsigned int scale = (unsigned int)pow(10.0, (int)log10((float)nbpatch) + 1);
2213
2214 void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ?
2217
2218 sub_name = prefix + "D";
2219 add("Average between pop Nei's D", sub_name, AGE, 0, 0,
2221
2222 if(pairwise) {
2223
2224 for(unsigned int i = 0; i < nbpatch -1; i++){
2225 for(unsigned int j = i+1; j < nbpatch; j++) {
2226 name = "Nei's D between pop" + tstring::int2str(i+1) + " and pop" + tstring::int2str(j+1);
2227 sub_name = prefix + "D.p" + tstring::int2str(i+1) + "p" + tstring::int2str(j+1);
2228 add(name, sub_name, AGE, i*scale + j, 0, 0, &TTNeutralGenesSH::getNeiGeneticDistance, 0, 0);
2229 }
2230 }
2231 }
2232}
double getNeiGeneticDistance(unsigned int i)
Definition: ttneutralgenes.h:501
void setOffsprgNeiGeneticDistance()
Definition: ttneutralgenes.h:499
void setAdltNeiGeneticDistance()
Definition: ttneutralgenes.h:498
double getMeanNeiGeneticDistance()
Definition: ttneutralgenes.h:506

References StatHandlerBase::_pop, StatHandler< SH >::add(), ADULTS, getMeanNeiGeneticDistance(), getNeiGeneticDistance(), Metapop::getPatchNbr(), tstring::int2str(), setAdltNeiGeneticDistance(), and setOffsprgNeiGeneticDistance().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgCoaBetween()

void TTNeutralGenesSH::setOffsprgCoaBetween ( )
inline
470{setCoaMatrix(OFFSx, 2);}

References OFFSx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgCoaMatrix()

void TTNeutralGenesSH::setOffsprgCoaMatrix ( )
inline
466{setCoaMatrix(OFFSx, 3);}

References OFFSx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgCoaWithin()

void TTNeutralGenesSH::setOffsprgCoaWithin ( )
inline
468{setCoaMatrix(OFFSx, 1);}

References OFFSx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgFstat()

void TTNeutralGenesSH::setOffsprgFstat ( )
inline

References OFFSPRG, and setFstat().

Referenced by setFstatRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgFstat2()

void TTNeutralGenesSH::setOffsprgFstat2 ( )
inline

References OFFSPRG, and setFstat2().

Referenced by setFstat2Recorders().

+ Here is the caller graph for this function:

◆ setOffsprgFstBetween()

void TTNeutralGenesSH::setOffsprgFstBetween ( )
inline

References OFFSPRG, and setFstMatrix().

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgFstMatrix()

void TTNeutralGenesSH::setOffsprgFstMatrix ( )
inline

References OFFSPRG, and setFstMatrix().

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgFstWithin()

void TTNeutralGenesSH::setOffsprgFstWithin ( )
inline

References OFFSPRG, and setFstMatrix().

Referenced by setFstMatrixRecorders().

+ Here is the caller graph for this function:

◆ setOffsprgNeiGeneticDistance()

void TTNeutralGenesSH::setOffsprgNeiGeneticDistance ( )
inline

References OFFSPRG, and setNeiGeneticDistance().

Referenced by setNeiGeneticDistanceRecorders().

+ Here is the caller graph for this function:

◆ setOffspringAlleleFreq()

void TTNeutralGenesSH::setOffspringAlleleFreq ( )
inline

References OFFSPRG, and setAlleleTables().

Referenced by setFreqRecorders().

+ Here is the caller graph for this function:

◆ setOffspringFstatWeirCockerham()

void TTNeutralGenesSH::setOffspringFstatWeirCockerham ( )
inline

References OFFSPRG, and setFstatWeirCockerham().

Referenced by setFstatWCRecorders().

+ Here is the caller graph for this function:

◆ setOffspringHeterozygosity()

void TTNeutralGenesSH::setOffspringHeterozygosity ( )
inline

References OFFSPRG, and setHeterozygosity().

Referenced by setFreqRecorders().

+ Here is the caller graph for this function:

◆ setSibCoa()

void TTNeutralGenesSH::setSibCoa ( Individual I1,
Individual I2 
)
351{
352 unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
353 //non sibs [3]
354 if((I1->getMotherID() != I2->getMotherID()) && (I1->getFatherID() != I2->getFatherID())) {
355 _sib_prop[3]++;
358 nb_locus);
359 }
360 //maternal half sibs [2]
361 else if((I1->getMotherID() == I2->getMotherID()) && (I1->getFatherID() != I2->getFatherID())) {
362 _sib_prop[2]++;
365 nb_locus);
366 }
367 //paternal half sibs [1]
368 else if((I1->getMotherID() != I2->getMotherID()) && (I1->getFatherID() == I2->getFatherID())) {
369 _sib_prop[1]++;
372 nb_locus);
373 }
374 //full sibs [0]
375 else if((I1->getMotherID() == I2->getMotherID()) && (I1->getFatherID() == I2->getFatherID())) {
376 _sib_prop[0]++;
379 nb_locus);
380 }
381}
unsigned long getMotherID()
Definition: individual.h:125
unsigned long getFatherID()
Definition: individual.h:124

References TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, _sib_coa, _sib_prop, Coancestry(), TProtoNeutralGenes::get_locus_num(), TTrait::get_sequence(), Individual::getFatherID(), Individual::getMotherID(), and Individual::getTrait().

Referenced by setSibStats().

+ Here is the caller graph for this function:

◆ setSibStats()

void TTNeutralGenesSH::setSibStats ( )
280{
281 unsigned int i,j,k;
282 unsigned int patchNbr = this->_pop->getPatchNbr(), Fsize, Msize;
283 Patch* current_patch;
284 Individual *I1, *I2;
285
286 for(i = 0;i < 4; ++i) {
287 _sib_prop[i] = 0.0;
288 _sib_coa[i] = 0.0;
289 }
290
291 for( i = 0; i < patchNbr; ++i) {
292
293 current_patch = _pop->getPatch(i);
294
295 if ( (Fsize = current_patch->size(FEM, OFFSx)) != 0) {
296
297 for(j = 0; j < Fsize -1; ++j) {
298
299 I1 = current_patch->get(FEM, OFFSx, j);
300
301 for(k = j+1; k < Fsize; ++k) {
302
303 I2 = current_patch->get(FEM, OFFSx, k);
304
305 setSibCoa(I1, I2);
306 }
307 }
308 }
309
310 if ( (Msize = current_patch->size(MAL, OFFSx)) != 0) {
311
312 for(j = 0; j < Msize -1; ++j) {
313
314 I1 = current_patch->get(MAL, OFFSx, j);
315
316 for(k = j+1; k < Msize; ++k) {
317
318 I2 = current_patch->get(MAL, OFFSx, k);
319
320 setSibCoa(I1, I2);
321 }
322 }
323 }
324
325 //male-female
326 for(j = 0; j < Msize; ++j) {
327
328 I1 = current_patch->get(MAL, OFFSx, j);
329
330 for(k = 0; k < Fsize; ++k) {
331
332 I2 = current_patch->get(FEM, OFFSx, k);
333
334 setSibCoa(I1, I2);
335 }
336 }
337
338 }
339
340 double tot = _sib_prop[0] + _sib_prop[1] + _sib_prop[2] + _sib_prop[3];
341
342 for(i = 0 ; i < 4; ++i) {
343 _sib_coa[i] = ( (_sib_prop[i] != 0) ? _sib_coa[i]/_sib_prop[i] : nanf("NULL"));
344 _sib_prop[i] /= tot;
345 }
346}
This class contains traits along with other individual information (sex, pedigree,...
Definition: individual.h:49
void setSibCoa(Individual *I1, Individual *I2)
Definition: stats_coa.cc:350

References StatHandlerBase::_pop, _sib_coa, _sib_prop, FEM, Patch::get(), Metapop::getPatch(), Metapop::getPatchNbr(), MAL, OFFSx, setSibCoa(), and Patch::size().

Referenced by setStatRecorders().

+ Here is the caller graph for this function:

◆ setStatRecorders()

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

Implements StatHandlerBase.

1743{
1744#ifdef _DEBUG_
1745 message("-TTNeutralGenesSH::setStatRecorders ");
1746#endif
1747 if(token.compare("coa") == 0) {
1748
1749 add("Wtn Patch Coancestry (offsprg)","off.theta",OFFSPRG,0,0,
1751 add("Btn Patch Coancestry (offsprg)","off.alpha",OFFSPRG,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
1752
1753 add("Wtn Patch Coancestry (adult)","adlt.theta",ADULTS,0,0,
1755 add("Btn Patch Coancestry (adult)","adlt.alpha",ADULTS,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
1756
1757 } else if(token.compare("adlt.coa") == 0) {
1758
1759 add("Wtn Patch Coancestry (adult)","adlt.theta",ADULTS,0,0,
1761 add("Btn Patch Coancestry (adult)","adlt.alpha",ADULTS,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
1762
1763 } else if(token.compare("off.coa") == 0) {
1764
1765 add("Wtn Patch Coancestry (offsprg)","off.theta",OFFSPRG,0,0,
1767 add("Btn Patch Coancestry (offsprg)","off.alpha",OFFSPRG,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
1768
1769 } else if(token.compare("adlt.coa.persex") == 0) {
1770
1771 add("Female Theta (adult)","adlt.theta",ADULTS,0,0,
1773 add("Female Theta (adult)","adlt.thetaFF",ADULTS,0,0,&TTNeutralGenesSH::getTheta_FF,0,0,0);
1774 add("Male Theta (adult)","adlt.thetaMM",ADULTS,0,0,&TTNeutralGenesSH::getTheta_MM,0,0,0);
1775 add("Female-Male Theta (adult)","adlt.thetaFM",ADULTS,0,0,&TTNeutralGenesSH::getTheta_FM,0,0,0);
1776
1777 } else if(token.compare("adlt.coa.within") == 0) {
1778
1779 add("Wtn Patch Coancestry (adult)","adlt.theta",ADULTS,0,0,
1781
1782 } else if(token.compare("off.coa.within") == 0) {
1783
1784 add("Wtn Patch Coancestry (offsprg)","off.theta",OFFSPRG,0,0,
1786
1787 } else if(token.compare("adlt.coa.between") == 0) {
1788
1789 add("Btn Patch Coancestry (adult)","adlt.alpha",ADULTS,0,0,
1791
1792 } else if(token.compare("off.coa.between") == 0) {
1793
1794 add("Btn Patch Coancestry (offsprg)","off.alpha",OFFSPRG,0,0,
1796
1797 } else if(token.compare("coa.matrix") == 0) {
1798
1801
1802 } else if(token.compare("off.coa.matrix") == 0) {
1803
1805
1806 } else if(token.compare("adlt.coa.matrix") == 0) {
1807
1809
1810 } else if(token.compare("coa.matrix.within") == 0) {
1811
1814
1815 } else if(token.compare("off.coa.matrix.within") == 0) {
1816
1818
1819 } else if(token.compare("adlt.coa.matrix.within") == 0) {
1820
1822
1823 } else if(token.compare("sibcoa") == 0) {
1824
1825 add("Proportion of full-sib offspring","prop.fsib",OFFSPRG,0,0,
1827 add("Proportion of paternal half-sib ","prop.phsib",OFFSPRG,1,0,0,&TTNeutralGenesSH::getSibProportions,0,0);
1828 add("Proportion of maternal half-sib ","prop.mhsib",OFFSPRG,2,0,0, &TTNeutralGenesSH::getSibProportions,0,0);
1829 add("Proportion of non-sib offspring ","prop.nsib",OFFSPRG,3,0,0,&TTNeutralGenesSH::getSibProportions,0,0);
1830 add("Coancestry of full-sib offspring","coa.fsib",OFFSPRG,0,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
1831 add("Coancestry of paternal half-sib ","coa.phsib",OFFSPRG,1,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
1832 add("Coancestry of maternal half-sib ","coa.mhsib",OFFSPRG,2,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
1833 add("Coancestry of non-sib offspring ","coa.nsib",OFFSPRG,3,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
1834
1835 } else if (token == "ntrl.freq") {
1836
1839
1840 } else if (token == "off.ntrl.freq") {
1841
1843
1844 } else if (token == "adlt.ntrl.freq") {
1845
1847
1848 // } else if (token == "ntrl.freq.patch") {
1849 //
1850 // setFreqRecordersPerPatch(ALL);
1851 //
1852 // } else if (token == "off.ntrl.freq.patch") {
1853 //
1854 // setFreqRecordersPerPatch(OFFSPRG);
1855 //
1856 // } else if (token == "adlt.ntrl.freq.patch") {
1857 //
1858 // setFreqRecordersPerPatch(ADULTS);
1859
1860 } else if(token.compare("offsprgfstat") == 0 || token.compare("off.fstat") == 0) {
1861
1863
1864 } else if(token.compare("adultfstat") == 0 || token.compare("adlt.fstat") == 0) {
1865
1867
1868 } else if(token.compare("fstat") == 0) {
1869
1871
1872 } else if(token == "off.fstat2") {
1873
1875
1876 } else if(token == "adlt.fstat2") {
1877
1879
1880 } else if(token == "fstat2") {
1881
1883
1884 } else if(token.compare("fstWC") == 0 || token.compare("fstatWC") == 0) {
1885
1887
1888 } else if(token.compare("off.fstWC") == 0 || token.compare("off.fstatWC") == 0) {
1889
1891
1892 } else if(token.compare("adlt.fstWC") == 0 || token.compare("adlt.fstatWC") == 0) {
1893
1895
1896 } else if(token.compare("weighted.fst") == 0) {
1897
1900
1901 } else if(token.compare("off.weighted.fst") == 0) {
1902
1904
1905 } else if(token.compare("adlt.weighted.fst") == 0) {
1906
1908
1909 } else if(token.compare("weighted.fst.matrix") == 0) {
1910
1913
1914 } else if(token.compare("off.weighted.fst.matrix") == 0) {
1915
1917
1918 } else if(token.compare("adlt.weighted.fst.matrix") == 0) {
1919
1921
1922 } else if(token.compare("weighted.fst.within") == 0) {
1923
1926
1927 } else if(token.compare("off.weighted.fst.within") == 0) {
1928
1930
1931 } else if(token.compare("adlt.weighted.fst.within") == 0) {
1932
1934
1935 } else if(token.compare("adlt.NeiDistance") == 0) {
1936
1938
1939 } else if(token.compare("off.NeiDistance") == 0) {
1940
1942
1943 } else if(token.compare("NeiDistance") == 0) {
1944
1947
1948 } else if(token.compare("mean.NeiDistance") == 0) {
1949
1952
1953 } else if(token.compare("adlt.mean.NeiDistance") == 0) {
1954
1956
1957 } else if(token.compare("off.mean.NeiDistance") == 0) {
1958
1960
1961 } else if(token.compare("Dxy") == 0) {
1962
1963 setDxyRecorders(ADULTS, false);
1964 setDxyRecorders(OFFSPRG, false);
1965
1966 } else if(token.compare("off.Dxy") == 0) {
1967
1968 setDxyRecorders(OFFSPRG, false);
1969
1970 } else if(token.compare("adlt.Dxy") == 0) {
1971
1972 setDxyRecorders(ADULTS, false);
1973
1974 } else if(token.compare("Dxy.patch") == 0) {
1975
1976 setDxyRecorders(ADULTS, true);
1977 setDxyRecorders(OFFSPRG, true);
1978
1979 } else if(token.compare("off.Dxy.patch") == 0) {
1980
1981 setDxyRecorders(OFFSPRG, true);
1982
1983 } else if(token.compare("adlt.Dxy.patch") == 0) {
1984
1985 setDxyRecorders(ADULTS, true);
1986
1987 } else
1988 return false;
1989
1990 return true;
1991}
double getSibCoaMeans(unsigned int i)
Definition: ttneutralgenes.h:493
double getSibProportions(unsigned int i)
Definition: ttneutralgenes.h:492
void setSibStats()
Definition: stats_coa.cc:279
void setDxyRecorders(age_t AGE, bool patchwise)
Definition: ttneutralgenes.cc:2237
void setCoaMatrixRecorders(age_t AGE, unsigned char dim)
Definition: ttneutralgenes.cc:1995
void setFstatRecorders(age_t AGE)
Definition: ttneutralgenes.cc:2083
void setFreqRecorders(age_t AGE)
Definition: ttneutralgenes.cc:2051
double getTheta_FF()
Gives the mean within females coancestry coefficient.
Definition: ttneutralgenes.h:485
double getTheta_MM()
Gives the mean within males coancestry coefficient.
Definition: ttneutralgenes.h:487
double getTheta_FM()
Gives the mean between males and females coancestry coefficient.
Definition: ttneutralgenes.h:489
void setFstatWCRecorders(age_t AGE)
Definition: ttneutralgenes.cc:2139
void setFstMatrixRecorders(age_t AGE, unsigned char dim)
Definition: ttneutralgenes.cc:2159
void setNeiGeneticDistanceRecorders(age_t AGE, bool pairwise)
Definition: ttneutralgenes.cc:2207
void setAdults_Theta()
Definition: stats_coa.cc:203
void setFstat2Recorders(age_t AGE)
Definition: ttneutralgenes.cc:2117
void message(const char *message,...)
Definition: output.cc:40
#define ALL
All ages age class flag.
Definition: types.h:56

References StatHandler< SH >::add(), ADULTS, ALL, getMeanAlpha(), getMeanTheta(), getSibCoaMeans(), getSibProportions(), getTheta_FF(), getTheta_FM(), getTheta_MM(), message(), OFFSPRG, setAdults_Theta(), setAdultsCoaBetween(), setAdultsCoaMatrix(), setAdultsCoaWithin(), setCoaMatrixRecorders(), setDxyRecorders(), setFreqRecorders(), setFstat2Recorders(), setFstatRecorders(), setFstatWCRecorders(), setFstMatrixRecorders(), setNeiGeneticDistanceRecorders(), setOffsprgCoaBetween(), setOffsprgCoaMatrix(), setOffsprgCoaWithin(), and setSibStats().

Member Data Documentation

◆ _alleleCountTable

DataTable< unsigned int > TTNeutralGenesSH::_alleleCountTable
private

◆ _alleleFreqTable

◆ _coa_matrix

TMatrix* TTNeutralGenesSH::_coa_matrix
private

◆ _D

TMatrix* TTNeutralGenesSH::_D
private

◆ _fis

double TTNeutralGenesSH::_fis
private

Referenced by getFis(), setFstat(), and setFstat2().

◆ _fis_WC

double TTNeutralGenesSH::_fis_WC
private

◆ _fis_WC_loc

double * TTNeutralGenesSH::_fis_WC_loc
private

◆ _fit

double TTNeutralGenesSH::_fit
private

Referenced by getFit(), setFstat(), and setFstat2().

◆ _fit_WC

double TTNeutralGenesSH::_fit_WC
private

◆ _fit_WC_loc

double * TTNeutralGenesSH::_fit_WC_loc
private

◆ _fix_loc_global

double TTNeutralGenesSH::_fix_loc_global
private

◆ _fix_loc_local

double TTNeutralGenesSH::_fix_loc_local
private

◆ _fst

double TTNeutralGenesSH::_fst
private

Referenced by getFst(), setFstat(), and setFstat2().

◆ _fst_matrix

TMatrix* TTNeutralGenesSH::_fst_matrix
private

Pairwise Fst matrix.

Referenced by getFst_ij(), setFstMatrix(), and ~TTNeutralGenesSH().

◆ _fst_W1

double TTNeutralGenesSH::_fst_W1
private

◆ _fst_W2

double TTNeutralGenesSH::_fst_W2
private

◆ _fst_WC

double TTNeutralGenesSH::_fst_WC
private

Weir & Cockerham (1984) F-stat estimates.

Referenced by getFstWC(), setFstatWeirCockerham(), and setFstatWeirCockerham_MS().

◆ _fst_WC_loc

double* TTNeutralGenesSH::_fst_WC_loc
private

Per-locus F-stats (Weir&Cockerham).

Referenced by setFstatWeirCockerham_MS(), and ~TTNeutralGenesSH().

◆ _fst_WH

double TTNeutralGenesSH::_fst_WH
private

Weir & Hill (2002) F-stat estimates.

Referenced by getWeightedFst(), and setFstMatrix().

◆ _globalAlleleFreq

◆ _heteroTable

◆ _ho

double TTNeutralGenesSH::_ho
private

F-statistics.

Referenced by getHo(), setFstat(), and setFstat2().

◆ _hs

double TTNeutralGenesSH::_hs
private

Referenced by getHs(), and setFstat().

◆ _hsnei

double TTNeutralGenesSH::_hsnei
private

Referenced by getHsnei(), setFstat(), and setFstat2().

◆ _ht

double TTNeutralGenesSH::_ht
private

Referenced by getHt(), and setFstat().

◆ _htnei

double TTNeutralGenesSH::_htnei
private

Referenced by getHtnei(), setFstat(), and setFstat2().

◆ _mean_alpha

double TTNeutralGenesSH::_mean_alpha
private

Referenced by getMeanAlpha(), and setCoaMatrix().

◆ _mean_theta

double TTNeutralGenesSH::_mean_theta
private

◆ _meanD

double TTNeutralGenesSH::_meanD
private

◆ _nb_all_global

double TTNeutralGenesSH::_nb_all_global
private

◆ _nb_all_local

double TTNeutralGenesSH::_nb_all_local
private

Referenced by getNbAllLocal(), and setLociDivCounter().

◆ _sib_coa

double TTNeutralGenesSH::_sib_coa[4]
private

◆ _sib_prop

double TTNeutralGenesSH::_sib_prop[4]
private

Kinship classes proportions.

Referenced by getSibProportions(), setSibCoa(), and setSibStats().

◆ _table_set_age

unsigned int TTNeutralGenesSH::_table_set_age
private

Referenced by allocateTables(), and setAlleleTables().

◆ _table_set_gen

unsigned int TTNeutralGenesSH::_table_set_gen
private

Referenced by allocateTables(), and setAlleleTables().

◆ _table_set_repl

unsigned int TTNeutralGenesSH::_table_set_repl
private

Referenced by allocateTables(), and setAlleleTables().

◆ Theta_FF

double TTNeutralGenesSH::Theta_FF
private

Referenced by getTheta_FF(), and setAdults_Theta().

◆ Theta_FM

double TTNeutralGenesSH::Theta_FM
private

Referenced by getTheta_FM(), and setAdults_Theta().

◆ Theta_MM

double TTNeutralGenesSH::Theta_MM
private

Referenced by getTheta_MM(), and setAdults_Theta().


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