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

TTQuanti_diallelic_no_pleio_epistasis : single or multiple non-pleiotropic traits, di-allelic. More...

#include <ttquanti_epistasis.h>

+ Inheritance diagram for TTQuanti_diallelic_no_pleio_epistasis:
+ Collaboration diagram for TTQuanti_diallelic_no_pleio_epistasis:

Public Member Functions

 TTQuanti_diallelic_no_pleio_epistasis ()
 
 TTQuanti_diallelic_no_pleio_epistasis (const TTQuanti_diallelic_no_pleio_epistasis &TT)
 
virtual ~TTQuanti_diallelic_no_pleio_epistasis ()
 
virtual void init_sequence ()
 
virtual TTQuanti_diallelic_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_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_no_pleio_epistasis : single or multiple non-pleiotropic traits, di-allelic.

Constructor & Destructor Documentation

◆ TTQuanti_diallelic_no_pleio_epistasis() [1/2]

TTQuanti_diallelic_no_pleio_epistasis::TTQuanti_diallelic_no_pleio_epistasis ( )
inline
132 : TTQuanti_diallelic() {}
TTQuanti_diallelic()
Definition: ttquanti.h:285

Referenced by clone().

◆ TTQuanti_diallelic_no_pleio_epistasis() [2/2]

TTQuanti_diallelic_no_pleio_epistasis::TTQuanti_diallelic_no_pleio_epistasis ( const TTQuanti_diallelic_no_pleio_epistasis TT)
inline
134  :
135  TTQuanti_diallelic(TT) {}

◆ ~TTQuanti_diallelic_no_pleio_epistasis()

virtual TTQuanti_diallelic_no_pleio_epistasis::~TTQuanti_diallelic_no_pleio_epistasis ( )
inlinevirtual
137 {}

Member Function Documentation

◆ clone()

virtual TTQuanti_diallelic_no_pleio_epistasis* TTQuanti_diallelic_no_pleio_epistasis::clone ( )
inlinevirtual

Implements TTrait.

141 {return new TTQuanti_diallelic_no_pleio_epistasis(*this);}
TTQuanti_diallelic_no_pleio_epistasis()
Definition: ttquanti_epistasis.h:132

References TTQuanti_diallelic_no_pleio_epistasis().

◆ copy_sequence_1locus()

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

Implements TTQuanti.

705 {
706  const unsigned char *orig = (const unsigned char*)parent->get_sequence()[chromosome];
707 
708  _sequence[SEX][at] = orig[at];
709 }
unsigned char ** _sequence
Definition: ttquanti.h:322
virtual void ** get_sequence() const =0
sequence accessor.

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

◆ copy_sequence_block()

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

689 {
690  const unsigned char *orig = (const unsigned char*)(parent->get_sequence()[chromosome]);
691  unsigned char *seq = _sequence[SEX];
692 
693  assert(to_locus >= from_locus);
694 
695  size_t block_size = (to_locus - from_locus) * _myProto->get_locus_byte_size();
696 
697  memcpy(&seq[from_locus], &orig[from_locus], block_size);
698 
699 }
size_t get_locus_byte_size()
Definition: ttquanti.h:429
TProtoQuanti * _myProto
Definition: ttquanti.h:107

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

◆ get_additive_genotype()

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

Implements TTQuanti.

622 {
623  double genotype = 0;
624  unsigned int L = _myProto->get_num_locus(trait); //number of loci affecting this trait
625  unsigned int pos = _myProto->get_locus_seq_pos(0, trait); //starting position, all loci contiguous on the map
626  vector<double> genotypes(L,0);
627 
628  for(unsigned int j = 0; j < L; ++j, ++pos) {
629  genotypes[j] = _myProto->get_seq_diallele_value( pos, _sequence[0][pos] )
630  + _myProto->get_seq_diallele_value( pos, _sequence[1][pos] );
631  genotype += genotypes[j];
632  }
633 
634  return genotype + get_epistatic_genotype(genotypes);
635 }
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:664

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

◆ get_dominant_genotype()

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

Implements TTQuanti.

641 {
642  double genotype = 0;
643  unsigned int L = _myProto->get_num_locus(trait); //number of loci affecting this trait
644  unsigned int pos = _myProto->get_locus_seq_pos(0, trait); //starting position, all loci contiguous on the map
645  unsigned int locID;
646  vector<double> genotypes(L,0);
647 
648  for(unsigned int j = 0; j < L; ++j, ++pos) {
649  locID = _myProto->get_locus_ID(j, trait);
650  genotypes[j] = _myProto->get_genotype_with_dominance(
651  _myProto->get_seq_diallele_value( pos, _sequence[0][pos] ) ,
652  _myProto->get_seq_diallele_value( pos, _sequence[1][pos] ),
653  locID,
654  trait);
655  genotype += genotypes[j];
656  }
657 
658 
659  return genotype + get_epistatic_genotype(genotypes);
660 }
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_diallelic::_sequence, get_epistatic_genotype(), TProtoQuanti::get_genotype_with_dominance(), TProtoQuanti::get_locus_ID(), TProtoQuanti::get_locus_seq_pos(), TProtoQuanti::get_num_locus(), and TProtoQuanti::get_seq_diallele_value().

◆ get_epistatic_genotype()

double TTQuanti_diallelic_no_pleio_epistasis::get_epistatic_genotype ( const vector< double > &  genotypes) const
inline
665 {
666  double genotype = 0;
667 
668  const TMatrix& epi_loc_idx = _myProto->get_epi_coef_index();
669  const TMatrix& epi_coef = _myProto->get_epi_coefs();
670  unsigned int pos0, pos1;
671 
672  // now add the epistatic effects
673  for(unsigned int i = 0; i < epi_coef.length(); ++i) {
674 
675  // Note _epistatic_coefs_indices has dimensions [_nb_epi_coefs; 2]
676  pos0 = epi_loc_idx.get(i,0);
677  pos1 = epi_loc_idx.get(i,1);
678 
679  genotype += epi_coef.get(0,i) * double(genotypes[pos0]) * double(genotypes[pos1]);
680  }
681 
682  return genotype;
683 }
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_no_pleio_epistasis::init_sequence ( )
inlinevirtual

Implements TTrait.

714 {
715  //options:
716  //0: no variation, all allele are set to allele 0
717  //1: max variation, randomly set alleles at each locus
718  //2: no variation and initialize with opposite effect alleles at alternating loci.
719 
720 
721  if (_myProto->get_doInitMutation() == 0){
722 
723  memset(_sequence[0], 0, _myProto->get_seq_length());
724  memset(_sequence[1], 0, _myProto->get_seq_length());
725 
726  } else if (_myProto->get_doInitMutation() == 1) {
727 
728  for(unsigned int i = 0; i < _myProto->get_num_locus(); ++i) {
729 
730  // randomly set the allelic values
731  _sequence[0][ i ] = RAND::RandBool();
732  _sequence[1][ i ] = RAND::RandBool();
733 
734  }
735 
736  } else {
737  //option 2; polarized di-allelic values; even locus set to +a, odd locus set to -a
738 
739  for(unsigned int i = 0; i < _myProto->get_num_locus(); i++) { // process trait-wise
740  // no need to set to 0 since it was initialized that way
741  if (i % 2 == 0){
742  _sequence[0][ i ] = 0; //_myProto->get_allele_value(i, 0); //set with the positive allele value
743  _sequence[1][ i ] = 0; //_myProto->get_allele_value(i, 0);
744  }
745  else { // odd locus
746  _sequence[1][ i ] = 1; // = _myProto->get_allele_value(i, 1); //set with the alternate/negative allele value
747  _sequence[0][ i ] = 1; // = _myProto->get_allele_value(i, 1);
748  }
749  }
750  }
751 
752 }
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_num_locus(), TProtoQuanti::get_seq_length(), and RAND::RandBool().

◆ show_up()

void TTQuanti_diallelic_no_pleio_epistasis::show_up ( )
virtual

Implements TTrait.

757 {
758  message("\
759  Trait's type: QUANTI (diallelic, no pleiotropy with epistasis)\n\
760  traits: %i\n\
761  loci: %i\n\
762  seq length: %i\n",_myProto->get_num_traits(),_myProto->get_num_locus()
764 
765  for(unsigned int i = 0; i < _myProto->get_num_traits(); i++)
766  message("phenotype %i: %f\n",i+1,_phenotypes[i]);
767 
768  message("genotype:");
769 
770  for(unsigned int i = 0; i < _myProto->get_num_traits(); ++i) {
771 
772  message("\ntrait %i (%i loci):\nloci: ",i+1, _myProto->get_num_locus(i));
773 
774  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
775  message("%i, ", _myProto->get_locus_ID(j, i));
776  }
777 
778  message("\n[0]: ");
779 
780  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
781  message("%.3f,",(double)_sequence[0][ _myProto->get_locus_seq_pos(j, i) ]);
782  }
783 
784  message("\n[1]: ");
785 
786  for(unsigned int j = 0; j < _myProto->get_num_locus(i); ++j) {
787  message("%.3f,",(double)_sequence[1][ _myProto->get_locus_seq_pos(j, i) ]);
788  }
789  }
790 
791  message("\n");
792 
793 }
unsigned int get_num_traits()
Definition: ttquanti.h:423
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_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