Loading [MathJax]/extensions/ams.js
ALMaSS  1.2 (after EcoStack, March 2024)
The Animal, Landscape and Man Simulation System
All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GeneticMaterial Class Reference

Class for the genetic material optionally carried by animals in ALMaSS. More...

#include <GeneticMaterial.h>

Public Member Functions

 GeneticMaterial ()
 
void ReadFrequencies ()
 
void SetAllele (int pos, uint32 value, int Chromosome)
 
uint32 GetAllele (int pos, int Chromosome)
 
void PrintChromosome (char *C, int Chromosome)
 
int HomozygosityCount ()
 
int HeterozygosityCount ()
 
void Recombine (GeneticMaterial *Gen21, GeneticMaterial *Gene2)
 
void Initiation (AlleleFreq *Al)
 
float ScoreReproduction ()
 
float ScoreHQThreshold ()
 
void SetGeneticFlag ()
 
void SetDirectFlag ()
 
void UnsetGeneticFlag ()
 
void UnsetDirectFlag ()
 
uint32 GetGeneticFlag ()
 
uint32 GetDirectFlag ()
 
void Mutation_1 ()
 
void Mutation_1ab ()
 
void Mutation_2 ()
 
void Mutation_3 ()
 
void Mutation_4 ()
 

Protected Attributes

uint32 Chromosome [6]
 

Detailed Description

Class for the genetic material optionally carried by animals in ALMaSS.

Constructor & Destructor Documentation

◆ GeneticMaterial()

GeneticMaterial::GeneticMaterial ( )
352  {
353  // ensure zeros in all loci
354  for ( int i = 0; i < 6; i++ ) Chromosome[ i ] = 0;
355 }

References Chromosome.

Member Function Documentation

◆ GetAllele()

uint32 GeneticMaterial::GetAllele ( int  pos,
int  Chromosome 
)
230  {
231  uint32 value;
232  // Get the right chromosome
233  // if Chromo==0 then 0-2, else 3-5
234  Chromo *= 3; // 0 or 3
235  if ( locus < 16 ) {
236  // Shift it so the locus is in the last two bits
237  // Does it twice because 32 bis coding for 16 loci
238  value = Chromosome[ Chromo ] >> locus;
239  value = ( value >> locus ) & 0x03;
240  } else {
241  Chromo++; // 1 or 4
242  locus -= 16; // Now 0 to 16
243  if ( locus >= 8 ) {
244  Chromo++; // 2 or 5
245  locus -= 8;
246  }
247  value = Chromosome[ Chromo ] >> ( locus * 4 );
248  value = value & 0x0f;
249  }
250  return value;
251 }

References Chromosome.

Referenced by GetDirectFlag(), GetGeneticFlag(), HeterozygosityCount(), HomozygosityCount(), Mutation_2(), Mutation_3(), Mutation_4(), PrintChromosome(), Recombine(), Vole_Base::SupplyAllele(), and Vole_Base::SupplyMyAllele().

◆ GetDirectFlag()

uint32 GeneticMaterial::GetDirectFlag ( )
188  {
189  return GetAllele(0,1);
190 }

References GetAllele().

Referenced by Vole_Base::GetDirectFlag(), and Vole_JuvenileMale::st_BecomeSubAdult().

◆ GetGeneticFlag()

uint32 GeneticMaterial::GetGeneticFlag ( )
184  {
185  return GetAllele(0,0);
186 }

References GetAllele().

Referenced by Vole_Base::GetGeneticFlag(), and Vole_JuvenileMale::st_BecomeSubAdult().

◆ HeterozygosityCount()

int GeneticMaterial::HeterozygosityCount ( )
325  {
326  int heterozyg = 0;
327  for ( int i = 0; i < 32; i++ ) {
328  if ( GetAllele( i, 0 ) != GetAllele( i, 1 ) ) heterozyg++;
329  }
330  return heterozyg;
331 }

References GetAllele().

Referenced by Vole_Base::SupplyHeteroZyg().

◆ HomozygosityCount()

int GeneticMaterial::HomozygosityCount ( )
314  {
315  // OK OK there is an easy way to do this by calling HeterozygosityCount and
316  // subtracting this from 32, but just is case that little bit of saved time is useful:
317  int homozyg=0;
318  for ( int i = 0; i < 32; i++ ) {
319  if ( GetAllele( i, 0 ) == GetAllele( i, 1 ) ) homozyg++;
320  }
321  return homozyg;
322 }

References GetAllele().

Referenced by Vole_Base::SupplyHomoZyg().

◆ Initiation()

void GeneticMaterial::Initiation ( AlleleFreq Al)

The method called to intialise genes on initiation of the simulation.
Gene frequencies are based on an external text file input read in on construction.

363  {
364  uint32 value;
365  for ( int l = 0; l < 32; l++ ) {
366  //if ( l < 16 ) c = 0; else if ( l < 24 ) c = 1; else c = 2;
367 
368  int chance = g_random_fnc( 1000 );
369  uint32 index = 0;
370  while ( chance > Al->SupplyAN( l, index ) ) {
371  index++;
372  }
373  value = index;
374  // set the value
375  SetAllele( l, value, 0 );
376  chance = g_random_fnc( 1000 );
377  index = 0;
378  while ( chance > Al->SupplyAN( l, index ) ) {
379  index++;
380  }
381  value = index;
382  // set the value
383  SetAllele( l, value, 1 );
384  }
385 }

References g_random_fnc(), SetAllele(), and AlleleFreq::SupplyAN().

Referenced by Vole_Population_Manager::Init().

◆ Mutation_1()

void GeneticMaterial::Mutation_1 ( )

random allele choice

442 {
443  for ( int i = 0; i < 16; i++ ) {
444  if ( g_rand_uni_fnc() < MutationChance ) // one chance in Mutation Chance
445  {
446  SetAllele( i, g_random_fnc( 4 ), g_random_fnc( 2 ) );
447  }
448  }
449  for ( int i = 16; i < 32; i++ ) {
450  if ( g_rand_uni_fnc() < MutationChance ) // one chance in Mutation Chance
451  {
452  SetAllele( i, g_random_fnc( 16 ), g_random_fnc( 2 ) );
453  }
454  }
455 }

References g_rand_uni_fnc(), g_random_fnc(), MutationChance, and SetAllele().

◆ Mutation_1ab()

void GeneticMaterial::Mutation_1ab ( )

random allele choice a & b only

463 {
464  // Only used when all loci have only two alleles!
465  for ( int i = 0; i < 32; i++ ) {
466  if ( g_rand_uni_fnc() < MutationChance ) // one chance in Mutation Chance
467  {
468  SetAllele( i, g_random_fnc( 2 ), g_random_fnc( 2 ) );
469  }
470  }
471 }

References g_rand_uni_fnc(), g_random_fnc(), MutationChance, and SetAllele().

◆ Mutation_2()

void GeneticMaterial::Mutation_2 ( )

Move one allele +/-

479 {
480  for ( int i = 0; i < 16; i++ ) {
481  if ( g_rand_uni_fnc() < MutationChance ) // one chance in Mutation Chance
482  {
483  int strand = g_random_fnc( 2 );
484  int allele = GetAllele( i, strand );
485  if ( g_random_fnc( 2 ) == 1 ) allele++; else allele--;
486  if ( allele == -1 ) allele = 3; else if ( allele == 4 ) allele = 0;
487  SetAllele( i, allele, strand );
488  }
489  }
490  for ( int i = 16; i < 32; i++ ) {
491  if ( g_rand_uni_fnc() < MutationChance ) // one chance in Mutation Chance
492  {
493  int strand = g_random_fnc( 2 );
494  int allele = GetAllele( i, strand );
495  if ( g_random_fnc( 2 ) == 1 ) allele++; else allele--;
496  if ( allele == -1 ) allele = 15; else if ( allele == 16 ) allele = 0;
497  SetAllele( i, allele, strand );
498  }
499  }
500 }

References g_rand_uni_fnc(), g_random_fnc(), GetAllele(), MutationChance, and SetAllele().

◆ Mutation_3()

void GeneticMaterial::Mutation_3 ( )

switch a<->b & c<->d

508 {
509  // NB Only works for the first 16 loci
510  for ( int i = 0; i < 16; i++ ) {
511  if ( g_rand_uni_fnc() < MutationChance ) // one chance in Mutation Chance
512  {
513  int strand = g_random_fnc( 2 );
514  int allele = GetAllele( i, strand );
515  switch ( allele ) {
516  case 0:
517  allele = 1;
518  break;
519  case 1:
520  allele = 0;
521  break;
522  case 2:
523  allele = 3;
524  break;
525  case 3:
526  allele = 2;
527  break;
528  }
529  SetAllele( i, allele, strand );
530  }
531  }
532 }

References g_rand_uni_fnc(), g_random_fnc(), GetAllele(), MutationChance, and SetAllele().

◆ Mutation_4()

void GeneticMaterial::Mutation_4 ( )

Specially mutation of only the first locus with two options 0/1

536 {
538  if (g_rand_uni_fnc() < MutationChance) // one chance in Mutation Chance
539  {
540  int strand = g_random_fnc(2);
541  int allele = GetAllele(0, strand);
542  allele = (allele + 1) && 1; // only 1 or 0 is allowed
543  }
544 }

References g_rand_uni_fnc(), g_random_fnc(), GetAllele(), and MutationChance.

◆ PrintChromosome()

void GeneticMaterial::PrintChromosome ( char *  C,
int  Chromosome 
)
255  {
256  for ( int i = 0; i < 16; i++ ) {
257  uint32 allele = GetAllele( i, Chromo );
258  switch ( allele ) {
259  case 0:
260  C[ i ] = 'a';
261  break;
262  case 1:
263  C[ i ] = 'b';
264  break;
265  case 2:
266  C[ i ] = 'c';
267  break;
268  case 3:
269  C[ i ] = 'd';
270  break;
271  case 4:
272  C[ i ] = 'e';
273  break;
274  case 5:
275  C[ i ] = 'f';
276  break;
277  case 6:
278  C[ i ] = 'g';
279  break;
280  case 7:
281  C[ i ] = 'h';
282  break;
283  case 8:
284  C[ i ] = 'i';
285  break;
286  case 9:
287  C[ i ] = 'j';
288  break;
289  case 10:
290  C[ i ] = 'k';
291  break;
292  case 11:
293  C[ i ] = 'l';
294  break;
295  case 12:
296  C[ i ] = 'm';
297  break;
298  case 13:
299  C[ i ] = 'n';
300  break;
301  case 14:
302  C[ i ] = 'o';
303  break;
304  case 15:
305  C[ i ] = 'p';
306  break;
307  }
308  }
309  C[ 16 ] = 0;
310 }

References GetAllele().

◆ ReadFrequencies()

void GeneticMaterial::ReadFrequencies ( )

◆ Recombine()

void GeneticMaterial::Recombine ( GeneticMaterial Gen21,
GeneticMaterial Gene2 
)
335  {
336  for ( int i = 0; i < 32; i++ ) {
337  // For each locus
338  // Choose which chromosome for each parent
339  int g0 = g_random_fnc( 2 );
340  int g1 = g_random_fnc( 2 );
341  // get the two alleles
342  uint32 a0 = Gene1->GetAllele( i, g0 );
343  uint32 a1 = Gene2->GetAllele( i, g1 );
344  // put a0 into chromo0 & a1 to chromo1 & vice versa
345  SetAllele( i, a0, 0 );
346  SetAllele( i, a1, 1 );
347  }
348 }

References g_random_fnc(), GetAllele(), and SetAllele().

◆ ScoreHQThreshold()

float GeneticMaterial::ScoreHQThreshold ( )

This function can be used to alter fitness based on associated genetic codes. These are only used in population genetic research, e.g. to create hybrid zones.

411  {
412  return 1.0;
413  // OLD CODE OUTDATED
414  /* Ditte's Simulation Version uint32 allele0a = GetAllele(0,0); // loci 0 uint32 allele1a = GetAllele(1,0); // loci 1
415  uint32 allele0b = GetAllele(0,1); // loci 0 uint32 allele1b = GetAllele(1,1); // loci 1
416  // Initial rules are that if 0a and 1a are 0 & 2 or 2 & 0 then OK (a,c)(c,a)
417  // likewise they may be 1 & 3 or 3 & 1 (b,d) (d,b) // any other combination is bad // Same for loci 1 bool IsOK0=false;
418  switch(allele0a) { case 0: if (allele1a==2) IsOK0=true; break; case 1: if (allele1a==3) IsOK0=true; break; case 2:
419  if (allele1a==0) IsOK0=true; break; case 3: if (allele1a==1) IsOK0=true; break; default: assert(NULL); break; }
420  bool IsOK1=false; switch(allele0b) { case 0: if (allele1b==2) IsOK1=true; break; case 1: if (allele1b==3) IsOK1=true; break;
421  case 2: if (allele1b==0) IsOK1=true; break; case 3: if (allele1b==1) IsOK1=true; break; default: assert(NULL); break; }
422  // determine the effect of the genetics // In the simple case it is good or bad
423  if (IsOK0 && IsOK1) return 1.0; else return 0.9;
424 
425  // Lar Bach's version uint32 allele2a = GetAllele(2,0); // loci 0 uint32 allele2b = GetAllele(2,1); // loci 0
426  float result= -0.5; switch (allele2a) { case 0: case 2: case 3: break; default: result+=0.5; } switch (allele2b) { case 0:
427  case 2: case 3: break; default: result+=0.5; break; } return result;
428  // returns -0.5 if homozygous aa, 0.5 if homozygous bb & het=0
429 
430  */
431 }

◆ ScoreReproduction()

float GeneticMaterial::ScoreReproduction ( )

This function can be used to alter reproductive effects based on genetic codes. These are only used in population genetic research.

391  {
392  return 1.0;
393  /* OLD CODE OUTDATED uint32 allele0a = GetAllele(0,0); // loci 0 uint32 allele1a = GetAllele(1,0); // loci 1
394  uint32 allele0b = GetAllele(0,1); // loci 0 uint32 allele1b = GetAllele(1,1); // loci 1
395  // Initial rules are that locus 0 and 1 are ac or bd then OK // likewise they may be ca or db
396  // any other combination is bad // Same for loci 1 bool IsOK0=false; switch(allele0a) { case 0: //a
397  if (allele1a==2) IsOK0=true; break; case 1: //b if (allele1a==3) IsOK0=true; break; case 2: //c
398  if (allele1a==0) IsOK0=true; break; case 3: //d if (allele1a==1) IsOK0=true; break; default:
399  FILE* errfile=fopen("GeneticErrorFile.Txt","w"); fprintf(errfile,"Unknown Allele Number\n"); fclose(errfile); exit(10);
400  break; } bool IsOK1=false; switch(allele0b) { case 0: if (allele1b==2) IsOK1=true; break; case 1:
401  if (allele1b==3) IsOK1=true; break; case 2: if (allele1b==0) IsOK1=true; break; case 3: if (allele1b==1) IsOK1=true; break;
402  default: FILE* errfile=fopen("GeneticErrorFile.Txt","w"); fprintf(errfile,"Unknown Allele Number\n"); fclose(errfile);
403  exit(11); break; } // determine the effect of the genetics // In the simple case it is good or bad
404  if (IsOK0 && IsOK1) return 1.0; else return 0.05; */
405 }

◆ SetAllele()

void GeneticMaterial::SetAllele ( int  pos,
uint32  value,
int  Chromosome 
)
193  {
194  Chromo*=3; // now 0 or 3
195  uint32 mask;
196  if (locus<16) {
197  // Get the right chromosome
198  // Create the mask
199  // Does it twice because 32 bits coding for 16 loci
200  mask = 0x03 << locus;
201  mask = mask << locus;
202  // just to make make sure it is 0-3
203  value = value & 0x03;
204  // create the value mask
205  value = value << locus;
206  value = value << locus;
207  // clear the locus
208  Chromosome[ Chromo ] &= ~mask;
209  // write the value
210  Chromosome[ Chromo ] |= value;
211  } else {
212  Chromo++; // now 1 or 4
213  locus-=16;
214  if (locus>=8) {
215  Chromo++;
216  locus-=8;
217  }
218  mask = 0x0F << (locus*4);
219  value = value & 0x0f; // make sure there was no extra stuff added!
220  // create the value mask
221  value = value << (locus*4);
222  Chromosome[ Chromo ] &= ~mask;
223  // write the value
224  Chromosome[ Chromo ] |= value;
225  }
226 }

References Chromosome.

Referenced by Initiation(), Mutation_1(), Mutation_1ab(), Mutation_2(), Mutation_3(), Recombine(), SetDirectFlag(), SetGeneticFlag(), UnsetDirectFlag(), and UnsetGeneticFlag().

◆ SetDirectFlag()

void GeneticMaterial::SetDirectFlag ( )
170  {
171  SetAllele(0,1,1);
172 }

References SetAllele().

Referenced by Vole_Base::SetDirectFlag().

◆ SetGeneticFlag()

void GeneticMaterial::SetGeneticFlag ( )
166  {
167  SetAllele(0,1,0);
168 }

References SetAllele().

Referenced by Vole_Base::SetGeneticFlag().

◆ UnsetDirectFlag()

void GeneticMaterial::UnsetDirectFlag ( )
179  {
180  SetAllele(0,0,1);
181 }

References SetAllele().

Referenced by Vole_Base::UnsetDirectFlag().

◆ UnsetGeneticFlag()

void GeneticMaterial::UnsetGeneticFlag ( )
175  {
176  SetAllele(0,0,0);
177 }

References SetAllele().

Referenced by Vole_Base::UnsetGeneticFlag().

Member Data Documentation

◆ Chromosome

uint32 GeneticMaterial::Chromosome[6]
protected

The documentation for this class was generated from the following files:
g_rand_uni_fnc
double g_rand_uni_fnc()
Definition: ALMaSS_Random.cpp:56
uint32
unsigned int uint32
Definition: ALMaSS_Setup.h:34
GeneticMaterial::Chromosome
uint32 Chromosome[6]
Definition: GeneticMaterial.h:97
GeneticMaterial::SetAllele
void SetAllele(int pos, uint32 value, int Chromosome)
Definition: GeneticMaterial.cpp:193
MutationChance
double MutationChance
Definition: GeneticMaterial.cpp:49
GeneticMaterial::GetAllele
uint32 GetAllele(int pos, int Chromosome)
Definition: GeneticMaterial.cpp:230
AlleleFreq::SupplyAN
int SupplyAN(int loc, int al)
Definition: GeneticMaterial.h:71
g_random_fnc
int g_random_fnc(const int a_range)
Definition: ALMaSS_Random.cpp:74