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

TTQuanti_diallelic_full_pleio : pleiotropic di-allelic loci, max PD = 2. More...

#include <ttquanti_epistasis.h>

+ Inheritance diagram for TTQuanti_diallelic_full_pleio_epistasis:
+ Collaboration diagram for TTQuanti_diallelic_full_pleio_epistasis:

Public Member Functions

 TTQuanti_diallelic_full_pleio_epistasis ()
 
 TTQuanti_diallelic_full_pleio_epistasis (const TTQuanti_diallelic_full_pleio_epistasis &TT)
 
virtual ~TTQuanti_diallelic_full_pleio_epistasis ()
 
virtual void init_sequence ()
 
virtual TTQuanti_diallelic_full_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_diallelic
 TTQuanti_diallelic ()
 
 TTQuanti_diallelic (const TTQuanti &T)
 
virtual ~TTQuanti_diallelic ()
 
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)
 
void inherit_free (const TTrait *mother, const TTrait *father)
 
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 locus, int allele) const
 
virtual void set_allele_value (unsigned int locus, unsigned int allele, double value)
 
virtual TTQuanti_diallelicoperator= (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)
 
- 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_diallelic
unsigned char ** _sequence
 
- Protected Attributes inherited from TTQuanti
double * _phenotypes
 
double * _genotypic_values
 
TProtoQuanti_myProto
 

Detailed Description

TTQuanti_diallelic_full_pleio : pleiotropic di-allelic loci, max PD = 2.

Constructor & Destructor Documentation

◆ TTQuanti_diallelic_full_pleio_epistasis() [1/2]

TTQuanti_diallelic_full_pleio_epistasis::TTQuanti_diallelic_full_pleio_epistasis ( )
inline
74 : TTQuanti_diallelic() {}
TTQuanti_diallelic()
Definition: ttquanti.h:285

Referenced by clone().

◆ TTQuanti_diallelic_full_pleio_epistasis() [2/2]

TTQuanti_diallelic_full_pleio_epistasis::TTQuanti_diallelic_full_pleio_epistasis ( const TTQuanti_diallelic_full_pleio_epistasis TT)
inline
76  :
77  TTQuanti_diallelic(TT) {}

◆ ~TTQuanti_diallelic_full_pleio_epistasis()

virtual TTQuanti_diallelic_full_pleio_epistasis::~TTQuanti_diallelic_full_pleio_epistasis ( )
inlinevirtual
79 {}

Member Function Documentation

◆ clone()

virtual TTQuanti_diallelic_full_pleio_epistasis* TTQuanti_diallelic_full_pleio_epistasis::clone ( )
inlinevirtual

Implements TTrait.

83 {return new TTQuanti_diallelic_full_pleio_epistasis(*this);}
TTQuanti_diallelic_full_pleio_epistasis()
Definition: ttquanti_epistasis.h:74

References TTQuanti_diallelic_full_pleio_epistasis().

◆ copy_sequence_1locus()

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

Implements TTQuanti.

342 {
343  const unsigned char *orig = (const unsigned char*)parent->get_sequence()[chromosome];
344  unsigned char *seq = _sequence[SEX];
345  unsigned int from_pos = at * _myProto->get_num_traits();
346 
347  memcpy(&seq[from_pos], &orig[from_pos], _myProto->get_locus_byte_size());
348 }
size_t get_locus_byte_size()
Definition: ttquanti.h:429
unsigned int get_num_traits()
Definition: ttquanti.h:423
unsigned char ** _sequence
Definition: ttquanti.h:322
TProtoQuanti * _myProto
Definition: ttquanti.h:107
virtual void ** get_sequence() const =0
sequence accessor.

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

◆ copy_sequence_block()

void TTQuanti_diallelic_full_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.

326 {
327  assert(to_locus >= from_locus);
328 
329  const unsigned char* orig = (const unsigned char*)parent->get_sequence()[chromosome];
330  unsigned char* seq = _sequence[SEX];
331  unsigned int from_pos = from_locus * _myProto->get_num_traits();
332 
333  size_t block_size = (to_locus - from_locus) * _myProto->get_locus_byte_size();
334 
335  memcpy(&seq[from_pos], &orig[from_pos], block_size);
336 }

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

◆ get_additive_genotype()

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

Implements TTQuanti.

255 {
256  double genotype = 0;
257  unsigned int L = _myProto->get_num_locus(); //number of loci
258  unsigned int pos = _myProto->get_locus_seq_pos(0, trait); //starting position, all loci contiguous on the map
259  vector<double> genotypes(L,0);
260 
261  for(unsigned int j = 0; j < L; ++j) {
262  genotypes[j] = _myProto->get_seq_diallele_value( pos, _sequence[0][pos] ) + _myProto->get_seq_diallele_value( pos, _sequence[1][pos] );
263  pos += _myProto->get_num_traits(); // stride size on the genome (all loci affect each trait)
264  genotype += genotypes[j];
265  }
266 
267  return genotype + get_epistatic_genotype(genotypes);
268 }
double get_seq_diallele_value(unsigned int position, unsigned int allele)
Definition: ttquanti.h:438
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:296

References TTQuanti::_myProto, TTQuanti_diallelic::_sequence, get_epistatic_genotype(), TProtoQuanti::get_locus_seq_pos(), TProtoQuanti::get_num_locus(), TProtoQuanti::get_num_traits(), and TProtoQuanti::get_seq_diallele_value().

◆ get_dominant_genotype()

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

Implements TTQuanti.

274 {
275  double genotype = 0;
276  unsigned int L = _myProto->get_num_locus(); //number of loci
277  unsigned int pos = _myProto->get_locus_seq_pos(0, trait); //starting position
278  vector<double> genotypes(L,0);
279 
280 
281  for(unsigned int j = 0; j < L; ++j) {
282  genotypes[j] = _myProto->get_genotype_with_dominance(
283  _myProto->get_seq_diallele_value( pos, _sequence[0][pos] ),
284  _myProto->get_seq_diallele_value( pos, _sequence[1][pos] ),
285  j,
286  trait);
287  pos += _myProto->get_num_traits(); // stride size on the genome (all loci affect each trait)
288  genotype += genotypes[j];
289  }
290 
291  return genotype + get_epistatic_genotype(genotypes);
292 }
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_diallelic::_sequence, get_epistatic_genotype(), TProtoQuanti::get_genotype_with_dominance(), TProtoQuanti::get_locus_seq_pos(), TProtoQuanti::get_num_locus(), TProtoQuanti::get_num_traits(), and TProtoQuanti::get_seq_diallele_value().

◆ get_epistatic_genotype()

double TTQuanti_diallelic_full_pleio_epistasis::get_epistatic_genotype ( const vector< double > &  genotypes) const
inline
297 {
298  double genotype = 0;
299 
300  const TMatrix& epi_loc_idx = _myProto->get_epi_coef_index();
301  const TMatrix& epi_coef = _myProto->get_epi_coefs();
302  unsigned int pos0, pos1;
303 
304  // now add the epistatic effects
305  for(unsigned int i = 0; i < epi_coef.length() ; ++i) {
306 
307  // Note _epistatic_coefs_indices has dimensions [_nb_epi_coefs; 2]
308  pos0 = epi_loc_idx.get(i,0);
309  pos1 = epi_loc_idx.get(i,1);
310 
311 // genotype += epi_coef.get(0,i)
312 // * (_sequence[0][pos0] + _sequence[1][pos0])
313 // * (_sequence[0][pos1] + _sequence[1][pos1]);
314 
315  genotype += epi_coef.get(0,i) * double(genotypes[pos0]) * double(genotypes[pos1]);
316 
317  }
318 
319  return genotype;
320 }
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_diallelic_full_pleio_epistasis::init_sequence ( )
inlinevirtual

Implements TTrait.

353 {
354  //options:
355  //0: no variation, init value = one of each allele, mean value is average genotype
356  //1: max variation, randomly set alleles at each locus
357  //2: no variation and initialize with opposite effect alleles at alternating loci.
358 
359 
360  if (_myProto->get_doInitMutation() == 0){
361 
362  memset(_sequence[0], 0, _myProto->get_seq_length());
363  memset(_sequence[1], 0, _myProto->get_seq_length());
364 
365  } else if (_myProto->get_doInitMutation() == 1) {
366 
367  for(unsigned int i = 0; i < _myProto->get_seq_length(); ++i) {
368 
369  // randomly set the allelic values
370  _sequence[0][ i ] = RAND::RandBool();
371  _sequence[1][ i ] = RAND::RandBool();
372 
373  }
374 
375  } else {
376  //option 2; polarized di-allelic values; even locus set to +a, odd locus set to -a
377 
378  for(unsigned int t = 0; t < _myProto->get_num_traits(); ++t) {
379  for(unsigned int i = 0; i < _myProto->get_num_locus(t); i++) { // process trait-wise
380  // even locus
381  if (i % 2 == 0){
382  _sequence[0][ _myProto->get_locus_seq_pos(i, t) ] = 0; //_myProto->get_allele_value(i, 0); //set with the positive allele value
383  _sequence[1][ _myProto->get_locus_seq_pos(i, t) ] = 0; //_myProto->get_allele_value(i, 0);
384  }
385  else { // odd locus
386  _sequence[1][ _myProto->get_locus_seq_pos(i, t) ] = 1; // = _myProto->get_allele_value(i, 1); //set with the alternate/negative allele value
387  _sequence[0][ _myProto->get_locus_seq_pos(i, t) ] = 1; // = _myProto->get_allele_value(i, 1);
388  }
389  }
390 
391  }
392  }
393 }
static bool RandBool()
Returns a random boolean.
Definition: Uniform.h:165
unsigned int get_seq_length()
Definition: ttquanti.h:427
unsigned int get_doInitMutation()
Definition: ttquanti.h:449

References TTQuanti::_myProto, TTQuanti_diallelic::_sequence, TProtoQuanti::get_doInitMutation(), TProtoQuanti::get_locus_seq_pos(), TProtoQuanti::get_num_locus(), TProtoQuanti::get_num_traits(), TProtoQuanti::get_seq_length(), and RAND::RandBool().

◆ show_up()

void TTQuanti_diallelic_full_pleio_epistasis::show_up ( )
virtual

Implements TTrait.

398 {
399  message("\
400  Trait's type: QUANTI (diallelic, full pleiotropy with epistasis)\n\
401  traits: %i\n\
402  loci: %i\n\
403  seq length: %i\n",_myProto->get_num_traits(),_myProto->get_num_locus()
405 
406  for(unsigned int i = 0; i < _myProto->get_num_traits(); i++)
407  message("phenotype %i: %f\n",i+1,_phenotypes[i]);
408 
409  message("_sequence:");
410 
411  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
412 
413  message("\nt%i (%i loci):\n[0]: ",i+1, _myProto->get_num_locus(i));
414 
415  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
416  message("%.3f,",(double)_sequence[0][ _myProto->get_locus_seq_pos(j, i) ]);
417  }
418 
419  message("\n[1]: ");
420 
421  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
422  message("%.3f,",(double)_sequence[1][ _myProto->get_locus_seq_pos(j, i) ]);
423  }
424  }
425 
426  message("\n");
427 
428 }
double * _phenotypes
Definition: ttquanti.h:105
void message(const char *message,...)
Definition: output.cc:40

References TTQuanti::_myProto, TTQuanti::_phenotypes, TTQuanti_diallelic::_sequence, 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