alpha_decay_rate_module Module Reference

This module contains subroutines to add alpha decays. More...

Functions/Subroutines

subroutine, public init_alpha_decay_rates ()
 Initialize alpha decay rates. More...
 
subroutine, private create_verbose_output_file (alpha_decay_rates, nalpha)
 Creates a file with the decays. More...
 
character(150) function, private reaction_string (reac)
 Return a string to represent a given reaction. More...
 
subroutine, public merge_alpha_decays (rrate_array, rrate_length)
 Merge alpha decays into the larger rate array. More...
 
subroutine, private unify_rate_array (rrate_array, alpha_dec_rate_array, merged_array, nrrate, nalpha, ntot)
 Merge reaclib rates and alpha decays into one array. More...
 
subroutine, private count_reactions (n)
 Count the amount of possible alpha decays. More...
 
subroutine, private read_reactions (alpha_rate_array, n)
 Read the alpha decays. More...
 

Variables

type(reactionrate_type), dimension(:), allocatable, private alpha_dec_rate
 Array storing the reaction rates. More...
 
integer, private nalpha_dec
 Number of alpha decays. More...
 
logical, private include_alpha_decays
 Flag whether alpha decays should be included or not. More...
 
character(len=4), dimension(:), allocatable, private src_ignore
 
integer, private src_ignore_length
 
character(15), parameter, private fileformat ="(A5,3X,1pE11.5)"
 

Detailed Description

This module contains subroutines to add alpha decays.

The alpha decays given by an external file. The default file in WinNet uses a file that has been calculated using the Viola-Seaborg formula. Hereby, a parameterization of the alpha-decay half life of Dong & Ren 2005 was used. Via the Viola-Seaborg formula it is possible to calculate alpha-decay half lifes by knowing only the Qvalue of the alpha decay. Within this model, reaclib rates can get replaced or additional alpha-decay rates can be added.

Author
Moritz Reichert
Date
10.03.23

Function/Subroutine Documentation

◆ count_reactions()

subroutine, private alpha_decay_rate_module::count_reactions ( integer, intent(out)  n)
private

Count the amount of possible alpha decays.

Check which isotopes between alpha_decay_zmin and alpha_decay_zmax should have an alpha decay.

Author
M. Reichert
Date
10.03.23
Parameters
[out]nNumber of possible alpha decays

Definition at line 372 of file alpha_decay_rate_module.f90.

Here is the call graph for this function:

◆ create_verbose_output_file()

subroutine, private alpha_decay_rate_module::create_verbose_output_file ( type(reactionrate_type), dimension(:), intent(in)  alpha_decay_rates,
integer, intent(in)  nalpha 
)
private

Creates a file with the decays.

This file is created for verbose_levels greater than 2. An example could look like

  Alpha decay rates that additionally have been added
 po181 =>   he4 + pb177  T_1/2: 1.49E-09 s;  a0: 2.00E+01  Qa: 1.02E+01
 po182 =>   he4 + pb178  T_1/2: 3.89E-09 s;  a0: 1.90E+01  Qa: 9.56E+00
 ..
 
Author
M. Reichert
Date
10.03.23
Parameters
[in]alpha_decay_ratesArray storing the reaction rates
[in]nalphaNumber of alpha decays

Definition at line 99 of file alpha_decay_rate_module.f90.

Here is the call graph for this function:

◆ init_alpha_decay_rates()

subroutine, public alpha_decay_rate_module::init_alpha_decay_rates

Initialize alpha decay rates.

This subroutine counts and reads the rates into the array beta_decays.

Author
M. Reichert
Date
25.01.21

Definition at line 53 of file alpha_decay_rate_module.f90.

Here is the call graph for this function:

◆ merge_alpha_decays()

subroutine, public alpha_decay_rate_module::merge_alpha_decays ( type(reactionrate_type), dimension(:), intent(inout), allocatable  rrate_array,
integer, intent(inout)  rrate_length 
)

Merge alpha decays into the larger rate array.

The return value of this routine will be a larger rate array and a new length

Author
M. Reichert
Date
25.01.21
Parameters
[in,out]rrate_arrayLarge rate array, containing all reactions
[in,out]rrate_lengthlength of rrate_array

Definition at line 203 of file alpha_decay_rate_module.f90.

Here is the call graph for this function:

◆ reaction_string()

character(150) function, private alpha_decay_rate_module::reaction_string ( type(reactionrate_type), intent(in)  reac)
private

Return a string to represent a given reaction.

This routine is useful for error messages in case a rate is not working as intented and is used for verbose output only.

Example

For rrate(i) being the alpha-decay of po181

a = reaction_string(rrate(i))

will return a="po181 => he4 + pb177 T_1/2: 1.49E-09 s; a0: 2.00E+01 Qa: 1.02E+01".

Author
Moritz Reichert

Definition at line 137 of file alpha_decay_rate_module.f90.

Here is the call graph for this function:

◆ read_reactions()

subroutine, private alpha_decay_rate_module::read_reactions ( type(reactionrate_type), dimension(:), intent(out), allocatable  alpha_rate_array,
integer, intent(in)  n 
)
private

Read the alpha decays.

All isotopes between alpha_decay_zmin and alpha_decay_zmax should have an alpha decay. This subroutine reads the file alpha_decay_file and creates the alpha decay rate array. An example of the file could look like:

 sb103   7.84058e+10
 sb104   4.29655e+05
 sb105   1.37511e+10
 sb106   1.03780e+17
 te103   8.30806e+02
 ...
 
Author
M. Reichert
Date
14.03.23
Parameters
[in]nLength of the rates

Definition at line 438 of file alpha_decay_rate_module.f90.

Here is the call graph for this function:

◆ unify_rate_array()

subroutine, private alpha_decay_rate_module::unify_rate_array ( type(reactionrate_type), dimension(nrrate), intent(inout)  rrate_array,
type(reactionrate_type), dimension(nalpha), intent(inout)  alpha_dec_rate_array,
type(reactionrate_type), dimension(:), intent(out), allocatable  merged_array,
integer, intent(in)  nrrate,
integer, intent(in)  nalpha,
integer, intent(out)  ntot 
)
private

Merge reaclib rates and alpha decays into one array.

Hereby the rates can be taken as addition or also replace the reaclib rates.

Author
M. Reichert
Date
25.01.21
Parameters
[in,out]rrate_arrayReaclib rates
[in,out]alpha_dec_rate_arrayAlpha decay rates
[out]merged_arrayCombined, merged array
[in]nrrateLength of rrate_array
[in]nalphaLength of alpha_dec_rate_array
[out]ntotLength of merged_array

Definition at line 272 of file alpha_decay_rate_module.f90.

Here is the call graph for this function:

Variable Documentation

◆ alpha_dec_rate

type(reactionrate_type), dimension(:), allocatable, private alpha_decay_rate_module::alpha_dec_rate
private

Array storing the reaction rates.

Definition at line 23 of file alpha_decay_rate_module.f90.

◆ fileformat

character(15), parameter, private alpha_decay_rate_module::fileformat ="(A5,3X,1pE11.5)"
private

Definition at line 31 of file alpha_decay_rate_module.f90.

◆ include_alpha_decays

logical, private alpha_decay_rate_module::include_alpha_decays
private

Flag whether alpha decays should be included or not.

Definition at line 25 of file alpha_decay_rate_module.f90.

◆ nalpha_dec

integer, private alpha_decay_rate_module::nalpha_dec
private

Number of alpha decays.

Definition at line 24 of file alpha_decay_rate_module.f90.

◆ src_ignore

character(len=4), dimension(:), allocatable, private alpha_decay_rate_module::src_ignore
private

Definition at line 28 of file alpha_decay_rate_module.f90.

◆ src_ignore_length

integer, private alpha_decay_rate_module::src_ignore_length
private

Definition at line 29 of file alpha_decay_rate_module.f90.