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

TTQuanti_continuous_no_pleio : multiple non-pleiotropic traits with epistasis. More...

#include <ttquanti_epistasis.h>

+ Inheritance diagram for TTQuanti_continuous_no_pleio_epistasis:
+ Collaboration diagram for TTQuanti_continuous_no_pleio_epistasis:

Public Member Functions

 TTQuanti_continuous_no_pleio_epistasis ()
 
 TTQuanti_continuous_no_pleio_epistasis (const TTQuanti_continuous_no_pleio_epistasis &TT)
 
virtual ~TTQuanti_continuous_no_pleio_epistasis ()
 
virtual void init_sequence ()
 
virtual TTQuanti_continuous_no_pleio_epistasisclone ()
 
virtual void show_up ()
 
virtual double get_additive_genotype (const unsigned int trait) const
 
virtual double get_dominant_genotype (const unsigned int trait) const
 
virtual void copy_sequence_block (sex_t SEX, unsigned int chromosome, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
 
virtual void copy_sequence_1locus (sex_t SEX, unsigned int chromosome, unsigned int at, const TTQuanti *parent)
 
double get_epistatic_genotype (const vector< double > &genotypes) const
 
- Public Member Functions inherited from TTQuanti_continuous
 TTQuanti_continuous ()
 
 TTQuanti_continuous (const TTQuanti &T)
 
virtual ~TTQuanti_continuous ()
 
virtual void reset ()
 
virtual void init ()
 
virtual void set_sequence (void **seq)
 
virtual void ** get_sequence () const
 
virtual unsigned int get_allele (int loc, int all) const
 
virtual double get_allele_value (int loc, int all) const
 
virtual void set_allele_value (unsigned int locus, unsigned int allele, double value)
 
virtual TTQuanti_continuousoperator= (const TTrait &T)
 
virtual bool operator== (const TTrait &T)
 
virtual bool operator!= (const TTrait &T)
 
virtual void store_data (BinaryStorageBuffer *saver)
 
virtual bool retrieve_data (BinaryStorageBuffer *reader)
 
virtual double get_full_genotype (unsigned int trait)
 
virtual void set_allele (int locus, int allele, double value)
 
virtual void mutate_add (unsigned int position, unsigned int allele, double value)
 
virtual void mutate_inplace (unsigned int position, unsigned int allele, double value)
 
virtual bool get_allele_bit (unsigned int position, unsigned int allele) const
 
virtual void set_allele_bit (unsigned int position, unsigned int allele, bool value)
 
- Public Member Functions inherited from TTQuanti
 TTQuanti ()
 
 TTQuanti (const TTQuanti &T)
 
virtual ~TTQuanti ()
 
virtual trait_t get_type () const
 
virtual void mutate ()
 
virtual void inherit (const TTrait *mother, const TTrait *father)
 
virtual void * set_trait (void *value)
 
virtual void set_value ()
 
virtual void * getValue () const
 
void set_proto (TProtoQuanti *proto)
 
TProtoQuantiget_proto ()
 
double get_phenotype (unsigned int trai)
 
void set_phenotype (unsigned int trait, double value)
 
void reset_phenotype_to_genotypic_value ()
 
- Public Member Functions inherited from TTrait
virtual ~TTrait ()
 
- Public Member Functions inherited from StorableComponent
virtual ~StorableComponent ()
 

Additional Inherited Members

- Protected Attributes inherited from TTQuanti_continuous
double ** _sequence
 
- Protected Attributes inherited from TTQuanti
double * _phenotypes
 
double * _genotypic_values
 
TProtoQuanti_myProto
 

Detailed Description

TTQuanti_continuous_no_pleio : multiple non-pleiotropic traits with epistasis.

Constructor & Destructor Documentation

◆ TTQuanti_continuous_no_pleio_epistasis() [1/2]

TTQuanti_continuous_no_pleio_epistasis::TTQuanti_continuous_no_pleio_epistasis ( )
inline
103 : TTQuanti_continuous() {}
TTQuanti_continuous()
Definition: ttquanti.h:118

Referenced by clone().

◆ TTQuanti_continuous_no_pleio_epistasis() [2/2]

TTQuanti_continuous_no_pleio_epistasis::TTQuanti_continuous_no_pleio_epistasis ( const TTQuanti_continuous_no_pleio_epistasis TT)
inline
105  :
106  TTQuanti_continuous(TT) {}

◆ ~TTQuanti_continuous_no_pleio_epistasis()

virtual TTQuanti_continuous_no_pleio_epistasis::~TTQuanti_continuous_no_pleio_epistasis ( )
inlinevirtual
108 {}

Member Function Documentation

◆ clone()

virtual TTQuanti_continuous_no_pleio_epistasis* TTQuanti_continuous_no_pleio_epistasis::clone ( )
inlinevirtual

Implements TTrait.

112 {return new TTQuanti_continuous_no_pleio_epistasis(*this);}
TTQuanti_continuous_no_pleio_epistasis()
Definition: ttquanti_epistasis.h:103

References TTQuanti_continuous_no_pleio_epistasis().

◆ copy_sequence_1locus()

void TTQuanti_continuous_no_pleio_epistasis::copy_sequence_1locus ( sex_t  SEX,
unsigned int  chromosome,
unsigned int  at,
const TTQuanti parent 
)
inlinevirtual

Implements TTQuanti.

519 {
520  const double *orig = (const double*)parent->get_sequence()[chromosome];
521 
522  _sequence[SEX][at] = orig[at];
523 }
double ** _sequence
Definition: ttquanti.h:163
virtual void ** get_sequence() const =0
sequence accessor.

References TTQuanti_continuous::_sequence, and TTrait::get_sequence().

◆ copy_sequence_block()

void TTQuanti_continuous_no_pleio_epistasis::copy_sequence_block ( sex_t  SEX,
unsigned int  chromosome,
unsigned int  from_locus,
unsigned int  to_locus,
const TTQuanti parent 
)
inlinevirtual

Implements TTQuanti.

503 {
504  assert(to_locus >= from_locus);
505 
506  const double *orig = (const double*)parent->get_sequence()[chromosome];
507  double *seq = _sequence[SEX];
508 
509  size_t block_size = (to_locus - from_locus) * _myProto->get_locus_byte_size();
510 
511  memcpy(&seq[from_locus], &orig[from_locus], block_size);
512 
513 }
size_t get_locus_byte_size()
Definition: ttquanti.h:429
TProtoQuanti * _myProto
Definition: ttquanti.h:107

References TTQuanti::_myProto, TTQuanti_continuous::_sequence, TProtoQuanti::get_locus_byte_size(), and TTrait::get_sequence().

◆ get_additive_genotype()

double TTQuanti_continuous_no_pleio_epistasis::get_additive_genotype ( const unsigned int  trait) const
inlinevirtual

Implements TTQuanti.

437 {
438  double genotype = 0;
439  unsigned int L = _myProto->get_num_locus(trait); //number of loci affecting this trait
440  unsigned int pos = _myProto->get_locus_seq_pos(0, trait); //starting position, all loci contiguous on the map
441  vector<double> genotypes(L,0);
442 
443  for(unsigned int j = 0; j < L; ++j, ++pos) {
444  genotypes[j] = _sequence[0][pos] + _sequence[1][pos];
445  genotype += genotypes[j];
446  }
447 
448  return genotype + get_epistatic_genotype(genotypes);
449 }
unsigned int get_locus_seq_pos(unsigned int loc, unsigned int trait)
Definition: ttquanti.h:453
unsigned int get_num_locus()
Definition: ttquanti.h:424
double get_epistatic_genotype(const vector< double > &genotypes) const
Definition: ttquanti_epistasis.cc:473

References TTQuanti::_myProto, TTQuanti_continuous::_sequence, get_epistatic_genotype(), TProtoQuanti::get_locus_seq_pos(), and TProtoQuanti::get_num_locus().

◆ get_dominant_genotype()

double TTQuanti_continuous_no_pleio_epistasis::get_dominant_genotype ( const unsigned int  trait) const
inlinevirtual

Implements TTQuanti.

454 {
455  double genotype = 0;
456  unsigned int L = _myProto->get_num_locus(trait); //number of loci affecting this trait
457  unsigned int pos = _myProto->get_locus_seq_pos(0, trait); //starting position, all loci contiguous on the map
458  unsigned int locID;
459  vector<double> genotypes(L,0);
460 
461  for(unsigned int j = 0; j < L; ++j, ++pos) {
462  locID = _myProto->get_locus_ID(j, trait);
463  genotypes[j] = _myProto->get_genotype_with_dominance(_sequence[0][pos] , _sequence[1][pos],
464  locID, trait);
465  genotype += genotypes[j];
466  }
467 
468  return genotype + get_epistatic_genotype(genotypes);
469 }
unsigned int get_locus_ID(unsigned int locus, unsigned int trait)
Definition: ttquanti.h:455
double get_genotype_with_dominance(const double a1, const double a2, const unsigned int locus, const unsigned int trait)
Definition: ttquanti.cc:2539

References TTQuanti::_myProto, TTQuanti_continuous::_sequence, get_epistatic_genotype(), TProtoQuanti::get_genotype_with_dominance(), TProtoQuanti::get_locus_ID(), TProtoQuanti::get_locus_seq_pos(), and TProtoQuanti::get_num_locus().

◆ get_epistatic_genotype()

double TTQuanti_continuous_no_pleio_epistasis::get_epistatic_genotype ( const vector< double > &  genotypes) const
inline
474 {
475  double genotype = 0;
476 
477  const TMatrix& epi_loc_idx = _myProto->get_epi_coef_index();
478  const TMatrix& epi_coef = _myProto->get_epi_coefs();
479  unsigned int pos0, pos1;
480 
481  // now add the epistatic effects
482  for(unsigned int i = 0; i < epi_coef.length(); ++i) {
483 
484  // Note _epistatic_coefs_indices has dimensions [_nb_epi_coefs; 2]
485  pos0 = epi_loc_idx.get(i,0);
486  pos1 = epi_loc_idx.get(i,1);
487 
488 // genotype += epi_coef.get(0,i)
489 // * (_sequence[0][pos0] + _sequence[1][pos0])
490 // * (_sequence[0][pos1] + _sequence[1][pos1]);
491 
492  genotype += epi_coef.get(0,i) * double(genotypes[pos0]) * double(genotypes[pos1]);
493 
494  }
495 
496  return genotype;
497 }
A class to handle matrix in params, coerces matrix into a vector of same total size.
Definition: tmatrix.h:50
double get(unsigned int i, unsigned int j) const
Accessor to element at row i and column j.
Definition: tmatrix.h:193
unsigned int length() const
Returns the number of elements in the matrix.
Definition: tmatrix.h:218
const TMatrix & get_epi_coefs() const
Definition: ttquanti.h:568
const TMatrix & get_epi_coef_index() const
Definition: ttquanti.h:569

References TTQuanti::_myProto, TMatrix::get(), TProtoQuanti::get_epi_coef_index(), TProtoQuanti::get_epi_coefs(), and TMatrix::length().

Referenced by get_additive_genotype(), and get_dominant_genotype().

◆ init_sequence()

void TTQuanti_continuous_no_pleio_epistasis::init_sequence ( )
inlinevirtual

Implements TTrait.

528 {
529  // cout << "\nStart of TTQuanti_var_pleio::init_sequence!\t\n";
530 
531  //options:
532  //0: no variation, init value = (trait value)/(2*_num_locus)
533  //1: init value = (trait value)/(2*_num_locus) + 1 mutation/locus
534  //2: init value = (trait value)/(2*_num_locus) + random deviate ~N(0, _init_variance/2*_num_locus)
535 
536  //Note: the initial values may have been set individually by LCE_quanti
537  // it wouldn't make sense then to store the init values in the prototype
538  // because init values are set patch-specific by LCE_quanti
539  // Problem: without LCE_quanti, the original init values must not change from one replicate
540  // to the other, so we use a local variable
541 
542 
543  double myinit, initSD;
544  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i){ // process trait-wise
545 
546  myinit = _myProto->get_init_value(i) / (2 * _myProto->get_num_locus(i)); //divide by num loci affecting that trait
547  initSD = sqrt(_myProto->get_init_variance(i)/(2*_myProto->get_num_locus(i)));
548 
549  for(unsigned int j = 0, pos = _myProto->get_locus_seq_pos(0, i); j < _myProto->get_num_locus(i); ++j, ++pos){
550  //cout << "\t_init_value after after: " << myinit << endl;
551  _sequence[0][ pos ] = myinit;
552  _sequence[1][ pos ] = myinit;
553  }
554 
555  if(_myProto->get_doInitMutation() == 2) {
556  // add a random deviate from a Normal distribution
557  for(unsigned int j = 0, pos = _myProto->get_locus_seq_pos(0, i); j < _myProto->get_num_locus(i); ++j, ++pos){
558  _sequence[0][ pos ] += RAND::Gaussian(initSD);
559  _sequence[1][ pos ] += RAND::Gaussian(initSD);
560  }
561  }
562  }
563 
564  //add random mutations to allele values
565  if(_myProto->get_doInitMutation() == 1) {
566 
567  for(unsigned int i = 0; i < _myProto->get_num_locus(); i++) {
568 
570 
571  }
572  }
573 
574 }
static double Gaussian(double sigma)
Definition: Uniform.h:264
static bool RandBool()
Returns a random boolean.
Definition: Uniform.h:165
double * getMutationEffects(unsigned int loc)
Definition: ttquanti.h:533
double get_init_value(unsigned int i)
Definition: ttquanti.h:447
double get_init_variance(unsigned int i)
Definition: ttquanti.h:448
unsigned int get_num_traits()
Definition: ttquanti.h:423
unsigned int get_doInitMutation()
Definition: ttquanti.h:449

References TTQuanti::_myProto, TTQuanti_continuous::_sequence, RAND::Gaussian(), TProtoQuanti::get_doInitMutation(), TProtoQuanti::get_init_value(), TProtoQuanti::get_init_variance(), TProtoQuanti::get_locus_seq_pos(), TProtoQuanti::get_num_locus(), TProtoQuanti::get_num_traits(), TProtoQuanti::getMutationEffects(), and RAND::RandBool().

◆ show_up()

void TTQuanti_continuous_no_pleio_epistasis::show_up ( )
virtual

Implements TTrait.

579 {
580  message("\
581  Trait's type: QUANTI (continuous, no pleiotropy with epistasis)\n\
582  traits: %i\n\
583  loci: %i\n\
584  seq length: %i\n",_myProto->get_num_traits(),_myProto->get_num_locus()
586 
587  for(unsigned int i = 0; i < _myProto->get_num_traits(); i++)
588  message("phenotype %i: %f\n",i+1,_phenotypes[i]);
589 
590  message("genotype:");
591 
592  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
593 
594  message("\ntrait %i (%i loci):\nloci: ",i+1, _myProto->get_num_locus(i));
595 
596  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
597  message("%i, ", _myProto->get_locus_ID(j, i));
598  }
599 
600  message("\n[0]: ");
601 
602  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
603  message("%.3f,",_sequence[0][ _myProto->get_locus_seq_pos(j, i) ]);
604  }
605 
606  message("\n[1]: ");
607 
608  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
609  message("%.3f,",_sequence[1][ _myProto->get_locus_seq_pos(j, i) ]);
610  }
611  }
612 
613  message("\n");
614 
615 }
unsigned int get_seq_length()
Definition: ttquanti.h:427
double * _phenotypes
Definition: ttquanti.h:105
void message(const char *message,...)
Definition: output.cc:40

References TTQuanti::_myProto, TTQuanti::_phenotypes, TTQuanti_continuous::_sequence, TProtoQuanti::get_locus_ID(), TProtoQuanti::get_locus_seq_pos(), TProtoQuanti::get_num_locus(), TProtoQuanti::get_num_traits(), TProtoQuanti::get_seq_length(), and message().


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