Nemo  2.4.0
Simulate forward-in-time genetic evolution in a spatially explicit, individual-based stochastic simulator
ttquanti.h
Go to the documentation of this file.
1 /*
2 *
3 * @file ttquanti.h
4 * Nemo2
5 *
6 * Copyright (C) 2006-2026 The Authors
7 *
8 * This file is part of Nemo
9 *
10 * Nemo is free software; you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation; either version 2 of the License, or
13 * (at your option) any later version.
14 *
15 * Nemo is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with this program; if not, write to the Free Software
22 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23 *
24 * created on @date 14.11.2005
25 *
26 * @author fred
27 */
28 
29 #ifndef TTQUANTI_H
30 #define TTQUANTI_H
31 
32 #include <cmath>
33 #include <vector>
34 #include <set>
35 #include "ttrait_with_map.h"
36 #include "bitstring.h"
37 #include "filehandler.h"
38 #include "stathandler.h"
39 #include "metapop.h"
40 #include "datatable.h"
41 #include "Uniform.h"
42 
43 #ifdef HAS_GSL
44 #include <gsl/gsl_vector.h>
45 #include <gsl/gsl_matrix.h>
46 #include <gsl/gsl_eigen.h>
47 #endif
48 
49 class TTQuantiSH;
50 class TTQuantiFH;
51 class TTQFreqExtractor;
52 class TTQOhtaStats;
53 class TProtoQuanti;
54 
55 // ------------------------------------------------------------------------------
59 // ------------------------------------------------------------------------------
60 class TTQuanti : public TTrait {
61 
62 public:
64  {
65 // message("CONSTRUCTOR::TTQuanti::CONSTRUCTOR\n");
66  }
67 
69 // { message("COPY CONSTRUCTOR::TTQuanti::COPY CONSTRUCTOR\n"); }
70 
71  virtual ~TTQuanti () {}
72 
73  virtual trait_t get_type () const {return QUANT;}
74  virtual void mutate ();
75  virtual void inherit (const TTrait* mother, const TTrait* father);
76  virtual void* set_trait (void* value) {return value;}
77  virtual void set_value ();
78  virtual void* getValue () const {return _phenotypes;}
79 
80  // TTQuanti interface:
81  virtual double get_additive_genotype (const unsigned int trait) const = 0;
82  virtual double get_dominant_genotype (const unsigned int trait) const = 0;
83  virtual double get_full_genotype (unsigned int trait) = 0;
84  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_pos,
85  unsigned int to_pos, const TTQuanti *parent) = 0;
86  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
87  const TTQuanti *parent) = 0;
88  virtual void mutate_add (unsigned int position, unsigned int allele, double value) = 0;
89  virtual void mutate_inplace (unsigned int position, unsigned int allele, double value) = 0;
90 
91  // interface to di-allelic architecture
92  virtual bool get_allele_bit (unsigned int position, unsigned int allele) const = 0;
93  virtual void set_allele_bit (unsigned int position, unsigned int allele, bool value) = 0;
94 
95  //accessors
96  void set_proto (TProtoQuanti* proto) {_myProto = proto;}
98  double get_phenotype (unsigned int trai);
99  void set_phenotype (unsigned int trait, double value);
101 
102 protected:
103 
104  double *_phenotypes;
107 };
108 // ------------------------------------------------------------------------------
112 // ------------------------------------------------------------------------------
114 
115 public:
116 
118  {
119 // message("CONSTRUCTOR::TTQuanti_continuous::CONSTRUCTOR\n");
120  }
121 
123  {
124 // message("COPY CONSTRUCTOR::TTQuanti_continuous::COPY CONSTRUCTOR\n");
125  }
126 
127  virtual ~TTQuanti_continuous () {
128 // message("DESTRUCTOR::TTQuanti_continuous::DESTRUCTOR\n");
129  reset();
130  }
131 
132  //implements TTrait:
133  virtual void reset ();
134  virtual void init ();
135  virtual void set_sequence (void** seq);
136  virtual void** get_sequence () const {return (void**)_sequence;}
137  virtual unsigned int get_allele (int loc, int all) const ;
138  virtual double get_allele_value (int loc, int all) const;
139  virtual void set_allele_value (unsigned int locus, unsigned int allele, double value);
140 
141  virtual TTQuanti_continuous& operator= (const TTrait& T);
142  virtual bool operator== (const TTrait& T);
143  virtual bool operator!= (const TTrait& T);
144 
145  //implements StorableComponent:
146  virtual void store_data ( BinaryStorageBuffer* saver );
147  virtual bool retrieve_data ( BinaryStorageBuffer* reader );
148 
149  //implements TTQuanti:
150  virtual double get_full_genotype (unsigned int trait);
151  virtual void set_allele (int locus, int allele, double value) {_sequence[allele][locus] = value;}
152  virtual void mutate_add (unsigned int position, unsigned int allele, double value)
153  {_sequence[allele][position] += value;}
154  virtual void mutate_inplace (unsigned int position, unsigned int allele, double value)
155  {_sequence[allele][position] = value;}
156  virtual bool get_allele_bit (unsigned int position, unsigned int allele) const {return false;}
157  virtual void set_allele_bit (unsigned int position, unsigned int allele, bool value) {}
158 
159 
160 protected:
161 
162  double **_sequence;
163 };
164 // Deriving specialized classes
165 // ------------------------------------------------------------------------------
169 // ------------------------------------------------------------------------------
171 
172 public:
174 
176 
178 
179  //TTtrait implementation:
180  virtual void init_sequence ();
182  virtual void show_up();
183 
184  //TTQuanti implementation:
185  virtual double get_additive_genotype (const unsigned int trait) const;
186  virtual double get_dominant_genotype (const unsigned int trait) const;
187  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_locus,
188  unsigned int to_locus, const TTQuanti *parent);
189  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
190  const TTQuanti *parent);
191 };
192 // ------------------------------------------------------------------------------
196 // ------------------------------------------------------------------------------
198 
199 public:
201 
203 
205 
206  unsigned int get_sequence_block_size (unsigned int from, unsigned int to);
207 
208  //TTtrait implementation:
209  virtual void init_sequence ();
211  virtual void show_up();
212 
213  //TTQuanti implementation:
214  virtual double get_additive_genotype (const unsigned int trait) const;
215  virtual double get_dominant_genotype (const unsigned int trait) const;
216  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_locus,
217  unsigned int to_locus, const TTQuanti *parent);
218  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
219  const TTQuanti *parent);
220 };
221 // ------------------------------------------------------------------------------
225 // ------------------------------------------------------------------------------
227 
228 public:
230 
232 
234 
235  //TTtrait implementation:
236  virtual void init_sequence ();
238  virtual void show_up();
239 
240  //TTQuanti implementation:
241  virtual double get_additive_genotype (const unsigned int trait) const;
242  virtual double get_dominant_genotype (const unsigned int trait) const;
243  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_locus,
244  unsigned int to_locus, const TTQuanti *parent);
245  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
246  const TTQuanti *parent);
247 };
248 // ------------------------------------------------------------------------------
252 // ------------------------------------------------------------------------------
254 
255 public:
257 
259 
261 
262  //TTtrait implementation:
263  virtual void init_sequence ();
265  virtual void show_up();
266 
267  //TTQuanti implementation:
268  virtual double get_additive_genotype (const unsigned int trait) const;
269  virtual double get_dominant_genotype (const unsigned int trait) const;
270  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_locus,
271  unsigned int to_locus, const TTQuanti *parent);
272  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
273  const TTQuanti *parent);
274 };
275 // ------------------------------------------------------------------------------
279 // ------------------------------------------------------------------------------
280 class TTQuanti_diallelic : public TTQuanti {
281 
282 public:
283 
285 
287 
288  virtual ~TTQuanti_diallelic () {reset();}
289 
290  virtual bool get_allele_bit (unsigned int position, unsigned int allele) const;
291  virtual void set_allele_bit (unsigned int position, unsigned int allele, bool value);
292  void inherit_free (const TTrait* mother, const TTrait* father);
293 
294  //implements TTrait:
295  virtual void reset ();
296  virtual void init ();
297  virtual void set_sequence (void** seq);
298  virtual void** get_sequence () const {return (void**)_sequence;}
299  virtual unsigned int get_allele (int loc, int all) const ;
300  virtual double get_allele_value (int locus, int allele) const;
301  virtual void set_allele_value (unsigned int locus, unsigned int allele, double value);
302 
303  virtual TTQuanti_diallelic& operator= (const TTrait& T);
304  virtual bool operator== (const TTrait& T);
305  virtual bool operator!= (const TTrait& T);
306 
307  //implements StorableComponent:
308  virtual void store_data ( BinaryStorageBuffer* saver );
309  virtual bool retrieve_data ( BinaryStorageBuffer* reader );
310 
311  //implements TTQuanti:
312  virtual double get_full_genotype (unsigned int trait);
313  virtual void set_allele (int locus, int allele, double value) {_sequence[allele][locus] = bool(value);}
314  virtual void mutate_add (unsigned int position, unsigned int allele, double value)
315  {_sequence[allele][position] = bool(value);}
316  virtual void mutate_inplace (unsigned int position, unsigned int allele, double value)
317  {_sequence[allele][position] ^= 1;}
318 
319 protected:
320 
321  unsigned char** _sequence;
322 };
323 // ------------------------------------------------------------------------------
327 // ------------------------------------------------------------------------------
329 
330 public:
332 
334 
336 
337  //TTtrait implementation:
338  virtual void init_sequence ();
340  virtual void show_up();
341 
342  //TTQuanti implementation:
343  virtual double get_additive_genotype (const unsigned int trait) const;
344  virtual double get_dominant_genotype (const unsigned int trait) const;
345  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_locus,
346  unsigned int to_locus, const TTQuanti *parent);
347  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
348  const TTQuanti *parent);
349 };
350 // ------------------------------------------------------------------------------
354 // ------------------------------------------------------------------------------
356 
357 public:
359 
361 
363 
364  //TTtrait implementation:
365  virtual void init_sequence ();
367  virtual void show_up();
368 
369  //TTQuanti implementation:
370  virtual double get_additive_genotype (const unsigned int trait) const;
371  virtual double get_dominant_genotype (const unsigned int trait) const;
372  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_locus,
373  unsigned int to_locus, const TTQuanti *parent);
374  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
375  const TTQuanti *parent);
376 };
377 // ------------------------------------------------------------------------------
381 // ------------------------------------------------------------------------------
383 
384 public:
386 
388 
390 
391  unsigned int get_sequence_block_size (unsigned int from, unsigned int to);
392 
393  //TTtrait implementation:
394  virtual void init_sequence ();
396  virtual void show_up();
397 
398  //TTQuanti implementation:
399  virtual double get_additive_genotype (const unsigned int trait) const;
400  virtual double get_dominant_genotype (const unsigned int trait) const;
401  virtual void copy_sequence_block (sex_t SEX, unsigned int strand, unsigned int from_locus,
402  unsigned int to_locus, const TTQuanti *parent);
403 
404  virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at,
405  const TTQuanti *parent);
406 };
407 
408 
409 // ------------------------------------------------------------------------------
413 // ------------------------------------------------------------------------------
414 class TProtoQuanti : public TTProtoWithMap {
415 
416 public:
417 
418  TProtoQuanti ();
419  TProtoQuanti (const TProtoQuanti& T);
420  virtual ~TProtoQuanti ();
421 
422  unsigned int get_num_traits() {return _num_traits;}
423  unsigned int get_num_locus() {return _num_locus;}
424  unsigned int get_num_locus (unsigned int trait) {return _trait_table[trait].size();}
425  unsigned int get_pleiotropy_type() {return _pleio_type;}
426  unsigned int get_seq_length () {return _seq_length;}
429  vector<double> get_env_var () {return _eVariance;}
430  vector<double> get_heritability () {return _h2;}
431  unsigned int get_h2_setTime () {return _h2_setTime;}
432  bool get_h2_isBroad () {return _h2_isBroad;}
433  unsigned int get_allele_model () {return _allele_model;}
435  double get_diallele_value (unsigned int locus, unsigned int allele)
436  {return _allele_value.get(locus, allele);}
437  double get_seq_diallele_value (unsigned int position, unsigned int allele)
438  {return _sequence_diallele_values[allele][position];}
441  double get_equal_val_0 () {return _equal_val_0;}
442  double get_equal_val_1 () {return _equal_val_1;}
444  const bitstring& get_trait_mask (unsigned int trait) const
445  {return _trait_masks[trait];}
446  double get_init_value (unsigned int i) {return _init_value.get(0,i); }
447  double get_init_variance (unsigned int i) {return _init_variance.get(0,i);}
448  unsigned int get_doInitMutation () {return _doInitMutation;}
450 
452  unsigned int get_locus_seq_pos (unsigned int loc, unsigned int trait)
453  {return _trait_table[trait][loc];}
454  unsigned int get_locus_ID (unsigned int locus, unsigned int trait)
455  {return _trait_locus_table[trait][locus];}
456  unsigned int get_locus_PD (unsigned int locus) {return _locus_table[locus][1];}
457  unsigned int get_locus_start_pos (unsigned int locus) {return _locus_table[locus][0];}
458  unsigned int get_sequence_block_size (unsigned int from, unsigned int to) {
459  return (to < _num_locus ? _locus_table[to][0] : _seq_length)
460  - _locus_table[from][0];
461  }
462 
463  void set_eVarianceSD (unsigned int trait, double SD);
464  void set_init_values (const double *values, unsigned int nval);
465  void set_trait_value_func_ptr (bool withVe);
466  double set_trait_value_VE (const TTQuanti* ind, const unsigned int trait);
467  double set_trait_value_noVE (const TTQuanti* ind, const unsigned int trait);
468 
469  double set_genotype_value_additive (const TTQuanti* ind, const unsigned int trait);
470  double set_genotype_value_dominance (const TTQuanti* ind, const unsigned int trait);
471 
472  int get_allele_position (const unsigned int locus, const unsigned int trait);
473  double get_allele_value (const TTQuanti* ind, const unsigned int allele, const unsigned int locus, const unsigned int trait);
474  double get_genotype_with_dominance (const double a1, const double a2, const unsigned int locus, const unsigned int trait);
475  double get_genotype_dominance_h (double a1, double a2, double h);
476  double get_genotype_dominance_k (double a1, double a2, double k);
477  unsigned int get_dominance_model () {return _dominance_model;}
478  double get_dominance (unsigned int locus, unsigned int trait)
479  {return _dominance_effects.get(trait, locus); }
481 
482  double get_trait_mutation_variance (unsigned int trait);
483 
484  double get_genotypic_value (const TTQuanti* ind, const unsigned int trait) {
485  return (this->*_set_genotype_func_ptr)(ind, trait);
486  }
487 
488  double get_phenotypic_value (const TTQuanti* ind, const unsigned int trait) {
489  return (this->*_set_trait_value_func_ptr)(ind, trait);
490  }
491 
492 
493  // Parameter setters
494  bool setHeritabilityParams ();
495  bool setInitialValuesParams ();
496  bool setGeneticMapParams ();
497  unsigned int setAlleleModel ();
503  bool setMutationCorrelation ();
504  bool setDominanceParameters ();
505 
507 
510  bool readMatrixFromQuantiMutationMatrix (vector<vector<double>>& varmat);
514  void allocate_gsl_mutation_matrix_space (unsigned int num_locus);
516  gsl_matrix* set_gsl_mutation_matrix_from_sigma (unsigned int loc, unsigned int pleio_deg);
517  gsl_matrix* set_gsl_mutation_matrix (unsigned int pleio_deg, const vector<double>& varcov);
518  void set_mutation_matrix_decomposition (unsigned int loc, unsigned int pleio_deg);
519 
520  // mutation generation functions
521  // first set for continuum-of-allele model with a single mutation matrix across all loci:
522  double* getMutationEffectMultivariateGaussian (unsigned int loc);
523  double* getMutationEffectBivariateGaussian (unsigned int loc);
524  double* getMutationEffectUnivariateGaussian (unsigned int loc);
525  // different implementation for locus-specific mutation matrices:
526  double* getMutationEffectMultivariateGaussianLocSpec (unsigned int loc);
527  double* getMutationEffectBivariateGaussianLocSpec (unsigned int loc);
528  double* getMutationEffectUnivariateGaussianLocSpec (unsigned int loc);
529 
530 // double* getMutationEffectUnivariateDiallelic (unsigned int loc);
531  double* getMutationEffectBivariateDiallelic (unsigned int loc);
532  double* getMutationEffects (unsigned int loc) {
533  return (this->* _getMutationValues)(loc);
534  }
535  double* getMutationEffectsVarPleio (unsigned int loc) {
536  return (this->* _getMutationValuesVarPleio[loc]) (loc);
537  }
538  void inherit (sex_t SEX, TTQuanti* ind, const TTQuanti* parent);
539  void inherit_free (sex_t SEX, TTQuanti* ind, const TTQuanti* parent);
540  void inherit_low (sex_t SEX, TTQuanti* ind, const TTQuanti* parent);
541 
542  unsigned int get_num_mutations () {
543  return (unsigned int)RAND::Binomial(_mutation_rate, _2L);
544  }
545  void mutate (TTQuanti* ind);
546  void mutate_nill (TTQuanti* ind);
547  void mutate_full_pleio (TTQuanti* ind);
548  void mutate_var_pleio (TTQuanti* ind);
549  void mutate_no_pleio (TTQuanti* ind);
552  void mutate_inplace_no_pleio (TTQuanti* ind);
553 
554  void mutate_diallelic_pleio (TTQuanti* ind);
557 
558  // EPISTASIS
559  bool setEpistasisParameters ();
560 // double get_epistatic_mean () {return _epistatic_mean;}
561 // double get_epistatic_variance () {return _epistatic_variance;}
562 // double get_epistatic_prop_large () {return _epistatic_prop_large;}
563 // double get_epistatic_large_mean () {return _epistatic_large_mean;}
564 // double get_epistatic_large_variance () {return _epistatic_large_variance;}
565  unsigned int get_num_epi_coefs () {return _num_epi_coefs;}
566  bool do_epistasis () {return _epistasis;}
569 
572  virtual void reset () {TTProtoWithMap::reset();}
573  virtual TTQuanti* hatch ();
574  virtual TProtoQuanti* clone () {return new TProtoQuanti(*this);}
575  virtual trait_t get_type () const {return QUANT;}
577  virtual int get_phenotype_dimension () {return _num_traits;}
579  virtual int get_allele_number () {return (_allele_model > 2 ? 3: 2);} //3 just means more than di-allelic
581  virtual int get_locus_number () {return _num_locus;}
585  virtual void store_data ( BinaryStorageBuffer* saver )
586  {saver->store(&_seq_length,sizeof(int));}
587 
588  virtual bool retrieve_data ( BinaryStorageBuffer* reader )
589  {reader->read(&_seq_length,sizeof(int));return true;}
593  virtual bool setParameters();
594  virtual void loadFileServices ( FileServices* loader );
595  virtual void loadStatServices ( StatServices* loader );
596  virtual bool resetParameterFromSource (std::string param, SimComponent* cmpt) {return false;}
598 
599 protected:
600 
602  unsigned int _num_locus;
604  unsigned int _2L;
606  unsigned int _num_traits;
608  unsigned int _seq_length;
609  unsigned int _allele_model;
611 
612  //mutations:
615 
616  //Pleiotropy type: 0 = "no pleio", 1 = "full", 2 = "variable"
617  unsigned int _pleio_type;
618 
619 private:
620 
621  TMatrix _allele_value; // num locus x 2
623  bool _equal_effects; // true when all loci share the same (val_0, val_1)
624  double _equal_val_0; // shared ancestral allele value
625  double _equal_val_1; // shared derived allele value
626  double _equal_val_diff; // val_1 - val_0, precomputed
627  vector<bitstring> _trait_masks; // one per trait, bit i set iff position i belongs to that trait
628  TMatrix _init_value; // initial value of each trait
629  TMatrix _init_variance; // initial variance of each trait
630 // double* _phenotypes; // num traits
631 
632  //mutations:
633  TMatrix _mutation_correlation; // 1 x n_locus, one correlation per locus, applied to all trait pairs
634  TMatrix _mutation_sigma; // n_locus x n_trait, variance per locus, used in the bivariate case, and to initialize the m-matrix
635  gsl_matrix **_gsl_mutation_matrix; // one mutation matrix per locus
636  gsl_matrix **_evect;
637  gsl_vector **_eval;
638  gsl_vector **_effects_multivar;
639  gsl_vector **_ws;
640  double _effects_bivar[2];
641  unsigned int _doInitMutation;
642  bool _mutationVarianceIsLocusSpecific; // a single mutation matrix for all loci
643  bool _mutationEffectIsFixedDiAllele; // bitstring sequence not storing effects
644 
645  // dominance
646  unsigned int _dominance_model;
649 
650  //recombination:
654 
655  //heritability
656  vector<double> _eVariance;
657  vector<double> _h2;
658  unsigned int _h2_setTime;
660 
661  //variable pleiotropy
662 
665 
667  vector< vector<unsigned int> > _trait_table;
668 
670  vector< vector<unsigned int> > _trait_locus_table;
671 
675  vector< vector<unsigned int> > _locus_table;
676 
677 
678  // EPISTASIS
680  unsigned int _num_epi_coefs;
683 
684 
685  // function pointers
686 
689 
692 
694  double* (TProtoQuanti::* _getMutationValues) (unsigned int);
695 
697  vector< double* (TProtoQuanti::* ) (unsigned int) > _getMutationValuesVarPleio;
698 
700  double (TProtoQuanti::* _getGenotypeWithDominance) (double, double, double);
701 
702 
706  double (TProtoQuanti::* _set_trait_value_func_ptr) (const TTQuanti*, const unsigned int);
707 
708 
712  double (TProtoQuanti::* _set_genotype_func_ptr) (const TTQuanti*, const unsigned int);
713 
714  friend class TTQuanti;
715  friend class TTQuanti_continuous;
720 
726 };
727 
728 // ------------------------------------------------------------------------------
732 // ------------------------------------------------------------------------------
733 class TTQuantiSH : public TraitStatHandler<TProtoQuanti, TTQuantiSH> {
734 
735  double *_meanP, *_meanG, *_Va, *_Vg, *_Vb, *_Vp, *_covar;
736  double **_pmeanP, **_pmeanG, **_pVa, **_pVp, **_pcovar, **_peigval, **_peigvect;
737  bool _eVar;
738 
740 
741  gsl_matrix *_G, *_evec;
742  gsl_vector *_eval;
743  gsl_eigen_symmv_workspace *_ws;
744 
747 
748 public:
749 
752  _meanP(0), _meanG(0), _Va(0), _Vg(0), _Vb(0), _Vp(0), _covar(0),
753  _pmeanP(0), _pmeanG(0), _pVa(0), _pVp(0), _pcovar(0), _peigval(0), _peigvect(0),
754  _eVar(0), _num_locus(0), _num_trait(0),_patchNbr(0), _G(0), _evec(0),_eval(0),_ws(0),
755  _table_set_gen(999999), _table_set_age(999999), _table_set_repl(999999)
756  {}
757 
758  virtual ~TTQuantiSH() {resetPtrs();}
759 
760  void resetPtrs();
761 
762  virtual void init ( );
763 
764  virtual bool setStatRecorders (std::string& token);
765 
766  void addQuanti (age_t AGE);
767  void addQuantiPerPatch (age_t AGE);
768  void addAvgPerPatch (age_t AGE);
769  void addGenotPerPatch (age_t AGE);
770  void addVarPerPatch (age_t AGE);
771  void addCovarPerPatch (age_t AGE);
772  void addEigenPerPatch (age_t AGE);
773  void addEigenValuesPerPatch (age_t AGE);
774  void addEigenVect1PerPatch (age_t AGE);
775  void addSkewPerPatch (age_t AGE);
776 
777  void setDataTables (age_t AGE);
780  void setStats (age_t AGE);
781  double getMeanGenot (unsigned int i) {return _meanG[i];}
782  double getMeanPhenot (unsigned int i) {return _meanP[i];}
783  double getVa (unsigned int i) {return _Va[i];}
784  double getVg (unsigned int i) {return _Vg[i];}
785  double getVb (unsigned int i) {return _Vb[i];}
786  double getVp (unsigned int i) {return _Vp[i];}
787  double getQst (unsigned int i) {return _Vb[i]/(_Vb[i]+2*_Va[i]);}
788  double getCovar (unsigned int i) {return _covar[i];}
789 
790  double getMeanGenotPerPatch (unsigned int i, unsigned int p) {return _pmeanG[i][p];}
791  double getMeanPhenotPerPatch (unsigned int i, unsigned int p) {return _pmeanP[i][p];}
792  double getVaPerPatch (unsigned int i, unsigned int p) {return _pVa[i][p];}
793  double getVpPerPatch (unsigned int i, unsigned int p) {return _pVp[i][p];}
794  double getEigenValuePerPatch (unsigned int i, unsigned int p) {return _peigval[i][p];}
795  double getCovarPerPatch (unsigned int p, unsigned int i) {return _pcovar[p][i];}
796  double getEigenVectorEltPerPatch (unsigned int p, unsigned int v) {return _peigvect[p][v];}
797  double getSkewPerPatch (unsigned int i, unsigned int p);
798 
799  vector<double> getSNPalleleFreqInPatch (Patch* patch, const age_idx AGE);
800  vector<double> getVaWithDominance (Patch* curPop, const age_idx AGE);
801  vector<double> getVaNoDominance (Patch* curPop, const age_idx AGE);
802 
803 };
804 // ------------------------------------------------------------------------------
808 // ------------------------------------------------------------------------------
809 class TTQuantiFH : public TraitFileHandler<TProtoQuanti> {
810 
813 
814 public:
816  _has_genetic_map(0) {}
817  virtual ~TTQuantiFH(){}
818 
819  void setOutputOption (string opt);
820 
821  virtual void FHwrite ();
822  void write_TABLE ();
823  void print(ofstream& FH, Patch* patch, sex_t SEX, age_idx Ax, unsigned int print_gene, bool print_genotype, bool print_additive_genotype);
824  void write_PLINK ();
825  void print_PLINK_PED(ofstream& FH, age_idx Ax, Patch *patch);
826  void print_PLINK_FAM(ofstream& FH, age_idx Ax, Patch *patch);
827 
828  virtual void FHread (string& filename);
829 };
830 // ------------------------------------------------------------------------------
834 // ------------------------------------------------------------------------------
835 class TTQFreqExtractor : public TraitFileHandler<TProtoQuanti> {
836 
837  vector< string > _records;
838 
839 public:
841 
842  void resetTable();
843 
844  virtual ~TTQFreqExtractor () {}
845  virtual void FHwrite ();
846  virtual void FHread (string& filename) {}
847 };
848 
849 // ------------------------------------------------------------------------------
853 // ------------------------------------------------------------------------------
854 class TTQOhtaStats : public TraitFileHandler<TProtoQuanti> {
855 
857 
858 public:
860 
861  virtual ~TTQOhtaStats () {}
862  virtual void FHwrite ();
863  virtual void FHread (string& filename) {}
864 };
865 #endif
866 
Nemo2.
A class to store any kind of data in a char buffer before unloading it in a binary data file.
Definition: binarystoragebuffer.h:43
void read(void *out, unsigned int nb_bytes)
Definition: binarystoragebuffer.h:220
void store(void *stream, unsigned int nb_bytes)
Definition: binarystoragebuffer.cc:37
A class to manage the files associated with each components of the simulation.
Definition: fileservices.h:51
Second class in the metapopulation design structure, between the Metapop and Individual classes.
Definition: metapop.h:431
static double Binomial(double p, unsigned int n)
Definition: Uniform.h:489
Interface to all basic components of a simulation (traits, life cycle events, pop,...
Definition: simcomponent.h:44
The Service class used to manage the StatHandler objects.
Definition: statservices.h:49
A class to handle matrix in params, coerces matrix into a vector of same total size.
Definition: tmatrix.h:49
double get(unsigned int i, unsigned int j) const
Accessor to element at row i and column j.
Definition: tmatrix.h:192
TProtoQuanti.
Definition: ttquanti.h:414
bool setMutationSigmaFromQuantiMutationVariance()
Definition: ttquanti.cc:1921
bool _equal_effects
Definition: ttquanti.h:623
vector< double > _eVariance
Definition: ttquanti.h:656
double * getMutationEffectMultivariateGaussian(unsigned int loc)
Definition: ttquanti.cc:2389
unsigned int setAlleleModel()
Definition: ttquanti.cc:804
double set_genotype_value_dominance(const TTQuanti *ind, const unsigned int trait)
Definition: ttquanti.cc:2531
bool do_epistasis()
Definition: ttquanti.h:566
unsigned int _allele_model
Definition: ttquanti.h:609
TMatrix _dominance_effects
Definition: ttquanti.h:647
bool _epistasis
Definition: ttquanti.h:679
size_t _sizeofLocusType
Definition: ttquanti.h:653
void set_mutation_matrix_decomposition(unsigned int loc, unsigned int pleio_deg)
Definition: ttquanti.cc:2258
unsigned int _doInitMutation
Definition: ttquanti.h:641
unsigned int get_locus_PD(unsigned int locus)
Definition: ttquanti.h:456
void mutate_inplace_var_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2787
double set_genotype_value_additive(const TTQuanti *ind, const unsigned int trait)
Definition: ttquanti.cc:2524
void setTraitAndLocusTables_full_pleio()
Definition: ttquanti.cc:1272
virtual ~TProtoQuanti()
Definition: ttquanti.cc:223
double set_trait_value_VE(const TTQuanti *ind, const unsigned int trait)
Definition: ttquanti.cc:2510
virtual int get_phenotype_dimension()
Returns the dimension of the phenotype of the trait (size of the array accessed with TTrait::getValue...
Definition: ttquanti.h:577
void mutate_diallelic_var_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2846
void mutate_var_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2716
bool setContinuousMutationModel_var_pleio()
Definition: ttquanti.cc:1668
unsigned int _dominance_model
Definition: ttquanti.h:646
bool _equal_dom_coeff
Definition: ttquanti.h:648
double get_equal_val_diff()
Definition: ttquanti.h:443
TMatrix _mutation_correlation
Definition: ttquanti.h:633
double * getMutationEffectUnivariateGaussian(unsigned int loc)
Definition: ttquanti.cc:2432
gsl_matrix * set_gsl_mutation_matrix_from_sigma(unsigned int loc, unsigned int pleio_deg)
Definition: ttquanti.cc:2222
bool setMutationModel_var_pleio()
Definition: ttquanti.cc:947
double get_allele_value(const TTQuanti *ind, const unsigned int allele, const unsigned int locus, const unsigned int trait)
Definition: ttquanti.cc:2486
void(TProtoQuanti::* _inherit_fun_ptr)(sex_t, TTQuanti *, const TTQuanti *)
Pointer to inheritance functions: either inherit_free() (r=0.5), or inherit_low() (r<0....
Definition: ttquanti.h:688
bool setMutationModel_no_pleio()
Definition: ttquanti.cc:1143
bool _h2_isBroad
Definition: ttquanti.h:659
const TMatrix & get_diallele_values()
Definition: ttquanti.h:434
const TMatrix & get_pleio_matrix()
Definition: ttquanti.h:451
TTQuantiFH * _reader
Definition: ttquanti.h:723
bool setDominanceParameters()
Definition: ttquanti.cc:2310
unsigned int _h2_setTime
Definition: ttquanti.h:658
unsigned int _seq_length
Total number of allelic values stored in individual sequences, no trait x no locus.
Definition: ttquanti.h:608
double _mutation_rate
Definition: ttquanti.h:613
gsl_matrix ** _evect
Definition: ttquanti.h:636
void mutate_nill(TTQuanti *ind)
Definition: ttquanti.cc:2662
bool setContinuousMutationModel_no_pleio()
Definition: ttquanti.cc:1825
gsl_vector ** _ws
Definition: ttquanti.h:639
TMatrix _mutation_sigma
Definition: ttquanti.h:634
double get_dominance(unsigned int locus, unsigned int trait)
Definition: ttquanti.h:478
bool setInitialValuesParams()
Definition: ttquanti.cc:593
bool _mutationVarianceIsLocusSpecific
Definition: ttquanti.h:642
const TMatrix & get_epi_coefs() const
Definition: ttquanti.h:567
bool readMatrixFromQuantiMutationMatrix(vector< vector< double >> &varmat)
Definition: ttquanti.cc:2057
TMatrix _allele_value
Definition: ttquanti.h:621
bool setGeneticMapParams()
Definition: ttquanti.cc:648
double * getMutationEffects(unsigned int loc)
Definition: ttquanti.h:532
double * getMutationEffectBivariateGaussian(unsigned int loc)
Definition: ttquanti.cc:2413
void set_trait_value_func_ptr(bool withVe)
Definition: ttquanti.cc:678
TTQuantiSH * get_stater()
Definition: ttquanti.h:449
vector< double > _h2
Definition: ttquanti.h:657
unsigned int _num_epi_coefs
Definition: ttquanti.h:680
size_t get_locus_byte_size()
Definition: ttquanti.h:428
void mutate_diallelic_no_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2877
double get_genotypic_value(const TTQuanti *ind, const unsigned int trait)
Definition: ttquanti.h:484
unsigned int _num_traits
Number of traits.
Definition: ttquanti.h:606
double * _sequence_diallele_values[2]
Definition: ttquanti.h:622
double get_init_value(unsigned int i)
Definition: ttquanti.h:446
void mutate_inplace_full_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2762
size_t get_size_locus_type()
Definition: ttquanti.h:427
unsigned int get_locus_ID(unsigned int locus, unsigned int trait)
Definition: ttquanti.h:454
TTQuantiSH * _stats
Definition: ttquanti.h:721
virtual TTQuanti * hatch()
Definition: ttquanti.cc:2897
void set_init_values(const double *values, unsigned int nval)
Definition: ttquanti.cc:2500
void(TProtoQuanti::* _mutation_func_ptr)(TTQuanti *)
Pointer to mutation function, which depends on allele on model (HC, noHC, diallelic)
Definition: ttquanti.h:691
gsl_vector ** _effects_multivar
Definition: ttquanti.h:638
double(TProtoQuanti::* _set_genotype_func_ptr)(const TTQuanti *, const unsigned int)
Pointer to functions get_genotype_value_additive() or get_genotype_value_dominance() computing the ge...
Definition: ttquanti.h:712
void inherit(sex_t SEX, TTQuanti *ind, const TTQuanti *parent)
Definition: ttquanti.cc:2593
TMatrix _pleio_matx
Pleiotropy matrix provided in input (num locu X num trait).
Definition: ttquanti.h:664
unsigned int get_num_epi_coefs()
Definition: ttquanti.h:565
string _diallele_datatype
Definition: ttquanti.h:610
const bitstring & get_trait_mask(unsigned int trait) const
Definition: ttquanti.h:444
double get_seq_diallele_value(unsigned int position, unsigned int allele)
Definition: ttquanti.h:437
double _genomic_mutation_rate
Definition: ttquanti.h:614
gsl_matrix ** _gsl_mutation_matrix
Definition: ttquanti.h:635
bool has_equal_effects()
Definition: ttquanti.h:440
double * getMutationEffectBivariateDiallelic(unsigned int loc)
Definition: ttquanti.cc:2456
virtual bool retrieve_data(BinaryStorageBuffer *reader)
Definition: ttquanti.h:588
double get_diallele_value(unsigned int locus, unsigned int allele)
Definition: ttquanti.h:435
vector< vector< unsigned int > > _trait_locus_table
Table storing the locus id of each locus affecting each trait (num trait X (variable length/trait)).
Definition: ttquanti.h:670
unsigned int get_allele_model()
Definition: ttquanti.h:433
double get_init_variance(unsigned int i)
Definition: ttquanti.h:447
bool * _all_chooser
Definition: ttquanti.h:651
bool setMutationSigmaFromQuantiMutationVariance_no_pleio()
Definition: ttquanti.cc:1991
double(TProtoQuanti::* _getGenotypeWithDominance)(double, double, double)
Pointer to either dominance_h() or dominance_k() function computing the genotypic value with dominanc...
Definition: ttquanti.h:700
double _equal_val_diff
Definition: ttquanti.h:626
void allocate_gsl_mutation_matrix_space(unsigned int num_locus)
Definition: ttquanti.cc:2116
TMatrix _epistatic_coefs_matrix
Definition: ttquanti.h:681
vector< bitstring > _trait_masks
Definition: ttquanti.h:627
double _equal_val_0
Definition: ttquanti.h:624
TMatrix _init_variance
Definition: ttquanti.h:629
void inherit_low(sex_t SEX, TTQuanti *ind, const TTQuanti *parent)
Definition: ttquanti.cc:2600
virtual int get_locus_number()
Returns the number of locus.
Definition: ttquanti.h:581
void mutate(TTQuanti *ind)
Definition: ttquanti.cc:2655
unsigned int get_h2_setTime()
Definition: ttquanti.h:431
const TMatrix & get_dominance_effects()
Definition: ttquanti.h:480
virtual TProtoQuanti * clone()
Definition: ttquanti.h:574
void mutate_inplace_no_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2742
double * getMutationEffectUnivariateGaussianLocSpec(unsigned int loc)
Definition: ttquanti.cc:2440
unsigned int get_locus_seq_pos(unsigned int loc, unsigned int trait)
Definition: ttquanti.h:452
void mutate_diallelic_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2811
vector< vector< unsigned int > > _trait_table
Trait table, (num trait X (variable length/trait)), holds, for each trait, the array position of caus...
Definition: ttquanti.h:667
void setTraitAndLocusTables_no_pleio(TMatrix &mat)
Definition: ttquanti.cc:1304
gsl_matrix * set_gsl_mutation_matrix(unsigned int pleio_deg, const vector< double > &varcov)
Definition: ttquanti.cc:2185
double * getMutationEffectBivariateGaussianLocSpec(unsigned int loc)
Definition: ttquanti.cc:2423
vector< vector< unsigned int > > _locus_table
Locus table, num_locus x 2, first column holds the start position of the alleles of each locus in the...
Definition: ttquanti.h:675
unsigned int get_num_mutations()
Definition: ttquanti.h:542
vector< double > get_heritability()
Definition: ttquanti.h:430
const TMatrix & get_epi_coef_index() const
Definition: ttquanti.h:568
size_t _locusByteSize
Definition: ttquanti.h:652
TMatrix _init_value
Definition: ttquanti.h:628
double get_equal_val_0()
Definition: ttquanti.h:441
double set_trait_value_noVE(const TTQuanti *ind, const unsigned int trait)
Definition: ttquanti.cc:2517
unsigned int get_sequence_block_size(unsigned int from, unsigned int to)
Definition: ttquanti.h:458
unsigned int get_num_locus()
Definition: ttquanti.h:423
virtual void loadStatServices(StatServices *loader)
Definition: ttquanti.cc:3090
unsigned int get_seq_length()
Definition: ttquanti.h:426
bool setEpistasisParameters()
Definition: ttquanti.cc:688
double(TProtoQuanti::* _set_trait_value_func_ptr)(const TTQuanti *, const unsigned int)
Pointer to either set_trait_value_VE() or set_trait_value_noVE() to compute phenotypic values.
Definition: ttquanti.h:706
double get_phenotypic_value(const TTQuanti *ind, const unsigned int trait)
Definition: ttquanti.h:488
unsigned int get_pleiotropy_type()
Definition: ttquanti.h:425
bool setHeritabilityParams()
Definition: ttquanti.cc:437
virtual trait_t get_type() const
Definition: ttquanti.h:575
double *(TProtoQuanti::* _getMutationValues)(unsigned int)
Pointer to mutation allele value function, which depends on allele model and number of traits affecte...
Definition: ttquanti.h:694
double get_genotype_dominance_k(double a1, double a2, double k)
Definition: ttquanti.cc:2572
virtual bool resetParameterFromSource(std::string param, SimComponent *cmpt)
Definition: ttquanti.h:596
double get_genotype_dominance_h(double a1, double a2, double h)
Definition: ttquanti.cc:2556
double _effects_bivar[2]
Definition: ttquanti.h:640
virtual void reset()
Definition: ttquanti.h:572
double * getMutationEffectMultivariateGaussianLocSpec(unsigned int loc)
Definition: ttquanti.cc:2401
bool setDiallelicMutationModel()
Definition: ttquanti.cc:1347
unsigned int get_num_traits()
Definition: ttquanti.h:422
double get_equal_val_1()
Definition: ttquanti.h:442
gsl_vector ** _eval
Definition: ttquanti.h:637
unsigned int get_doInitMutation()
Definition: ttquanti.h:448
TMatrix _epistatic_coefs_indices
Definition: ttquanti.h:682
TTQuantiFH * _writer
Definition: ttquanti.h:722
virtual bool setParameters()
Definition: ttquanti.cc:248
bool has_equal_domCoeff()
Definition: ttquanti.h:439
bool setMutationCorrelation()
Definition: ttquanti.cc:1863
bool get_h2_isBroad()
Definition: ttquanti.h:432
virtual void loadFileServices(FileServices *loader)
Definition: ttquanti.cc:2976
void mutate_full_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2690
virtual void store_data(BinaryStorageBuffer *saver)
Definition: ttquanti.h:585
bool _mutationEffectIsFixedDiAllele
Definition: ttquanti.h:643
TTQOhtaStats * _ohtaStats
Definition: ttquanti.h:725
double get_trait_mutation_variance(unsigned int trait)
Definition: ttquanti.cc:2579
virtual int get_allele_number()
Returns the number of allele per locus.
Definition: ttquanti.h:579
unsigned int _2L
Diploid locus size, to save on useless operations during mutation.
Definition: ttquanti.h:604
unsigned int get_dominance_model()
Definition: ttquanti.h:477
bool setContinuousMutationModel_full_pleio()
Definition: ttquanti.cc:1501
double * getMutationEffectsVarPleio(unsigned int loc)
Definition: ttquanti.h:535
TProtoQuanti()
Definition: ttquanti.cc:63
bool setMutationModel_full_pleio()
Definition: ttquanti.cc:868
unsigned int _num_locus
Total number of loci, for all traits.
Definition: ttquanti.h:602
unsigned int get_num_locus(unsigned int trait)
Definition: ttquanti.h:424
unsigned int get_locus_start_pos(unsigned int locus)
Definition: ttquanti.h:457
double _equal_val_1
Definition: ttquanti.h:625
void set_eVarianceSD(unsigned int trait, double SD)
Definition: ttquanti.cc:2379
double get_genotype_with_dominance(const double a1, const double a2, const unsigned int locus, const unsigned int trait)
Definition: ttquanti.cc:2538
int get_allele_position(const unsigned int locus, const unsigned int trait)
Definition: ttquanti.cc:2467
void inherit_free(sex_t SEX, TTQuanti *ind, const TTQuanti *parent)
Definition: ttquanti.cc:2638
vector< double > get_env_var()
Definition: ttquanti.h:429
unsigned int _pleio_type
Definition: ttquanti.h:617
void mutate_no_pleio(TTQuanti *ind)
Definition: ttquanti.cc:2669
void deallocate_gsl_mutation_matrix_space()
Definition: ttquanti.cc:2136
TTQFreqExtractor * _freqExtractor
Definition: ttquanti.h:724
vector< double *(TProtoQuanti::*)(unsigned int) > _getMutationValuesVarPleio
Collection of pointers to mutation functions, which generate allele values in dependence of pleiotrop...
Definition: ttquanti.h:697
TTProtoWithMap.
Definition: ttrait_with_map.h:183
virtual void reset()
Definition: ttrait_with_map.cc:637
TTQFreqExtractor.
Definition: ttquanti.h:835
virtual ~TTQFreqExtractor()
Definition: ttquanti.h:844
TTQFreqExtractor(TProtoQuanti *T)
Definition: ttquanti.h:840
vector< string > _records
Definition: ttquanti.h:837
void resetTable()
Definition: ttquanti.cc:6319
virtual void FHread(string &filename)
Definition: ttquanti.h:846
virtual void FHwrite()
Definition: ttquanti.cc:6326
TTQOhtaStats.
Definition: ttquanti.h:854
TMatrix _pairwiseCombs
Definition: ttquanti.h:856
virtual void FHread(string &filename)
Definition: ttquanti.h:863
TTQOhtaStats(TProtoQuanti *T)
Definition: ttquanti.h:859
virtual ~TTQOhtaStats()
Definition: ttquanti.h:861
virtual void FHwrite()
Definition: ttquanti.cc:6465
TTQuantiFH.
Definition: ttquanti.h:809
void write_PLINK()
Definition: ttquanti.cc:5884
void print_PLINK_FAM(ofstream &FH, age_idx Ax, Patch *patch)
Definition: ttquanti.cc:6116
TTQuantiFH(TProtoQuanti *T)
Definition: ttquanti.h:815
virtual void FHread(string &filename)
Definition: ttquanti.cc:6191
void print_PLINK_PED(ofstream &FH, age_idx Ax, Patch *patch)
Definition: ttquanti.cc:6019
void write_TABLE()
Definition: ttquanti.cc:5720
virtual void FHwrite()
Definition: ttquanti.cc:5698
virtual ~TTQuantiFH()
Definition: ttquanti.h:817
string _output_option
Definition: ttquanti.h:811
bool _has_genetic_map
Definition: ttquanti.h:812
void setOutputOption(string opt)
Definition: ttquanti.cc:5676
void print(ofstream &FH, Patch *patch, sex_t SEX, age_idx Ax, unsigned int print_gene, bool print_genotype, bool print_additive_genotype)
Definition: ttquanti.cc:5815
TTQuantiSH.
Definition: ttquanti.h:733
double * _Vp
Definition: ttquanti.h:735
double getMeanGenotPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:790
void addQuanti(age_t AGE)
Definition: ttquanti.cc:4737
void addAvgPerPatch(age_t AGE)
Definition: ttquanti.cc:4822
double ** _pVa
Definition: ttquanti.h:736
double getEigenVectorEltPerPatch(unsigned int p, unsigned int v)
Definition: ttquanti.h:796
vector< double > getVaNoDominance(Patch *curPop, const age_idx AGE)
Definition: ttquanti.cc:5389
double * _Vb
Definition: ttquanti.h:735
unsigned int _num_locus
Definition: ttquanti.h:739
void addCovarPerPatch(age_t AGE)
Definition: ttquanti.cc:4945
double * _meanP
Definition: ttquanti.h:735
void addEigenVect1PerPatch(age_t AGE)
Definition: ttquanti.cc:5064
void addVarPerPatch(age_t AGE)
Definition: ttquanti.cc:4892
TTQuantiSH(TProtoQuanti *TP)
Definition: ttquanti.h:750
virtual ~TTQuantiSH()
Definition: ttquanti.h:758
double ** _pVp
Definition: ttquanti.h:736
double getVa(unsigned int i)
Definition: ttquanti.h:783
bool _eVar
Definition: ttquanti.h:737
gsl_matrix * _evec
Definition: ttquanti.h:741
double getEigenValuePerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:794
void addGenotPerPatch(age_t AGE)
Definition: ttquanti.cc:4857
unsigned int _table_set_age
Definition: ttquanti.h:746
void addEigenValuesPerPatch(age_t AGE)
Definition: ttquanti.cc:5030
double getVb(unsigned int i)
Definition: ttquanti.h:785
double getSkewPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.cc:5139
double ** _pmeanP
Definition: ttquanti.h:736
double * _meanG
Definition: ttquanti.h:735
double getMeanPhenotPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:791
void setStats(age_t AGE)
Definition: ttquanti.cc:5225
void setAdultStats()
Definition: ttquanti.h:778
double getVg(unsigned int i)
Definition: ttquanti.h:784
gsl_vector * _eval
Definition: ttquanti.h:742
double getMeanPhenot(unsigned int i)
Definition: ttquanti.h:782
double * _covar
Definition: ttquanti.h:735
DataTable< double > _phenoTable
Definition: ttquanti.h:745
double getCovarPerPatch(unsigned int p, unsigned int i)
Definition: ttquanti.h:795
double getQst(unsigned int i)
Definition: ttquanti.h:787
void addEigenPerPatch(age_t AGE)
Definition: ttquanti.cc:4986
vector< double > getVaWithDominance(Patch *curPop, const age_idx AGE)
computation of the additive genetic variance from the average excess of each allele exact under rando...
Definition: ttquanti.cc:5488
double getVp(unsigned int i)
Definition: ttquanti.h:786
gsl_matrix * _G
Definition: ttquanti.h:741
double getMeanGenot(unsigned int i)
Definition: ttquanti.h:781
double ** _peigvect
Definition: ttquanti.h:736
unsigned int _table_set_repl
Definition: ttquanti.h:746
double getVpPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:793
gsl_eigen_symmv_workspace * _ws
Definition: ttquanti.h:743
unsigned int _table_set_gen
Definition: ttquanti.h:746
double getCovar(unsigned int i)
Definition: ttquanti.h:788
vector< double > getSNPalleleFreqInPatch(Patch *patch, const age_idx AGE)
Definition: ttquanti.cc:5346
virtual bool setStatRecorders(std::string &token)
Definition: ttquanti.cc:4675
DataTable< double > _genoTable
Definition: ttquanti.h:745
unsigned int _num_trait
Definition: ttquanti.h:739
double * _Va
Definition: ttquanti.h:735
void resetPtrs()
Definition: ttquanti.cc:4544
double * _Vg
Definition: ttquanti.h:735
void setDataTables(age_t AGE)
Definition: ttquanti.cc:5154
void addQuantiPerPatch(age_t AGE)
Definition: ttquanti.cc:4806
double ** _peigval
Definition: ttquanti.h:736
unsigned int _patchNbr
Definition: ttquanti.h:739
double ** _pcovar
Definition: ttquanti.h:736
double ** _pmeanG
Definition: ttquanti.h:736
double getVaPerPatch(unsigned int i, unsigned int p)
Definition: ttquanti.h:792
void addSkewPerPatch(age_t AGE)
Definition: ttquanti.cc:5103
void setOffsprgStats()
Definition: ttquanti.h:779
virtual void init()
Definition: ttquanti.cc:4602
TTQuanti_continuous_full_pleio : universal pleiotropy.
Definition: ttquanti.h:170
TTQuanti_continuous_full_pleio(const TTQuanti_continuous_full_pleio &TT)
Definition: ttquanti.h:175
virtual void show_up()
Definition: ttquanti.cc:3451
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
Definition: ttquanti.cc:3362
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
Definition: ttquanti.cc:3345
virtual ~TTQuanti_continuous_full_pleio()
Definition: ttquanti.h:177
virtual TTQuanti_continuous_full_pleio * clone()
Definition: ttquanti.h:181
TTQuanti_continuous_full_pleio()
Definition: ttquanti.h:173
virtual double get_dominant_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3330
virtual double get_additive_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3313
virtual void init_sequence()
Definition: ttquanti.cc:3374
TTQuanti_continuous_no_pleio : multiple non-pleiotropic traits.
Definition: ttquanti.h:226
virtual double get_dominant_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3685
virtual ~TTQuanti_continuous_no_pleio()
Definition: ttquanti.h:233
TTQuanti_continuous_no_pleio()
Definition: ttquanti.h:229
virtual TTQuanti_continuous_no_pleio * clone()
Definition: ttquanti.h:237
TTQuanti_continuous_no_pleio(const TTQuanti_continuous_no_pleio &TT)
Definition: ttquanti.h:231
virtual void init_sequence()
Definition: ttquanti.cc:3754
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
Definition: ttquanti.cc:3744
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
Definition: ttquanti.cc:3728
virtual double get_additive_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3664
virtual void show_up()
Definition: ttquanti.cc:3808
TTQuanti_continuous_single : simple implementation for a single quantitative trait with continuous al...
Definition: ttquanti.h:253
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
virtual TTQuanti_continuous_single * clone()
Definition: ttquanti.h:264
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
virtual ~TTQuanti_continuous_single()
Definition: ttquanti.h:260
virtual double get_dominant_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3866
virtual void init_sequence()
virtual double get_additive_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3853
TTQuanti_continuous_single(const TTQuanti_continuous_single &TT)
Definition: ttquanti.h:258
TTQuanti_continuous_single()
Definition: ttquanti.h:256
TTQuanti_continuous_var_pleio : variable pleiotropy.
Definition: ttquanti.h:197
unsigned int get_sequence_block_size(unsigned int from, unsigned int to)
Definition: ttquanti.cc:3524
virtual double get_additive_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3490
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
Definition: ttquanti.cc:3547
TTQuanti_continuous_var_pleio()
Definition: ttquanti.h:200
virtual TTQuanti_continuous_var_pleio * clone()
Definition: ttquanti.h:210
virtual void show_up()
Definition: ttquanti.cc:3619
virtual double get_dominant_genotype(const unsigned int trait) const
Definition: ttquanti.cc:3507
TTQuanti_continuous_var_pleio(const TTQuanti_continuous_var_pleio &TT)
Definition: ttquanti.h:202
virtual ~TTQuanti_continuous_var_pleio()
Definition: ttquanti.h:204
virtual void init_sequence()
Definition: ttquanti.cc:3560
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
Definition: ttquanti.cc:3532
TTQuanti_continuous.
Definition: ttquanti.h:113
virtual bool retrieve_data(BinaryStorageBuffer *reader)
Definition: ttquanti.cc:3290
virtual void set_allele_value(unsigned int locus, unsigned int allele, double value)
Definition: ttquanti.cc:3230
TTQuanti_continuous(const TTQuanti &T)
Definition: ttquanti.h:122
TTQuanti_continuous()
Definition: ttquanti.h:117
virtual void set_allele(int locus, int allele, double value)
Definition: ttquanti.h:151
virtual void store_data(BinaryStorageBuffer *saver)
Definition: ttquanti.cc:3282
virtual bool operator!=(const TTrait &T)
Definition: ttquanti.cc:3204
virtual unsigned int get_allele(int loc, int all) const
Definition: ttquanti.cc:3214
virtual void ** get_sequence() const
Definition: ttquanti.h:136
virtual void set_allele_bit(unsigned int position, unsigned int allele, bool value)
Definition: ttquanti.h:157
virtual ~TTQuanti_continuous()
Definition: ttquanti.h:127
virtual void set_sequence(void **seq)
Definition: ttquanti.cc:3299
virtual void mutate_add(unsigned int position, unsigned int allele, double value)
Definition: ttquanti.h:152
virtual void init()
Definition: ttquanti.cc:3238
virtual void mutate_inplace(unsigned int position, unsigned int allele, double value)
Definition: ttquanti.h:154
virtual TTQuanti_continuous & operator=(const TTrait &T)
Definition: ttquanti.cc:3166
virtual double get_full_genotype(unsigned int trait)
Definition: ttquanti.cc:3275
virtual double get_allele_value(int loc, int all) const
Definition: ttquanti.cc:3222
virtual bool operator==(const TTrait &T)
Definition: ttquanti.cc:3187
virtual void reset()
Definition: ttquanti.cc:3258
double ** _sequence
Definition: ttquanti.h:162
virtual bool get_allele_bit(unsigned int position, unsigned int allele) const
Definition: ttquanti.h:156
TTQuanti_diallelic_full_pleio : pleiotropic di-allelic loci, max PD = 2.
Definition: ttquanti.h:355
virtual TTQuanti_diallelic_full_pleio * clone()
Definition: ttquanti.h:366
virtual double get_dominant_genotype(const unsigned int trait) const
Definition: ttquanti.cc:4243
virtual void show_up()
Definition: ttquanti.cc:4337
TTQuanti_diallelic_full_pleio()
Definition: ttquanti.h:358
virtual ~TTQuanti_diallelic_full_pleio()
Definition: ttquanti.h:362
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
Definition: ttquanti.cc:4264
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
Definition: ttquanti.cc:4280
virtual double get_additive_genotype(const unsigned int trait) const
Definition: ttquanti.cc:4225
virtual void init_sequence()
Definition: ttquanti.cc:4292
TTQuanti_diallelic_full_pleio(const TTQuanti_diallelic_full_pleio &TT)
Definition: ttquanti.h:360
TTQuanti_diallelic_no_pleio : single or multiple non-pleiotropic traits, di-allelic.
Definition: ttquanti.h:328
TTQuanti_diallelic_no_pleio()
Definition: ttquanti.h:331
virtual ~TTQuanti_diallelic_no_pleio()
Definition: ttquanti.h:335
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
Definition: ttquanti.cc:4111
virtual TTQuanti_diallelic_no_pleio * clone()
Definition: ttquanti.h:339
virtual double get_dominant_genotype(const unsigned int trait) const
Definition: ttquanti.cc:4090
virtual double get_additive_genotype(const unsigned int trait) const
Definition: ttquanti.cc:4062
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
Definition: ttquanti.cc:4127
virtual void init_sequence()
Definition: ttquanti.cc:4137
TTQuanti_diallelic_no_pleio(const TTQuanti_diallelic_no_pleio &TT)
Definition: ttquanti.h:333
virtual void show_up()
Definition: ttquanti.cc:4182
TTQuanti_diallelic_var_pleio : variable pleiotropic di-allelic loci, max PD = 2.
Definition: ttquanti.h:382
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
Definition: ttquanti.cc:4422
TTQuanti_diallelic_var_pleio(const TTQuanti_diallelic_var_pleio &TT)
Definition: ttquanti.h:387
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
Definition: ttquanti.cc:4439
TTQuanti_diallelic_var_pleio()
Definition: ttquanti.h:385
virtual TTQuanti_diallelic_var_pleio * clone()
Definition: ttquanti.h:395
virtual ~TTQuanti_diallelic_var_pleio()
Definition: ttquanti.h:389
virtual double get_additive_genotype(const unsigned int trait) const
Definition: ttquanti.cc:4375
virtual void show_up()
Definition: ttquanti.cc:4497
unsigned int get_sequence_block_size(unsigned int from, unsigned int to)
Definition: ttquanti.cc:4414
virtual double get_dominant_genotype(const unsigned int trait) const
Definition: ttquanti.cc:4393
virtual void init_sequence()
Definition: ttquanti.cc:4453
TTQuanti_diallelic.
Definition: ttquanti.h:280
virtual void ** get_sequence() const
Definition: ttquanti.h:298
virtual void mutate_inplace(unsigned int position, unsigned int allele, double value)
Definition: ttquanti.h:316
virtual unsigned int get_allele(int loc, int all) const
Definition: ttquanti.cc:3988
TTQuanti_diallelic()
Definition: ttquanti.h:284
TTQuanti_diallelic(const TTQuanti &T)
Definition: ttquanti.h:286
virtual TTQuanti_diallelic & operator=(const TTrait &T)
Definition: ttquanti.cc:3921
virtual void init()
Definition: ttquanti.cc:3883
virtual bool operator==(const TTrait &T)
Definition: ttquanti.cc:3942
virtual void mutate_add(unsigned int position, unsigned int allele, double value)
Definition: ttquanti.h:314
virtual void reset()
Definition: ttquanti.cc:3903
unsigned char ** _sequence
Definition: ttquanti.h:321
virtual double get_full_genotype(unsigned int trait)
Definition: ttquanti.cc:4025
virtual ~TTQuanti_diallelic()
Definition: ttquanti.h:288
virtual void store_data(BinaryStorageBuffer *saver)
Definition: ttquanti.cc:4032
virtual bool retrieve_data(BinaryStorageBuffer *reader)
Definition: ttquanti.cc:4040
virtual void set_allele(int locus, int allele, double value)
Definition: ttquanti.h:313
virtual double get_allele_value(int locus, int allele) const
Definition: ttquanti.cc:3997
virtual void set_allele_value(unsigned int locus, unsigned int allele, double value)
Definition: ttquanti.cc:4008
virtual bool get_allele_bit(unsigned int position, unsigned int allele) const
Definition: ttquanti.cc:3969
virtual void set_allele_bit(unsigned int position, unsigned int allele, bool value)
Definition: ttquanti.cc:3978
virtual void set_sequence(void **seq)
Definition: ttquanti.cc:4050
void inherit_free(const TTrait *mother, const TTrait *father)
virtual bool operator!=(const TTrait &T)
Definition: ttquanti.cc:3959
TTQuanti.
Definition: ttquanti.h:60
virtual void mutate_add(unsigned int position, unsigned int allele, double value)=0
TProtoQuanti * _myProto
Definition: ttquanti.h:106
virtual void copy_sequence_1locus(sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)=0
virtual ~TTQuanti()
Definition: ttquanti.h:71
virtual void set_value()
Definition: ttquanti.cc:3124
virtual void set_allele_bit(unsigned int position, unsigned int allele, bool value)=0
TProtoQuanti * get_proto()
Definition: ttquanti.h:97
virtual bool get_allele_bit(unsigned int position, unsigned int allele) const =0
virtual double get_dominant_genotype(const unsigned int trait) const =0
void set_phenotype(unsigned int trait, double value)
Definition: ttquanti.cc:3153
virtual double get_additive_genotype(const unsigned int trait) const =0
TTQuanti()
Definition: ttquanti.h:63
virtual void * getValue() const
Definition: ttquanti.h:78
virtual void * set_trait(void *value)
Definition: ttquanti.h:76
TTQuanti(const TTQuanti &T)
Definition: ttquanti.h:68
void set_proto(TProtoQuanti *proto)
Definition: ttquanti.h:96
virtual void mutate_inplace(unsigned int position, unsigned int allele, double value)=0
virtual void inherit(const TTrait *mother, const TTrait *father)
Definition: ttquanti.cc:3115
double get_phenotype(unsigned int trai)
Definition: ttquanti.cc:3144
double * _phenotypes
Definition: ttquanti.h:104
virtual double get_full_genotype(unsigned int trait)=0
void reset_phenotype_to_genotypic_value()
Definition: ttquanti.cc:3136
virtual void mutate()
Definition: ttquanti.cc:3108
double * _genotypic_values
Definition: ttquanti.h:105
virtual void copy_sequence_block(sex_t SEX, unsigned int strand, unsigned int from_pos, unsigned int to_pos, const TTQuanti *parent)=0
virtual trait_t get_type() const
Definition: ttquanti.h:73
Interface for all trait types, declares all basic trait operations.
Definition: ttrait.h:45
Template class for the trait's FileHandler.
Definition: filehandler.h:220
Template class for the trait's StatHandler.
Definition: stathandler.h:167
Non-template and faster implementation of std::bitset.
Definition: bitstring.h:56
Nemo2.
#define QUANT
Definition: types.h:71
std::string trait_t
Trait types.
Definition: types.h:62
sex_t
Sex types, males are always 0 and females 1!!
Definition: types.h:35
unsigned int age_t
Age class flags.
Definition: types.h:45
#define ADULTS
Adults age class flag (breeders).
Definition: types.h:53
#define OFFSPRG
Offspring age class flag.
Definition: types.h:49
age_idx
Array index of the age classes in the patch sizes and containers arrays.
Definition: types.h:40

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

Locations of visitors to this page
Catalogued on GSR