detailed_balance Module Reference

Module that deals with inverse reaction rates. More...

Functions/Subroutines

subroutine, public init_inverse_rates ()
 Initialize everything for the inverse rates. More...
 
subroutine, public merge_inverse_rates (rrate_array, rrate_length)
 Delete and create new inverse rates. More...
 
subroutine, private create_inverse_rate (forward_rate, inverse_rate)
 Routine to create an inverse rate, based on the forward rate. More...
 
subroutine, public calculate_detailed_balance (temp, rrate_array, inverse_rate, rate)
 Calculate the inverse rate from the forward one via detailed balance. More...
 
subroutine, private write_reac_rate (rate_array, nrea, rate_array_old, nrea_old)
 Verbose routine to output one inverse reaction of chapter 4-9. More...
 
subroutine, private write_reac_verbose_out ()
 Write the verbose output of the reaction rates. More...
 

Variables

integer, private ninv
 Number of inverse rates. More...
 
integer, private ninv_old
 Number of old inverse rates. More...
 
character(len=4), dimension(:), allocatable, private src_ignore
 
integer, private src_ignore_length
 
character(len=4), dimension(:), allocatable, private src_q_reacl
 
integer, private src_q_reacl_length
 
character(len=4), dimension(:), allocatable, private src_q_winvn
 
integer, private src_q_winvn_length
 
 checklist
 
 program
 
 testdir
 
 logfile
 
 arcfile
 

Detailed Description

Module that deals with inverse reaction rates.

This module is only used when parameter_class::use_detailed_balance is enabled. The detailed balance is calculated as follows:

\[ \sigma_\mathrm{inv} = N_A^{-n} \delta \left(\frac{\prod_i A_i}{\prod_j A_j}\right)^{1.5} \left(\frac{\prod_i g_i}{\prod_j g_j}\right) \left( \frac{m_u k_B T}{2\pi \hbar^2}\right)^{1.5 \cdot n} e^{-Q/(k_B T)} \sigma_\mathrm{forw}, \]

where \(n\) is the number of reactants minus the number of products, \( N_A \) is Avogadros number, \(A_i\) is the mass number of reactants of the forward reaction, \(A_j\) the number of products of the forward reaction, \(g_i\) the spin factor \( 2J_i +1\) of reactants (j, for products), Q the q-value of the forward reaction, T the temperature and \( \sigma_\mathrm{forw} \) the forward reaction.

The necessity of calculating detailed balance within the network is given by the tabulated rates that may not easily obey the detailed balance principle otherwise. Furthermore, there exist inconsistencies between the reaclib Q-value and the Q-value calculated by the winvn (version downloaded 24.07.22). The Q-value given in the reaclib are mostly used for calculating the inverse rates. However there seem to exist also some other inconsistencies. An example of an inconsistent q-value in reaclib is given by the reaction: n + k55 -> k56 other inconsistencies are included in, e.g., n + ge89 -> ge90 or p + co52 -> ni53 which deviates by a (almost) constant factor of 4 indicating a different spin that was used in reaclib to calculate the inverse reaction.

Reaction calculated with detailed balance and taken from reaclib (25.07.22)
Reaction calculated with detailed balance and taken from reaclib (25.07.22)

This problem was also topic of Lippuner & Roberts 2017.

The essence of a discussion with H. Schatz is that reaclib reactions can be based on other mass models/data than the winvn table. Every reaclib reaction can have its own set of mass tables. It might be better to calculate inverse reaction yourself in the network as you will be more consistent, especially close to equilibrium values such as NSE. However, this goes to the cost of an inconsistency between the inverse and forward rate as the forward rate involves other mass models for the calculation. For (n,gamma) reactions that are based on Hauser-Feshbach, Hendrik suggested that one could use a scaling relation of the forward rate based on the Q-value differences.

Author
M. Reichert
Date
22.07.22

Function/Subroutine Documentation

◆ calculate_detailed_balance()

subroutine, public detailed_balance::calculate_detailed_balance ( real(r_kind), intent(in)  temp,
type(reactionrate_type), dimension(nreac), intent(in)  rrate_array,
type(reactionrate_type), intent(in)  inverse_rate,
real(r_kind), intent(out)  rate 
)

Calculate the inverse rate from the forward one via detailed balance.

This subroutine calculates the inverse rate via detailed balance. Detailed balance gives a direct relation between forward and inverse reaction.

Note
This subroutine also sets a threshold for inverse rates in case the Q-value for a photodissociation is positive. In this case the rates can diverge towards lower temperatures. This can mess up convergence of the system due to precision errors. For the default FRDM winvn/reaclib q-values this is not a problem, however, when using other mass models or when varying the masses it can become a problem.
See also
Fowler, Caughlan, and Zimmerman (1967)
Author
M. Reichert
Date
22.07.22

Edited:

  • 18.11.22; M.R: implemented more restrictive cutoff for the rates in case of positive q-values
Parameters
[in]rrate_arrayLarge rate array, containing all reactions
[in]inverse_rateinverse reaction rate to evaluate
[in]tempTemperature [GK]
[out]rateValue for the reaction rate

Definition at line 449 of file detailed_balance.f90.

Here is the call graph for this function:

◆ create_inverse_rate()

subroutine, private detailed_balance::create_inverse_rate ( type(reactionrate_type), intent(in)  forward_rate,
type(reactionrate_type), intent(out)  inverse_rate 
)
private

Routine to create an inverse rate, based on the forward rate.

Writes the correct nuclei in the parts array of each reaction. Furthermore, as the param() array of the rates is not used, it is redefined to contain:

  • 1: Index to forward rate in rate array
  • 2: #Products - #Reactants.
  • 3: \( 1.5 \cdot \log(A_{r1} \cdot ...) / (A_{p1} \cdot ...) \), with A_r being the mass number of the reactants and A_p of the products.
  • 4: \( 1.5 \cdot \log(M_u/(2.0\cdot \pi \cdot \hbar^2)) \).
  • 5: \( - \log((2J_{r1}+1) \cdot ...) / ((2J_{p1}+1) \cdot ...) \), with J_r being the spin of the reactants and J_p of the products. . Additionally, the Q-value is calculated from the mass excess.
Author
M. Reichert
Date
22.07.22
Parameters
[in]forward_rateTemporary reaction rate
[out]inverse_rateTemporary reaction rate

Definition at line 323 of file detailed_balance.f90.

Here is the call graph for this function:

◆ init_inverse_rates()

subroutine, public detailed_balance::init_inverse_rates

Initialize everything for the inverse rates.

This subroutine analyzes the parameters to determine which rates have to be replaced and which Q-value has to be used.

Author
M. Reichert
Date
04.08.22

Definition at line 94 of file detailed_balance.f90.

Here is the call graph for this function:

◆ merge_inverse_rates()

subroutine, public detailed_balance::merge_inverse_rates ( type(reactionrate_type), dimension(:), intent(inout), allocatable  rrate_array,
integer, intent(inout)  rrate_length 
)

Delete and create new inverse rates.

This routine first identifies included reverse rates and deletes them. Afterwards, new reverse rates are created and attached to the rate array at the end. The position of the rates at the end of the array is important as they will use the previously calculated forward rates. Weak reactions, fission, and chapter 7 reactions do not have any inverse reaction. This routine also stores the index of the forward reaction in param(1).

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

Definition at line 136 of file detailed_balance.f90.

Here is the call graph for this function:

◆ write_reac_rate()

subroutine, private detailed_balance::write_reac_rate ( type(reactionrate_type), dimension(nrea), intent(inout)  rate_array,
integer, intent(in)  nrea,
type(reactionrate_type), dimension(nrea), intent(inout)  rate_array_old,
integer, intent(in)  nrea_old 
)
private

Verbose routine to output one inverse reaction of chapter 4-9.

This subroutine outputs one forward, one calculated inverse, and the same inverse as given by the reaclib for chapter 4-9. All other chapters do not have forward rates so far. The calculated inverse and reaclib inverse should to large extends agree. There might be a small deviation due to different q-values used in the reaclib and in the winvn. This routine is only called when verbose_level is greater or equal two.

Note
For all chapter 7 reactions, the double counting factor seems to be missing in the inverse rate of Reaclib (31.7.2022).
Author
M. Reichert
Date
30.07.22
Parameters
[in]nrea_oldLength of rate arrays
[in,out]rate_arrayReaction rates with db
[in,out]rate_array_oldReaction rates from reacl.

Definition at line 559 of file detailed_balance.f90.

Here is the call graph for this function:

◆ write_reac_verbose_out()

subroutine, private detailed_balance::write_reac_verbose_out
private

Write the verbose output of the reaction rates.

The rates are always counted, for a certain verbose level they are also printed to the OUT file

Author
M. Reichert
Date
21.07.22

Definition at line 902 of file detailed_balance.f90.

Here is the call graph for this function:

Variable Documentation

◆ arcfile

detailed_balance.arcfile

Definition at line 9 of file detailed_balance.py.

◆ checklist

detailed_balance.checklist

Definition at line 3 of file detailed_balance.py.

◆ logfile

detailed_balance.logfile

Definition at line 8 of file detailed_balance.py.

◆ ninv

integer, private detailed_balance::ninv
private

Number of inverse rates.

Definition at line 63 of file detailed_balance.f90.

◆ ninv_old

integer, private detailed_balance::ninv_old
private

Number of old inverse rates.

Definition at line 64 of file detailed_balance.f90.

◆ program

detailed_balance.program

Definition at line 6 of file detailed_balance.py.

◆ src_ignore

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

Definition at line 67 of file detailed_balance.f90.

◆ src_ignore_length

integer, private detailed_balance::src_ignore_length
private

Definition at line 68 of file detailed_balance.f90.

◆ src_q_reacl

character(len=4), dimension(:), allocatable, private detailed_balance::src_q_reacl
private

Definition at line 69 of file detailed_balance.f90.

◆ src_q_reacl_length

integer, private detailed_balance::src_q_reacl_length
private

Definition at line 70 of file detailed_balance.f90.

◆ src_q_winvn

character(len=4), dimension(:), allocatable, private detailed_balance::src_q_winvn
private

Definition at line 71 of file detailed_balance.f90.

◆ src_q_winvn_length

integer, private detailed_balance::src_q_winvn_length
private

Definition at line 72 of file detailed_balance.f90.

◆ testdir

detailed_balance.testdir

Definition at line 7 of file detailed_balance.py.