tw_rate_module Module Reference

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...
 

Detailed Description

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).

See also
Langanke & Martinez-Pinedo 2001, Fuller et al. 1985, Oda et al. 1994, Weak rate library

Function/Subroutine Documentation

◆ calculate_twr_rate()

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:

  • 26.07.22, MR: Created this subroutine from previously existing code in Jacobi_init
Warning
In case of parameter_class::iwformat equal to two and using timmes chemical potential (parameter_class::use_timmes_mue) mue has to be calculated prior to calling this subroutine.
Parameters
[in,out]rraterate instance
[in]tempTemperature [GK]
[in]rhoDensity [g/ccm]
[in]yeElectron fraction
[out]rat_calcrate value
nulosscalculate neutrino loss?

Definition at line 484 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ init_theoretical_weak_rates()

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.

Author
Moritz Reichert
Date
24.01.21

Definition at line 90 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ initialize_cubic_interp()

subroutine tw_rate_module::initialize_cubic_interp

Calculate derivatives of weak reactions.

These derivatives are necessary for the bicubic interpolation.

Author
: M. Reichert

Definition at line 121 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ merge_theoretical_weak_rates()

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.

Parameters
[in,out]rrate_arrayLarge rate array, containing all reactions
[in,out]rrate_lengthlength of rrate_array

Definition at line 208 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ output_binary_weak_reaction_data()

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.

Author
M. Reichert
Date
27.07.21

Definition at line 390 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ output_n_p()

subroutine tw_rate_module::output_n_p

Definition at line 1217 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ read_binary_weak_reaction_data()

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.

Author
M. Reichert
Date
27.07.21
Parameters
[in]pathPath to binary file

Definition at line 259 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ read_theoretical_weak_rates()

subroutine, private tw_rate_module::read_theoretical_weak_rates
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:

  • 24.01.21, MR - copied it into this separate subroutine

Definition at line 573 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ readweak_logft()

subroutine, private tw_rate_module::readweak_logft ( integer, intent(in)  sourcefile,
integer, intent(out)  cnt 
)
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
  ... 
Remarks
Some rates in the file will replace rates from the reaclib, changing them there may not have an effect as other rates are used then.
See also
parameter_class::temp_reload_exp_weak_rates, parameter_class::iwformat

Edited:

  • 12.01.14
  • 11.03.22 (MR: extended routine for more general grids)
Parameters
[in]sourcefilefile to read weak rates from
[out]cnttotal count of weak rates

Definition at line 661 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ reload_exp_weak_rates()

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)

Remarks
Some rates will replace rates from the reaclib, changing them there may not have an effect for low temperatures as other rates are used then.
See also
parameter_class::temp_reload_exp_weak_rates, parameter_class::iwformat weak_index, weak_inter, temp_grid_weak

Edited:

  • 27.01.21, M.R.
  • 23.03.23, M.R.

Definition at line 1116 of file tw_rate_module.f90.

◆ sort()

subroutine, private tw_rate_module::sort ( integer, intent(in)  cnt)
private

Sorts the entries in weak_rate in increasing order of the decaying nuclide.

Definition at line 1280 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ weak_index()

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.

See also
weak_inter, readweak_logft, parameter_class::iwinterp, inter_module::get_indice_2D

Edited:

  • MR : 19.01.2021 - Removed hard coded temperatures
  • MR : 11.03.2022 - Rates are now on two dimensional grid, rewrote this routine
  • MR : 06.07.2022 - Moved things to inter_module, now it is only an interface
Parameters
[in]ltempTemperature [GK]
[in]rhoDensity [g/ccm]
[in]yeElectron fraction
[in]wrtInput weak rate type

Definition at line 986 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ weak_inter()

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

Author
: M. Reichert
Date
: 14.03.22
Warning
The bicubic interpolation fails at values where -99 is neighboring. There, it produces artifacts. For mue, it is, however, okay.

Edited:

  • MR : 06.07.2022 - Moved things to interpolation module
Parameters
[in]lrhoLog( \( \rho \cdot Y_e \))
[in]ltempTemperature [GK]
[in]force_cubicWhether or not a cubic interpolation should be applied
rateIndex for rate. 1: ft_rate, 2: beta_rate, 3: mue, 4: nu_loss
[in]wrtInput weak rate type
Returns
Weak reaction rate, interpolated at the desired temp and dens

Definition at line 1015 of file tw_rate_module.f90.

Here is the call graph for this function:

◆ write_reac_verbose_out()

subroutine, private tw_rate_module::write_reac_verbose_out
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

Author
M. Reichert
Date
27.01.21

Definition at line 180 of file tw_rate_module.f90.

Here is the call graph for this function:

Variable Documentation

◆ beta_ident

integer, parameter, private tw_rate_module::beta_ident =2
private

Identifier for beta rate.

Definition at line 64 of file tw_rate_module.f90.

◆ ft_ident

integer, parameter, private tw_rate_module::ft_ident =3
private

Identifier for forward rate.

Definition at line 65 of file tw_rate_module.f90.

◆ mue

real(r_kind), public tw_rate_module::mue

chemical potential of the electron

Definition at line 62 of file tw_rate_module.f90.

◆ mue_ident

integer, parameter, private tw_rate_module::mue_ident =1
private

Identifier for chemical potential.

Definition at line 63 of file tw_rate_module.f90.

◆ n_ec

integer, private tw_rate_module::n_ec
private

Definition at line 68 of file tw_rate_module.f90.

◆ n_o

integer, private tw_rate_module::n_o
private

individual reaction types

Definition at line 68 of file tw_rate_module.f90.

◆ nu_loss_ident

integer, parameter, private tw_rate_module::nu_loss_ident =4
private

Identifier for the neutrino loss.

Definition at line 66 of file tw_rate_module.f90.

◆ nweak

integer, private tw_rate_module::nweak
private

number of weak reactions taken from Langanke&Martinez-pinedo

Definition at line 59 of file tw_rate_module.f90.

◆ rrate_weak

type(reactionrate_type), dimension(:), allocatable, private tw_rate_module::rrate_weak
private

Reaction rate representative in global_class::rrate.

Definition at line 60 of file tw_rate_module.f90.

◆ weak

logical, public tw_rate_module::weak

switch for weak rates

Definition at line 58 of file tw_rate_module.f90.

◆ weak_binary_name

character(len=*), parameter, private tw_rate_module::weak_binary_name ='weak_rates.windat'
private

Filename of binary file to save weak rates.

Definition at line 69 of file tw_rate_module.f90.

◆ weak_rate

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.

◆ wk_index

integer, dimension(2,2), private tw_rate_module::wk_index
private

Multi-index for the weak rates.

Definition at line 61 of file tw_rate_module.f90.