readini.f90
Go to the documentation of this file.
1 !> @file readini.f90
2 !!
3 !! The error file code for this file is ***W39***.
4 !!
5 !! Contains the module \ref readini
6 !!
7 
8 !> Subroutines for initialization and parameter parsing.
9 !!
10 !! The module is responsible for reading several input files
11 !! such as initial seeds via \ref read_seed_urs or the hydrodynamic trajectory
12 !! via \ref custom_read
13 !!
14 #include "macros.h"
15 module readini
17  use error_msg_class
21  implicit none
22 
23  ! Only make necessary routines public
29 
30  contains
31 
32 !>
33 !! Initialize the readini module and open files.
34 !!
35 !! Depending on user defined parameters
36 !! (e.g., \ref parameter_class::custom_snapshots)
37 !! this subroutine calls several subroutines to read
38 !! input files. After reading the files, it calls
39 !! \ref consistency_check to make a couple of checks and
40 !! raise an error if something is not correct.
41 !!
42 !! @author: Moritz Reichert
43 !!
44 !! \b Edited: 28.02.17
45 subroutine readini_init()
46  implicit none
47 
48  ! read trajectory
49  if (trim(trajectory_mode).EQ.'from_file') then
50  call custom_read()
51  end if
52 
53  ! Custom snapshots
56  end if
57 
58  ! Read the tracked nuclei
59  if ((track_nuclei_every .gt. 0) .or. (h_track_nuclei_every .gt. 0)) then
60  call read_track_nuclei()
61  end if
62 
63  ! Check for consistency
64  call consistency_check()
65 
66 
67 end subroutine readini_init
68 
69 
70 
71 !>
72 !! @brief Read in times for snapshots.
73 !!
74 !! These times should be given by a separate file
75 !! with the parameter \ref parameter_class::custom_snapshots.
76 !! The units in this file is days. The subroutines that are based on
77 !! the custom snapshots assume a sorted time array and it is therefore
78 !! sorted in the end of this subroutine. An example of the file could be:
79 !! \file{
80 !! 1.1574074074074073e-04
81 !! 0.0208333
82 !! 1 }
83 !! This would output a snapshot at 10 s, 30 min, and 1 day.
84 !!
85 !! @see timestep_module::restrict_timestep, analysis::output_iteration
86 !!
87 !! @author Moritz Reichert
88 !!
92  implicit none
93  integer :: c_snapshots !< File id
94  integer :: rstat !< Reading status
95  integer :: file_length !< Amount of lines of the file
96  integer :: i !< Loop variable
97  logical :: sorted !< For bubblesort
98  real(r_kind) :: helper !< Storage for sorting array
99 
100  ! Open the file
101  c_snapshots= open_infile(snapshot_file)
102 
103  ! Get the length of the file
104  file_length = 0
105  do
106  read(c_snapshots,*,iostat=rstat)
107  if (rstat .ne. 0)exit
108  file_length = file_length + 1
109  end do
110 
111  ! Raise error if file_length is 0
112  if (file_length .eq. 0) then
113  call raise_exception('No custom snapshots found in the file given with the parameter "snapshot_file".',&
114  "read_custom_snapshots",390027)
115  end if
116 
117  ! Allocate Array
118  allocate(snapshot_time(file_length),stat=rstat)
119  if (rstat /= 0) call raise_exception('Allocation of "snapshot_time" failed.',"read_custom_snapshots",&
120  390001)
121  ! Rewind the file to read it again
122  rewind(c_snapshots)
123 
124  ! Read the file again
125  do i=1,file_length
126  read(c_snapshots,*,iostat=rstat) snapshot_time(i)
127  if (rstat .ne. 0) call raise_exception("Unable to read custom snapshots."//new_line("A")//&
128  'Check the input file given with the parameter "snapshot_file".',&
129  "read_custom_snapshots",&
130  390003)
131  snapshot_time(i) = snapshot_time(i)*24.*60.*60. ! Convert Days -> Seconds
132  end do
133 
134  ! Store also length of the array
135  snapshot_amount = file_length
136  ! Bubblesort to ensure array is sorted
137  ! It is also the fastest way for already sorted arrays
138  sorted = .false.
139  do while (.not. sorted)
140  sorted = .true.
141  ! Sort the array. Use Bubblesort
142  do i=1,file_length-1
143  if (snapshot_time(i+1) .lt. snapshot_time(i)) then
144  ! Swap entries
145  helper = snapshot_time(i+1)
146  snapshot_time(i+1) = snapshot_time(i)
147  snapshot_time(i) = helper
148  sorted = .false.
149  end if
150  end do
151  end do
152 
153  ! Close the file again
154  call close_io_file(c_snapshots,snapshot_file)
155  return
156 
157 end subroutine read_custom_snapshots
158 
159 
160 
161 
162 
163 
164 !>
165 !! @brief Reads in nuclei that will be tracked and where the abundance will
166 !! be stored in a separate file.
167 !!
168 !! If parameter_class::track_nuclei_every is larger than 0,
169 !! a file will be created, containing the information of the abundances of the nuclei.
170 !! This is based on the input parameter \ref parameter_class::track_nuclei_file
171 !! that contains a path to a file with the format of the nuclei being the same
172 !! as the format of the sunet file. An example could be:
173 !! \file{
174 !! c13
175 !! ti44
176 !! ni56 }
177 !!
178 !! @note Verbose level greater than 2 will create the file _debug_track_nuclei.dat_
179 !! to debug this subroutine.
180 !!
181 !! @author Moritz Reichert
182 !!
183 !! \b Edited:
184 !! - 09.06.17
185 !! - 03.10.24 - MR
186 !! .
187 subroutine read_track_nuclei()
190  use benam_class, only: benam
192  implicit none
193  integer :: track_id !< File ID
194  integer :: istat !< Status of reading and allocation
195  integer :: i !< Loop variable
196  integer :: i_tmp !< temporary index of the loop
197  integer :: file_length !< Loop variable
198  integer :: ind_tmp !< Temporary index of the read nucleus
199  integer :: debug_track !< File ID of debug file
200  character*5 :: nam_tmp !< Temporary name
201 
202  info_entry("read_track_nuclei")
203 
204  ! Debugging
205  if (verbose_level .gt. 2) then
206  ! Open debug file and write a header for the subroutine variance (???)
207  debug_track = open_outfile('debug_track_nuclei.dat')
208  write(debug_track,'(A)') &
209  '# Debugging info for read_track_nuclei() from readini.f90'
210  write(debug_track,'("Opened File: ",A)') track_nuclei_file
211  end if
212 
213  ! Amount of nuclei to track
214  track_nuclei_nr = 0
215  ! File length can be larger, nuclei are not necessary contained in the list
216  file_length = 0
217 
218  ! Open the file
219  track_id= open_infile(track_nuclei_file)
220 
221  do ! determine the number tracked nuclei
222  read(track_id,'(a5)',iostat=istat) nam_tmp
223  if (istat .ne. 0)exit
224 
225  ind_tmp = benam(nam_tmp)
226  ! Check that the nucleus is also contained, otherwise don't track it
227  if (ind_tmp .ne. 0) track_nuclei_nr= track_nuclei_nr + 1
228  file_length = file_length + 1
229  end do
230 
231  ! Check if at least one was found. If not raise an error.
232  if (track_nuclei_nr .eq. 0) then
233  call raise_exception('No nuclei to track found in the track_nuclei_file. '//new_line("A")//&
234  'Check the file and the format of the file! '//new_line("A")//&
235  'Each line must have five characters and be right oriented, e.g., " ni56".',&
236  "read_track_nuclei", 390026)
237  end if
238 
239  ! Allocate the array that will contain the indices of the nuclei
240  allocate(track_nuclei_indices(track_nuclei_nr),stat=istat)
241  if (istat /= 0) call raise_exception('Allocation of "track_nuclei_indices" failed.',"read_track_nuclei",&
242  390001)
243 
244  ! start at the beginning of the file
245  rewind(track_id)
246 
247  ! Index of the array
248  i_tmp = 1
249  ! Read the nuclei and store the index
250  do i=1,file_length
251  read(track_id,'(a5)',iostat=istat)nam_tmp
252  ! Get the index of the nucleus
253  ind_tmp = benam(nam_tmp)
254 
255  if (verbose_level .gt. 2) then
256  !Track which nucleus is read
257  write(debug_track,'("Reading",A,7x,"Index: ",I10)') nam_tmp,ind_tmp
258  end if
259  ! write it to the array
260  if (ind_tmp .ne. 0) then
261  track_nuclei_indices(i_tmp) = ind_tmp
262  i_tmp = i_tmp +1
263  else
264  ! Say something if there is a nucleus that is not known
265  if (verbose_level .gt. 2) then
266  write(debug_track,'("Nucleus",2x,A,2x,"is not contained in network")') nam_tmp
267  end if
268  end if
269  end do
270 
271  ! Close everything again
272  call close_io_file(track_id,track_nuclei_file)
273 
274  if (verbose_level .gt. 2) then
275  write(debug_track,'("Tracking ",I10," nuclei")') track_nuclei_nr
276  call close_io_file(debug_track,'debug_track_nuclei.dat')
277  end if
278 
279  ! Output the amount of tracked nuclei
280  call write_data_to_std_out("Amount of tracked nuclei",int_to_str(track_nuclei_nr))
281 
282  info_exit("read_track_nuclei")
283 
284  return
285 end subroutine read_track_nuclei
286 
287 
288 
289 !>
290 !! @brief Reads in seed abundances from file \ref parameter_class::seed_file
291 !!
292 !! The input file can have a custom format specified by \ref parameter_class::seed_format.
293 !! For the seed_format "name X" The file would look like:
294 !! \file{
295 !! \# I'm a comment
296 !! n 2.000000e-01
297 !! p 1.000000e-01
298 !! c12 2.000000e-01
299 !! o16 3.000000e-01
300 !! ca40 1.500000e-01
301 !! tc98 0.500000e-01 }
302 !!
303 !! Note that empty lines or lines starting with "#" at the beginning of the file are skipped.
304 !! Another example for seed format "A Z X skip" looks like:
305 !! \file{
306 !! 1 0 2.000000e-01 unimportant
307 !! 1 1 1.000000e-01 unimportant
308 !! 12 6 2.000000e-01 unimportant
309 !! 16 8 3.000000e-01 unimportant
310 !! 40 20 1.500000e-01 unimportant
311 !! 98 48 0.500000e-01 unimportant }
312 !!
313 !! If the mass fractions do not sum up to one, this subroutine will rescale
314 !! them. In order to make this routine work, also
315 !! \ref parameter_class::read_initial_composition has to be set to "yes".
316 !!
317 !! @see parameter_class::seed_file, parameter_class::seed_format
318 !! @author Moritz Reichert
319 !!
320 !! @remark This subroutine replaced the previous "read_seed_urs" subroutine of
321 !! M. Ugliano.
322 !!
323 !! \b Edited:
324 !! - 01.06.22 - MR
325 !! .
326 subroutine read_seed(iniab)
328  use global_class, only: isotope,net_size
331  implicit none
332  ! The pass
333  real(r_kind),dimension(net_size),intent(out) :: iniab !< initial abundances
334  ! Internal variables
335  character(max_fname_len) :: current_col,col_name,col_unit !< Helper variables
336  character(max_fname_len) :: help_reader !< Helper variable
337  integer :: col_count !< Amount of columns
338  integer :: i,j !< Loop variable
339  integer :: a_i,z_i,n_i,x_i,y_i,name_i !< Column indices
340  integer :: seed_unit !< Seed file identifier
341  integer :: istat !< Status variable
342  integer :: skip_count !< Amount of lines skipped
343  logical :: skipping !< helper variable
344  integer :: l_file !< length of the file
345  real(r_kind) :: mafra_norm !< Mass fraction norm
346  integer,dimension(:),allocatable :: a,z,n !< Mass number, atomic number, neutron number
347  real(r_kind),dimension(:),allocatable :: y,x !< Abundance, mass fraction
348  character*5,dimension(:),allocatable :: name !< Name of nucleus
349  character(max_fname_len),dimension(:),allocatable :: row_read !< temporary array for one line
350 
351  info_entry("read_seed")
352 
353  ! Initialize all index variables
354  a_i = -1; z_i = -1; n_i = -1; x_i = -1; y_i = -1; name_i = -1
355 
356  ! Convert the custom string to lower case letters
357  do i = 1, max_fname_len
358  j = iachar(seed_format(i:i))
359  if (j>= iachar("A") .and. j<=iachar("Z") ) then
360  seed_format(i:i) = achar(iachar(seed_format(i:i))+32)
361  else
362  continue
363  end if
364  end do
365 
366  ! Get correct index from the columns and count them
367  current_col = ''
368  col_count = 0
369  do i = 2,max_fname_len
370  current_col = trim(adjustl(current_col))//seed_format(i-1:i-1)
371  if ((seed_format(i-1:i-1) .ne. ' ') .and. (seed_format(i:i) .eq. ' ')) then
372  ! Count the columns
373  col_count = col_count +1
374  ! Get the name of the column
375  col_name = ''
376  name_loop : do j=1,max_fname_len
377  if ((current_col(j:j) .eq. ':' ) .or. (current_col(j:j) .eq. ' ' )) exit name_loop
378  col_name = trim(adjustl(col_name))//current_col(j:j)
379  end do name_loop
380  ! Get the unit of the column
381  j=j+1
382  col_unit = ''
383  unit_loop : do j=j,max_fname_len
384  if ((current_col(j:j) .eq. ' ' )) exit unit_loop
385  col_unit = trim(adjustl(col_unit))//current_col(j:j)
386  end do unit_loop
387 
388  ! The column and everything is known here
389  ! Store the column index and units
390  select case (trim(adjustl(col_name)))
391  case('a')
392  a_i = col_count !index of the column (mass number)
393  case('z')
394  z_i = col_count !index of the column (proton number)
395  case('n')
396  n_i = col_count !index of the column (neutron number)
397  case('x')
398  x_i = col_count !index of the column (mass fraction)
399  case('y')
400  y_i = col_count !index of the column (abundance)
401  case('name')
402  name_i = col_count !index of the column (nucleus name)
403  case('skip')
404  continue
405  case default
406  call raise_exception('Problem when analyzing: "'//&
407  trim(adjustl(seed_format))//'". '//new_line('A')//&
408  'Can not read column "'//trim(adjustl(col_name))//'".'//&
409  new_line('A')//'Column type unknown.',"read_seeds",&
410  390021)
411  end select
412 
413  ! Clear the string for the current column
414  current_col = ''
415  end if
416  end do
417 
418  ! open the trajectory file
419  seed_unit= open_infile(seed_file)
420  ! Count how many lines to skip
421  l_file = 0
422  skip_count = 0
423  skipping = .true.
424  do ! determine the number of records
425  read(seed_unit,"(A)",iostat=istat)help_reader
426  if (istat .ne. 0)exit
427  help_reader = adjustl(help_reader) ! Move everything to the left
428 
429  ! Check if the line is to skip or not
430  if ((skipping .eqv. .true.) .and. (len_trim(help_reader) .eq. 0)) then
431  skip_count = skip_count+1
432  cycle
433  elseif ((skipping .eqv. .true.) .and. (help_reader(1:1) .eq. "#")) then
434  skip_count = skip_count+1
435  cycle
436  else
437  skipping = .false.
438  end if
439 
440  ! Check for empty lines in between
441  if ((skipping .eqv. .false.) .and. (len_trim(help_reader) .eq. 0)) then
442  if (verbose_level .gt. 1) write(*,*) "Warning, seed file ended unexpectedly."
443  exit
444  end if
445 
446  ! Count the amount of entries
447  l_file = l_file + 1
448  end do
449 
450  ! Allocate arrays
451  allocate(row_read(col_count),stat=istat)
452  if (istat /= 0) call raise_exception('allocation failed.',"read_seeds",390001)
453  allocate(a(l_file),z(l_file),n(l_file),y(l_file),x(l_file),name(l_file),stat=istat)
454  if (istat /= 0) call raise_exception('allocation failed.',"read_seeds",390001)
455 
456  ! Skip first rows
457  rewind(seed_unit)
458  do i=1,skip_count
459  read(seed_unit,*)
460  end do
461 
462  ! Initialize
463  a(:) = 0; z(:) = 0; n(:) = 0; y(:) = 0; x(:) = 0 ; name(:) = " "
464  ! Finally read the file
465  do i=1,l_file
466  read(seed_unit,*)row_read
467  if (a_i .ne. -1) a(i) = str_to_int(row_read(a_i))
468  if (z_i .ne. -1) z(i) = str_to_int(row_read(z_i))
469  if (n_i .ne. -1) n(i) = str_to_int(row_read(n_i))
470  if (y_i .ne. -1) y(i) = str_to_float(row_read(y_i))
471  if (x_i .ne. -1) x(i) = str_to_float(row_read(x_i))
472  if (name_i .ne. -1) name(i) = trim(adjustl(row_read(name_i)))
473  ! write(*,*)row_read(2)
474  end do
475 
476  ! Construct other entries
477  if ((a_i .eq. -1) .and. (z_i .ne. -1) .and. (n_i .ne. -1)) then
478  a(:) = z(:)+n(:)
479  elseif ((a_i .ne. -1) .and. (z_i .eq. -1) .and. (n_i .ne. -1)) then
480  z(:) = a(:)-n(:)
481  elseif ((a_i .ne. -1) .and. (z_i .ne. -1) .and. (n_i .eq. -1)) then
482  n(:) = a(:)-z(:)
483  elseif ((a_i .ne. -1) .and. (z_i .ne. -1) .and. (n_i .ne. -1)) then
484  continue
485  else
486  if (name_i .ne. -1) then
487  do j=1,l_file
488  iso_loop: do i=1,net_size
489  if(trim(adjustl(name(j))).eq.trim(adjustl(isotope(i)%name))) then
490  a(j) = isotope(i)%mass
491  z(j) = isotope(i)%p_nr
492  n(j) = isotope(i)%n_nr
493  exit iso_loop
494  end if
495  end do iso_loop
496  end do
497  else
498  call raise_exception("Not possible to read seeds. Columns missing!",&
499  "read_seeds",390022)
500  end if
501  end if
502 
503  ! Store abundances / mass fractions
504  if ((y_i .ne. -1) .and. (x_i .eq. -1)) then
505  x(:) = y(:)*a(:)
506  elseif ((y_i .eq. -1) .and. (x_i .ne. -1)) then
507  do i=1,l_file
508  ! Avoid division by zero, filter with sunet later
509  if (a(i) .eq. 0) a(i) = 1
510  y(i) = x(i)/a(i)
511  end do
512  elseif ((y_i .ne. -1) .and. (x_i .ne. -1)) then
513  continue
514  else
515  call raise_exception("Not possible to read seeds. Give either Y or X!",&
516  "read_seeds",390023)
517  end if
518  ! At this point, Y, X, A, Z, and N are known (or the code crashed)
519  ! Renormalize
520  mafra_norm = sum(x)
521 
522  ! Check if it is a reasonable norm
523  if (abs(mafra_norm-1d0) .gt. 0.1 ) then
524  call raise_exception("Norm of seed nuclei was terrible ("//num_to_str(abs(mafra_norm-1d0))//")! "//&
525  "Check your initial abundances or seed format!",&
526  "read_seeds",390025)
527  end if
528 
529  call write_data_to_std_out("Mass conservation pre norm",num_to_str(mafra_norm))
530  x(:) = x(:)/mafra_norm
531  y(:) = y(:)/mafra_norm
532  mafra_norm = sum(x)
533  call write_data_to_std_out("Mass conservation post norm",num_to_str(mafra_norm))
534 
535  ! Also ensure initial zero Y array
536  iniab(:) = 0
537  ! Find correct entries
538  do j=1, l_file
539  inner_loop: do i=1, net_size
540  if (name_i .ne. -1) then
541  ! Use the name as it is more unique
542  if(trim(adjustl(name(j))).eq.trim(adjustl(isotope(i)%name))) then
543  iniab(i) = y(j)
544  exit inner_loop
545  end if
546  else
547  ! Otherwise use the proton and mass number and hope for the best
548  if ((z(j) .eq. isotope(i)%p_nr) .and. (a(j) .eq. isotope(i)%mass)) then
549  iniab(i) = y(j)
550  exit inner_loop
551  end if
552  end if
553  end do inner_loop
554  end do
555 
556  ! Close the file again
557  call close_io_file(seed_unit,seed_file)
558  info_exit("read_seed")
559 
560 end subroutine read_seed
561 
562 
563 
564 !>
565 !! Function to convert from a given time unit to seconds.
566 !!
567 !! ### Example
568 !!~~~~~~~~~~~~~~.f90
569 !! b = time_unit_conversion("ms")
570 !!~~~~~~~~~~~~~~
571 !! b will be \f$ 10^{-3} \f$ as it is the conversion from ms to s.
572 !!
573 !! @note New time units can be defined here.
574 !!
575 !! @author Moritz Reichert
576 function time_unit_conversion(unit)
578  implicit none
579  character(len=*),intent(in) :: unit !< Input unit
580  real(r_kind) :: time_unit_conversion !< Conversion factor to obtain seconds
581 
582  select case (trim(adjustl(unit)))
583  case('')
585  case('s')
587  case('ms')
588  time_unit_conversion = 1.0d-3
589  case('h')
591  case('min')
593  case('yrs')
594  time_unit_conversion = 31536000
595  case default
596  call raise_exception('Problem when analyzing: "'//&
597  trim(adjustl(trajectory_format))//'". '//new_line('A')//'Time unit "'//&
598  trim(adjustl(unit))//'" unknown.',"time_unit_conversion",&
599  390004)
600  end select
601 end function time_unit_conversion
602 
603 
604 !>
605 !! Function to convert from a given temperature unit to GK
606 !!
607 !! ### Example
608 !!~~~~~~~~~~~~~~.f90
609 !! b = temp_unit_conversion("K")
610 !!~~~~~~~~~~~~~~
611 !! b will be \f$ 10^{-9} \f$ as it is the conversion from K to GK.
612 !!
613 !! @note New temperature units can be defined here.
614 !!
615 !! @author Moritz Reichert
616 function temp_unit_conversion(unit)
618  implicit none
619  character(len=*),intent(in) :: unit !< Input unit
620  real(r_kind) :: temp_unit_conversion !< Conversion factor to obtain GK
621 
622  select case (trim(adjustl(unit)))
623  case('')
625  case('-')
627  case('gk')
629  case('k')
630  temp_unit_conversion = 1.0d-9
631  case('mev')
632  temp_unit_conversion = 11.604
633  case('ev')
634  temp_unit_conversion = 11604.52*1d-9
635  case('kev')
636  temp_unit_conversion = 11604.52*1d-6
637  case default
638  call raise_exception('Problem when analyzing: "'//&
639  trim(adjustl(trajectory_format))//'". '//new_line('A')//'Temperature unit "'//&
640  trim(adjustl(unit))//'" unknown.',"temp_unit_conversion",&
641  390005)
642  end select
643 end function temp_unit_conversion
644 
645 
646 !>
647 !! Function to convert from a given density unit to g/ccm
648 !!
649 !! ### Example
650 !!~~~~~~~~~~~~~~.f90
651 !! b = dens_unit_conversion("g/m^3")
652 !!~~~~~~~~~~~~~~
653 !! b will be \f$ 10^{-6} \f$ as it is the conversion from g/m^3 to g/ccm.
654 !!
655 !! @note New density units can be defined here.
656 !!
657 !! @author Moritz Reichert
658 function dens_unit_conversion(unit)
660  implicit none
661  character(len=*),intent(in) :: unit !< Input unit
662  real(r_kind) :: dens_unit_conversion !< Conversion factor to obtain g/ccm
663 
664  select case (trim(adjustl(unit)))
665  case('')
667  case('-')
669  case('g/cm^3')
671  case('g/m^3')
672  dens_unit_conversion = 1d-6
673  case('kg/cm^3')
675  case('kg/m^3')
676  dens_unit_conversion = 1d-3
677  case default
678  call raise_exception('Problem when analyzing: "'//&
679  trim(adjustl(trajectory_format))//'". '//new_line('A')//'Density unit "'//&
680  trim(adjustl(unit))//'" unknown.',"dens_unit_conversion",&
681  390006)
682  end select
683 end function dens_unit_conversion
684 
685 
686 !>
687 !! Function to convert from a given distance unit to km
688 !!
689 !! ### Example
690 !!~~~~~~~~~~~~~~.f90
691 !! b = dist_unit_conversion("km")
692 !!~~~~~~~~~~~~~~
693 !! b will be \f$ 1 \f$ as it is the conversion from km to km.
694 !!
695 !! @note New distance units can be defined here.
696 !!
697 !! @author Moritz Reichert
698 function dist_unit_conversion(unit)
700  implicit none
701  character(len=*),intent(in) :: unit !< Input unit
702  real(r_kind) :: dist_unit_conversion !< Conversion factor to obtain km
703 
704  select case (trim(adjustl(unit)))
705  case('')
707  case('-')
709  case('km')
711  case('m')
712  dist_unit_conversion = 1e-3
713  case('cm')
714  dist_unit_conversion = 1e-5
715  case default
716  call raise_exception('Problem when analyzing: "'//&
717  trim(adjustl(trajectory_format))//'". '//new_line('A')//'Radius unit "'//&
718  trim(adjustl(unit))//'" unknown.',"dist_unit_conversion",&
719  390007)
720  end select
721 end function dist_unit_conversion
722 
723 
724 
725 !>
726 !! Function to convert from a given neutrino temperature unit to MeV
727 !!
728 !! ### Example
729 !!~~~~~~~~~~~~~~.f90
730 !! b = nutemp_unit_conversion("mev")
731 !!~~~~~~~~~~~~~~
732 !! b will be \f$ 1 \f$ as it is the conversion from Mev to Mev.
733 !!
734 !! @note New (anti-)neutrino temperature units can be defined here.
735 !!
736 !! @author Moritz Reichert
739  implicit none
740  character(len=*),intent(in) :: unit !> Input unit
741  real(r_kind) :: nutemp_unit_conversion !> Conversion factor to obtain MeV
742 
743  select case (trim(adjustl(unit)))
744  case('')
746  case('-')
748  case('mev')
750  case default
751  call raise_exception('Problem when analyzing: "'//&
752  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
753  'Neutrino temperature unit "'//&
754  trim(adjustl(unit))//'" unknown.',"nutemp_unit_conversion",&
755  390008)
756  end select
757 end function nutemp_unit_conversion
758 
759 
760 !>
761 !! Function to convert from a given neutrino energy unit to neutrino temperatures in MeV
762 !!
763 !! ### Example
764 !!~~~~~~~~~~~~~~.f90
765 !! b = nuenergy_unit_conversion("mev")
766 !!~~~~~~~~~~~~~~
767 !! b will be \f$ 3.15^{-1} \f$ as it is the conversion from neutrino energies to
768 !! neutrino temperatures (assuming average energies).
769 !!
770 !! @note New (anti-)neutrino energy units can be defined here.
771 !!
772 !! @author Moritz Reichert
775  implicit none
776  character(len=*),intent(in) :: unit !< Input unit
777  real(r_kind) :: nuenergy_unit_conversion !< Conversion factor to obtain temperatures in MeV
778 
779  select case (trim(adjustl(unit)))
780  case('')
781  nuenergy_unit_conversion = 1.0/3.151374374 ! Energy to temperature
782  case('-')
783  nuenergy_unit_conversion = 1.0/3.151374374
784  case('mev')
785  nuenergy_unit_conversion = 1.0/3.151374374
786  case default
787  call raise_exception('Problem when analyzing: "'//&
788  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
789  'Neutrino energy unit "'//&
790  trim(adjustl(unit))//'" unknown.',"nuenergy_unit_conversion",&
791  390009)
792  end select
793 end function nuenergy_unit_conversion
794 
795 
796 !>
797 !! Function to convert from a given neutrino luminosity unit to erg/s
798 !!
799 !! ### Example
800 !!~~~~~~~~~~~~~~.f90
801 !! b = nulumin_unit_conversion("erg/s")
802 !!~~~~~~~~~~~~~~
803 !! b will be \f$ 1 \f$ as it is the conversion from erg/s to erg/s.
804 !!
805 !! @note New (anti-)neutrino luminosity units can be defined here.
806 !!
807 !! @author Moritz Reichert
810  implicit none
811  character(len=*),intent(in) :: unit !< Input unit
812  real(r_kind) :: nulumin_unit_conversion !< Conversion factor to obtain erg/s
813 
814  select case (trim(adjustl(unit)))
815  case('')
817  case('-')
819  case('erg/s')
821  case default
822  call raise_exception('Problem when analyzing: "'//&
823  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
824  'Neutrino luminosity unit "'//&
825  trim(adjustl(unit))//'" unknown.',"nulumin_unit_conversion",&
826  390010)
827  end select
828 end function nulumin_unit_conversion
829 
830 
831 !>
832 !! Reads the trajectory file.
833 !!
834 !! The columns and structure of the file is determined by the user
835 !! defined parameter \ref parameter_class::trajectory_format.
836 !! The default trajectory format is given by: <b> "time temp dens rad ye" </b>.
837 !! This format will assume a file that looks like:
838 !! \file{
839 !! \#time[s], T9[GK], density[g/cm^3], R[km], Ye
840 !! \#--------------------------------------------
841 !! 0.0000e+00 1.0965e+01 8.7096e+12 4.9273e+01 0.03488
842 !! 2.5000e-05 9.9312e+00 7.4131e+12 5.0361e+01 0.03488
843 !! ...}
844 !! A more complex format, <b> "time:ms x y z log_dens log_temp:K ye" </b>, would assume a file like:
845 !! \file{
846 !! \#t [ms] x[km] y[km] z[km] log(rho[cgs]) log(T[K]) Ye
847 !! \#----------------------------------------------------------------------------------
848 !! 0.0000E+00 0.1051E+02 0.3774E+01 -0.4538E+01 0.1364E+02 0.1016E+02 0.1604E-01
849 !! 0.2521E-01 0.1046E+02 0.4724E+01 -0.4541E+01 0.1363E+02 0.1060E+02 0.1604E-01
850 !! ...}
851 !! Neutrino luminosities and energies are also read here. A format could look like:
852 !! <b> "time temp dens rad ye lnue lanue tnue tanue" </b>. \n
853 !! Possible column names are:
854 !! <pre style="line-height:.3">
855 !! - time : The time (s) of the trajectory
856 !! - temp : The temperature (GK) of the trajectory
857 !! - dens : The density (g/ccm) of the trajectory
858 !! - rad : The radius (km) of the trajectory
859 !! - x,y,z : x, y, and z coordinates (km) of the trajectory
860 !! - ye : The electron fraction
861 !! - lanue : Anti-electron neutrino luminosities (erg/s)
862 !! - lnue : Electron neutrino luminosities (erg/s)
863 !! - tanue : Anti-electron neutrino temperatures (MeV)
864 !! - tnue : Electron neutrino temperatures (MeV)
865 !! - eanue : Anti-electron neutrino energies (MeV)
866 !! - enue : Electron neutrino energies (MeV)
867 !! - lanux : Anti-neutrino luminosities of heavy neutrinos (erg/s)
868 !! - lnux : Neutrino luminosities of heavy neutrinos (erg/s)
869 !! - tanux : Anti-neutrino temperatures of heavy neutrinos (MeV)
870 !! - tnux : Neutrino temperatures of heavy neutrinos (MeV)
871 !! - eanux : Anti-neutrino energies of heavy neutrinos (MeV)
872 !! - enux : Neutrino energies of heavy neutrinos (MeV)
873 !! - skip : skip a column
874 !! .
875 !! </pre>
876 !!
877 !! @note Lines starting with "#" or blank lines are skipped
878 !!
879 !! @author Moritz Reichert
880 !!
881 !! \b Created: 17.12.20
882 subroutine custom_read()
885  use nuflux_class
887  use hydro_trajectory
888  implicit none
889  !
890  integer :: traj_unit, i,j, istat
891  integer :: col_count
892  integer :: time_i,temp_i,dens_i,rad_i,ye_i !< index of the columns
893  integer :: x_i,y_i,z_i !< index for x, y, and z if r not given
894  real(r_kind) :: conv_time,conv_temp,conv_dens !< unit conversion factors
895  real(r_kind) :: conv_rad,conv_ye !< unit conversion factors
896  real(r_kind) :: conv_x,conv_y,conv_z !< unit conversion factors
897  real(r_kind) :: x_tmp,y_tmp,z_tmp !< Temporary storage of x,y, and z to convert them
898  integer :: nuetemp_i,anuetemp_i !< index of the neutrino columns
899  integer :: lnue_i,lanue_i !< index of the neutrino columns
900  real(r_kind) :: conv_tnue,conv_tanue !< unit conversion factors
901  real(r_kind) :: conv_lnue,conv_lanue !< unit conversion factors
902  integer :: nuxtemp_i,anuxtemp_i !< index of the neutrino columns
903  integer :: lnux_i,lanux_i !< index of the neutrino columns
904  real(r_kind) :: conv_tnux,conv_tanux !< unit conversion factors
905  real(r_kind) :: conv_lnux,conv_lanux !< unit conversion factors
906  character(max_fname_len) :: current_col,col_name,col_unit
907  real(r_kind),dimension(:),allocatable :: row_read !< temporary array for one line
908  real(r_kind) :: tny !< tiny value, preventing zeros.
909  logical :: log_switch !< flag to indicate if one columns is logarithmic
910  logical :: ls_time,ls_temp,ls_dens,ls_rad !< which column is it logarithmic?
911  logical :: ls_x,ls_y,ls_z,ls_ye !< which column is it logarithmic?
912  logical :: ls_tnue,ls_tanue,ls_lnue,ls_lanue!< which column is it logarithmic?
913  logical :: ls_tnux,ls_tanux,ls_lnux,ls_lanux!< which column is it logarithmic?
914  integer :: skip_count !< Check how many lines to skip
915  logical :: skipping !< Helper variable
916  character(max_fname_len) :: help_reader !< Helper variable
917 
918  info_entry("custom_read")
919 
920  tny = tiny(tny)
921 
922  ! open the trajectory file
923  traj_unit= open_infile(trajectory_file)
924 
925  ! Convert the custom string to lower case letters
926  do i = 1, max_fname_len
927  j = iachar(trajectory_format(i:i))
928  if (j>= iachar("A") .and. j<=iachar("Z") ) then
929  trajectory_format(i:i) = achar(iachar(trajectory_format(i:i))+32)
930  else
931  continue
932  end if
933  end do
934 
935 
936  ! Count how many lines to skip and how many are there
937  zsteps=0
938  skip_count = 0
939  skipping = .true.
940  do ! determine the number of records
941  read(traj_unit,"(A)",iostat=istat)help_reader
942  ! Go out, file ended
943  if (istat .ne. 0)exit
944  ! Move string to the left
945  help_reader = adjustl(help_reader)
946  ! Check if the line is to skip or not
947  if ((skipping .eqv. .true.) .and. (len_trim(help_reader) .eq. 0)) then
948  skip_count = skip_count+1
949  cycle
950  elseif ((skipping .eqv. .true.) .and. (help_reader(1:1) .eq. "#")) then
951  skip_count = skip_count+1
952  cycle
953  else
954  skipping = .false.
955  end if
956  zsteps= zsteps + 1
957  end do
958 
959  ! allocate arrays
960  allocate(ztime(zsteps), ztemp(zsteps), zdens(zsteps), &
961  zrad(zsteps) , zye(zsteps) , stat=istat)
962  if (istat /= 0) call raise_exception('allocation failed.',"custom_read",390001)
963 
964  ! Allocate neutrino arrays
965  if (((nuflag.ge.1)) .and. (neutrino_mode .eq. 'from_file')) then
966  allocate(tnue(zsteps) ,tnuebar(zsteps) , &
968  tnux(zsteps) ,tnuxbar(zsteps) , &
970  stat=istat)
971  if (istat /= 0) call raise_exception('allocation failed.',"custom_read",390001)
972  end if
973 
974  ! Set default values
975  time_i = 0
976  temp_i = 0
977  dens_i = 0
978  rad_i = 0
979  x_i = 0
980  y_i = 0
981  z_i = 0
982  ye_i = 0
983  ! Same for neutrino quantities (electron neutrinos)
984  nuetemp_i = 0; anuetemp_i = 0
985  lnue_i = 0; lanue_i = 0
986  ! Same for the heavy neutrinos
987  nuxtemp_i = 0; anuxtemp_i = 0
988  lnux_i = 0; lanux_i = 0
989 
990  ! Get correct index from the columns and count them
991  current_col = ''
992  col_count = 0
993  do i = 2,max_fname_len
994  current_col = trim(adjustl(current_col))//trajectory_format(i-1:i-1)
995  if ((trajectory_format(i-1:i-1) .ne. ' ') .and. (trajectory_format(i:i) .eq. ' ')) then
996  ! Count the columns
997  col_count = col_count +1
998  ! Get the name of the column
999  col_name = ''
1000  name_loop : do j=1,max_fname_len
1001  if ((current_col(j:j) .eq. ':' ) .or. (current_col(j:j) .eq. ' ' )) exit name_loop
1002  col_name = trim(adjustl(col_name))//current_col(j:j)
1003  end do name_loop
1004  ! Get the unit of the column
1005  j=j+1
1006  col_unit = ''
1007  unit_loop : do j=j,max_fname_len
1008  if ((current_col(j:j) .eq. ' ' )) exit unit_loop
1009  col_unit = trim(adjustl(col_unit))//current_col(j:j)
1010  end do unit_loop
1011 
1012  ! Check if the column is already the logarithm
1013  select case (col_name(1:4))
1014  case('log_')
1015  log_switch = .true.
1016  case default
1017  log_switch = .false.
1018  end select
1019 
1020  ! Remove the "log_" from the string to identify the column
1021  if (log_switch) then
1022  col_name = col_name(5:)
1023  end if
1024 
1025  ! Store the column index and units
1026  select case (trim(adjustl(col_name)))
1027  case('time')
1028  time_i = col_count !index of the column
1029  ls_time=log_switch !check if it should be logarithmic
1030  ! Get the conversion factor to seconds
1031  conv_time = time_unit_conversion(trim(adjustl(col_unit)))
1032  case('temp')
1033  temp_i = col_count !index of the column
1034  ls_temp=log_switch !check if it should be logarithmic
1035  ! Get the conversion factor to GK
1036  conv_temp = temp_unit_conversion(trim(adjustl(col_unit)))
1037  case('dens')
1038  dens_i = col_count !index of the column
1039  ls_dens=log_switch !check if it should be logarithmic
1040  ! Get the conversion factor to g/ccm
1041  conv_dens = dens_unit_conversion(trim(adjustl(col_unit)))
1042  case('rad')
1043  rad_i = col_count !index of the column
1044  ls_rad=log_switch !check if it should be logarithmic
1045  ! Get the conversion factor to km
1046  conv_rad = dist_unit_conversion(trim(adjustl(col_unit)))
1047  case('x')
1048  x_i = col_count !index of the column
1049  ls_x=log_switch !check if it should be logarithmic
1050  ! Get the conversion factor to km
1051  conv_x = dist_unit_conversion(trim(adjustl(col_unit)))
1052  case('y')
1053  y_i = col_count !index of the column
1054  ls_y=log_switch !check if it should be logarithmic
1055  ! Get the conversion factor to km
1056  conv_y = dist_unit_conversion(trim(adjustl(col_unit)))
1057  case('z')
1058  z_i = col_count !index of the column
1059  ls_z=log_switch !check if it should be logarithmic
1060  ! Get the conversion factor to km
1061  conv_z = dist_unit_conversion(trim(adjustl(col_unit)))
1062  case('ye')
1063  ye_i = col_count
1064  ls_ye = log_switch !check if it should be logarithmic, for ye not very probable
1065  ! Ye does not have a unit
1066  conv_ye = 1.0
1067  case('tnue')
1068  nuetemp_i = col_count !index of the column
1069  ls_tnue = log_switch !check if it should be logarithmic
1070  ! Get the conversion factor to MeV
1071  conv_tnue = nutemp_unit_conversion(trim(adjustl(col_unit)))
1072  case('tanue')
1073  anuetemp_i = col_count !index of the column
1074  ls_tanue = log_switch !check if it should be logarithmic
1075  ! Get the conversion factor to MeV
1076  conv_tanue = nutemp_unit_conversion(trim(adjustl(col_unit)))
1077  case('enue')
1078  nuetemp_i = col_count
1079  ls_tnue = log_switch !check if it should be logarithmic
1080  ! Get the conversion factor to MeV and neutrino temperatures
1081  conv_tnue = nuenergy_unit_conversion(trim(adjustl(col_unit)))
1082  case('eanue')
1083  anuetemp_i = col_count
1084  ls_tanue = log_switch !check if it should be logarithmic
1085  ! Get the conversion factor to MeV and neutrino temperatures
1086  conv_tanue = nuenergy_unit_conversion(trim(adjustl(col_unit)))
1087  case('lnue')
1088  lnue_i = col_count
1089  ls_lnue= log_switch !check if it should be logarithmic
1090  ! Get the conversion factor to erg/s
1091  conv_lnue= nulumin_unit_conversion(trim(adjustl(col_unit)))
1092  case('lanue')
1093  lanue_i = col_count
1094  ls_lanue= log_switch !check if it should be logarithmic
1095  ! Get the conversion factor to erg/s
1096  conv_lanue= nulumin_unit_conversion(trim(adjustl(col_unit)))
1097  case("tnux")
1098  nuxtemp_i = col_count !index of the column
1099  ls_tnux = log_switch !check if it should be logarithmic
1100  ! Get the conversion factor to MeV
1101  conv_tnux = nutemp_unit_conversion(trim(adjustl(col_unit)))
1102  case("tanux")
1103  anuxtemp_i = col_count !index of the column
1104  ls_tanux = log_switch !check if it should be logarithmic
1105  ! Get the conversion factor to MeV
1106  conv_tanux = nutemp_unit_conversion(trim(adjustl(col_unit)))
1107  case("enux")
1108  nuxtemp_i = col_count
1109  ls_tnux = log_switch !check if it should be logarithmic
1110  ! Get the conversion factor to MeV and neutrino temperatures
1111  conv_tnux = nuenergy_unit_conversion(trim(adjustl(col_unit)))
1112  case("eanux")
1113  anuxtemp_i = col_count
1114  ls_tanux = log_switch !check if it should be logarithmic
1115  ! Get the conversion factor to MeV and neutrino temperatures
1116  conv_tanux = nuenergy_unit_conversion(trim(adjustl(col_unit)))
1117  case("lnux")
1118  lnux_i = col_count
1119  ls_lnux= log_switch !check if it should be logarithmic
1120  ! Get the conversion factor to erg/s
1121  conv_lnux= nulumin_unit_conversion(trim(adjustl(col_unit)))
1122  case("lanux")
1123  lanux_i = col_count
1124  ls_lanux= log_switch !check if it should be logarithmic
1125  ! Get the conversion factor to erg/s
1126  conv_lanux= nulumin_unit_conversion(trim(adjustl(col_unit)))
1127  case('skip')
1128  continue
1129  case default
1130  call raise_exception('Problem when analyzing: "'//&
1131  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
1132  'Can not read column "'//trim(adjustl(col_name))//'".'//&
1133  new_line('A')//'Column type unknown.',"custom_read",&
1134  390011)
1135  end select
1136  ! Clear the string for the current column
1137  current_col = ''
1138  end if
1139  end do
1140 
1141  ! Check if all necessary hydro-columns were there
1142  if ((time_i .eq. 0) .or. (temp_i .eq. 0) .or. (dens_i .eq. 0)) then
1143  call raise_exception('Problem when analyzing: "'//&
1144  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
1145  'Necessary column is missing.',"custom_read",&
1146  390012)
1147  end if
1148 
1149  ! Ye is only necessary if no seed file is given
1150  if ((ye_i .eq. 0) .and. (.not. read_initial_composition)) then
1151  call raise_exception('Problem when analyzing: "'//&
1152  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
1153  'Ye was missing.',"custom_read",&
1154  390012)
1155  end if
1156 
1157  ! Radius is necessary for neutrinos
1158  if (((rad_i .eq. 0) .and. (x_i .eq. 0)) .and. (nuflag.ge.1) ) then
1159  call raise_exception('Problem when analyzing: "'//&
1160  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
1161  'Radius was missing.',"custom_read",&
1162  390012)
1163  end if
1164 
1165  ! Check if all neutrino-columns were there
1166  if (((nuflag.ge.1)) .and. (neutrino_mode .eq. 'from_file')) then
1167  if ((nuetemp_i .eq. 0) .or. (anuetemp_i .eq. 0) .or. (lnue_i .eq. 0) .or. &
1168  (lanue_i .eq. 0)) then
1169  call raise_exception('Problem when analyzing: "'//&
1170  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
1171  'Necessary column with neutrino information is missing.',"custom_read",&
1172  390013)
1173  end if
1174  end if
1175 
1176  ! Check that the columns for the heavies neutrinos are there in case of neutral
1177  ! current neutrino reactions
1178  if (((nuflag .eq. 4) .or. (nuflag .eq. 5)) .and. (neutrino_mode .eq. 'from_file')) then
1179  if ((nuxtemp_i .eq. 0) .or. (anuxtemp_i .eq. 0) .or. (lnux_i .eq. 0) .or. &
1180  (lanux_i .eq. 0)) then
1181  call raise_exception('Problem when analyzing: "'//&
1182  trim(adjustl(trajectory_format))//'". '//new_line('A')//&
1183  'Necessary column with (heavy) neutrino information is missing.',"custom_read",&
1184  390013)
1185  end if
1186  end if
1187 
1188 
1189  ! read the file
1190  rewind(traj_unit)
1191  ! Skip the first rows again
1192  do i=1,skip_count
1193  read(traj_unit,*)
1194  end do
1195  ! allocate array to read one row
1196  allocate(row_read(col_count),stat=istat)
1197  if (istat /= 0) call raise_exception('Allocation failed.',"custom_read",390001)
1198 
1199  ! Read the trajectory and convert the units
1200  do i=1,zsteps
1201  read(traj_unit,*,iostat=istat)row_read
1202  ! Check that everything is okay
1203  if (istat .ne. 0) then
1204  call raise_exception("Could not read trajectory file."//new_line('A')//&
1205  "Line "//int_to_str(skip_count+i)//" raised a problem.",&
1206  "custom_read",390024)
1207  end if
1208 
1209  ! Time
1210  if (ls_time) ztime(i) = 10**(row_read(time_i))
1211  if (.not. ls_time) ztime(i) = row_read(time_i)
1212  ztime(i) = ztime(i)*conv_time
1213 
1214  ! Density
1215  if (ls_dens) zdens(i) = 10**(row_read(dens_i))
1216  if (.not. ls_dens) zdens(i) = row_read(dens_i)
1217  zdens(i) = zdens(i)*conv_dens
1218 
1219  ! Temperature
1220  if (ls_temp) ztemp(i) = 10**(row_read(temp_i))
1221  if (.not. ls_temp) ztemp(i) = row_read(temp_i)
1222  ztemp(i) = ztemp(i)*conv_temp
1223 
1224  ! Electron fraction
1225  if (ye_i .ne. 0) then
1226  if (ls_ye) zye(i) = 10**(row_read(ye_i))
1227  if (.not. ls_ye) zye(i) = row_read(ye_i)
1228  zye(i) = zye(i)*conv_ye
1229  else
1230  zye(i) = 0.5 ! Fill with dummy value if ye will be taken from a seed file
1231  end if
1232 
1233  ! Radius is special in the case that x, y, or z was given instead of r
1234  if (rad_i .ne. 0) then
1235  if (ls_rad) zrad(i) = 10**(row_read(rad_i))
1236  if (.not. ls_rad) zrad(i) = row_read(rad_i)
1237  zrad(i) = zrad(i)*conv_rad
1238  elseif ((rad_i .eq. 0) .and. (x_i .ne. 0) .and. (y_i .eq. 0) .and. (z_i .eq. 0)) then
1239  if (ls_x) x_tmp = 10**(row_read(x_i))
1240  if (.not. ls_x) x_tmp = row_read(x_i)
1241  zrad(i) = x_tmp*conv_x
1242  elseif ((rad_i .eq. 0) .and. (x_i .ne. 0) .and. (y_i .ne. 0) .and. (z_i .eq. 0)) then
1243  ! check if it was converted by the logarithm
1244  if (ls_x) x_tmp = 10**(row_read(x_i))
1245  if (.not. ls_x) x_tmp = row_read(x_i)
1246  if (ls_y) y_tmp = 10**(row_read(y_i))
1247  if (.not. ls_y) y_tmp = row_read(y_i)
1248  ! Calculate radius
1249  zrad(i) = sqrt((x_tmp*conv_x)**2+(y_tmp*conv_y)**2)
1250  elseif ((rad_i .eq. 0) .and. (x_i .ne. 0) .and. (y_i .ne. 0) .and. (z_i .ne. 0)) then
1251  ! check if it was converted by the logarithm
1252  if (ls_x) x_tmp = 10**(row_read(x_i))
1253  if (.not. ls_x) x_tmp = row_read(x_i)
1254  if (ls_y) y_tmp = 10**(row_read(y_i))
1255  if (.not. ls_y) y_tmp = row_read(y_i)
1256  if (ls_z) z_tmp = 10**(row_read(z_i))
1257  if (.not. ls_z) z_tmp = row_read(z_i)
1258  zrad(i) = sqrt((x_tmp*conv_x)**2+(y_tmp*conv_y)**2+(z_tmp*conv_z)**2)
1259  else
1260  ! The radius is not always necessary. E.g., if there is no extrapolation
1261  ! and no neutrinos.
1262  zrad(i) = 1
1263  ! call raise_exception('Not possible to read the radius, specify r or x!.',"custom_read",&
1264  ! 390014)
1265  end if
1266  ! Read neutrino quantities
1267  if ((nuflag.ge.1) .and. (neutrino_mode .eq. 'from_file')) then
1268 
1269  ! neutrino temperature
1270  if (ls_tnue) tnue(i) = 10**row_read(nuetemp_i)
1271  if (.not. ls_tnue) tnue(i) = row_read(nuetemp_i)
1272  tnue(i) = tnue(i)*conv_tnue
1273 
1274  ! anti-neutrino temperature
1275  if (ls_tanue) tnuebar(i) = 10**row_read(anuetemp_i)
1276  if (.not. ls_tanue) tnuebar(i) = row_read(anuetemp_i)
1277  tnuebar(i) = tnuebar(i)*conv_tanue
1278 
1279  ! neutrino luminosity
1280  if (ls_lnue) nlume(i) = 10**row_read(lnue_i)
1281  if (.not. ls_lnue) nlume(i) = row_read(lnue_i)
1282  nlume(i) = nlume(i)*conv_lnue
1283 
1284  if (ls_lanue) nlumebar(i) = 10**row_read(lanue_i)
1285  if (.not. ls_lanue) nlumebar(i) = row_read(lanue_i)
1286  nlumebar(i) = nlumebar(i)*conv_lanue
1287 
1288  ! Avoid small and negative numbers
1289  if (nlume(i).le.0.d0) nlume(i) = tny
1290  if (nlumebar(i).le.0.d0) nlumebar(i) = tny
1291  if (tnue(i).le.0.d0) tnue(i) = tny
1292  if (tnuebar(i).le.0.d0) tnuebar(i) = tny
1293  end if
1294  ! Do the same for the heavy neutrinos
1295  if ( ((nuflag .eq. 4) .or. (nuflag .eq. 5)) .and. (neutrino_mode .eq. 'from_file')) then
1296  ! Neutrino temperature of heavy neutrinos
1297  if (ls_tnux) tnux(i) = 10**row_read(nuxtemp_i)
1298  if (.not. ls_tnux) tnux(i) = row_read(nuxtemp_i)
1299  tnux(i) = tnux(i)*conv_tnux
1300 
1301  ! Neutrino temperature of heavy anti neutrinos
1302  if (ls_tanux) tnuxbar(i) = 10**row_read(anuxtemp_i)
1303  if (.not. ls_tanux) tnuxbar(i) = row_read(anuxtemp_i)
1304  tnuxbar(i) = tnuxbar(i)*conv_tanux
1305 
1306  ! Neutrino luminosity of heavy neutrinos
1307  if (ls_lnux) nlumx(i) = 10**row_read(lnux_i)
1308  if (.not. ls_lnux) nlumx(i) = row_read(lnux_i)
1309  nlumx(i) = nlumx(i)*conv_lnux
1310 
1311  ! Neutrino luminosity of heavy anti neutrinos
1312  if (ls_lanux) nlumxbar(i) = 10**row_read(lanux_i)
1313  if (.not. ls_lanux) nlumxbar(i) = row_read(lanux_i)
1314  nlumxbar(i) = nlumxbar(i)*conv_lanux
1315 
1316  ! Avoid small and negative numbers
1317  if (nlumx(i).le.0.d0) nlumx(i) = tny
1318  if (tnux(i).le.0.d0) tnux(i) = tny
1319  if (tnuxbar(i).le.0.d0) tnuxbar(i) = tny
1320  if (nlumxbar(i).le.0.d0) nlumxbar(i) = tny
1321  end if
1322  end do
1323 
1324  ! log representation
1325  ztemp = dlog(ztemp)
1326  zdens = dlog(zdens)
1327  zrad = dlog(zrad)
1328 
1329  ! Same for neutrino quantities if necessary
1330  if ((nuflag.ge.1) .and. (neutrino_mode .eq. 'from_file')) then
1331  nlume = dlog(nlume)
1332  nlumebar = dlog(nlumebar)
1333  tnue = dlog(tnue)
1334  tnuebar = dlog(tnuebar)
1335  end if
1336  ! Same for heavy neutrinos if necessary
1337  if ( ((nuflag .eq. 4) .or. (nuflag .eq. 5)) .and. &
1338  (neutrino_mode .eq. 'from_file')) then
1339  nlumx = dlog(nlumx)
1340  nlumxbar = dlog(nlumxbar)
1341  tnux = dlog(tnux)
1342  tnuxbar = dlog(tnuxbar)
1343  end if
1344 
1345  ! close file
1346  call close_io_file(traj_unit,trajectory_file)
1347 
1348  info_exit("custom_read")
1349 
1350 end subroutine custom_read
1351 
1352 
1353 !>
1354 !! Make tests to ensure the trajectory is okay.
1355 !!
1356 !! Otherwise complain. Tests are, for example, the occurence of
1357 !! negative radii at some point or check if the electron fraction
1358 !! is out of bounds (\f$ 0 \le y_e \le 1 \f$).
1359 !!
1360 !! @author Moritz Reichert
1361 subroutine consistency_check()
1362  use hydro_trajectory
1363  implicit none
1364  integer :: i !< Loop variable
1365 
1366  ! Bunch of consistency checks for "from_file" mode
1367  if (trim(trajectory_mode).EQ.'from_file') then
1368 
1369 
1370  ! Negative time
1371  if (ztime(zsteps) .lt. 0) then
1372  call raise_exception("Last point of trajectory had negative time ("//&
1373  trim(adjustl((num_to_str(ztime(zsteps)))))//&
1374  ")!"//new_line('A')//'Please check your input trajectory.',&
1375  "consistency_check",&
1376  390015)
1377  ! Electron fraction out of bounds
1378  elseif ((maxval(zye) .gt. 1) .or. (minval(zye) .lt. 0)) then
1379  call raise_exception("The electron fraction was out of bounds (0 < Ye < 1)! "//&
1380  new_line('A')//'Please check your input trajectory.',&
1381  "consistency_check",&
1382  390016)
1383  end if
1384 
1385 
1386  ! Check negative entries
1387  do i=1, zsteps
1388  ! Negative radius
1389  if (isnan(dexp(zrad(i)))) then
1390  call raise_exception('Problem with input radius: "'//&
1391  'The radius was below zero."'//new_line('A')//&
1392  'Please check your input trajectory.',"consistency_check",&
1393  390017)
1394  ! Negative temperature
1395  elseif (isnan(dexp(ztemp(i)))) then
1396  call raise_exception('Problem with input temperature: "'//&
1397  'The temperature was below zero".'//new_line('A')//&
1398  'Please check your input trajectory.',"consistency_check",&
1399  390018)
1400  ! Negative density
1401  elseif (isnan(dexp(zdens(i)))) then
1402  call raise_exception('Problem with input density: "'//&
1403  'The density was below zero".'//new_line('A')//&
1404  'Please check your input trajectory.',"consistency_check",&
1405  390019)
1406  end if
1407  end do
1408 
1409  ! Check for increasing time
1410  do i=1, zsteps-1
1411  if (ztime(i) .ge. ztime(i+1)) call raise_exception("Time in the trajectory was not monotonically increasing."//&
1412  new_line('A')//"Check your trajectory!","consistency_check",&
1413  390020)
1414  end do
1415  end if
1416 
1417 
1418 end subroutine consistency_check
1419 
1420 
1421 
1422 !>
1423 !! Readini cleanup subroutine
1424 !!
1425 !! At the moment, the subroutine does nothing.
1426 !!
1427 !! @author Moritz Reichert
1428 !!
1429 !! \b Edited:
1430 !! - MR 28.02.2017
1431 !! - OK 18.06.2017
1432 !! .
1433 subroutine readini_finalize()
1435  implicit none
1436  ! not much here really
1437 end subroutine readini_finalize
1438 
1439 
1440 
1441 end module readini
parameter_class::h_track_nuclei_every
integer h_track_nuclei_every
frequency of track nuclei output in hdf5 format
Definition: parameter_class.f90:74
nuflux_class::nlumxbar
real(r_kind), dimension(:), allocatable, public nlumxbar
(anti-)mu and tau neutrino number luminosities from trajectory
Definition: nuflux_class.f90:92
parameter_class::track_nuclei_every
integer track_nuclei_every
frequency of track nuclei output
Definition: parameter_class.f90:73
nuflux_class::tnuxbar
real(r_kind), dimension(:), allocatable, public tnuxbar
(anti-)mu and tau neutrino temperatures from trajectory
Definition: nuflux_class.f90:90
readini::consistency_check
subroutine, private consistency_check()
Make tests to ensure the trajectory is okay.
Definition: readini.f90:1362
analysis::snapshot_amount
integer snapshot_amount
Amount of custom snapshots.
Definition: analysis.f90:31
nuflux_class::tnuebar
real(r_kind), dimension(:), allocatable, public tnuebar
(anti-)electron neutrino temperatures from trajectory
Definition: nuflux_class.f90:89
readini::nuenergy_unit_conversion
real(r_kind) function, private nuenergy_unit_conversion(unit)
Function to convert from a given neutrino energy unit to neutrino temperatures in MeV.
Definition: readini.f90:774
analysis::snapshot_time
real, dimension(:), allocatable snapshot_time
point in time for custom snapshot
Definition: analysis.f90:30
nuflux_class::nlumebar
real(r_kind), dimension(:), allocatable, public nlumebar
(anti-)electron neutrino number luminosities from trajectory
Definition: nuflux_class.f90:91
readini::readini_finalize
subroutine, public readini_finalize()
Readini cleanup subroutine.
Definition: readini.f90:1434
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
readini::time_unit_conversion
real(r_kind) function, private time_unit_conversion(unit)
Function to convert from a given time unit to seconds.
Definition: readini.f90:577
readini::temp_unit_conversion
real(r_kind) function, private temp_unit_conversion(unit)
Function to convert from a given temperature unit to GK.
Definition: readini.f90:617
nuflux_class::tnue
real(r_kind), dimension(:), allocatable, public tnue
Definition: nuflux_class.f90:89
readini::dist_unit_conversion
real(r_kind) function, private dist_unit_conversion(unit)
Function to convert from a given distance unit to km.
Definition: readini.f90:699
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
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
readini::read_track_nuclei
subroutine, private read_track_nuclei()
Reads in nuclei that will be tracked and where the abundance will be stored in a separate file.
Definition: readini.f90:188
readini
Subroutines for initialization and parameter parsing.
Definition: readini.f90:15
error_msg_class::str_to_float
real(r_kind) function, public str_to_float(input_string)
Converts a string to a float.
Definition: error_msg_class.f90:179
readini::dens_unit_conversion
real(r_kind) function, private dens_unit_conversion(unit)
Function to convert from a given density unit to g/ccm.
Definition: readini.f90:659
parameter_class::neutrino_mode
character(max_fname_len) neutrino_mode
Similar to trajectory mode.
Definition: parameter_class.f90:185
parameter_class::max_fname_len
integer, parameter, public max_fname_len
maximum length of filenames
Definition: parameter_class.f90:58
parameter_class::trajectory_file
character(max_fname_len) trajectory_file
name of trajectory data file
Definition: parameter_class.f90:152
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
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
nuflux_class::nlume
real(r_kind), dimension(:), allocatable, public nlume
Definition: nuflux_class.f90:91
hydro_trajectory::ztime
real(r_kind), dimension(:), allocatable, public ztime
time information from trajectory
Definition: hydro_trajectory.f90:19
nuflux_class
Contains variables and parameters related to neutrino fluxes.
Definition: nuflux_class.f90:14
parameter_class::trajectory_format
character(max_fname_len) trajectory_format
Format to read the trajectory.
Definition: parameter_class.f90:184
hydro_trajectory::zrad
real(r_kind), dimension(:), allocatable, public zrad
radii from trajectory
Definition: hydro_trajectory.f90:23
parameter_class::read_initial_composition
logical read_initial_composition
specify whether initial distribution of abundances should be read from file
Definition: parameter_class.f90:91
readini::read_seed
subroutine, public read_seed(iniab)
Reads in seed abundances from file parameter_class::seed_file.
Definition: readini.f90:327
readini::readini_init
subroutine, public readini_init()
Initialize the readini module and open files.
Definition: readini.f90:46
parameter_class::seed_file
character(max_fname_len) seed_file
name of file with initial seeds for trajectory run
Definition: parameter_class.f90:153
global_class::net_size
integer, public net_size
total number of isotopes (network size)
Definition: global_class.f90:93
parameter_class::track_nuclei_file
character(max_fname_len) track_nuclei_file
File of nuclei to track. Gives an output similar to mainout.dat.
Definition: parameter_class.f90:181
file_handling_class
Provide some basic file-handling routines.
Definition: file_handling_class.f90:18
error_msg_class::str_to_int
integer function, public str_to_int(input_string)
Converts a string to an integer.
Definition: error_msg_class.f90:155
parameter_class::trajectory_mode
character *20 trajectory_mode
determines how trajectory is calculated (from_file|analytic|expand)
Definition: parameter_class.f90:90
error_msg_class::num_to_str
character(:) function, allocatable, public num_to_str(num)
Converts a given real to a string with format "(1pE10.2)".
Definition: error_msg_class.f90:205
analysis
Contains various routines for analysis and diagnostic output.
Definition: analysis.f90:9
readini::nulumin_unit_conversion
real(r_kind) function, private nulumin_unit_conversion(unit)
Function to convert from a given neutrino luminosity unit to erg/s.
Definition: readini.f90:809
global_class::track_nuclei_nr
integer, public track_nuclei_nr
amount of tracked nuclei
Definition: global_class.f90:38
r_kind
#define r_kind
Definition: macros.h:46
parameter_class::h_custom_snapshots
logical h_custom_snapshots
Same, but in hdf5 format.
Definition: parameter_class.f90:108
parameter_class::seed_format
character(max_fname_len) seed_format
Seed format.
Definition: parameter_class.f90:154
readini::read_custom_snapshots
subroutine, private read_custom_snapshots()
Read in times for snapshots.
Definition: readini.f90:90
hydro_trajectory::zsteps
integer zsteps
number of timesteps in the hydro trajectory
Definition: hydro_trajectory.f90:17
hydro_trajectory
Contains arrays representing thermodynamic conditions from hydro trajectory file.
Definition: hydro_trajectory.f90:14
hydro_trajectory::zdens
real(r_kind), dimension(:), allocatable, public zdens
density information from trajectory
Definition: hydro_trajectory.f90:21
readini::nutemp_unit_conversion
real(r_kind) function, private nutemp_unit_conversion(unit)
Function to convert from a given neutrino temperature unit to MeV.
Definition: readini.f90:738
parameter_class::snapshot_file
character(max_fname_len) snapshot_file
File that contains days, where a snapshot should be made.
Definition: parameter_class.f90:177
file_handling_class::open_infile
integer function, public open_infile(file_name)
Same for reading (input file)
Definition: file_handling_class.f90:126
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
hydro_trajectory::zye
real(r_kind), dimension(:), allocatable, public zye
electron fraction information from trajectory
Definition: hydro_trajectory.f90:22
global_class
Contains types and objects shared between multiple modules.
Definition: global_class.f90:19
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
hydro_trajectory::ztemp
real(r_kind), dimension(:), allocatable, public ztemp
temperature information from trajectory
Definition: hydro_trajectory.f90:20
readini::custom_read
subroutine, private custom_read()
Reads the trajectory file.
Definition: readini.f90:883
nuflux_class::nlumx
real(r_kind), dimension(:), allocatable, public nlumx
Definition: nuflux_class.f90:92
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
parameter_class::nuflag
integer nuflag
defines type of neutrino reactions used
Definition: parameter_class.f90:85
parameter_class
Contains all runtime parameters as well as phys and math constants.
Definition: parameter_class.f90:24
nuflux_class::tnux
real(r_kind), dimension(:), allocatable, public tnux
Definition: nuflux_class.f90:90
benam_class
Subroutines needed to initialise the network.
Definition: benam_class.f90:11