fission_rate_module.f90
Go to the documentation of this file.
1 !> @file fission_rate_module.f90
2 !!
3 !! The error file code for this file is ***W19***.
4 !! Contains the module \ref fission_rate_module
5 
6 !> Module to deal with fission reactions
7 !!
8 !! This module contains subroutines to read fission reactions
9 !! and fission fragment distributions. Furthermore, it counts the
10 !! amount of reactions. Many subroutines here were originally implemented
11 !! by M. Eichler and C. Winteler.
12 !!
13 !! @remark All subroutines have often been existed at other places in
14 !! the network and possibly as parts of other subroutines.
15 #include "macros.h"
17  implicit none
18 
19 
20  !> fission rate type, designed to save fragment distribution at initialization
21  type,public :: fissionrate_type
22  integer :: fissnuc_index !< index of parent nucleus
23  integer :: channels !< number of fragment pairs for each fission rate
24  integer :: mode !< fission mode; 1:n-induced fission 2: spontaneous and beta-delayed fission
25  integer :: released_n !< Amount of released neutrons for beta-delayed fission
26  integer :: dimens !< dimension of fissparts array
27  real(r_kind) :: cached !< computed rate
28  real(r_kind) :: averageq !< average Q value, weighted over all fragment channels
29  character(4) :: src !< Source of the reaction
30  integer :: reac_type !< "rrt_sf": spontaneous fission, "rrt_bf": beta-delayed fission, "rrt_nf": neutron-induced fission
31  integer,dimension(:),allocatable :: fissparts !< corresponds to rrate()%parts() array
32  integer,dimension(:,:),allocatable :: cscf_ind !< cscf_ind (i,j) gives the position of the entry \f$\frac{d\dot{Y}_{j}}{dY_{i}}\f$ in the cscf data array
33  real(r_kind),dimension(9) :: param !< rate parameters a0-a9, as given in reaclib file
34  real(r_kind),dimension(:),allocatable :: q_value !< Q value of each fission channel
35  real(r_kind),dimension(:),allocatable :: channelprob !< probability for each fragment pair
36  real(r_kind),dimension(:),allocatable :: ch_amount !< corresponds to rrate(i)%ch_amount(j); for fragment i ch_amount(i) = +1 * channelprob(i)
37  end type fissionrate_type
38 
39  type(fissionrate_type),dimension(:),allocatable,public :: fissrate !< Array storing fission reactions
40 
41  real(r_kind),dimension(:,:),allocatable,private :: beta_delayed_fiss_probs
42  integer,private :: amount_cols !< Amount of columns in the beta-delayed fission file
43 
44 
45  character(len=*), private, parameter :: fiss_binary_name='fiss_rates.windat' !< Name of the binary file containing the fission rates
46 
47 
48  type,private :: fragtype
49  integer :: zp !< Z of fissioning nucleus
50  integer :: ap !< A of fissioning nucleus
51  integer :: neutrons !< number of neutrons emitted
52  integer :: nr_frags !< Number of fragments
53  real(r_kind) :: yn !< Yield for neutrons
54  integer,dimension(:),allocatable :: net_idx !< index of fragment in network
55  integer,dimension(:),allocatable :: z !< Z of fragment
56  integer,dimension(:),allocatable :: a !< A of fragment
57  real(r_kind),dimension(:),allocatable :: y !< Yields of fragment
58  end type fragtype
59  type(fragtype),dimension(:),allocatable,private :: fissfrags_n_induced !< Array storing fragment distributions of neutron-induced fission from file
60  type(fragtype),dimension(:),allocatable,private :: fissfrags_spontaneous !< Array storing fragment distributions of spontaneous fission from file
61  type(fragtype),dimension(:),allocatable,private :: fissfrags_beta_delayed !< Array storing fragment distributions of beta-delayed fission from file
62 
63 
64 
65 
66  integer, public, parameter :: fiss_neglect=5 !< Amount of fission fragments not to be neglected
67  ! in the jacobian in case of vanishing parent
68  integer, public :: nfiss !< Amount of fission rates
69  integer, private :: nfiss_spont !< Amount of spontaneous fission rates
70  integer, private :: nfiss_n_ind !< Amount of neutron induced fission rates
71  integer, private :: nfiss_bdel !< Amount of beta-delayed fission rates
72  integer, private :: nufiss !< total number of fission neutrons
73 
74  integer, private :: n_nf,n_bf,n_sf !< Amount of individual reactions
75  !
76  ! Public and private fields and methods of the module
77  !
78  public:: &
80  private:: &
87 contains
88 
89 
90  !> Initialize the fission reactions
91  !!
92  !! This subroutine counts and reads fission reactions.
93  !! After calling it \ref fissrate will be filled as well as the integer
94  !! \ref nfiss will store the amount of fission reactions.
95  !!
96  !! @author M. Reichert
97  !! @date 25.01.21
98  subroutine init_fission_rates()
102  implicit none
103  integer :: alloc_stat !< Allocation status flag
104 
105  n_nf=0; n_bf=0; n_sf=0; ! Set counters to 0
106  nfiss=0
107 
108  if (fissflag.ne.0) then
109  if (use_prepared_network) then
111  else
112  !-- Count the amount of fission rates
113  call count_fission_rates()
114 
115  !-- Allocate the fission rate arrays
116  allocate(fissrate(nfiss),stat=alloc_stat)
117  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate" failed.',&
118  "init_fission_rates",190001)
119 
120  !-- Read the fission rates into both arrays
121  call read_fission_rates()
122  call add_fission_fragments()
123 
124  !-- Reorder the fragments, having the ones with the largest
125  !! ch_amount at the beginning.
126  call reorder_fragments()
127 
128  end if
129  !-- Write an overview of the reactions
131  end if
132 
133 
134  end subroutine init_fission_rates
135 
136 
137 
138  !> Reorder the fission fragments
139  !!
140  !! Subroutine to reorder the fragments, having the ones with the largest
141  !! ch_amount at the beginning. This is benificial as we can ignore
142  !! some of the entries in the jacobian in case the parent abundance
143  !! is zero.
144  !!
145  !! @author M. Reichert
146  !! @date 22.02.23
147  subroutine reorder_fragments()
149  implicit none
150  integer :: i !< loop variable
151  integer,dimension(:),allocatable :: indices
152 
153  info_entry("reorder_fragments")
154 
155  do i=1,nfiss
156  ! Dont do this for negligible amount of fragments
157  if (fissrate(i)%dimens .le. fiss_neglect+2) cycle
158 
159  ! Allocate indice array
160  if (allocated(indices)) deallocate(indices)
161  allocate(indices(fissrate(i)%dimens-2))
162 
163  ! Order descending according to ch_amount,
164  ! start with the 3rd entry, since the first two are the fissioning nuclei
165  call bubblesort(0,fissrate(i)%dimens-2,fissrate(i)%ch_amount(3:),indices)
166  ! Also reorder the fissparts array
167  call reorder_int(fissrate(i)%fissparts(3:),indices,fissrate(i)%dimens-2)
168 
169  end do
170 
171  info_exit("reorder_fragments")
172 
173  end subroutine reorder_fragments
174 
175 
176  !> Write the amount of individual reactions to the out
177  !!
178  !! The rates are always counted, for a certain verbose level they
179  !! are also printed to the OUT file
180  !!
181  !! @author M. Reichert
182  !! @date 27.01.21
185  implicit none
186  character(len=7) :: tmp !< temporary character for pretty output
187 
188  if (verbose_level .ge. 1) then
189  call write_data_to_std_out("Amount fission rates",int_to_str(nfiss))
190  call write_data_to_std_out("Total number of fission neutrons",int_to_str(nufiss))
191  elseif (verbose_level .ge. 2) then
192  if (nfiss .gt. 0) write(*,"(A)") ""
193  if (nfiss .gt. 0) write(*,"(A)") " Fission rates: "
194  if (nfiss .gt. 0) write(*,"(A)") " |----------------------|"
195  tmp = int_to_str(nfiss)
196  if (nfiss .gt. 0) write(*,"(A)") " | Total :"//adjustr(tmp)//" |"
197  tmp = int_to_str(n_nf)
198  if (n_nf .gt. 0) write(*,"(A)") " | n-induced :"//adjustr(tmp)//" |"
199  tmp = int_to_str(n_bf)
200  if (n_bf .gt. 0) write(*,"(A)") " | b-delayed :"//adjustr(tmp)//" |"
201  tmp = int_to_str(n_sf)
202  if (n_sf .gt. 0) write(*,"(A)") " | spontaneous :"//adjustr(tmp)//" |"
203  if (nfiss .gt. 0) write(*,"(A)") " |----------------------|"
204  if (nfiss .gt. 0) write(*,"(A)") ""
205  end if
206 
207  end subroutine write_reac_verbose_out
208 
209 
210 
211  !> Read the fission reactions and fragment distributions from a binary file
212  !!
213  !! This subroutine reads the fission reactions and fragment distributions
214  !! from a binary, unformatted file. The file is assumed to be created
215  !! beforehand. It is only read when use_prepared_network is set to true.
216  !!
217  !! @author M. Reichert
218  !! @date 21.07.23
223  implicit none
224  character(len=*), intent(in) :: path !< Path to folder with binary files
225  integer :: i !< Loop variable
226  integer :: file_id !< File id
227  integer :: alloc_stat !< Allocation status
228  logical :: is_allocated !< Logical to check if array is allocated
229 
230  if (fissflag .eq. 0) return
231 
232  ! Open an unformatted file
233  file_id = open_unformatted_infile(trim(adjustl(path))//trim(adjustl(fiss_binary_name)))
234 
235  read(file_id) nfiss
236  read(file_id) nufiss
237  read(file_id) n_nf
238  read(file_id) n_bf
239  read(file_id) n_sf
240 
241  ! Allocate the fission rate array
242  allocate(fissrate(nfiss),stat=alloc_stat)
243  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate" failed.',&
244  "read_binary_fission_reaction_data",190001)
245 
246  do i=1,nfiss
247  ! Write all fission rates
248  read(file_id) fissrate(i)%fissnuc_index
249  read(file_id) fissrate(i)%channels
250  read(file_id) fissrate(i)%mode
251  read(file_id) fissrate(i)%src
252  read(file_id) fissrate(i)%dimens
253  read(file_id) fissrate(i)%cached
254  read(file_id) fissrate(i)%averageQ
255  read(file_id) fissrate(i)%reac_type
256 
257  read(file_id) is_allocated
258  if (is_allocated) then
259  ! Allocate fissparts, cscf_ind, q_value, channelprob, and ch_amount array
260  allocate(fissrate(i)%fissparts(fissrate(i)%dimens),stat=alloc_stat)
261  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate(i)%fissparts" failed.',&
262  "read_binary_fission_reaction_data",190001)
263  read(file_id) fissrate(i)%fissparts ! has dimension dimens
264  end if
265 
266  read(file_id) is_allocated
267  if (is_allocated) then
268  allocate(fissrate(i)%cscf_ind(fissrate(i)%dimens,fissrate(i)%dimens),stat=alloc_stat)
269  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate(i)%cscf_ind" failed.',&
270  "read_binary_fission_reaction_data",190001)
271  read(file_id) fissrate(i)%cscf_ind ! has dimension (dimens,dimens)
272  end if
273 
274  read(file_id) fissrate(i)%param ! has dimension 9
275 
276  read(file_id) is_allocated
277  if (is_allocated) then
278  allocate(fissrate(i)%q_value(fissrate(i)%channels),stat=alloc_stat)
279  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate(i)%q_value" failed.',&
280  "read_binary_fission_reaction_data",190001)
281  read(file_id) fissrate(i)%q_value ! has dimension channels
282  end if
283 
284  read(file_id) is_allocated
285  if (is_allocated) then
286  allocate(fissrate(i)%channelprob(fissrate(i)%channels),stat=alloc_stat)
287  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate(i)%channelprob" failed.',&
288  "read_binary_fission_reaction_data",190001)
289  read(file_id) fissrate(i)%channelprob ! has dimension channels
290  end if
291 
292  read(file_id) is_allocated
293  if (is_allocated) then
294  allocate(fissrate(i)%ch_amount(fissrate(i)%dimens),stat=alloc_stat)
295  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate(i)%ch_amount" failed.',&
296  "read_binary_fission_reaction_data",190001)
297  read(file_id) fissrate(i)%ch_amount ! has dimension dimens
298  end if
299 
300 
301  end do
302 
303  close(file_id)
304 
305 
306  end subroutine read_binary_fission_reaction_data
307 
308 
309  !> Save the fission data to a unformatted binary file
310  !!
311  !! This subroutine saves the fission data to a unformatted binary file.
312  !!
313  !! @author M. Reichert
314  !! @date 21.07.23
317  use parameter_class, only: fissflag
318  implicit none
319  character(len=*), intent(in) :: path
320  integer :: i
321  integer :: file_id
322 
323  if (fissflag .eq. 0) return
324  ! Open an unformatted file
325  file_id = open_unformatted_outfile(trim(adjustl(path))//trim(adjustl(fiss_binary_name)))
326  ! Write the fission rates to binary file
327  write(file_id) nfiss
328  write(file_id) nufiss
329  write(file_id) n_nf
330  write(file_id) n_bf
331  write(file_id) n_sf
332 
333 
334  do i=1,nfiss
335  ! Write all fission rates
336  write(file_id) fissrate(i)%fissnuc_index
337  write(file_id) fissrate(i)%channels
338  write(file_id) fissrate(i)%mode
339  write(file_id) fissrate(i)%src
340  write(file_id) fissrate(i)%dimens
341  write(file_id) fissrate(i)%cached
342  write(file_id) fissrate(i)%averageQ
343  write(file_id) fissrate(i)%reac_type
344  write(file_id) allocated(fissrate(i)%fissparts)
345  if (allocated(fissrate(i)%fissparts)) then
346  write(file_id) fissrate(i)%fissparts ! has dimension dimens
347  end if
348  write(file_id) allocated(fissrate(i)%cscf_ind)
349  if (allocated(fissrate(i)%cscf_ind)) then
350  write(file_id) fissrate(i)%cscf_ind ! has dimension (dimens,dimens)
351  end if
352  write(file_id) fissrate(i)%param ! has dimension 9
353  write(file_id) allocated(fissrate(i)%q_value)
354  if (allocated(fissrate(i)%q_value)) then
355  write(file_id) fissrate(i)%q_value ! has dimension channels
356  end if
357  write(file_id) allocated(fissrate(i)%channelprob)
358  if (allocated(fissrate(i)%channelprob)) then
359  write(file_id) fissrate(i)%channelprob ! has dimension channels
360  end if
361  write(file_id) allocated(fissrate(i)%ch_amount)
362  if (allocated(fissrate(i)%ch_amount)) then
363  write(file_id) fissrate(i)%ch_amount ! has dimension dimens
364  end if
365  end do
366 
367  ! Close the file again
368  close(file_id)
369 
371 
372 
373  !> Merge fission rates with larger array.
374  !!
375  !! This subroutine will merge the fission rates.
376  !! Since fission rates are treated separately, there is not really much
377  !! to do. However, in case of beta-delayed fission, the half lifes
378  !! of the beta decays have to be modified to include the fission rates.
379  !!
380  !! @note The other fission rate array \ref fissrate will exist independently
381  !! and are used to calculate the correct equations for the
382  !! fission fragments in \ref jacobian_class::jacobi_init,
383  !! \ref jacobian_class::abchange.
384  !!
385  !! @author M. Reichert
386  !! @date 25.01.21
387  subroutine merge_fission_rates(rrate_array,rrate_length,fiss_count)
391  implicit none
392  type(reactionrate_type),dimension(:),allocatable,intent(inout) :: rrate_array !< Large rate array, containing all reactions
393  integer,intent(inout) :: rrate_length !< length of rrate_array
394  integer,intent(out) :: fiss_count !< Amount of fission rates
395  integer :: new_length !< New length of rrate_array
396 
397  !-- Only do something if there are fission rates
398  if (allocated(fissrate) .and. (.not. use_prepared_network)) then
399  ! New length of the array
400  new_length = rrate_length
401  if (nfiss .ne. 0) then
402  ! Prepare something for the case that fission_format_beta_delayed is
403  ! Set to probability format. In this case, the halflifes of the beta decays
404  ! have to be modified.
405  if (fission_format_beta_delayed .eq. 3) then
406  call modify_halflifes(rrate_array,rrate_length,fissrate,nfiss)
407  ! The length may have gotten modified
408  new_length = rrate_length
409  end if
410  !-- Output the new length
411  rrate_length = new_length
412  end if
413  end if
414  fiss_count = nfiss
415  end subroutine merge_fission_rates
416 
417 
418  !> Modifies half lifes of beta decays to include fission
419  !!
420  !! In case the beta delayed fission rates are given as a probability,
421  !! the beta decays are rescaled to include the fission rates without
422  !! changing the total half life.
423  !! The fission rates are given by:
424  !! \f[
425  !! \lambda_{\text{fiss}} =P_{\text{fiss}} * \lambda_{\beta}
426  !! \f]
427  !! where \f$P_{\text{fiss}}\f$ is the probability of fission and
428  !! \f$\lambda_{\beta}\f$ is the beta decay rate as given initially
429  !! in the Reaclib. The beta decay rate is then modified to
430  !! \f[
431  !! \lambda_{\beta, new} = \left(1-P_{\text{fiss}}\right) * \lambda_{\beta}
432  !! \f]
433  !! There are special cases to consider, namely \f$P_{\text{fiss}}=1\f$ in
434  !! which case the reaclib rate has to be removed completely and
435  !! \f$\lambda_{\beta}\f = 0$ in which case the fission rate has to be removed.
436  !! Note also that fission can have multiple channels, which are by emitting
437  !! different amount of neutrons.
438  !! The fission rates in form of probabilities are given, e.g., in
439  !! [Panov et al. 2005](https://ui.adsabs.harvard.edu/abs/2005NuPhA.747..633P/abstract) or
440  !! [Mumpower et al. 2018](https://ui.adsabs.harvard.edu/abs/2018ApJ...869...14M/abstract).
441  !!
442  !! @author M. Reichert
443  !! @date 30.05.24
444  subroutine modify_halflifes(rrate,rrate_length,fissrate_in,nfiss_in)
449  implicit none
450  type(reactionrate_type),dimension(:),allocatable,intent(inout) :: rrate !< Large rate array, containing all reactions
451  integer,intent(inout) :: rrate_length !< length of rrate_array
452  type(fissionrate_type),dimension(:),allocatable,intent(inout) :: fissrate_in !< Fission rate array
453  integer,intent(inout) :: nfiss_in !< Amount of fission rates
454  integer :: i !< Loop variable
455  real(r_kind),dimension(net_size) :: lambdas !< Lambdas of the beta decay
456  integer :: j !< Loop variable
457  real(r_kind) :: tot_fiss_prob !< Total fission probability
458  integer :: parent_idx !< Index of fissioning nucleus
459  integer :: px_idx !< Channel of the fission rate
460  real(r_kind) :: px !< Probability of the fission rate
461  logical,dimension(rrate_length) :: rrate_mask !< Mask for the rrate array
462  logical,dimension(nfiss_in) :: fissrate_mask !< Mask for the fissrate array
463  type(reactionrate_type),dimension(:),allocatable :: rrate_tmp !< Temporary array to store the rates
464  type(fissionrate_type),dimension(:),allocatable :: fissrate_tmp !< Temporary array to store the fission rates
465  integer :: new_rrate_length !< New length of the rrate array
466  integer :: new_nfiss_in !< New amount of fission rates
467  integer :: alloc_stat !< Allocation status
468  real(r_kind) :: lambda_rate !< Rate of the beta decay
469 
470  info_entry("modify_halflifes")
471 
472  lambdas(:) = 0d0
473 
474  fissrate_mask(:) = .true.
475 
476  ! Ensure that there are reactions in rrate
477  if (allocated(rrate)) then
478  rrate_mask(:) = .true.
479  ! Modify beta decays in rrate array and get the half lifes of the beta decays
480  do i=1,rrate_length
481 
482  ! Ignore things that are definitely not beta decays
483  ! and also theoretical weak rates
484  if ((rrate(i)%reac_src .eq. rrs_twr) .or. &
485  (rrate(i)%reac_src .eq. rrs_nu) .or. &
486  (rrate(i)%reac_src .eq. rrs_fiss) .or. &
487  (rrate(i)%reac_src .eq. rrs_aext) .or. &
488  (rrate(i)%reac_src .eq. rrs_detb)) cycle
489 
490  if (rrate(i)%reac_type.eq.rrt_betm) then
491  part_loop: do j=1,get_nr_reactants(rrate(i)%group)
492  if ((rrate(i)%parts(j).ne.ineu) .and. (rrate(i)%parts(j).ne.ipro)) then
493  ! Get the total fission fraction over all channels (P0n, P1n, ...)
494  tot_fiss_prob = sum(beta_delayed_fiss_probs(rrate(i)%parts(j),:))
495 
496  ! Take care of reaclib type reactions
497  if ((rrate(i)%reac_src .eq. rrs_reacl) .or. &
498  (rrate(i)%reac_src .eq. rrs_wext)) then
499  lambda_rate = dexp(rrate(i)%param(1))
500  lambdas(rrate(i)%parts(j)) = lambdas(rrate(i)%parts(j))+lambda_rate
501  if ((1d0-tot_fiss_prob) .ne. 0d0) then
502  rrate(i)%param(1) = rrate(i)%param(1)+dlog(1d0-tot_fiss_prob)
503  else
504  ! Remove the rate, it will become a fission rate
505  rrate_mask(i) = .false.
506  end if
507  ! Take care of tabulated rates
508  else if (rrate(i)%reac_src .eq. rrs_tabl) then
509 
510  ! The rate should be constant so just take 1GK
511  call tabulated_index(1d0)
512  call calculate_tab_rate(rrate(i),1d0,lambda_rate)
513  ! Save it
514  lambdas(rrate(i)%parts(j)) = lambdas(rrate(i)%parts(j))+lambda_rate
515  if ((1d0-tot_fiss_prob) .ne. 0d0) then
516  call multiply_tab_rate_by_factor(rrate(i),1d0-tot_fiss_prob)
517  else
518  ! Remove the rate, it will become a fission rate
519  rrate_mask(i) = .false.
520  end if
521  else
522  call raise_exception('Unknown reac_src in modify_halflifes.',&
523  "modify_halflifes",190014)
524  end if
525  ! Go out of the loop
526  exit part_loop
527  end if
528  end do part_loop
529  end if
530  end do
531 
532  ! Reduce the rrate array to the correct size
533  new_rrate_length = count(rrate_mask)
534  if (new_rrate_length .ne. rrate_length) then
535 
536  ! Say something if the verbose level is high enough
537  if (verbose_level .ge. 2) then
538  write(*,*) "Reducing the rate array to the correct size by removing "//&
539  int_to_str(rrate_length-new_rrate_length)//&
540  " beta-decays and putting it into fission rates."
541  end if
542 
543  ! Allocate a temporary array to store the content of rrate
544  allocate(rrate_tmp(new_rrate_length),stat=alloc_stat)
545  if ( alloc_stat /= 0) call raise_exception('Allocation of "rrate_tmp" failed.',&
546  "modify_halflifes",190001)
547 
548  rrate_tmp(1:new_rrate_length) = pack(rrate,rrate_mask)
549  deallocate(rrate,stat=alloc_stat)
550  if ( alloc_stat /= 0) call raise_exception('Deallocation of "rrate" failed.',&
551  "modify_halflifes",190002)
552 
553  allocate(rrate(new_rrate_length),stat=alloc_stat)
554  if ( alloc_stat /= 0) call raise_exception('Allocation of "rrate" failed.',&
555  "modify_halflifes",190001)
556 
557  rrate(1:new_rrate_length) = rrate_tmp(1:new_rrate_length)
558  deallocate(rrate_tmp,stat=alloc_stat)
559  if ( alloc_stat /= 0) call raise_exception('Deallocation of "rrate_tmp" failed.',&
560  "modify_halflifes",190002)
561  rrate_length = new_rrate_length
562  end if
563  end if
564 
565  ! Make the correct lambdas in the fission array
566  do i=1,size(fissrate_in)
567  if (fissrate_in(i)%reac_type.eq.rrt_bf) then
568  ! Beta-delayed fission
569  parent_idx = fissrate_in(i)%fissnuc_index
570  px_idx = int(fissrate_in(i)%released_n+1d0)
571  px = beta_delayed_fiss_probs(parent_idx,px_idx)
572  if (lambdas(parent_idx) .eq. 0d0) then
573  ! Remove the rate, the parent is not beta-decaying
574  fissrate_mask(i) = .false.
575  else
576  ! Fissrate
577  fissrate_in(i)%param(:) = 0d0
578  fissrate_in(i)%param(1) = dlog(px * lambdas(parent_idx))
579  end if
580  end if
581  end do
582 
583  ! Reduce the fission rate array to the correct size
584  new_nfiss_in = count(fissrate_mask)
585  if (new_nfiss_in .ne. nfiss_in) then
586  ! Say something if the verbose level is high enough
587  if (verbose_level .ge. 2) then
588  write(*,*) "Reducing the fission rate array to the correct size by removing "//&
589  int_to_str(nfiss_in-new_nfiss_in)//" beta-delayed fission rates."
590  end if
591  ! Allocate a temporary array to store the content of fissrate
592  allocate(fissrate_tmp(new_nfiss_in),stat=alloc_stat)
593  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_tmp" failed.',&
594  "modify_halflifes",190001)
595  j=1
596  do i=1,nfiss_in
597  if (fissrate_mask(i)) then
598  ! Allocate all the necessary stuff in fissrate_tmp
599  allocate(fissrate_tmp(j)%fissparts(fissrate_in(i)%dimens),stat=alloc_stat)
600  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_tmp(j)%fissparts" failed.',&
601  "modify_halflifes",190001)
602  allocate(fissrate_tmp(j)%cscf_ind(fissrate_in(i)%dimens,fissrate_in(i)%dimens),stat=alloc_stat)
603  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_tmp(j)%cscf_ind" failed.',&
604  "modify_halflifes",190001)
605  allocate(fissrate_tmp(j)%q_value(fissrate_in(i)%channels),stat=alloc_stat)
606  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_tmp(j)%q_value" failed.',&
607  "modify_halflifes",190001)
608  allocate(fissrate_tmp(j)%channelprob(fissrate_in(i)%channels),stat=alloc_stat)
609  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_tmp(j)%channelprob" failed.',&
610  "modify_halflifes",190001)
611  allocate(fissrate_tmp(j)%ch_amount(fissrate_in(i)%dimens),stat=alloc_stat)
612  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_tmp(j)%ch_amount" failed.',&
613  "modify_halflifes",190001)
614  ! Finally copy... A "pack" would be easier, but we need to ensure the allocation of the arrays inside the types
615  fissrate_tmp(j) = fissrate_in(i)
616  j = j+1
617  end if
618  end do
619 
620  ! fissrate_tmp(1:new_nfiss_in) = pack(fissrate_in,fissrate_mask)
621  deallocate(fissrate_in,stat=alloc_stat)
622  if ( alloc_stat /= 0) call raise_exception('Deallocation of "fissrate_in" failed.',&
623  "modify_halflifes",190002)
624  allocate(fissrate_in(new_nfiss_in),stat=alloc_stat)
625  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_in" or failed.',&
626  "modify_halflifes",190001)
627 
628  ! Copy back from temporary array
629  do i=1,new_nfiss_in
630  allocate(fissrate_in(i)%fissparts(fissrate_tmp(i)%dimens),stat=alloc_stat)
631  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_in(i)%fissparts" failed.',&
632  "modify_halflifes",190001)
633  allocate(fissrate_in(i)%cscf_ind(fissrate_tmp(i)%dimens,fissrate_tmp(i)%dimens),stat=alloc_stat)
634  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_in(i)%cscf_ind" failed.',&
635  "modify_halflifes",190001)
636  allocate(fissrate_in(i)%q_value(fissrate_tmp(i)%channels),stat=alloc_stat)
637  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_in(i)%q_value" failed.',&
638  "modify_halflifes",190001)
639  allocate(fissrate_in(i)%channelprob(fissrate_tmp(i)%channels),stat=alloc_stat)
640  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_in(i)%channelprob" failed.',&
641  "modify_halflifes",190001)
642  allocate(fissrate_in(i)%ch_amount(fissrate_tmp(i)%dimens),stat=alloc_stat)
643  if ( alloc_stat /= 0) call raise_exception('Allocation of "fissrate_in(i)%ch_amount" failed.',&
644  "modify_halflifes",190001)
645  fissrate_in(i) = fissrate_tmp(i)
646  end do
647 
648 
649  ! fissrate_in(1:new_nfiss_in) = fissrate_tmp(1:new_nfiss_in)
650  deallocate(fissrate_tmp,stat=alloc_stat)
651  if ( alloc_stat /= 0) call raise_exception('Deallocation of "fissrate_tmp" or failed.',&
652  "modify_halflifes",190002)
653  nfiss_in = new_nfiss_in
654  end if
655 
656  ! Deallocate the beta-delayed fission probabilities
657  ! as everything is now in the fission rate array
658  deallocate(beta_delayed_fiss_probs,stat=alloc_stat)
659  if ( alloc_stat /= 0) call raise_exception('Deallocation of "beta_delayed_fiss_probs" failed.',&
660  "modify_halflifes",190002)
661 
662  info_exit("modify_halflifes")
663  end subroutine modify_halflifes
664 
665 
666 
667  !> Count the amount of fission rates
668  !!
669  !! This subroutine counts the amount of fission rates and
670  !! stores the result in \ref nfiss.
671  !! In case that a new fission type will be implemented,
672  !! the reaction rates should be counted here.
673  !! The different formats for the different types of fission reactions are:
674  !! <table>
675  !! <caption id="multi_row">Fission formats</caption>
676  !! <tr><th> Type <th> Value <th> Description
677  !! <tr><td> Spontaneous <td> 0 <td> No rates are read
678  !! <tr><td> Spontaneous <td> 1 <td> Reaclib format
679  !! <tr><td> Spontaneous <td> 2 <td> Half life format
680  !! <tr><td> n-induced <td> 0 <td> No rates are read
681  !! <tr><td> n-induced <td> 1 <td> Reaclib format
682  !! <tr><td> \f$\beta\f$-delayed <td> 0 <td> No rates are read
683  !! <tr><td> \f$\beta\f$-delayed <td> 1 <td> Reaclib format
684  !! <tr><td> \f$\beta\f$-delayed <td> 2 <td> Half life format
685  !! <tr><td> \f$\beta\f$-delayed <td> 3 <td> Probability format
686  !! </table>
687  !!
688  !! @see \ref read_fission_rates
689  !!
690  !! @author M. Reichert
691  !! @date 25.01.21
692  subroutine count_fission_rates()
696  use global_class, only: net_size
698  implicit none
699  integer :: alloc_stat !< Allocation status
700 
701  info_entry("count_fission_rates")
702 
703  ! Spontaneous fission rates
704  select case(fission_format_spontaneous)
705  case(0)
706  nfiss_spont = 0
707  case(1)
709  case(2)
711  case default
712  call raise_exception("Fission format for spontaneous fission rates not implemented, got '"//&
713  int_to_str(fission_format_spontaneous)//"'. Change the "//&
714  "'fission_format_spontaneous' parameter.",&
715  "count_fission_rates",190010)
716  end select
717 
718  ! Neutron induced fission rates
719  select case(fission_format_n_induced)
720  case(0)
721  nfiss_n_ind = 0
722  case(1)
724  case default
725  call raise_exception("Fission format for neutron-induced fission rates not implemented, got '"//&
726  int_to_str(fission_format_n_induced)//"'. Change the "//&
727  "'fission_format_n_induced' parameter.",&
728  "count_fission_rates",190010)
729  end select
730 
731  ! Beta-delayed fission rates
732  select case(fission_format_beta_delayed)
733  case(0)
734  nfiss_bdel = 0
735  case(1)
737  case(2)
739  case(3)
741  ! Initialize an array to store the beta-delayed fission probabilities
742  allocate(beta_delayed_fiss_probs(net_size,amount_cols),stat=alloc_stat)
743  if ( alloc_stat /= 0) call raise_exception('Allocation of "beta_delayed_fiss_probs" failed.',&
744  "count_fission_rates",190001)
745  beta_delayed_fiss_probs(:,:) = 0d0
746  case default
747  call raise_exception("Fission format for beta-delayed fission rates not implemented, got '"//&
748  int_to_str(fission_format_beta_delayed)//"'. Change the "//&
749  "'fission_format_beta_delayed' parameter.",&
750  "count_fission_rates",190010)
751  end select
752 
753  ! Calculate the total amount of fission rates
755 
756  info_exit("count_fission_rates")
757  end subroutine count_fission_rates
758 
759 
760  !> Count the amount of fission rates in reaclib format
761  !!
762  !! This subroutine counts the amount of fission rates in reaclib format
763  !! and stores the result in count_rates. This is necessary to allocate
764  !! the fission array with the correct length. The subroutine can be called
765  !! for any kind of fission (spontaneous, n-induced, beta-delayed).
766  !!
767  !! \b Edited:
768  !! - 30.05.2024 (MR): Made this its own subroutine to add more file formats
769  !! .
770  subroutine count_fission_rates_reaclib_format(fission_rate_file,count_rates)
772  use benam_class, only: benam
773  use format_class
774  implicit none
775  character(len=*), intent(in) :: fission_rate_file !< File containing fission rates
776  integer, intent(out) :: count_rates !< Amount of fission rates
777  integer :: j !< Loop variable
778  integer :: read_stat !< Read status
779  integer :: fissionlib !< File id of fisionrates
780  character(5), dimension(6) :: parts !< Participating nuclei names
781  integer, dimension(6) :: parts_index !< Participating nuclei indices
782  real(r_kind), dimension(9) :: params !< Reaclib fitting parameter
783  character(4) :: src !< Reaclib source label
784  character(1) :: res !< Reaclib weak flag label
785  character(1) :: rev !< Reaclib reverse label
786  real(r_kind) :: q !< Reaclib Q-Value
787  integer :: grp !< Reaclib chapter
788  integer :: ifiss !< reaction rate counter
789 
790  ! Counter for fission rates
791  ifiss = 0
792 
793  ! Open fission file
794  fissionlib = open_infile(fission_rate_file)
795 
796  ! Count the rates
797  fission: do
798  read(fissionlib,my_format(1), iostat = read_stat) &
799  grp, parts(1:6), src, res, rev, q
800  if (read_stat /= 0) exit
801  read(fissionlib,"(4E13.6)") params(1:4)
802  read(fissionlib,"(5E13.6)") params(5:9)
803  if (grp.ne.0) cycle fission
804  parts_index = 0
805  indices: do j=1,6
806  if (parts(j) .eq. ' ') exit indices
807  parts_index(j) = benam(parts(j))
808  ! parts_index(i) = 0 if nuclide i is not in network
809  if (parts_index(j) .eq. 0) cycle fission
810  end do indices
811  ifiss = ifiss + 1
812  end do fission
813 
814  ! Close the fission file again
815  close(fissionlib)
816 
817  count_rates = ifiss
819 
820 
821  !> Count the amount of fission rates in half-life format
822  !!
823  !! This subroutine counts the amount of fission rates in half-life format
824  !! and stores the result in count_rates. This is necessary to allocate
825  !! the fission array with the correct length. The subroutine can be called
826  !! for spontaneous and beta-delayed fission.
827  !! An example of an entry can be:
828  !! \file{
829  !! th198 4.760787e-09 kh20
830  !! th199 2.808520e-08 kh20
831  !! th200 2.841403e-08 kh20
832  !! ...}
833  !! where the first columns gives the nucleus, the second the halflife in seconds
834  !! and the third the source. Note that in this format, no beta-delayed neutrons
835  !! can be specified.
836  !!
837  !! @author M. Reichert
838  !! @date 30.05.24
839  subroutine count_fission_rates_halflife_format(fission_rate_file,count_rates)
841  use benam_class, only: benam
842  use format_class
843  implicit none
844  character(len=*), intent(in) :: fission_rate_file !< File containing fission rates
845  integer, intent(out) :: count_rates !< Amount of fission rates
846  integer :: j !< Loop variable
847  integer :: read_stat !< Read status
848  integer :: fissionlib !< File id of fisionrates
849  character(5), dimension(6) :: parts !< Participating nuclei names
850  integer, dimension(6) :: parts_index !< Participating nuclei indices
851  integer :: ifiss !< reaction rate counter
852  real(r_kind), dimension(9) :: params !< Reaclib fitting parameter
853  character(4) :: src !< Reaclib source label
854 
855  ! Counter for fission rates
856  ifiss = 0
857 
858  ! Open fission file
859  fissionlib = open_infile(fission_rate_file)
860 
861  ! Count the rates
862  fission_hl: do
863  ! Set parts to empty strings
864  parts(1:6) = ' '
865  params(:) = 0d0
866 
867  ! An example line could be: pa231 3.446020e+23 kh20
868  read(fissionlib, *, iostat = read_stat) &
869  parts(1), params(1), src
870  if (read_stat /= 0) exit
871 
872  parts_index = 0
873  indices_hl: do j=1,6
874  if (parts(j) .eq. ' ') exit indices_hl
875  parts_index(j) = benam(parts(j))
876  ! parts_index(i) = 0 if nuclide i is not in network
877  if (parts_index(j) .eq. 0) cycle fission_hl
878  end do indices_hl
879  ifiss = ifiss + 1
880  end do fission_hl
881 
882  ! Close the fission file again
883  close(fissionlib)
884 
885  ! Store the amount of fission rates
886  count_rates = ifiss
888 
889 
890  !> Count the amount of fission rates in probability format
891  !!
892  !! This subroutine counts the amount of fission rates in probability format
893  !! and stores the result in count_rates. This is necessary to allocate
894  !! the fission array with the correct length. The subroutine can be called
895  !! for beta-delayed fission.
896  !! An example of an entry can be:
897  !! \file{
898  !! pa291 mp18
899  !! 0.000 0.008 0.000 0.015 0.000 0.018 0.000 0.025 0.000 0.003
900  !! pa292 mp18
901  !! 0.000 0.001 0.002 0.000 0.004 0.000 0.006 0.000 0.018 0.000
902  !! pa293 mp18
903  !! 0.239 0.001 0.000 0.009 0.000 0.010 0.000 0.014 0.000 0.009
904  !! ...}
905  !! where the first line of an entry gives nucleus name and the source. The second line
906  !! gives the probabilities of the fission channels (i.e., P0nf, P1nf, ...).
907  !! Note that there is a maximum number of columns (i.e., neutron emission) set to 30 here.
908  !! Throughout the file, the amount of possible channels have to stay the same, but they
909  !! can be different for two different files. For example, for Panov et al. 2005 the file
910  !! could look like:
911  !! \file{
912  !! u259 pa05
913  !! 0.990
914  !! u260 pa05
915  !! 0.020
916  !! u261 pa05
917  !! 0.990
918  !! ...}.
919  !!
920  !! @author M. Reichert
921  !! @date 30.05.24
922  subroutine count_fission_rates_probability_format(fission_rate_file,count_rates,nr_cols)
924  use benam_class, only: benam
926  use format_class
927  implicit none
928  character(len=*), intent(in) :: fission_rate_file !< File containing fission rates
929  integer, intent(out) :: count_rates !< Amount of fission rates
930  integer, intent(out) :: nr_cols !< Number of columns in the file
931  integer,parameter :: maxcols = 30!< Maximum number of columns
932  integer :: j !< Loop variable
933  integer :: read_stat !< Read status
934  integer :: fissionlib !< File id of fisionrates
935  character(5) :: parent !< Participating nuclei names
936  integer :: idx_parent !< Participating nuclei indices
937  integer :: ifiss !< reaction rate counter
938  real(r_kind), dimension(maxcols):: buff !< Reaclib fitting parameter
939  character(500) :: dummy !< Reaclib source label
940  character(4) :: src !< Reaclib source label
941  integer :: fiss_probs !< Number of fission probabilities not equal zero
942 
943  ! Counter for fission rates
944  ifiss = 0
945 
946  ! Open fission file
947  fissionlib = open_infile(fission_rate_file)
948 
949  ! Count the rates
950  fission_prob: do
951  buff(:) = -1d0
952 
953  ! An example entry could be: pa250 pa05
954  ! 0.010
955  read(fissionlib, * , iostat = read_stat) parent, src
956 
957  if (read_stat /= 0) exit
958 
959  ! Count columns, maximum allowed ones are 30
960  if (ifiss .eq. 0) then
961  read(fissionlib,'(A)', iostat = read_stat) dummy
962  do j=1,maxcols
963  read(dummy,*,iostat=read_stat) buff(1:j)
964  if (read_stat .ne. 0) exit
965  end do
966  ! Count the number of columns
967  nr_cols = count(buff .ne. -1d0)
968  else
969  read(fissionlib,*, iostat = read_stat) buff(1:nr_cols)
970  if (read_stat .ne. 0) call raise_exception('Reading the fission probabilities failed.',&
971  "count_fission_rates_probability_format",190012)
972  end if
973 
974  ! Count the number of fission probabilities not equal zero
975  fiss_probs = count((buff .ne. -1d0) .and. (buff .ne. 0d0))
976 
977  ! Check if the parent is included in the network
978  idx_parent = benam(parent)
979 
980  if (idx_parent .eq. 0) cycle fission_prob
981 
982  ! Count the amount of fission rates that have a
983  ! probability not equal zero
984  ifiss = ifiss + fiss_probs
985  end do fission_prob
986 
987  ! Close the fission file again
988  close(fissionlib)
989 
990  ! Store the amount of fission rates
991  count_rates = ifiss
993 
994 
995  !> Initializes the rates with fragments
996  !!
997  !! This routine initializes the fission rates with the specified fragments.
998  !! The fragments can differ for the types of fission.
999  !!
1000  !!
1001  !!
1002  !!
1003  !! \b Edited:
1004  !! - 31.05.2024 (MR): Created this subroutine from code parts contained in other routines
1005  !! .
1007  use parameter_class, only: fissflag, &
1013  use global_class, only: isotope
1016  implicit none
1017  integer :: neutfission=0 !< File IDs for abla fiss. distr. file
1018  integer :: betafission=0 !< File IDs for abla fiss. distr. file
1019  integer :: fissindex !< Index for current fission rate
1020  integer :: astat !< Allocation status
1021  real(r_kind) :: q !< Reaclib Q-Value
1022  integer :: nfrac_distr !< Number of fragment distributions in file
1023 
1024 
1025  !-- Open the abla files, containing probabilities of the fiss. fragments
1026  if (fissflag .eq. 3) then
1028  ! Create the same for beta-delayed fission without explicitely reading the file again
1029  nfrac_distr = size(fissfrags_n_induced)
1030  allocate(fissfrags_beta_delayed(nfrac_distr),stat=astat)
1031  if (astat /= 0) call raise_exception('Allocation of "fissfrags_bdel" failed.',&
1032  "add_fission_fragments",190001)
1034  else if (fissflag .eq. 4) then
1035  ! Read the fission fragment files if necessary
1036  ! n_induced
1037  if ((fission_frag_n_induced .eq. 3) .and. (fission_format_n_induced .ne. 0)) then
1039  end if
1040  ! beta-delayed
1041  if ((fission_frag_beta_delayed .eq. 3) .and. (fission_format_beta_delayed .ne. 0)) then
1043  end if
1044  ! spontaneous
1045  if ((fission_frag_spontaneous .eq. 3) .and. (fission_format_spontaneous .ne. 0)) then
1047  end if
1048  else if (fissflag .eq. 5) then
1049  neutfission= open_infile(nfission_file)
1050  betafission= open_infile(bfission_file)
1051  end if
1052 
1053  ! Loop through reactions and add fission fragments
1054  do fissindex=1,nfiss
1055 
1056  if (fissrate(fissindex)%reac_type.eq.rrt_nf) then ! neutron-induced fission
1057  select case (fissflag)
1058  case(1)
1059  call fiss_dist(fissrate(fissindex))
1060  case(2)
1061  call kodtakdist(fissrate(fissindex))
1062  case(3)
1063  call file_fiss_frag(fissrate(fissindex),2)
1064  case(4)
1065  select case (fission_frag_n_induced)
1066  case(1)
1067  call fiss_dist(fissrate(fissindex))
1068  case(2)
1069  call kodtakdist(fissrate(fissindex))
1070  case(3)
1072  case default
1073  call raise_exception("Fission fragment flag ("//trim(adjustl(int_to_str(fission_frag_n_induced)))&
1074  //") not implemented yet. "//&
1075  "Set it to a supported value!",&
1076  "add_fission_fragments",&
1077  190003)
1078  end select
1079  case(5)
1080  call abla_nfiss(fissindex,isotope(fissrate(fissindex)%fissnuc_index)%mass, &
1081  isotope(fissrate(fissindex)%fissnuc_index)%p_nr,neutfission,q)
1082  case default
1083  call raise_exception("Fission flag ("//trim(adjustl(int_to_str(fissflag)))&
1084  //") not implemented yet. "//&
1085  "Set it to a supported value!",&
1086  "add_fission_fragments",&
1087  190003)
1088  end select
1089  else if (fissrate(fissindex)%reac_type.eq.rrt_sf) then ! spontaneous fission
1090  select case (fissflag)
1091  case(1)
1092  call fiss_dist(fissrate(fissindex))
1093  case(2)
1094  call kodtakdist(fissrate(fissindex))
1095  case(3)
1096  call kodtakdist(fissrate(fissindex))
1097  case(4)
1098  select case (fission_frag_spontaneous)
1099  case(1)
1100  call fiss_dist(fissrate(fissindex))
1101  case(2)
1102  call kodtakdist(fissrate(fissindex))
1103  case(3)
1105  case default
1106  call raise_exception("Fission fragment flag ("//trim(adjustl(int_to_str(fission_frag_spontaneous)))&
1107  //") not implemented yet. "//&
1108  "Set it to a supported value!",&
1109  "add_fission_fragments",&
1110  190003)
1111  end select
1112  case(5)
1113  call abla_betafiss(fissindex,isotope(fissrate(fissindex)%fissnuc_index)%mass, &
1114  isotope(fissrate(fissindex)%fissnuc_index)%p_nr,betafission,q)
1115  case default
1116  call raise_exception("Fission flag ("//trim(adjustl(int_to_str(fissflag)))&
1117  //") not implemented yet. "//&
1118  "Set it to a supported value!",&
1119  "add_fission_fragments",&
1120  190003)
1121  end select
1122  else if (fissrate(fissindex)%reac_type.eq.rrt_bf) then ! beta-delayed fission
1123  select case (fissflag)
1124  case(1)
1125  call fiss_dist(fissrate(fissindex))
1126  case(2)
1127  call kodtakdist(fissrate(fissindex))
1128  case(3)
1129  call file_fiss_frag(fissrate(fissindex),2)
1130  case(4)
1131  select case (fission_frag_beta_delayed)
1132  case(1)
1133  call fiss_dist(fissrate(fissindex))
1134  case(2)
1135  call kodtakdist(fissrate(fissindex))
1136  case(3)
1138  case default
1139  call raise_exception("Fission fragment flag ("//trim(adjustl(int_to_str(fission_frag_beta_delayed)))&
1140  //") not implemented yet. "//&
1141  "Set it to a supported value!",&
1142  "add_fission_fragments",&
1143  190003)
1144  end select
1145  case(5)
1146  call abla_betafiss(fissindex,isotope(fissrate(fissindex)%fissnuc_index)%mass, &
1147  isotope(fissrate(fissindex)%fissnuc_index)%p_nr,betafission,q)
1148  case default
1149  call raise_exception("Fission flag ("//trim(adjustl(int_to_str(fissflag)))&
1150  //") not implemented yet. "//&
1151  "Set it to a supported value!",&
1152  "add_fission_fragments",&
1153  190003)
1154  end select
1155  end if
1156 
1157 
1158  end do
1159 
1160  ! Clean up
1161  ! Close the neutron and betafission files again
1162  if (neutfission.gt.0) call close_io_file(neutfission,nfission_file)
1163  if (betafission.gt.0) call close_io_file(betafission,bfission_file)
1164  ! Deallocate the fission fragment arrays in case they were read from file
1165  if (allocated(fissfrags_beta_delayed)) then
1166  deallocate(fissfrags_beta_delayed,stat=astat)
1167  if (astat.ne.0) call raise_exception("Could not deallocate fissfrags_beta_delayed array.",&
1168  "add_fission_fragments",190002)
1169  end if
1170  if (allocated(fissfrags_n_induced)) then
1171  deallocate(fissfrags_n_induced,stat=astat)
1172  if (astat.ne.0) call raise_exception("Could not deallocate fissfrags_n_induced array.",&
1173  "add_fission_fragments",190002)
1174  end if
1175  if (allocated(fissfrags_spontaneous)) then
1176  deallocate(fissfrags_spontaneous,stat=astat)
1177  if (astat.ne.0) call raise_exception("Could not deallocate fissfrags_spontaneous array.",&
1178  "add_fission_fragments",190002)
1179  end if
1180  end subroutine add_fission_fragments
1181 
1182 
1183 
1184 
1185  !> Read the fission rates, splitted into different types of fission
1186  !!
1187  !! This routine serves as interface between the reading routines and
1188  !! the different fission formats and types. In case a new type of
1189  !! fission will be implemented, the reaction rates should be read here.
1190  !! The different formats for the different types of fission reactions are:
1191  !! <table>
1192  !! <caption id="multi_row">Fission formats</caption>
1193  !! <tr><th> Type <th> Value <th> Description
1194  !! <tr><td> Spontaneous <td> 0 <td> No rates are read
1195  !! <tr><td> Spontaneous <td> 1 <td> Reaclib format
1196  !! <tr><td> Spontaneous <td> 2 <td> Half life format
1197  !! <tr><td> n-induced <td> 0 <td> No rates are read
1198  !! <tr><td> n-induced <td> 1 <td> Reaclib format
1199  !! <tr><td> \f$\beta\f$-delayed <td> 0 <td> No rates are read
1200  !! <tr><td> \f$\beta\f$-delayed <td> 1 <td> Reaclib format
1201  !! <tr><td> \f$\beta\f$-delayed <td> 2 <td> Half life format
1202  !! <tr><td> \f$\beta\f$-delayed <td> 3 <td> Probability format
1203  !! </table>
1204  !!
1205  !! @author M. Reichert
1206  !! @date 31.05.2024
1207  subroutine read_fission_rates()
1212  implicit none
1213  integer :: fisscount !< Amount of fission rates
1214 
1215  ! Keep track of fission rates
1216  fisscount = 1
1217 
1218  ! Read the spontaneous fission rates
1219  select case (fission_format_spontaneous)
1220  case(0)
1221  ! Dont read anything
1222  continue
1223  case(1)
1224  ! Read the spontaneous fission rates in reaclib format
1226  case(2)
1227  ! Read the spontaneous fission rates in halflife format
1229  case default
1230  call raise_exception("Fission format for spontaneous fission rates not implemented, got '"//&
1231  int_to_str(fission_format_spontaneous)//"'. Change the "//&
1232  "'fission_format_spontaneous' parameter.",&
1233  "read_fission_rates",190010)
1234  end select
1235 
1236  ! Read the neutron-induced fission rates
1237  select case (fission_format_n_induced)
1238  case(0)
1239  ! Dont read anything
1240  continue
1241  case(1)
1243  case default
1244  call raise_exception("Fission format for neutron-induced fission rates not implemented, got '"//&
1245  int_to_str(fission_format_n_induced)//"'. Change the "//&
1246  "'fission_format_n_induced' parameter.",&
1247  "read_fission_rates",190010)
1248  end select
1249 
1250  ! Read the beta-delayed fission rates
1251  select case (fission_format_beta_delayed)
1252  case(0)
1253  ! Dont read anything
1254  continue
1255  case(1)
1256  ! Read the beta-delayed fission rates in reaclib format
1258  case(2)
1259  ! Read the beta-delayed fission rates in halflife format
1261  case(3)
1262  ! Read the beta-delayed fission rates in probability format
1264  case default
1265  call raise_exception("Fission format for beta-delayed fission rates not implemented, got '"//&
1266  int_to_str(fission_format_beta_delayed)//"'. Change the "//&
1267  "'fission_format_beta_delayed' parameter.",&
1268  "read_fission_rates",190010)
1269  end select
1270 
1271  end subroutine read_fission_rates
1272 
1273 
1274  !> Reads fission rates in Reaclib format
1275  !!
1276  !! This subroutine reads the fission rates that are specified in the file
1277  !! parameter_class::fission_rates. Fission rates are implemented with the same
1278  !! fits and a similar format as usual reaclib rates (see reaclib_rate_module::read_reaclib).
1279  !! An example looks like:
1280  !! \file{
1281  !! th256 ms99w 0.00000E+00
1282  !! 1.988512E-02 0.000000E+00 0.000000E+00 0.000000E+00
1283  !! 0.000000E+00 0.000000E+00 0.000000E+00
1284  !! th257 ms99w 0.00000E+00
1285  !! -1.483947E+00 0.000000E+00 0.000000E+00 0.000000E+00
1286  !! 0.000000E+00 0.000000E+00 0.000000E+00
1287  !! ...}
1288  !! The chapter structure is the same as in reaclib, but only the parent nucleus
1289  !! has an entry. The products are given by the respective fission fragment
1290  !! distribution.
1291  !!
1292  !! @returns \ref fissrate filled with fission rates
1293  !! @author M. Eichler
1294  !!
1295  !! \b Edited:
1296  !! - 25.01.21 (MR): Moved from init_network to this new independent subroutine
1297  !! - 27.05.24 (MR): Moved initializing fragments to independent subroutine
1298  !! -
1299  subroutine read_fission_rates_reaclib_format(fission_path,reac_type,start_idx)
1300  use parameter_class, only: fissflag
1301  use benam_class, only: benam
1302  use global_class, only: isotope, ineu
1304  use format_class
1305  use error_msg_class
1306  implicit none
1307  character(len=*), intent(in) :: fission_path !< Path to fission rates file
1308  integer, intent(in) :: reac_type !< Reaction type
1309  integer, intent(inout) :: start_idx !< Start index for fission rates
1310  integer :: read_stat !< Read status
1311  integer :: fissionlib !< File id of fisionrates
1312  character(5), dimension(6) :: parts !< Participating nuclei names
1313  integer, dimension(6) :: parts_index !< Participating nuclei indices
1314  real(r_kind), dimension(9) :: params !< Reaclib fitting parameter
1315  character(4) :: src !< Reaclib source label
1316  character(1) :: res !< Reaclib weak flag label
1317  character(1) :: rev !< Reaclib reverse label
1318  real(r_kind) :: q !< Reaclib Q-Value
1319  integer :: grp !< Reaclib chapter
1320  integer :: group_index !< Storage for the current chapter
1321  integer :: fissindex !< Index for current fission rate
1322  integer :: j !< Loop variable
1323 
1324 
1325 
1326  !-- Open the fission library file
1327  fissionlib = open_infile(fission_path)
1328 
1329  ! make sure counting index for reaction rates continues for fission rates
1330  fissindex = start_idx
1331  fission_loop: do
1332  read(fissionlib,my_format(1), iostat = read_stat) &
1333  grp, parts(1:6), src, res, rev, q
1334  if (read_stat /= 0) exit fission_loop
1335  read(fissionlib,"(4E13.6)") params(1:4)
1336  read(fissionlib,"(5E13.6)") params(5:9)
1337 
1338  select case (grp)
1339  case (1:11)
1340  group_index = grp
1341  cycle fission_loop
1342  case default
1343  parts_index = 0
1344  inner_fission: do j=1,6
1345  if (parts(j) .eq. ' ') exit inner_fission
1346  parts_index(j) = benam(parts(j))
1347  ! parts_index(i) = 0 if nuclide i is not in network
1348  if (parts_index(j) .eq. 0) cycle fission_loop
1349  end do inner_fission
1350  end select
1351 
1352  ! Keep track for debugging
1353  if (reac_type .eq. rrt_nf) n_nf=n_nf+1 !< count neutron induced fission
1354  if (reac_type .eq. rrt_sf) n_sf=n_sf+1 !< count spontaneous fission
1355  if (reac_type .eq. rrt_bf) n_bf=n_bf+1 !< count beta delayed fission
1356 
1357  ! Set the amount of released neutrons per default to 0
1358  fissrate(fissindex)%released_n = 0
1359 
1360  if (reac_type.eq.rrt_nf) then ! neutron-induced fission
1361  fissrate(fissindex)%reac_type = reac_type
1362  fissrate(fissindex)%param = params
1363  ! Save fission nucleus index
1364  if (parts_index(1) .eq. ineu) then
1365  fissrate(fissindex)%fissnuc_index = parts_index(2)
1366  else
1367  fissrate(fissindex)%fissnuc_index = parts_index(1)
1368  end if
1369 
1370  else if (reac_type.eq.rrt_sf) then ! spontaneous fission
1371  fissrate(fissindex)%reac_type = reac_type
1372  fissrate(fissindex)%param = params
1373  ! Save fission nucleus index
1374  fissrate(fissindex)%fissnuc_index = parts_index(1)
1375 
1376  else if (reac_type.eq.rrt_bf) then ! beta-delayed fission
1377  ! Check if there should be neutrons released
1378  fissrate(fissindex)%released_n = count(parts_index(:) .eq. ineu)
1379  fissrate(fissindex)%reac_type = reac_type
1380  fissrate(fissindex)%param = params
1381  ! Save fission nucleus index
1382  fissrate(fissindex)%fissnuc_index = parts_index(1)
1383 
1384 
1385  end if
1386  fissrate(fissindex)%src = src
1387 
1388 
1389  fissindex = fissindex + 1
1390  end do fission_loop
1391 
1392  ! Return the number of fission rates that have been added
1393  start_idx = fissindex
1394 
1395  ! Close the fission file again
1396  close(fissionlib)
1397 
1398  end subroutine read_fission_rates_reaclib_format
1399 
1400 
1401 
1402  !> Read the fission rates in probability format
1403  !!
1404  !! This subroutine reads fission rates in probability format
1405  !! and stores the result in count_rates. The subroutine can be called
1406  !! for beta-delayed fission and will throw an error for other types of fission.
1407  !! An example of an entry can be:
1408  !! \file{
1409  !! pa291 mp18
1410  !! 0.000 0.008 0.000 0.015 0.000 0.018 0.000 0.025 0.000 0.003
1411  !! pa292 mp18
1412  !! 0.000 0.001 0.002 0.000 0.004 0.000 0.006 0.000 0.018 0.000
1413  !! pa293 mp18
1414  !! 0.239 0.001 0.000 0.009 0.000 0.010 0.000 0.014 0.000 0.009
1415  !! ...}
1416  !! where the first line of an entry gives nucleus name and the source. The second line
1417  !! gives the probabilities of the fission channels (i.e., P0nf, P1nf, ...).
1418  !! Note that there is a maximum number of columns (i.e., neutron emission) set to 30 here.
1419  !! Throughout the file, the amount of possible channels have to stay the same, but they
1420  !! can be different for two different files. For example, for Panov et al. 2005 the file
1421  !! could look like:
1422  !! \file{
1423  !! u259 pa05
1424  !! 0.990
1425  !! u260 pa05
1426  !! 0.020
1427  !! u261 pa05
1428  !! 0.990
1429  !! ...}.
1430  !!
1431  !! @author M. Reichert
1432  !! @date 31.05.24
1433  subroutine read_fission_rates_probability_format(fission_path,reac_type,start_idx)
1434  use benam_class, only: benam
1435  use global_class, only: isotope, ineu, net_names
1436  use error_msg_class, only: raise_exception
1438  use format_class
1439  implicit none
1440  character(len=*), intent(in) :: fission_path!< Path to fission rates file
1441  integer, intent(in) :: reac_type !< Reaction type
1442  integer, intent(inout) :: start_idx !< Start index for fission rates
1443  ! Internal variables
1444  integer :: read_stat !< Read status
1445  integer :: fissionlib !< File id of fisionrates
1446  character(5), dimension(6) :: parts !< Participating nuclei names
1447  integer, dimension(6) :: parts_index !< Participating nuclei indices
1448  real(r_kind), dimension(9) :: params !< Reaclib fitting parameter
1449  character(4) :: src !< Reaclib source label
1450  integer :: group_index !< Storage for the current chapter
1451  integer :: fissindex !< Index for current fission rate
1452  integer :: j !< Loop variable
1453  integer :: idx_parent !< Index of the parent nucleus
1454  character(5) :: parent !< Parent nucleus name
1455  real(r_kind),dimension(amount_cols):: pnf !< Probability of beta-delayed neutron emission fission
1456 
1457  ! Counter for fission rates
1458  fissindex = start_idx
1459 
1460  ! Open fission file
1461  fissionlib = open_infile(fission_path)
1462 
1463  ! Count the rates
1464  fission_prob: do
1465  ! An example entry could be: pa250 pa05
1466  ! 0.010
1467  read(fissionlib,*, iostat = read_stat) parent, src
1468  if (read_stat /= 0) exit
1469  read(fissionlib,*, iostat = read_stat) pnf
1470  ! Throw error if the file ends with a parent and a source
1471  if (read_stat /= 0) call raise_exception("Could not read fission probabilities.",&
1472  "read_fission_rates_probability_format",&
1473  190012)
1474 
1475 
1476  ! Check if the parent is included in the network
1477  idx_parent = benam(parent)
1478  if (idx_parent .eq. 0) cycle fission_prob
1479 
1480 
1481  parts_index(:) = 0
1482  parts_index(1) = idx_parent
1483 
1484  do j=1,amount_cols
1485  if (pnf(j) .eq. 0d0) cycle
1486 
1487  params(:) = 0d0
1488  params(1) = pnf(j) ! This will be changed once the reactions are merged
1489  params(2) = j-1 ! Gives the amount of released neutrons
1490 
1491  if (reac_type.eq.rrt_bf) then ! beta-delayed fission
1492  group_index = 2
1493  fissrate(fissindex)%reac_type = reac_type
1494  fissrate(fissindex)%param = params
1495  ! Save fission nucleus index
1496  fissrate(fissindex)%fissnuc_index = parts_index(1)
1497  fissrate(fissindex)%released_n = j-1
1498  else
1499  call raise_exception("Fission type not implemented for probability format!",&
1500  "read_fission_rates_probability_format",&
1501  190011)
1502  end if
1503  ! Put also the source there
1504  fissrate(fissindex)%src = src
1505 
1506  ! Also save it in a separate array that is more easier to
1507  ! access by index
1508  beta_delayed_fiss_probs(idx_parent,j) = pnf(j)
1509 
1510  ! Keep track for debugging
1511  if (reac_type .eq. rrt_bf) n_bf=n_bf+1 ! count beta delayed fission
1512 
1513  fissindex = fissindex + 1
1514  end do
1515  ! Check input
1516  if (sum(beta_delayed_fiss_probs(idx_parent,:)) .gt. 1d0) then
1517  call raise_exception("Sum of beta-delayed fission probabilities for "//net_names(idx_parent)//&
1518  " is larger than one.", "read_fission_rates_probability_format",&
1519  190017)
1520  end if
1521  end do fission_prob
1522 
1523  ! Close the fission file again
1524  close(fissionlib)
1525 
1527 
1528 
1529  !> Read fission rates in half-life format
1530  !!
1531  !! This subroutine reads fission rates in half-life format.
1532  !! The subroutine can be called for spontaneous and beta-delayed fission
1533  !! and throws an error otherwise.
1534  !! An example of an entry can be:
1535  !! \file{
1536  !! th198 4.760787e-09 kh20
1537  !! th199 2.808520e-08 kh20
1538  !! th200 2.841403e-08 kh20
1539  !! ...}
1540  !! where the first columns gives the nucleus, the second the halflife in seconds
1541  !! and the third the source. Note that in this format, no beta-delayed neutrons
1542  !! can be specified.
1543  !!
1544  !! @author M. Reichert
1545  !! @date 30.05.24
1546  subroutine read_fission_rates_halflife_format(fission_path,reac_type,start_idx)
1547  use benam_class, only: benam
1548  use global_class, only: isotope, ineu
1549  use error_msg_class, only: raise_exception
1551  use format_class
1552  implicit none
1553  character(len=*), intent(in) :: fission_path !< Path to fission rates file
1554  integer, intent(in) :: reac_type !< Reaction type
1555  integer, intent(inout) :: start_idx !< Start index for fission rates
1556  integer :: read_stat !< Read status
1557  integer :: fissionlib !< File id of fisionrates
1558  character(5), dimension(6) :: parts !< Participating nuclei names
1559  integer, dimension(6) :: parts_index !< Participating nuclei indices
1560  real(r_kind), dimension(9) :: params !< Reaclib fitting parameter
1561  character(4) :: src !< Reaclib source label
1562  integer :: group_index !< Storage for the current chapter
1563  integer :: fissindex !< Index for current fission rate
1564  integer :: j !< Loop variable
1565 
1566 
1567 
1568  !-- Open the fission library file
1569  fissionlib = open_infile(fission_path)
1570 
1571  ! make sure counting index for reaction rates continues for fission rates
1572  fissindex = start_idx
1573  fission_loop: do
1574  ! Set parts to empty strings
1575  parts(1:6) = ' '
1576  params(:) = 0d0
1577  read(fissionlib, *, iostat = read_stat) &
1578  parts(1), params(1), src
1579  if (read_stat /= 0) exit fission_loop
1580 
1581  params(1) = dlog(dlog(2d0)/params(1)) ! Convert half-life to reaclib a0
1582  parts_index = 0
1583  inner_fission: do j=1,6
1584  if (parts(j) .eq. ' ') exit inner_fission
1585  parts_index(j) = benam(parts(j))
1586  ! parts_index(i) = 0 if nuclide i is not in network
1587  if (parts_index(j) .eq. 0) cycle fission_loop
1588  end do inner_fission
1589 
1590  ! Keep track for debugging
1591  if (reac_type .eq. rrt_nf) n_nf=n_nf+1 ! count neutron induced fission
1592  if (reac_type .eq. rrt_sf) n_sf=n_sf+1 ! count spontaneous fission
1593  if (reac_type .eq. rrt_bf) n_bf=n_bf+1 ! count beta delayed fission
1594 
1595  ! Per default set number of released neutrons to zero (important for bdel fission)
1596  fissrate(fissindex)%released_n = 0
1597  if (reac_type .eq.rrt_sf) then ! spontaneous fission
1598  fissrate(fissindex)%reac_type = reac_type
1599  fissrate(fissindex)%param = params
1600  ! Save fission nucleus index
1601  fissrate(fissindex)%fissnuc_index = parts_index(1)
1602 
1603  else if (reac_type .eq. rrt_bf) then ! beta-delayed fission
1604  fissrate(fissindex)%reac_type = reac_type
1605  fissrate(fissindex)%param = params
1606  ! Save fission nucleus index
1607  fissrate(fissindex)%fissnuc_index = parts_index(1)
1608  else
1609  call raise_exception("Fission type not implemented for half life format!",&
1610  "read_fission_rates_halflife_format",&
1611  190011)
1612  end if
1613 
1614  ! Set some reaclib chapter, note that n induced fission should
1615  ! not have a halflife format
1616  group_index = 2
1617 
1618  ! Put the source to the rate
1619  fissrate(fissindex)%src = src
1620 
1621  fissindex = fissindex + 1
1622  end do fission_loop
1623 
1624  ! Return the number of fission rates that have been added
1625  start_idx = fissindex
1626 
1627  ! Close the fission file again
1628  close(fissionlib)
1629 
1630  end subroutine read_fission_rates_halflife_format
1631 
1632 
1633 
1634 
1635 
1636 
1637  !> Determines fission fragment mass distribution as described in Panov et al. 2001.
1638  !!
1639  !! This routine is called for \ref parameter_class::fissflag = 1 and fills
1640  !! the array \ref fissrates with indices of the fragment nuclei.
1641  !!
1642  !! @see [Panov et al., Nuc. Phys. A688 2001](https://ui.adsabs.harvard.edu/abs/2001NuPhA.688..587P/abstract)
1643  !! @author C. Winteler
1644  subroutine fiss_dist(fissrate_inout)
1645  use global_class, only: isotope, ineu
1646  use benam_class, only: minmax, findaz
1648  implicit none
1649  type(fissionrate_type),intent(inout) :: fissrate_inout !< Temporary fission rate
1650  ! Internal variables
1651  integer :: mass !< mass of fissioning nucleus
1652  integer :: pnr !< proton number of fissioning nucleus
1653  integer :: af,zf !< mass and proton number of "compound" nucleus
1654  integer :: a1,z1 !< mass and proton number of fragment 1
1655  integer :: a2,z2 !< mass and proton number of fragment 2
1656  integer :: nemiss !< number of fission neutrons
1657  integer :: fiss_mode !< specifies fission mode
1658  real(r_kind) :: qval !< Qvalue of fission reaction
1659  integer :: alloc_stat!< Allocation status
1660 
1661  info_entry("fiss_dist")
1662 
1663  ! Counter for additional neutrons in case of
1664  ! fragment beyond dripline.
1665  nemiss = 0
1666 
1667  ! Get the mass and proton number
1668  mass = isotope(fissrate_inout%fissnuc_index)%mass
1669  pnr = isotope(fissrate_inout%fissnuc_index)%p_nr
1670 
1671  !neutron induced fission
1672  if (fissrate_inout%reac_type .eq. rrt_nf) then
1673  fiss_mode = 1
1674  af = mass+1
1675  zf = pnr
1676  !spontaneous fission
1677  else if (fissrate_inout%reac_type .eq. rrt_sf) then
1678  fiss_mode = 2
1679  af = mass
1680  zf = pnr
1681  !beta delayed fission
1682  else if (fissrate_inout%reac_type .eq. rrt_bf) then
1683  fiss_mode = 3
1684  ! Take care of beta-delayed neutron emission fission
1685  ! The amount of neutrons is stored in the second parameter
1686  ! of the fission rate
1687  nemiss = fissrate_inout%released_n
1688  af = mass-nemiss
1689  zf = pnr + 1
1690  end if
1691 
1692  select case (af)
1693  case(255:265) ! -> symmetric fission
1694  a2 = nint(af/2.d0)
1695  z2 = nint(zf/2.d0)
1696  a1 = af - a2
1697  z1 = zf - z2
1698  case default ! -> asymmetric fission
1699  a2 = 130
1700  z2 = nint(52.01 - (zf - 80)/1.d1)
1701  ! a2 = max(a2,af-a2) !modification to prevent very n-rich
1702  ! z2 = max(z2,zf-z2) !fragments
1703  a1 = af - a2
1704  z1 = zf - z2
1705  end select
1706 
1707  if (a1.gt.minmax(z1,2)) then
1708  nemiss = a1 - minmax(z1,2)
1709  a1 = minmax(z1,2)
1710  else if (a1.lt.minmax(z1,1)) then
1711  if (verbose_level .ge. 2) then
1712  print *, 'a1 error:', z1, a1, zf, af
1713  print *, 'adjusting...'
1714  end if
1715  a1 = minmax(z1,1)
1716  a2 = af - a1
1717  end if
1718 
1719  if (a2.gt.minmax(z2,2)) then
1720  nemiss = nemiss + a2 - minmax(z2,2)
1721  a2 = minmax(z2,2)
1722 
1723  else if (a2.lt.minmax(z2,1)) then
1724  if (verbose_level .ge. 2) then
1725  print *, 'a2 error:', z2, a2, zf, af
1726  print *, 'adjusting...'
1727  end if
1728  a2 = minmax(z2,1)
1729  a1 = af - a2
1730  end if
1731  nufiss = nufiss + nemiss
1732 
1733  fissrate_inout%fissnuc_index = findaz(mass,pnr)
1734  fissrate_inout%channels = 1
1735 
1736  allocate(fissrate_inout%channelprob(1),fissrate_inout%q_value(1),stat=alloc_stat)
1737  if (alloc_stat.ne.0) call raise_exception("Could not allocate memory for channelprob and q_value.",&
1738  "fiss_dist",190001)
1739 
1740  fissrate_inout%channelprob(1) = 1.d0
1741  fissrate_inout%mode = fiss_mode
1742 
1743  if (fiss_mode.eq.1) then
1744  fissrate_inout%dimens = 4
1745  allocate(fissrate_inout%fissparts(fissrate_inout%dimens), &
1746  fissrate_inout%ch_amount(fissrate_inout%dimens),stat=alloc_stat)
1747  if (alloc_stat.ne.0) call raise_exception("Could not allocate 'fissparts' and 'ch_amount'.",&
1748  "fiss_dist",190001)
1749  fissrate_inout%fissparts(1) = ineu
1750  fissrate_inout%ch_amount(1) = float(nemiss) - 1.d0
1751  fissrate_inout%fissparts(2) = fissrate_inout%fissnuc_index
1752  fissrate_inout%ch_amount(2) = -1.d0
1753  fissrate_inout%fissparts(3) = findaz(a1,z1)
1754  fissrate_inout%ch_amount(3) = 1.d0
1755  fissrate_inout%fissparts(4) = findaz(a2,z2)
1756  fissrate_inout%ch_amount(4) = 1.d0
1757  fissrate_inout%q_value(1) = (1-nemiss) * isotope(ineu)%mass_exc + &
1758  isotope(fissrate_inout%fissnuc_index)%mass_exc - &
1759  isotope(fissrate_inout%fissparts(3))%mass_exc - &
1760  isotope(fissrate_inout%fissparts(4))%mass_exc
1761  else ! fiss_mode 2 and 3
1762  if (nemiss .eq. 0) then ! no fission neutrons
1763  fissrate_inout%dimens = 3
1764  allocate(fissrate_inout%fissparts(fissrate_inout%dimens),&
1765  fissrate_inout%ch_amount(fissrate_inout%dimens),stat=alloc_stat)
1766  if (alloc_stat.ne.0) call raise_exception("Could not allocate 'fissparts' and 'ch_amount'.",&
1767  "fiss_dist",190001)
1768  else ! fission neutrons are produced
1769  fissrate_inout%dimens = 4
1770  allocate(fissrate_inout%fissparts(fissrate_inout%dimens),&
1771  fissrate_inout%ch_amount(fissrate_inout%dimens),stat=alloc_stat)
1772  if (alloc_stat.ne.0) call raise_exception("Could not allocate 'fissparts' and 'ch_amount'.",&
1773  "fiss_dist",190001)
1774  fissrate_inout%fissparts(4) = ineu
1775  fissrate_inout%ch_amount(4) = float(nemiss)
1776  end if
1777  fissrate_inout%fissparts(1) = fissrate_inout%fissnuc_index
1778  fissrate_inout%ch_amount(1) = -1.d0
1779  fissrate_inout%fissparts(2) = findaz(a1,z1)
1780  fissrate_inout%ch_amount(2) = 1.d0
1781  fissrate_inout%fissparts(3) = findaz(a2,z2)
1782  fissrate_inout%ch_amount(3) = 1.d0
1783 
1784  if (fissrate_inout%fissparts(2) .eq. -1) then
1785  call raise_exception("Fragment was not found in nucleus library, does your library have 'holes'? "//&
1786  "(A="//int_to_str(a1)//", Z="//int_to_str(z1)//")",&
1787  "fiss_dist",190013)
1788  end if
1789  if (fissrate_inout%fissparts(3) .eq. -1) then
1790  call raise_exception("Fragment was not found in nucleus library, does your library have 'holes'? "//&
1791  "(A="//int_to_str(a2)//", Z="//int_to_str(z2)//")",&
1792  "fiss_dist",190013)
1793  end if
1794 
1795  fissrate_inout%q_value(1) = isotope(fissrate_inout%fissparts(1))%mass_exc - &
1796  isotope(fissrate_inout%fissparts(2))%mass_exc - &
1797  isotope(fissrate_inout%fissparts(3))%mass_exc - &
1798  nemiss * isotope(ineu)%mass_exc
1799  end if
1800 
1801  qval = fissrate_inout%q_value(1)
1802  fissrate_inout%averageQ = qval
1803  info_exit("fiss_dist")
1804 
1805  return
1806 
1807  end subroutine fiss_dist
1808 
1809 
1810 
1811  !>
1812  !! Kodama-Takahashi distribution
1813  !!
1814  !! Subroutine that calculates fission fragment distribution as
1815  !! described in Kodama & Takahashi, Nuc. Phys. Sec.A, Vol. 239,
1816  !! Issue 3, p. 489-510, 1975.
1817  !!
1818  !! @see [Kodama & Takahashi 1975](https://ui.adsabs.harvard.edu/abs/1975NuPhA.239..489K/abstract)
1819  !!
1820  !! \b Edited:
1821  !! - 11.01.14
1822  !! - 05.12.17 (ME)
1823  !! - 30.05.2024 (MR): Fixed bug related to neutron emission if fragment was not in Sunet
1824  !! .
1825  !! @author C. Winteler
1826  subroutine kodtakdist(fissrate_inout)
1827  use parameter_class, only: unit
1828  use global_class, only: isotope,ineu
1829  use benam_class, only: minmax,findaz
1831  implicit none
1832  type(fissionrate_type),intent(inout) :: fissrate_inout !< Temporary fission rate
1833  ! Internal variables
1834  integer :: mass !< Mass number of parent nucleus
1835  integer :: pnr !< Proton number of parent nucleus
1836  integer :: nf !< number of fission channels
1837  integer :: fiss_mode !< Fission mode (1: n-induced, 2: spontaneous, 3: beta-delayed fission)
1838  integer :: af,zf !< mass and proton number of fissioning nucleus
1839  integer :: a,z !< mass and proton number of fragment
1840  integer :: a1,a2,z1,z2
1841  real(r_kind) :: paz,ptot
1842  real(r_kind) :: za,al,ah,cz,ca
1843  real(r_kind),parameter :: lim = 1.d-6 !< Limit for fragments to be taken into account
1844  integer :: nemiss !< Neutrons missing in case of fragment beyond dripline for one fragment
1845  integer :: alloc_stat !< Allocation status
1846  integer :: i !< Loop variable
1847  integer :: dimens !< Dimension of the parts in fission rate (amount fragments + 1)
1848  integer :: neutronflag !< Flag to check whether any fragment is beyond dripline
1849  real(r_kind) :: nemiss_total !< Total number of neutrons emitted (sum of nemiss)
1850  integer :: additional_neutrons !< Neutrons that may be emitted by beta decay
1851  real(r_kind) :: qval !< Q-value of the fission reaction
1852 
1853  info_entry("kodtakdist")
1854 
1855  ! Set neutron flag to zero
1856  neutronflag = 0
1857  additional_neutrons = 0d0
1858 
1859  ! Get the mass and proton number
1860  mass = isotope(fissrate_inout%fissnuc_index)%mass
1861  pnr = isotope(fissrate_inout%fissnuc_index)%p_nr
1862 
1863  !neutron induced fission
1864  if (fissrate_inout%reac_type .eq. rrt_nf) then
1865  fiss_mode = 1
1866  af = mass+1
1867  zf = pnr
1868  !spontaneous fission
1869  else if (fissrate_inout%reac_type .eq. rrt_sf) then
1870  fiss_mode = 2
1871  af = mass
1872  zf = pnr
1873  !beta delayed fission
1874  else if (fissrate_inout%reac_type .eq. rrt_bf) then
1875  fiss_mode = 3
1876  ! Take care of beta-delayed neutron emission fission
1877  ! The amount of neutrons is stored in the second parameter
1878  ! of the fission rate
1879  additional_neutrons = fissrate_inout%released_n
1880  af = mass-additional_neutrons
1881  zf = pnr + 1
1882  end if
1883 
1884  al = 0.85d0*af - 104.98d0
1885  ah = 0.15d0*af + 103.87d0
1886  cz = 0.8d0
1887  ca = 78.d0
1888 
1889  nf = 0
1890  ptot = 0.d0
1891 
1892  do a = 1,af
1893  za = zf*(a+0.6d0)/af
1894  do z=1,zf
1895  paz = 0.5d0*dexp(-((z-za)**2)/cz) * &
1896  ((dexp(-((a-al)**2)/ca))+(dexp(-((a-ah)**2)/ca)))/ &
1897  (unit%pi*dsqrt(cz*ca))
1898  if(paz.ge.lim) then
1899  nf = nf+1
1900  ptot = ptot + paz
1901  a1 = a
1902  z1 = z
1903  a2 = af - a1
1904  z2 = zf - z1
1905 
1906  ! check if any neutrons are emitted (needed to determine fissrate()%dimens)
1907  if (a1.gt.minmax(z1,2)) then
1908  neutronflag = 1
1909  else if (a1.lt.minmax(z1,1)) then
1910  if (verbose_level .ge. 2) then
1911  print *, 'Warning in kodtakdist, fragment was lighter than included nuclei'
1912  print *, 'Fragment was: ', a1, z1
1913  end if
1914  end if
1915 
1916  if (a2.gt.minmax(z2,2)) then
1917  neutronflag = 1
1918  else if (a2.lt.minmax(z2,1)) then
1919  if (verbose_level .ge. 2) then
1920  print *, 'Warning in kodtakdist, fragment was lighter than included nuclei'
1921  print *, 'Fragment was: ', a2, z2
1922  end if
1923  end if
1924 
1925  end if
1926  end do
1927  end do
1928 
1929  fissrate_inout%channels = nf
1930  fissrate_inout%fissnuc_index = findaz(mass,pnr)
1931  fissrate_inout%mode = fiss_mode
1932 
1933  ! Allocate and complain in case of failure
1934  allocate(fissrate_inout%channelprob(nf),fissrate_inout%q_value(nf), &
1935  stat=alloc_stat)
1936  if (alloc_stat .ne. 0) call raise_exception("Allocation of fission rate prob. arrays failed.",&
1937  "kodtakdist",190001)
1938 
1939 
1940  select case (fiss_mode)
1941  case(1)
1942  dimens = 2 * fissrate_inout%channels + 2
1943  fissrate_inout%dimens = dimens
1944  allocate(fissrate_inout%fissparts(dimens),fissrate_inout%ch_amount(dimens),&
1945  stat=alloc_stat)
1946  if (alloc_stat .ne. 0) call raise_exception("Allocation of fission rate parts failed.",&
1947  "kodtakdist",190001)
1948 
1949  fissrate_inout%fissparts(1) = ineu
1950  fissrate_inout%fissparts(2) = fissrate_inout%fissnuc_index
1951  fissrate_inout%ch_amount(2) = -1.d0
1952  case(2:)
1953  dimens = 2 * fissrate_inout%channels + 1
1954  if (neutronflag.eq.1) dimens = dimens+1
1955  fissrate_inout%dimens = dimens
1956  allocate(fissrate_inout%fissparts(dimens),fissrate_inout%ch_amount(dimens),&
1957  stat=alloc_stat)
1958  if (alloc_stat .ne. 0) call raise_exception("Allocation of fission rate parts failed.",&
1959  "kodtakdist",190001)
1960  fissrate_inout%ch_amount(1) = -1.d0
1961  fissrate_inout%fissparts(1) = fissrate_inout%fissnuc_index
1962  end select
1963 
1964  ! Keep track of total number of neutrons emitted
1965  nemiss_total = float(additional_neutrons)
1966  !---- start writing the fission rates to rrate
1967  i = 0
1968  do z = 1,zf
1969  massloop: do a = 1,af
1970  za = zf*(a+0.6d0)/af
1971  paz = 0.5d0*dexp(-((z-za)**2)/cz) * &
1972  ((dexp(-((a-al)**2)/ca))+(dexp(-((a-ah)**2)/ca)))/ &
1973  (unit%pi*dsqrt(cz*ca))
1974  if(paz.lt.lim) cycle massloop
1975  i = i+1
1976  paz = paz/ptot ! normalise probabilities to 1
1977  a1 = a
1978  z1 = z
1979  a2 = af - a1
1980  z2 = zf - z1
1981 
1982  nemiss = 0
1983 
1984  if (a1.gt.minmax(z1,2)) then
1985  nemiss = a1 - minmax(z1,2)
1986  a1 = minmax(z1,2)
1987 
1988  else if (a1.lt.minmax(z1,1)) then
1989  if (verbose_level .ge. 2) then
1990  print *, 'Warning in kodtakdist, fragment was lighter than included nuclei'
1991  print *, 'Fragment was: ', a1, z1
1992  end if
1993  end if
1994 
1995  if (a2.gt.minmax(z2,2)) then
1996  nemiss = nemiss + a2 - minmax(z2,2)
1997  a2 = minmax(z2,2)
1998  else if (a2.lt.minmax(z2,1)) then
1999  ! TODO: Maybe borrow neutrons here?
2000  if (verbose_level .ge. 2) then
2001  print *, 'Warning in kodtakdist, fragment was lighter than included nuclei'
2002  print *, 'Fragment was: ', a2, z2
2003  end if
2004  end if
2005 
2006  ! Count the total amount of neutrons emitted
2007  nemiss_total = nemiss_total+ float(nemiss)*paz
2008  ! fill in fissrate()%fissparts(); same as rrate()%parts(1:6), but with individual array sizes
2009  select case(fiss_mode)
2010  case(1)
2011  fissrate_inout%channelprob(i) = paz
2012  fissrate_inout%fissparts(i+2) = findaz(a1,z1)
2013  ! Check that the nucleus was really found!
2014  if (fissrate_inout%fissparts(i+2) .eq. -1) then
2015  call raise_exception("Fragment was not found in nucleus library, does your library have 'holes'? "//&
2016  "(A="//int_to_str(a1)//", Z="//int_to_str(z1)//")",&
2017  "kodtakdist",190013)
2018  end if
2019 
2020  fissrate_inout%ch_amount(i+2) = fissrate_inout%channelprob(i) ! rate at which fragment is produced per destroyed parent nucleus
2021  fissrate_inout%fissparts(i+2+fissrate_inout%channels) = findaz(a2,z2)
2022  ! Check that the nucleus was really found!
2023  if (fissrate_inout%fissparts(i+2+fissrate_inout%channels) .eq. -1) then
2024  call raise_exception("Fragment was not found in nucleus library, does your library have 'holes'? "//&
2025  "(A="//int_to_str(a2)//", Z="//int_to_str(z2)//")",&
2026  "kodtakdist",190013)
2027  end if
2028  fissrate_inout%ch_amount(i+2+fissrate_inout%channels) = fissrate_inout%channelprob(i)
2029  fissrate_inout%q_value(i) = isotope(fissrate_inout%fissnuc_index)%mass_exc - &
2030  (nemiss - 1) * isotope(ineu)%mass_exc - & ! (nemiss - 1) to account for one neutron that is destroyed
2031  isotope(fissrate_inout%fissparts(i+2))%mass_exc - &
2032  isotope(fissrate_inout%fissparts(i+2+fissrate_inout%channels))%mass_exc
2033  case(2:)
2034  fissrate_inout%channelprob(i) = paz
2035  fissrate_inout%fissparts(i+1) = findaz(a1,z1)
2036  ! Check that the nucleus was really found!
2037  if (fissrate_inout%fissparts(i+1) .eq. -1) then
2038  call raise_exception("Fragment was not found in nucleus library, does your library have 'holes'?",&
2039  "kodtakdist",190013)
2040  end if
2041 
2042  fissrate_inout%ch_amount(i+1) = fissrate_inout%channelprob(i) ! rate at which fragment is produced per destroyed parent nucleus
2043  fissrate_inout%fissparts(i+1+fissrate_inout%channels) = findaz(a2,z2)
2044  ! Check that the nucleus was really found!
2045  if (fissrate_inout%fissparts(i+1+fissrate_inout%channels) .eq. -1) then
2046  call raise_exception("Fragment was not found in nucleus library, does your library have 'holes'?",&
2047  "kodtakdist",190013)
2048  end if
2049 
2050  fissrate_inout%ch_amount(i+1+fissrate_inout%channels) = fissrate_inout%channelprob(i)
2051  if (neutronflag.eq.1) then
2052  fissrate_inout%fissparts(dimens) = ineu
2053  fissrate_inout%ch_amount(dimens) = real(nemiss)
2054  end if
2055  fissrate_inout%q_value(i) = isotope(fissrate_inout%fissnuc_index)%mass_exc - &
2056  nemiss * isotope(ineu)%mass_exc - &
2057  isotope(fissrate_inout%fissparts(i+1))%mass_exc - &
2058  isotope(fissrate_inout%fissparts(i+1+fissrate_inout%channels))%mass_exc
2059  end select
2060  end do massloop
2061  end do
2062 
2063  ! Set the neutrons correctly
2064  select case(fiss_mode)
2065  case(1)
2066  fissrate_inout%ch_amount(1) = nemiss_total - 1.d0
2067  case(2:)
2068  if (neutronflag.eq.1) then
2069  fissrate_inout%fissparts(dimens) = ineu
2070  fissrate_inout%ch_amount(dimens) = nemiss_total
2071  end if
2072  end select
2073 
2074  qval = 0.d0
2075  do i=1,fissrate_inout%dimens
2076  qval = qval - isotope(fissrate_inout%fissparts(i))%mass_exc*fissrate_inout%ch_amount(i) ! weighted average Q-value
2077  end do
2078  fissrate_inout%averageQ = qval
2079 
2080  nufiss = nufiss + int(nemiss_total)
2081 
2082  info_exit("kodtakdist")
2083  return
2084 
2085  end subroutine kodtakdist
2086 
2087 
2088 
2089  !> Read the fission distribution from a file
2090  !!
2091  !! The fission distribution is read from a file and stored in the
2092  !! fragment_array. The default file contains the fission distribution
2093  !! according to [Mumpower et al. (2020)](https://ui.adsabs.harvard.edu/abs/2020PhRvC.101e4607M/abstract).
2094  !! However, the format has been modified compared to the original publication.
2095  !! An example of the file looks like:
2096  !! \file{
2097  !! 93 220 574
2098  !! 0 0 0
2099  !! 019 048 1.96607e-06
2100  !! 020 048 2.24289e-05
2101  !! 020 049 1.88041e-05
2102  !! 020 050 7.65944e-06
2103  !! }
2104  !! The first line contains the Z, A and number of the fissioning nucleus.
2105  !! The next line is a dummy line. Then followed by Z, A and the yield Y(Z,A).
2106  !! Note that \f$\sum Y(Z,A) \cdot A \f$ should be the mass number of the
2107  !! parent nucleus. Furthermore, \f$\sum Y(Z,A) = 2\f$ in case 2 fragments are created.
2108  !! A large part of the subroutine deals with the problem what to do if the fragment
2109  !! is not included in the network. In case of a nucleus heavier than the ones contained in
2110  !! the network, the fragment is split into a lighter isotope plus neutrons.
2111  !! For lighter nuclei where neutrons would be necessary to create a nucleus that is in the network,
2112  !! the neutrons are "borrowed" from the heaviest fragment. This latter case should almost never happen.
2113  !! If the remapping fails for some reason, an error will be raised.
2114  !!
2115  !! @see [Mumpower et al. (2020)](https://ui.adsabs.harvard.edu/abs/2020PhRvC.101e4607M/abstract).
2116  !!
2117  !! @author M. Reichert
2118  !! @date 14.02.23
2119  subroutine read_fiss_fragment_file(file,fragment_array)
2121  use benam_class, only: findaz, minmax
2122  use error_msg_class, only: raise_exception
2123  implicit none
2124 
2125  character(len=*), intent(in) :: file !< File path to fission fragment file
2126  type(fragtype), dimension(:), allocatable, intent(out) :: fragment_array !< Array of fission fragments
2127 
2128  integer :: file_id !< Integer id of the file
2129  integer :: frag_counter !< Count the number of distributions
2130  integer :: zf !< Z fission nucleus
2131  integer :: af !< A fission nucleus
2132  integer :: nof !< Number of fragments
2133  integer :: i,j,m !< Loop variable
2134  integer :: read_stat !< Status variable
2135  integer :: astat !< Allocation status
2136  integer :: nem !< Number of emitted neutrons
2137  integer :: bn !< Number of borrowed neutrons
2138  real(r_kind) :: Ynem !< Yield of emitted neutrons
2139  real(r_kind) :: helper !< Yield of emitted neutrons
2140  integer :: mina,maxa !< Minimum and maximum A
2141 
2142 
2143  info_entry("read_fiss_fragment_file")
2144 
2145  ! Open the file
2146  file_id =open_infile(file)
2147 
2148  ! Count the number of relevant
2149  ! fragment distributions
2150  frag_counter = 0
2151  do
2152  read(file_id,*,iostat=read_stat) zf, af, nof
2153  if (read_stat .ne. 0) exit ! end of file
2154  read(file_id,*)
2155 
2156  do i=1,nof
2157  read(file_id,*)
2158  end do
2159 
2160  if (((findaz(af,zf) .ne. -1) .or. &
2161  (findaz(af-1,zf) .ne. -1) .or. &
2162  (findaz(af,zf-1) .ne. -1))) then
2163  frag_counter = frag_counter + 1
2164  end if
2165  end do
2166 
2167  ! Allocate the distribtions
2168  allocate(fragment_array(frag_counter),stat=astat)
2169  if (astat .ne. 0) then
2170  call raise_exception("Could not allocate 'fragment_array'.", &
2171  "read_fiss_fragment_file",190001)
2172  end if
2173 
2174  ! Read the file again
2175  rewind(file_id)
2176 
2177  fragloop: do i = 1,frag_counter
2178  read(file_id,*,iostat=read_stat) zf, af, nof
2179  read(file_id,*)
2180 
2181 
2182  ! Nucleus not included nor the decay of
2183  ! a nucleus in the network, nor reachable
2184  ! via a neutron capture
2185  if (.not. ((findaz(af,zf) .ne. -1) .or. &
2186  (findaz(af-1,zf) .ne. -1) .or. &
2187  (findaz(af,zf-1) .ne. -1))) then
2188  do j=1,nof
2189  read(file_id,*)
2190  end do
2191  cycle fragloop
2192  end if
2193 
2194 
2195  ! Allocate the fragments
2196  allocate(fragment_array(i)%net_idx(nof),&
2197  fragment_array(i)%Z(nof),&
2198  fragment_array(i)%A(nof),&
2199  fragment_array(i)%Y(nof),&
2200  stat=astat)
2201  if (astat .ne. 0) then
2202  call raise_exception("Could not allocate fragment properties.", &
2203  "read_fiss_fragment_file",190001)
2204  end if
2205 
2206  ! Store the properties
2207  fragment_array(i)%nr_frags = nof
2208  fragment_array(i)%Zp = zf
2209  fragment_array(i)%Ap = af
2210 
2211  ! Count the neutrons emitted
2212  nem = 0
2213  ynem = 0
2214 
2215  ! Read the fragments
2216  do j=1,nof
2217  read(file_id,*) fragment_array(i)%Z(j),fragment_array(i)%A(j),fragment_array(i)%Y(j)
2218  fragment_array(i)%net_idx(j) =findaz(fragment_array(i)%A(j),fragment_array(i)%Z(j))
2219 
2220  ! What to do when the fragment is not in the network?
2221  ! Ideally make neutrons until the fragment is in the network
2222  if (fragment_array(i)%net_idx(j) .eq. -1) then
2223  ! Tell the world that this is done
2224  if (verbose_level .ge. 2) then
2225  write(*,*) "Fragment ",fragment_array(i)%Z(j),fragment_array(i)%A(j),&
2226  " not in network. Splitting into lighter isotope and neutrons."
2227  end if
2228 
2229  mina = minmax(fragment_array(i)%Z(j),1)
2230  maxa = minmax(fragment_array(i)%Z(j),2)
2231  if (fragment_array(i)%A(j) .gt. maxa) then
2232  ! Rest goes to neutrons
2233  nem = nem + (fragment_array(i)%A(j)-maxa)
2234  ! Yield of neutrons is the amount scaled with the
2235  ! yield of the not included fragment
2236  ! Y * A should be conserved
2237  ynem = ynem + fragment_array(i)%Y(j)*(fragment_array(i)%A(j)-maxa)
2238  fragment_array(i)%Y(j) = fragment_array(i)%Y(j) ! THis stays conserved
2239  fragment_array(i)%A(j) = maxa
2240  fragment_array(i)%net_idx(j) = findaz(maxa,fragment_array(i)%Z(j))
2241  end if
2242  end if
2243  end do
2244 
2245  ! Check how many neutrons are avaiable and give it to the missing
2246  ! fragments
2247  do j=1,nof
2248  if (fragment_array(i)%net_idx(j) .eq. -1) then
2249 
2250  ! Say something if verbose
2251  if (verbose_level .ge. 2) then
2252  write(*,*) "Fission fragment Z=",fragment_array(i)%Z(j), &
2253  "A=",fragment_array(i)%A(j), &
2254  "Y=",fragment_array(i)%Y(j), &
2255  "not in network, trying to remap."
2256  end if
2257 
2258  ! Needed neutrons are:
2259  mina = minmax(fragment_array(i)%Z(j),1)
2260  bn = (mina-fragment_array(i)%A(j)) ! Amount of missing neutrons
2261 
2262  ! Borrow neutrons, take heaviest fragment that has enough
2263  ! abundance
2264  do m=nof,1,-1
2265  if ((fragment_array(i)%net_idx(m) .ne. -1) .and. &
2266  ((bn*fragment_array(i)%Y(j) / fragment_array(i)%A(m)) .lt. fragment_array(i)%Y(m))) then
2267  exit
2268  end if
2269  end do
2270 
2271  ! None found? Too bad, make an error
2272  if (m .eq. 0) then
2273  call raise_exception("Could not remap fission fragments, "// &
2274  "consider to include more nuclides in the network. "//&
2275  new_line('A')//&
2276  "Nucleus Z="//trim(int_to_str(fragment_array(i)%Z(j)))// &
2277  ", A="//trim(int_to_str(fragment_array(i)%A(j)))// &
2278  " was not included. "//new_line('A')//&
2279  "Need "//int_to_str(bn)//" neutrons.",&
2280  "read_fiss_fragment_file",&
2281  190008)
2282  end if
2283 
2284  ! Otherwise borrow the neutrons from the heavier fragment
2285  fragment_array(i)%Y(m) = fragment_array(i)%Y(m) - (bn*fragment_array(i)%Y(j) / fragment_array(i)%A(m))
2286  fragment_array(i)%Y(j) = fragment_array(i)%Y(j)
2287  ! Give the fragment a new identity
2288  fragment_array(i)%A(j) = mina
2289  fragment_array(i)%net_idx(j) = findaz(mina,fragment_array(i)%Z(j))
2290  end if
2291  end do
2292 
2293  ! Proof the yields, this should be Af!
2294  helper = ynem
2295  do j=1,nof
2296  helper = helper+fragment_array(i)%A(j)*fragment_array(i)%Y(j)
2297  end do
2298  if (abs(helper - float(fragment_array(i)%Ap)) .gt. 1e-2) then
2299  call raise_exception("Yield of fission fragments is not conserved. "// &
2300  "This happened for fragment distribution of "// &
2301  "Nucleus Z="//trim(int_to_str(fragment_array(i)%Zp))// &
2302  ", A="//trim(int_to_str(fragment_array(i)%Ap))// &
2303  ". Ensure that the sum is 2.","read_fiss_fragment_file",&
2304  190009)
2305  end if
2306 
2307  fragment_array(i)%neutrons = nem
2308  fragment_array(i)%Yn = ynem
2309  nufiss = nufiss + ynem
2310  end do fragloop
2311 
2312  ! Close the file
2313  close(file_id)
2314 
2315  info_exit("read_fiss_fragment_file")
2316 
2317  end subroutine read_fiss_fragment_file
2318 
2319 
2320  !> Fill the rates with the correct fragments
2321  !!
2322  !! This subroutine fills the rates with the correct fragments
2323  !! from [Mumpower et al. (2020)](https://ui.adsabs.harvard.edu/abs/2020PhRvC.101e4607M/abstract).
2324  !! This subroutine makes use of the fragtype struct.
2325  !! If no fragment is found, the fragment distribution is set to the
2326  !! one given by missing_frags.
2327  !! For beta-delayed fission the fissioning nucleus is the one of (Z+1, A)
2328  !! for neutron induced fission the nucleus with (Z, A+1), and for spontanous
2329  !! fission (not used within fissflag = 3) it is (Z,A).
2330  !!
2331  !! @author M. Reichert
2332  !! @date 14.02.2023
2333  subroutine file_fiss_frag(fissrate_inout,missing_frags)
2334  use global_class, only: isotope,ineu
2335  use benam_class, only: minmax,findaz
2337  implicit none
2338  type(fissionrate_type),intent(inout) :: fissrate_inout !< Temporary fission rate
2339  integer,optional, intent(in) :: missing_frags!< Distribution to use if not in the file (1:panov, 2:kodama)
2340  ! Internal variables
2341  integer :: mass !< Mass number of reacting nucleus
2342  integer :: pnr !< Proton number of reacting nucleus
2343  real(r_kind) :: qval !< Q-value of the reaction
2344  integer :: nfrac_distr !< Number of fragment distributions
2345  integer :: i !< Loop variable
2346  integer :: astat !< Allocation status variable
2347  integer :: ind_parent !< Index of reacting nucleus
2348  integer :: afiss !< A of fissioning nucleus
2349  integer :: zfiss !< Z of fissioning nucleus
2350  type(fragtype) :: fissfrags_tmp!< Temporary fragment distribution
2351  logical :: found !< Helper variable flag if the fragment
2352  ! distribution was found
2353  integer :: additional_neutrons !< Neutrons to be added for beta-delayed fission
2354  integer :: missing_frags_internal !< Internal variable for missing fragments
2355 
2356  info_entry("file_fiss_frag")
2357 
2358  ! Set the default missing fragments
2359  if (.not. present(missing_frags)) then
2360  missing_frags_internal = 0
2361  else
2362  missing_frags_internal = missing_frags
2363  end if
2364 
2365  ! Get the mass and proton number
2366  mass = isotope(fissrate_inout%fissnuc_index)%mass
2367  pnr = isotope(fissrate_inout%fissnuc_index)%p_nr
2368 
2369  ! Loop through the fission fragments
2370  if (fissrate_inout%reac_type .eq. rrt_bf) then
2371  nfrac_distr = size(fissfrags_beta_delayed)
2372  else if (fissrate_inout%reac_type .eq. rrt_nf) then
2373  nfrac_distr = size(fissfrags_n_induced)
2374  else if (fissrate_inout%reac_type .eq. rrt_sf) then
2375  nfrac_distr = size(fissfrags_spontaneous)
2376  end if
2377 
2378  ! Assume first that there are no additional neutrons
2379  additional_neutrons = fissrate_inout%released_n
2380 
2381  ! Check which nuclei is fissioning
2382  zfiss = pnr
2383  afiss = mass - additional_neutrons
2384 
2385  ! Z-1 for beta delayed fission
2386  if (fissrate_inout%reac_type .eq. rrt_bf) then
2387  zfiss = zfiss+1
2388  ! A+1 for neutron induced fission
2389  elseif (fissrate_inout%reac_type .eq. rrt_nf) then
2390  afiss = afiss+1
2391  end if
2392  ! Don't change anything for spontaneous fission
2393 
2394  ! Find relevant indices
2395  ind_parent = findaz(mass,pnr)
2396  fissrate_inout%fissnuc_index = ind_parent
2397  found = .false.
2398  do i=1,nfrac_distr
2399 
2400  if (fissrate_inout%reac_type .eq. rrt_bf) then
2401  fissfrags_tmp = fissfrags_beta_delayed(i)
2402  else if (fissrate_inout%reac_type .eq. rrt_nf) then
2403  fissfrags_tmp = fissfrags_n_induced(i)
2404  else if (fissrate_inout%reac_type .eq. rrt_sf) then
2405  fissfrags_tmp = fissfrags_spontaneous(i)
2406  end if
2407 
2408  if ((fissfrags_tmp%Zp .eq. zfiss) .and. (fissfrags_tmp%Ap .eq. afiss)) then
2409  found = .true.
2410  exit
2411  end if
2412  end do
2413 
2414  ! Not found, check what to do
2415  if (.not. found) then
2416 
2417  ! Raise an error when there is no alternative distribution specified
2418  if (missing_frags_internal .eq. 0) then
2419  call raise_exception("No fragment distribution found for Z="//int_to_str(zfiss)//", "//&
2420  "A="//int_to_str(afiss)//" in fragment file. Either add it to the fission fragment file or set "//&
2421  "'fission_frag_missing' to 1 (Panov et al. 2001) or 2 (Kodama & Takahashi 1975)",&
2422  "file_fiss_frag", 190016)
2423  end if
2424 
2425  ! Output this information
2426  if (verbose_level .ge.2) then
2427  write(*,*) "No fragment distribution found for Z=",zfiss,", A=",afiss
2428  if (missing_frags_internal .eq. 1) then
2429  write(*,*) "Using Panov 2001 instead."
2430  else if (missing_frags_internal .eq. 2) then
2431  write(*,*) "Using Kodama & Takahashi 1975 instead."
2432  end if
2433  end if
2434 
2435  if (missing_frags_internal .eq. 1) then
2436  ! Get panov distribution
2437  call fiss_dist(fissrate_inout)
2438  else if (missing_frags_internal .eq. 2) then
2439  ! Get kodama distribution
2440  call kodtakdist(fissrate_inout)
2441  else
2442  call raise_exception("Unknown missing fragment distribution. "//&
2443  "Got "//int_to_str(missing_frags_internal)//".",&
2444  "file_fiss_frag",190015)
2445  end if
2446  else
2447  ! Found, do mumpower
2448  if (fissrate_inout%reac_type .eq. rrt_bf) then
2449  fissfrags_tmp = fissfrags_beta_delayed(i)
2450  else if (fissrate_inout%reac_type .eq. rrt_nf) then
2451  fissfrags_tmp = fissfrags_n_induced(i)
2452  else if (fissrate_inout%reac_type .eq. rrt_sf) then
2453  fissfrags_tmp = fissfrags_spontaneous(i)
2454  end if
2455 
2456  ! Allocate the parts, the first two are the fissioning nucleus and neutrons
2457  allocate(fissrate_inout%fissparts(fissfrags_tmp%nr_frags+2),&
2458  fissrate_inout%ch_amount(fissfrags_tmp%nr_frags+2),&
2459  stat=astat)
2460  ! Complain if not possible
2461  if (astat .ne. 0) then
2462  call raise_exception("Could not allocate memory for fission fragments.",&
2463  "file_fiss_frag",190001)
2464  end if
2465 
2466  ! Set up the number of participants,
2467  ! Number of fragments + parent nucleus + neutrons
2468  fissrate_inout%dimens = fissfrags_tmp%nr_frags+2
2469 
2470  ! Use the correct order of the reacting nuclei
2471  if (.not. (fissrate_inout%reac_type .eq. rrt_nf)) then
2472  fissrate_inout%fissparts(1) = ind_parent
2473  fissrate_inout%fissparts(2) = ineu
2474  fissrate_inout%ch_amount(1) = -1d0
2475  fissrate_inout%ch_amount(2) = fissfrags_tmp%Yn+float(additional_neutrons)
2476  else
2477  fissrate_inout%fissparts(1) = ineu
2478  fissrate_inout%fissparts(2) = ind_parent
2479  fissrate_inout%ch_amount(1) = fissfrags_tmp%Yn-1d0+float(additional_neutrons)
2480  fissrate_inout%ch_amount(2) = -1d0
2481  end if
2482 
2483 
2484  do i=1,fissfrags_tmp%nr_frags
2485  fissrate_inout%fissparts(i+2) = fissfrags_tmp%net_idx(i)
2486  fissrate_inout%fissparts(i+2) = fissfrags_tmp%net_idx(i)
2487  fissrate_inout%ch_amount(i+2) = fissfrags_tmp%Y(i)
2488  end do
2489  end if
2490 
2491  ! Calculate the Q-value
2492  qval = 0.d0
2493  qvalue: do i=1,fissrate_inout%dimens
2494  qval = qval - isotope(fissrate_inout%fissparts(i))%mass_exc*fissrate_inout%ch_amount(i) ! weighted average Q-value
2495  end do qvalue
2496  fissrate_inout%averageQ = qval
2497 
2498  info_exit("file_fiss_frag")
2499 
2500  end subroutine file_fiss_frag
2501 
2502 
2503 
2504  !> Calculates the (neutron-induced) fission fragment mass distribution according to
2505  !! the ABLA07 model: Kelic, Ricciardi, & Schmidt (arXiv:0906.4193)
2506  !!
2507  !! @see [Kelic et al. 2009](https://ui.adsabs.harvard.edu/abs/2009arXiv0906.4193K/abstract)
2508  subroutine abla_nfiss(pos,mass,pnr,neutfission,qval)
2509  use global_class, only: isotope,ineu
2510  use benam_class, only: minmax,findaz
2511  use error_msg_class, only: raise_exception
2512 
2513  integer,intent(in) :: mass,pnr
2514  integer, intent(in) :: pos
2515  real(r_kind),intent(out) :: qval
2516  integer :: fmass,fpnr
2517  integer :: read_stat
2518  integer :: alloc_stat
2519  integer :: af,zf
2520  integer,dimension(:),allocatable :: a1,z1
2521  integer,dimension(:),allocatable :: a2,z2
2522  integer,intent(in) :: neutfission
2523  integer :: i,nof, nemiss
2524  integer :: dimens
2525  real(r_kind) :: nrate,npre,nafter
2526  real(r_kind),dimension(:),allocatable :: paz
2527  real(r_kind) :: psum
2528  real(r_kind) :: n_mult
2529  character(2) :: dummy
2530 
2531  info_entry("abla_nfiss")
2532 
2533  nemiss = 0
2534  do
2535  read(neutfission,"(I3,1X,I3,5X,I3)",iostat=read_stat) zf, af, nof
2536  if (read_stat /= 0) then !< this condition can happen and the code should be allowed to continue (stop -> exit)
2537  ! print*, 'fission reaction not found in abla_nfiss: ', fmass, fpnr
2538  exit ! parent nucleus not found in ABLA table -> continue with old Panov model (below, corresponds to fissflag=1 for this nucleus)
2539  endif
2540 
2541  read(neutfission,*,iostat=read_stat) nrate, npre, nafter
2542  if (read_stat /= 0) then
2543  call raise_exception("Could not read fission file.","abla_nfiss",190004)
2544  endif
2545 
2546  if ((af.eq.mass).and.(zf.eq.pnr).and.(nof.ne.0)) then
2547 
2548  n_mult = npre + nafter
2549  nemiss = int(dnint(n_mult))
2550  nufiss = nufiss + nemiss
2551 
2552  allocate(a1(nof),stat=alloc_stat)
2553  allocate(a2(nof),stat=alloc_stat)
2554  allocate(z1(nof),stat=alloc_stat)
2555  allocate(z2(nof),stat=alloc_stat)
2556  allocate(paz(nof),stat=alloc_stat)
2557 
2558  fissrate(pos)%channels = nof
2559  fissrate(pos)%fissnuc_index = findaz(af,zf)
2560  fissrate(pos)%mode = 1
2561 
2562  allocate(fissrate(pos)%channelprob(nof))
2563  allocate(fissrate(pos)%q_value(nof))
2564  dimens = 2 * fissrate(pos)%channels + 2
2565  fissrate(pos)%dimens = dimens
2566  allocate(fissrate(pos)%fissparts(dimens))
2567  allocate(fissrate(pos)%ch_amount(dimens))
2568  fissrate(pos)%fissparts(1) = ineu
2569  fissrate(pos)%ch_amount(1) = float(nemiss) - 1.d0
2570  fissrate(pos)%fissparts(2) = fissrate(pos)%fissnuc_index
2571  fissrate(pos)%ch_amount(2) = -1.d0
2572 
2573 
2574  psum = 0.d0
2575  inner: do i=1,nof
2576  read(neutfission,*,iostat=read_stat) &
2577  & a1(i),z1(i),paz(i)
2578  if (read_stat /= 0) then
2579  call raise_exception("Could not read fission file.","abla_nfiss",190004)
2580  endif
2581  ! check if fragment 1 is in the network
2582  if (a1(i).gt.minmax(z1(i),2)) then
2583  nemiss = nemiss + a1(i) - minmax(z1(i),2)
2584  a1(i) = minmax(z1(i),2)
2585  else if (a1(i).lt.minmax(z1(i),1)) then
2586  print *, 'a1 error:'
2587  end if
2588  a2(i) = mass + 1 - a1(i) - nemiss
2589  z2(i) = pnr - z1(i)
2590  psum = psum + paz(i)
2591 
2592  ! check if fragment 2 is in the network
2593  if (a2(i).gt.minmax(z2(i),2)) then
2594  nemiss = nemiss + a2(i) - minmax(z2(i),2)
2595  a2(i) = minmax(z2(i),2)
2596  else if (a2(i).lt.minmax(z2(i),1)) then
2597  print *, 'a2 error:', a2(i), z2(i), a1(i), z1(i), mass, pnr
2598  end if
2599  end do inner
2600 
2601  again: do i=1,nof
2602  paz(i) = paz(i)/psum ! normalize probabilities
2603  fissrate(pos)%channelprob(i) = paz(i)
2604 
2605  !---------rate for neutron-induced fission
2606  fissrate(pos)%fissparts(i+2) = findaz(a1(i),z1(i))
2607  fissrate(pos)%ch_amount(i+2) = fissrate(pos)%channelprob(i) ! rate at which fragment is produced per destroyed parent nucleus
2608  fissrate(pos)%fissparts(i+2+fissrate(pos)%channels) = findaz(a2(i),z2(i))
2609  fissrate(pos)%ch_amount(i+2+fissrate(pos)%channels) = fissrate(pos)%channelprob(i)
2610  fissrate(pos)%q_value(i) = (1-nemiss) * isotope(ineu)%mass_exc + & ! one neutron is destroyed, a number equal to nemiss are produced
2611  isotope(fissrate(pos)%fissnuc_index)%mass_exc - &
2612  isotope(fissrate(pos)%fissparts(i+2))%mass_exc - &
2613  isotope(fissrate(pos)%fissparts(i+2+fissrate(pos)%channels))%mass_exc
2614 
2615  end do again
2616 
2617  qval = 0.d0
2618  qvalue: do i=1,fissrate(pos)%channels
2619  qval = qval + fissrate(pos)%q_value(i)*fissrate(pos)%channelprob(i) ! weighted average Q-value
2620  end do qvalue
2621  fissrate(pos)%averageQ = qval
2622 
2623 
2624  deallocate(a1)
2625  deallocate(a2)
2626  deallocate(z1)
2627  deallocate(z2)
2628  deallocate(paz)
2629 
2630  rewind(neutfission,iostat=read_stat)
2631  if (read_stat /= 0) then
2632  call raise_exception("Could not rewind fission file.",&
2633  "abla_nfiss",190005)
2634  endif
2635  return
2636 
2637  else ! read nof lines to advance to next parent nucleus in the table
2638 
2639  do i=1,nof
2640  read(neutfission,*,iostat=read_stat) dummy
2641  if (read_stat /= 0) then
2642  call raise_exception("Could not read fission file.",&
2643  "abla_nfiss",190004)
2644  endif
2645  end do
2646 
2647  end if
2648 
2649  end do
2650 
2651  !--------old Panov model for nuclei not found in ABLA fission table---------------
2652  allocate(a1(1),stat=alloc_stat)
2653  allocate(a2(1),stat=alloc_stat)
2654  allocate(z1(1),stat=alloc_stat)
2655  allocate(z2(1),stat=alloc_stat)
2656 
2657  fmass = mass + 1
2658  fpnr = pnr
2659  select case (fmass)
2660  case(:240) ! -> symmetric fission
2661  a2(1) = nint(fmass/2.d0)
2662  z2(1) = nint(fpnr/2.d0)
2663  a1(1) = fmass - a2(1)
2664  z1(1) = fpnr - z2(1)
2665  case(255:265) ! -> symmetric fission
2666  a2(1) = nint(fmass/2.d0)
2667  z2(1) = nint(fpnr/2.d0)
2668  a1(1) = fmass - a2(1)
2669  z1(1) = fpnr - z2(1)
2670  case default ! -> asymmetric fission
2671  a2(1) = 130
2672  z2(1) = nint(52.01 - (fpnr - 80)/1.d1)
2673  ! a2 = max(a2,af-a2) !modification to prevent very n-rich
2674  ! z2 = max(z2,zf-z2) !fragments
2675  a1(1) = fmass - a2(1)
2676  z1(1) = fpnr - z2(1)
2677  end select
2678 
2679  nemiss = 0
2680 
2681  if (a1(1).gt.minmax(z1(1),2)) then
2682  nemiss = a1(1) - minmax(z1(1),2)
2683  a1(1) = minmax(z1(1),2)
2684  else if (a1(1).lt.minmax(z1(1),1)) then
2685  print *, 'a1 error:'
2686  end if
2687 
2688  if (a2(1).gt.minmax(z2(1),2)) then
2689  nemiss = nemiss + a2(1) - minmax(z2(1),2)
2690  a2(1) = minmax(z2(1),2)
2691  else if (a2(1).lt.minmax(z2(1),1)) then
2692  print *, 'a2 error:'
2693  end if
2694  nufiss = nufiss + nemiss
2695 
2696  fissrate(pos)%dimens = 4
2697  allocate(fissrate(pos)%fissparts(fissrate(pos)%dimens))
2698  allocate(fissrate(pos)%ch_amount(fissrate(pos)%dimens))
2699  allocate(fissrate(pos)%channelprob(1))
2700  allocate(fissrate(pos)%q_value(1))
2701  fissrate(pos)%channels = 1
2702  fissrate(pos)%fissnuc_index = findaz(mass,pnr)
2703  fissrate(pos)%mode = 1
2704  fissrate(pos)%channelprob(1) = 1.d0
2705 
2706  fissrate(pos)%fissparts(1) = ineu
2707  fissrate(pos)%ch_amount(1) = float(nemiss) - 1.d0
2708  fissrate(pos)%fissparts(2) = fissrate(pos)%fissnuc_index
2709  fissrate(pos)%ch_amount(2) = -1.d0
2710 
2711  fissrate(pos)%fissparts(3) = findaz(a1(1),z1(1))
2712  fissrate(pos)%ch_amount(3) = 1.d0
2713  fissrate(pos)%fissparts(4) = findaz(a2(1),z2(1))
2714  fissrate(pos)%ch_amount(4) = 1.d0
2715 
2716  fissrate(pos)%averageQ = (1-nemiss) * isotope(ineu)%mass_exc + &
2717  isotope(fissrate(pos)%fissparts(2))%mass_exc - &
2718  isotope(fissrate(pos)%fissparts(3))%mass_exc - &
2719  isotope(fissrate(pos)%fissparts(4))%mass_exc
2720 
2721  deallocate(a1)
2722  deallocate(a2)
2723  deallocate(z1)
2724  deallocate(z2)
2725 
2726 
2727  rewind(neutfission,iostat=read_stat)
2728  if (read_stat /= 0) then
2729  call raise_exception("Could not rewind fission file.",&
2730  "abla_nfiss",190005)
2731  endif
2732 
2733  info_exit("abla_nfiss")
2734 
2735  end subroutine abla_nfiss
2736 
2737 
2738  !> Calculates the (beta-delayed and spontaneous) fission fragment mass distribution according to
2739  !! the ABLA07 model: Kelic, Ricciardi, & Schmidt (arXiv:0906.4193)
2740  !!
2741  !! @see [Kelic et al. 2009](https://ui.adsabs.harvard.edu/abs/2009arXiv0906.4193K/abstract)
2742  subroutine abla_betafiss(pos,mass,pnr,betafission,qval)
2743  use global_class, only: isotope,ineu
2744  use benam_class, only: minmax,findaz
2745  use error_msg_class, only: raise_exception
2746  integer,intent(in) :: mass,pnr
2747  integer, intent(in) :: pos
2748  real(r_kind),intent(out) :: qval
2749  integer :: nemiss,fmass,fpnr
2750  integer :: read_stat
2751  integer :: alloc_stat
2752  integer :: af,zf
2753  integer,dimension(:),allocatable :: a1,z1
2754  integer,dimension(:),allocatable :: a2,z2
2755  integer,intent(in) :: betafission
2756  integer :: i,nof
2757  integer :: dimens
2758  real(r_kind) :: nrate,npre,nafter
2759  real(r_kind),dimension(:),allocatable :: paz
2760  real(r_kind) :: psum
2761  real(r_kind) :: n_mult
2762  character(2) :: dummy
2763 
2764  info_entry("abla_betafiss")
2765 
2766  nemiss = 0
2767  do
2768  read(betafission,*, iostat = read_stat) zf, af, nof
2769 
2770  if (read_stat /= 0) exit
2771  read(betafission,*,iostat=read_stat) nrate, npre, nafter
2772  if (read_stat /= 0) then
2773  call raise_exception("Could not read fission file.",&
2774  "abla_betafiss",190006)
2775  endif
2776  if ((af.eq.mass).and.(zf.eq.pnr).and.(nof.ne.0)) then
2777 
2778  n_mult = npre + nafter
2779  nemiss = int(dnint(n_mult))
2780  nufiss = nufiss + nemiss
2781 
2782  allocate(a1(nof),stat=alloc_stat)
2783  allocate(a2(nof),stat=alloc_stat)
2784  allocate(z1(nof),stat=alloc_stat)
2785  allocate(z2(nof),stat=alloc_stat)
2786  allocate(paz(nof),stat=alloc_stat)
2787 
2788  fissrate(pos)%channels = nof
2789  fissrate(pos)%fissnuc_index = findaz(mass,pnr)
2790  if (fissrate(pos)%reac_type.eq. rrt_sf) then
2791  fissrate(pos)%mode = 2 ! spontaneous fission
2792  else if (fissrate(pos)%reac_type.eq. rrt_bf) then
2793  fissrate(pos)%mode = 3 ! beta-delayed fission
2794  end if
2795 
2796  allocate(fissrate(pos)%channelprob(nof))
2797  allocate(fissrate(pos)%q_value(nof))
2798  dimens = 2 * fissrate(pos)%channels + 1
2799  if (nemiss.gt.0) dimens = dimens + 1 ! neutrons which would otherwise not be part of the reaction
2800  fissrate(pos)%dimens = dimens
2801  allocate(fissrate(pos)%fissparts(dimens))
2802  allocate(fissrate(pos)%ch_amount(dimens))
2803  fissrate(pos)%fissparts(1) = fissrate(pos)%fissnuc_index
2804  fissrate(pos)%ch_amount(1) = -1.d0
2805 
2806  if (nemiss.gt.0) then
2807  fissrate(pos)%fissparts(dimens) = ineu
2808  fissrate(pos)%ch_amount(dimens) = float(nemiss)
2809  end if
2810 
2811  psum = 0.d0
2812  inner: do i=1,nof
2813  read(betafission,"(I3,1X,I3,2X,F6.4)",iostat=read_stat) a1(i), z1(i), paz(i)
2814  if (read_stat /= 0) then
2815  call raise_exception("Could not read fission file.",&
2816  "abla_betafiss",190006)
2817  endif
2818 
2819  ! check if fragment 1 is in the network
2820  if (a1(i).gt.minmax(z1(i),2)) then
2821  nemiss = nemiss + a1(i) - minmax(z1(i),2)
2822  a1(i) = minmax(z1(i),2)
2823  else if (a1(i).lt.minmax(z1(i),1)) then
2824  print *, 'a1 error:'
2825  end if
2826 
2827  a2(i) = mass - a1(i) - nemiss
2828  z2(i) = pnr - z1(i)
2829  if (fissrate(pos)%mode.eq.3) z2(i) = z2(i) + 1 ! beta-delayed fission
2830  psum = psum + paz(i)
2831 
2832  ! check if fragment 2 is in the network
2833  if (a2(i).gt.minmax(z2(i),2)) then
2834  nemiss = nemiss + a2(i) - minmax(z2(i),2)
2835  a2(i) = minmax(z2(i),2)
2836  else if (a2(i).lt.minmax(z2(i),1)) then
2837  print *, 'a2 error:',a2(i),minmax(z2(i),1)
2838  end if
2839  end do inner
2840 
2841  again: do i=1,nof
2842  paz(i) = paz(i)/psum ! normalize probabilities
2843 
2844  fissrate(pos)%channelprob(i) = paz(i)
2845 
2846  !---------rate for spontaneous and beta-delayed fission
2847  fissrate(pos)%fissparts(i+1) = findaz(a1(i),z1(i)) ! first fragment
2848  fissrate(pos)%ch_amount(i+1) = fissrate(pos)%channelprob(i) ! rate at which fragment is produced per destroyed parent nucleus
2849  fissrate(pos)%fissparts(i+1+fissrate(pos)%channels) = findaz(a2(i),z2(i)) ! second fragment
2850  fissrate(pos)%ch_amount(i+1+fissrate(pos)%channels) = fissrate(pos)%channelprob(i)
2851  fissrate(pos)%q_value(i) = isotope(fissrate(pos)%fissnuc_index)%mass_exc - &
2852  nemiss * isotope(ineu)%mass_exc - &
2853  isotope(fissrate(pos)%fissparts(i+1))%mass_exc - &
2854  isotope(fissrate(pos)%fissparts(i+1+fissrate(pos)%channels))%mass_exc
2855  end do again
2856 
2857 
2858  qval = 0.d0
2859  qvalue: do i=1,fissrate(pos)%channels
2860  qval = qval + fissrate(pos)%q_value(i)*fissrate(pos)%channelprob(i) ! weighted average Q-value
2861  end do qvalue
2862  fissrate(pos)%averageQ = qval
2863 
2864  deallocate(a1)
2865  deallocate(a2)
2866  deallocate(z1)
2867  deallocate(z2)
2868  deallocate(paz)
2869 
2870  rewind(betafission,iostat=read_stat)
2871  if (read_stat /= 0) then
2872  call raise_exception("Could not rewind fission file.",&
2873  "abla_betafiss",190007)
2874  endif
2875  return
2876 
2877  else
2878 
2879  do i=1,nof
2880  read(betafission,*,iostat=read_stat) dummy
2881  if (read_stat /= 0) then
2882  call raise_exception("Could not read fission file.",&
2883  "abla_betafiss",190007)
2884  endif
2885  end do
2886 
2887  end if
2888 
2889  end do
2890 
2891  !--------old Panov model for nuclei not found in ABLA fission file---------------
2892  allocate(a1(1),stat=alloc_stat)
2893  allocate(a2(1),stat=alloc_stat)
2894  allocate(z1(1),stat=alloc_stat)
2895  allocate(z2(1),stat=alloc_stat)
2896 
2897  fmass = mass
2898  if (fissrate(pos)%reac_type.eq. rrt_sf) then
2899  fpnr = pnr
2900  fissrate(pos)%mode = 2 ! spontaneous fission
2901  else if (fissrate(pos)%reac_type.eq. rrt_bf) then
2902  fpnr = pnr + 1
2903  fissrate(pos)%mode = 3 ! beta-delayed fission
2904  end if
2905 
2906  select case (fmass)
2907  case(:240) ! -> symmetric fission
2908  a2(1) = nint(fmass/2.d0)
2909  z2(1) = nint(fpnr/2.d0)
2910  a1(1) = fmass - a2(1)
2911  z1(1) = fpnr - z2(1)
2912  case(255:265) ! -> symmetric fission
2913  a2(1) = nint(fmass/2.d0)
2914  z2(1) = nint(fpnr/2.d0)
2915  a1(1) = fmass - a2(1)
2916  z1(1) = fpnr - z2(1)
2917  case default ! -> asymmetric fission
2918  a2(1) = 130
2919  z2(1) = nint(52.01 - (fpnr - 80)/1.d1)
2920  ! a2 = max(a2,af-a2) !modification to prevent very n-rich
2921  ! z2 = max(z2,zf-z2) !fragments
2922  a1(1) = fmass - a2(1)
2923  z1(1) = fpnr - z2(1)
2924  end select
2925 
2926  nemiss = 0
2927 
2928  if (a1(1).gt.minmax(z1(1),2)) then
2929  nemiss = a1(1) - minmax(z1(1),2)
2930  a1(1) = minmax(z1(1),2)
2931  else if (a1(1).lt.minmax(z1(1),1)) then
2932  print *, 'a1 error:'
2933  end if
2934 
2935  if (a2(1).gt.minmax(z2(1),2)) then
2936  nemiss = nemiss + a2(1) - minmax(z2(1),2)
2937  a2(1) = minmax(z2(1),2)
2938  else if (a2(1).lt.minmax(z2(1),1)) then
2939  print *, 'a2 error:',a2(1),minmax(z2(1),1)
2940  end if
2941  nufiss = nufiss + nemiss
2942 
2943  if (nemiss .eq. 0) then
2944  fissrate(pos)%dimens = 3
2945  else
2946  fissrate(pos)%dimens = 4
2947  end if
2948 
2949  allocate(fissrate(pos)%fissparts(fissrate(pos)%dimens))
2950  allocate(fissrate(pos)%ch_amount(fissrate(pos)%dimens))
2951  allocate(fissrate(pos)%channelprob(1))
2952  allocate(fissrate(pos)%q_value(1))
2953  fissrate(pos)%channels = 1
2954  fissrate(pos)%fissnuc_index = findaz(mass,pnr)
2955  fissrate(pos)%channelprob(1) = 1.d0
2956 
2957  fissrate(pos)%fissparts(1) = fissrate(pos)%fissnuc_index
2958  fissrate(pos)%ch_amount(1) = -1.d0
2959 
2960  fissrate(pos)%fissparts(2) = findaz(a1(1),z1(1))
2961  fissrate(pos)%ch_amount(2) = 1.d0
2962  fissrate(pos)%fissparts(3) = findaz(a2(1),z2(1))
2963  fissrate(pos)%ch_amount(3) = 1.d0
2964  if (nemiss.gt.0) then
2965  fissrate(pos)%fissparts(4) = ineu
2966  fissrate(pos)%ch_amount(4) = float(nemiss)
2967  end if
2968 
2969  fissrate(pos)%averageQ = isotope(fissrate(pos)%fissparts(1))%mass_exc - &
2970  isotope(fissrate(pos)%fissparts(2))%mass_exc - &
2971  isotope(fissrate(pos)%fissparts(3))%mass_exc - &
2972  nemiss * isotope(ineu)%mass_exc
2973 
2974  deallocate(a1)
2975  deallocate(a2)
2976  deallocate(z1)
2977  deallocate(z2)
2978 
2979 
2980  rewind(betafission,iostat=read_stat)
2981  if (read_stat /= 0) then
2982  call raise_exception("Could not rewind fission file.",&
2983  "abla_betafiss",190007)
2984  endif
2985 
2986  info_exit("abla_betafiss")
2987 
2988  end subroutine abla_betafiss
2989 
2990 
2991 
2992 end module fission_rate_module
parameter_class::sfission_file
character(max_fname_len) sfission_file
Fission table for spontaneous fission.
Definition: parameter_class.f90:180
parameter_class::bfission_file
character(max_fname_len) bfission_file
Fission table for beta-delayed fission.
Definition: parameter_class.f90:178
fission_rate_module::count_fission_rates
subroutine, private count_fission_rates()
Count the amount of fission rates.
Definition: fission_rate_module.f90:693
fission_rate_module::n_bf
integer, private n_bf
Definition: fission_rate_module.f90:74
parameter_class::use_prepared_network
logical use_prepared_network
Use a prepared folder with all necessary data in binary format.
Definition: parameter_class.f90:94
global_class::net_names
character *5, dimension(:), allocatable, public net_names
list of isotopes contained in the network
Definition: global_class.f90:95
rrt_nf
#define rrt_nf
Definition: macros.h:78
global_class::ipro
integer, public ipro
index of alphas, neutrons and protons
Definition: global_class.f90:94
global_class::reactionrate_type
reaction rate type
Definition: global_class.f90:44
fission_rate_module::beta_delayed_fiss_probs
real(r_kind), dimension(:,:), allocatable, private beta_delayed_fiss_probs
Definition: fission_rate_module.f90:41
fission_rate_module::fragtype
Definition: fission_rate_module.f90:48
fission_rate_module::nfiss_spont
integer, private nfiss_spont
Amount of spontaneous fission rates.
Definition: fission_rate_module.f90:69
benam_class::findaz
integer function, public findaz(ai, zi)
Finds an isotope index in the network table, given its A and Z. In case it was not found,...
Definition: benam_class.f90:887
fission_rate_module::init_fission_rates
subroutine, public init_fission_rates()
Initialize the fission reactions.
Definition: fission_rate_module.f90:99
fission_rate_module::n_sf
integer, private n_sf
Amount of individual reactions.
Definition: fission_rate_module.f90:74
parameter_class::fission_format_spontaneous
integer fission_format_spontaneous
Format of spontaneous fission rates (0: Off, 1: Reaclib, 2: Halflifes)
Definition: parameter_class.f90:164
fission_rate_module::read_fiss_fragment_file
subroutine read_fiss_fragment_file(file, fragment_array)
Definition: fission_rate_module.f90:2120
fission_rate_module::read_fission_rates_probability_format
subroutine, private read_fission_rates_probability_format(fission_path, reac_type, start_idx)
Definition: fission_rate_module.f90:1434
fission_rate_module::count_fission_rates_reaclib_format
subroutine, private count_fission_rates_reaclib_format(fission_rate_file, count_rates)
Count the amount of fission rates in reaclib format.
Definition: fission_rate_module.f90:771
rrs_reacl
#define rrs_reacl
Definition: macros.h:49
error_msg_class
Error handling routines.
Definition: error_msg_class.f90:16
benam_class::benam
integer function, public benam(name)
Returns the index number of isotope 'name'.
Definition: benam_class.f90:251
fission_rate_module
Module to deal with fission reactions.
Definition: fission_rate_module.f90:16
fission_rate_module::abla_betafiss
subroutine, private abla_betafiss(pos, mass, pnr, betafission, qval)
Calculates the (beta-delayed and spontaneous) fission fragment mass distribution according to the ABL...
Definition: fission_rate_module.f90:2743
error_msg_class::int_to_str
character(:) function, allocatable, public int_to_str(num)
Converts a given integer to a string.
Definition: error_msg_class.f90:224
rrs_aext
#define rrs_aext
Definition: macros.h:56
fission_rate_module::fissfrags_beta_delayed
type(fragtype), dimension(:), allocatable, private fissfrags_beta_delayed
Array storing fragment distributions of beta-delayed fission from file.
Definition: fission_rate_module.f90:61
error_msg_class::raise_exception
subroutine, public raise_exception(msg, sub, error_code)
Raise a exception with a given error message.
Definition: error_msg_class.f90:245
fission_rate_module::count_fission_rates_halflife_format
subroutine, private count_fission_rates_halflife_format(fission_rate_file, count_rates)
Definition: fission_rate_module.f90:840
fission_rate_module::nfiss
integer, public nfiss
Amount of fission rates.
Definition: fission_rate_module.f90:68
fission_rate_module::read_fission_rates_halflife_format
subroutine, private read_fission_rates_halflife_format(fission_path, reac_type, start_idx)
Definition: fission_rate_module.f90:1547
fission_rate_module::nfiss_bdel
integer, private nfiss_bdel
Amount of beta-delayed fission rates.
Definition: fission_rate_module.f90:71
parameter_class::max_fname_len
integer, parameter, public max_fname_len
maximum length of filenames
Definition: parameter_class.f90:58
global_class::ineu
integer, public ineu
Definition: global_class.f90:94
mergesort_module::reorder_int
subroutine, public reorder_int(arr1, arr2, length)
Takes an array and an array with indices to reorder the first array.
Definition: mergesort_module.f90:932
global_class::isotope
type(isotope_type), dimension(:), allocatable, public isotope
all nuclides used in the network
Definition: global_class.f90:34
error_msg_class::write_data_to_std_out
subroutine, public write_data_to_std_out(str_msg, value_str, unit)
Write data to the standard output (usually OUT)
Definition: error_msg_class.f90:122
format_class::my_format
character(60), dimension(100) my_format
Definition: format_class.f90:19
fission_rate_module::nufiss
integer, private nufiss
total number of fission neutrons
Definition: fission_rate_module.f90:72
tabulated_rate_module::multiply_tab_rate_by_factor
subroutine, public multiply_tab_rate_by_factor(rrate, factor)
Multiply a tabulated rate by a factor.
Definition: tabulated_rate_module.f90:212
parameter_class::fission_rates_spontaneous
character(max_fname_len) fission_rates_spontaneous
reaction library for spontaneous fission rates
Definition: parameter_class.f90:159
nucstuff_class
Module with helper functions such as rate tabulation, partition functions, and theoretical weak rates...
Definition: nucstuff_class.f90:11
benam_class::minmax
integer, dimension(:,:), allocatable minmax
TODO: add description.
Definition: benam_class.f90:14
fission_rate_module::fiss_dist
subroutine, private fiss_dist(fissrate_inout)
Determines fission fragment mass distribution as described in Panov et al. 2001.
Definition: fission_rate_module.f90:1645
rrs_wext
#define rrs_wext
Definition: macros.h:55
parameter_class::fissflag
integer fissflag
defines type of fission fragment distribution used
Definition: parameter_class.f90:86
fission_rate_module::read_binary_fission_reaction_data
subroutine, private read_binary_fission_reaction_data(path)
Read the fission reactions and fragment distributions from a binary file.
Definition: fission_rate_module.f90:220
rrs_twr
#define rrs_twr
Definition: macros.h:52
fission_rate_module::fissfrags_spontaneous
type(fragtype), dimension(:), allocatable, private fissfrags_spontaneous
Array storing fragment distributions of spontaneous fission from file.
Definition: fission_rate_module.f90:60
fission_rate_module::write_reac_verbose_out
subroutine, private write_reac_verbose_out()
Write the amount of individual reactions to the out.
Definition: fission_rate_module.f90:184
fission_rate_module::nfiss_n_ind
integer, private nfiss_n_ind
Amount of neutron induced fission rates.
Definition: fission_rate_module.f90:70
global_class::net_size
integer, public net_size
total number of isotopes (network size)
Definition: global_class.f90:93
parameter_class::unit
type(unit_type), public unit
constants and unit conversion factors
Definition: parameter_class.f90:54
rrs_nu
#define rrs_nu
Definition: macros.h:51
fission_rate_module::fissrate
type(fissionrate_type), dimension(:), allocatable, public fissrate
Array storing fission reactions.
Definition: fission_rate_module.f90:39
fission_rate_module::abla_nfiss
subroutine, private abla_nfiss(pos, mass, pnr, neutfission, qval)
Calculates the (neutron-induced) fission fragment mass distribution according to the ABLA07 model: Ke...
Definition: fission_rate_module.f90:2509
fission_rate_module::add_fission_fragments
subroutine add_fission_fragments()
Initializes the rates with fragments.
Definition: fission_rate_module.f90:1007
file_handling_class
Provide some basic file-handling routines.
Definition: file_handling_class.f90:18
parameter_class::fission_frag_n_induced
integer fission_frag_n_induced
Fragment distribution of n-induced fission rates in case of custom fragments (fissflag=4)
Definition: parameter_class.f90:166
fission_rate_module::kodtakdist
subroutine, private kodtakdist(fissrate_inout)
Kodama-Takahashi distribution.
Definition: fission_rate_module.f90:1827
parameter_class::fission_frag_spontaneous
integer fission_frag_spontaneous
Fragment distribution of spontaneous fission rates in case of custom fragments (fissflag=4)
Definition: parameter_class.f90:165
parameter_class::nfission_file
character(max_fname_len) nfission_file
Fission table for neutron-induced fission.
Definition: parameter_class.f90:179
fission_rate_module::modify_halflifes
subroutine, private modify_halflifes(rrate, rrate_length, fissrate_in, nfiss_in)
Modifies half lifes of beta decays to include fission.
Definition: fission_rate_module.f90:445
fission_rate_module::output_binary_fission_reaction_data
subroutine, public output_binary_fission_reaction_data(path)
Save the fission data to a unformatted binary file.
Definition: fission_rate_module.f90:316
r_kind
#define r_kind
Definition: macros.h:46
tabulated_rate_module::tabulated_index
subroutine, public tabulated_index(temp)
Set tab_index for a given temperature.
Definition: tabulated_rate_module.f90:537
fission_rate_module::fiss_binary_name
character(len= *), parameter, private fiss_binary_name
Name of the binary file containing the fission rates.
Definition: fission_rate_module.f90:45
parameter_class::fission_format_beta_delayed
integer fission_format_beta_delayed
Format of beta-delayed fission rates (0: Off, 1: Reaclib, 2: Halflifes, 3: Probability)
Definition: parameter_class.f90:162
fission_rate_module::fissfrags_n_induced
type(fragtype), dimension(:), allocatable, private fissfrags_n_induced
Array storing fragment distributions of neutron-induced fission from file.
Definition: fission_rate_module.f90:59
fission_rate_module::fissionrate_type
fission rate type, designed to save fragment distribution at initialization
Definition: fission_rate_module.f90:21
fission_rate_module::fiss_neglect
integer, parameter, public fiss_neglect
Amount of fission fragments not to be neglected.
Definition: fission_rate_module.f90:66
fission_rate_module::n_nf
integer, private n_nf
Definition: fission_rate_module.f90:74
rrt_betm
#define rrt_betm
Definition: macros.h:59
nucstuff_class::get_nr_reactants
integer function get_nr_reactants(group)
Get number of reactants of the reaction based on the reaclib chapter.
Definition: nucstuff_class.f90:301
parameter_class::fission_rates_n_induced
character(max_fname_len) fission_rates_n_induced
reaction library for neutron induced fission rates
Definition: parameter_class.f90:160
rrs_fiss
#define rrs_fiss
Definition: macros.h:54
file_handling_class::open_infile
integer function, public open_infile(file_name)
Same for reading (input file)
Definition: file_handling_class.f90:126
tabulated_rate_module
This module contains everything for the tabulated rates that can replace reaclib rates.
Definition: tabulated_rate_module.f90:28
fission_rate_module::reorder_fragments
subroutine, private reorder_fragments()
Reorder the fission fragments.
Definition: fission_rate_module.f90:148
global_class
Contains types and objects shared between multiple modules.
Definition: global_class.f90:19
fission_rate_module::file_fiss_frag
subroutine file_fiss_frag(fissrate_inout, missing_frags)
Fill the rates with the correct fragments.
Definition: fission_rate_module.f90:2334
parameter_class::fission_format_n_induced
integer fission_format_n_induced
Format of neutron-induced fission rates (0: Off, 1: Reaclib)
Definition: parameter_class.f90:163
tabulated_rate_module::calculate_tab_rate
subroutine, public calculate_tab_rate(rrate, temp, rat_calc)
Calculates the tabulated rate.
Definition: tabulated_rate_module.f90:191
fission_rate_module::merge_fission_rates
subroutine, public merge_fission_rates(rrate_array, rrate_length, fiss_count)
Merge fission rates with larger array.
Definition: fission_rate_module.f90:388
file_handling_class::close_io_file
subroutine, public close_io_file(unit_no, file_name)
Close an external file.
Definition: file_handling_class.f90:144
rrt_bf
#define rrt_bf
Definition: macros.h:79
file_handling_class::open_unformatted_outfile
integer function, public open_unformatted_outfile(file_name)
Shorthand for opening a new unformatted file for writing (output file)
Definition: file_handling_class.f90:74
rrt_sf
#define rrt_sf
Definition: macros.h:80
format_class
Define custom format statements used to read in major data files.
Definition: format_class.f90:16
fission_rate_module::count_fission_rates_probability_format
subroutine, private count_fission_rates_probability_format(fission_rate_file, count_rates, nr_cols)
Definition: fission_rate_module.f90:923
file_handling_class::open_unformatted_infile
integer function open_unformatted_infile(file_name)
Open an unformatted file for reading.
Definition: file_handling_class.f90:89
rrs_detb
#define rrs_detb
Definition: macros.h:53
parameter_class::fission_frag_missing
integer fission_frag_missing
Fragment distribution in case of missing fragments in case of custom fragments (fissflag=4)
Definition: parameter_class.f90:168
mergesort_module::bubblesort
subroutine, public bubblesort(mode, length, arr1, arr2)
Bubblesort of array.
Definition: mergesort_module.f90:724
mergesort_module
Module mergesort_module for merging arrays of rates.
Definition: mergesort_module.f90:16
fission_rate_module::read_fission_rates
subroutine, private read_fission_rates()
Read the fission rates, splitted into different types of fission.
Definition: fission_rate_module.f90:1208
parameter_class::prepared_network_path
character(max_fname_len) prepared_network_path
Prepared network folder.
Definition: parameter_class.f90:198
rrs_tabl
#define rrs_tabl
Definition: macros.h:50
fission_rate_module::read_fission_rates_reaclib_format
subroutine, private read_fission_rates_reaclib_format(fission_path, reac_type, start_idx)
Definition: fission_rate_module.f90:1300
parameter_class::fission_rates_beta_delayed
character(max_fname_len) fission_rates_beta_delayed
reaction library for beta delayed fission rates
Definition: parameter_class.f90:161
fission_rate_module::amount_cols
integer, private amount_cols
Amount of columns in the beta-delayed fission file.
Definition: fission_rate_module.f90:42
parameter_class
Contains all runtime parameters as well as phys and math constants.
Definition: parameter_class.f90:24
parameter_class::fission_frag_beta_delayed
integer fission_frag_beta_delayed
Fragment distribution of beta-delayed fission rates in case of custom fragments (fissflag=4)
Definition: parameter_class.f90:167
benam_class
Subroutines needed to initialise the network.
Definition: benam_class.f90:11