25 character (len=*),
parameter,
private ::
modname =
"(THE_HORMONES)"
32 real(
srp) :: growhorm_level
35 real(
srp) :: thyroid_level
40 real(
srp) :: adrenaline_level
44 real(
srp) :: cortisol_level
54 real(
srp) :: testosterone_level
57 real(
srp) :: estrogen_level
59 real(
srp) :: testosterone_baseline
61 real(
srp) :: estrogen_baseline
63 real(
srp),
dimension(HISTORY_SIZE_AGENT_PROP) :: testosterone_history
64 real(
srp),
dimension(HISTORY_SIZE_AGENT_PROP) :: estrogen_history
119 procedure,
public :: testosterone_base_get => &
129 procedure,
public :: reproductive_factor => &
149 class(
hormones),
intent(inout) :: this
153 call this%trait_init(this%growhorm_level, &
157 call this%trait_init(this%thyroid_level, &
161 call this%trait_init(this%adrenaline_level, &
165 call this%trait_init(this%cortisol_level, &
169 call this%trait_init(this%testosterone_level, &
174 call this%trait_init(this%estrogen_level, &
180 this%testosterone_baseline = this%testosterone_level
181 this%estrogen_baseline = this%estrogen_level
184 call this%hormone_history_clean()
187 call add_to_history(this%testosterone_history, this%testosterone_level)
196 class(
hormones),
intent(inout) :: this
198 this%testosterone_history =
missing
199 this%estrogen_history =
missing
206 class(
hormones),
intent(inout) :: this
209 call add_to_history(this%testosterone_history, this%testosterone_level)
224 real(
srp) :: reprfact
226 if ( this%is_male() )
then
227 reprfact = this%testosterone_level
229 reprfact = this%estrogen_level
249 real(
srp) :: value_get
251 value_get = this%growhorm_level
258 class(
hormones),
intent(inout) :: this
261 real(
srp),
intent(in) :: value_set
263 this%growhorm_level = value_set
273 real(
srp) :: value_get
275 value_get = this%thyroid_level
282 class(
hormones),
intent(inout) :: this
285 real(
srp),
intent(in) :: value_set
287 this%thyroid_level = value_set
297 real(
srp) :: value_get
299 value_get = this%adrenaline_level
306 class(
hormones),
intent(inout) :: this
309 real(
srp),
intent(in) :: value_set
311 this%adrenaline_level = value_set
321 real(
srp) :: value_get
323 value_get = this%cortisol_level
330 class(
hormones),
intent(inout) :: this
333 real(
srp),
intent(in) :: value_set
335 this%cortisol_level = value_set
345 real(
srp) :: value_get
347 value_get = this%testosterone_level
354 class(
hormones),
intent(inout) :: this
357 real(
srp),
intent(in) :: value_set
362 logical,
optional,
intent(in) :: update_history
364 this%testosterone_level = value_set
366 if (
present(update_history))
then
367 if (.NOT. update_history)
return
379 real(
srp) :: value_get
381 value_get = this%estrogen_level
388 class(
hormones),
intent(inout) :: this
391 real(
srp),
intent(in) :: value_set
396 logical,
optional,
intent(in) :: update_history
398 this%estrogen_level = value_set
400 if (
present(update_history))
then
401 if (.NOT. update_history)
return
414 real(
srp) :: value_get
416 value_get = this%testosterone_baseline
426 real(
srp) :: value_get
428 value_get = this%estrogen_baseline
Simple history stack function, add to the end of the stack. We need only to add components on top (en...
COMMONDATA – definitions of global constants and procedures.
character(len= *), parameter, private modname
MODNAME always refers to the name of the current module for use by the LOGGER function LOG_DBG....
real(srp), parameter, public adrenaline_init
Genotype to phenotype gamma2gene initialisation value for adrenaline
logical, dimension(max_nalleles, n_chromosomes), parameter, public thyroid_genotype_phenotype
Genotype x Phenotype matrix for thyroid.
logical, dimension(max_nalleles, n_chromosomes), parameter, public estrogen_genotype_phenotype
Genotype x Phenotype matrix for ESTROGEN.
real(srp), parameter, public testosterone_gerror_cv
Genotype to phenotype gamma2gene Gaussian error parameter. This is really the coefficient of variatio...
real(srp), parameter, public thyroid_gerror_cv
Genotype to phenotype gamma2gene Gaussian error parameter. This is really the coefficient of variatio...
integer, parameter, public srp
Definition of the standard real type precision (SRP).
logical, dimension(max_nalleles, n_chromosomes), parameter, public adrenaline_genotype_phenotype
Genotype x Phenotype matrix for adrenaline
logical, dimension(max_nalleles, n_chromosomes), parameter, public cortisol_genotype_phenotype
Genotype x Phenotype matrix for cortisol.
real(srp), parameter, public adrenaline_gerror_cv
Genotype to phenotype gamma2gene Gaussian error parameter. This is really the coefficient of variatio...
real(srp), parameter, public growhorm_init
Genotype to phenotype gamma2gene initialisation value for growth hormone
real(srp), parameter, public cortisol_init
Genotype to phenotype gamma2gene initialisation value for cortisol
real(srp), parameter, public missing
Numerical code for missing and invalid real type values.
logical, dimension(max_nalleles, n_chromosomes), parameter, public testosterone_genotype_phenotype
Genotype x Phenotype matrix for testosterone.
logical, dimension(max_nalleles, n_chromosomes), parameter, public growhorm_genotype_phenotype
Genotype x Phenotype matrix for growth hormone.
real(srp), parameter, public estrogen_init
Genotype to phenotype gamma2gene initialisation value for estrogen
real(srp), parameter, public thyroid_init
Genotype to phenotype gamma2gene initialisation value for thyroid
real(srp), parameter, public testosterone_init
Genotype to phenotype gamma2gene initialisation value for testosterone
real(srp), parameter, public growhorm_gerror_cv
Genotype to phenotype gamma2gene Gaussian error parameter. This is really the coefficient of variatio...
real(srp), parameter, public cortisol_gerror_cv
Genotype to phenotype gamma2gene Gaussian error parameter. This is really the coefficient of variatio...
real(srp), parameter, public estrogen_gerror_cv
Genotype to phenotype gamma2gene Gaussian error parameter. This is really the coefficient of variatio...
Definition the genetic architecture of the agent.
Definition the hormonal architecture of the agent.
elemental real(srp) function thyroid_get_level(this)
Get the value of thyroid.
elemental real(srp) function growhorm_get_level(this)
Get the value of growth hormone.
subroutine hormones_init_genotype(this)
Initialise hormone levels based on the genome value. Two alleles are selected at random and input int...
elemental real(srp) function estrogen_get_level(this)
Get the value of estrogen.
elemental real(srp) function cortisol_get_level(this)
Get the value of cortisol.
elemental real(srp) function adrenaline_get_level(this)
Get the value of adrenaline.
elemental subroutine estrogen_set_level(this, value_set, update_history)
Set the value of estrogen.
elemental real(srp) function estrogen_baseline_get_level(this)
Get the value of estrogen baseline.
elemental real(srp) function testosteron_baseline_get_level(this)
Get the value of testosterone baseline.
elemental real(srp) function testosterone_get_level(this)
Get the value of testosterone.
elemental subroutine cortisol_set_level(this, value_set)
Set the value of cortisol.
elemental subroutine thyroid_set_level(this, value_set)
Set the value of thyroid.
elemental subroutine growhorm_set_level(this, value_set)
Set the value of growth hormone.
elemental subroutine testosterone_set_level(this, value_set, update_history)
Set the value of testosterone.
elemental real(srp) function hormones_reproductive_factor_calc(this)
Calculate the reproductive factor. Reproductive factor is defined as the current level of the_hormone...
elemental subroutine adrenaline_set_level(this, value_set)
Set the value of adrenaline.
elemental subroutine hormones_update_history(this)
Update the sex steroid hormones history stack from the current level.
elemental subroutine hormones_clean_history_stack(this)
Clean the history stack of hormones: testosterone and estrogen histories are set to MISSING.
This type describes parameters of the individual agent's genome The genome is an array of allocatable...
This type adds hormonal architecture extending the genome object.