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

TTQuanti_continuous_full_pleio : universal pleiotropy. More...

#include <ttquanti.h>

+ Inheritance diagram for TTQuanti_continuous_full_pleio:
+ Collaboration diagram for TTQuanti_continuous_full_pleio:

Public Member Functions

 TTQuanti_continuous_full_pleio ()
 
 TTQuanti_continuous_full_pleio (const TTQuanti_continuous_full_pleio &TT)
 
virtual ~TTQuanti_continuous_full_pleio ()
 
virtual void init_sequence ()
 
virtual TTQuanti_continuous_full_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_full_pleio : universal pleiotropy.

Constructor & Destructor Documentation

◆ TTQuanti_continuous_full_pleio() [1/2]

TTQuanti_continuous_full_pleio::TTQuanti_continuous_full_pleio ( )
inline
169 : TTQuanti_continuous() {}
TTQuanti_continuous()
Definition: ttquanti.h:114

Referenced by clone().

◆ TTQuanti_continuous_full_pleio() [2/2]

TTQuanti_continuous_full_pleio::TTQuanti_continuous_full_pleio ( const TTQuanti_continuous_full_pleio TT)
inline
171 : TTQuanti_continuous(TT) {}

◆ ~TTQuanti_continuous_full_pleio()

virtual TTQuanti_continuous_full_pleio::~TTQuanti_continuous_full_pleio ( )
inlinevirtual
173 {}

Member Function Documentation

◆ clone()

virtual TTQuanti_continuous_full_pleio* TTQuanti_continuous_full_pleio::clone ( )
inlinevirtual

Implements TTrait.

177 {return new TTQuanti_continuous_full_pleio(*this);}
TTQuanti_continuous_full_pleio()
Definition: ttquanti.h:169

References TTQuanti_continuous_full_pleio().

◆ copy_sequence_1locus()

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

Implements TTQuanti.

3311 {
3312  const double *orig = (const double*)parent->get_sequence()[strand];
3313  double *seq = _sequence[SEX];
3314  unsigned int from_pos = at * _myProto->get_num_traits();
3315 
3316  memcpy(&seq[from_pos], &orig[from_pos], _myProto->get_locus_byte_size());
3317 }
size_t get_locus_byte_size()
Definition: ttquanti.h:423
unsigned int get_num_traits()
Definition: ttquanti.h:417
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_byte_size(), TProtoQuanti::get_num_traits(), and TTrait::get_sequence().

◆ copy_sequence_block()

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

Implements TTQuanti.

3294 {
3295  const double *orig = (const double*)parent->get_sequence()[strand];
3296  double *seq = _sequence[SEX];
3297  unsigned int from_pos = from_locus * _myProto->get_num_traits();
3298 
3299  assert(to_locus >= from_locus);
3300 
3301  size_t block_size = (to_locus - from_locus) * _myProto->get_locus_byte_size();
3302 
3303  memcpy(&seq[from_pos], &orig[from_pos], block_size);
3304 
3305 }

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

◆ get_additive_genotype()

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

Implements TTQuanti.

3260 {
3261  unsigned int loc = trait; // initial position in the sequence, in the first locus
3262  double genotype = 0;
3263 
3264 
3265  for(unsigned int j = 0; j < _myProto->get_num_locus(); ++j) {
3266  // alleles are ordered locus-trait: seq[0]:l0t0,l0t1,l0t2, ..., l1t0,l1t1,l1t2, ...
3267  genotype += _sequence[0][loc] + _sequence[1][loc];
3268  loc += _myProto->get_num_traits(); // stride size on the genome, between loci
3269  }
3270 
3271  return genotype;
3272 // return _myProto->set_genotype_value_additive_full_pleio((const double**)get_sequence(), trait);
3273 }
unsigned int get_num_locus()
Definition: ttquanti.h:418

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

◆ get_dominant_genotype()

double TTQuanti_continuous_full_pleio::get_dominant_genotype ( const unsigned int  trait) const
virtual

Implements TTQuanti.

3278 {
3279  unsigned int loc = trait;
3280  double genotype = 0;
3281 
3282  for(unsigned int j = 0; j < _myProto->get_num_locus(); ++j) {
3283  genotype += _myProto->get_genotype_with_dominance( _sequence[0][loc] , _sequence[1][loc], j, trait);
3284  loc += _myProto->get_num_traits(); // stride size on the genome, between loci
3285  }
3286 
3287  return genotype;
3288 }
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_num_locus(), and TProtoQuanti::get_num_traits().

◆ init_sequence()

void TTQuanti_continuous_full_pleio::init_sequence ( )
inlinevirtual

Implements TTrait.

3322 {
3323 // cout<<" TTQuanti_continuous_full_pleio::init_sequence\n";
3324 
3325  //options:
3326  //0: no variation, init value = (trait value)/(2*_num_locus)
3327  //1: init value = (trait value)/(2*_num_locus) + 1 mutation/locus
3328  //2: init value = (trait value)/(2*_num_locus) + random deviate ~N(0, _init_variance/2*_num_locus)
3329 
3330  //Note: the initial values may have been set individually by LCE_quanti
3331  // it wouldn't make sense then to store the init values in the prototype
3332  // because init values are set patch-specific by LCE_quanti
3333  // Problem: without LCE_quanti, the original init values must not change from one replicate
3334  // to the other, so we use a local variable
3335 
3336  double my_init[_myProto->get_num_traits()];
3337 
3338  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++)
3339  my_init[j] = _myProto->get_init_value(j) / (2*_myProto->get_num_locus());
3340 
3341  unsigned int pos = 0;
3342 
3343  for(unsigned int j = 0; j < _myProto->get_num_traits(); ++j) {
3344  pos = j;
3345  //set the allele values from the init value, mutations may be added later
3346  for(unsigned int i = 0; i < _myProto->get_num_locus(); ++i) {
3347  _sequence[0][pos] = my_init[j];
3348  _sequence[1][pos] = my_init[j];
3349  pos += _myProto->get_num_traits();
3350  }
3351  }
3352 
3353  //add random effects to allele values if needed
3354  if(_myProto->get_doInitMutation() == 1) {
3355 
3356  double *mut;
3357 
3358  for(unsigned int i = 0; i < _myProto->get_num_locus(); i++) {
3359 
3360  pos = i * _myProto->get_num_traits();
3361 
3362  // change only one of the two allele at the locus, randomly
3363  mut = (_myProto->*_myProto->_getMutationValues)(i);
3364 
3365  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++) {
3366  _sequence[RAND::RandBool()][pos] += mut[j];
3367  pos++;
3368  }
3369  } // end for_num_locus
3370  } // end if_doInit == 1
3371 
3372  // random deviate to initial additive effect
3373  if(_myProto->get_doInitMutation() == 2) {
3374 
3375  double initSD[_myProto->get_num_traits()];
3376 
3377  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++) {
3378  initSD[j] = sqrt(_myProto->get_init_variance(j) / (2*_myProto->get_num_locus()) );
3379  }
3380 
3381  for(unsigned int i = 0; i < _myProto->get_num_locus(); i++) {
3382 
3383  pos = i * _myProto->get_num_traits();
3384 
3385  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++) {
3386  _sequence[0][pos] += RAND::Gaussian(initSD[j]);
3387  _sequence[1][pos] += RAND::Gaussian(initSD[j]);
3388  pos++;
3389  }
3390  } // end for_num_locus
3391  } // end if_doInit == 2
3392 
3393 // cout<<" TTQuanti_continuous_full_pleio::init_sequence::end\n";
3394 }
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
double *(TProtoQuanti::* _getMutationValues)(unsigned int)
Pointer to mutation allele value function, which depends on allele model and number of traits affecte...
Definition: ttquanti.h:663
unsigned int get_doInitMutation()
Definition: ttquanti.h:436

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

◆ show_up()

void TTQuanti_continuous_full_pleio::show_up ( )
virtual

Implements TTrait.

3399 {
3400  message("\
3401  Trait's type: QUANTI (continuous, full pleiotropy)\n\
3402  traits: %i\n\
3403  loci: %i\n\
3404  seq length: %i\n",_myProto->get_num_traits(),_myProto->get_num_locus()
3405  ,_myProto->get_seq_length());
3406 
3407  for(unsigned int i = 0; i < _myProto->get_num_traits(); i++)
3408  message("phenotype %i: %f\n",i+1,_phenotypes[i]);
3409 
3410  message("_sequence: \n0:");
3411  unsigned int loc;
3412  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
3413  message("\nt%i:",i+1);
3414  for(unsigned int j = 0; (j < _myProto->get_num_locus()); ++j) {
3415  loc = j * _myProto->get_num_traits() + i;
3416  message("%.3f,",_sequence[0][loc]);
3417  }
3418  }
3419  message("\n1:");
3420  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
3421  message("\nt%i:",i+1);
3422  for(unsigned int j = 0; (j < _myProto->get_num_locus()); ++j) {
3423  loc = j * _myProto->get_num_traits() + i;
3424  message("%.3f,",_sequence[1][loc]);
3425  }
3426  }
3427  message("\n");
3428 
3429 }
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_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