The AHA Model  Revision: 12809
Reference implementation 04 (HEDG02_04)
m_indiv.f90
Go to the documentation of this file.
1 !> @file m_indiv.f90
2 !! Definition of the individual agent in the AHA Model.
3 !! @author Sergey Budaev <sergey.budaev@uib.no>
4 !! @author Jarl Giske <jarl.giske@uib.no>
5 !! @date 2016-2017
6 
7 !-------------------------------------------------------------------------------
8 ! $Id: m_indiv.f90 8329 2019-05-24 06:34:57Z sbu062 $
9 !-------------------------------------------------------------------------------
10 
11 !-------------------------------------------------------------------------------
12 !> @brief An umbrella module that collects all the components of the
13 !! individual agent.
14 !> @section the_individual_module THE_INDIVIDUAL module
15 !> Define the individual fish object and its properties and methods/functions
17 
18  use commondata ! Global definitions of the model objects
19  use the_genome ! ... and genome
20  use the_neurobio
21  use the_behaviour
22 
23  implicit none
24 
25  character (len=*), parameter, private :: modname = "(THE_INDIVIDUAL)"
26 
27  !> @brief This type describes parameters of the individual agent
28  !! @details `INDIVIDUAL_AGENT` extends the neurobiological architecture type
29  !! by adding the alive flag and explicit fitness.
30  !! Fitness can be either set or calculated or assessed
31  !! approximately. In models, it can be used in calculations
32  !! or used just for data output.
33  !! NOTE: procedure, final may not be implemented in all
34  !! compilers. E.g. gfortran issues this error prior to 4.9:
35  !! Error: Finalization at (1) is not yet implemented
36  !! @note `INDIVIDUAL_AGENT` is an umbrella object for the individual
37  type, public, extends(architecture_neuro) :: individual_agent
38  !> fitness is a fixed criterion for the evolution at the starting
39  !! 'pre-evolution' phase of the genetic algorithm.
40  !! @warning Note that fitness here is actually an "antifitness", the
41  !! higher its value, the **worse** fitting is the agent.
42  integer :: fitness
43  !> Individually specific mortality risk for various undefined or defined
44  !! reasons, such as immune status. Note that this risk is not linked to
45  !! the predation or habitat-specific factors. The latter two causes are
46  !! treated separately.
47  real(srp) :: ind_mortality
48  contains
49  private
50  !> Generate a random agent from the genotype.
51  !! See `the_individual::individual_init_random()`.
52  procedure, public :: init => individual_init_random
53  !> Generate a new *empty* agent.
54  !! See `the_individual::individual_create_zero()`.
55  procedure, public :: create_ind => individual_create_zero
56  !> Set the individual to be **dead**. This method overrides the
57  !! the_genome::individual_genome::dies() method, nullifying all
58  !! reproductive and neurobiological and behavioural objects.
59  !! However, this function does not deallocate the individual
60  !! agent object, this may be a separate destructor function.
61  !! The `dies` method is implemented at the following levels
62  !! of the agent object hierarchy (upper overrides the lower level):
63  !! - the_genome::individual_genome::dies();
64  !! - the_neurobio::appraisal::dies();
65  !! - the_neurobio::gos_global::dies();
66  !! - the_individual::individual_agent::dies().
67  !! .
68  !! See `the_individual::individual_agent_set_dead()`.
69  procedure, public :: dies => individual_agent_set_dead
70  !> Get the individually-specific mortality risk for the agent.
71  !! See `the_individual::individual_get_risk_mortality_individual()`.
72  procedure, public :: get_mortality => &
74  !> Calculate fitness for the pre-evolution phase of the genetic algorithm.
75  !! Pre-evolution is based on selection for a simple criterion without
76  !! explicit reproduction etc. The criterion for selection at this phase
77  !! is set by the integer the_individual::individual_agent::fitness
78  !! component. See `the_individual::individual_preevol_fitness_calc()`.
79  procedure, public :: fitness_calc => individual_preevol_fitness_calc
80  end type individual_agent
81 
82  ! we never use the actual proc names, only object-bound
83  private :: individual_init_random
84 
85 contains ! ........ implementation of procedures for this level ................
86 
87  !-----------------------------------------------------------------------------
88  !> @brief Generate a random agent from the genotype.
89  !! @details This subroutine is used to initialise the individual agent with
90  !! random data from the genotype. Its use is the most natural for
91  !! the initialisation of a population of agents. The init procedure
92  !! calls the init functions for the lower order layers of the class
93  !! hierarchy (genome, hormones, neurobio), and sets values.
94  subroutine individual_init_random(this, exclude_genome)
95  class(individual_agent), intent(inout) :: this
96  !> @param[in] exclude_genome is a logical flag to exclude initialisation
97  !! of random genome. If absent, assumed FALSE.
98  !! @note Only the first generation is initialised with random
99  !! genome (exclude_genome = FALSE), all subsequent
100  !! generations use pre-existing genomes from the
101  !! ancestors (exclude_genome = TRUE) randomised by
102  !! crossover in the genetic algorithm.
103  logical, optional, intent(in) :: exclude_genome
104 
105  !> We first initialise all the components of the agent down the class
106  !! hierarchy: from individual genome to the neurobiological architecture.
107  if (present(exclude_genome)) then
108  if (exclude_genome .eqv. .false.) then
109  call this%init_genome() ! initialise **random** genome
110  end if
111  else
112  call this%init_genome() ! initialise **random** genome
113  end if
114 
115  !> Clean the spatial location history stack of the agent.
116  call this%spatial_history_clean()
117 
118  call this%init_hormones() !> initialise hormone objects **from genome**
119  call this%init_condition() !> initialise condition **from genome**
120  call this%init_reproduction() !> initialise empty reproduction objects
121  call this%init_neurobio() !> initialise empty neuro objects
122 
123  !> Finally, we bring the agent to life by setting alive boolean flag
124  call this%lives()
125 
126  !> Set the individually specific mortality risk initially as a Gaussian
127  !! variable with mean commondata::individual_mortality_risk_def and CV
128  !! commondata::individual_mortality_risk_cv. There is also a restriction
129  !! that the risk of mortality should never be smaller than
130  !! commondata::zero.
131  this%ind_mortality = &
132  max( zero, rnorm( individual_mortality_risk_def, &
135 
136  !> Calculate the initial value of fitness for pre-evolution stage.
137  call this%fitness_calc()
138 
139  end subroutine individual_init_random
140 
141  !-----------------------------------------------------------------------------
142  !> @brief Generate a new empty agent.
143  !! @details This subroutine is used to create a new empty individual agent.
144  !! It should be used to make newborn agents that inherit traits
145  !! from their parents.
146  !! @warning This procedure is not used so far and is candidate for being
147  !! deprecated. Not clear if necessary.
148  subroutine individual_create_zero(this)
149  class(individual_agent), intent(inout) :: this
150 
151  !> Clean the spatial location history stack of the agent.
152  call this%spatial_history_clean()
153 
154  call this%create_genome() !> create **empty** genome
155 
156  !> Create empty hormone objects. No genome-based initialisation is done.
157  call this%hormone_history_clean()
158  call this%growhorm_set(missing)
159  call this%thyroid_set(missing)
160  call this%adrenaline_set(missing)
161  call this%cortisol_set(missing)
162  call this%testosterone_set(missing, .false.)
163  call this%estrogen_set(missing, .false.)
164 
165  !> Empty condition objects. Clean history stack.
166  ! @warning: energy_current, control_unselected and smr that are set
167  ! from the genome are not initialised so far due to absence of
168  ! "setter" methods: it is not clear if this procedure is ever
169  ! necessary at all (may be obsolete).
170  call this%body_history_clean()
171  call this%set_mass(missing, .false.)
172  call this%set_length(missing, .false.)
173 
174  !> Initialise reproduction objects and neurobiology to a zero state.
175  call this%init_reproduction()
176  call this%init_neurobio()
177 
178  !> Finally, we bring the agent to life by setting alive boolean flag
179  call this%lives()
180 
181  this%ind_mortality = missing
182 
183  !> Calculate the initial value of fitness for pre-evolution stage.
184  this%fitness = missing
185 
186  end subroutine individual_create_zero
187 
188  !-----------------------------------------------------------------------------
189  !> Set the individual to be **dead**. Note that this function does not
190  !! deallocate the individual agent object, this may be a separate destructor
191  !! function.
192  !!
193  !! The `dies` method is implemented at the following levels
194  !! of the agent object hierarchy (upper overrides the lower level):
195  !! - the_genome::individual_genome::dies();
196  !! - the_neurobio::appraisal::dies();
197  !! - the_neurobio::gos_global::dies();
198  !! - the_individual::individual_agent::dies().
199  !! .
200  !! @note This method overrides the the_genome::individual_genome::dies()
201  !! method, nullifying all reproductive and neurobiological and
202  !! behavioural objects.
203  elemental subroutine individual_agent_set_dead(this)
204  class(individual_agent), intent(inout) :: this
205 
206  call this%set_dead() !> - Set the agent "dead";
207  call this%init_reproduction() !> - emptify reproduction objects;
208  call this%init_neurobio() !> - emptify neurobiological objects.
209  !> .
210  end subroutine individual_agent_set_dead
211 
212  !-----------------------------------------------------------------------------
213  !> Finalisation procedure. Note that finalisation of objects may
214  !! not yet be implemented in the compiler. Therefore this subroutine
215  !! is not used so far, just a stub.
216  subroutine kill_destroy(this)
217  class(individual_agent), intent(inout) :: this
218 
219  call this%dies()
220 
221  ! call deallocation etc when final procedures are implemented
222  ! so far it is only a stub.
223 
224  end subroutine kill_destroy
225 
226  !-----------------------------------------------------------------------------
227  !> Get the individually-specific mortality risk for the agent.
228  elemental function individual_get_risk_mortality_individual(this) &
229  result(get_val)
230  class(individual_agent), intent(in) :: this
231  !> @return Individually-specific mortality risk for the agent.
232  real(srp) :: get_val
233 
234  get_val = this%ind_mortality
235 
237 
238  !-----------------------------------------------------------------------------
239  !> Calculate fitness for the pre-evolution phase of the genetic algorithm.
240  !! Pre-evolution is based on selection for a simple criterion without
241  !! explicit reproduction etc. The criterion for selection at this phase
242  !! is set by the integer the_individual::individual_agent::fitness component.
243  !! @warning Note that fitness here is actually an "antifitness", the
244  !! higher its value, the **worse** fitting is the agent.
245  elemental subroutine individual_preevol_fitness_calc(this)
246  class(individual_agent), intent(inout) :: this
247 
248  !> `INT_FITNESS_DEAD` is the fitness ascribed to the dead agent, it should
249  !! be a very large value greater than for any alive.
250  integer, parameter :: int_fitness_dead = ga_fitness_dead
251 
252  !> First check if the agent is dead. If so, give very high value that
253  !! is never selected.
254  if ( this%is_dead() ) then
255  this%fitness = int_fitness_dead
256  return
257  end if
258 
259  !> Now, the reverse of fitness can be calculated by various methods.
260  !! Specific calculation functions are implemented within this function.
261  !! So far the following routines were implemented:
262  !! - ::fitness_energy_reprfact_mass
263  !! - ::fitness_energy_reprfact())
264  !! - ::fitness_mass_incr_abs()
265  !! - ::fitness_birth_mass_ratio()
266  !! - ::fitness_stomach_mass_ratio()
267  !! - ::fitness_stomach_mass_abs()
268  !! - ::fitness_food_mass_sum()
269  !! - ::fitness_mass_incr_ratio()
270  !! - ::fitness_reprod_factor()
271  !! .
272  this%fitness = fitness_energy_reprfact_mass()
273 
274  contains !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
275 
276  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
277  !> Fitness is calculated as the ratio of the birth mass to the current
278  !! mass. This value is weighted by the multiplier (1000) to get a fairly
279  !! large integer (so decimals are unimportant) and also weighted by the
280  !! number of food items eaten and the number of offspring produced.
281  !! @note This procedure is internal to
282  !! the_individual::individual_agent::fitness_calc().
283  elemental function fitness_birth_mass_ratio() result (fitness_out)
284  !> @return Fitness value.
285  integer :: fitness_out
286  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
287  !! scaling for the initial the_individual::individual_agent::fitness.
288  real(srp), parameter :: int_multiplier_fitness = 1000.0_srp
289  !> `INT_WEIGHT_FEEDING` is an integer weight given to the any non-zero
290  !! successful feeding.
291  integer, parameter :: int_weight_feeding = 80
292  !> `INT_WEIGHT_OFFSPRING` is an integer weight given to the any
293  !! non-zero successful reproductions
294  integer, parameter :: int_weight_offspring = 400
295  !> Initial fitness is the ratio of the birth mass to the current mass
296  !! weighted by the `INT_MULTIPLIER_FITNESS`.
297  fitness_out = nint( int_multiplier_fitness * &
298  ( this%get_mass_birth() / this%get_mass() ) )
299  !> If the agent successfully caught and eaten any number of food items,
300  !! its fitness is divided by the number of food items eaten weighted by
301  !! the `INT_WEIGHT_FEEDING` parameter.
302  if (this%n_eaten_indicator>0) fitness_out = fitness_out / &
303  ( this%n_eaten_indicator * int_weight_feeding )
304  !> If the agent has successfully done any reproduction, its fitness is
305  !! divided by the number of offspring weighted by `INT_WEIGHT_OFFSPRING`.
306  ! @note wxMaxima quick plot:
307  ! `wxplot3d ( (40/m*1000)/(o*10) , [m,30,100], [o,1,10]);`
308  ! Gnuplot code:
309  ! @verbatim
310  ! set xrange [30:100]
311  ! set yrange [1:100]
312  ! set xlabel "mass"
313  ! set xlabel "offspring"
314  ! splot (40/x*1000)/(y*10) with lines
315  ! @endverbatim
316  if (this%get_offspring()>0) fitness_out = fitness_out / &
317  ( this%get_offspring() * int_weight_offspring )
318  end function fitness_birth_mass_ratio
319 
320  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
321  !> Fitness is calculated as the ratio of the body mass to the stomach
322  !! content
323  !! @note This procedure is internal to
324  !! the_individual::individual_agent::fitness_calc().
325  elemental function fitness_stomach_mass_ratio() result (fitness_out)
326  !> @return Fitness value.
327  integer :: fitness_out
328 
329  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
330  !! scaling for the initial the_individual::individual_agent::fitness.
331  real(srp), parameter :: int_multiplier_fitness = 10.0_srp
332 
333  if (is_near_zero(this%get_stom_content(), tolerance_low_def_srp)) then
334  fitness_out = 2000
335  return
336  end if
337 
338  fitness_out = nint( int_multiplier_fitness * &
339  this%get_mass() / this%get_stom_content() )
340 
341  end function fitness_stomach_mass_ratio
342 
343  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
344  !> Fitness is calculated as the absolute stomach content.
345  !! @note This procedure is internal to
346  !! the_individual::individual_agent::fitness_calc().
347  elemental function fitness_stomach_mass_abs() result (fitness_out)
348  !> @return Fitness value.
349  integer :: fitness_out
350 
351  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
352  !! scaling for the initial the_individual::individual_agent::fitness.
353  real(srp), parameter :: int_multiplier_fitness = 100.0_srp
354 
355  if (is_near_zero(this%get_stom_content(), tolerance_low_def_srp)) then
356  fitness_out = 2000
357  return
358  end if
359 
360  fitness_out = nint( int_multiplier_fitness * &
361  1.0_srp / this%get_stom_content() )
362 
363  end function fitness_stomach_mass_abs
364 
365  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
366  !> Fitness is calculated as the cumulative mass of all food objects eaten.
367  !! @note This procedure is internal to
368  !! the_individual::individual_agent::fitness_calc().
369  elemental function fitness_food_mass_sum() result (fitness_out)
370  !> @return Fitness value.
371  integer :: fitness_out
372 
373  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
374  !! scaling for the initial the_individual::individual_agent::fitness.
375  real(srp), parameter :: int_multiplier_fitness = 10000.0_srp
376 
377  if (is_near_zero(this%mass_eaten_indicator, tolerance_low_def_srp)) then
378  fitness_out = 1000
379  return
380  end if
381 
382  fitness_out = nint( int_multiplier_fitness * &
383  1.0_srp / this%mass_eaten_indicator )
384 
385  end function fitness_food_mass_sum
386 
387  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
388  !> Fitness is calculated as the mass increment in units of birth mass.
389  !! @note This procedure is internal to
390  !! the_individual::individual_agent::fitness_calc().
391  elemental function fitness_mass_incr_ratio() result (fitness_out)
392  !> @return Fitness value.
393  integer :: fitness_out
394 
395  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
396  !! scaling for the initial the_individual::individual_agent::fitness.
397  real(srp), parameter :: int_multiplier_fitness = 10.0_srp
398 
399  if (this%get_mass() <= this%get_mass_birth()) then
400  fitness_out = 1000 - this%n_eaten_indicator
401  return
402  end if
403 
404  fitness_out = nint( int_multiplier_fitness * &
405  this%get_mass_birth() / (this%get_mass() - this%get_mass_birth()) )
406 
407  end function fitness_mass_incr_ratio
408 
409  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
410  !> Fitness is calculated as absolute mass increment from the birth mass.
411  !! @note This procedure is internal to
412  !! the_individual::individual_agent::fitness_calc().
413  elemental function fitness_mass_incr_abs() result (fitness_out)
414  !> @return Fitness value.
415  integer :: fitness_out
416 
417  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
418  !! scaling for the initial the_individual::individual_agent::fitness.
419  real(srp), parameter :: int_multiplier_fitness = 10.0_srp
420 
421  if (this%get_mass() <= this%get_mass_birth()) then
422  fitness_out = 1000 - this%n_eaten_indicator
423  return
424  end if
425 
426  fitness_out = nint( int_multiplier_fitness * &
427  1.0_srp / (this%get_mass() - this%get_mass_birth()) )
428 
429  end function fitness_mass_incr_abs
430 
431  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
432  !> Fitness as the reproductive factor.
433  !! @note This procedure is internal to
434  !! the_individual::individual_agent::fitness_calc().
435  elemental function fitness_reprod_factor() result (fitness_out)
436  !> @return Fitness value.
437  integer :: fitness_out
438 
439  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
440  !! scaling for the initial the_individual::individual_agent::fitness.
441  real, parameter :: int_multiplier_fitness = 100.0_srp
442 
443  if (this%get_mass() <= this%get_mass_birth()) then
444  fitness_out = 1000 - within( this%n_eaten_indicator, 0, 100 )
445  return
446  end if
447 
448  fitness_out =nint( int_multiplier_fitness * &
449  1.0_srp / this%reproductive_factor() )
450 
451  end function fitness_reprod_factor
452 
453  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
454  !> Fitness as a weighted combination of energy and reproductive factor.
455  !! @note This procedure is internal to
456  !! the_individual::individual_agent::fitness_calc().
457  elemental function fitness_energy_reprfact() result (fitness_out)
458  !> @return Fitness value.
459  integer :: fitness_out
460  real(hrp) :: fitness_real
461 
462  !> `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
463  !! scaling for the initial the_individual::individual_agent::fitness.
464  real(hrp), parameter :: int_multiplier_fitness = 100000.0_hrp
465 
466  !> `FITNESS_WEIGHT_ENERGY` is the weighting factor for energy,
467  !! respective weighting for the reproductive factor is then
468  !! `1 - FITNESS_WEIGHT_ENERGY`.
469  real(srp), parameter :: fitness_weight_energy = 0.9_srp
470 
471  fitness_real = real( this%get_energy() * fitness_weight_energy + &
472  this%reproductive_factor() * &
473  (1.0_srp - fitness_weight_energy ), hrp )
474 
475  fitness_out = &
476  nint( int_multiplier_fitness / fitness_real )
477 
478  !> Fitness of an alive agent cannot be smaller than in a dead agent.
479  fitness_out = within( fitness_out, 0, int_fitness_dead )
480 
481  end function fitness_energy_reprfact
482 
483  !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
484  !> Fitness as a weighted combination of mass increment, energy, and
485  !! reproductive factor.
486  !! @note This procedure is internal to
487  !! the_individual::individual_agent::fitness_calc().
488  elemental function fitness_energy_reprfact_mass() result (fitness_out)
489  !> @return Fitness value.
490  integer :: fitness_out
491  real(hrp) :: fitness_real
492 
493  !> - `INT_MULTIPLIER_FITNESS` is a multiplier to set the appropriate
494  !! scaling for the initial
495  !! the_individual::individual_agent::fitness.
496  real(hrp), parameter :: int_multiplier_fitness = 10000.0_hrp
497 
498  !> - `FITNESS_WEIGHT_MASINC` is the weighting factor for
499  !! relative weight increment.
500  !! - `FITNESS_WEIGHT_ENERGY` is the weighting factor for energy,
501  !! .
502  !! @note their sum must not exceed 1.0
503  real(srp), parameter :: fitness_weight_masinc = 0.30_srp, &
504  fitness_weight_energy = 0.60_srp
505 
506 
507  associate( bm => this%body_mass_birth, &
508  m=>this%body_mass, &
509  e => this%energy_current, &
510  r => this%reproductive_factor(), &
511  z1 => int_multiplier_fitness, &
512  zm => fitness_weight_masinc, &
513  ze => fitness_weight_energy )
514 
515  !> If mass increment is negative, fitness is given a big value
516  !! **1000,000**,
517  !! @note Note that this should be smaller than the fitness of a
518  !! dead agent defined in commondata::ga_fitness_dead .
519  if (m <= bm) then
520  fitness_out = z1 * 100.0_hrp - z1 * (m - bm) / bm
521  return
522  end if
523 
524  !> Then fitness function is calculated as a weighted sum of
525  !! relative weight increment W1 (M-Mb)/Mb + W2 E + W3 Rf
526  fitness_real = &
527  real( ((m - bm) / bm) * zm + e * ze + r * (1.0_srp - ze -zm )&
528  , hrp )
529 
530  fitness_out = &
531  nint( z1 / fitness_real )
532 
533  end associate
534 
535  !> Fitness of an alive agent is always limited to lay within the
536  !! range from 0 to fitness of a dead agent.
537  fitness_out = within( fitness_out, 0, int_fitness_dead )
538 
539  end function fitness_energy_reprfact_mass
540 
541  end subroutine individual_preevol_fitness_calc
542 
543 end module the_individual
Checks if a real number is near 0.0. Thus function can be used for comparing two real values like thi...
Definition: m_common.f90:5385
Force a value within the range set by the vmin and vmax dummy parameter values.
Definition: m_common.f90:5350
elemental integer function fitness_mass_incr_ratio()
Fitness is calculated as the mass increment in units of birth mass.
Definition: m_indiv.f90:392
elemental integer function fitness_stomach_mass_ratio()
Fitness is calculated as the ratio of the body mass to the stomach content.
Definition: m_indiv.f90:326
elemental integer function fitness_birth_mass_ratio()
Fitness is calculated as the ratio of the birth mass to the current mass. This value is weighted by t...
Definition: m_indiv.f90:284
elemental integer function fitness_food_mass_sum()
Fitness is calculated as the cumulative mass of all food objects eaten.
Definition: m_indiv.f90:370
elemental integer function fitness_energy_reprfact()
Fitness as a weighted combination of energy and reproductive factor.
Definition: m_indiv.f90:458
elemental integer function fitness_stomach_mass_abs()
Fitness is calculated as the absolute stomach content.
Definition: m_indiv.f90:348
elemental integer function fitness_reprod_factor()
Fitness as the reproductive factor.
Definition: m_indiv.f90:436
elemental integer function fitness_energy_reprfact_mass()
Fitness as a weighted combination of mass increment, energy, and reproductive factor.
Definition: m_indiv.f90:489
elemental integer function fitness_mass_incr_abs()
Fitness is calculated as absolute mass increment from the birth mass.
Definition: m_indiv.f90:414
COMMONDATA – definitions of global constants and procedures.
Definition: m_common.f90:1497
real(srp), parameter, public individual_mortality_risk_cv
The coefficient of variation for Gaussian stochastic individually-specific mortality risk of the agen...
Definition: m_common.f90:2345
real(srp), parameter, public tolerance_low_def_srp
Default value of low tolerance (high precision). This is the standard commondata::srp precision....
Definition: m_common.f90:1671
integer, parameter, public ga_fitness_dead
Fitness value ascribed to dead agent in pre-evol. See the_individual::individual_agent::fitness_calc(...
Definition: m_common.f90:5194
integer, parameter, public srp
Definition of the standard real type precision (SRP).
Definition: m_common.f90:1551
elemental real(srp) function cv2variance(cv, mean)
Calculate the variance from the coefficient of variation.
Definition: m_common.f90:6956
real(srp), parameter, public missing
Numerical code for missing and invalid real type values.
Definition: m_common.f90:1699
integer, parameter, public hrp
Definition of the high real precision (HRP). This real type kind is used in pieces where a higher lev...
Definition: m_common.f90:1556
real(srp), parameter, public zero
Some parameters should never be zero or below. In such cases they could be set to some smallest disti...
Definition: m_common.f90:1644
real(srp), parameter, public individual_mortality_risk_def
Default individually-specific mortality risk. It can increase or decrease depending on various factor...
Definition: m_common.f90:2341
logical, parameter, public false
Definition: m_common.f90:1632
Definition of high level behavioural architecture.
Definition: m_behav.f90:17
character(len= *), parameter, private modname
Definition: m_behav.f90:26
Definition the genetic architecture of the agent.
Definition: m_genome.f90:16
An umbrella module that collects all the components of the individual agent.
Definition: m_indiv.f90:16
subroutine individual_create_zero(this)
Generate a new empty agent.
Definition: m_indiv.f90:149
subroutine, private individual_init_random(this, exclude_genome)
Generate a random agent from the genotype.
Definition: m_indiv.f90:95
elemental real(srp) function individual_get_risk_mortality_individual(this)
Get the individually-specific mortality risk for the agent.
Definition: m_indiv.f90:230
elemental subroutine individual_agent_set_dead(this)
Set the individual to be dead. Note that this function does not deallocate the individual agent objec...
Definition: m_indiv.f90:204
subroutine kill_destroy(this)
Finalisation procedure. Note that finalisation of objects may not yet be implemented in the compiler....
Definition: m_indiv.f90:217
elemental subroutine individual_preevol_fitness_calc(this)
Calculate fitness for the pre-evolution phase of the genetic algorithm. Pre-evolution is based on sel...
Definition: m_indiv.f90:246
Definition of the decision making and behavioural the architecture.
Definition: m_neuro.f90:17
This type is an "umbrella" for all the lower-level classes.
Definition: m_behav.f90:627
This type describes parameters of the individual agent.
Definition: m_indiv.f90:37