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

TTQuanti_continuous_var_pleio : variable pleiotropy. More...

#include <ttquanti.h>

+ Inheritance diagram for TTQuanti_continuous_var_pleio:
+ Collaboration diagram for TTQuanti_continuous_var_pleio:

Public Member Functions

 TTQuanti_continuous_var_pleio ()
 
 TTQuanti_continuous_var_pleio (const TTQuanti_continuous_var_pleio &TT)
 
virtual ~TTQuanti_continuous_var_pleio ()
 
unsigned int get_sequence_block_size (unsigned int from, unsigned int to)
 
virtual void init_sequence ()
 
virtual TTQuanti_continuous_var_pleioclone ()
 
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 strand, unsigned int from_locus, unsigned int to_locus, const TTQuanti *parent)
 
virtual void copy_sequence_1locus (sex_t SEX, unsigned int strand, unsigned int at, const TTQuanti *parent)
 
- 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 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)
 
- 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
 
TProtoQuanti_myProto
 

Detailed Description

TTQuanti_continuous_var_pleio : variable pleiotropy.

Constructor & Destructor Documentation

◆ TTQuanti_continuous_var_pleio() [1/2]

TTQuanti_continuous_var_pleio::TTQuanti_continuous_var_pleio ( )
inline
196 : TTQuanti_continuous() {}
TTQuanti_continuous()
Definition: ttquanti.h:114

Referenced by clone().

◆ TTQuanti_continuous_var_pleio() [2/2]

TTQuanti_continuous_var_pleio::TTQuanti_continuous_var_pleio ( const TTQuanti_continuous_var_pleio TT)
inline
198 : TTQuanti_continuous(TT) {}

◆ ~TTQuanti_continuous_var_pleio()

virtual TTQuanti_continuous_var_pleio::~TTQuanti_continuous_var_pleio ( )
inlinevirtual
200 {}

Member Function Documentation

◆ clone()

virtual TTQuanti_continuous_var_pleio* TTQuanti_continuous_var_pleio::clone ( )
inlinevirtual

Implements TTrait.

206 {return new TTQuanti_continuous_var_pleio(*this);}
TTQuanti_continuous_var_pleio()
Definition: ttquanti.h:196

References TTQuanti_continuous_var_pleio().

◆ copy_sequence_1locus()

void TTQuanti_continuous_var_pleio::copy_sequence_1locus ( sex_t  SEX,
unsigned int  strand,
unsigned int  at,
const TTQuanti parent 
)
inlinevirtual

Implements TTQuanti.

3503 {
3504  const double *orig = (const double*)parent->get_sequence()[strand];
3505  double *seq = _sequence[SEX];
3506  unsigned int from_pos = _myProto->get_locus_start_pos(at);
3507  size_t size = _myProto->get_locus_PD(at) * _myProto->get_size_locus_type();
3508 
3509  memcpy(&seq[from_pos], &orig[from_pos], size);
3510 }
unsigned int get_locus_PD(unsigned int locus)
Definition: ttquanti.h:444
size_t get_size_locus_type()
Definition: ttquanti.h:422
unsigned int get_locus_start_pos(unsigned int locus)
Definition: ttquanti.h:445
double ** _sequence
Definition: ttquanti.h:158
TProtoQuanti * _myProto
Definition: ttquanti.h:103
virtual void ** get_sequence() const =0
sequence accessor.

References TTQuanti::_myProto, TTQuanti_continuous::_sequence, TProtoQuanti::get_locus_PD(), TProtoQuanti::get_locus_start_pos(), TTrait::get_sequence(), and TProtoQuanti::get_size_locus_type().

◆ copy_sequence_block()

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

Implements TTQuanti.

3488 {
3489  const double *orig = (const double*)parent->get_sequence()[strand];
3490  double *seq = _sequence[SEX];
3491  unsigned int from_pos = _myProto->get_locus_start_pos(from_locus);
3492 
3493  size_t block_size = get_sequence_block_size(from_locus, to_locus);
3494 
3495  memcpy(&seq[from_pos], &orig[from_pos], block_size);
3496 
3497 }
unsigned int get_sequence_block_size(unsigned int from, unsigned int to)
Definition: ttquanti.cc:3471

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

◆ get_additive_genotype()

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

Implements TTQuanti.

3438 {
3439  double genotype = 0;
3440  unsigned int L = _myProto->get_num_locus(trait), pos;
3441 
3442  for(unsigned int j = 0; j < L; ++j) { // each trait can be affected by a different number of loci
3443  // alleles are ordered locus-trait: seq[0]:l0t0,l0t1,l0t2, ..., l1t0,l1t1,l1t2, ...
3444  pos = _myProto->get_locus_seq_pos(j, trait);
3445  genotype += _sequence[0][pos] + _sequence[1][pos];
3446  }
3447 
3448  return genotype;
3449 }
unsigned int get_locus_seq_pos(unsigned int loc, unsigned int trait)
Definition: ttquanti.h:440
unsigned int get_num_locus()
Definition: ttquanti.h:418

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

◆ get_dominant_genotype()

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

Implements TTQuanti.

3455 {
3456  double genotype = 0;
3457  unsigned int L = _myProto->get_num_locus(trait), pos, locID;
3458 
3459  for(unsigned int j = 0; j < L; ++j) {
3460  pos = _myProto->get_locus_seq_pos(j, trait);
3461  locID = _myProto->get_locus_ID(j, trait);
3462  genotype += _myProto->get_genotype_with_dominance(_sequence[0][pos] ,
3463  _sequence[1][pos], locID, trait);
3464  }
3465 
3466  return genotype;
3467 }
unsigned int get_locus_ID(unsigned int locus, unsigned int trait)
Definition: ttquanti.h:442
double get_genotype_with_dominance(const double a1, const double a2, const unsigned int locus, const unsigned int trait)
Definition: ttquanti.cc:2504

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

◆ get_sequence_block_size()

unsigned int TTQuanti_continuous_var_pleio::get_sequence_block_size ( unsigned int  from,
unsigned int  to 
)
inline
3472 {
3473  unsigned int size = 0;
3474  assert(to > from);
3475 
3476  for (unsigned int i = from; i < to; ++i) {
3477 
3478  size += _myProto->get_locus_PD(i); // count number of allelic values to copy
3479 
3480  }
3481  return size * _myProto->get_size_locus_type(); // multiply with number of bytes used to store allele values
3482 }

References TTQuanti::_myProto, TProtoQuanti::get_locus_PD(), and TProtoQuanti::get_size_locus_type().

Referenced by copy_sequence_block().

◆ init_sequence()

void TTQuanti_continuous_var_pleio::init_sequence ( )
inlinevirtual

Implements TTrait.

3515 {
3516 // cout << "TTQuanti_var_pleio::init_sequence\n";
3517 
3518  //options:
3519  //0: no variation, init value = (trait value)/(2*_num_locus)
3520  //1: init value = (trait value)/(2*_num_locus) + 1 mutation/locus
3521  //2: init value = (trait value)/(2*_num_locus) + random deviate ~N(0, _init_variance/2*_num_locus)
3522 
3523  //Note: the initial values may have been set individually by LCE_quanti
3524  // it wouldn't make sense then to store the init values in the prototype
3525  // because init values are set patch-specific by LCE_quanti
3526  // Problem: without LCE_quanti, the original init values must not change from one replicate
3527  // to the other, so we use a local variable
3528 
3529 
3530  //set the allele values from the trait value and the pleiotropic degree of each locus
3531  double myinit, initSD;
3532 
3533  for(unsigned int t = 0; t < _myProto->get_num_traits(); ++t){ // process trait-wise
3534 
3535  myinit = _myProto->get_init_value(t) / (2*_myProto->get_num_locus(t)); //(2 * ttable[i].size()); //divide by num loci affecting that trait
3536 
3537  initSD = sqrt(_myProto->get_init_variance(t)/(2*_myProto->get_num_locus(t)));
3538 
3539  for(unsigned int i = 0; i < _myProto->get_num_locus(t); ++i){
3540  _sequence[0][ _myProto->get_locus_seq_pos(i, t) ] = myinit;
3541  _sequence[1][ _myProto->get_locus_seq_pos(i, t) ] = myinit; //ttable[i][j]
3542  }
3543  //add a random deviate to initial additive effect
3544  if(_myProto->get_doInitMutation() == 2) {
3545  for(unsigned int i = 0; i < _myProto->get_num_locus(t); i++) {
3546  _sequence[0][ _myProto->get_locus_seq_pos(i, t) ] += RAND::Gaussian(initSD);
3547  _sequence[1][ _myProto->get_locus_seq_pos(i, t) ] += RAND::Gaussian(initSD);
3548  } // end for_num_locus
3549  } // end if_doInit == 2
3550 
3551  }
3552  //add random effects to allele values
3553  if(_myProto->get_doInitMutation() == 1) {
3554 
3555  double *mut;
3556 
3557  for(unsigned int i = 0; i < _myProto->get_num_locus(); i++) {
3558 
3559  mut = _myProto->getMutationEffectsVarPleio(i); // returns as many values as the locus' PD
3560 
3561  for(unsigned int j = 0; j < _myProto->get_locus_PD(i); j++){
3562  //mutate only one of the alleles, randomly
3563  _sequence[RAND::RandBool()][ _myProto->get_locus_start_pos(i) + j ] += mut[j];
3564  }
3565  }
3566  } // end init model 1
3567 
3568 }
static double Gaussian(double sigma)
Definition: Uniform.h:261
static bool RandBool()
Returns a random boolean.
Definition: Uniform.h:162
double get_init_value(unsigned int i)
Definition: ttquanti.h:434
double get_init_variance(unsigned int i)
Definition: ttquanti.h:435
unsigned int get_num_traits()
Definition: ttquanti.h:417
unsigned int get_doInitMutation()
Definition: ttquanti.h:436
double * getMutationEffectsVarPleio(unsigned int loc)
Definition: ttquanti.h:519

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

◆ show_up()

void TTQuanti_continuous_var_pleio::show_up ( )
virtual

Implements TTrait.

3574 {
3575  message("\
3576  Trait's type: QUANTI (continuous, variable pleiotropy)\n\
3577  traits: %i\n\
3578  loci: %i\n\
3579  seq length: %i\n",_myProto->get_num_traits(),_myProto->get_num_locus()
3580  ,_myProto->get_seq_length());
3581 
3582  for(unsigned int i = 0; i < _myProto->get_num_traits(); i++)
3583  message("phenotype %i: %f\n",i+1,_phenotypes[i]);
3584 
3585  message("genotype:");
3586 
3587  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
3588 
3589  message("\ntrait %i (%i loci):\nloci: ",i+1, _myProto->get_num_locus(i));
3590 
3591  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
3592  message("%i, ", _myProto->get_locus_ID(j, i));
3593  }
3594 
3595  message("\n[0]: ");
3596 
3597  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
3598  message("%.3f,",_sequence[0][ _myProto->get_locus_seq_pos(j, i) ]);
3599  }
3600 
3601  message("\n[1]: ");
3602 
3603  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
3604  message("%.3f,",_sequence[1][ _myProto->get_locus_seq_pos(j, i) ]);
3605  }
3606  }
3607 
3608  message("\n");
3609 
3610 }
unsigned int get_seq_length()
Definition: ttquanti.h:421
double * _phenotypes
Definition: ttquanti.h:102
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