beta_decay_rate_module Module Reference

This module contains subroutines to include external beta decays. More...

Data Types

type  beta_pn_type
 Type for storing the beta decay data with Pn channels. More...
 

Functions/Subroutines

subroutine, public init_ext_beta_rates ()
 Initialize external beta decay rates. More...
 
subroutine, public merge_beta_decays (rrate_array, rrate_length)
 Merge external beta decays into the larger rate array. More...
 
subroutine ignore_reactions (rrate_array, rrate_length)
 Subroutine to remove weak rates from the beta decay rate array. More...
 
subroutine, private count_reactions (sourcefile)
 Count the amount of external beta decays. More...
 
subroutine, private remove_weak_rates (rrate_array, length_rate_array)
 Remove beta decays from the rate array. More...
 
subroutine create_rrate_array (rrate_beta, rrate_beta_length)
 Convert beta_pn_type to reactionrate_type. More...
 
subroutine, private read_beta_decays (sourcefile)
 Reads beta decays from a separate file. More...
 

Variables

type(beta_pn_type), dimension(:), allocatable, private beta_pn
 Array storing the reaction rates. More...
 
integer, private nbeta_pn
 Number of beta decays. More...
 
type(reactionrate_type), dimension(:), allocatable, public beta_decays
 array containing external beta decays More...
 
integer, public nbeta
 total number of external beta decays More...
 
logical, private ext_decays
 Flag if external beta decays are used. More...
 
character(len=4), dimension(:), allocatable, private src_ignore
 
integer, private src_ignore_length
 

Detailed Description

This module contains subroutines to include external beta decays.

The external beta decays are given in a different format compared to reaclib. They tabulate the halflife of a nucleus together with Pn probabilities (beta delayed neutron emission). The file includes up to P10n (10 beta delayed neutrons), which can not be calculated in the default reaclib format (There P2n is max. in chapter 3, or P3n in the new reaclib format in chapter 11). An example entry looks like:

...
 ca73    6.200000e-04
 0.0000   0.0100   0.0200   0.0100   0.2800   0.0800   0.4000   0.0700   0.1100   0.0100   0.0100
 ...

Some tables, however, tabulate until P10n, for example, Moeller et al. 2019

Author
Moritz Reichert
Date
24.01.21

Function/Subroutine Documentation

◆ count_reactions()

subroutine, private beta_decay_rate_module::count_reactions ( integer, intent(in)  sourcefile)
private

Count the amount of external beta decays.

After this routine has been called, nbeta will contain the amount of external beta decay rates.

Author
M. Reichert
Date
25.01.21
Parameters
[in]sourcefilefile id to read beta decays

Definition at line 258 of file beta_decay_rate_module.f90.

Here is the call graph for this function:

◆ create_rrate_array()

subroutine beta_decay_rate_module::create_rrate_array ( type(reactionrate_type), dimension(:), intent(out), allocatable  rrate_beta,
integer, intent(out)  rrate_beta_length 
)

Convert beta_pn_type to reactionrate_type.

Creates a reactionrate_type array from the beta_pn_type array. Rates with Pn of 0 will not result in an extra rate.

Author
Moritz Reichert
Date
28.02.2023
Parameters
[out]rrate_betaReaction rates in reactionrate_type format
[out]rrate_beta_lengthLength of rrate_beta

Definition at line 420 of file beta_decay_rate_module.f90.

Here is the call graph for this function:

◆ ignore_reactions()

subroutine beta_decay_rate_module::ignore_reactions ( type(reactionrate_type), dimension(:), allocatable  rrate_array,
integer, intent(in)  rrate_length 
)

Subroutine to remove weak rates from the beta decay rate array.

This subroutine ensures that parameter_class::beta_decay_src_ignore is working. It removes all reactions from the array beta_pn that are listed in parameter_class::beta_decay_src_ignore.

Author
M. Reichert

Definition at line 170 of file beta_decay_rate_module.f90.

Here is the call graph for this function:

◆ init_ext_beta_rates()

subroutine, public beta_decay_rate_module::init_ext_beta_rates

Initialize external beta decay rates.

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

Author
M. Reichert
Date
25.01.21

Definition at line 73 of file beta_decay_rate_module.f90.

Here is the call graph for this function:

◆ merge_beta_decays()

subroutine, public beta_decay_rate_module::merge_beta_decays ( type(reactionrate_type), dimension(:), intent(inout), allocatable  rrate_array,
integer, intent(inout)  rrate_length 
)

Merge external beta 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 123 of file beta_decay_rate_module.f90.

Here is the call graph for this function:

◆ read_beta_decays()

subroutine, private beta_decay_rate_module::read_beta_decays ( integer, intent(in)  sourcefile)
private

Reads beta decays from a separate file.

This routine is thought to read beta decays in a different format than reaclib. An advantage of this is that we can implement more than the P2n channel. The new format allows beta delayed neutrons up to 10 neutrons that become more important close to stability. An example input file could look like:

 ca73    6.200000e-04
 0.0000   0.0100   0.0200   0.0100   0.2800   0.0800   0.4000   0.0700   0.1100   0.0100   0.0100 

which will replace the decay of ca73 with a half life of 6.2e-04s and the listet Pn probabilities. Additionally, the average Q-Value and the average neutrino Q-Value can be given optionally in the first line. In this case an entry could look like: ca73 1.268195e-03 2.443690e+01 9.034672e+00 0.0030 0.0000 0.0010 0.0000 0.0020 0.9950 0.0000 0.0000 0.0000 0.0000 0.0000 } where the first line gives the name of the parent nucleus, the half life, the average Q-Value and the average neutrino Q-Value.

Author
Moritz Reichert

Edited:

  • 20.05.19
Parameters
[in]sourcefilefile id to read beta decays

Definition at line 532 of file beta_decay_rate_module.f90.

Here is the call graph for this function:

◆ remove_weak_rates()

subroutine, private beta_decay_rate_module::remove_weak_rates ( type(reactionrate_type), dimension(:), allocatable  rrate_array,
integer  length_rate_array 
)
private

Remove beta decays from the rate array.

This routine removes beta decays in case that they exist in another format ( see parameter_class::use_beta_decay_file ).

Author
Moritz Reichert

Definition at line 313 of file beta_decay_rate_module.f90.

Here is the call graph for this function:

Variable Documentation

◆ beta_decays

type(reactionrate_type), dimension(:), allocatable, public beta_decay_rate_module::beta_decays

array containing external beta decays

Definition at line 49 of file beta_decay_rate_module.f90.

◆ beta_pn

type(beta_pn_type), dimension(:), allocatable, private beta_decay_rate_module::beta_pn
private

Array storing the reaction rates.

Definition at line 44 of file beta_decay_rate_module.f90.

◆ ext_decays

logical, private beta_decay_rate_module::ext_decays
private

Flag if external beta decays are used.

Definition at line 51 of file beta_decay_rate_module.f90.

◆ nbeta

integer, public beta_decay_rate_module::nbeta

total number of external beta decays

Definition at line 50 of file beta_decay_rate_module.f90.

◆ nbeta_pn

integer, private beta_decay_rate_module::nbeta_pn
private

Number of beta decays.

Definition at line 45 of file beta_decay_rate_module.f90.

◆ src_ignore

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

Definition at line 53 of file beta_decay_rate_module.f90.

◆ src_ignore_length

integer, private beta_decay_rate_module::src_ignore_length
private

Definition at line 54 of file beta_decay_rate_module.f90.