The AHA Model
Revision: 12809
Reference implementation 04 (HEDG02_04)
|
Definition the genetic architecture of the agent. More...
Data Types | |
type | gene |
This describes an individual gene object. See the genome structure for as general description and gene for details. More... | |
type | chromosome |
This type describes the chromosome object. Chromosome consists of an array of alleles and a descriptive string label. See "the genome structure" for as general description and "chromosome" for details. More... | |
type | individual_genome |
This type describes parameters of the individual agent's genome The genome is an array of allocatable the_genome::chromosome objects, different kinds of agents may have different genomes with different number of chromosomes. See "the genome structure" for as general description and "genome" for details. More... | |
Functions/Subroutines | |
elemental subroutine | chromosome_sort_rank_id (this) |
Sort GENE objects within the CHROMOSOME by their rank_id . The two subroutines qsort and qs_partition_rank_id are a variant of the recursive quick-sort algorithm adapted for sorting integer components of the the CHROMOSOME object. More... | |
subroutine | allele_init_random (this) |
Initialises allele with a random integer. Note that we do not set the labels for the alleles here during the random initialisation. More... | |
elemental subroutine | allele_create_zero (this) |
Create allele with zero value. We don't set labels for alleles here. More... | |
subroutine | allele_label_init_random (this) |
The (pair of) alleles here are assigned random string labels Not sure if that is necessary for any application though. More... | |
elemental subroutine | allele_label_set (this, label) |
Set labels for the alleles. The subroutine parameter is array of labels. More... | |
elemental character(len=label_length) function | allele_label_get (this) |
Get the i-th allele label. More... | |
elemental subroutine | allele_value_set (this, set_value, nr) |
Set a single value of the allele additive component. More... | |
pure subroutine | alleles_value_vector_set (this, values) |
Set values of the alleles as a vector, i.e. sets the whole gene values. More... | |
elemental integer function | allele_value_get (this, nr) |
Get the value of the allele. More... | |
pure subroutine | allele_values_vector_get (this, values) |
Get the vector of all values of the alleles, i.e. gets the gene values. More... | |
elemental subroutine | allele_rank_id_set (this, value_id) |
subroutine | allele_mutate_random (this, prob) |
Introduce a random point mutation to a random allele component. More... | |
subroutine | allele_mutate_random_batch (this, prob) |
Introduce a random mutation of the whole set of additive allele components. More... | |
subroutine | chromosome_init_allocate_random (this, length, label) |
This subroutine initialises the chromosome with, and allocates, random alleles, sets one of them randomly dominant and optionally defines the chromosome label. More... | |
subroutine | chromosome_create_allocate_zero (this, length, label) |
Init a new chromosome, zero, non-random. More... | |
elemental subroutine | chromosome_recalculate_rank_ids (this) |
This subroutine recalculates rank_id indices for consecutive gene objects within the chromosome. This may be necessary after reordering by random relocation mutation. More... | |
subroutine | chromosome_mutate_relocate_swap_random (this, prob) |
Mutate within the same chromosome, relocate a gene (unit of alleles) to a different random position within the same chromosome, the misplaced gene moves to the relocated gene position, so they are just swap. More... | |
subroutine | chromosome_mutate_relocate_shift_random (this, prob) |
Mutate within the same chromosome, relocate a gene (unit of alleles) to a different random position within the same chromosome, shifting all other genes within the chromosome down one position. This works as follows: first, we randomly determine the gene to relocate, assign it a new random rank_id. Then re-sort the chromosome according to the new ranks with qsort with the qs_partition_rank_id backend. More... | |
subroutine | genome_init_random (this, label) |
Initialise the genome at random, and set sex as determined by the sex determination locus. More... | |
subroutine | genome_create_zero (this) |
Create a new empty genome, and set sex as determined by the sex determination locus. Genome values are from parents using inherit functions. More... | |
subroutine | genome_label_set (this, label) |
Label genome. If label is not provided, make a random string. More... | |
elemental character(len=label_length) function | genome_label_get (this) |
Accessor function to get the genome label. The label is a kind of a (random) text string name of the genome and the individual agent. More... | |
subroutine | genome_sex_determine_init (this) |
Sex determination initialisation subroutine. More... | |
elemental logical function | genome_get_sex_is_male (this) |
Get the logical sex ID of the genome object component. More... | |
elemental logical function | genome_get_sex_is_female (this) |
Get the logical sex ID of the genome object component. More... | |
elemental character(len=label_length) function | genome_get_sex_label (this) |
Get the descriptive sex label: male or female. More... | |
elemental subroutine | genome_individual_set_alive (this) |
Set the individual to be alive, normally this function is used after init or birth. More... | |
elemental subroutine | genome_individual_set_dead (this) |
Set the individual to be dead. Note that this function does not deallocate the individual agent object, this may be a separate destructor function. More... | |
elemental logical function | genome_individual_check_alive (this) |
Check if the individual is alive. More... | |
elemental logical function | genome_individual_check_dead (this) |
Check if the individual is dead (the opposite of is_alive ). More... | |
subroutine | genome_individual_recombine_homol_full_rand_alleles (this, mother, father, exchange_ratio) |
Internal genetic recombination backend, exchange individual alleles between homologous chromosomes in mother and father genomes to form the this (offspring) genome. Fully random recombination. More... | |
subroutine | genome_individual_recombine_homol_part_rand_alleles (this, mother, father, exchange_ratio) |
Internal genetic recombination backend, exchange individual alleles between homologous chromosomes in mother and father genomes to form the this (offspring) genome. Partially random recombination, identical across the homologous chromosomes. More... | |
subroutine | genome_individual_crossover_homol_fix (this, mother, father, pattern_matrix) |
Internal fixed genetic crossover backend, exchange blocks of alleles between homologous chromosomes in mother and father genomes to form the this (offspring) genome. More... | |
subroutine | genome_mutate_wrapper (this, p_point, p_set, p_swap, p_shift) |
Perform a probabilistic random mutation(s) on the individual genome. This is a high level wrapper to build mutations from various components. More... | |
Neuronal response functions | |
There are two separate functions that produce a trait value from the genotype. The procedure | |
subroutine | trait_init_genotype_gamma2gene (this, this_trait, g_p_matrix, init_val, gerror_cv, label) |
Initialise an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy. More... | |
subroutine | trait_set_genotype_gamma2gene (this, this_trait, g_p_matrix, init_val, gerror_cv) |
Set an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy. More... | |
The genotype to phenotype functions based on fixed linear | |
transformation. | |
subroutine | trait_init_linear_sum_additive_comps_2genes_r (this, this_trait, g_p_matrix, phenotype_min, phenotype_max, label) |
Initialise an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy. More... | |
subroutine | trait_set_linear_sum_additive_comps_2genes_r (this, this_trait, g_p_matrix, phenotype_min, phenotype_max) |
Set an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy. More... | |
Variables | |
character(len= *), parameter, private | modname = "(THE_GENOME)" |
Definition the genetic architecture of the agent.
This module defines the genetic architecture objects of the agent. See The genome structure for an overview.
elemental subroutine the_genome::chromosome_sort_rank_id | ( | class(chromosome), intent(inout) | this | ) |
Sort GENE objects within the CHROMOSOME by their rank_id
. The two subroutines qsort
and qs_partition_rank_id
are a variant of the recursive quick-sort algorithm adapted for sorting integer components of the the CHROMOSOME
object.
Definition at line 323 of file m_genome.f90.
subroutine the_genome::allele_init_random | ( | class(gene), intent(inout) | this | ) |
Initialises allele with a random integer. Note that we do not set the labels for the alleles here during the random initialisation.
The allele components are initialised by a random integers within the range ALLELERANGE_MIN
and ALLELERANGE_MAX
parameetr values.
Definition at line 410 of file m_genome.f90.
elemental subroutine the_genome::allele_create_zero | ( | class(gene), intent(inout) | this | ) |
Create allele with zero value. We don't set labels for alleles here.
allele_value
array as it has fixed shape.Definition at line 424 of file m_genome.f90.
subroutine the_genome::allele_label_init_random | ( | class(gene), intent(inout) | this | ) |
The (pair of) alleles here are assigned random string labels Not sure if that is necessary for any application though.
Definition at line 436 of file m_genome.f90.
elemental subroutine the_genome::allele_label_set | ( | class(gene), intent(inout) | this, |
character(len=*), intent(in) | label | ||
) |
Set labels for the alleles. The subroutine parameter is array of labels.
[in] | label | label, provides an array of labels to set for the allele. |
Definition at line 445 of file m_genome.f90.
elemental character(len=label_length) function the_genome::allele_label_get | ( | class(gene), intent(in) | this | ) |
Get the i-th allele label.
Definition at line 456 of file m_genome.f90.
elemental subroutine the_genome::allele_value_set | ( | class(gene), intent(inout) | this, |
integer, intent(in) | set_value, | ||
integer, intent(in) | nr | ||
) |
Set a single value of the allele additive component.
[in] | set_value | value, provides the value to set for the allele and the allele number. |
[in] | nr | number, provides the number of the allele component to set. |
Definition at line 467 of file m_genome.f90.
pure subroutine the_genome::alleles_value_vector_set | ( | class(gene), intent(inout) | this, |
integer, dimension(additive_comps), intent(in) | values | ||
) |
Set values of the alleles as a vector, i.e. sets the whole gene values.
[in] | values | values, provides vector of values to set for the alleles. |
Definition at line 482 of file m_genome.f90.
elemental integer function the_genome::allele_value_get | ( | class(gene), intent(in) | this, |
integer, intent(in) | nr | ||
) |
Get the value of the allele.
[in] | nr | number, provides the number of the allele component to set. |
nr
's allele. Definition at line 493 of file m_genome.f90.
pure subroutine the_genome::allele_values_vector_get | ( | class(gene), intent(in) | this, |
integer, dimension(additive_comps), intent(out) | values | ||
) |
Get the vector of all values of the alleles, i.e. gets the gene values.
[out] | values | values, Gets the vector of the values for the alleles. |
Definition at line 508 of file m_genome.f90.
elemental subroutine the_genome::allele_rank_id_set | ( | class(gene), intent(inout) | this, |
integer, intent(in) | value_id | ||
) |
[in] | value_id | rank_id, set this value to the allele rank_id . |
Definition at line 518 of file m_genome.f90.
subroutine the_genome::allele_mutate_random | ( | class(gene), intent(inout) | this, |
real(srp), intent(in), optional | prob | ||
) |
Introduce a random point mutation to a random allele component.
[in] | prob | prob optional probability of mutation, if absent, the default value commondata::mutationrate_point is used. |
Do mutate if a random value is smaller than the commondata::mutationrate parameter constant.
First, determine which of the alleles components gets mutation.
Second, change the value of this allele component to an random integer.
Definition at line 530 of file m_genome.f90.
subroutine the_genome::allele_mutate_random_batch | ( | class(gene), intent(inout) | this, |
real(srp), intent(in), optional | prob | ||
) |
Introduce a random mutation of the whole set of additive allele components.
[in] | prob | prob optional probability of mutation, if absent, the default value commondata::mutationrate_batch is used. |
This mutation just re-init the whole allele set as random.
Definition at line 562 of file m_genome.f90.
subroutine the_genome::chromosome_init_allocate_random | ( | class(chromosome), intent(inout) | this, |
integer, intent(in) | length, | ||
character(len=*), intent(in), optional | label | ||
) |
This subroutine initialises the chromosome with, and allocates, random alleles, sets one of them randomly dominant and optionally defines the chromosome label.
[in] | length,sets | the length of the chromosome object, N of alleles. |
[in] | label,sets | the label for the chromosome object, optional. |
We set the chromosome label if such a parameter is provided, or a random string if not.
First, set the chromosome length using the procedure parameter length
.
Then, allocate the array of the allele objects with this length.
Initialise all the alleles within this chromosome.
Specifically, initialise the allele.
Set initial rank_id ID of the allele.
Finally, set the label for the alleles within this chromosome. Labels can be set random using this function (disabled):
But in this implementation we construct the label for the allele from the chromosome label and the allele number;
Long chromosome labels are trimmed at right to fit the allele number.
Definition at line 596 of file m_genome.f90.
subroutine the_genome::chromosome_create_allocate_zero | ( | class(chromosome), intent(inout) | this, |
integer, intent(in) | length, | ||
character(len=*), intent(in), optional | label | ||
) |
Init a new chromosome, zero, non-random.
Set the chromosome label if provided as parameter, or random string if not.
First, set the chromosome length using the procedure parameter length
.
Then, we allocate the array of the allele objects with this length.
Initialise all the alleles within this chromosome.
do concurrent
construct is used here.Definition at line 647 of file m_genome.f90.
elemental subroutine the_genome::chromosome_recalculate_rank_ids | ( | class(chromosome), intent(inout) | this | ) |
This subroutine recalculates rank_id indices for consecutive gene objects within the chromosome. This may be necessary after reordering by random relocation mutation.
Definition at line 697 of file m_genome.f90.
subroutine the_genome::chromosome_mutate_relocate_swap_random | ( | class(chromosome), intent(inout) | this, |
real(srp), intent(in), optional | prob | ||
) |
Mutate within the same chromosome, relocate a gene (unit of alleles) to a different random position within the same chromosome, the misplaced gene moves to the relocated gene position, so they are just swap.
[in] | prob | prob optional probability of mutation, if absent, the default value commondata::relocation_swap_rate is used. |
Do mutate if a random value is smaller than the commondata::relocation_swap_rate parameter constant value.
If yes, randomly determine the gene (gene_move
) that initiates the mutation move within the chromosome.
Randomly determine the gene that will be swapped with the gene_move
.
Then, cycle through the alleles and select new random allele if it happens to coincide with gene_move
.
After this, do the gene swap, and gene rank_id ID swap.
Definition at line 714 of file m_genome.f90.
subroutine the_genome::chromosome_mutate_relocate_shift_random | ( | class(chromosome), intent(inout) | this, |
real(srp), intent(in), optional | prob | ||
) |
Mutate within the same chromosome, relocate a gene (unit of alleles) to a different random position within the same chromosome, shifting all other genes within the chromosome down one position. This works as follows: first, we randomly determine the gene to relocate, assign it a new random rank_id. Then re-sort the chromosome according to the new ranks with qsort
with the qs_partition_rank_id
backend.
[in] | prob | prob optional probability of mutation, if absent, the default value commondata::relocation_shift_rate is used. |
Do mutate if a random value is smaller than the commondata::relocation_shift_rate parameter constant value.
If yes, randomly determine the gene that initiates move within the chromosome.
Randomly determine the new position of this gene.
Then, cycle through alleles and select new random allele if it happens to coincide with gene_move
.
After the cycle, adjust rank_id's of the alleles.
Then re-sort the allele objects by their updated rank_id's.
Finally, recalculate rank_id's so they are again ordered.
Definition at line 772 of file m_genome.f90.
subroutine the_genome::genome_init_random | ( | class(individual_genome), intent(inout) | this, |
character(len=*), intent(in), optional | label | ||
) |
Initialise the genome at random, and set sex as determined by the sex determination locus.
First, create spatial moving object component of the individual genome. But we do not yet position the genome object.
Allocate the genome object, it must have genome_size
chromosomes and CHROMOSOME_PLOIDY
homologs.
Now cycle over all the chromosomes and homologs and initialise each of them.
On exit from the cycle, set the genome label if provided, or random string if not.
Then, determine the sex of the genome by the genome sex determination locus taking into account the sex ratio.
Definition at line 827 of file m_genome.f90.
subroutine the_genome::genome_create_zero | ( | class(individual_genome), intent(inout) | this | ) |
Create a new empty genome, and set sex as determined by the sex determination locus. Genome values are from parents using inherit functions.
First of all, Create spatial moving object component of the individual genome.
Allocate the genome object, it must have genome_size
chromosomes and CHROMOSOME_PLOIDY
homologs
Now cycle over all the chromosomes and homologs and initialise each of them with empty genes (zero).
Initialise the label the genome with a random string.
Determine the sex of the genome by the genome sex determination locus taking into account the sex ratio.
Definition at line 874 of file m_genome.f90.
subroutine the_genome::genome_label_set | ( | class(individual_genome), intent(inout) | this, |
character(len=*), intent(in), optional | label | ||
) |
Label genome. If label is not provided, make a random string.
Set the genome label if provided, or random string if not.
Definition at line 912 of file m_genome.f90.
elemental character(len=label_length) function the_genome::genome_label_get | ( | class(individual_genome), intent(in) | this | ) |
Accessor function to get the genome label. The label is a kind of a (random) text string name of the genome and the individual agent.
Definition at line 935 of file m_genome.f90.
subroutine the_genome::genome_sex_determine_init | ( | class(individual_genome), intent(inout) | this | ) |
Sex determination initialisation subroutine.
Determine the genome's sex, sex is set by a logical identifier, sex_is_male
TRUE is male. Sex is calculated from the genome and based on the average values of the sex determination alleles in homologous chromosomes, rescaled to 0:1. This rescaled value is then compared with the sex ratio parameter.
The implementation is based on genotype x phenotype matrix (logical type): commondata::sex_genotype_phenotype.
First, initialise the average sex locus sum across the homologous chromosomes.
Then loop across homologs, chromosomes and alleles until the value of SEX_GENOTYPE_PHENOTYPE
gets TRUE. This means it is the sex locus.
Upon exit from the loop, check if the average sex locus across all homologous chromosomes and additive allele components, scaled to 0:1 is less than the SEX_RATIO
, the subject becomes the male genotype.
Otherwise, the subject becomes the female.
Definition at line 955 of file m_genome.f90.
elemental logical function the_genome::genome_get_sex_is_male | ( | class(individual_genome), intent(in) | this | ) |
Get the logical sex ID of the genome object component.
Definition at line 1016 of file m_genome.f90.
elemental logical function the_genome::genome_get_sex_is_female | ( | class(individual_genome), intent(in) | this | ) |
Get the logical sex ID of the genome object component.
Definition at line 1028 of file m_genome.f90.
elemental character(len=label_length) function the_genome::genome_get_sex_label | ( | class(individual_genome), intent(in) | this | ) |
Get the descriptive sex label: male or female.
Definition at line 1044 of file m_genome.f90.
elemental subroutine the_genome::genome_individual_set_alive | ( | class(individual_genome), intent(inout) | this | ) |
Set the individual to be alive, normally this function is used after init or birth.
Definition at line 1057 of file m_genome.f90.
elemental subroutine the_genome::genome_individual_set_dead | ( | class(individual_genome), intent(inout) | this | ) |
Set the individual to be dead. Note that this function does not deallocate the individual agent object, this may be a separate destructor function.
The dies
method is implemented at the following levels of the agent object hierarchy (upper overrides the lower level):
Definition at line 1076 of file m_genome.f90.
elemental logical function the_genome::genome_individual_check_alive | ( | class(individual_genome), intent(in) | this | ) |
Check if the individual is alive.
Definition at line 1085 of file m_genome.f90.
elemental logical function the_genome::genome_individual_check_dead | ( | class(individual_genome), intent(in) | this | ) |
Check if the individual is dead (the opposite of is_alive
).
Definition at line 1095 of file m_genome.f90.
subroutine the_genome::genome_individual_recombine_homol_full_rand_alleles | ( | class(individual_genome), intent(inout) | this, |
class(individual_genome), intent(in) | mother, | ||
class(individual_genome), intent(in) | father, | ||
real(srp), intent(in), optional | exchange_ratio | ||
) |
Internal genetic recombination backend, exchange individual alleles between homologous chromosomes in mother and father genomes to form the this
(offspring) genome. Fully random recombination.
[in] | mother | mother The mother object the_genome::individual_genome class. |
[in] | father | father The father object the_genome::individual_genome class. |
acomp_vect_mother and acomp_vect_father are the vectors of additive allele components that are obtained from the mother and the father.
Check optional exchange_ratio
parameter that defines the ratio of the alleles that are inherited from the mother. If absent, get the default value from the commondata::genome_recombination_ratio_mother parameter.
Loop through the chromosomes (CHROMOSOMES
block over i
) and their homologues (HOMOLOGS
block over j
) and set the genetic make up of the this
object from the genes of the mother
and the father
objects.
The outline of the main loop is:
do concurrent
loops here due to non-pure random call of the RAND_R4()
function in the innermost loop (ALLELES
block over k
).First, get the number of alleles for each of the chromosomes: n_alleles
.
Then, loop over the n_alleles
alleles (ALLELES
block over k
):
Randomly select if this specific allele (k) is copied from the mother or is a subject to (random) recombination and gets the values from the father. This is determined stochastically using the exchange_ratio
ratio dummy argument or the commondata::genome_recombination_ratio_mother
parameter as the probability to get the allele from the mother.
Finally, set the vector of additive allele components of the this
agent from the mother's vector.
Finally, set the vector of additive allele components of the this
agent from the father's vector.
Definition at line 1116 of file m_genome.f90.
subroutine the_genome::genome_individual_recombine_homol_part_rand_alleles | ( | class(individual_genome), intent(inout) | this, |
class(individual_genome), intent(in) | mother, | ||
class(individual_genome), intent(in) | father, | ||
real(srp), intent(in), optional | exchange_ratio | ||
) |
Internal genetic recombination backend, exchange individual alleles between homologous chromosomes in mother and father genomes to form the this
(offspring) genome. Partially random recombination, identical across the homologous chromosomes.
[in] | mother | mother The mother object the_genome::individual_genome class. |
[in] | father | father The father object the_genome::individual_genome class. |
acomp_vect_mother and acomp_vect_father are the vectors of additive allele components that are obtained from the mother and the father.
Check optional exchange_ratio
parameter that defines the ratio of the alleles that are inherited from the mother. If absent, get the default value from the commondata::genome_recombination_ratio_mother parameter.
Loop through the chromosomes (CHROMOSOMES
block over i
) and their and set the genetic make up of the this
object from the genes of the mother
and the father
objects.
The outline of the main loop is:
do concurrent
loops here due to non-pure random call of the RAND_R4()
function in the innermost loop (ALLELES
block over k
).First, get the number of alleles for each of the chromosomes: n_alleles
.
Then, loop over the n_alleles
alleles (ALLELES
block over k
):
Randomly select if this specific allele (k) is copied from the mother or is a subject to (random) recombination and gets the values from the father. This is determined stochastically using the exchange_ratio
ratio dummy argument or the commondata::genome_recombination_ratio_mother
parameter as the probability to get the allele from the mother.
Finally, loop through the homologues of the i
th chromosome (HOMOLOGS
block over j
) and set the same allele across all homologues the same way (transferred from mother or the father). Thus, all homologues have identical non-random recombination pattern, either from the mother or from the father.
Finally, set the vector of additive allele components of the this
agent from the mother's vector.
Finally, set the vector of additive allele components of the this
agent from the father's vector.
Definition at line 1227 of file m_genome.f90.
subroutine the_genome::genome_individual_crossover_homol_fix | ( | class(individual_genome), intent(inout) | this, |
class(individual_genome), intent(in) | mother, | ||
class(individual_genome), intent(in) | father, | ||
logical, dimension(max_nalleles,n_chromosomes), intent(in), optional | pattern_matrix | ||
) |
Internal fixed genetic crossover backend, exchange blocks of alleles between homologous chromosomes in mother and father genomes to form the this
(offspring) genome.
[in] | mother | mother The mother object the_genome::individual_genome class. |
[in] | father | father The father object the_genome::individual_genome class. |
[in] | pattern_matrix | pattern_matrix is an optional boolean pattern matrix that determines the pattern of fixed chromosome crossover. For each chromosome, the alleles that are marked with the TRUE (YES) values are inherited from the mother whereas those marked FALSE (NO) are inherited from the father. see commondata::genome_crossover_fixed_mother. |
acomp_vect_mother and acomp_vect_father are the vectors of additive allele components that are obtained from the mother and the father.
Check optional pattern_matrix
parameter matrix. If it is not provided, use the default commondata::genome_crossover_fixed_mother.
Loop through the chromosomes (CHROMOSOMES
block over i
) and their homologues (HOMOLOGS
block over j
) and set the genetic make up of the this
object from the genes of the mother
and the father
objects.
First, get the number of alleles for each of the chromosomes: n_alleles
.
Then, loop over the n_alleles
alleles (ALLELES
block over k
):
Check the boolean pattern matrix pattern_matrix_loc
value, if it is TRUE, the allele value is copied from the mother, otherwise from the father.
Finally, set the vector of additive allele components of the this
agent from the mother's vector.
Finally, set the vector of additive allele components of the this
agent from the father's vector.
Definition at line 1347 of file m_genome.f90.
subroutine the_genome::trait_init_genotype_gamma2gene | ( | class(individual_genome), intent(inout) | this, |
real(srp), intent(out) | this_trait, | ||
logical, dimension(:,:), intent(in) | g_p_matrix, | ||
real(srp), intent(in) | init_val, | ||
real(srp), intent(in), optional | gerror_cv, | ||
character(len=*), intent(in) | label | ||
) |
Initialise an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy.
Create a new trait object from the genotype according to the boolean genotype x phenotype matrix g_p_matrix
by commondata::gamma2gene()
function and the init_val
baseline value. It also introduces an initial Gaussian variance with the coefficient of variation=gerror_cv
. Normally, g_p_matrix
, init_val
and gerror_cv
are parameters set in commondata
for each specific trait. For example, for thyroid hormone (the_hormones::hormones
class level up) they are:
this
, the agent as it modifies the label. Therefore, it should be mainly used for initialisation of agent traits. commondata::chromosome_ploidy = 1
, in such a case there is no need to select random homologous chromosome and it is impossible to set the two parameters of the commondata::gamma2gene
function (there is only a single chromosome). [out] | this_trait | this_trait the component of the individual agent object hierarchy that we are initialising in. |
[in] | g_p_matrix | g_p_matrix matrix genotype x phenotype, the matrix setting the correspondence between the genome and the phenotype. |
[in] | init_val | init_val start value, trait initialisation baseline value parameter, e.g. for thyroid it is commondata::thyroid_init . |
[in] | gerror_cv | gerror_cv error value, trait Gaussian error value, e.g. for thyroid it is commondata::thyroid_gerror_cv , if absent, we don't introduce random error and the initial trait values are deterministic by the genome. |
[in] | label | label a label for the allele locus that sets the phenotypic object. |
Set trait values from the genome using the g_p_matrix
matrix. We first need to select two chromosomes from the available set (normally two, but possibly more) for input into commondata::gamma2gene() function. We do it random.
Then, cycle through and select a different chromosomes (i.e. cycle if happens to coincide with the first). (Use do while chromosome1 = chromosome2
construct.)
As a result there are two distinct chromosomes k1
and k2
. This unique chromosomes selection part of the code precludes the use of haploid genome.
Definition at line 1471 of file m_genome.f90.
subroutine the_genome::trait_set_genotype_gamma2gene | ( | class(individual_genome), intent(in) | this, |
real(srp), intent(out) | this_trait, | ||
logical, dimension(:,:), intent(in) | g_p_matrix, | ||
real(srp), intent(in) | init_val, | ||
real(srp), intent(in), optional | gerror_cv | ||
) |
Set an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy.
This is almost the same as the_genome::trait_init_genotype_gamma2gene()
, but does not modify the this
object (it has intent [in]). Therefore it should be used for setting such traits as behavioural expectancies.
this
, the agent as it does not modify the object (e.g. does not set label). Therefore, it can be used in assessing the subjective expectancies. commondata::chromosome_ploidy = 1
, in such a case there is no need to select random homologous chromosome and it is impossible to set the two parameters of the commondata::gamma2gene()
function (there is only a single chromosome). [out] | this_trait | this_trait the component of the individual agent object hierarchy that we are initialising. |
[in] | g_p_matrix | matrix genotype x phenotype, the matrix setting the correspondence between the genome and the phenotype. |
[in] | init_val | init_val start value, trait initialisation baseline value parameter, e.g. for thyroid it is commondata::thyroid_init . |
[in] | gerror_cv | gerror_cv error value, trait Gaussian error value, e.g. for thyroid it is commondata::thyroid_gerror_cv , if absent, we don't introduce random error and the initial trait values are deterministic by the genome. |
Set trait values from the genome using the g_p_matrix
matrix. We first need to select two chromosomes from the available set (normally two, but possibly more) for input into commondata::gamma2gene() function. We do it random.
Then, we cycle through and select a different chromosomes (i.e. cycle if happens to coincide with the first.
As a result there are two distinct chromosomes k1
and k2
. This unique chromosomes selection part of the code precludes the use of haploid genome.
commondata::gamma2gene()
. Definition at line 1561 of file m_genome.f90.
subroutine the_genome::trait_init_linear_sum_additive_comps_2genes_r | ( | class(individual_genome), intent(inout) | this, |
real(srp), intent(out) | this_trait, | ||
logical, dimension(:,:), intent(in) | g_p_matrix, | ||
real(srp), intent(in) | phenotype_min, | ||
real(srp), intent(in) | phenotype_max, | ||
character(len=*), intent(in) | label | ||
) |
Initialise an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy.
Create a new trait object from the genotype according to the boolean genotype x phenotype matrix g_p_matrix
by a simple linear scaling transformation function. Normally, g_p_matrix
, is a parameters set in commondata for each specific trait.
[out] | this_trait | this_trait the component of the individual agent object hierarchy that we are initialising in. |
[in] | g_p_matrix | g_p_matrix matrix genotype x phenotype, the matrix setting the correspondence between the genome and the phenotype. |
[in] | phenotype_min | init_val start value, trait initialisation baseline value parameter, e.g. for thyroid it is commondata::thyroid_init . |
[in] | label | label a label for the allele locus that sets the phenotypic object. |
Set trait values from the genome using the g_p_matrix
matrix. We first need to select two chromosomes from the available set (normally two, but possibly more) for input into gamma2gene
. We do it random. Assign the first chromosome.
Also assign the second chromosome.
Then, cycle through and select a different chromosomes (i.e. cycle if happens to coincide with the first). (Use do while chromosome1 = chromosome2` construct.)
As a result there are two distinct chromosomes k1
and k2
. This unique chromosomes selection part of the code precludes the use of haploid genome.
phenotype_min
and phenotype_max
parameters.Definition at line 1646 of file m_genome.f90.
subroutine the_genome::trait_set_linear_sum_additive_comps_2genes_r | ( | class(individual_genome), intent(in) | this, |
real(srp), intent(out) | this_trait, | ||
logical, dimension(:,:), intent(in) | g_p_matrix, | ||
real(srp), intent(in) | phenotype_min, | ||
real(srp), intent(in) | phenotype_max | ||
) |
Set an individual trait of the agent that depends on the genotype. This can be any trait upwards in the class hierarchy.
This is almost the same as the_genome::trait_init_linear_sum_additive_comps_2genes_r(), but does not modify the this
object (it has intent [in]). Therefore it should be used for setting such traits as behavioural expectancies. Create a new trait object from the genotype according to the boolean genotype x phenotype matrix g_p_matrix
by a simple linear scaling transformation function. Normally, g_p_matrix
, is a parameters set in commondata for each specific trait.
[out] | this_trait | this_trait the component of the individual agent object hierarchy that we are initialising in. |
[in] | g_p_matrix | g_p_matrix matrix genotype x phenotype, the matrix setting the correspondence between the genome and the phenotype. |
[in] | phenotype_min | init_val start value, trait initialisation baseline value parameter, e.g. for thyroid it is commondata::thyroid_init . |
Set trait values from the genome using the g_p_matrix
matrix. We first need to select two chromosomes from the available set (normally two, but possibly more) for input into gamma2gene
. We do it random. Assign the first chromosome.
Also assign the second chromosome.
Then, cycle through and select a different chromosomes (i.e. cycle if happens to coincide with the first). (Use do while chromosome1 = chromosome2` construct.)
As a result there are two distinct chromosomes k1
and k2
. This unique chromosomes selection part of the code precludes the use of haploid genome.
phenotype_min
and phenotype_max
parameters.Definition at line 1742 of file m_genome.f90.
subroutine the_genome::genome_mutate_wrapper | ( | class(individual_genome), intent(inout) | this, |
real(srp), intent(in), optional | p_point, | ||
real(srp), intent(in), optional | p_set, | ||
real(srp), intent(in), optional | p_swap, | ||
real(srp), intent(in), optional | p_shift | ||
) |
Perform a probabilistic random mutation(s) on the individual genome. This is a high level wrapper to build mutations from various components.
[in] | p_point | p_point optional probability of mutation, if absent, the default value commondata::mutationrate_point is used. |
[in] | p_set | p_set optional probability of mutation, if absent, the default value commondata::mutationrate_batch is used. |
[in] | p_swap | p_swap optional probability of mutation, if absent, the default value commondata::mutationrate_batch is used. |
[in] | p_shift | p_shift optional probability of mutation, if absent, the default value commondata::mutationrate_batch is used. |
Each of the mutations below is a random process that occurs with a specific probability that is set by its respective mutation rate parameter.
Call for a point mutation on a whole random allele (batch of allele components) of a randomly chosen chromosome using the the_genome::gene::mutate_set() method.
Definition at line 1820 of file m_genome.f90.
|
private |
Definition at line 27 of file m_genome.f90.