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

TTQuanti_continuous_full_pleio : universal pleiotropy. More...

#include <ttquanti_epistasis.h>

+ Inheritance diagram for TTQuanti_continuous_full_pleio_epistasis:
+ Collaboration diagram for TTQuanti_continuous_full_pleio_epistasis:

Public Member Functions

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

Constructor & Destructor Documentation

◆ TTQuanti_continuous_full_pleio_epistasis() [1/2]

TTQuanti_continuous_full_pleio_epistasis::TTQuanti_continuous_full_pleio_epistasis ( )
inline
TTQuanti_continuous()
Definition: ttquanti.h:118

Referenced by clone().

◆ TTQuanti_continuous_full_pleio_epistasis() [2/2]

TTQuanti_continuous_full_pleio_epistasis::TTQuanti_continuous_full_pleio_epistasis ( const TTQuanti_continuous_full_pleio_epistasis TT)
inline
46  :
47  TTQuanti_continuous(TT) {}

◆ ~TTQuanti_continuous_full_pleio_epistasis()

virtual TTQuanti_continuous_full_pleio_epistasis::~TTQuanti_continuous_full_pleio_epistasis ( )
inlinevirtual
49 {}

Member Function Documentation

◆ clone()

virtual TTQuanti_continuous_full_pleio_epistasis* TTQuanti_continuous_full_pleio_epistasis::clone ( )
inlinevirtual

Implements TTrait.

TTQuanti_continuous_full_pleio_epistasis()
Definition: ttquanti_epistasis.h:44

References TTQuanti_continuous_full_pleio_epistasis().

◆ copy_sequence_1locus()

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

Implements TTQuanti.

130 {
131  const double *orig = (const double*)parent->get_sequence()[chromosome];
132  double *seq = _sequence[SEX];
133  unsigned int from_pos = at * _myProto->get_num_traits();
134 
135  memcpy(&seq[from_pos], &orig[from_pos], _myProto->get_locus_byte_size());
136 }
size_t get_locus_byte_size()
Definition: ttquanti.h:429
unsigned int get_num_traits()
Definition: ttquanti.h:423
double ** _sequence
Definition: ttquanti.h:163
TProtoQuanti * _myProto
Definition: ttquanti.h:107
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_epistasis::copy_sequence_block ( sex_t  SEX,
unsigned int  chromosome,
unsigned int  from_locus,
unsigned int  to_locus,
const TTQuanti parent 
)
inlinevirtual

Implements TTQuanti.

113 {
114  const double *orig = (const double*)parent->get_sequence()[chromosome];
115  double *seq = _sequence[SEX];
116  unsigned int from_pos = from_locus * _myProto->get_num_traits();
117 
118  assert(to_locus >= from_locus);
119 
120  size_t block_size = (to_locus - from_locus) * _myProto->get_locus_byte_size();
121 
122  memcpy(&seq[from_pos], &orig[from_pos], block_size);
123 
124 }

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_epistasis::get_additive_genotype ( const unsigned int  trait) const
inlinevirtual

Implements TTQuanti.

46 {
47  unsigned int loc = trait; // initial position in the sequence, in the first locus
48  unsigned int L = _myProto->get_num_locus(); // same for all traits
49  vector<double> genotypes(L,0);
50  double genotype = 0;
51 
52 
53  for(unsigned int j = 0; j < L; ++j) {
54  // alleles are ordered locus-trait: seq[0]:l0t0,l0t1,l0t2, ..., l1t0,l1t1,l1t2, ...
55  genotypes[j] = _sequence[0][loc] + _sequence[1][loc];
56  loc += _myProto->get_num_traits(); // stride size on the genome, between loci
57  genotype += genotypes[j];
58  }
59 
60  return genotype + get_epistatic_genotype(genotypes);
61 }
unsigned int get_num_locus()
Definition: ttquanti.h:424
double get_epistatic_genotype(const vector< double > &genotypes) const
Definition: ttquanti_epistasis.cc:83

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

◆ get_dominant_genotype()

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

Implements TTQuanti.

66 {
67  unsigned int loc = trait; // position of first allele in the sequence at first locus
68  unsigned int L = _myProto->get_num_locus(); // same for all traits
69  vector<double> genotypes(L,0);
70  double genotype = 0;
71 
72  for(unsigned int j = 0; j < L; ++j) {
73  genotypes[j] = _myProto->get_genotype_with_dominance( _sequence[0][loc] , _sequence[1][loc], j, trait);
74  loc += _myProto->get_num_traits(); // stride size on the genome, between loci
75  genotype += genotypes[j];
76  }
77 
78  return genotype + get_epistatic_genotype(genotypes);
79 }
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_num_locus(), and TProtoQuanti::get_num_traits().

◆ get_epistatic_genotype()

double TTQuanti_continuous_full_pleio_epistasis::get_epistatic_genotype ( const vector< double > &  genotypes) const
inline
84 {
85  double genotype = 0;
86 
87  const TMatrix& epi_loc_idx = _myProto->get_epi_coef_index();
88  const TMatrix& epi_coef = _myProto->get_epi_coefs();
89  unsigned int pos0, pos1;
90 
91  // now add the epistatic effects
92  for(unsigned int i = 0; i < epi_coef.length(); ++i) {
93 
94  // Note _epistatic_coefs_indices has dimensions [_nb_epi_coefs; 2]
95  pos0 = epi_loc_idx.get(i,0);
96  pos1 = epi_loc_idx.get(i,1);
97 
98 // genotype += epi_coef.get(0,i)
99 // * (_sequence[0][pos0] + _sequence[1][pos0])
100 // * (_sequence[0][pos1] + _sequence[1][pos1]);
101 
102  genotype += epi_coef.get(0,i) * genotypes[pos0] * genotypes[pos1];
103 
104  }
105 
106  return genotype;
107 }
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_full_pleio_epistasis::init_sequence ( )
inlinevirtual

Implements TTrait.

141 {
142 // cout<<" TTQuanti_continuous_full_pleio_epistasis::init_sequence\n";
143 
144  //options:
145  //0: no variation, init value = (trait value)/(2*_num_locus)
146  //1: init value = (trait value)/(2*_num_locus) + 1 mutation/locus
147  //2: only di-allelic
148 
149  //Note: the initial values may have been set individually by LCE_quanti
150  // it wouldn't make sense then to store the init values in the prototype
151  // because init values are set patch-specific by LCE_quanti
152  // Problem: without LCE_quanti, the original init values must not change from one replicate
153  // to the other, so we use a local variable
154 
155  double my_init[_myProto->get_num_traits()];
156 
157  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++)
158  my_init[j] = _myProto->get_init_value(j) / (2*_myProto->get_num_locus());
159 
160  unsigned int pos = 0;
161 
162  for(unsigned int j = 0; j < _myProto->get_num_traits(); ++j) {
163  pos = j;
164  //set the allele values from the init value, mutations may be added later
165  for(unsigned int i = 0; i < _myProto->get_num_locus(); ++i) {
166  _sequence[0][pos] = my_init[j];
167  _sequence[1][pos] = my_init[j];
168  pos += _myProto->get_num_traits();
169  }
170  }
171 
172  //add random effects to allele values if needed
173  if(_myProto->get_doInitMutation() == 1) {
174 
175  double *mut;
176 
177  for(unsigned int i = 0; i < _myProto->get_num_locus(); i++) {
178 
179  pos = i * _myProto->get_num_traits();
180 
181  // change only one of the two allele at the locus, randomly
182  mut = _myProto->getMutationEffects(i);
183 
184  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++) {
185  _sequence[RAND::RandBool()][pos] += mut[j];
186  pos++;
187  }
188  } // end for_num_locus
189  } // end if_doInit == 1
190 
191  // random deviate to initial additive effect
192  if(_myProto->get_doInitMutation() == 2) {
193 
194  double initSD[_myProto->get_num_traits()];
195 
196  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++) {
197  initSD[j] = sqrt(_myProto->get_init_variance(j) / (2*_myProto->get_num_locus()) );
198  }
199 
200  for(unsigned int i = 0; i < _myProto->get_num_locus(); i++) {
201 
202  pos = i * _myProto->get_num_traits();
203 
204  for(unsigned int j = 0; j < _myProto->get_num_traits(); j++) {
205  _sequence[0][pos] += RAND::Gaussian(initSD[j]);
206  _sequence[1][pos] += RAND::Gaussian(initSD[j]);
207  pos++;
208  }
209  } // end for_num_locus
210  } // end if_doInit == 2
211 
212 // cout<<" TTQuanti_continuous_full_pleio_epistasis::init_sequence::end\n";
213 }
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_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_num_locus(), TProtoQuanti::get_num_traits(), TProtoQuanti::getMutationEffects(), and RAND::RandBool().

◆ show_up()

void TTQuanti_continuous_full_pleio_epistasis::show_up ( )
virtual

Implements TTrait.

218 {
219  message("\
220  Trait's type: QUANTI (continuous, full pleiotropy with epistasis)\n\
221  traits: %i\n\
222  loci: %i\n\
223  seq length: %i\n",_myProto->get_num_traits(),_myProto->get_num_locus()
225 
226  for(unsigned int i = 0; i < _myProto->get_num_traits(); i++)
227  message("phenotype %i: %f\n",i+1,_phenotypes[i]);
228 
229  message("_sequence: \n0:");
230  unsigned int loc;
231  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
232  message("\nt%i:",i+1);
233  for(unsigned int j = 0; (j < _myProto->get_num_locus()); ++j) {
234  loc = j * _myProto->get_num_traits() + i;
235  message("%.3f,",_sequence[0][loc]);
236  }
237  }
238  message("\n1:");
239  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
240  message("\nt%i:",i+1);
241  for(unsigned int j = 0; (j < _myProto->get_num_locus()); ++j) {
242  loc = j * _myProto->get_num_traits() + i;
243  message("%.3f,",_sequence[1][loc]);
244  }
245  }
246  message("\n");
247 
248 }
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_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