analysis.f90
Go to the documentation of this file.
1 !> @file analysis.f90
2 !!
3 !! The error file code for this file is ***W10***.
4 !! @brief Module \ref analysis with various analysis subroutines
5 !!
6 
7 !> Contains various routines for analysis and diagnostic output
8 #include "macros.h"
9 module analysis
10 
12  use global_class, only: net_size, ihe4, ineu, ipro
14 
15  implicit none
16 
17  real(r_kind),dimension(:),allocatable :: sn !< array of n-separation energies
18  real(r_kind) :: sn_ave !< average n-separation energy
19  integer :: mainout_unit !< file descriptor for the main analysis output file
20  integer :: nrdiag_unit !< diagnostic output inside the Newton-Rhapson loop
21  integer :: track_unit !< file id for tracked nuclei
22  integer :: toplist_unit !< file id for top nuclear energy contributors
23  integer :: nu_loss_gain_unit !< file id for neutrino loss/gain
24  integer :: cl_start, cl_end, cl_rate, cl_cmax !< system clock variables
25  integer,private :: tsfile !< timescale file unit
26  integer,private :: nuc_heat_file !< nuclear heating file unit
27  real(r_kind) :: tsnp, t9snp, rosnp !< for snapshot printing triggers
28  integer :: snapcnt !< snapshot count
29  integer :: flowcnt !< flow printing counter
30  real,dimension(:),allocatable :: snapshot_time !< point in time for custom snapshot
31  integer :: snapshot_amount !< Amount of custom snapshots
32  integer,parameter :: toplist_amount=10 !< Amount of reactions in the toplist
33 
34  !
35  ! Public and private fields and methods of the module
36  !
37  public:: &
40 
41  private:: &
45 
46 contains
47 
48 !>
49 !! Open files, write headers, allocate space etc.
50 !!
51 subroutine analysis_init()
52  use parameter_class
54  use error_msg_class, only: int_to_str
55  use benam_class, only: findaz
57 #ifdef USE_HDF5
58  use hdf5_module, only: init_hdf5_module
59 #endif
60  implicit none
61 
62  ! neutron separation energies
63  real(r_kind) :: val !< aux variable for reading Sn file
64  integer :: ai, ni, zi, snfile, fstat !< temporary variables needed for reading the Sn file
65  character :: nl, buf(30) !< newline character, buffer
66  integer :: i,k
67  character*10,dimension(track_nuclei_nr) :: tmp_nucleinames !< Name of the nuclei that will be tracked
68  character*1500 :: head_toplist
69  character*5 :: help_char
70  info_entry("analysis_init")
71 
72 #ifdef USE_HDF5
73  ! Only initialize if it is really possible to execute it
74  ! (i.e., correct compiler h5fc)
75  call init_hdf5_module()
76 #endif
77 
78  cl_cmax= 100000
79  snapcnt= 1
80  flowcnt= 1
81  nl = new_line('A')
82 
83  ! main output file: create, write header
84  if (mainout_every .gt. 0) then
85  mainout_unit= open_outfile("mainout.dat")
86  write(mainout_unit,'(3A)') &
87  '# 1:iteration 2:time[s], 3:T[GK], 4:rho[g/cm3], 5:Ye '//nl &
88  , '# 6:R[km], 7:Y_n, 8:Y_p, 9:Y_alpha, 10:Y_lights '//nl &
89  , '# 11:Y_heavies 12:<Z> 13:<A> 14:entropy [kB/baryon] (15:Sn [MeV])'
90  endif
91 
92  ! read neutron separation energies
93  if (calc_nsep_energy) then
95 
96  allocate(sn(net_size))
97  sn(1:net_size) = 0d0
98  read(snfile, *) buf ! skip the first line
99  do
100  read(snfile, *, iostat=fstat) zi, ni, ai, val
101  if (fstat /= 0) exit
102  k = findaz(ai, zi)
103  if (k>0) then
104  sn(k) = val
105  endif
106  enddo
107  call close_io_file(snfile,nsep_energies_file)
108  endif
109 
110  ! Open the timescale file and write a header
111  if (timescales_every .gt. 0) then
112  tsfile = open_outfile('timescales.dat')
113  write(tsfile, '((A1,(A22,3x),*(A23,3x)))') &
114  "#","time[s]","Temperature [GK]","tau_ga [s]","tau_ag [s]","tau_ng [s]","tau_gn [s]", &
115  "tau_pg [s]", "tau_gp [s]", "tau_np [s]", "tau_pn [s]", &
116  "tau_an [s]", "tau_na [s]", "tau_ap [s]", "tau_pa [s]", "tau_beta [s]", "tau_alpha [s]", &
117  "tau_nfiss [s]", "tau_sfiss [s]", "tau_bfiss [s]"
118  end if
119 
120  ! Open the track_nuclei file and write a header
121  if (track_nuclei_every .gt. 0) then
122  ! Convert indices to names
123  do i=1, track_nuclei_nr
124  tmp_nucleinames(i) = "Y( "//net_names(track_nuclei_indices(i))//" )"
125  end do
126  ! Write the header
127  track_unit = open_outfile('tracked_nuclei.dat')
128  write(track_unit, '((A,20X,*(A,18X)))') "# time[s]",tmp_nucleinames(:)
129  end if
130 
131  ! Write the header for the toplist
132  if (top_engen_every .gt. 0) then
133  toplist_unit= open_outfile("toplist.dat")
134 
135  ! Give some hint on what is what
136  write (toplist_unit,"('#',31X,A,434X,A,413X,'|',5X,A,185X,A)") "Top nuclear energy producers for each reaction given in erg/g/s", &
137  "Top nuclear energy consumers for each reaction given in erg/g/s", &
138  "Top nuclear energy producers for each isotope given in erg/g/s", &
139  "Top nuclear energy consumers for each isotope given in erg/g/s"
140 
141  ! Construct the rather complicated header of the file
142  head_toplist ="# Time[s] Etot[erg/g/s] Reaction_1"
143  do i=2,toplist_amount
144  head_toplist = trim(adjustl(head_toplist))//" Reaction_"//trim(adjustl(int_to_str(i)))
145  end do
146  head_toplist = trim(adjustl(head_toplist))//" Ereac_1"
147  do i=2,toplist_amount
148  head_toplist = trim(adjustl(head_toplist))//" Ereac_"//trim(adjustl(int_to_str(i)))
149  end do
150  head_toplist = trim(adjustl(head_toplist))//" Reaction_1"
151  do i=2,toplist_amount
152  head_toplist = trim(adjustl(head_toplist))//" Reaction_"//trim(adjustl(int_to_str(i)))
153  end do
154  head_toplist = trim(adjustl(head_toplist))//" Ereac_1"
155  do i=2,toplist_amount
156  head_toplist = trim(adjustl(head_toplist))//" Ereac_"//trim(adjustl(int_to_str(i)))
157  end do
158  head_toplist = trim(adjustl(head_toplist))//" | Iso1"
159  do i=2,toplist_amount
160  help_char = "Iso"//trim(adjustl(int_to_str(i)))
161  head_toplist = trim(adjustl(head_toplist))//" "//adjustr(help_char)
162  end do
163  head_toplist = trim(adjustl(head_toplist))//" Eiso_1"
164  do i=2,toplist_amount
165  head_toplist = trim(adjustl(head_toplist))//" Eiso_"//trim(adjustl(int_to_str(i)))
166  end do
167  head_toplist = trim(adjustl(head_toplist))//" Iso1"
168  do i=2,toplist_amount
169  help_char = "Iso"//trim(adjustl(int_to_str(i)))
170  head_toplist = trim(adjustl(head_toplist))//" "//adjustr(help_char)
171  end do
172  head_toplist = trim(adjustl(head_toplist))//" Eiso_1"
173  do i=2,toplist_amount
174  head_toplist = trim(adjustl(head_toplist))//" Eiso_"//trim(adjustl(int_to_str(i)))
175  end do
176  ! Finally write it out
177  write (toplist_unit,"(A)") trim(adjustl(head_toplist))
178  end if
179 
180 
181  if (nrdiag_every.gt.0) then
182  nrdiag_unit= open_outfile("nrdiag.dat")
183  write(nrdiag_unit,'(11A)') &
184  '# Newton-Raphson diagnostic output'//nl &
185  , '# 1 : global iteration count (cnt)'//nl &
186  , '# 2 : global time [s]'//nl &
187  , '# 3 : temperature T9 [GK]'//nl &
188  , '# 4 : adapt stepsize loop counter (k)'//nl &
189  , '# 5 : NR loop counter (nr_count)'//nl &
190  , '# 6 : global timestep [s]'//nl &
191  , '# 7 : mass conservation error'//nl &
192  , '# 8 : maximal abundance change (eps)'//nl &
193  , '# 9 : most rapidly evolving isotope (epsl)'//nl &
194  , '# 10: abundance of the isotope epsl, y(epsl)'
195  end if
196 
197  ! Open the energy generation file and write a header
198  if (engen_every .gt. 0) then
199  nuc_heat_file = open_outfile('generated_energy.dat')
200  ! Header
201  write(nuc_heat_file, '((A,5X))') "# Generated Energy for each nuclear reaction type given in erg/g/s"
202  write(nuc_heat_file, '(17(A23,3X))') &
203  "# time[s]", "Engen(Total)", "S_src","Engen(weak)","Engen(n,g)",&
204  "Engen(p,g)","Engen(a,g)","Engen(n,p)","Engen(n,a)", &
205  "Engen(p,a)", "Engen(fiss)"
206  end if
207 
208 
209 
210  if ((nu_loss_every .gt. 0)) then
211  nu_loss_gain_unit = open_outfile('nu_loss_gain.dat')
212  write(nu_loss_gain_unit, '(A)') "# File containing the neutrino loss/gain for the nuclear heating."
213  write(nu_loss_gain_unit, '(A)') "# Negative values mean that energy is added to the system."
214  write(nu_loss_gain_unit, '(A)') "# Neutrino energies are given in Erg/g/s."
215  write(nu_loss_gain_unit, '((A1,(A18,3x),*(A19,3x)))') &
216  "#","time[s]","T[GK]", "rho[g/cm3]", "R[km]", "nu_total", "nu_beta",&
217  "nu_nuheat", "nu_thermal"
218  end if
219 
220 
221  info_exit("analysis_init")
222 
223 end subroutine analysis_init
224 
225 
226 !>
227 !! Output first step
228 !!
229 !! This routine calls output_iteration at the first step
230 !!
231 !! @see output_final_step
232 subroutine output_initial_step(ctime,T9,rho_b,entropy,Rkm,Y,pf)
233  implicit none
234  real(r_kind),intent(in) :: ctime !< actual time [s]
235  real(r_kind),intent(in) :: t9 !< initial temperature [GK]
236  real(r_kind),intent(in) :: rho_b !< initial density [gcc]
237  real(r_kind),intent(in) :: entropy !< entropy [kB/nucleon]
238  real(r_kind),intent(in) :: rkm !< Radius [km]
239  real(r_kind),dimension(net_size),intent(in) :: y !< array of abundances (Y_i)
240  real(r_kind),dimension(0:net_size),intent(inout) :: pf !< partition functions
241 
242  info_entry('output_initial_step')
243  call output_iteration(0,ctime,t9,rho_b,entropy,rkm,y,pf)
244  info_exit('output_initial_step')
245 
246 end subroutine output_initial_step
247 
248 
249 !>
250 !! Output final step
251 !!
252 !! This routine calls output_iteration at the last step
253 !!
254 !! @see output_final_step
255 !!
256 !! \b Edited: OK 26.08.2017
257 subroutine output_final_step(cnt,ctime,T9,rho_b,entropy,Rkm,Y,pf)
260  implicit none
261  integer,intent(in) :: cnt !< current iteration counter
262  real(r_kind),intent(in) :: ctime !< actual time [s]
263  real(r_kind),intent(in) :: t9 !< initial temperature [GK]
264  real(r_kind),intent(in) :: rho_b !< initial density [gcc]
265  real(r_kind),intent(in) :: entropy !< entropy [kB/nucleon]
266  real(r_kind),intent(in) :: rkm !< Radius [km]
267  real(r_kind),dimension(net_size),intent(in) :: y !< array of abundances (Y_i)
268  real(r_kind),dimension(0:net_size),intent(inout) :: pf !< partition functions
269 
270  info_entry('output_final_step')
271 
272  ! hack output frequency parameters to output the final step out-of-order
273  if (out_every.gt.0) out_every= 1
274  if (mainout_every.gt.0) mainout_every= 1
275  if (snapshot_every.gt.0) snapshot_every= 1
276  if (flow_every.gt.0) flow_every= 1
277  if (timescales_every.gt.0) timescales_every= 1
279  print '("---------------")'
280  call output_iteration(cnt,ctime,t9,rho_b,entropy,rkm,y,pf)
281  print '("============================")'
282  select case(termination_criterion)
283  case(0)
284  print'(A)',"End of trajectory file reached."
285  case(1)
286  print'(A)',"Final time reached."
287  case(2)
288  print'(A)',"Final temperature reached."
289  case(3)
290  print'(A)',"Final density reached."
291  case(4)
292  print'(A)',"Freeze-out reached."
293  endselect
294 
295  info_exit('output_final_step')
296 
297 end subroutine output_final_step
298 
299 
300 !>
301 !! Periodic output of analysis data for current iteration
302 !!
303 !! Prints various outputs depending on the input parameter.
304 !!
305 !! \b Edited:
306 !! - MR : 19.01.21
307 !! .
308 subroutine output_iteration(cnt,ctime,T9,rho_b,entropy,Rkm,Y,pf)
316  use nucstuff_class, only: inter_partf
321  use global_class, only: isotope
322  use flow_module, only: flowcalc, flows
323 #ifdef USE_HDF5
324  use hdf5_module
325 #endif
326  implicit none
327 
328  integer,intent(in) :: cnt !< current iteration counter
329  real(r_kind),intent(in) :: ctime !< actual time [s]
330  real(r_kind),intent(in) :: t9 !< temperature [GK]
331  real(r_kind),intent(in) :: rho_b !< density [gcc]
332  real(r_kind),intent(in) :: entropy !< entropy [kB/nucleon]
333  real(r_kind),intent(in) :: rkm !< length scale [km]
334  real(r_kind),dimension(net_size),intent(in) :: y !< array of abundances (Y_i)
335  real(r_kind),dimension(0:net_size),intent(inout) :: pf !< partition functions
336  integer :: i !< Loop variable
337  real(r_kind) :: ysum,yasum,abar,ye
338 
339  info_entry("output_iteration")
340 
341  ! screen output
342  if (out_every.gt.0) then
343  if (mod(cnt,50*out_every)==0) then
344  print *
345  print '(A)',"#-- 1:it 2:time[s] 3:T9[GK] 4:rho[g/cm3] &
346  5:entropy[kB/baryon]"
347  endif
348  if (mod(cnt,out_every) .eq. 0) then
349  print '(I8,X,ES14.7,X,F10.7,X,ES14.7,X,F12.7)', &
350  cnt,ctime,t9,rho_b,entropy
351  endif
352 
353  if(mod(cnt,out_every).eq.0 .and. (heating_mode .eq. 1)) then
354  if (verbose_level.ge.2) then
355  call output_nuclear_heating(cnt,ctime)
356  end if
357  endif
358  endif
359 
360  ! This is only done if the correct compiler (h5fc) was used
361 #ifdef USE_HDF5
362  ! Store snapshots in hdf5
363  if (h_snapshot_every.gt.0) then
364  if (mod(cnt,h_snapshot_every).eq. 0) then
365  ! Ye = sum(isotope(1:net_size)%p_nr*Y(1:net_size))
366  call extend_snaps(ctime,y)
367  end if
368  end if
369 
370  ! Store mainout in hdf5
371  if (h_mainout_every.gt.0) then
372  if (mod(cnt,h_mainout_every).eq. 0) then
373  call extend_mainout(cnt,ctime,t9,rho_b,entropy,rkm,y)
374  end if
375  end if
376 
377  ! Snapshot output for specific times
378  if (h_custom_snapshots) then
379  ! Test if a snapshot should be done
380  h_make_snap: do i=1,snapshot_amount
381  ! Time is within certain tolerance? Do the snapshot!
382  if (abs(ctime-snapshot_time(i)) .le. 1e-20) then
383  call extend_snaps(ctime,y)
384  exit h_make_snap
385  end if
386  enddo h_make_snap
387  end if
388 
389  ! Output the tracked nuclei
390  if (h_track_nuclei_every .gt. 0) then
391  if (mod(cnt,h_track_nuclei_every).eq.0) then
392  call extend_track_nuclei(ctime,y)
393  end if
394  end if
395 
396  ! Output average timescales
397  if (h_timescales_every .gt. 0) then
398  if ((mod(cnt,h_timescales_every).eq.0) .and. (evolution_mode .ne. em_nse)) then
399  ! calculate averaged timescales of various types of reactions
400  ! calculate partition function first
401  call inter_partf (t9, pf)
402  call calc_av_timescales(ctime,t9,y,.true.)
403  end if
404  end if
405 
406  ! Output the flow
407  if ((h_flow_every .gt. 0) .and. (cnt .ne. 0)) then
408  if ((mod(cnt,h_flow_every).eq.0)) then
409  call flowcalc(y)
410  call extend_flow(ctime,t9,rho_b,flows,size(flows),y)
411  end if
412  end if
413 
414  ! Calculate and output generated energy
415  if (h_engen_every .gt. 0) then
416  if (mod(cnt,h_engen_every).eq.0) then
417  call calc_nuclear_heating(.true., .true., .false.)
418  end if
419  end if
420 
421  ! Calculate and output generated energy
422  if ((h_nu_loss_every .gt. 0)) then
423  if (mod(cnt,h_nu_loss_every).eq.0) then
424 
425  if ((heating_mode .eq. 0) .and. (use_thermal_nu_loss)) then
426  ysum= sum(y(1:net_size))
427  yasum= sum(y(1:net_size)*isotope(1:net_size)%mass)
428  ye = sum(isotope(1:net_size)%p_nr*y(1:net_size))
429  abar= yasum / ysum
430  call thermal_neutrinos(abar, ye, t9, rho_b, stepsize, qdot_th)
431  qdot_th = qdot_th * unit%hix ! erg/g -> MeV/baryon
432  else
433  qdot_th = 0d0
434  end if
435 
436  call extend_nu_loss_gain(ctime,t9,rho_b,rkm,qdot_th/unit%hix/stepsize,&
437  qdot_nu/unit%hix/stepsize,&
438  qdot_bet/unit%hix/stepsize)
439  end if
440  end if
441 
442 #endif
443 
444  ! main data log
445  if (mainout_every .gt.0) then
446  if (mod(cnt,mainout_every).eq.0) then
447  call output_mainout(cnt,ctime,t9,rho_b,entropy,rkm,y)
448  end if
449  end if
450 
451  ! Iterative Snapshot and flow output
452  if(snapshot_every.gt.0) then
453  if(mod(cnt,snapshot_every).eq.0) then
454  call output_snapshot(ctime,t9,rho_b,y)
455  end if
456  end if
457 
458  if((flow_every.gt.0) .and. (cnt .ne. 0)) then
459  if(mod(cnt,flow_every).eq.0) then
460  call output_flow(ctime,t9,rho_b,y)
461  end if
462  end if
463 
464  if ((nu_loss_every .gt. 0)) then
465  if (mod(cnt,nu_loss_every).eq.0) then
466 
467  if ((heating_mode .eq. 0) .and. (use_thermal_nu_loss)) then
468  ysum= sum(y(1:net_size))
469  yasum= sum(y(1:net_size)*isotope(1:net_size)%mass)
470  ye = sum(isotope(1:net_size)%p_nr*y(1:net_size))
471  abar= yasum / ysum
472  call thermal_neutrinos(abar, ye, t9, rho_b, stepsize, qdot_th)
473  qdot_th = qdot_th * unit%hix ! erg/g -> MeV/baryon
474  else
475  qdot_th = 0d0
476  end if
477 
478  call output_nu_loss(ctime,t9,rho_b,rkm)
479  end if
480  end if
481 
482 
483  ! Snapshot output for specific times
484  if (custom_snapshots) then
485  ! Test if a snapshot should be done
486  make_snap: do i=1,snapshot_amount
487  ! Time is within certain tolerance? Do the snapshot!
488  if (abs(ctime-snapshot_time(i)) .le. 1e-20) then
489  call output_snapshot(ctime,t9,rho_b,y)
490  exit make_snap
491  end if
492  enddo make_snap
493  end if
494 
495  if (timescales_every .gt. 0) then
496  if ((mod(cnt,timescales_every).eq.0) .and. (evolution_mode .ne. em_nse)) then
497  ! calculate averaged timescales of various types of reactions
498  ! calculate partition function first
499  call inter_partf (t9, pf)
500  call calc_av_timescales(ctime,t9,y)
501  end if
502  end if
503 
504  ! Output the tracked nuclei
505  if (track_nuclei_every .gt. 0) then
506  if (mod(cnt,track_nuclei_every).eq.0) then
507  call output_track_nuclei(ctime,y)
508  end if
509  end if
510 
511  ! Calculate generated energy
512  if (engen_every .gt. 0) then
513  if (mod(cnt,engen_every).eq.0) then
514  call calc_nuclear_heating(.false., .true., .false.)
515  end if
516  end if
517 
518  ! Calculate generated energy
519  if (top_engen_every .gt. 0) then
520  if (mod(cnt,top_engen_every).eq.0) then
521  call calc_nuclear_heating(.false., .false., .true.)
522  end if
523  end if
524 
525 
526  info_exit("output_iteration")
527 end subroutine output_iteration
528 
529 
530 !> Output the abundances of specific nuclei
531 !!
532 !! First column is always the time, followed by the abundances of
533 !! nuclei given in "track_nuclei_file". An example of the file looks like:
534 !! \l_file{
535 !! # time[s] Y( n ) Y( o24 ) Y( f24 )
536 !! 0.0000000000000000E+00 0.0000000000000000E+00 2.0833333333333332E-02 0.0000000000000000E+00
537 !! 1.9999999999999999E-20 8.5150206007805894E-19 2.0833333333333332E-02 1.8661570495808601E-21
538 !! 5.9999999999999994E-20 2.5545061802341771E-18 2.0833333333333332E-02 5.5984711487425806E-21
539 !! 1.3999999999999998E-19 5.9605144205464133E-18 2.0833333333333332E-02 1.3063099347066021E-20
540 !! ...}
541 !!
542 !! @author: Moritz Reichert
543 !!
544 !! \b Edited:
545 !! - 09.06.17
546 !! .
547 subroutine output_track_nuclei(ctime,Y)
549  implicit none
550  real(r_kind),intent(in) :: ctime !< Actual time
551  real(r_kind),dimension(net_size),intent(in) :: y !< array of abundances (Y_i)
552  real(r_kind),dimension(track_nuclei_nr) :: tmp_y !< Stores the abundance of the tracked nuclei
553  integer :: i !< Loop variable
554 
555  info_entry("output_track_nuclei")
556  ! Store the abundance temporary
557  do i=1, track_nuclei_nr
558  tmp_y(i) = y(track_nuclei_indices(i))
559  if (tmp_y(i) .lt. 1d-99) tmp_y(i) = 0.
560  end do
561  ! Output
562  write(track_unit,'(*(es23.16,5x))') ctime,tmp_y(:)
563  info_exit("output_track_nuclei")
564 
565 end subroutine output_track_nuclei
566 
567 
568 !> Output the energy that enters and leaves the system
569 !!
570 !! Outputs the energy that enters and leaves the system when
571 !! nuclear_heating is turned on.
572 !! This output is controlled by the parameter \ref nu_loss_every .
573 !! The output is written to the file "nu_loss_gain_file".
574 !! An example of the file could look like:
575 !! \file{
576 !! \# File containing the neutrino loss/gain for the nuclear heating.
577 !! \# Negative values mean that energy is added to the system.
578 !! \# Neutrino energies are given in MeV/baryon/s.
579 !! \# time[s] T[GK] rho[g/cm3] R[km] nu_total ...
580 !! 1.10385E+00 1.00000E+01 4.89334E+07 3.55626E+02 0.00000E+00 ...
581 !! 1.10385E+00 1.00000E+01 4.89334E+07 3.55626E+02 -5.05524E+00 ...
582 !! 1.10385E+00 1.00000E+01 4.89334E+07 3.55626E+02 -5.05524E+00 ...
583 !! 1.10385E+00 1.00000E+01 4.89334E+07 3.55626E+02 -5.05524E+00 ...
584 !! ...}
585 !! @author Moritz Reichert
586 !! @date 12.04.2023
587 subroutine output_nu_loss(ctime,temp,rho,Rkm)
588  use single_zone_vars,only: stepsize
590  use parameter_class, only: unit
591  implicit none
592  real(r_kind),intent(in) :: ctime !< actual time (s)
593  real(r_kind),intent(in) :: temp !< temperature (GK)
594  real(r_kind),intent(in) :: rho !< density (g/cm3)
595  real(r_kind),intent(in) :: rkm !< radius (km)
596  ! local variables
597  real(r_kind) :: qdot_tot
598 
599  qdot_tot = qdot_bet + qdot_nu + qdot_th
600 
601  write(nu_loss_gain_unit, '(*(es19.5,3x))') ctime, temp, rho, rkm, &
602  qdot_tot/unit%hix/stepsize, &
603  qdot_bet/unit%hix/stepsize, &
604  qdot_nu/unit%hix/stepsize, &
605  qdot_th/unit%hix/stepsize
606 
607 end subroutine output_nu_loss
608 
609 
610 !>
611 !! Output Mainout file
612 !!
613 !! This file contains all basic quantities such as neutron abundance Yn,
614 !! time, temperature, density, Abar, and much more. An example file
615 !! looks like:
616 !! \l_file{
617 !! # 1:iteration 2:time[s], 3:T[GK], 4:rho[g/cm3], 5:Ye
618 !! # 6:R[km], 7:Y_n, 8:Y_p, 9:Y_alpha, 10:Y_lights
619 !! # 11:Y_heavies 12:<Z> 13:<A> 14:entropy [kB/baryon] (15:Sn [MeV])
620 !! 0 0.0000000000000000E+00 1.0964999999999998E+01 8.7095999999999795E+12 3.4880000000001680E-02 4.9272999999999975E+01 8.7430471311642388E-01 1.8154340816658842E-15 1.8313910002116333E-13 8.4293473893924100E-09 1.6154109228777289E-03 3.9820982195819934E-02 1.1416565996507526E+00 9.4638587929828290E-03
621 !! 10 4.0940000000000001E-09 1.0965038071994945E+01 8.7093701186135752E+12 3.4880006806212963E-02 4.9273176232580205E+01 8.7430472419630612E-01 1.8157922458135375E-15 1.8317495847108614E-13 8.4302150496944692E-09 1.6154086944791908E-03 3.9820989563731306E-02 1.1416565881127740E+00 9.4639764072160116E-03
622 !! ...}
623 !!
624 !!
625 subroutine output_mainout(cnt,ctime,T9,rho_b,entropy,Rkm,Y)
627  use global_class, only: isotope
628  implicit none
629 
630  integer,intent(in) :: cnt !< current iteration counter
631  real(r_kind),intent(in) :: ctime !< current time
632  real(r_kind),intent(in) :: t9 !< temperature [GK]
633  real(r_kind),intent(in) :: rho_b !< density [gcc]
634  real(r_kind),intent(in) :: entropy !< entropy [kB/nucleon]
635  real(r_kind),intent(in) :: rkm !< length scale [km]
636  real(r_kind),dimension(net_size),intent(in) :: y !< array of abundances (Y_i)
637  !
638  real(r_kind) :: ye,yneu,ypro,yhe4,ylight,yheavies,ysum,yasum,yzsum,abar,zbar
639 
640  ! calculate neutron separation energy
641  if (calc_nsep_energy) then
643  endif
644 
645  ! output diagnostics
646  ye = sum(isotope(1:net_size)%p_nr*y(1:net_size))
647  ylight= max(1.0d-99,sum(y(ipro+1:ihe4-1)))
648  yheavies= 0.0
649  if (ihe4.ge.1) yheavies= max(1.0d-99,sum(y(ihe4+1:net_size)))
650 
651  ysum = sum(y(1:net_size))
652  yasum = sum(y(1:net_size)*isotope(1:net_size)%mass)
653  yzsum = sum(y(1:net_size)*isotope(1:net_size)%p_nr)
654  abar = yasum/ysum
655  zbar = yzsum/ysum
656  !zbar = yzsum/sum(Y(2:net_size)) !as neutrons have Z=0
657 
658  yneu= 0.0
659  if (ineu.ge.1) yneu= max(1.0d-99,y(ineu))
660 
661  ypro= 0.0
662  if (ipro.ge.1) ypro= max(1.0d-99,y(ipro))
663 
664  yhe4= 0.0
665  if (ihe4.ge.1) yhe4= max(1.0d-99,y(ihe4))
666 
667  if (calc_nsep_energy) then
668  write (mainout_unit,'(i8,X,100(es23.16,2X))') cnt,ctime,t9,rho_b,ye,rkm, &
669  yneu,ypro,yhe4,ylight,yheavies,zbar,abar,entropy,sn_ave
670  else
671  write (mainout_unit,'(i8,X,100(es23.16,2X))') cnt,ctime,t9,rho_b,ye,rkm, &
672  yneu,ypro,yhe4,ylight,yheavies,zbar,abar,entropy
673  endif
674 
675 end subroutine output_mainout
676 
677 
678 !>
679 !! Output Snapshot files
680 !!
681 !! Snapshots store the abundances and mass fractions of all nuclei
682 !! at certain times. An example file looks like:
683 !! \file{
684 !! time temp dens
685 !! 9.77198412010000E-02 9.99927332947922E+00 3.88380373937259E+06
686 !! nin zin y x
687 !! 1 0 5.72194135572689E-01 5.72194135572689E-01
688 !! 0 1 3.92194484636783E-01 3.92194484636783E-01
689 !! 1 1 7.83033243870819E-05 1.56606648774164E-04
690 !! ...}
691 !!
692 !! \b Edited:
693 !! - 03.04.18, M. Jacobi
694 !! .
695 subroutine output_snapshot(ctime,T9,rho_b,Y)
696  implicit none
697 
698  real(r_kind),intent(in) :: ctime !< current time
699  real(r_kind),intent(in) :: t9 !< temperature [GK]
700  real(r_kind),intent(in) :: rho_b !< density [gcc]
701  real(r_kind),dimension(net_size),intent(in):: y !< array of abundances Y_i
702 
703 
704  ! Output Snapshots
705  call printsnap (ctime,t9,rho_b,y)
706 
707 end subroutine output_snapshot
708 
709 
710 !>
711 !! Output flow files
712 !!
713 !! \b Edited:
714 !! - 19.02.20, M. Jacobi
715 subroutine output_flow(ctime,T9,rho_b,Y)
716  use flow_module, only: flowprint, flowcalc
717  implicit none
718 
719  real(r_kind),intent(in) :: ctime !< current time
720  real(r_kind),intent(in) :: T9 !< temperature [GK]
721  real(r_kind),intent(in) :: rho_b !< density [gcc]
722  real(r_kind),dimension(net_size),intent(in):: Y !< array of abundances Y_i
723 
724  call flowcalc(y)
725  call flowprint(ctime,t9,rho_b,y,flowcnt)
726  flowcnt = flowcnt + 1
727 end subroutine output_flow
728 
729 
730 !>
731 !! Output of the Newton-Raphson loop diagnostics
732 !!
733 !! The diagnostic output contains NR iterations, timesteps,
734 !! the nucleus that restricts the timestep and more.
735 !! An example file looks like:
736 !! \l_file{
737 !! # Newton-Raphson diagnostic output
738 !! # 1 : global iteration count (cnt)
739 !! # 2 : global time [s]
740 !! # 3 : temperature T9 [GK]
741 !! # 4 : adapt stepsize loop counter (k)
742 !! # 5 : NR loop counter (nr_count)
743 !! # 6 : global timestep [s]
744 !! # 7 : mass conservation error
745 !! # 8 : maximal abundance change (eps)
746 !! # 9 : most rapidly evolving isotope (epsl)
747 !! # 10: abundance of the isotope epsl, y(epsl)
748 !! 0 5.20740000000E-22 7.42450000000E+00 0 1 5.20740000000E-22 9.32587340685E-15 2.22044604925E-16 cl42 2.43899441583E-08
749 !! 0 5.20740000000E-22 7.42450000000E+00 0 2 5.20740000000E-22 9.32587340685E-15 2.22044604925E-16 s40 2.73320077046E-08
750 !! 1 1.56220000000E-21 7.42450000000E+00 0 1 1.04146000000E-21 9.32587340685E-15 2.22044604925E-16 he4 3.13811866442E-03
751 !! 1 1.56220000000E-21 7.42450000000E+00 0 2 1.04146000000E-21 9.32587340685E-15 2.22044604925E-16 mg30 1.57736344914E-10
752 !! ...}
753 subroutine output_nr_diagnostic(cnt,k,nr_c,ctime,T9,stepsize,mtot,Y_p,Y)
754  use parameter_class, only: nrdiag_every
755  use global_class, only: isotope
756  implicit none
757  integer,intent(in) :: cnt !< global iteration counter
758  integer,intent(in) :: k !< counter in adapt_stepsize loop
759  integer,intent(in) :: nr_c !< counter in the NR loop
760  real(r_kind),intent(in) :: ctime !< global current time
761  real(r_kind),intent(in) :: t9 !< current temperature [GK]
762  real(r_kind),intent(in) :: stepsize !< current step size
763  real(r_kind),intent(in) :: mtot !< total mass
764  real(r_kind),dimension(net_size) :: y_p,y !< old/new abundances
765  !
766  integer :: i
767  real(r_kind) :: eps,epst
768  integer :: epsl
769 
770  info_entry("output_nr_diagnostic")
771 
772  ! check periodicity condition
773  if (nrdiag_every.gt.0) then
774  if (mod(cnt,nrdiag_every) .ne. 0) return
775  else
776  return
777  endif
778 
779  ! which isotope has changed the most
780  eps = 0.d0
781  epsl= 1
782  do i=1,net_size
783  if(y_p(i).lt.1.d-10) cycle
784  if(dabs(y_p(i)-y(i)).eq.0.d0) cycle
785  epst = dabs(y(i)/y_p(i)-1.d0)
786  if (epst.gt.eps) then
787  eps = epst
788  epsl = i
789  end if
790  end do
791 
792  ! output convergence diagnostic
793  write(nrdiag_unit, &
794  '(i7,2(1x,es19.11),2(1x,i4),3(1x,es19.11),a5,x,es19.11)') &
795  cnt,ctime,t9,k,nr_c,stepsize,dabs(mtot-1d0),eps,isotope(epsl)%name, &
796  y(epsl)
797 
798 
799  info_exit("output_nr_diagnostic")
800 
801 end subroutine output_nr_diagnostic
802 
803 
804 !>
805 !! Print final statistics
806 !!
807 !! An example looks like:
808 !! \file{
809 !! Final time reached.
810 !! Final time: 2.000000000000E+17
811 !! Total number of iterations: 5334
812 !! Elapsed time [s]: 5.41572E+02
813 !! }
814 !!
815 !! \b Edited:
816 !! - 2023-02-06, M.R. took care of wrapping of the simulation time
817 !! that could lead to negative simulation times.
818 !! .
819 subroutine output_final_stats(cnt,ctime)
820  implicit none
821  integer,intent(in) :: cnt
822  real(r_kind),intent(in) :: ctime
823  !
824  real(r_kind) :: simtime
825 
826  info_entry('output_final_stats')
827 
828  ! calculate elapsed time, check wrapping at cl_cmax
829  if (cl_end .gt. cl_start) then
830  simtime= float(cl_end-cl_start)/float(cl_rate)
831  else
832  simtime= (float(cl_end)+float(cl_cmax)-float(cl_start))/float(cl_rate)
833  end if
834 
835  ! output final statistics
836  print '(A,ES19.12)', 'Final time: ', ctime
837  print '(A,I10)', 'Total number of iterations: ', cnt
838  print '(A,ES12.5)', 'Elapsed time [s]: ', simtime
839 
840  info_exit('output_final_stats')
841 
842 end subroutine output_final_stats
843 
844 
845 !>
846 !! Print current snapshot
847 !!
848 !! Snapshots store the abundances and mass fractions of all nuclei
849 !! at certain times. An example file looks like:
850 !! \file{
851 !! time temp dens
852 !! 9.77198412010000E-02 9.99927332947922E+00 3.88380373937259E+06
853 !! nin zin y x
854 !! 1 0 5.72194135572689E-01 5.72194135572689E-01
855 !! 0 1 3.92194484636783E-01 3.92194484636783E-01
856 !! 1 1 7.83033243870819E-05 1.56606648774164E-04
857 !! ...}
858 subroutine printsnap (t,T9,rho_b,Y)
859  use global_class, only: isotope
860  implicit none
861 
862  real(r_kind),intent(in) :: t,t9,rho_b
863  real(r_kind),dimension(net_size),intent(in) :: y
864  character (len=21) :: snapfile
865  integer :: snapunit, z
866  real(r_kind) :: xz
867 
868  info_entry("printsnap")
869 
870  write(snapfile,'("snaps/snapsh_",i4.4,".dat")') snapcnt
871  snapunit= open_outfile(snapfile)
872 
873  write(snapunit,'(3a8)')'time','temp','dens'
874  write(snapunit,'(3es23.14)')t,t9,rho_b
875  write(snapunit,'(4(a8))')'nin','zin','y','x'
876  do z=1,net_size
877  xz= y(z)*isotope(z)%mass
878  write(snapunit,'(2(i3,2x),2(es23.14))') &
879  isotope(z)%n_nr, isotope(z)%p_nr, max(1d-99,y(z)), max(1d-99,xz)
880  end do
881  call close_io_file(snapunit,snapfile)
882  snapcnt= snapcnt + 1
883  tsnp= t
884  t9snp= t9
885  rosnp=rho_b
886 
887  info_exit("printsnap")
888 
889 end subroutine printsnap
890 
891 
892 !>
893 !! Calculates average neutron separation energy
894 !!
895 !! This neutron separation energy is printed in the mainout when
896 !! parameter_class::calc_nsep_energy was set to "yes" and
897 !! a valid file for parameter_class::nsep_energies_file.
899  use global_class, only: isotope
900  implicit none
901 
902  real(r_kind),dimension(net_size),intent(in) :: y !< abundances
903  !
904  real(r_kind) :: ysum,yasum,yzsum,abar_ions,zbar_ions
905  integer :: i
906 
907  info_entry("calc_nseparation_energy")
908 
909  abar_ions = 0d0
910  zbar_ions = 0d0
911  do i=max(1,ihe4),net_size
912  abar_ions = abar_ions + isotope(i)%mass * y(i)
913  zbar_ions = zbar_ions + isotope(i)%p_nr * y(i)
914  !Sn_ave = Sn_ave + Sn(i) * Y(i)
915  enddo
916  ysum = sum(y(max(1,ihe4):net_size))
917  abar_ions = abar_ions / ysum
918  zbar_ions = zbar_ions / ysum
919  ysum = 0d0
920  sn_ave = 0d0
921  do i=max(1,ihe4),net_size
922  if (sn(i) /= 0) then
923  if ((isotope(i)%mass.ge.220).and.(isotope(i)%mass.le.260)) then
924  sn_ave = sn_ave + sn(i) * y(i)
925  ysum = ysum + y(i)
926  end if
927  endif
928  enddo
929  sn_ave = sn_ave / ysum
930  ysum = sum(y)
931  yasum = sum(y*isotope%mass)
932  yzsum = sum(y*isotope%p_nr)
933 
934  info_exit("calc_nseparation_energy")
935 
936 end subroutine calc_nseparation_energy
937 
938 
939 
940 !>
941 !! @brief: Calculate the generated energy for each class of reactions separately.
942 !!
943 !! This subroutine calculates the generated energy of the nuclear reactions
944 !! with the help of the Q-values. This energy is written to an output only.
945 !! It is not fed back to the temperature and has therefore
946 !! no impact on the final result. Also, depending on the input the
947 !! subroutine can calculate the top contributors to the energy generation.
948 !! These top contributors can be reactions and isotopes.
949 !!
950 !! @see \ref nuclear_heating
951 !!
952 !! @author: Moritz Reichert
953 !! @date 29.02.20
954 !!
955 !! \b Edited:
956 !! - 08.05.24: M.R. Moved the top list calculation from the nuclear heating to this subroutine.
957 !! .
958 subroutine calc_nuclear_heating(hdf5_mode, write_engen, write_toplist)
961  use nucstuff_class, only: isotope,inter_partf
964  use mergesort_module, only: quicksort
965  use benam_class, only: reaction_string
966  use hdf5_module
967  implicit none
968  ! Declare input variables
969  logical, intent(in) :: hdf5_mode !< whether to store as hdf5 or ascii
970  logical, intent(in) :: write_engen !< whether to write the energy generation
971  logical, intent(in) :: write_toplist!< whether to write the energy toplist
972  ! Local variables
973  type(reactionrate_type) :: rr_tmp !< auxiliary shorthand for reaction
974  type(fissionrate_type) :: fr_tmp !< auxiliary shorthand for fission reaction
975  character*4 :: src !< reaction source
976  integer :: i,j !< Loop variable
977  integer :: help_parts !< Helper variable
978  real(r_kind) :: rat !< storage of the cached reaction rate
979  ! Heating variables
980  real(r_kind) :: engen_q_value !< All reactions
981  real(r_kind) :: engen_q_weak !< weak reactions
982  real(r_kind) :: engen_q_ng,engen_q_gn !< gamma-n reactions
983  real(r_kind) :: engen_q_pg,engen_q_gp !< gamma-p reactions
984  real(r_kind) :: engen_q_ga,engen_q_ag !< gamma-alpha reactions
985  real(r_kind) :: engen_q_pa,engen_q_ap !< p-alpha reactions
986  real(r_kind) :: engen_q_pn,engen_q_np !< p-n reactions
987  real(r_kind) :: engen_q_na,engen_q_an !< n-alpha reactions
988  real(r_kind) :: engen_fiss !< fission reactions
989  real(r_kind) :: engen_tot_ncap, engen_tot_pcap, engen_tot_acap, &
990  engen_tot_npcap,engen_tot_nacap, engen_tot_pacap !< net energy
991  real(r_kind) :: engen_tot_exc !< total energy with mass excess formula
992  real(r_kind) :: entropy_src_term !< Source term of entropy * kBT
993  real(r_kind) :: mic2,gi,twopi_h2c2,mui!< Helper variables
994  real(r_kind) :: kbt !< Helper variables
995  real(r_kind),dimension(0:net_size)::pf!< partition functions
996  real(r_kind) :: etaele,etapos !< electron and positron chemical potental
997  real(r_kind) :: mue !< electron chemical potental [MeV]
998  ! Helpers
999  real(r_kind) :: tmp_energy !< Temporary energy per reaction
1000  real(r_kind),dimension(net_size) ::det_bet_engen !< Energy of decays per parent nucleus
1001  real(r_kind),dimension(net_size) ::det_ag_engen !< Energy of (a,g)+(g,a) per parent nucleus
1002  real(r_kind),dimension(net_size) ::det_pg_engen !< Energy of (p,g)+(g,p) per parent nucleus
1003  real(r_kind),dimension(net_size) ::det_ng_engen !< Energy of (n,g)+(g,n) per parent nucleus
1004  real(r_kind),dimension(net_size) ::det_pn_engen !< Energy of (p,n)+(n,p) per parent nucleus
1005  real(r_kind),dimension(net_size) ::det_ap_engen !< Energy of (a,p)+(p,a) per parent nucleus
1006  real(r_kind),dimension(net_size) ::det_an_engen !< Energy of (a,n)+(n,a) per parent nucleus
1007  real(r_kind),dimension(net_size) ::det_other_engen!< Energy of other reactions per parent nucleus
1008  real(r_kind),dimension(net_size) ::det_fiss_engen !< Energy of fission per parent nucleus
1009  real(r_kind),dimension(net_size) ::det_tot_engen !< Total energy per parent nucleus
1010  ! Toplist stuff
1011  integer, dimension(nreac) ::toplist_idx !< Toplist indices of reactions
1012  real(r_kind), dimension(nreac) ::toplist_energy !< Energy of reactions
1013  logical, dimension(nreac) ::toplist_mask !< Mask which energy has been calculated
1014  character*30 ::reac_dummy !< reaction source
1015  character*30, dimension(toplist_amount)::tplist_str !< Reduced list of top energies
1016  real(r_kind), dimension(toplist_amount)::tplist !< Reduced list of top energies
1017  character*30, dimension(toplist_amount)::btlist_str !< Reduced list of top energies
1018  real(r_kind), dimension(toplist_amount)::btlist !< Reduced list of bottom energies
1019  real(r_kind), dimension(net_size) ::toplist_energy_iso !< Reduced list of bottom energies
1020  integer, dimension(net_size) ::toplist_idx_iso !< Indices of isotopes
1021  character*5, dimension(toplist_amount) ::tplist_str_iso !< Reduced list of top energies
1022  real(r_kind), dimension(toplist_amount)::tplist_iso !< Reduced list of top energies
1023  character*5, dimension(toplist_amount) ::btlist_str_iso !< Reduced list of top energies
1024  real(r_kind), dimension(toplist_amount)::btlist_iso !< Reduced list of bottom energies
1025 
1026  info_entry("calc_nuclear_heating")
1027 
1028  ! Calculate the heating with help of the q-value
1029  engen_q_value = 0 ! Total energy [erg/g/s]
1030  engen_q_weak = 0 ! Energy generated by weak reactions [erg/g/s]
1031  engen_q_ng = 0 ! Energy generated by n-gamma reactions [erg/g/s]
1032  engen_q_gn = 0 ! Energy generated by gamma-n reactions [erg/g/s]
1033  engen_q_pg = 0 ! Energy generated by p-gamma reactions [erg/g/s]
1034  engen_q_gp = 0 ! Energy generated by gamma-p reactions [erg/g/s]
1035  engen_q_ga = 0 ! Energy generated by gamma-alpha reactions [erg/g/s]
1036  engen_q_ag = 0 ! Energy generated by alpha-gamma reactions [erg/g/s]
1037  engen_q_np = 0 ! Energy generated by n-p reactions [erg/g/s]
1038  engen_q_pn = 0 ! Energy generated by p-n reactions [erg/g/s]
1039  engen_q_na = 0 ! Energy generated by n-alpha reactions [erg/g/s]
1040  engen_q_an = 0 ! Energy generated by alpha-n reactions [erg/g/s]
1041  engen_q_pa = 0 ! Energy generated by p-alpha reactions [erg/g/s]
1042  engen_q_ap = 0 ! Energy generated by alpha-p reactions [erg/g/s]
1043  engen_fiss = 0 ! Energy generated by fission [erg/g/s]
1044 
1045  ! Initialize decay energy
1046  if ((h_engen_detailed) .and. (write_engen)) then
1047  det_bet_engen(:) = 0.0
1048  det_ag_engen(:) = 0.0
1049  det_pg_engen(:) = 0.0
1050  det_ng_engen(:) = 0.0
1051  det_pn_engen(:) = 0.0
1052  det_ap_engen(:) = 0.0
1053  det_an_engen(:) = 0.0
1054  det_other_engen(:) = 0.0
1055  det_fiss_engen(:) = 0.0
1056  det_tot_engen(:) = 0.0
1057  end if
1058 
1059  ! Initialize the energy for the top contributors
1060  ! This list has one entry per reaction
1061  if (write_toplist) then
1062  toplist_energy(:) = 0d0
1063  toplist_energy_iso(:) = 0d0
1064  toplist_mask(:) = .false.
1065  end if
1066 
1067  ! For NSE the reaction rates are not cached and the generated energy would be therefore weird
1068  ! Set it to zero for this case
1069  ! Loop through all reactions
1070  do i=1,nreac
1071  rr_tmp = rrate(i)
1072  src = rr_tmp%source
1073  rat = rr_tmp%cached
1074 
1075  if ((evolution_mode .eq. em_nse) .and. (.not. rr_tmp%is_weak)) cycle
1076 
1077  ! Do nothing if the rate is zero
1078  if (rat .eq. 0) cycle
1079 
1080  ! dYdt depends on the reaclib chapter and we therefore have to separate them
1081  select case (rr_tmp%group)
1082  case(1:3,11)! Reaclib chapter 1-3 and 11
1083  ! Energy per reaction
1084  tmp_energy = rat*y(rr_tmp%parts(1)) * rr_tmp%q_value
1085  ! Total energy
1086  engen_q_value = engen_q_value + tmp_energy
1087  ! Seperate the reaction types, weak reactions and photodissociations should be in chapter 1-3
1088 
1089  ! weak rates (only contained in chapter 1-3)
1090  if (rr_tmp%is_weak) then
1091  engen_q_weak = engen_q_weak + tmp_energy
1092  if ((h_engen_detailed) .and. (write_engen)) then
1093  det_bet_engen(rr_tmp%parts(1)) = det_bet_engen(rr_tmp%parts(1))+tmp_energy
1094  end if
1095  end if
1096  ! (gamma,n) reactions
1097  if(rr_tmp%reac_type.eq.rrt_gn) then
1098  engen_q_gn = engen_q_gn + tmp_energy
1099  ! Detailed output
1100  if ((h_engen_detailed) .and. (write_engen)) then
1101  if (rr_tmp%parts(2) .eq. ineu) then
1102  help_parts = rr_tmp%parts(3)
1103  else
1104  help_parts = rr_tmp%parts(2)
1105  end if
1106  det_ng_engen(help_parts) = det_ng_engen(help_parts)+tmp_energy
1107  end if
1108 
1109  endif
1110  ! (gamma,p) reactions
1111  if(rr_tmp%reac_type.eq.rrt_gp) then
1112  engen_q_gp = engen_q_gp + tmp_energy
1113  ! Detailed output
1114  if ((h_engen_detailed) .and. (write_engen)) then
1115  if (rr_tmp%parts(2) .eq. ipro) then
1116  help_parts = rr_tmp%parts(3)
1117  else
1118  help_parts = rr_tmp%parts(2)
1119  end if
1120  det_pg_engen(help_parts) = det_pg_engen(help_parts)+tmp_energy
1121  end if
1122  endif
1123  ! (gamma,alpha) reactions
1124  if(rr_tmp%reac_type.eq.rrt_ga) then
1125  engen_q_ga = engen_q_ga + tmp_energy
1126  ! Detailed output
1127  if ((h_engen_detailed) .and. (write_engen)) then
1128  if (rr_tmp%parts(2) .eq. ihe4) then
1129  help_parts = rr_tmp%parts(3)
1130  else
1131  help_parts = rr_tmp%parts(2)
1132  end if
1133  det_ag_engen(help_parts) = det_ag_engen(help_parts)+tmp_energy
1134  end if
1135  endif
1136 
1137 
1138  if(rr_tmp%reac_type.eq.rrt_o) then
1139  if ((h_engen_detailed) .and. (write_engen)) then
1140  det_other_engen(rr_tmp%parts(1)) = det_other_engen(rr_tmp%parts(1))+tmp_energy
1141  end if
1142  end if
1143 
1144  case(4:7)! Reaclib chapter 4-7
1145  ! Energy per reaction
1146  tmp_energy = rat*y(rr_tmp%parts(1)) * y(rr_tmp%parts(2)) &
1147  * rr_tmp%q_value
1148  ! Total energy
1149  engen_q_value = engen_q_value + tmp_energy
1150 
1151  ! Seperate the reaction types
1152  ! (n-gamma) reactions
1153  if(rr_tmp%reac_type.eq.rrt_ng) then
1154  engen_q_ng = engen_q_ng + tmp_energy
1155  ! Detailed output
1156  if ((h_engen_detailed) .and. (write_engen)) then
1157  if (rr_tmp%parts(2) .eq. ineu) then
1158  help_parts = rr_tmp%parts(1)
1159  else
1160  help_parts = rr_tmp%parts(2)
1161  end if
1162  det_ng_engen(help_parts) = det_ng_engen(help_parts)+tmp_energy
1163  end if
1164  endif
1165  ! (p-gamma) reactions
1166  if(rr_tmp%reac_type.eq.rrt_pg) then
1167  engen_q_pg = engen_q_pg + tmp_energy
1168  ! Detailed output
1169  if ((h_engen_detailed) .and. (write_engen)) then
1170  if (rr_tmp%parts(2) .eq. ipro) then
1171  help_parts = rr_tmp%parts(1)
1172  else
1173  help_parts = rr_tmp%parts(2)
1174  end if
1175  det_pg_engen(help_parts) = det_pg_engen(help_parts)+tmp_energy
1176  end if
1177  endif
1178  ! (alpha gamma) reactions
1179  if(rr_tmp%reac_type.eq.rrt_ag) then
1180  engen_q_ag = engen_q_ag + tmp_energy
1181  ! Detailed output
1182  if ((h_engen_detailed) .and. (write_engen)) then
1183  if (rr_tmp%parts(2) .eq. ihe4) then
1184  help_parts = rr_tmp%parts(1)
1185  else
1186  help_parts = rr_tmp%parts(2)
1187  end if
1188  det_ag_engen(help_parts) = det_ag_engen(help_parts)+tmp_energy
1189  end if
1190  endif
1191  ! (p,alpha) reactions
1192  if(rr_tmp%reac_type.eq.rrt_pa) then
1193  engen_q_pa = engen_q_pa + tmp_energy
1194  ! Detailed output
1195  if ((h_engen_detailed) .and. (write_engen)) then
1196  if (rr_tmp%parts(3) .eq. ihe4) then
1197  help_parts = rr_tmp%parts(4)
1198  else
1199  help_parts = rr_tmp%parts(3)
1200  end if
1201  det_ap_engen(help_parts) = det_ap_engen(help_parts)+tmp_energy
1202  end if
1203  endif
1204  ! (alpha,p) reactions
1205  if(rr_tmp%reac_type.eq.rrt_ap) then
1206  engen_q_ap = engen_q_ap + tmp_energy
1207  ! Detailed output
1208  if ((h_engen_detailed) .and. (write_engen)) then
1209  if (rr_tmp%parts(2) .eq. ihe4) then
1210  help_parts = rr_tmp%parts(1)
1211  else
1212  help_parts = rr_tmp%parts(2)
1213  end if
1214  det_ap_engen(help_parts) = det_ap_engen(help_parts)+tmp_energy
1215  end if
1216  endif
1217  ! (n,p) reactions
1218  if(rr_tmp%reac_type.eq.rrt_np) then
1219  engen_q_np = engen_q_np + tmp_energy
1220  ! Detailed output
1221  if ((h_engen_detailed) .and. (write_engen)) then
1222  if (rr_tmp%parts(3) .eq. ipro) then
1223  help_parts = rr_tmp%parts(4)
1224  else
1225  help_parts = rr_tmp%parts(3)
1226  end if
1227  det_pn_engen(help_parts) = det_pn_engen(help_parts)+tmp_energy
1228  end if
1229  endif
1230  ! (p,n) reactions
1231  if(rr_tmp%reac_type.eq.rrt_pn) then
1232  engen_q_pn = engen_q_pn + tmp_energy
1233  ! Detailed output
1234  if ((h_engen_detailed) .and. (write_engen)) then
1235  if (rr_tmp%parts(1) .eq. ipro) then
1236  help_parts = rr_tmp%parts(2)
1237  else
1238  help_parts = rr_tmp%parts(1)
1239  end if
1240  det_pn_engen(help_parts) = det_pn_engen(help_parts)+tmp_energy
1241  end if
1242  endif
1243  ! (alpha,n) reactions
1244  if(rr_tmp%reac_type.eq.rrt_an) then
1245  engen_q_an = engen_q_an + tmp_energy
1246  ! Detailed output
1247  if ((h_engen_detailed) .and. (write_engen)) then
1248  if (rr_tmp%parts(1) .eq. ihe4) then
1249  help_parts = rr_tmp%parts(2)
1250  else
1251  help_parts = rr_tmp%parts(1)
1252  end if
1253  det_an_engen(help_parts) = det_an_engen(help_parts)+tmp_energy
1254  end if
1255  endif
1256  ! (n,alpha) reactions
1257  if(rr_tmp%reac_type.eq.rrt_na) then
1258  engen_q_na = engen_q_na + tmp_energy
1259  ! Detailed output
1260  if ((h_engen_detailed) .and. (write_engen)) then
1261  if (rr_tmp%parts(3) .eq. ihe4) then
1262  help_parts = rr_tmp%parts(4)
1263  else
1264  help_parts = rr_tmp%parts(3)
1265  end if
1266  det_an_engen(help_parts) = det_an_engen(help_parts)+tmp_energy
1267  end if
1268  endif
1269 
1270  if(rr_tmp%reac_type.eq.rrt_o) then
1271  if ((h_engen_detailed) .and. (write_engen)) then
1272  det_other_engen(rr_tmp%parts(1)) = det_other_engen(rr_tmp%parts(1))+tmp_energy
1273  end if
1274  end if
1275 
1276  case(8:9)! Reaclib chapter 8-9
1277  ! Energy per reaction
1278  tmp_energy = rat*y(rr_tmp%parts(1)) * y(rr_tmp%parts(2)) &
1279  *y(rr_tmp%parts(3)) * rr_tmp%q_value
1280 
1281  ! Total energy
1282  engen_q_value = engen_q_value + tmp_energy
1283 
1284  if(rr_tmp%reac_type.eq.rrt_o) then
1285  if ((h_engen_detailed) .and. (write_engen)) then
1286  det_other_engen(rr_tmp%parts(1)) = det_other_engen(rr_tmp%parts(1))+tmp_energy
1287  end if
1288  end if
1289 
1290  case(10)! Reaclib chapter 8
1291  ! Energy per reaction
1292  tmp_energy = rat*y(rr_tmp%parts(1)) * y(rr_tmp%parts(2)) &
1293  *y(rr_tmp%parts(3)) * y(rr_tmp%parts(4))* rr_tmp%q_value
1294  ! Total energy
1295  engen_q_value = engen_q_value + tmp_energy
1296 
1297  if(rr_tmp%reac_type.eq.rrt_o) then
1298  if ((h_engen_detailed) .and. (write_engen)) then
1299  det_other_engen(rr_tmp%parts(1)) = det_other_engen(rr_tmp%parts(1))+tmp_energy
1300  end if
1301  end if
1302 
1303 
1304  case default ! This should not happen
1305  cycle
1306  end select
1307 
1308  if (write_toplist) then
1309  ! Store the energy of the top contributors
1310  ! This list has one entry per reaction and will be later sorted
1311  if (tmp_energy .ne. 0d0) then
1312  toplist_energy(i) = tmp_energy
1313  toplist_mask(i) = .true.
1314  end if
1315  end if
1316  end do
1317 
1318  ! Do the same for the fission reactions
1319  if (evolution_mode .ne. em_nse) then
1320  if (fissflag.ne.0) then
1321  fission:do i=1,nfiss
1322  fr_tmp = fissrate(i)
1323  rat = fr_tmp%cached
1324  if (rat .eq. 0) cycle
1325  if ((fr_tmp%mode.eq.2).or.(fr_tmp%mode.eq.3)) then ! spont. and b-delayed fission
1326  engen_q_value = engen_q_value + rat*y(fissrate(i)%fissparts(1)) * fissrate(i)%averageQ
1327  engen_fiss = engen_fiss + rat*y(fissrate(i)%fissparts(1)) * fissrate(i)%averageQ
1328  tmp_energy = rat*y(fissrate(i)%fissparts(1)) * fissrate(i)%averageQ
1329  else if (fr_tmp%mode.eq.1) then ! n-induced fission
1330  engen_q_value = engen_q_value + rat*y(fissrate(i)%fissparts(1))*y(fissrate(i)%fissparts(2)) * fissrate(i)%averageQ
1331  engen_fiss = engen_fiss + rat*y(fissrate(i)%fissparts(1))*y(fissrate(i)%fissparts(2)) * fissrate(i)%averageQ
1332  tmp_energy = rat*y(fissrate(i)%fissparts(1))*y(fissrate(i)%fissparts(2)) * fissrate(i)%averageQ
1333  end if
1334 
1335  if ((h_engen_detailed) .and. (write_engen)) then
1336  det_fiss_engen(fissrate(i)%fissparts(1)) = det_fiss_engen(fissrate(i)%fissparts(1))+tmp_energy
1337  end if
1338 
1339  end do fission
1340  end if
1341  end if
1342 
1343 
1344  engen_tot_exc = 0
1345  do i=1,net_size
1346  if (y(i).gt.1d-25) then
1347  tmp_energy = -(isotope(i)%mass_exc)*(y(i)-y_p(i))
1348  engen_tot_exc = engen_tot_exc + tmp_energy ! energy generated by species i [MeV/baryon]
1349 
1350  if (write_toplist) then
1351  toplist_energy_iso(i) = tmp_energy
1352  end if
1353 
1354  if ((h_engen_detailed) .and. (write_engen)) then
1355  det_tot_engen(i) = -(isotope(i)%mass_exc)*(y(i)-y_p(i))
1356  end if
1357  endif
1358  end do
1359  engen_tot_exc = engen_tot_exc/stepsize / unit%hix
1360 
1361  entropy_src_term = 0.d0
1362  call inter_partf (t9, pf)
1363 
1364  ! Take care that one does not go off-grid
1365  if (rhob .gt. 1d-10) then
1366  ! Get the electron chemical potential
1367  call chempot(t9*1d9,rhob,ye,etaele,etapos)
1368  mue = etaele * unit%k_mev * t9*1d9
1369  else
1370  mue = 0.d0
1371  end if
1372  mue = mue + unit%mass_e
1373 
1374  kbt = t9*1d9 * unit%k_mev ! 8.617332d-11 * T [MeV]
1375  do i=1,net_size
1376  if (y(i).gt.1d-25) then
1377  mic2 = isotope(i)%mass*unit%mass_u - isotope(i)%p_nr * unit%mass_e + isotope(i)%mass_exc ! [MeV]
1378  gi = 2d0*isotope(i)%spin + 1d0 ! statistical weight
1379  twopi_h2c2 = 2d0*unit%pi/(unit%h_mevs*unit%clight)**2 ! 2\pi/(hc)^2 [MeV^-2 cm^-2]
1380  mui = -kbt*log( gi*pf(i) / (rhob*y(i)*unit%n_a) * &
1381  (twopi_h2c2 * mic2 * kbt)**1.5d0 )
1382  ! entropy generated by species i [MeV/baryon]
1383  entropy_src_term = entropy_src_term -(mic2 + mui + isotope(i)%p_nr*mue)*(y(i)-y_p(i))
1384  endif
1385  end do
1386 
1387  entropy_src_term = entropy_src_term/ unit%hix/stepsize
1388  ! Convert from MeV/Baryon/s - > Erg/g/s
1389  engen_q_value = engen_q_value / unit%hix
1390  engen_q_weak = engen_q_weak / unit%hix
1391  engen_q_ng = engen_q_ng / unit%hix
1392  engen_q_gn = engen_q_gn / unit%hix
1393  engen_q_pg = engen_q_pg / unit%hix
1394  engen_q_gp = engen_q_gp / unit%hix
1395  engen_q_ga = engen_q_ga / unit%hix
1396  engen_q_ag = engen_q_ag / unit%hix
1397  engen_q_np = engen_q_np / unit%hix
1398  engen_q_pn = engen_q_pn / unit%hix
1399  engen_q_na = engen_q_na / unit%hix
1400  engen_q_an = engen_q_an / unit%hix
1401  engen_q_pa = engen_q_pa / unit%hix
1402  engen_q_ap = engen_q_ap / unit%hix
1403  engen_fiss = engen_fiss / unit%hix
1404  ! Detailed output
1405  if ((h_engen_detailed) .and. (write_engen)) then
1406  det_tot_engen = det_tot_engen/stepsize / unit%hix
1407  det_bet_engen = det_bet_engen / unit%hix
1408  det_ag_engen = det_ag_engen / unit%hix
1409  det_pg_engen = det_pg_engen / unit%hix
1410  det_ng_engen = det_ng_engen / unit%hix
1411  det_pn_engen = det_pn_engen / unit%hix
1412  det_ap_engen = det_ap_engen / unit%hix
1413  det_an_engen = det_an_engen / unit%hix
1414  det_other_engen = det_other_engen / unit%hix
1415  det_fiss_engen = det_fiss_engen / unit%hix
1416  end if
1417 
1418  !Calculate net energy
1419  engen_tot_ncap = engen_q_ng + engen_q_gn
1420  engen_tot_pcap = engen_q_pg + engen_q_gp
1421  engen_tot_acap = engen_q_ag + engen_q_ga
1422  engen_tot_npcap= engen_q_np + engen_q_pn
1423  engen_tot_nacap= engen_q_na + engen_q_an
1424  engen_tot_pacap= engen_q_pa + engen_q_ap
1425 
1426  if (write_toplist) then
1427  ! Store the energy of the top contributors
1428  ! This list has one entry per reaction and will be later sorted
1429  toplist_energy(:) = toplist_energy(:) / unit%hix
1430  toplist_energy_iso(:) = toplist_energy_iso(:) /stepsize/ unit%hix
1431 
1432  ! Sort the list
1433  call quicksort(0,nreac,toplist_energy,toplist_idx)
1434 
1435  if (toplist_energy(1) .lt. toplist_energy(2)) then
1436  call raise_exception("There was a bug in the toplist. Energy is not sorted.", &
1437  "calc_nuclear_heating")
1438  end if
1439 
1440  i=0
1441  j=1
1442  tplist_str(:) = "---------"
1443  tplist(:) = 0d0
1444  do while(i<10)
1445 
1446  if (j .gt. nreac) exit
1447 
1448  if (toplist_mask(toplist_idx(j))) then
1449  reac_dummy=trim(adjustl(reaction_string(rrate(toplist_idx(j)))))
1450  tplist_str(i+1) = reaction_string(rrate(toplist_idx(j)))
1451  tplist(i+1) = toplist_energy(j)
1452  i = i+1
1453  else
1454  continue
1455  end if
1456 
1457  j = j+1
1458  end do
1459 
1460  if (tplist(1) .lt. tplist(2)) then
1461  call raise_exception("There was a bug in the toplist. Energy is not sorted.", &
1462  "calc_nuclear_heating")
1463  end if
1464 
1465  i=0
1466  j=1
1467  btlist_str(:) = "---------"
1468  btlist(:) = 0d0
1469  do while(i<10)
1470 
1471  if (j .gt. nreac) exit
1472 
1473  if (toplist_mask(toplist_idx(nreac+1-j))) then
1474  reac_dummy=trim(adjustl(reaction_string(rrate(toplist_idx(nreac+1-j)))))
1475  btlist_str(i+1) = reaction_string(rrate(toplist_idx(nreac+1-j)))
1476  btlist(i+1) = toplist_energy(nreac+1-j)
1477  i = i+1
1478  else
1479  continue
1480  end if
1481 
1482  j = j+1
1483  end do
1484 
1485 
1486  ! Do the same for the isotopes
1487  call quicksort(0,net_size,toplist_energy_iso,toplist_idx_iso)
1488 
1489  i=0
1490  j=1
1491  tplist_str_iso(:) = "-----"
1492  tplist_iso(:) = 0d0
1493  do while(i<10)
1494 
1495  if (j .gt. net_size) exit
1496 
1497  if (toplist_energy_iso(j) .ne. 0d0) then
1498  tplist_str_iso(i+1) = net_names(toplist_idx_iso(j))
1499  tplist_iso(i+1) = toplist_energy_iso(j)
1500  i = i+1
1501  else
1502  continue
1503  end if
1504 
1505  j = j+1
1506  end do
1507 
1508  i=0
1509  j=1
1510  btlist_str_iso(:) = "-----"
1511  btlist_iso(:) = 0d0
1512  do while(i<10)
1513 
1514  if (j .gt. net_size) exit
1515 
1516  if (toplist_energy_iso(net_size+1-j) .ne. 0d0) then
1517  btlist_str_iso(i+1) = net_names(toplist_idx_iso(net_size+1-j))
1518  btlist_iso(i+1) = toplist_energy_iso(net_size+1-j)
1519  i = i+1
1520  else
1521  continue
1522  end if
1523 
1524  j = j+1
1525  end do
1526 
1527 
1528 
1529  write(toplist_unit,"(es12.4, 4X, es12.4, 4X, "//int_to_str(toplist_amount)//"(A,3X), 3X,"//&
1530  int_to_str(toplist_amount)//"(es12.4,3X), 4X, "//int_to_str(toplist_amount)//&
1531  "(A,3X), 2X,"//int_to_str(toplist_amount)//"(es12.4,3X),' | ',"//&
1532  int_to_str(toplist_amount)//"(A,3X), 3X,"//int_to_str(toplist_amount)//&
1533  "(es12.4,3X), 4X, "//int_to_str(toplist_amount)//&
1534  "(A,3X), 3X,"//int_to_str(toplist_amount)//"(es12.4,3X))" )&
1535  time,engen_tot_exc, tplist_str, tplist, btlist_str,btlist, tplist_str_iso, tplist_iso, &
1536  btlist_str_iso, btlist_iso
1537 
1538  end if
1539 
1540 
1541 
1542  if ((.not. hdf5_mode) .and. write_engen) then
1543  ! Write the output (detailed)
1544  !write(nuc_heat_file,'(*(es23.14,3x))') time,engen_q_value,engen_q_weak,&
1545  ! engen_q_ng,engen_q_gn,engen_tot_ncap,&
1546  ! engen_q_pg,engen_q_gp,engen_tot_pcap,&
1547  ! engen_q_ag,engen_q_ga,engen_tot_acap,&
1548  ! engen_q_np,engen_q_pn,engen_tot_npcap,&
1549  ! engen_q_na,engen_q_an,engen_tot_nacap,&
1550  ! engen_q_pa,engen_q_ap,engen_tot_pacap,&
1551  ! engen_fiss
1552  ! Write the output (more relaxed)
1553  write(nuc_heat_file,'(*(es23.14E3,3x))') time,engen_tot_exc,entropy_src_term,&
1554  engen_q_weak,&
1555  engen_tot_ncap,&
1556  engen_tot_pcap,&
1557  engen_tot_acap,&
1558  engen_tot_npcap,&
1559  engen_tot_nacap,&
1560  engen_tot_pacap,&
1561  engen_fiss
1562  end if
1563 
1564 #ifdef USE_HDF5
1565  if ((write_engen) .and. (hdf5_mode)) then
1566  call extend_engen(time,engen_tot_exc, entropy_src_term,engen_tot_ncap, &
1567  engen_tot_pcap, engen_tot_acap,engen_tot_npcap, &
1568  engen_tot_nacap, engen_tot_pacap, engen_q_weak, &
1569  engen_fiss, &
1570  det_bet_engen, det_ag_engen, det_pg_engen, &
1571  det_ng_engen, det_pn_engen, det_ap_engen, &
1572  det_an_engen, det_other_engen, det_fiss_engen, &
1573  det_tot_engen, &
1575  end if
1576 #endif
1577 
1578 
1579  info_exit("calc_nuclear_heating")
1580 
1581 end subroutine calc_nuclear_heating
1582 
1583 
1584 
1585 
1586 
1587 
1588 !>
1589 !! Calculate average timescales of different classes of reactions
1590 !!
1591 !! These timescales are calculated as in
1592 !! [Arcones et al. 2012](https://ui.adsabs.harvard.edu/abs/2012ApJ...750...18A/abstract)
1593 !! (Eq. 1-6). The (p,gamma) timescale is given by e.g.,
1594 !! \f[
1595 !! \tau_{p\gamma} = \left[ \frac{\rho Y_p}{Y_h}N_A \sum \langle \sigma \nu \rangle_{p,\gamma}(Z,A)Y(Z,A) \right]^{-1}
1596 !! \f]
1597 !!
1598 !! An example file looks like:
1599 !! \l_file{
1600 !! # time[s] Temperature [GK] tau_ng [s] tau_gn [s] tau_pg [s] tau_gp [s] tau_np [s] tau_pn [s] tau_an [s] tau_na [s] tau_beta [s] tau_alpha [s] tau_nfiss [s] tau_sfiss [s] tau_bfiss [s]
1601 !! 1.1130000000000000E+00 9.1821750000000009E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 1.2642767370327558E+01 2.9487624276191804E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00
1602 !! 1.1220000000000001E+00 8.2638489999999951E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00 1.0557999180036711E+02 2.5923655028366777E+01 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00
1603 !! 1.1247618491767180E+00 7.9904974735555152E+00 9.5690472538317728E-13 5.3150487750907188E-11 7.4419112944811684E-10 9.8479480155439862E-08 1.9792689873245490E-12 1.9444248037144314E-12 4.8686627517806918E-10 3.9657685448682009E-12 1.4600540103519526E+01 5.2634156341152767E+01 0.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00
1604 !! ...}
1605 subroutine calc_av_timescales(ctime,T9,Y,hdf5)
1609 #ifdef USE_HDF5
1611 #endif
1612  implicit none
1613  real(r_kind),intent(in) :: ctime !< current time [s]
1614  real(r_kind),intent(in) :: t9 !< Temperature [GK]
1615  logical,optional,intent(in) :: hdf5 !< Whether the output is written to hdf5 or not
1616  real(r_kind),dimension(net_size),intent(in) :: y !< abundances
1617 
1618  ! average timescales
1619  real(r_kind) :: tau_ng !< (n,gamma)
1620  real(r_kind) :: tau_gn !< (gamma,n)
1621  real(r_kind) :: tau_pg !< (p,gamma)
1622  real(r_kind) :: tau_gp !< (gamma,p)
1623  real(r_kind) :: tau_np !< (n,p)
1624  real(r_kind) :: tau_pn !< (p,n)
1625  real(r_kind) :: tau_an !< (alpha,n)
1626  real(r_kind) :: tau_na !< (n,alpha)
1627  real(r_kind) :: tau_ag !< (alpha,gamma)
1628  real(r_kind) :: tau_ga !< (gamma,alpha)
1629  real(r_kind) :: tau_pa !< (p,alpha)
1630  real(r_kind) :: tau_ap !< (alpha,p)
1631  !real(r_kind) :: tau_npro !< neutron production timescale
1632  !real(r_kind) :: tau_fdis !< photodissociation timescale
1633  real(r_kind) :: tau_beta !< beta-decay timescale
1634  real(r_kind) :: tau_alpha !< alpha-decay timescale
1635  real(r_kind) :: tau_nfiss !< n-induced fission timescale
1636  real(r_kind) :: tau_sfiss !< spontaneous fission timescale
1637  real(r_kind) :: tau_bfiss !< beta-delayed fission timescale
1638  integer :: nuc !< auxiliary var
1639  real(r_kind) :: rat !< Value of reaction rate
1640  real(r_kind) :: ng_sum,gn_sum,pg_sum,gp_sum,beta_sum
1641  real(r_kind) :: np_sum,pn_sum,an_sum,na_sum,ag_sum,ga_sum
1642  real(r_kind) :: ap_sum,pa_sum
1643  real(r_kind) :: alpha_sum, nfiss_sum, sfiss_sum, bfiss_sum
1644  real(r_kind) :: ysum !< Sum of abundances heavier than helium
1645  integer :: i !< Loop variable
1646  type(reactionrate_type) :: rr_tmp !< auxiliary shorthand for reaction
1647  type(fissionrate_type) :: fr_tmp !< auxiliary shorthand for reaction
1648  character*4 :: src !< reaction source
1649  real(r_kind) :: q_sm = 1.e-20 !< small quantity to add to denominator
1650 
1651  info_entry("calc_av_timescales")
1652 
1653 
1654  ! to avoid division by zero
1655  ng_sum = q_sm
1656  gn_sum = q_sm
1657  pg_sum = q_sm
1658  gp_sum = q_sm
1659  np_sum = q_sm
1660  pn_sum = q_sm
1661  an_sum = q_sm
1662  na_sum = q_sm
1663  ag_sum = q_sm
1664  ga_sum = q_sm
1665  ap_sum = q_sm
1666  pa_sum = q_sm
1667  beta_sum = q_sm
1668  alpha_sum = q_sm
1669  nfiss_sum = q_sm
1670  sfiss_sum = q_sm
1671  bfiss_sum = q_sm
1672 
1673 
1674  ysum = sum(y(max(1,ihe4):net_size))
1675 
1676  do i=1,nreac
1677  rr_tmp = rrate(i)
1678  src = rr_tmp%source
1679  rat = rr_tmp%cached
1680 
1681  ! beta-decays (not including beta-delayed fission)
1682  if ((rr_tmp%group.eq.1).or.(rr_tmp%group.eq.2).or.(rr_tmp%group.eq.3) .or.(rr_tmp%group.eq.11)) then
1683  if ((rr_tmp%reac_type .eq.rrt_betm).or.(rr_tmp%reac_type .eq.rrt_betp).and.(rr_tmp%source.ne.'ms99')) then
1684  ! only include significant reactions
1685  nuc = rr_tmp%parts(1)
1686  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1687  beta_sum = beta_sum + rat*y(nuc)
1688  end if
1689  endif
1690  endif
1691 
1692  ! alpha-decays
1693  if ((rr_tmp%reac_type.eq.rrt_alpd)) then
1694  nuc = rr_tmp%parts(1)
1695  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1696  alpha_sum = alpha_sum + rat*y(nuc)
1697  end if
1698  end if
1699 
1700  ! (a,p)
1701  if(rr_tmp%reac_type.eq.rrt_ap) then
1702  nuc = rr_tmp%parts(2)
1703  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1704  ap_sum = ap_sum + rat*y(nuc)*y(ihe4)
1705  end if
1706  endif
1707 
1708  ! (p,a)
1709  if(rr_tmp%reac_type.eq.rrt_pa) then
1710  nuc = rr_tmp%parts(2)
1711  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1712  pa_sum = pa_sum + rat*y(nuc)*y(ipro)
1713  end if
1714  endif
1715 
1716  ! (n,gamma)
1717  if(rr_tmp%reac_type.eq.rrt_ng) then
1718  nuc = rr_tmp%parts(2)
1719  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1720  ng_sum = ng_sum + rat*y(nuc)*y(ineu)
1721  end if
1722  endif
1723 
1724  ! (gamma,n)
1725  if(rr_tmp%reac_type.eq.rrt_gn) then
1726  nuc = rr_tmp%parts(1)
1727  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1728  gn_sum = gn_sum + rat*y(nuc)
1729  end if
1730  endif
1731 
1732  ! (p,gamma)
1733  if(rr_tmp%reac_type.eq.rrt_pg) then
1734  nuc = rr_tmp%parts(2)
1735  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1736  pg_sum = pg_sum + rat*y(nuc)*y(ipro)
1737  end if
1738  endif
1739 
1740  ! (gamma,p)
1741  if(rr_tmp%reac_type.eq.rrt_gp) then
1742  nuc = rr_tmp%parts(1)
1743  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1744  gp_sum = gp_sum + rat*y(nuc)
1745  end if
1746  endif
1747 
1748  ! (a,g)
1749  if(rr_tmp%reac_type.eq.rrt_ag) then
1750  nuc = rr_tmp%parts(2)
1751  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1752  ag_sum = ag_sum + rat*y(nuc)*y(ihe4)
1753  end if
1754  endif
1755 
1756  ! (g,a)
1757  if(rr_tmp%reac_type.eq.rrt_ga) then
1758  nuc = rr_tmp%parts(1)
1759  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1760  ga_sum = ga_sum + rat*y(nuc)
1761  end if
1762  endif
1763 
1764  ! (n,p)
1765  if(rr_tmp%reac_type.eq.rrt_np) then
1766  nuc = rr_tmp%parts(2)
1767  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1768  np_sum = np_sum + rat*y(nuc)*y(ineu)
1769  end if
1770  endif
1771 
1772  ! (p,n)
1773  if(rr_tmp%reac_type.eq.rrt_pn) then
1774  nuc = rr_tmp%parts(2)
1775  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1776  pn_sum = pn_sum + rat*y(nuc)*y(ipro)
1777  end if
1778  endif
1779 
1780  ! (a,n)
1781  if(rr_tmp%reac_type.eq.rrt_an) then
1782  nuc = rr_tmp%parts(2)
1783  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1784  an_sum = an_sum + rat*y(nuc)*y(ihe4)
1785  end if
1786  endif
1787 
1788  ! (n,a)
1789  if(rr_tmp%reac_type.eq.rrt_na) then
1790  nuc = rr_tmp%parts(2)
1791  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1792  na_sum = na_sum + rat*y(nuc)*y(ineu)
1793  end if
1794  endif
1795 
1796  enddo
1797 
1798  fission:do i=1,nfiss
1799  fr_tmp = fissrate(i)
1800  rat = fr_tmp%cached
1801  ! neutron-induced fission
1802  if(fr_tmp%reac_type.eq.rrt_nf) then
1803  nuc = fr_tmp%fissparts(2)
1804  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1805  nfiss_sum = nfiss_sum + rat*y(nuc)*y(ineu)
1806  end if
1807  end if
1808 
1809  ! spontaneous fission
1810  if(fr_tmp%reac_type.eq.rrt_sf) then
1811  nuc = fr_tmp%fissparts(1)
1812  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1813  sfiss_sum = sfiss_sum + rat*y(nuc)
1814  end if
1815  end if
1816 
1817  ! beta-delayed fission
1818  if(fr_tmp%reac_type.eq.rrt_bf) then
1819  nuc = fr_tmp%fissparts(1)
1820  if(rat.gt.1.d-20 .and. y(nuc).gt.1.d-20) then
1821  bfiss_sum = bfiss_sum + rat*y(nuc)
1822  end if
1823  end if
1824 
1825  end do fission
1826 
1827  tau_ng = ysum / ng_sum
1828  tau_gn = ysum / gn_sum
1829  tau_pg = ysum / pg_sum
1830  tau_gp = ysum / gp_sum
1831  tau_ag = ysum / ag_sum
1832  tau_ga = ysum / ga_sum
1833  tau_np = ysum / np_sum
1834  tau_pn = ysum / pn_sum
1835  tau_an = ysum / an_sum
1836  tau_na = ysum / na_sum
1837  tau_pa = ysum / pa_sum
1838  tau_ap = ysum / ap_sum
1839  tau_beta = ysum / beta_sum
1840  tau_alpha = ysum / alpha_sum
1841  tau_nfiss = ysum / nfiss_sum
1842  tau_sfiss = ysum / sfiss_sum
1843  tau_bfiss = ysum / bfiss_sum
1844 
1845  if (.not. present(hdf5)) then
1846  ! write it to an ascii file
1847  write(tsfile, '(*(es23.16,3x))') &
1848  ctime, t9, tau_ga,tau_ag, tau_ng, tau_gn, tau_pg, tau_gp, tau_np, &
1849  tau_pn, tau_an, tau_na,tau_pa,tau_ap, tau_beta, tau_alpha, &
1850  tau_nfiss, tau_sfiss, tau_bfiss
1851  else
1852 #ifdef USE_HDF5
1853  ! Write it to the hdf5
1854  call extend_av_timescales(ctime,tau_ga, tau_ag, tau_ng, tau_gn, tau_pg, &
1855  tau_gp, tau_np, tau_pn, tau_an, tau_na, &
1856  tau_ap, tau_pa, tau_beta, tau_alpha, &
1857  tau_nfiss, tau_sfiss, tau_bfiss)
1858 
1859 #else
1860  continue
1861 #endif
1862  endif
1863 
1864  info_exit('calc_av_timescales')
1865 
1866 end subroutine calc_av_timescales
1867 
1868 
1869 !>
1870 !! Compute final abundances and write final output
1871 !!
1872 !! There are three files written _finabsum.dat_
1873 !! containing abundances summed over equal A,
1874 !! _finab.dat_ containing all final abundances,
1875 !! _finabelem.dat_ containing all abundances
1876 !! summed over equal Z.
1877 subroutine finab(Y)
1878  use global_class, only: isotope, net_size
1880  implicit none
1881 
1882  real(r_kind),dimension(net_size),intent(in) :: y
1883  integer,dimension(net_size) :: list
1884  real(r_kind),dimension(:),allocatable :: xa,yz
1885  integer :: finfile,finsfile,finzfile
1886  integer :: mval,zmval
1887  integer :: i,k
1888  integer :: mass,protonnum
1889 
1890  info_entry('finab')
1891 
1892  finfile= open_outfile('finab.dat')
1893  write(finfile,'(A)') "# 1:A 2:Z 3:N 4:Yi 5:Xi"
1894 
1895  call finab_sort(net_size,list)
1896 
1897  do k=1,net_size
1898  i=list(k)
1899  if (y(i).le.1.d-20) cycle
1900  write(finfile,'(3(i5,2x),2es23.14)') isotope(i)%mass, &
1901  isotope(i)%p_nr,isotope(i)%n_nr,y(i), &
1902  y(i)*isotope(i)%mass
1903  end do
1904  call close_io_file(finfile,'finab.dat')
1905 
1906  mval = maxval(isotope%mass)
1907  allocate(xa(mval))
1908  zmval = maxval(isotope%p_nr)
1909  allocate(yz(zmval))
1910 
1911  xa = 0.d0
1912  yz = 0.d0
1913 
1914  do i=1,net_size
1915  mass=isotope(i)%mass
1916  protonnum = isotope(i)%p_nr
1917  if(y(i).gt.1.d-20) then
1918  xa(mass) = xa(mass) + y(i)*mass
1919  if (protonnum.gt.0) yz(protonnum) = yz(protonnum) + y(i) ! omit neutrons
1920  end if
1921  end do
1922 
1923  finsfile= open_outfile('finabsum.dat')
1924  write(finsfile,'(A)') "# 1:A 2:Y(A) 3:X(A)"
1925  do i=1,mval
1926  !i:Mass, XA(i)/i: Abundance, XA(i): Mass fraction
1927  write(finsfile,'(i5,2x,es23.14,2x,es23.14)')i, xa(i)/i ,xa(i)
1928  end do
1929  finzfile= open_outfile('finabelem.dat')
1930  write(finzfile,'(A)') "# 1:Z 2:Y(Z)"
1931  do i=1,zmval
1932  write(finzfile,'(i5,2x,es23.14)')i, yz(i)
1933  end do
1934 
1935  call close_io_file(finsfile,'finabsum.dat')
1936 
1937  info_exit('finab')
1938 end subroutine finab
1939 
1940 
1941 !>
1942 !! auxiliary sorting subroutine
1943 !!
1944 subroutine finab_sort(nsize,list)
1945  use global_class, only:isotope
1947 
1948  implicit none
1949 
1950  type :: sorttype
1951  integer :: a
1952  integer :: z
1953  logical :: chk
1954  end type sorttype
1955 
1956  integer,intent(in) :: nsize
1957  integer,dimension(nsize),intent(out) :: list
1958  type(sorttype),dimension(nsize) :: tmp
1959  integer :: i,j
1960  integer :: amin,zmin
1961  integer :: nxt
1962 
1963  info_entry('finab_sort')
1964 
1965  amin = 10000
1966  zmin = 10000
1967 
1968  do i=1,nsize
1969  tmp(i)%a = isotope(i)%mass
1970  tmp(i)%z = isotope(i)%p_nr
1971  tmp(i)%chk = .false.
1972  end do
1973 
1974  do i=1,nsize
1975  do j=1,nsize
1976  if(tmp(j)%chk) cycle
1977  select case (tmp(j)%a-amin)
1978  case(:-1)
1979  nxt = j
1980  amin = tmp(j)%a
1981  zmin = tmp(j)%z
1982  case(0)
1983  if (tmp(j)%z.lt.zmin) then
1984  nxt = j
1985  amin = tmp(j)%a
1986  zmin = tmp(j)%z
1987  end if
1988  case default
1989  cycle
1990  end select
1991  end do
1992  list(i) = nxt
1993  tmp(nxt)%chk = .true.
1994  amin=10000
1995  zmin=10000
1996  end do
1997 
1998  info_exit('finab_sort')
1999 
2000  return
2001 
2002 end subroutine finab_sort
2003 
2004 
2005 !>
2006 !! Close files at the end of the calculation
2007 !!
2008 !! @author: Moritz Reichert
2009 !!
2010 !! \b Edited:
2011 !! - 01.03.17
2012 !! - 25.01.21
2013 !! .
2014 subroutine analysis_finalize()
2017  implicit none
2018 
2019  info_entry('analysis_finalize')
2020 
2021  ! Close the timescale file
2022  if (timescales_every .gt. 0) then
2023  call close_io_file(tsfile,'timescales.dat')
2024  end if
2025 
2026  ! Close the Newton-Raphson diagnostic file
2027  if (nrdiag_every.gt.0) then
2028  call close_io_file(nrdiag_unit,"nrdiag.dat")
2029  end if
2030 
2031  ! Close track nuclei file
2032  if (track_nuclei_every .gt. 0) then
2033  call close_io_file(track_unit,"tracked_nuclei.dat")
2034  end if
2035 
2036  ! Cleanup the toplist
2037  if (top_engen_every .gt. 0) then
2038  call close_io_file(toplist_unit,"toplist.dat")
2039  end if
2040 
2041  info_exit('analysis_finalize')
2042 
2043 end subroutine analysis_finalize
2044 
2045 
2046 end module analysis
parameter_class::h_track_nuclei_every
integer h_track_nuclei_every
frequency of track nuclei output in hdf5 format
Definition: parameter_class.f90:74
single_zone_vars::y
real(r_kind), dimension(:), allocatable y
Definition: single_zone_vars.f90:21
parameter_class::engen_every
integer engen_every
Energy generation output frequency.
Definition: parameter_class.f90:68
nuclear_heating::qdot_bet
real(r_kind), save, public qdot_bet
Total energy radiated away by.
Definition: nuclear_heating.f90:21
analysis::output_mainout
subroutine, private output_mainout(cnt, ctime, T9, rho_b, entropy, Rkm, Y)
Output Mainout file.
Definition: analysis.f90:626
hdf5_module
Module that contains necessary subroutines to write abundances to an hdf5 file.
Definition: hdf5_module.f90:24
hdf5_module::write_engen
subroutine, private write_engen(include_decay)
Write energy generations and time into the hdf5 file.
Definition: hdf5_module.f90:722
parameter_class::out_every
integer out_every
– runtime parameters set when reading parameter file
Definition: parameter_class.f90:61
analysis::tsfile
integer, private tsfile
timescale file unit
Definition: analysis.f90:25
parameter_class::track_nuclei_every
integer track_nuclei_every
frequency of track nuclei output
Definition: parameter_class.f90:73
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
analysis::snapshot_amount
integer snapshot_amount
Amount of custom snapshots.
Definition: analysis.f90:31
ls_timmes_eos_module
Definition: ls_timmes_eos_module.f:2
global_class::rrate
type(reactionrate_type), dimension(:), allocatable, public rrate
array containing all reaction rates used in the network
Definition: global_class.f90:65
analysis::cl_cmax
integer cl_cmax
system clock variables
Definition: analysis.f90:24
analysis::sn_ave
real(r_kind) sn_ave
average n-separation energy
Definition: analysis.f90:18
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
hdf5_module::extend_snaps
subroutine, public extend_snaps(time, Y_array)
Write abundances and time into the hdf5 file.
Definition: hdf5_module.f90:1866
analysis::cl_start
integer cl_start
Definition: analysis.f90:24
parameter_class::h_timescales_every
integer h_timescales_every
timescales output frequency in hdf5 format
Definition: parameter_class.f90:67
single_zone_vars::rhob
real(r_kind) rhob
Definition: single_zone_vars.f90:17
hdf5_module::extend_engen
subroutine, public extend_engen(time, engen_tot, heat, engen_ng_gn, engen_pg_gp, engen_ag_ga, engen_np_pn, engen_an_na, engen_ap_pa, engen_beta, engen_fiss, det_bet_engen, det_ag_engen, det_pg_engen, det_ng_engen, det_pn_engen, det_ap_engen, det_an_engen, det_other_engen, det_fiss_engen, det_tot_engen, include_decay)
Subroutine to write the energy generation per reaction type.
Definition: hdf5_module.f90:665
analysis::printsnap
subroutine, private printsnap(t, T9, rho_b, Y)
Print current snapshot.
Definition: analysis.f90:859
analysis::snapshot_time
real, dimension(:), allocatable snapshot_time
point in time for custom snapshot
Definition: analysis.f90:30
parameter_class::flow_every
integer flow_every
flow output frequency
Definition: parameter_class.f90:64
single_zone_vars::t9
real(r_kind) t9
Definition: single_zone_vars.f90:15
rrt_an
#define rrt_an
Definition: macros.h:72
nuclear_heating
Takes care of the temperature and entropy updates due to the nuclear energy release.
Definition: nuclear_heating.f90:11
analysis::output_iteration
subroutine, public output_iteration(cnt, ctime, T9, rho_b, entropy, Rkm, Y, pf)
Periodic output of analysis data for current iteration.
Definition: analysis.f90:309
analysis::output_final_step
subroutine, public output_final_step(cnt, ctime, T9, rho_b, entropy, Rkm, Y, pf)
Output final step.
Definition: analysis.f90:258
analysis::nrdiag_unit
integer nrdiag_unit
diagnostic output inside the Newton-Rhapson loop
Definition: analysis.f90:20
analysis::output_track_nuclei
subroutine, private output_track_nuclei(ctime, Y)
Output the abundances of specific nuclei.
Definition: analysis.f90:548
parameter_class::h_flow_every
integer h_flow_every
flow output frequency in hdf5 format
Definition: parameter_class.f90:65
error_msg_class
Error handling routines.
Definition: error_msg_class.f90:16
global_class::ihe4
integer, public ihe4
Definition: global_class.f90:94
analysis::toplist_amount
integer, parameter toplist_amount
Amount of reactions in the toplist.
Definition: analysis.f90:32
analysis::output_nu_loss
subroutine, private output_nu_loss(ctime, temp, rho, Rkm)
Output the energy that enters and leaves the system.
Definition: analysis.f90:588
parameter_class::termination_criterion
integer termination_criterion
condition to terminate the simulation ([0]=trajectory_file, 1=final_time, 2=final_temp,...
Definition: parameter_class.f90:127
rrt_np
#define rrt_np
Definition: macros.h:73
benam_class::reaction_string
character(50) function, public reaction_string(reac)
Return a string to represent a given reaction.
Definition: benam_class.f90:119
analysis::rosnp
real(r_kind) rosnp
for snapshot printing triggers
Definition: analysis.f90:27
analysis::finab_sort
subroutine, private finab_sort(nsize, list)
auxiliary sorting subroutine
Definition: analysis.f90:1945
analysis::finab
subroutine, public finab(Y)
Compute final abundances and write final output.
Definition: analysis.f90:1878
fission_rate_module
Module to deal with fission reactions.
Definition: fission_rate_module.f90:16
rrt_gp
#define rrt_gp
Definition: macros.h:70
analysis::output_final_stats
subroutine, public output_final_stats(cnt, ctime)
Print final statistics.
Definition: analysis.f90:820
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
single_zone_vars::evolution_mode
integer evolution_mode
NSE, network hot/cold, etc.
Definition: single_zone_vars.f90:23
flow_module::flows
type(flow_vector), dimension(:), allocatable, public flows
Definition: flow_module.f90:25
parameter_class::timescales_every
integer timescales_every
timescales output frequency
Definition: parameter_class.f90:66
thermal_neutrino_module::thermal_neutrinos
subroutine, public thermal_neutrinos(abar, Ye, temp, den, timestep, neutrino_loss)
Calculates the neutrino emissivity.
Definition: thermal_neutrino_module.f90:109
flow_module::flowprint
subroutine, public flowprint(t, t9, dens, abu, cnt)
Output reaction flows to a file.
Definition: flow_module.f90:206
parameter_class::nsep_energies_file
character(max_fname_len) nsep_energies_file
neutron separation energies
Definition: parameter_class.f90:173
analysis::t9snp
real(r_kind) t9snp
Definition: analysis.f90:27
hdf5_module::init_hdf5_module
subroutine, public init_hdf5_module()
Initialize the hdf5 module.
Definition: hdf5_module.f90:293
parameter_class::top_engen_every
integer top_engen_every
frequency of energy generators toplist
Definition: parameter_class.f90:75
fission_rate_module::nfiss
integer, public nfiss
Amount of fission rates.
Definition: fission_rate_module.f90:68
nuclear_heating::qdot_nu
real(r_kind), save, public qdot_nu
Total energy added to the system by neutrino.
Definition: nuclear_heating.f90:22
analysis::calc_nuclear_heating
subroutine, private calc_nuclear_heating(hdf5_mode, write_engen, write_toplist)
: Calculate the generated energy for each class of reactions separately.
Definition: analysis.f90:959
global_class::ineu
integer, public ineu
Definition: global_class.f90:94
rrt_pn
#define rrt_pn
Definition: macros.h:74
global_class::track_nuclei_indices
integer, dimension(:), allocatable, public track_nuclei_indices
Variables related to tracking individual nuclei.
Definition: global_class.f90:37
global_class::isotope
type(isotope_type), dimension(:), allocatable, public isotope
all nuclides used in the network
Definition: global_class.f90:34
rrt_ga
#define rrt_ga
Definition: macros.h:68
nucstuff_class
Module with helper functions such as rate tabulation, partition functions, and theoretical weak rates...
Definition: nucstuff_class.f90:11
nuclear_heating::qdot_th
real(r_kind), save, public qdot_th
Total energy lost by thermal neutrinos.
Definition: nuclear_heating.f90:23
rrt_alpd
#define rrt_alpd
Definition: macros.h:62
rrt_betp
#define rrt_betp
Definition: macros.h:60
rrt_na
#define rrt_na
Definition: macros.h:71
parameter_class::fissflag
integer fissflag
defines type of fission fragment distribution used
Definition: parameter_class.f90:86
analysis::mainout_unit
integer mainout_unit
file descriptor for the main analysis output file
Definition: analysis.f90:19
rrt_pg
#define rrt_pg
Definition: macros.h:69
hdf5_module::extend_flow
subroutine, public extend_flow(time, temp, dens, flow, n_flows, Y)
Extend the flow in the file With every call a new group is created and the flows are written in this ...
Definition: hdf5_module.f90:762
parameter_class::h_nu_loss_every
integer h_nu_loss_every
Output neutrino loss and gain in hdf5 format.
Definition: parameter_class.f90:101
rrt_ap
#define rrt_ap
Definition: macros.h:76
parameter_class::h_engen_detailed
logical h_engen_detailed
Output the energy per parent nucleus and reaction type.
Definition: parameter_class.f90:122
analysis::toplist_unit
integer toplist_unit
file id for top nuclear energy contributors
Definition: analysis.f90:22
flow_module::flowcalc
subroutine, public flowcalc(Y)
Flow calculation from jacobian. It is calculated with the help of the Jacobian.
Definition: flow_module.f90:110
thermal_neutrino_module
The thermal_neutrino_module serves as interface to the neutrino emission routines from the sneut5....
Definition: thermal_neutrino_module.f90:19
analysis::cl_end
integer cl_end
Definition: analysis.f90:24
hdf5_module::extend_nu_loss_gain
subroutine, public extend_nu_loss_gain(time, temp, dens, rad, the, heat, bet)
Write neutrino gain and loss.
Definition: hdf5_module.f90:415
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
fission_rate_module::fissrate
type(fissionrate_type), dimension(:), allocatable, public fissrate
Array storing fission reactions.
Definition: fission_rate_module.f90:39
nuclear_heating::output_nuclear_heating
subroutine, public output_nuclear_heating(cnt, ctime)
Output data related to nuclear heating.
Definition: nuclear_heating.f90:757
analysis::analysis_init
subroutine, public analysis_init()
Open files, write headers, allocate space etc.
Definition: analysis.f90:52
file_handling_class
Provide some basic file-handling routines.
Definition: file_handling_class.f90:18
parameter_class::heating_mode
integer heating_mode
Mode for heating: 0 - no heating, 1 - heating using an entropy equation, 2 - heating from the energy ...
Definition: parameter_class.f90:148
parameter_class::nrdiag_every
integer nrdiag_every
frequency of NR loop diagnostic messages (wrt iteration counter)
Definition: parameter_class.f90:70
analysis
Contains various routines for analysis and diagnostic output.
Definition: analysis.f90:9
analysis::tsnp
real(r_kind) tsnp
Definition: analysis.f90:27
single_zone_vars::time
real(r_kind) time
Definition: single_zone_vars.f90:13
global_class::track_nuclei_nr
integer, public track_nuclei_nr
amount of tracked nuclei
Definition: global_class.f90:38
single_zone_vars::stepsize
real(r_kind) stepsize
Stepsize.
Definition: single_zone_vars.f90:14
r_kind
#define r_kind
Definition: macros.h:46
analysis::calc_nseparation_energy
subroutine, private calc_nseparation_energy(Y)
Calculates average neutron separation energy.
Definition: analysis.f90:899
parameter_class::h_mainout_every
integer h_mainout_every
HDF5 output frequency of the mainout.
Definition: parameter_class.f90:72
rrt_gn
#define rrt_gn
Definition: macros.h:66
parameter_class::h_custom_snapshots
logical h_custom_snapshots
Same, but in hdf5 format.
Definition: parameter_class.f90:108
chempot
subroutine chempot(temp, den, ye, etaele, etapos)
Given a temperature temp [K], density den [g/cm**3], and a composition characterized by y_e,...
Definition: chem_pot.f90:25
analysis::track_unit
integer track_unit
file id for tracked nuclei
Definition: analysis.f90:21
fission_rate_module::fissionrate_type
fission rate type, designed to save fragment distribution at initialization
Definition: fission_rate_module.f90:21
flow_module
Provides subroutines to calculate reaction flows.
Definition: flow_module.f90:18
rrt_pa
#define rrt_pa
Definition: macros.h:75
analysis::output_flow
subroutine output_flow(ctime, T9, rho_b, Y)
Output flow files.
Definition: analysis.f90:716
parameter_class::use_thermal_nu_loss
logical use_thermal_nu_loss
Whether to include thermal neutrino loss or not.
Definition: parameter_class.f90:99
analysis::output_snapshot
subroutine, private output_snapshot(ctime, T9, rho_b, Y)
Output Snapshot files.
Definition: analysis.f90:696
analysis::analysis_finalize
subroutine, public analysis_finalize()
Close files at the end of the calculation.
Definition: analysis.f90:2015
rrt_betm
#define rrt_betm
Definition: macros.h:59
single_zone_vars::y_p
real(r_kind), dimension(:), allocatable y_p
Abundances.
Definition: single_zone_vars.f90:21
rrt_o
#define rrt_o
Definition: macros.h:81
parameter_class::nu_loss_every
integer nu_loss_every
Output neutrino loss and gain.
Definition: parameter_class.f90:100
hdf5_module::extend_track_nuclei
subroutine, public extend_track_nuclei(time, Y_array)
Write tracked abundances and time into the hdf5 file.
Definition: hdf5_module.f90:1807
analysis::calc_av_timescales
subroutine, private calc_av_timescales(ctime, T9, Y, hdf5)
Calculate average timescales of different classes of reactions.
Definition: analysis.f90:1606
file_handling_class::open_outfile
integer function, public open_outfile(file_name)
Shorthand for opening a new file for writing (output file)
Definition: file_handling_class.f90:108
nucstuff_class::inter_partf
subroutine, public inter_partf(temp, interpol)
Calculates partition function.
Definition: nucstuff_class.f90:132
hdf5_module::extend_av_timescales
subroutine, public extend_av_timescales(time, tau_ga, tau_ag, tau_ng, tau_gn, tau_pg, tau_gp, tau_np, tau_pn, tau_an, tau_na, tau_ap, tau_pa, tau_beta, tau_alpha, tau_nfiss, tau_sfiss, tau_bfiss)
Subroutine to write the average timescales (1/s) into the hdf5.
Definition: hdf5_module.f90:1753
global_class
Contains types and objects shared between multiple modules.
Definition: global_class.f90:19
parameter_class::h_snapshot_every
integer h_snapshot_every
snapshot output frequency in hdf5 format
Definition: parameter_class.f90:63
single_zone_vars::ye
real(r_kind) ye
Definition: single_zone_vars.f90:19
global_class::nreac
integer, public nreac
total number of reactions
Definition: global_class.f90:98
mergesort_module::quicksort
recursive subroutine, public quicksort(mode, length, arr1, arr2, is_recursion)
Quicksort of array.
Definition: mergesort_module.f90:810
analysis::output_nr_diagnostic
subroutine, public output_nr_diagnostic(cnt, k, nr_c, ctime, T9, stepsize, mtot, Y_p, Y)
Output of the Newton-Raphson loop diagnostics.
Definition: analysis.f90:754
rrt_ag
#define rrt_ag
Definition: macros.h:67
hdf5_module::extend_mainout
subroutine, public extend_mainout(cnt, time, temp, dens, entr, rad, Y)
Write the mainout data to the hdf5.
Definition: hdf5_module.f90:1917
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
parameter_class::mainout_every
integer mainout_every
frequency of mainout output
Definition: parameter_class.f90:71
rrt_sf
#define rrt_sf
Definition: macros.h:80
analysis::sn
real(r_kind), dimension(:), allocatable sn
array of n-separation energies
Definition: analysis.f90:17
parameter_class::custom_snapshots
logical custom_snapshots
If true, a file must be provided with numbers in days. Snapshots will be created for these points in ...
Definition: parameter_class.f90:107
rrt_ng
#define rrt_ng
Definition: macros.h:65
analysis::cl_rate
integer cl_rate
Definition: analysis.f90:24
analysis::nu_loss_gain_unit
integer nu_loss_gain_unit
file id for neutrino loss/gain
Definition: analysis.f90:23
mergesort_module
Module mergesort_module for merging arrays of rates.
Definition: mergesort_module.f90:16
single_zone_vars
Simulation variables for a single zone model.
Definition: single_zone_vars.f90:11
parameter_class::h_engen_every
integer h_engen_every
Energy generation output frequency in hdf5 format.
Definition: parameter_class.f90:69
analysis::output_initial_step
subroutine, public output_initial_step(ctime, T9, rho_b, entropy, Rkm, Y, pf)
Output first step.
Definition: analysis.f90:233
analysis::flowcnt
integer flowcnt
flow printing counter
Definition: analysis.f90:29
parameter_class::snapshot_every
integer snapshot_every
snapshot output frequency
Definition: parameter_class.f90:62
analysis::snapcnt
integer snapcnt
snapshot count
Definition: analysis.f90:28
parameter_class::calc_nsep_energy
logical calc_nsep_energy
calculate neutron separation energy?
Definition: parameter_class.f90:121
parameter_class
Contains all runtime parameters as well as phys and math constants.
Definition: parameter_class.f90:24
analysis::nuc_heat_file
integer, private nuc_heat_file
nuclear heating file unit
Definition: analysis.f90:26
benam_class
Subroutines needed to initialise the network.
Definition: benam_class.f90:11