Nemo  2.4.0b
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 (const TTrait *ind1, const TTrait *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 ~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
 
bool _is_diallelic_bitstring
 
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
361  _fit_WC_loc(0), _fst_matrix(0), _D(0)
362  { }
double * _fis_WC_loc
Definition: ttneutralgenes.h:346
TMatrix * _fst_matrix
Pairwise Fst matrix.
Definition: ttneutralgenes.h:350
unsigned int _table_set_gen
Definition: ttneutralgenes.h:327
TMatrix * _D
Definition: ttneutralgenes.h:353
double * _fst_WC_loc
Per-locus F-stats (Weir&Cockerham).
Definition: ttneutralgenes.h:346
double * _fit_WC_loc
Definition: ttneutralgenes.h:346
TMatrix * _coa_matrix
Definition: ttneutralgenes.h:332
bool _is_diallelic_bitstring
Definition: ttneutralgenes.h:328
unsigned int _table_set_age
Definition: ttneutralgenes.h:327
unsigned int _table_set_repl
Definition: ttneutralgenes.h:327

◆ ~TTNeutralGenesSH()

virtual TTNeutralGenesSH::~TTNeutralGenesSH ( )
inlinevirtual
365  {
366  if(_coa_matrix != NULL) delete _coa_matrix;
367  if(_fst_matrix != NULL) delete _fst_matrix;
368  if(_fst_WC_loc) delete[]_fst_WC_loc;
369  if(_fis_WC_loc) delete[]_fis_WC_loc;
370  if(_fit_WC_loc) delete[]_fit_WC_loc;
371  if(_D != NULL) delete _D;
372  }

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 
)
46 {
47  unsigned int nb_patch = _pop->getPatchNbr();
48  unsigned int **sizes;
49 
50  sizes = new unsigned int * [nb_patch];
51 
52  for(unsigned int i = 0; i < nb_patch; ++i) {
53  sizes[i] = new unsigned int [loci];
54  for(unsigned int j = 0; j < loci; ++j)
55  sizes[i][j] = all;
56  }
57 
58  _alleleCountTable.allocate(nb_patch, loci, sizes);
59 
60  _alleleFreqTable.allocate(nb_patch, loci, sizes);
61 
62  _globalAlleleFreq.reset(loci, all);
63 
64  for(unsigned int i = 0; i < nb_patch; ++i)
65  delete [] sizes[i];
66  delete [] sizes;
67 
68  //reset the time info, to force recalculation after allocate
70  _table_set_gen = 0;
71  _table_set_repl = 0;
72 }
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 and all elements to 0.
Definition: tmatrix.h:161
DataTable< unsigned int > _alleleCountTable
Definition: ttneutralgenes.h:323
TMatrix _globalAlleleFreq
Definition: ttneutralgenes.h:326
DataTable< double > _alleleFreqTable
Definition: ttneutralgenes.h:324
#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().

◆ Coancestry()

double TTNeutralGenesSH::Coancestry ( const TTrait ind1,
const TTrait 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 
47  for (unsigned int k = 0; k < nb_locus; ++k)
48  p += (ind1->get_allele(k,0)==ind2->get_allele(k,0))
49  + (ind1->get_allele(k,0)==ind2->get_allele(k,1))
50  + (ind1->get_allele(k,1)==ind2->get_allele(k,0))
51  + (ind1->get_allele(k,1)==ind2->get_allele(k,1));
52 
53  return (double)p/(4.0*nb_locus);
54 }
virtual unsigned int get_allele(int loc, int all) const =0
Called to read the allele identity at a locus.

References TTrait::get_allele().

◆ getAlleleCountTable()

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

References _alleleCountTable.

Referenced by TTNeutralGenesFH::write_varcompWC().

◆ getAlleleFreqTable()

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

Accessor to the table of allele frequencies, per patch.

415 {return &_alleleFreqTable;}

References _alleleFreqTable.

Referenced by TTNeutralGenesFH::write_varcompWC().

◆ 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.
523  {
524  unsigned int scale = (unsigned int)pow( 10.0, (int)log10((float)_coa_matrix->getNbCols()) + 1 );
525  return _coa_matrix->get(i/scale, i%scale);
526  }
unsigned int getNbCols() const
Gives the number of columns.
Definition: tmatrix.h:215
double get(unsigned int i, unsigned int j) const
Accessor to element at row i and column j.
Definition: tmatrix.h:193

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

Referenced by setCoaMatrixRecorders().

◆ getDxy()

double TTNeutralGenesSH::getDxy ( unsigned int  age_class)
1270 {
1271  double D = 0;
1272  unsigned int p = 0;
1273 
1274  age_t AGE = (static_cast<age_idx>(age_class) == OFFSx ? OFFSPRG : ADULTS);
1275 
1276  for (unsigned int p1 = 0; p1 < _pop->getPatchNbr(); ++p1) {
1277 
1278  if( !_pop->size( AGE , p1) ) continue;
1279 
1280  for (unsigned int p2 = p1 + 1; p2 < _pop->getPatchNbr(); ++p2) {
1281 
1282  if( !_pop->size( AGE , p2) ) continue;
1283 
1284  D += getDxyPerPatch(static_cast<age_idx>(age_class), p1, p2);
1285 
1286  p++;
1287  }
1288  }
1289 
1290  return (p != 0 ? D/p : nanf("NULL"));
1291 }
unsigned int size()
Get the total number of individuals present in the population, all sex and age classes together.
Definition: metapop.h:312
double getDxyPerPatch(age_idx age, unsigned int patch1, unsigned patch2)
Definition: stats_fstat.cc:1295
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().

◆ getDxyAdultPerPatch()

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

References ADLTx, and getDxyPerPatch().

Referenced by setDxyRecorders().

◆ getDxyOffspringPerPatch()

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

References getDxyPerPatch(), and OFFSx.

Referenced by setDxyRecorders().

◆ getDxyPerPatch()

double TTNeutralGenesSH::getDxyPerPatch ( age_idx  age,
unsigned int  patch1,
unsigned  patch2 
)
1296 {
1297  double D = 0;
1298  Patch* patch_1 = _pop->getPatchPtr(patch1);
1299  Patch* patch_2 = _pop->getPatchPtr(patch2);
1300  unsigned int size_1 = patch_1->size(FEM, age);
1301  unsigned int size_2 = patch_2->size(FEM, age);
1302  unsigned int N = 0;//4*size_1 * size_2;// (2N)^2, tot num of haplotype comparisons
1303 
1304 
1305  unsigned int num_loc = _SHLinkedTrait->get_locus_num();
1306 
1307  TTrait *trait_1, *trait_2;
1308 
1309  //females:
1310  for (unsigned int i = 0; i < size_1; ++i) {
1311 
1312  trait_1 = patch_1->get(FEM, age, i)->getTrait(_SHLinkedTraitIndex);
1313 
1314  for (unsigned int j = 0; j < size_2; ++j) {
1315 
1316  trait_2 = patch_2->get(FEM, age, j)->getTrait(_SHLinkedTraitIndex);
1317 
1318  for (unsigned int l = 0; l < num_loc; ++l) {
1319  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 0));
1320  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 1));
1321  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 1));
1322  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 0));
1323  }
1324  N += 4;
1325  }
1326 
1327  for (unsigned int j = 0; j < patch_2->size(MAL, age); ++j) {
1328 
1329  trait_2 = patch_2->get(MAL, age, j)->getTrait(_SHLinkedTraitIndex);
1330 
1331  for (unsigned int l = 0; l < num_loc; ++l) {
1332  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 0));
1333  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 1));
1334  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 1));
1335  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 0));
1336  }
1337  N += 4;
1338 
1339  }
1340  }
1341 
1342  //same for the males:
1343  size_1 = patch_1->size(MAL, age);
1344  size_2 = patch_2->size(MAL, age);
1345 
1346 // N += 4 * size_1 * size_2;
1347 
1348 // if( N == 0) return 0; //needless to go further
1349 
1350  for (unsigned int i = 0; i < size_1; ++i) {
1351 
1352  trait_1 = patch_1->get(MAL, age, i)->getTrait(_SHLinkedTraitIndex);
1353 
1354  for (unsigned int j = 0; j < size_2; ++j) {
1355 
1356  trait_2 = patch_2->get(MAL, age, j)->getTrait(_SHLinkedTraitIndex);
1357 
1358  for (unsigned int l = 0; l < num_loc; ++l) {
1359  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 0));
1360  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 1));
1361  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 1));
1362  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 0));
1363  }
1364  N += 4;
1365  }
1366 
1367  for (unsigned int j = 0; j < patch_2->size(FEM, age); ++j) {
1368 
1369  trait_2 = patch_2->get(FEM, age, j)->getTrait(_SHLinkedTraitIndex);
1370 
1371  for (unsigned int l = 0; l < num_loc; ++l) {
1372  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 0));
1373  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 1));
1374  D += (trait_1->get_allele(l, 0) != trait_2->get_allele(l, 1));
1375  D += (trait_1->get_allele(l, 1) != trait_2->get_allele(l, 0));
1376  }
1377  N += 4;
1378 
1379  }
1380  }
1381 
1382  if( N == 0) return 0;
1383 
1384  return D/N;
1385 }
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:432
unsigned int size(age_t AGE)
Returns the size of the container of the appropriate age class(es) for both sexes.
Definition: metapop.h:498
Individual * get(sex_t SEX, age_idx AGE, unsigned int at)
Returns a pointer to the individual sitting at the index passed.
Definition: metapop.h:534
unsigned int get_locus_num()
Definition: ttneutralgenes.h:200
Interface for all trait types, declares all basic trait operations.
Definition: ttrait.h:46
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(), TTrait::get_allele(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatchPtr(), Individual::getTrait(), MAL, and Patch::size().

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

◆ getFis()

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

References _fis.

Referenced by setFstatRecorders().

◆ getFisWC()

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

References _fis_WC.

Referenced by setFstatWCRecorders().

◆ getFit()

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

References _fit.

Referenced by setFstatRecorders().

◆ getFitWC()

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

References _fit_WC.

Referenced by setFstatWCRecorders().

◆ getFixLocGlobal()

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

References _fix_loc_global.

Referenced by setFstatRecorders().

◆ getFixLocLocal()

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

References _fix_loc_local.

Referenced by setFstatRecorders().

◆ getFst()

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

References _fst.

Referenced by setFstat2Recorders(), and setFstatRecorders().

◆ getFst_ij()

double TTNeutralGenesSH::getFst_ij ( unsigned int  i)
inline

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

445  {
446  unsigned int scale = (unsigned int)pow( 10.0, (int)log10((float)_fst_matrix->getNbCols()) + 1 );
447  return _fst_matrix->get(i/scale, i%scale);
448  }

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

Referenced by setFstMatrixRecorders().

◆ getFstWC()

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

References _fst_WC.

Referenced by setFstatWCRecorders().

◆ getGlobalAlleleFreq()

double TTNeutralGenesSH::getGlobalAlleleFreq ( unsigned int  loc,
unsigned int  all 
)
inline
396  {
397  return _globalAlleleFreq.get(loc, all);
398  }

References _globalAlleleFreq, and TMatrix::get().

Referenced by setFreqRecorders().

◆ getGlobalFreqs()

TMatrix* TTNeutralGenesSH::getGlobalFreqs ( )
inline

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

422 {return &_globalAlleleFreq;}

References _globalAlleleFreq.

Referenced by TTNeutralGenesFH::write_varcompWC().

◆ getHeteroTable()

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

References _heteroTable.

Referenced by TTNeutralGenesFH::write_varcompWC().

◆ getHeterozygosity()

double TTNeutralGenesSH::getHeterozygosity ( unsigned int  loc)
inline
400  {
401  double het = 0;
402  for(unsigned int i = 0; i < _heteroTable.getNumGroups(); ++i )
403  het += _heteroTable.get(i, loc, 0);
404  return het/_heteroTable.getNumGroups(); //mean per patch heterozygosity
405  }
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().

◆ getHo()

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

References _ho.

Referenced by setFstat2Recorders(), and setFstatRecorders().

◆ getHs()

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

References _hs.

◆ getHsnei()

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

References _hsnei.

Referenced by setFstat2Recorders(), and setFstatRecorders().

◆ getHt()

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

References _ht.

◆ getHtnei()

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

References _htnei.

Referenced by setFstat2Recorders(), and setFstatRecorders().

◆ getMeanAlpha()

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

References _mean_alpha.

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ getMeanNeiGeneticDistance()

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

References _meanD.

Referenced by setNeiGeneticDistanceRecorders().

◆ getMeanTheta()

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

References _mean_theta.

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ getNbAllGlobal()

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

References _nb_all_global.

Referenced by setFstatRecorders().

◆ getNbAllLocal()

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

References _nb_all_local.

Referenced by setFstatRecorders().

◆ getNeiGeneticDistance()

double TTNeutralGenesSH::getNeiGeneticDistance ( unsigned int  i)
inline
547  {
548  unsigned int scale = (unsigned int)pow( 10.0, (int)log10((float)_D->getNbCols()) + 1 );
549  return _D->get(i/scale,i%scale);
550  }

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

Referenced by setNeiGeneticDistanceRecorders().

◆ getSibCoaMeans()

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

References _sib_coa.

Referenced by setStatRecorders().

◆ getSibProportions()

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

References _sib_prop.

Referenced by setStatRecorders().

◆ getTheta_FF()

double TTNeutralGenesSH::getTheta_FF ( )
inline

Gives the mean within females coancestry coefficient.

530 {return Theta_FF;}
double Theta_FF
Definition: ttneutralgenes.h:330

References Theta_FF.

Referenced by setStatRecorders().

◆ getTheta_FM()

double TTNeutralGenesSH::getTheta_FM ( )
inline

Gives the mean between males and females coancestry coefficient.

534 {return Theta_FM;}
double Theta_FM
Definition: ttneutralgenes.h:330

References Theta_FM.

Referenced by setStatRecorders().

◆ getTheta_MM()

double TTNeutralGenesSH::getTheta_MM ( )
inline

Gives the mean within males coancestry coefficient.

532 {return Theta_MM;}
double Theta_MM
Definition: ttneutralgenes.h:330

References Theta_MM.

Referenced by setStatRecorders().

◆ getWeightedFst()

double TTNeutralGenesSH::getWeightedFst ( )
inline

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

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

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

References _fst_WH.

Referenced by setFstMatrixRecorders().

◆ init()

void TTNeutralGenesSH::init ( )
virtual

Reimplemented from StatHandlerBase.

2068 {
2070 
2072 
2074 }
virtual void init()
Definition: stathandler.cc:39
unsigned int get_allele_num()
Definition: ttneutralgenes.h:201
void allocateTables(unsigned int loci, unsigned int all)
Definition: stats_fstat.cc:45

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

◆ setAdltNeiGeneticDistance()

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

References ADULTS, and setNeiGeneticDistance().

Referenced by setNeiGeneticDistanceRecorders().

◆ setAdultAlleleFreq()

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

References ADULTS, and setAlleleTables().

Referenced by setFreqRecorders().

◆ setAdultHeterozygosity()

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

References ADULTS, and setHeterozygosity().

Referenced by setFreqRecorders().

◆ setAdults_Theta()

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

References ADLTx, FEM, Patch::get(), Individual::getTrait(), MAL, and Patch::size().

Referenced by setStatRecorders().

◆ setAdultsCoaBetween()

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

References ADLTx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ setAdultsCoaMatrix()

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

References ADLTx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ setAdultsCoaWithin()

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

References ADLTx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ setAdultsFstat()

void TTNeutralGenesSH::setAdultsFstat ( )
inline
457 {setFstat(ADULTS);}
void setFstat(age_t AGE)
Computes the F-statistics following Nei & Chesser (1983).
Definition: stats_fstat.cc:282

References ADULTS, and setFstat().

Referenced by setFstatRecorders().

◆ setAdultsFstat2()

void TTNeutralGenesSH::setAdultsFstat2 ( )
inline
472 {setFstat2(ADULTS);}
void setFstat2(age_t AGE)
New version of Nei & Chesser.
Definition: stats_fstat.cc:523

References ADULTS, and setFstat2().

Referenced by setFstat2Recorders().

◆ 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:830

References ADULTS, and setFstatWeirCockerham().

Referenced by setFstatWCRecorders().

◆ setAdultsFstBetween()

void TTNeutralGenesSH::setAdultsFstBetween ( )
inline
436 {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:724

References ADULTS, and setFstMatrix().

Referenced by setFstMatrixRecorders().

◆ setAdultsFstMatrix()

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

References ADULTS, and setFstMatrix().

Referenced by setFstMatrixRecorders().

◆ setAdultsFstWithin()

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

References ADULTS, and setFstMatrix().

Referenced by setFstMatrixRecorders().

◆ setAlleleTables()

void TTNeutralGenesSH::setAlleleTables ( age_t  AGE)
77 {
78 
79  if(_table_set_age == AGE
82  return;
83 
84  unsigned int nb_locus = _SHLinkedTrait->get_locus_num(), nb_allele = _SHLinkedTrait->get_allele_num();
85  unsigned int patch_size, patchNbr = _pop->getPatchNbr();
86  age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
87  Patch* crnt_patch;
88  TTrait* trait;
89 
90  //check if the population size has changed:
91  if(_alleleCountTable.getNumGroups() != patchNbr)
93 
95 
97 
98  unsigned int* count1 = new unsigned int[nb_locus]();
99 
100  for(unsigned int k = 0; k < patchNbr; ++k) {
101 
102  crnt_patch = _pop->getPatch(k);
103  memset(count1, 0, nb_locus * sizeof(unsigned int));
104 
105  for (int sx = 0; sx < 2; ++sx) {
106 
107  for(unsigned int i = 0, size = crnt_patch->size(sex_t(sx), age_pos); i < size; ++i) {
108 
109  const TTNeutralGenes_bitstring* bs = static_cast<const TTNeutralGenes_bitstring*>(
110  crnt_patch->get(sex_t(sx), age_pos, i)->getTrait(_SHLinkedTraitIndex));
111 
112  for (int chr = 0; chr < 2; ++chr) {
113  const bitstring& seq = bs->get_bit_sequence(chr);
114  size_t nwords = seq.nb_words();
115  for (size_t w = 0; w < nwords; ++w) {
116  bitstring::_ul word = *seq.getword_atIdx(w);
117  while (word) {
118  count1[w * BITS_PER_WORD + __builtin_ctzl(word)]++;
119  word &= word - 1;
120  }
121  }
122  }
123  }
124  }
125 
126  for (unsigned int j = 0; j < nb_locus; ++j) {
127  _alleleCountTable.set(k, j, 1, count1[j]);
128  _alleleCountTable.set(k, j, 0,
129  crnt_patch->size(age_pos) * 2 - count1[j]);
130  }
131  }
132 
133  delete[] count1;
134 
135  } else {
136 
137  //counting the copies of each allele present in each patch:
138  for(unsigned int k = 0; k < patchNbr; ++k) {
139 
140  crnt_patch = _pop->getPatch(k);
141 
142  for (int sx = 0; sx < 2; ++sx) {
143 
144  for(unsigned int i = 0, size = crnt_patch->size(sex_t(sx), age_pos); i < size; ++i) {
145 
146  trait = crnt_patch->get(sex_t(sx), age_pos, i)->getTrait(_SHLinkedTraitIndex);
147 
148  for (unsigned int j = 0; j < nb_locus; ++j)
149  _alleleCountTable.increment(k, j, trait->get_allele(j, 0));
150 
151  for (unsigned int j = 0; j < nb_locus; ++j)
152  _alleleCountTable.increment(k, j, trait->get_allele(j, 1));
153 
154  }
155  }
156  }
157  }
158 
159  //allelic frequencies:
160  for (unsigned int i = 0; i < patchNbr; ++i) {
161 
162  patch_size = _pop->getPatch(i)->size(age_pos) * 2;
163 
164  if (patch_size) {
165 
166  for (unsigned int l = 0; l < nb_locus; ++l)
167  for (unsigned int u = 0; u < nb_allele; ++u)
168  _alleleFreqTable.set(i, l, u, (double)_alleleCountTable.get(i, l, u) / (double)patch_size);
169 
170  } else {
171  for (unsigned int l = 0; l < nb_locus; ++l)
172  for (unsigned int u = 0; u < nb_allele; ++u)
173  _alleleFreqTable.set(i, l, u, 0 );
174  }
175  }
176 
177  //population-wide allelic frequencies:
178  unsigned int tot_size = _pop->size(AGE) * 2;
179 
181 
182  for(unsigned int i = 0; i < patchNbr; i++)
183  for (unsigned int l = 0; l < nb_locus; ++l)
184  for (unsigned int u = 0; u < nb_allele; ++u)
186 
187  for (unsigned int l = 0; l < nb_locus; ++l)
188  for (unsigned int u = 0; u < nb_allele; ++u)
189  _globalAlleleFreq.divide(l, u, tot_size);
190 
191  _table_set_age = AGE;
192 
194 
196 }
#define BITS_PER_WORD
Definition: bitstring.h:41
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:295
unsigned int getCurrentGeneration()
Definition: metapop.h:296
void assign(double val)
Assigns a value to all element of the matrix.
Definition: tmatrix.h:155
void divide(unsigned int i, unsigned int j, double value)
Divide an element of the matrix by a value.
Definition: tmatrix.h:289
void plus(unsigned int i, unsigned int j, double value)
Adds a value to an element of the matrix.
Definition: tmatrix.h:256
TTNeutralGenes_bitstring : diallelic neutral loci encoded as bitstrings.
Definition: ttneutralgenes_bitstring.h:40
const bitstring & get_bit_sequence(bool chromosome) const
Definition: ttneutralgenes_bitstring.h:52
Non-template and faster implementation of std::bitset.
Definition: bitstring.h:57
unsigned long _ul
Definition: bitstring.h:61
_ul * getword_atIdx(size_t index) const
Definition: bitstring.h:159
size_t nb_words() const
Definition: bitstring.h:164
sex_t
Sex types, males are always 0 and females 1!!
Definition: types.h:36

References _alleleCountTable, _alleleFreqTable, _globalAlleleFreq, _is_diallelic_bitstring, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, _table_set_age, _table_set_gen, _table_set_repl, ADLTx, ADULTS, allocateTables(), TMatrix::assign(), BITS_PER_WORD, TMatrix::divide(), Patch::get(), DataTable< T >::get(), TTrait::get_allele(), TProtoNeutralGenes::get_allele_num(), TTNeutralGenes_bitstring::get_bit_sequence(), TProtoNeutralGenes::get_locus_num(), Metapop::getCurrentGeneration(), Metapop::getCurrentReplicate(), DataTable< T >::getNumGroups(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), bitstring::getword_atIdx(), DataTable< T >::increment(), DataTable< T >::init(), bitstring::nb_words(), 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().

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

References FEM, Patch::get(), Individual::getTrait(), MAL, and Patch::size().

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

◆ setCoaMatrixRecorders()

void TTNeutralGenesSH::setCoaMatrixRecorders ( age_t  AGE,
unsigned char  dim 
)
2332 {
2333  std::ostringstream name, sub_name;
2334 
2335  void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ?
2342 
2343  const char *prefix = (AGE == ADULTS ? "adlt." : "off.");
2344 
2345  unsigned int nbpatch = _pop->getPatchNbr();
2346  unsigned int scale = (unsigned int)pow(10.0, (int)log10((float)nbpatch) + 1);
2347 
2348  if(dim & 1) {
2349  name<<"Wtn Patch Coancestry ("<<prefix<<")";
2350  sub_name<< prefix << "theta";
2351  add(name.str(), sub_name.str(), AGE,0, 0, &TTNeutralGenesSH::getMeanTheta, 0, 0, setter);
2352  name.str("");
2353  sub_name.str("");
2354  }
2355  if(dim & 2){
2356  name<<"Btn Patch Coancestry ("<<prefix<<")";
2357  sub_name<< prefix << "alpha";
2358  add(name.str(), sub_name.str(), AGE, 0, 0, &TTNeutralGenesSH::getMeanAlpha, 0, 0, (!(dim&1) ? setter : 0) );
2359  name.str("");
2360  sub_name.str("");
2361  }
2362  if(dim & 1) {
2363  for(unsigned int i = 0; i < nbpatch; ++i) {
2364  name<<"coancestry "<<i+1<<"."<<i+1;
2365  sub_name<< prefix << "coa" << i+1 << "." << i+1;
2366  add(name.str(), sub_name.str(), AGE, i*scale + i, 0, 0, &TTNeutralGenesSH::getCoa, 0, 0);
2367  name.str("");
2368  sub_name.str("");
2369  }
2370  }
2371  if(dim & 2){
2372  for(unsigned int i = 0; i < nbpatch; ++i) {
2373  for(unsigned int j = i+1; j < nbpatch; ++j) {
2374  name<<"coancestry "<<i+1<< "." <<j+1;
2375  sub_name<< prefix << "coa" << i+1 << "." << j+1;
2376  add(name.str(), sub_name.str(), AGE, i*scale + j, 0, 0, &TTNeutralGenesSH::getCoa, 0, 0);
2377  name.str("");
2378  sub_name.str("");
2379  }
2380  }
2381  }
2382 
2383 }
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:321
void setAdultsCoaWithin()
Definition: ttneutralgenes.h:512
void setOffsprgCoaWithin()
Definition: ttneutralgenes.h:513
double getCoa(unsigned int i)
Gets the given coancestry coefficient from the coancestry matrix.
Definition: ttneutralgenes.h:522
void setOffsprgCoaBetween()
Definition: ttneutralgenes.h:515
double getMeanTheta()
Definition: ttneutralgenes.h:527
double getMeanAlpha()
Definition: ttneutralgenes.h:528
void setOffsprgCoaMatrix()
Definition: ttneutralgenes.h:511
void setAdultsCoaBetween()
Definition: ttneutralgenes.h:514
void setAdultsCoaMatrix()
Definition: ttneutralgenes.h:510

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

Referenced by setStatRecorders().

◆ setDxyRecorders()

void TTNeutralGenesSH::setDxyRecorders ( age_t  AGE,
bool  patchwise 
)
2574 {
2575  string prefix = (AGE == OFFSPRG ? "off." : "adlt.");
2576  age_idx age = (AGE == OFFSPRG ? OFFSx : ADLTx );
2577  string name = "Sequence divergence - Dxy";
2578  string sub_name = "Dxy";
2579 
2580 
2581  if (!patchwise) {
2582 
2583  add(name, prefix + sub_name, AGE, (unsigned int)age, 0, 0, &TTNeutralGenesSH::getDxy, 0, 0);
2584 
2585  } else {
2586 
2587  for (unsigned int p1 = 0; p1 < _pop->getPatchNbr(); ++p1) {
2588  for (unsigned int p2 = p1 + 1; p2 < _pop->getPatchNbr(); ++p2) {
2589  add(name, prefix + sub_name + ".p" + tstring::int2str(p1+1) + "p" + tstring::int2str(p2+1),
2590  AGE, p1, p2,
2593  0);
2594  }
2595  }
2596 
2597  }
2598 }
double getDxyAdultPerPatch(unsigned int patch1, unsigned patch2)
Definition: ttneutralgenes.h:557
double getDxy(unsigned int age_class)
Definition: stats_fstat.cc:1269
double getDxyOffspringPerPatch(unsigned int patch1, unsigned patch2)
Definition: ttneutralgenes.h:556
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().

◆ setFreqRecorders()

void TTNeutralGenesSH::setFreqRecorders ( age_t  AGE)
2388 {
2389  string prefix = (AGE == OFFSPRG ? "off." : "adlt.");
2390 
2391  void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ? &TTNeutralGenesSH::setAdultAlleleFreq :
2393 
2394  unsigned int nb_allele = _SHLinkedTrait->get_allele_num();
2395  unsigned int nb_locus = _SHLinkedTrait->get_locus_num();
2396 
2397  for (unsigned int l = 0; l < nb_locus; ++l) {
2398  for (unsigned int u = 0; u < nb_allele-1; ++u) {
2399  add("", prefix + "ntrl.l" + tstring::int2str(l+1) + ".a" + tstring::int2str(u+1), AGE, l, u,
2400  0, 0, &TTNeutralGenesSH::getGlobalAlleleFreq, setter);
2401  }
2402  }
2403 
2404  setter = (AGE == ADULTS ? &TTNeutralGenesSH::setAdultHeterozygosity :
2406 
2407  for (unsigned int l = 0; l < nb_locus; ++l) {
2408  add("", prefix + "ntrl.l" + tstring::int2str(l+1) + ".Het", AGE, l, 0, 0,
2410  // add("", prefix + "ntrl.l" + tstring::int2str(l) + ".Hom", AGE, l, 0, 0,
2411  // &TTNeutralGenesSH::getHomozygosity, 0, 0);
2412  }
2413 
2414 }
void setAdultHeterozygosity()
Definition: ttneutralgenes.h:393
double getGlobalAlleleFreq(unsigned int loc, unsigned int all)
Definition: ttneutralgenes.h:396
double getHeterozygosity(unsigned int loc)
Definition: ttneutralgenes.h:400
void setAdultAlleleFreq()
Definition: ttneutralgenes.h:390
void setOffspringHeterozygosity()
Definition: ttneutralgenes.h:394
void setOffspringAlleleFreq()
Definition: ttneutralgenes.h:391

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

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

1117 {
1118  unsigned int patchNbr = this->_pop->getPatchNbr();
1119  unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
1120  unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
1121 
1122  assert(N == patchNbr);
1123  assert(L == nb_locus);
1124 
1126 
1127  double *pop_weights = new double[patchNbr];
1128  double *pop_sizes = new double[patchNbr];
1129  double *denominator = new double[nb_locus];
1130  double sum_weights = 0;
1131 
1132  double tot_size = _pop->size(ADULTS) * 2; //diploids!
1133 
1134  for(unsigned int i = 0; i < patchNbr; ++i) {
1135  pop_sizes[i] = _pop->size(ADULTS, i) * 2;
1136  pop_weights[i] = pop_sizes[i] - (pop_sizes[i] * pop_sizes[i] / tot_size); //n_ic in Weir & Hill 2002
1137  sum_weights += pop_weights[i];
1138  for(unsigned int j = 0; j < nb_locus; ++j)
1139  array[i][j] = nanf("NULL");
1140  }
1141 
1142  for(unsigned int l = 0; l < nb_locus; ++l)
1143  denominator[l] = 0;
1144 
1145  double p, pq, var;
1146 
1147  for (unsigned int i = 0; i < patchNbr; ++i) {
1148 
1149  if( !pop_sizes[i] ) continue;
1150 
1151  for (unsigned int l = 0; l < nb_locus; ++l) {
1152 
1153  array[i][l] = 0;
1154 
1155  for (unsigned int u = 0; u < nb_allele; ++u) {
1156 
1157  p = _alleleFreqTable.get(i, l, u); //p_liu
1158 
1159  pq = p * (1 - p);
1160 
1161  var = p - _globalAlleleFreq.get(l, u); //(p_liu - pbar_u)
1162 
1163  var *= var;
1164 
1165  array[i][l] += pq * pop_sizes[i] / (pop_sizes[i] -1);
1166 
1167  denominator[l] += pop_sizes[i] * var + pop_weights[i] * pq;
1168 
1169  } // end for allele
1170  }// end for locus
1171  }//end for pop
1172 
1173 
1174  for (unsigned int i = 0; i < patchNbr; ++i) {
1175  if( !pop_sizes[i] ) continue;
1176  for (unsigned int l = 0; l < nb_locus; ++l)
1177  array[i][l] = 1 - ( array[i][l] * sum_weights / denominator[l]);
1178  }
1179 
1180  delete [] pop_weights;
1181  delete [] pop_sizes;
1182  delete [] denominator;
1183 }

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

◆ setFstat()

void TTNeutralGenesSH::setFstat ( age_t  AGE)

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

283 {
284  double harmonic = 0, nbpatch = 0, nbind;
285  unsigned int patchNbr = _pop->getPatchNbr();
286  age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
287 
288  //check if the population size has changed:
289  if(_alleleCountTable.getNumGroups() != patchNbr)
291 
292  setAlleleTables(AGE);
293 
294  setLociDivCounter(AGE);
295 
296  //harmonic mean of patch sizes:
297  for (unsigned int i = 0; i < patchNbr; ++i){
298 
299  nbind = _pop->size(AGE, i);
300  if(nbind != 0){
301  nbpatch++;
302  harmonic += 1.0/nbind;
303  }
304  }
305 
306  harmonic = nbpatch / harmonic;
307 
308  _ho = setHo(age_pos);
309  _hs = setHs(age_pos);
310  _ht = setHt(age_pos);
311  //Nei's corrections:
312  _hsnei = (nbpatch != 0 ? harmonic/(harmonic-1.0)*(_hs-(_ho/(2.0*harmonic))) : nanf("NULL") );
313  _htnei = (nbpatch != 0 ? _ht + (_hsnei/(harmonic*nbpatch))
314  -(_ho/(2.0*harmonic*nbpatch)) : nanf("NULL") );
315 
316  _fis = ( _hsnei ? 1.0-(_ho/_hsnei) : nanf("NULL") );
317  _fit = ( _htnei ? 1.0-(_ho/_htnei) : nanf("NULL") );
318  _fst = ( _htnei ? 1.0-(_hsnei/_htnei) : nanf("NULL") );
319 }
double setHo(age_idx age_pos)
Definition: stats_fstat.cc:391
double setHt(age_idx age_pos)
Definition: stats_fstat.cc:496
void setLociDivCounter(age_t AGE)
Sets the allelic diversity counters.
Definition: stats_fstat.cc:323
double setHs(age_idx age_pos)
Definition: stats_fstat.cc:451

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

◆ setFstat2()

void TTNeutralGenesSH::setFstat2 ( age_t  AGE)

New version of Nei & Chesser.

524 {
525  double harmonic = 0, nbpatch = 0, nbind;
526  unsigned int patchNbr = _pop->getPatchNbr();
527  age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
528  unsigned int nloc = _SHLinkedTrait->get_locus_num();
529 
530  //check if the population size has changed:
531  if(_alleleCountTable.getNumGroups() != patchNbr)
533 
534  setAlleleTables(AGE);
535 
536  setLociDivCounter(AGE);
537 
538  //harmonic mean of patch sizes:
539  for (unsigned int i = 0; i < patchNbr; ++i){
540 
541  nbind = _pop->size(AGE, i);
542  if(nbind != 0){
543  nbpatch++;
544  harmonic += 1.0/nbind;
545  }
546  }
547 
548  harmonic = nbpatch / harmonic;
549 
550  deque<double> ho = setHo2(age_pos);
551  deque<double> hs = setHs2(age_pos);
552  deque<double> ht = setHt2(age_pos);
553 
554  _fis = _fit = _fst = 0;
555 
556  double hsnei, htnei;
557 
558  for (unsigned int l = 0; l < nloc; ++l) {
559 
560  //Nei's corrections:
561  hsnei = (nbpatch != 0 ? harmonic/(harmonic-1.0)*(hs[l]-(ho[l]/(2.0*harmonic))) : nanf("NULL") );
562  htnei = (nbpatch != 0 ? ht[l] + (hsnei/(harmonic*nbpatch))
563  -(ho[l]/(2.0*harmonic*nbpatch)) : 0 );
564  _hsnei += hsnei;
565  _htnei += htnei;
566 // _hsnei += hs[l];
567 // _htnei += ht[l];
568  _ho += ho[l];
569 // _fis = ( _hsnei ? 1.0-(_ho/_hsnei) : nanf("NULL") );
570 // _fit = ( _htnei ? 1.0-(_ho/_htnei) : nanf("NULL") );
571  if(ht[l] == 0) nloc--; //would mean locus is fixed
572  _fst += ( ht[l] != 0 ? 1.0-(hsnei/htnei) : 0);
573 
574 // if(ht[l] == 0) nloc--;
575 //
576 // _fst += (ht[l] != 0 ? 1 - (hs[l]/ht[l]) : 0);
577  }
578  _hsnei /= nloc;
579  _htnei /= nloc;
580  _ho /= nloc;
581  _fst /= nloc;
582 }
deque< double > setHo2(age_idx age_pos)
Definition: stats_fstat.cc:586
deque< double > setHt2(age_idx age_pos)
Definition: stats_fstat.cc:698
deque< double > setHs2(age_idx age_pos)
Definition: stats_fstat.cc:648

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

◆ setFstat2Recorders()

void TTNeutralGenesSH::setFstat2Recorders ( age_t  AGE)
2454 {
2455  if(AGE & OFFSPRG) {
2456  add("Ho (offsprg)","off.ho2",OFFSPRG,0,0,&TTNeutralGenesSH::getHo,0,0,
2458  add("Hs-Nei (offsprg)","off.hsnei2",OFFSPRG,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2459  add("Ht-Nei (offsprg)","off.htnei2",OFFSPRG,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2460  add("Fst (offsprg)","off.fst2",OFFSPRG,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2461  }
2462 
2463  if(AGE & ADULTS) {
2464  add("Ho (adult)","adlt.ho2",ADULTS,0,0,&TTNeutralGenesSH::getHo,0,0,
2466  add("Hs-Nei (adult)","adlt.hsnei2",ADULTS,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2467  add("Ht-Nei (adult)","adlt.htnei2",ADULTS,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2468  add("Fst (adult)","adlt.fst2",ADULTS,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2469  }
2470 }
void setAdultsFstat2()
Definition: ttneutralgenes.h:472
double getHsnei()
Definition: ttneutralgenes.h:461
double getHo()
Definition: ttneutralgenes.h:463
double getHtnei()
Definition: ttneutralgenes.h:462
void setOffsprgFstat2()
Definition: ttneutralgenes.h:471
double getFst()
Definition: ttneutralgenes.h:466

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

Referenced by setStatRecorders().

◆ setFstatRecorders()

void TTNeutralGenesSH::setFstatRecorders ( age_t  AGE)
2420 {
2421  if(AGE & OFFSPRG) {
2422  add("Nbr of Alleles per Locus - local (offsprg)","off.allnbp",OFFSPRG,0,0,
2424  add("Nbr of Alleles per Locus - global (offsprg)","off.allnb",OFFSPRG,0,0,&TTNeutralGenesSH::getNbAllGlobal,0,0,0);
2425  add("Nbr of Fixed Loci per Patch (offsprg)","off.fixlocp",OFFSPRG,0,0,&TTNeutralGenesSH::getFixLocLocal,0,0,0);
2426  add("Nbr of Fixed Loci in the Pop (offsprg)","off.fixloc",OFFSPRG,0,0,&TTNeutralGenesSH::getFixLocGlobal,0,0,0);
2427  add("Ho (offsprg)","off.ho",OFFSPRG,0,0,&TTNeutralGenesSH::getHo,0,0,0);
2428  add("Hs-Nei (offsprg)","off.hsnei",OFFSPRG,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2429  add("Ht-Nei (offsprg)","off.htnei",OFFSPRG,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2430  add("Fis (offsprg)","off.fis",OFFSPRG,0,0,&TTNeutralGenesSH::getFis,0,0,0);
2431  add("Fst (offsprg)","off.fst",OFFSPRG,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2432  add("Fit (offsprg)","off.fit",OFFSPRG,0,0,&TTNeutralGenesSH::getFit,0,0,0);
2433  }
2434 
2435  if(AGE & ADULTS) {
2436  add("Nbr of Alleles per Locus - local (adult)","adlt.allnbp",ADULTS,0,0,
2438  add("Nbr of Alleles per Locus - global (adult)","adlt.allnb",ADULTS,0,0,&TTNeutralGenesSH::getNbAllGlobal,0,0,0);
2439  add("Nbr of Fixed Loci per Patch (adult)","adlt.fixlocp",ADULTS,0,0,&TTNeutralGenesSH::getFixLocLocal,0,0,0);
2440  add("Nbr of Fixed Loci in the Pop (adult)","adlt.fixloc",ADULTS,0,0,&TTNeutralGenesSH::getFixLocGlobal,0,0,0);
2441  add("Ho (adult)","adlt.ho",ADULTS,0,0,&TTNeutralGenesSH::getHo,0,0,0);
2442  add("Hs-Nei (adult)","adlt.hsnei",ADULTS,0,0,&TTNeutralGenesSH::getHsnei,0,0,0);
2443  add("Ht-Nei (adult)","adlt.htnei",ADULTS,0,0,&TTNeutralGenesSH::getHtnei,0,0,0);
2444  add("Fis (adult)","adlt.fis",ADULTS,0,0,&TTNeutralGenesSH::getFis,0,0,0);
2445  add("Fst (adult)","adlt.fst",ADULTS,0,0,&TTNeutralGenesSH::getFst,0,0,0);
2446  add("Fit (adult)","adlt.fit",ADULTS,0,0,&TTNeutralGenesSH::getFit,0,0,0);
2447  }
2448 }
double getFixLocGlobal()
Definition: ttneutralgenes.h:490
double getFit()
Definition: ttneutralgenes.h:468
void setAdultsFstat()
Definition: ttneutralgenes.h:457
double getNbAllGlobal()
Definition: ttneutralgenes.h:488
double getFis()
Definition: ttneutralgenes.h:467
double getNbAllLocal()
Definition: ttneutralgenes.h:487
double getFixLocLocal()
Definition: ttneutralgenes.h:489
void setOffsprgFstat()
Definition: ttneutralgenes.h:456

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

Referenced by setStatRecorders().

◆ setFstatWCRecorders()

void TTNeutralGenesSH::setFstatWCRecorders ( age_t  AGE)
2476 {
2477  if(AGE & OFFSPRG) {
2478 
2479  add("Fis Weir & Cockerham","off.fis.WC",OFFSPRG,0,0,&TTNeutralGenesSH::getFisWC,0,0,&TTNeutralGenesSH::setOffspringFstatWeirCockerham);
2480  add("Fst Weir & Cockerham","off.fst.WC",OFFSPRG,0,0,&TTNeutralGenesSH::getFstWC,0,0,0);
2481  add("Fit Weir & Cockerham","off.fit.WC",OFFSPRG,0,0,&TTNeutralGenesSH::getFitWC,0,0,0);
2482  }
2483 
2484  if(AGE & ADULTS) {
2485 
2486  add("Fis Weir & Cockerham","adlt.fis.WC",ADULTS,0,0,&TTNeutralGenesSH::getFisWC,0,0,&TTNeutralGenesSH::setAdultsFstatWeirCockerham);
2487  add("Fst Weir & Cockerham","adlt.fst.WC",ADULTS,0,0,&TTNeutralGenesSH::getFstWC,0,0,0);
2488  add("Fit Weir & Cockerham","adlt.fit.WC",ADULTS,0,0,&TTNeutralGenesSH::getFitWC,0,0,0);
2489 
2490  }
2491 }
void setAdultsFstatWeirCockerham()
Definition: ttneutralgenes.h:481
void setOffspringFstatWeirCockerham()
Definition: ttneutralgenes.h:480
double getFitWC()
Definition: ttneutralgenes.h:484
double getFstWC()
Definition: ttneutralgenes.h:482
double getFisWC()
Definition: ttneutralgenes.h:483

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

Referenced by setStatRecorders().

◆ setFstatWeirCockerham()

void TTNeutralGenesSH::setFstatWeirCockerham ( age_t  AGE)

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

831 {
832  unsigned int patchNbr = this->_pop->getPatchNbr();
833  unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
834  unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
835  //age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
836 
837  setAlleleTables(AGE);
838  setHeteroTable(AGE);
839 
840  //init
841  double *pop_sizes = new double [patchNbr];
842  double tot_size, inv_ntot;
843  double sum_weights = 0;
844  double nbar, nc, inv_nbar;
845  unsigned int extantPs = 0;
846 
847  tot_size = _pop->size(AGE);
848 
849  for(unsigned int i = 0; i < patchNbr; i++) {
850  pop_sizes[i] = _pop->size(AGE, i);
851  if(pop_sizes[i]) {
852  extantPs++;
853  sum_weights += (pop_sizes[i] * pop_sizes[i] / tot_size);
854  }
855  }
856  nbar = tot_size/extantPs;
857  nc = (tot_size - sum_weights)/(extantPs-1);
858  inv_nbar = 1/(nbar - 1);
859  inv_ntot = 1/tot_size;
860 
861  double var;
862  double s2, pbar, hbar;
863  double s2_denom = 1.0/((extantPs-1)*nbar),
864  r = (double)(extantPs-1)/extantPs,
865  hbar_factor=(2*nbar-1)/(4*nbar);
866  double a = 0, b = 0, c = 0, x;
867 
868  for (unsigned int l = 0; l < nb_locus; ++l) {
869 
870  for (unsigned int u = 0; u < nb_allele; ++u) {
871 
872  s2 = pbar = hbar = 0;
873 
874  for (unsigned int i = 0; i < patchNbr; ++i) {
875 
876  var = _alleleFreqTable.get(i, l, u) - _globalAlleleFreq.get(l, u); //(p_liu - pbar_u)^2
877 
878  var *= var;
879 
880  s2 += var * pop_sizes[i];
881 
882  hbar += _heteroTable.get(i, l, u);
883 
884  }//end for pop
885 
886  s2 *= s2_denom;
887  pbar = _globalAlleleFreq.get(l, u);
888  hbar *= inv_ntot;
889 
890  x = pbar * (1 - pbar) - r * s2;
891  a += s2 - inv_nbar*( x - 0.25 * hbar);
892  b += x - hbar_factor * hbar;
893  c += hbar;
894 
895 
896  } // end for allele
897 
898  }// end for locus
899 
900  a *= nbar/nc;
901  b *= nbar/(nbar - 1);
902  c *= 0.5;
903 
904  _fst_WC = a / (a + b + c);
905  _fit_WC = (a + b) / (a + b + c);
906  _fis_WC = b / (b + c);
907 
908  delete [] pop_sizes;
909 }
void setHeteroTable(age_t AGE)
Definition: stats_fstat.cc:200

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

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

914 {
917  unsigned int patchNbr = this->_pop->getPatchNbr();
918  unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
919  unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
920  //age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
921 
922  setAlleleTables(AGE);
923  setHeteroTable(AGE);
924 
925  //init
926  double *pop_sizes = new double [patchNbr];
927  double tot_size;
928  double sum_weights = 0;
929  double nc;
930  unsigned int extantPs = 0;
931 
932  tot_size = _pop->size(AGE);
933 
934  for(unsigned int i = 0; i < patchNbr; i++) {
935  pop_sizes[i] = _pop->size(AGE, i);
936  if(pop_sizes[i]) {
937  extantPs++;
938  sum_weights += (pop_sizes[i] * pop_sizes[i] / tot_size);
939  }
940  }
941 
942  nc = (tot_size - sum_weights)/(extantPs-1);
943 
944 // unsigned int np = extantPs;
945  unsigned int npl = extantPs; //all loci typed in all patches
946 
947  //p = _alleleFreqTable
948  //pb = _globalAlleleFreq
949 
950  unsigned int *alploc = new unsigned int [nb_locus];
951 
952  unsigned int **alploc_table = new unsigned int* [nb_locus];
953 
954  for(unsigned int i = 0; i < nb_locus; ++i)
955  alploc_table[i] = new unsigned int[nb_allele];
956 
957  unsigned int tot_num_allele = 0;
958 
959  for(unsigned int l = 0; l < nb_locus; ++l){
960 
961  alploc[l] = 0;
962 
963  for(unsigned int cnt, a = 0; a < nb_allele; ++a) {
964 
965  cnt=0;
966 
967  for(unsigned int i = 0; i < patchNbr; i++) {
968 
969  cnt += _alleleCountTable.get(i,l,a);
970 
971  }
972  alploc_table[l][a] = (cnt != 0);
973  alploc[l] += (cnt != 0);
974  }
975 
976  tot_num_allele += alploc[l];
977  }
978 
979  //n, and nal are given by pop_sizes, same num ind typed at all loci in each patch
980  //nc is the same for each locus
981  //nt is given by tot_size, same tot num of ind typed for all loci
982 
983  //SSG: het/2 for each allele
984  double *SSG = new double[tot_num_allele];
985  double *SSP = new double[tot_num_allele];
986  double *SSi = new double[tot_num_allele];
987 
988  unsigned int all_cntr = 0;
989 
990  double het, freq, var;
991 
992  for(unsigned int l = 0; l < nb_locus; ++l) {
993 
994  for(unsigned int a = 0; a < nb_allele & all_cntr < tot_num_allele; ++a) {
995 
996  if(alploc_table[l][a] == 0) continue; //do not consider alleles not present in the pop
997 
998  SSG[all_cntr] = 0;
999  SSi[all_cntr] = 0;
1000  SSP[all_cntr] = 0;
1001 
1002  for(unsigned int p = 0; p < patchNbr; ++p){
1003 
1004  if(!_pop->size(AGE, p)) continue; //skip empty patches
1005 
1006  het = _heteroTable.get(p, l, a);
1007 
1008  freq = _alleleFreqTable.get(p, l, a);
1009 
1010  var = freq - _globalAlleleFreq.get(l, a); //(p_liu - pbar_u)^2
1011 
1012  var *= var;
1013 
1014  SSG[all_cntr] += het;
1015 
1016  SSi[all_cntr] += 2*pop_sizes[p]*freq*(1-freq) - het/2;
1017 
1018  SSP[all_cntr] += 2*pop_sizes[p]*var;
1019  }
1020 
1021  all_cntr++;
1022  }
1023 
1024  }
1025 
1026 
1027  assert(all_cntr == tot_num_allele);
1028 
1029  double *MSG = new double[tot_num_allele];
1030  double *MSP = new double[tot_num_allele];
1031  double *MSI = new double[tot_num_allele];
1032  double *sigw = new double[tot_num_allele];
1033  double *siga = new double[tot_num_allele];
1034  double *sigb = new double[tot_num_allele];
1035 
1036 // double *FST_pal = new double[tot_num_allele];
1037 // double *FIS_pal = new double[tot_num_allele];
1038 
1039  double SIGA = 0, SIGB = 0, SIGW = 0;
1040 
1041  for(unsigned int i = 0; i < tot_num_allele; ++i){
1042 
1043  MSG[i] = SSG[i] / (2 * tot_size);
1044  sigw[i] = MSG[i]; //wasted!
1045 
1046  MSP[i] = SSP[i] / (npl-1);
1047 
1048  MSI[i] = SSi[i]/ (tot_size - npl);
1049 
1050  sigb[i] = 0.5*(MSI[i] - MSG[i]);
1051 
1052  siga[i] = (MSP[i] - MSI[i])/(2*nc);
1053 
1054 // FST_pal[i] = siga[i]/(siga[i]+sigb[i]+sigw[i]);
1055 // FIS_pal[i] = sigb[i]/(sigb[i]+sigw[i]);
1056 
1057  SIGA += siga[i];
1058  SIGB += sigb[i];
1059  SIGW += sigw[i];
1060  }
1061 
1062  //per locus stats:
1063  if(_fst_WC_loc) delete [] _fst_WC_loc; _fst_WC_loc = new double [nb_locus];
1064  if(_fis_WC_loc) delete [] _fis_WC_loc; _fis_WC_loc = new double [nb_locus];
1065  if(_fit_WC_loc) delete [] _fit_WC_loc; _fit_WC_loc = new double [nb_locus];
1066 
1067  double lsiga, lsigb, lsigw;
1068 
1069 // cout<<" computing sigma per locus\n";
1070 
1071  for(unsigned int allcntr = 0, i = 0; i < nb_locus; ++i) {
1072 
1073  lsiga = lsigb = lsigw = 0;
1074 
1075  for(unsigned int l = 0; l < alploc[i]; ++l) {
1076 
1077  lsiga += siga[allcntr];
1078  lsigb += sigb[allcntr];
1079  lsigw += sigw[allcntr];
1080 
1081  allcntr++;
1082 
1083  }
1084 
1085  _fst_WC_loc[i] = lsiga /(lsiga + lsigb + lsigw);
1086  _fis_WC_loc[i] = lsigb /(lsigb + lsigw);
1087  _fit_WC_loc[i] = (lsiga +lsigb) /(lsiga + lsigb + lsigw);
1088 
1089  }
1090 
1091 
1092  // Total F-stats
1093  _fst_WC = SIGA / (SIGA + SIGB + SIGW);
1094  _fit_WC = (SIGA + SIGB) / (SIGA + SIGB + SIGW);
1095  _fis_WC = SIGB / (SIGB + SIGW);
1096 
1097  delete[]pop_sizes;
1098  delete[]alploc;
1099  for(unsigned int i = 0; i < nb_locus; ++i)
1100  delete[]alploc_table[i];
1101  delete[]alploc_table;
1102  delete[]SSG;
1103  delete[]SSi;
1104  delete[]SSP;
1105  delete[]MSG;
1106  delete[]MSI;
1107  delete[]MSP;
1108  delete[]sigw;
1109  delete[]siga;
1110  delete[]sigb;
1111 
1112 }

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
725 {
726  unsigned int patchNbr = this->_pop->getPatchNbr();
727  unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
728  unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
729 
730  setAlleleTables(AGE);
731 
732  if(_fst_matrix == NULL)
733 
734  _fst_matrix = new TMatrix(patchNbr, patchNbr);
735 
736  else if( _fst_matrix->length() != patchNbr * patchNbr)
737 
738  _fst_matrix->reset(patchNbr, patchNbr);
739 
740  _fst_matrix->assign(nanf("NULL"));
741 
742  //init
743  double *pop_weights = new double[patchNbr];
744  double *pop_sizes = new double[patchNbr];
745  double **numerator = new double*[patchNbr];
746  for(unsigned int i = 0; i < patchNbr; i++) numerator[i] = new double [patchNbr];
747  double tot_size;
748  double numerator_W = 0;
749  double denominator = 0;
750  double sum_weights = 0;
751 
752  tot_size = _pop->size(AGE) * 2;
753 
754  for(unsigned int i = 0; i < patchNbr; ++i) {
755  pop_sizes[i] = _pop->size(AGE, i) * 2;
756  pop_weights[i] = pop_sizes[i] - (pop_sizes[i] * pop_sizes[i] / tot_size); //n_ic in Weir & Hill 2002
757  sum_weights += pop_weights[i];
758  for(unsigned int j = 0; j < patchNbr; j++)
759  numerator[i][j] = 0;
760  }
761 
762  double p, pq, var, num;
763 
764  for (unsigned int i = 0; i < patchNbr; ++i) {
765 
766  if( !pop_sizes[i] ) continue;
767 
768  for (unsigned int l = 0; l < nb_locus; ++l) {
769 
770  for (unsigned int u = 0; u < nb_allele; ++u) {
771 
772  p = _alleleFreqTable.get(i, l, u); //p_liu
773 
774  pq = p * (1 - p);
775 
776  var = p - _globalAlleleFreq.get(l, u); //(p_liu - pbar_u)^2
777 
778  var *= var;
779 
780  num = pq * pop_sizes[i] / (pop_sizes[i] -1);
781 
782  numerator[i][i] += num;
783 
784  numerator_W += num * pop_sizes[i]; //see equ. 9, Weir & Hill 2002
785 
786  denominator += pop_sizes[i] * var + pop_weights[i] * pq; //common denominator
787 
788  } // end for allele
789  }// end for locus
790  }//end for pop
791 
792  for (unsigned int i = 0; i < patchNbr; ++i) {
793  if( !pop_sizes[i] ) continue;
794  _fst_matrix->set(i, i, 1 - (numerator[i][i] * sum_weights / denominator) );
795  }
796  _fst_WH = 1 - ((numerator_W * sum_weights) / (denominator * tot_size)); //equ. 9 Weir & Hill 2002
797 
798  //pairwise Fst:
799  if(dim & 2) {
800  double pi, pj;
801  for (unsigned int l = 0; l < nb_locus; ++l)
802  for (unsigned int u = 0; u < nb_allele; ++u)
803  for (unsigned int i = 0; i < patchNbr - 1; ++i) {
804  if( !pop_sizes[i] ) continue;
805  for (unsigned int j = i + 1; j < patchNbr; ++j) {
806  if( !pop_sizes[j] ) continue;
807  pi = _alleleFreqTable.get(i, l, u);
808  pj = _alleleFreqTable.get(j, l, u);
809  numerator[i][j] += pi * (1 - pj) + pj * (1 - pi); //equ. 7 of Weir & Hill 2002
810  }
811  }
812  for (unsigned int i = 0; i < patchNbr - 1; ++i){
813  if( !pop_sizes[i] ) continue;
814  for (unsigned int j = i + 1; j < patchNbr; ++j){
815  if( !pop_sizes[j] ) continue;
816  _fst_matrix->set(i, j, 1 - ( (numerator[i][j] * sum_weights) / (2 * denominator)) );
817  }
818  }
819  }
820 
821 
822  delete [] pop_weights;
823  delete [] pop_sizes;
824  for(unsigned int i = 0; i < patchNbr; i++) delete [] numerator[i];
825  delete [] numerator;
826 }
void set(unsigned int i, unsigned int j, double val)
Sets element at row i and column j to value val.
Definition: tmatrix.h:103

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

◆ setFstMatrixRecorders()

void TTNeutralGenesSH::setFstMatrixRecorders ( age_t  AGE,
unsigned char  dim 
)
2496 {
2497  std::ostringstream name, sub_name;
2498 
2499  void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ?
2506 
2507  const char *prefix = (AGE == ADULTS ? "adlt." : "off.");
2508 
2509  unsigned int nbpatch = _pop->getPatchNbr();
2510  unsigned int scale = (unsigned int)pow(10.0, (int)log10((float)nbpatch) + 1);
2511 
2512  name<<"Weir&Hill weighted Fst ("<<prefix<<")";
2513  sub_name<< prefix << "fst.WH";
2514  add(name.str(), sub_name.str(), AGE, 0, 0, &TTNeutralGenesSH::getWeightedFst, 0, 0, setter);
2515  name.str("");
2516  sub_name.str("");
2517 
2518  if(dim & 1) {
2519  for(unsigned int i = 0; i < nbpatch; ++i) {
2520  name<<"Weighted Fst "<<i+1<<"."<<i+1;
2521  sub_name<< prefix << "fst" << i+1 << "." << i+1;
2522  add(name.str(), sub_name.str(), AGE, i*scale + i, 0, 0, &TTNeutralGenesSH::getFst_ij, 0, 0);
2523  name.str("");
2524  sub_name.str("");
2525  }
2526  }
2527  if(dim & 2){
2528  for(unsigned int i = 0; i < nbpatch; ++i) {
2529  for(unsigned int j = i+1; j < nbpatch; ++j) {
2530  name<<"Weighted Fst "<<i+1<< "." <<j+1;
2531  sub_name<< prefix << "fst" << i+1 << "." << j+1;
2532  add(name.str(), sub_name.str(), AGE, i*scale + j, 0, 0, &TTNeutralGenesSH::getFst_ij, 0, 0);
2533  name.str("");
2534  sub_name.str("");
2535  }
2536  }
2537  }
2538 
2539 }
double getWeightedFst()
Returns the weighted Fst using Weir & Hill (2002) method.
Definition: ttneutralgenes.h:442
void setAdultsFstBetween()
Definition: ttneutralgenes.h:436
void setOffsprgFstBetween()
Definition: ttneutralgenes.h:439
void setAdultsFstMatrix()
Definition: ttneutralgenes.h:434
void setOffsprgFstWithin()
Definition: ttneutralgenes.h:438
double getFst_ij(unsigned int i)
Accessor to the Fst matrix as set by setFstMatrix().
Definition: ttneutralgenes.h:444
void setAdultsFstWithin()
Definition: ttneutralgenes.h:435
void setOffsprgFstMatrix()
Definition: ttneutralgenes.h:437

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

Referenced by setStatRecorders().

◆ setHeteroTable()

void TTNeutralGenesSH::setHeteroTable ( age_t  AGE)
201 {
202  unsigned int patchNbr = this->_pop->getPatchNbr();
203  unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
204  unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
205  Patch* patch;
206  unsigned int isHetero, a0, a1;
207  TTrait* trait;
208  age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
209 
210  unsigned int **sizes;
211 
212  sizes = new unsigned int * [patchNbr];
213 
214  for(unsigned int i = 0; i < patchNbr; ++i) {
215  sizes[i] = new unsigned int [nb_locus];
216  for(unsigned int j = 0; j < nb_locus; ++j)
217  sizes[i][j] = nb_allele;
218  }
219 
220  _heteroTable.update(patchNbr, nb_locus, sizes);
221  _heteroTable.init(0);
222 
223 
224  for (unsigned int i = 0; i < patchNbr; ++i) {
225 
226  patch = _pop->getPatch(i);
227 
228  if(!patch->size(age_pos)) continue;
229 
230  for (int sx = 0; sx < 2; ++sx) {
231 
232  for(unsigned int j = 0, size = patch->size(sex_t(sx), age_pos); j < size; ++j) {
233 
234  trait = patch->get(sex_t(sx), age_pos, j)->getTrait(_SHLinkedTraitIndex);
235 
236  for (unsigned int l = 0; l < nb_locus; ++l) {
237  a0 = trait->get_allele(l, 0);
238  a1 = trait->get_allele(l, 1);
239  isHetero = (a0 != a1);
240  _heteroTable.plus(i, l, a0, isHetero);
241  _heteroTable.plus(i, l, a1, isHetero);
242  }
243  }
244  }
245  }
246  for(unsigned int i = 0; i < patchNbr; ++i) delete [] sizes[i];
247  delete [] sizes;
248 
249 }
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, Patch::get(), TTrait::get_allele(), TProtoNeutralGenes::get_allele_num(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatch(), Metapop::getPatchNbr(), DataTable< T >::init(), OFFSx, DataTable< T >::plus(), Patch::size(), and DataTable< T >::update().

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

◆ setHeterozygosity()

void TTNeutralGenesSH::setHeterozygosity ( age_t  AGE)
254 {
255  unsigned int patchNbr = this->_pop->getPatchNbr();
256  unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
257  unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
258  Patch* patch;
259  age_idx age_pos = (AGE == ADULTS ? ADLTx : OFFSx);
260 
261  setHeteroTable(AGE);
262 
263  for (unsigned int i = 0, size; i < patchNbr; ++i) {
264 
265  patch = _pop->getPatch(i);
266  size = patch->size(age_pos);
267 
268  if(!size) continue;
269 
270  for (unsigned int l = 0; l < nb_locus; ++l) {
271  for (unsigned int u = 0; u < nb_allele; ++u) {
272  _heteroTable.divide(i, l , u, size);
273  }
274  }
275  }
276 
277 
278 }
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().

◆ setHo()

double TTNeutralGenesSH::setHo ( age_idx  age_pos)
392 {
393  unsigned int nloc = _SHLinkedTrait->get_locus_num(), nbpatch = _pop->getPatchNbr(), psize;
394  double indloc = 0, hetero = 0;
395  Patch* current_patch;
396 
398 
399  for (unsigned int i = 0; i < nbpatch; ++i) {
400  current_patch = _pop->getPatch(i);
401 
402  psize = current_patch->size(FEM, age_pos);
403  indloc += psize;
404  for(unsigned int j = 0; j < psize; ++j) {
405  const TTNeutralGenes_bitstring* bs = static_cast<const TTNeutralGenes_bitstring*>(
406  current_patch->get(FEM, age_pos, j)->getTrait(_SHLinkedTraitIndex));
407  hetero += bs->get_bit_sequence(0).count_xor(bs->get_bit_sequence(1));
408  }
409 
410  psize = current_patch->size(MAL, age_pos);
411  indloc += psize;
412  for(unsigned int j = 0; j < psize; ++j) {
413  const TTNeutralGenes_bitstring* bs = static_cast<const TTNeutralGenes_bitstring*>(
414  current_patch->get(MAL, age_pos, j)->getTrait(_SHLinkedTraitIndex));
415  hetero += bs->get_bit_sequence(0).count_xor(bs->get_bit_sequence(1));
416  }
417  }
418 
419  } else {
420 
421  TTrait* trait;
422 
423  for (unsigned int i = 0; i < nbpatch; ++i) {
424  current_patch = _pop->getPatch(i);
425  psize = current_patch->size(FEM, age_pos);
426  indloc += psize;
427  for(unsigned int j = 0; j < psize; ++j) {
428  trait = current_patch->get(FEM, age_pos, j)->getTrait(_SHLinkedTraitIndex);
429  for (unsigned int k = 0; k < nloc; ++k) hetero += (trait->get_allele(k, 0) != trait->get_allele(k, 1));
430  }
431  }
432  for (unsigned int i = 0; i < nbpatch; ++i) {
433  current_patch = _pop->getPatch(i);
434  psize = current_patch->size(MAL, age_pos);
435  indloc += psize;
436  for(unsigned int j = 0; j < psize; ++j) {
437  trait = current_patch->get(MAL, age_pos, j)->getTrait(_SHLinkedTraitIndex);
438  for (unsigned int k = 0; k < nloc; ++k) hetero += (trait->get_allele(k, 0) != trait->get_allele(k, 1));
439  }
440  }
441  }
442 
443  indloc *= nloc;
444 
445  return (indloc != 0 ? hetero/indloc : 0.0);
446 }
unsigned int count_xor(const bitstring &other) const
Fused XOR popcount: count set bits in (this XOR other).
Definition: bitstring.h:309

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

Referenced by setFstat().

◆ setHo2()

deque< double > TTNeutralGenesSH::setHo2 ( age_idx  age_pos)
587 {
588  unsigned int nloc = _SHLinkedTrait->get_locus_num(), nbpatch = _pop->getPatchNbr();
589  unsigned int psize;
590  deque<double> hetero(nloc,0);
591  Patch* current_patch;
592 
594 
595  for (unsigned int i = 0; i < nbpatch; ++i) {
596  current_patch = _pop->getPatch(i);
597 
598  for (int sx = 0; sx < 2; ++sx) {
599  sex_t sex = (sx == 0) ? FEM : MAL;
600  for(unsigned int j = 0, size = current_patch->size(sex, age_pos); j < size; ++j) {
601  const TTNeutralGenes_bitstring* bs = static_cast<const TTNeutralGenes_bitstring*>(
602  current_patch->get(sex, age_pos, j)->getTrait(_SHLinkedTraitIndex));
603  const bitstring& s0 = bs->get_bit_sequence(0);
604  const bitstring& s1 = bs->get_bit_sequence(1);
605  size_t nwords = s0.nb_words();
606  for (size_t w = 0; w < nwords; ++w) {
607  bitstring::_ul xword = *s0.getword_atIdx(w) ^ *s1.getword_atIdx(w);
608  while (xword) {
609  hetero[w * BITS_PER_WORD + __builtin_ctzl(xword)]++;
610  xword &= xword - 1;
611  }
612  }
613  }
614  }
615  }
616 
617  } else {
618 
619  TTrait* trait;
620 
621  for (unsigned int i = 0; i < nbpatch; ++i) {
622  current_patch = _pop->getPatch(i);
623 
624  for(unsigned int j = 0, fsize = current_patch->size(FEM, age_pos); j < fsize; ++j) {
625  trait = current_patch->get(FEM, age_pos, j)->getTrait(_SHLinkedTraitIndex);
626  for (unsigned int k = 0; k < nloc; ++k) hetero[k] += (trait->get_allele(k, 0) != trait->get_allele(k, 1));
627  }
628 
629  for(unsigned int j = 0, msize = current_patch->size(MAL, age_pos); j < msize; ++j) {
630  trait = current_patch->get(MAL, age_pos, j)->getTrait(_SHLinkedTraitIndex);
631  for (unsigned int k = 0; k < nloc; ++k) hetero[k] += (trait->get_allele(k, 0) != trait->get_allele(k, 1));
632  }
633  }
634  }
635 
636  psize = _pop->size(age_pos);
637 
638  if(psize != 0)
639  for (unsigned int k = 0; k < nloc; ++k)
640  hetero[k] /= psize;
641 
642  return hetero;
643 }

References _is_diallelic_bitstring, StatHandlerBase::_pop, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTrait, TraitStatHandler< TProtoNeutralGenes, TTNeutralGenesSH >::_SHLinkedTraitIndex, BITS_PER_WORD, FEM, Patch::get(), TTrait::get_allele(), TTNeutralGenes_bitstring::get_bit_sequence(), TProtoNeutralGenes::get_locus_num(), Metapop::getPatch(), Metapop::getPatchNbr(), Individual::getTrait(), bitstring::getword_atIdx(), MAL, bitstring::nb_words(), Metapop::size(), and Patch::size().

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

◆ setHs()

double TTNeutralGenesSH::setHs ( age_idx  age_pos)
452 {
453  unsigned int i, k, x, nb_patch=0;
454  double hs = 0;
455  unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
456  unsigned int patchNbr = _pop->getPatchNbr();
457  deque<double>genehet(nbLoc,1);
458  deque<double>hspop(patchNbr,0);
459  double freq;
460  Patch* current_patch;
461 
462  //compute the expected heterozygosity for each patch and return average:
463  for (i = 0; i < patchNbr; ++i) {
464 
465  current_patch = _pop->getPatch(i);
466 
467  genehet.assign(nbLoc, 1);
468 
469  if(current_patch->size(age_pos) != 0) {
470 
471  nb_patch++;
472 
473  for (k = 0; k < nbLoc; ++k) {
474 
475  for (x = 0; x < nbAll; ++x) {
476 
477  freq = _alleleFreqTable.get(i, k, x);
478 
479  freq *= freq; //squared frequencies (expected _homozygosity)
480 
481  genehet[k] -= freq; //1 - sum of p2 = expected heterozygosity
482  }
483  //expected heterozygosity:
484  hspop[i] += genehet[k];
485  }
486  }//end_if size!=0
487 
488  hs += hspop[i];
489  }//end patch loop
490 
491  return (nb_patch !=0 ? hs/(nbLoc*nb_patch) : 0.0);
492 }
void assign(sex_t SEX, age_idx AGE, size_t 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().

◆ setHs2()

deque< double > TTNeutralGenesSH::setHs2 ( age_idx  age_pos)
649 {
650  unsigned int i, k, x, nb_patch = 0;
651  unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
652  unsigned int patchNbr = _pop->getPatchNbr();
653  deque<double>genehet(nbLoc,1);
654  deque<double>hsloc(nbLoc,0);
655  double freq;
656  Patch* current_patch;
657 
658  for (i = 0; i < patchNbr; ++i) {
659  if(_pop->size(age_pos, i) != 0)
660  nb_patch++;
661  }
662 
663  //compute the expected heterozygosity for each patch and return average:
664  for (k = 0; k < nbLoc; ++k) {
665 
666  for (i = 0; i < patchNbr; ++i) {
667 
668  current_patch = _pop->getPatch(i);
669 
670  genehet[k] = 1;
671 
672  if(current_patch->size(age_pos) != 0) {
673 
674  for (x = 0; x < nbAll; ++x) {
675 
676  freq = _alleleFreqTable.get(i, k, x);
677 
678  freq *= freq; //squared frequencies (expected _homozygosity)
679 
680  genehet[k] -= freq; //1 - sum of p2 = expected heterozygosity
681  }
682 
683  }//end_if size!=0
684 
685  hsloc[k] += genehet[k];
686 
687  }//end patch loop
688 
689  if(nb_patch !=0) hsloc[k] /= nb_patch;
690 
691  }//end loci loop
692 
693  return hsloc;
694 }

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

◆ setHt()

double TTNeutralGenesSH::setHt ( age_idx  age_pos)
497 {
498  double ht = 0;
499  unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
500  deque<double> genehet(nbLoc, 1);
501  double freq;
502 
503  //get global allele frequencies per locus
504  for (unsigned int l = 0; l < nbLoc; ++l) {
505  for (unsigned int u = 0; u < nbAll; ++u) {
506 
507  freq = _globalAlleleFreq.get(l, u);
508 
509  freq *= freq; //squared frequencies
510 
511  genehet[l] -= freq; //1 - sum of p2 = expected heterozygosity
512  }
513 
514  ht += genehet[l];
515  }
516 
517  return (ht/nbLoc);
518 }

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

Referenced by setFstat().

◆ setHt2()

deque< double > TTNeutralGenesSH::setHt2 ( age_idx  age_pos)
699 {
700  unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
701  deque<double> genehet(nbLoc, 1);
702  double freq;
703 
704  //get global allele frequencies per locus
705  for (unsigned int l = 0; l < nbLoc; ++l) {
706 
707  for (unsigned int u = 0; u < nbAll; ++u) {
708 
709  freq = _globalAlleleFreq.get(l, u);
710 
711  freq *= freq; //squared frequencies
712 
713  genehet[l] -= freq; //1 - sum of p2 = expected heterozygosity
714  }
715  }
716 
717  return genehet;
718 }

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

Referenced by setFstat2().

◆ setLociDivCounter()

void TTNeutralGenesSH::setLociDivCounter ( age_t  AGE)

Sets the allelic diversity counters.

324 {
325  unsigned int i, j, k;
326  unsigned int nbLoc = _SHLinkedTrait->get_locus_num(), nbAll = _SHLinkedTrait->get_allele_num();
327  unsigned int nbpatch = 0, patchNbr = _pop->getPatchNbr();
328  double patch_mean, pop_mean = 0;
329  bool **pop_div;
330 
331  //number of alleles per locus, Patch and pop counters:
332  pop_div = new bool * [nbLoc];
333 
334  for(i = 0; i < nbLoc; ++i) {
335  pop_div[i] = new bool [nbAll];
336  for(j = 0; j < nbAll;++j)
337  pop_div[i][j] = 0;
338  }
339 
340  for(k = 0; k < patchNbr; ++k)
341  nbpatch += (_pop->size(AGE, k) != 0);
342 
343  for(k = 0; k < patchNbr; ++k) {
344 
345  patch_mean = 0;
346 
347  for(i = 0; i < nbLoc; ++i)
348  for(j = 0; j < nbAll;++j) {
349  patch_mean += (_alleleCountTable.get(k,i,j) != 0);
350  pop_div[i][j] |= (_alleleCountTable.get(k,i,j) != 0);
351  }
352  //add mean nb of alleles per locus for Patch k to the pop mean
353  pop_mean += patch_mean/nbLoc;
354  }
355 
356  _nb_all_local = (nbpatch ? pop_mean/nbpatch : nanf("NULL"));
357 
358  _nb_all_global = 0;
359 
360  for(i = 0; i < nbLoc; ++i)
361  for(j = 0; j < nbAll;++j)
362  _nb_all_global += pop_div[i][j];
363 
364  _nb_all_global /= nbLoc;
365 
366  for(i = 0; i < nbLoc; ++i)
367  delete [] pop_div[i];
368  delete [] pop_div;
369 
370  //number of fixed loci, local and global counters:
371  _fix_loc_local = 0;
372 
373  for (k = 0; k < patchNbr; ++k)
374  for (i = 0; i < nbLoc; ++i)
375  for (j = 0; j < nbAll; ++j)
376  _fix_loc_local += ( _alleleFreqTable.get(k, i, j) == 1 );
377 
378  _fix_loc_local /= nbpatch;
379 
380  _fix_loc_global = 0;
381 
382  //globally:
383  for (i = 0; i < nbLoc; ++i)
384  for (j = 0; j < nbAll; ++j)
385  _fix_loc_global += ( _globalAlleleFreq.get(i, j) == 1 );
386 
387 }

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

◆ setNeiGeneticDistance()

void TTNeutralGenesSH::setNeiGeneticDistance ( age_t  AGE)
1189 {
1190  //see: Nei, M. 1978. Genetics 89:583-590
1191  unsigned int patchNbr = this->_pop->getPatchNbr();
1192  unsigned int nb_locus = this->_SHLinkedTrait->get_locus_num();
1193  unsigned int nb_allele = this->_SHLinkedTrait->get_allele_num();
1194 
1195  setAlleleTables(AGE);
1196 
1197  if(_D == NULL)
1198 
1199  _D = new TMatrix(patchNbr, patchNbr);
1200 
1201  else if( _D->length() != patchNbr*patchNbr )
1202 
1203  _D->reset(patchNbr, patchNbr);
1204 
1205  _D->assign(nanf("NULL"));
1206 
1207  double num, denom, p1, p2, sum;
1208  double *sum_square = new double [patchNbr];
1209  double *sample_size = new double [patchNbr];
1210 
1211  for(unsigned int i = 0; i < patchNbr; ++i) {
1212 
1213  sample_size[i] = 2 * _pop->size(AGE, i);
1214 
1215  sum_square[i] = 0;
1216 
1217  if( !sample_size[i] ) continue;
1218 
1219  for(unsigned int l = 0; l < nb_locus; ++l) {
1220 
1221  sum = 0;
1222 
1223  for(unsigned int u = 0; u < nb_allele ; ++u) {
1224  p1 = _alleleFreqTable.get(i, l, u);
1225  sum += p1*p1;
1226  }
1227 
1228  sum_square[i] += (sample_size[i] * sum -1) / (sample_size[i] -1); //unbiased estimate, equ. 6 in Nei 1978 (Genetics)
1229  }
1230  }
1231 
1232  unsigned int pairs = 0;
1233  double Dpair;
1234  _meanD = 0;
1235 
1236  for(unsigned int i = 0; i < patchNbr-1; ++i) {
1237 
1238  if( !sample_size[i] ) continue;
1239 
1240  for(unsigned int j = i + 1; j < patchNbr; ++j) {
1241 
1242  if( !sample_size[j] ) continue;
1243 
1244  num = 0;
1245  pairs++;
1246 
1247  for(unsigned int l = 0; l < nb_locus; ++l) {
1248  for(unsigned int u = 0; u < nb_allele ; ++u) {
1249  p1 = _alleleFreqTable.get(i, l, u);
1250  p2 = _alleleFreqTable.get(j, l, u);
1251  num += p1*p2;
1252  }
1253  }
1254 
1255  denom = sqrt(sum_square[i] * sum_square[j]);
1256  Dpair = -log(num/denom);
1257  _meanD += Dpair;
1258  _D->set(i, j, Dpair);
1259  }
1260  }
1261  _meanD /= pairs;
1262 
1263  delete [] sum_square;
1264  delete [] sample_size;
1265 }

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

◆ setNeiGeneticDistanceRecorders()

void TTNeutralGenesSH::setNeiGeneticDistanceRecorders ( age_t  AGE,
bool  pairwise 
)
2544 {
2545  string name, sub_name;
2546  unsigned int nbpatch = _pop->getPatchNbr();
2547  string prefix = (AGE == ADULTS ? "adlt." : "off.");
2548  unsigned int scale = (unsigned int)pow(10.0, (int)log10((float)nbpatch) + 1);
2549 
2550  void (TTNeutralGenesSH::* setter) () = (AGE == ADULTS ?
2553 
2554  sub_name = prefix + "D";
2555  add("Average between pop Nei's D", sub_name, AGE, 0, 0,
2557 
2558  if(pairwise) {
2559 
2560  for(unsigned int i = 0; i < nbpatch -1; i++){
2561  for(unsigned int j = i+1; j < nbpatch; j++) {
2562  name = "Nei's D between pop" + tstring::int2str(i+1) + " and pop" + tstring::int2str(j+1);
2563  sub_name = prefix + "D.p" + tstring::int2str(i+1) + "p" + tstring::int2str(j+1);
2564  add(name, sub_name, AGE, i*scale + j, 0, 0, &TTNeutralGenesSH::getNeiGeneticDistance, 0, 0);
2565  }
2566  }
2567  }
2568 }
double getNeiGeneticDistance(unsigned int i)
Definition: ttneutralgenes.h:546
void setOffsprgNeiGeneticDistance()
Definition: ttneutralgenes.h:544
void setAdltNeiGeneticDistance()
Definition: ttneutralgenes.h:543
double getMeanNeiGeneticDistance()
Definition: ttneutralgenes.h:551

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

Referenced by setStatRecorders().

◆ setOffsprgCoaBetween()

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

References OFFSx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ setOffsprgCoaMatrix()

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

References OFFSx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ setOffsprgCoaWithin()

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

References OFFSx, and setCoaMatrix().

Referenced by setCoaMatrixRecorders(), and setStatRecorders().

◆ setOffsprgFstat()

void TTNeutralGenesSH::setOffsprgFstat ( )
inline
456 {setFstat(OFFSPRG);}

References OFFSPRG, and setFstat().

Referenced by setFstatRecorders().

◆ setOffsprgFstat2()

void TTNeutralGenesSH::setOffsprgFstat2 ( )
inline
471 {setFstat2(OFFSPRG);}

References OFFSPRG, and setFstat2().

Referenced by setFstat2Recorders().

◆ setOffsprgFstBetween()

void TTNeutralGenesSH::setOffsprgFstBetween ( )
inline
439 {setFstMatrix(OFFSPRG, 2);}

References OFFSPRG, and setFstMatrix().

Referenced by setFstMatrixRecorders().

◆ setOffsprgFstMatrix()

void TTNeutralGenesSH::setOffsprgFstMatrix ( )
inline
437 {setFstMatrix(OFFSPRG, 3);}

References OFFSPRG, and setFstMatrix().

Referenced by setFstMatrixRecorders().

◆ setOffsprgFstWithin()

void TTNeutralGenesSH::setOffsprgFstWithin ( )
inline
438 {setFstMatrix(OFFSPRG, 1);}

References OFFSPRG, and setFstMatrix().

Referenced by setFstMatrixRecorders().

◆ setOffsprgNeiGeneticDistance()

void TTNeutralGenesSH::setOffsprgNeiGeneticDistance ( )
inline

◆ setOffspringAlleleFreq()

void TTNeutralGenesSH::setOffspringAlleleFreq ( )
inline

References OFFSPRG, and setAlleleTables().

Referenced by setFreqRecorders().

◆ setOffspringFstatWeirCockerham()

void TTNeutralGenesSH::setOffspringFstatWeirCockerham ( )
inline

◆ setOffspringHeterozygosity()

void TTNeutralGenesSH::setOffspringHeterozygosity ( )
inline

References OFFSPRG, and setHeterozygosity().

Referenced by setFreqRecorders().

◆ setSibCoa()

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

References Individual::getFatherID(), Individual::getMotherID(), and Individual::getTrait().

◆ setSibStats()

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

References FEM, Patch::get(), MAL, OFFSx, and Patch::size().

Referenced by setStatRecorders().

◆ setStatRecorders()

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

Implements StatHandlerBase.

2079 {
2080 #ifdef _DEBUG_
2081  message("-TTNeutralGenesSH::setStatRecorders ");
2082 #endif
2083  if(token == "coa") {
2084 
2085  add("Wtn Patch Coancestry (offsprg)","off.theta",OFFSPRG,0,0,
2087  add("Btn Patch Coancestry (offsprg)","off.alpha",OFFSPRG,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
2088 
2089  add("Wtn Patch Coancestry (adult)","adlt.theta",ADULTS,0,0,
2091  add("Btn Patch Coancestry (adult)","adlt.alpha",ADULTS,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
2092 
2093  } else if(token == "adlt.coa") {
2094 
2095  add("Wtn Patch Coancestry (adult)","adlt.theta",ADULTS,0,0,
2097  add("Btn Patch Coancestry (adult)","adlt.alpha",ADULTS,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
2098 
2099  } else if(token == "off.coa") {
2100 
2101  add("Wtn Patch Coancestry (offsprg)","off.theta",OFFSPRG,0,0,
2103  add("Btn Patch Coancestry (offsprg)","off.alpha",OFFSPRG,0,0,&TTNeutralGenesSH::getMeanAlpha,0,0,0);
2104 
2105  } else if(token == "adlt.coa.persex") {
2106 
2107  add("Female Theta (adult)","adlt.theta",ADULTS,0,0,
2109  add("Female Theta (adult)","adlt.thetaFF",ADULTS,0,0,&TTNeutralGenesSH::getTheta_FF,0,0,0);
2110  add("Male Theta (adult)","adlt.thetaMM",ADULTS,0,0,&TTNeutralGenesSH::getTheta_MM,0,0,0);
2111  add("Female-Male Theta (adult)","adlt.thetaFM",ADULTS,0,0,&TTNeutralGenesSH::getTheta_FM,0,0,0);
2112 
2113  } else if(token == "adlt.coa.within") {
2114 
2115  add("Wtn Patch Coancestry (adult)","adlt.theta",ADULTS,0,0,
2117 
2118  } else if(token == "off.coa.within") {
2119 
2120  add("Wtn Patch Coancestry (offsprg)","off.theta",OFFSPRG,0,0,
2122 
2123  } else if(token == "adlt.coa.between") {
2124 
2125  add("Btn Patch Coancestry (adult)","adlt.alpha",ADULTS,0,0,
2127 
2128  } else if(token == "off.coa.between") {
2129 
2130  add("Btn Patch Coancestry (offsprg)","off.alpha",OFFSPRG,0,0,
2132 
2133  } else if(token == "coa.matrix") {
2134 
2137 
2138  } else if(token == "off.coa.matrix") {
2139 
2141 
2142  } else if(token == "adlt.coa.matrix") {
2143 
2145 
2146  } else if(token == "coa.matrix.within") {
2147 
2150 
2151  } else if(token == "off.coa.matrix.within") {
2152 
2154 
2155  } else if(token == "adlt.coa.matrix.within") {
2156 
2158 
2159  } else if(token == "sibcoa") {
2160 
2161  add("Proportion of full-sib offspring","prop.fsib",OFFSPRG,0,0,
2163  add("Proportion of paternal half-sib ","prop.phsib",OFFSPRG,1,0,0,&TTNeutralGenesSH::getSibProportions,0,0);
2164  add("Proportion of maternal half-sib ","prop.mhsib",OFFSPRG,2,0,0, &TTNeutralGenesSH::getSibProportions,0,0);
2165  add("Proportion of non-sib offspring ","prop.nsib",OFFSPRG,3,0,0,&TTNeutralGenesSH::getSibProportions,0,0);
2166  add("Coancestry of full-sib offspring","coa.fsib",OFFSPRG,0,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
2167  add("Coancestry of paternal half-sib ","coa.phsib",OFFSPRG,1,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
2168  add("Coancestry of maternal half-sib ","coa.mhsib",OFFSPRG,2,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
2169  add("Coancestry of non-sib offspring ","coa.nsib",OFFSPRG,3,0,0,&TTNeutralGenesSH::getSibCoaMeans,0,0);
2170 
2171  } else if (token == "ntrl.freq") {
2172 
2175 
2176  } else if (token == "off.ntrl.freq") {
2177 
2179 
2180  } else if (token == "adlt.ntrl.freq") {
2181 
2183 
2184  // } else if (token == "ntrl.freq.patch") {
2185  //
2186  // setFreqRecordersPerPatch(ALL);
2187  //
2188  // } else if (token == "off.ntrl.freq.patch") {
2189  //
2190  // setFreqRecordersPerPatch(OFFSPRG);
2191  //
2192  // } else if (token == "adlt.ntrl.freq.patch") {
2193  //
2194  // setFreqRecordersPerPatch(ADULTS);
2195 
2196  } else if(token == "offsprgfstat" || token == "off.fstat") {
2197 
2199 
2200  } else if(token == "adultfstat" || token == "adlt.fstat") {
2201 
2203 
2204  } else if(token == "fstat") {
2205 
2207 
2208  } else if(token == "off.fstat2") {
2209 
2211 
2212  } else if(token == "adlt.fstat2") {
2213 
2215 
2216  } else if(token == "fstat2") {
2217 
2219 
2220  } else if(token == "fstWC" || token == "fstatWC") {
2221 
2223 
2224  } else if(token == "off.fstWC" || token == "off.fstatWC") {
2225 
2227 
2228  } else if(token == "adlt.fstWC" || token == "adlt.fstatWC") {
2229 
2231 
2232  } else if(token == "weighted.fst") {
2233 
2236 
2237  } else if(token == "off.weighted.fst") {
2238 
2240 
2241  } else if(token == "adlt.weighted.fst") {
2242 
2244 
2245  } else if(token == "weighted.fst.matrix") {
2246 
2249 
2250  } else if(token == "off.weighted.fst.matrix") {
2251 
2253 
2254  } else if(token == "adlt.weighted.fst.matrix") {
2255 
2257 
2258  } else if(token == "weighted.fst.within") {
2259 
2262 
2263  } else if(token == "off.weighted.fst.within") {
2264 
2266 
2267  } else if(token == "adlt.weighted.fst.within") {
2268 
2270 
2271  } else if(token == "adlt.NeiDistance") {
2272 
2274 
2275  } else if(token == "off.NeiDistance") {
2276 
2278 
2279  } else if(token == "NeiDistance") {
2280 
2283 
2284  } else if(token == "mean.NeiDistance") {
2285 
2288 
2289  } else if(token == "adlt.mean.NeiDistance") {
2290 
2292 
2293  } else if(token == "off.mean.NeiDistance") {
2294 
2296 
2297  } else if(token == "Dxy") {
2298 
2299  setDxyRecorders(ADULTS, false);
2300  setDxyRecorders(OFFSPRG, false);
2301 
2302  } else if(token == "off.Dxy") {
2303 
2304  setDxyRecorders(OFFSPRG, false);
2305 
2306  } else if(token == "adlt.Dxy") {
2307 
2308  setDxyRecorders(ADULTS, false);
2309 
2310  } else if(token == "Dxy.patch") {
2311 
2312  setDxyRecorders(ADULTS, true);
2313  setDxyRecorders(OFFSPRG, true);
2314 
2315  } else if(token == "off.Dxy.patch") {
2316 
2317  setDxyRecorders(OFFSPRG, true);
2318 
2319  } else if(token == "adlt.Dxy.patch") {
2320 
2321  setDxyRecorders(ADULTS, true);
2322 
2323  } else
2324  return false;
2325 
2326  return true;
2327 }
double getSibCoaMeans(unsigned int i)
Definition: ttneutralgenes.h:538
double getSibProportions(unsigned int i)
Definition: ttneutralgenes.h:537
void setSibStats()
Definition: stats_coa.cc:278
void setDxyRecorders(age_t AGE, bool patchwise)
Definition: ttneutralgenes.cc:2573
void setCoaMatrixRecorders(age_t AGE, unsigned char dim)
Definition: ttneutralgenes.cc:2331
void setFstatRecorders(age_t AGE)
Definition: ttneutralgenes.cc:2419
void setFreqRecorders(age_t AGE)
Definition: ttneutralgenes.cc:2387
double getTheta_FF()
Gives the mean within females coancestry coefficient.
Definition: ttneutralgenes.h:530
double getTheta_MM()
Gives the mean within males coancestry coefficient.
Definition: ttneutralgenes.h:532
double getTheta_FM()
Gives the mean between males and females coancestry coefficient.
Definition: ttneutralgenes.h:534
void setFstatWCRecorders(age_t AGE)
Definition: ttneutralgenes.cc:2475
void setFstMatrixRecorders(age_t AGE, unsigned char dim)
Definition: ttneutralgenes.cc:2495
void setNeiGeneticDistanceRecorders(age_t AGE, bool pairwise)
Definition: ttneutralgenes.cc:2543
void setAdults_Theta()
Definition: stats_coa.cc:202
void setFstat2Recorders(age_t AGE)
Definition: ttneutralgenes.cc:2453
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

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

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

◆ _is_diallelic_bitstring

bool TTNeutralGenesSH::_is_diallelic_bitstring
private

Referenced by init(), setAlleleTables(), setHo(), and setHo2().

◆ _mean_alpha

double TTNeutralGenesSH::_mean_alpha
private

Referenced by getMeanAlpha().

◆ _mean_theta

double TTNeutralGenesSH::_mean_theta
private

Referenced by getMeanTheta().

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

Referenced by getSibCoaMeans().

◆ _sib_prop

double TTNeutralGenesSH::_sib_prop[4]
private

Kinship classes proportions.

Referenced by getSibProportions().

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

◆ Theta_FM

double TTNeutralGenesSH::Theta_FM
private

Referenced by getTheta_FM().

◆ Theta_MM

double TTNeutralGenesSH::Theta_MM
private

Referenced by getTheta_MM().


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

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

Locations of visitors to this page
Catalogued on GSR