![]() |
This module contains everything for the theoretical weak rates that are added into the rate array. More...
Data Types | |
type | weakrate_type |
data fields for weak rates given in Langanke&Martinez-Pinedo 2001 More... | |
Functions/Subroutines | |
subroutine, public | init_theoretical_weak_rates () |
Initialize theoretical weak rates. More... | |
subroutine | initialize_cubic_interp |
Calculate derivatives of weak reactions. More... | |
subroutine, private | write_reac_verbose_out () |
Write the amount of individual reactions to the out. More... | |
subroutine | merge_theoretical_weak_rates (rrate_array, rrate_length) |
Merge theoretical weak rates into larger rate array. More... | |
subroutine | read_binary_weak_reaction_data (path) |
Read the theoretical weak rates from a unformatted binary file. More... | |
subroutine, public | output_binary_weak_reaction_data (path) |
Save the theoretical weak rates to a unformatted binary file. More... | |
subroutine, public | calculate_twr_rate (rrate, temp, rho, Ye, rat_calc, nuloss) |
Calculate the theoretical weak rate. More... | |
subroutine, private | read_theoretical_weak_rates () |
Reads and counts the amount of theoretical weak rates. More... | |
subroutine, private | readweak_logft (sourcefile, cnt) |
Reads weak reaction rates in log format. More... | |
subroutine, public | weak_index (ltemp, rho, Ye, wrt) |
Set wk_index for given T, \(\rho\) and \(Y_e\) The index is written into module variable wk_index. More... | |
real(r_kind) function, public | weak_inter (rate, ltemp, lrho, wrt, force_cubic) |
Interpolation interface. More... | |
subroutine, public | reload_exp_weak_rates () |
This routine is used to replace theoretical weak rates with experimental ones it's called when T9 < parameter_class::temp_reload_exp_weak_rates (default 1d-2) More... | |
subroutine | output_n_p |
subroutine, private | sort (cnt) |
Sorts the entries in weak_rate in increasing order of the decaying nuclide. More... | |
Variables | |
type(weakrate_type), dimension(:), allocatable, public | weak_rate |
array containing the weak reactions More... | |
logical, public | weak |
switch for weak rates More... | |
integer, private | nweak |
number of weak reactions taken from Langanke&Martinez-pinedo More... | |
type(reactionrate_type), dimension(:), allocatable, private | rrate_weak |
Reaction rate representative in global_class::rrate. More... | |
integer, dimension(2, 2), private | wk_index |
Multi-index for the weak rates. More... | |
real(r_kind), public | mue |
chemical potential of the electron More... | |
integer, parameter, private | mue_ident =1 |
Identifier for chemical potential. More... | |
integer, parameter, private | beta_ident =2 |
Identifier for beta rate. More... | |
integer, parameter, private | ft_ident =3 |
Identifier for forward rate. More... | |
integer, parameter, private | nu_loss_ident =4 |
Identifier for the neutrino loss. More... | |
integer, private | n_ec |
integer, private | n_o |
individual reaction types More... | |
character(len= *), parameter, private | weak_binary_name ='weak_rates.windat' |
Filename of binary file to save weak rates. More... | |
This module contains everything for the theoretical weak rates that are added into the rate array.
The theoretical weak rates in the current form were implemented by M. Ugliano. They contain electron captures and beta decays. Rates that are contained in the weak rate file and in reaclib get replaced by the ones contained in weak rates for small temperatures (see parameter_class::temp_reload_exp_weak_rates).
subroutine, public tw_rate_module::calculate_twr_rate | ( | type(reactionrate_type), intent(inout) | rrate, |
real(r_kind), intent(in) | temp, | ||
real(r_kind), intent(in) | rho, | ||
real(r_kind), intent(in) | Ye, | ||
real(r_kind), intent(out) | rat_calc, | ||
logical, optional | nuloss | ||
) |
Calculate the theoretical weak rate.
For parameter_class::iwformat equal to 1, this routine interpolates the rate on the temperature - density grid. Otherwise, for parameter_class::iwformat equal 2, this subroutine calculates the effective phase space integral to translate the tabulated log <ft> values to actual rates. For this also the chemical potential has to be known.
Edited:
[in,out] | rrate | rate instance |
[in] | temp | Temperature [GK] |
[in] | rho | Density [g/ccm] |
[in] | ye | Electron fraction |
[out] | rat_calc | rate value |
nuloss | calculate neutrino loss? |
Definition at line 484 of file tw_rate_module.f90.
subroutine, public tw_rate_module::init_theoretical_weak_rates |
Initialize theoretical weak rates.
This subroutine set the flag weak and calls the necessary subroutines to read the rates.
Definition at line 90 of file tw_rate_module.f90.
subroutine tw_rate_module::initialize_cubic_interp |
Calculate derivatives of weak reactions.
These derivatives are necessary for the bicubic interpolation.
Definition at line 121 of file tw_rate_module.f90.
subroutine tw_rate_module::merge_theoretical_weak_rates | ( | type(reactionrate_type), dimension(:), intent(inout), allocatable | rrate_array, |
integer, intent(inout) | rrate_length | ||
) |
Merge theoretical weak rates into larger rate array.
This merge is partly done instantaneously and partly when the temperature reaches parameter_class::reload_exp_weak_rates.
[in,out] | rrate_array | Large rate array, containing all reactions |
[in,out] | rrate_length | length of rrate_array |
Definition at line 208 of file tw_rate_module.f90.
subroutine, public tw_rate_module::output_binary_weak_reaction_data | ( | character(len=*), intent(in) | path | ) |
Save the theoretical weak rates to a unformatted binary file.
Definition at line 390 of file tw_rate_module.f90.
subroutine tw_rate_module::output_n_p |
subroutine tw_rate_module::read_binary_weak_reaction_data | ( | character(len=*), intent(in) | path | ) |
Read the theoretical weak rates from a unformatted binary file.
In case the binary file is read, no other data has to be read.
[in] | path | Path to binary file |
Definition at line 259 of file tw_rate_module.f90.
|
private |
Reads and counts the amount of theoretical weak rates.
Furthermore, this routine creates the rrate_weak array that is merged into the total rate array rrate.
Edited:
Definition at line 573 of file tw_rate_module.f90.
|
private |
Reads weak reaction rates in log format.
This subroutine is called when parameter_class::iwformat was set to 1. It assumes weak reaction rates that are logarithmic. It will fill the global_class::weak_rate array with rates that are later indicated by the label " ffn". An example of the file that is read here is given by:
neg. daughter Sc45 z=21 n=24 a=45 Q= 0.7678 pos. daughter Ca45 z=20 n=25 a=45 Q= -0.7678 +++ Sc45 --> Ca45 +++ --- Ca45 --> Sc45 --- t9 lrho uf lbeta+ lfte- lnu lbeta- lfte+ lnubar 0.01 1.0 -0.003 -99.999 -99.999 -99.999 -7.302 -99.999 -0.746 0.10 1.0 -0.058 -91.517 5.796 -1.577 -7.302 7.565 -0.746 0.20 1.0 -0.134 -49.637 5.792 -1.266 -7.302 7.083 -0.746 ...
Edited:
[in] | sourcefile | file to read weak rates from |
[out] | cnt | total count of weak rates |
Definition at line 661 of file tw_rate_module.f90.
subroutine, public tw_rate_module::reload_exp_weak_rates |
This routine is used to replace theoretical weak rates with experimental ones it's called when T9 < parameter_class::temp_reload_exp_weak_rates (default 1d-2)
Edited:
Definition at line 1116 of file tw_rate_module.f90.
|
private |
Sorts the entries in weak_rate in increasing order of the decaying nuclide.
Definition at line 1280 of file tw_rate_module.f90.
subroutine, public tw_rate_module::weak_index | ( | real(r_kind), intent(in) | ltemp, |
real(r_kind), intent(in) | rho, | ||
real(r_kind), intent(in) | Ye, | ||
type(weakrate_type), intent(in) | wrt | ||
) |
Set wk_index for given T, \(\rho\) and \(Y_e\) The index is written into module variable wk_index.
After the subroutine has been called the correct indices for a two dimensional interpolation are stored in wk_index.
Edited:
[in] | ltemp | Temperature [GK] |
[in] | rho | Density [g/ccm] |
[in] | ye | Electron fraction |
[in] | wrt | Input weak rate type |
Definition at line 986 of file tw_rate_module.f90.
real(r_kind) function, public tw_rate_module::weak_inter | ( | integer | rate, |
real(r_kind), intent(in) | ltemp, | ||
real(r_kind), intent(in) | lrho, | ||
type(weakrate_type), intent(in) | wrt, | ||
logical, intent(in), optional | force_cubic | ||
) |
Interpolation interface.
This function interpolates the weak rates either bilinear or bicubic
Edited:
[in] | lrho | Log( \( \rho \cdot Y_e \)) |
[in] | ltemp | Temperature [GK] |
[in] | force_cubic | Whether or not a cubic interpolation should be applied |
rate | Index for rate. 1: ft_rate, 2: beta_rate, 3: mue, 4: nu_loss | |
[in] | wrt | Input weak rate type |
Definition at line 1015 of file tw_rate_module.f90.
|
private |
Write the amount of individual reactions to the out.
The rates are always counted, for a certain verbose level they are also printed to the OUT file
Definition at line 180 of file tw_rate_module.f90.
|
private |
Identifier for beta rate.
Definition at line 64 of file tw_rate_module.f90.
|
private |
Identifier for forward rate.
Definition at line 65 of file tw_rate_module.f90.
real(r_kind), public tw_rate_module::mue |
chemical potential of the electron
Definition at line 62 of file tw_rate_module.f90.
|
private |
Identifier for chemical potential.
Definition at line 63 of file tw_rate_module.f90.
|
private |
Definition at line 68 of file tw_rate_module.f90.
|
private |
individual reaction types
Definition at line 68 of file tw_rate_module.f90.
|
private |
Identifier for the neutrino loss.
Definition at line 66 of file tw_rate_module.f90.
|
private |
number of weak reactions taken from Langanke&Martinez-pinedo
Definition at line 59 of file tw_rate_module.f90.
|
private |
Reaction rate representative in global_class::rrate.
Definition at line 60 of file tw_rate_module.f90.
logical, public tw_rate_module::weak |
switch for weak rates
Definition at line 58 of file tw_rate_module.f90.
|
private |
Filename of binary file to save weak rates.
Definition at line 69 of file tw_rate_module.f90.
type(weakrate_type), dimension(:), allocatable, public tw_rate_module::weak_rate |
array containing the weak reactions
Definition at line 56 of file tw_rate_module.f90.
|
private |
Multi-index for the weak rates.
Definition at line 61 of file tw_rate_module.f90.