readini Module Reference

Subroutines for initialization and parameter parsing. More...

Functions/Subroutines

subroutine, public readini_init ()
 Initialize the readini module and open files. More...
 
subroutine, private read_custom_snapshots ()
 Read in times for snapshots. More...
 
subroutine, private read_track_nuclei ()
 Reads in nuclei that will be tracked and where the abundance will be stored in a separate file. More...
 
subroutine, public read_seed (iniab)
 Reads in seed abundances from file parameter_class::seed_file. More...
 
real(r_kind) function, private time_unit_conversion (unit)
 Function to convert from a given time unit to seconds. More...
 
real(r_kind) function, private temp_unit_conversion (unit)
 Function to convert from a given temperature unit to GK. More...
 
real(r_kind) function, private dens_unit_conversion (unit)
 Function to convert from a given density unit to g/ccm. More...
 
real(r_kind) function, private dist_unit_conversion (unit)
 Function to convert from a given distance unit to km. More...
 
real(r_kind) function, private nutemp_unit_conversion (unit)
 Function to convert from a given neutrino temperature unit to MeV. More...
 
real(r_kind) function, private nuenergy_unit_conversion (unit)
 Function to convert from a given neutrino energy unit to neutrino temperatures in MeV. More...
 
real(r_kind) function, private nulumin_unit_conversion (unit)
 Function to convert from a given neutrino luminosity unit to erg/s. More...
 
subroutine, private custom_read ()
 Reads the trajectory file. More...
 
subroutine, private consistency_check ()
 Make tests to ensure the trajectory is okay. More...
 
subroutine, public readini_finalize ()
 Readini cleanup subroutine. More...
 

Detailed Description

Subroutines for initialization and parameter parsing.

The module is responsible for reading several input files such as initial seeds via read_seed_urs or the hydrodynamic trajectory via custom_read

Function/Subroutine Documentation

◆ consistency_check()

subroutine, private readini::consistency_check
private

Make tests to ensure the trajectory is okay.

Otherwise complain. Tests are, for example, the occurence of negative radii at some point or check if the electron fraction is out of bounds ( \( 0 \le y_e \le 1 \)).

Author
Moritz Reichert

Definition at line 1342 of file readini.f90.

Here is the call graph for this function:

◆ custom_read()

subroutine, private readini::custom_read
private

Reads the trajectory file.

The columns and structure of the file is determined by the user defined parameter parameter_class::trajectory_format. The default trajectory format is given by: "time temp dens rad ye" . This format will assume a file that looks like:

 #time[s], T9[GK], density[g/cm^3], R[km], Ye
 #--------------------------------------------
 0.0000e+00 1.0965e+01 8.7096e+12 4.9273e+01 0.03488
 2.5000e-05 9.9312e+00 7.4131e+12 5.0361e+01 0.03488
 ...

A more complex format, "time:ms x y z log_dens log_temp:K ye" , would assume a file like:

 #t [ms]       x[km]       y[km]      z[km]      log(rho[cgs]) log(T[K])   Ye
 #----------------------------------------------------------------------------------
  0.0000E+00  0.1051E+02  0.3774E+01 -0.4538E+01  0.1364E+02  0.1016E+02  0.1604E-01
  0.2521E-01  0.1046E+02  0.4724E+01 -0.4541E+01  0.1363E+02  0.1060E+02  0.1604E-01
 ...

Neutrino luminosities and energies are also read here. A format could look like: "time temp dens rad ye lnue lanue tnue tanue" .
Possible column names are:

  • time : The time (s) of the trajectory
  • temp : The temperature (GK) of the trajectory
  • dens : The density (g/ccm) of the trajectory
  • rad : The radius (km) of the trajectory
  • x,y,z : x, y, and z coordinates (km) of the trajectory
  • ye : The electron fraction
  • lanue : Anti-electron neutrino luminosities (erg/s)
  • lnue : Electron neutrino luminosities (erg/s)
  • tanue : Anti-electron neutrino temperatures (MeV)
  • tnue : Electron neutrino temperatures (MeV)
  • eanue : Anti-electron neutrino energies (MeV)
  • enue : Electron neutrino energies (MeV)
  • lanux : Anti-neutrino luminosities of heavy neutrinos (erg/s)
  • lnux : Neutrino luminosities of heavy neutrinos (erg/s)
  • tanux : Anti-neutrino temperatures of heavy neutrinos (MeV)
  • tnux : Neutrino temperatures of heavy neutrinos (MeV)
  • eanux : Anti-neutrino energies of heavy neutrinos (MeV)
  • enux : Neutrino energies of heavy neutrinos (MeV)
  • skip : skip a column
Note
Lines starting with "#" or blank lines are skipped
Author
Moritz Reichert

Created: 17.12.20

Definition at line 863 of file readini.f90.

Here is the call graph for this function:

◆ dens_unit_conversion()

real(r_kind) function, private readini::dens_unit_conversion ( character(len=*), intent(in)  unit)
private

Function to convert from a given density unit to g/ccm.

Example

b = dens_unit_conversion("g/m^3")

b will be \( 10^{-6} \) as it is the conversion from g/m^3 to g/ccm.

Note
New density units can be defined here.
Author
Moritz Reichert
Parameters
[in]unitInput unit
Returns
Conversion factor to obtain g/ccm

Definition at line 639 of file readini.f90.

Here is the call graph for this function:

◆ dist_unit_conversion()

real(r_kind) function, private readini::dist_unit_conversion ( character(len=*), intent(in)  unit)
private

Function to convert from a given distance unit to km.

Example

b = dist_unit_conversion("km")

b will be \( 1 \) as it is the conversion from km to km.

Note
New distance units can be defined here.
Author
Moritz Reichert
Parameters
[in]unitInput unit
Returns
Conversion factor to obtain km

Definition at line 679 of file readini.f90.

Here is the call graph for this function:

◆ nuenergy_unit_conversion()

real(r_kind) function, private readini::nuenergy_unit_conversion ( character(len=*), intent(in)  unit)
private

Function to convert from a given neutrino energy unit to neutrino temperatures in MeV.

Example

b = nuenergy_unit_conversion("mev")

b will be \( 3.15^{-1} \) as it is the conversion from neutrino energies to neutrino temperatures (assuming average energies).

Note
New (anti-)neutrino energy units can be defined here.
Author
Moritz Reichert
Parameters
[in]unitInput unit
Returns
Conversion factor to obtain temperatures in MeV

Definition at line 754 of file readini.f90.

Here is the call graph for this function:

◆ nulumin_unit_conversion()

real(r_kind) function, private readini::nulumin_unit_conversion ( character(len=*), intent(in)  unit)
private

Function to convert from a given neutrino luminosity unit to erg/s.

Example

b = nulumin_unit_conversion("erg/s")

b will be \( 1 \) as it is the conversion from erg/s to erg/s.

Note
New (anti-)neutrino luminosity units can be defined here.
Author
Moritz Reichert
Parameters
[in]unitInput unit
Returns
Conversion factor to obtain erg/s

Definition at line 789 of file readini.f90.

Here is the call graph for this function:

◆ nutemp_unit_conversion()

real(r_kind) function, private readini::nutemp_unit_conversion ( character(len=*), intent(in)  unit)
private

Function to convert from a given neutrino temperature unit to MeV.

Example

b = nutemp_unit_conversion("mev")

b will be \( 1 \) as it is the conversion from Mev to Mev.

Note
New (anti-)neutrino temperature units can be defined here.
Author
Moritz Reichert

Definition at line 718 of file readini.f90.

Here is the call graph for this function:

◆ read_custom_snapshots()

subroutine, private readini::read_custom_snapshots
private

Read in times for snapshots.

These times should be given by a separate file with the parameter parameter_class::custom_snapshots. The units in this file is days. The subroutines that are based on the custom snapshots assume a sorted time array and it is therefore sorted in the end of this subroutine. An example of the file could be:

 1.1574074074074073e-04
 0.0208333
 1  

This would output a snapshot at 10 s, 30 min, and 1 day.

See also
timestep_module::restrict_timestep, analysis::output_iteration
Author
Moritz Reichert

Definition at line 89 of file readini.f90.

Here is the call graph for this function:

◆ read_seed()

subroutine, public readini::read_seed ( real(r_kind), dimension(net_size), intent(out)  iniab)

Reads in seed abundances from file parameter_class::seed_file.

The input file can have a custom format specified by parameter_class::seed_format. For the seed_format "name X" The file would look like:

 # I'm a comment
    n       2.000000e-01
    p       1.000000e-01
  c12       2.000000e-01
  o16       3.000000e-01
 ca40       1.500000e-01
 tc98       0.500000e-01 

Note that empty lines or lines starting with "#" at the beginning of the file are skipped. Another example for seed format "A Z X skip" looks like:

    1   0    2.000000e-01  unimportant
    1   1    1.000000e-01  unimportant
   12   6    2.000000e-01  unimportant
   16   8    3.000000e-01  unimportant
   40  20    1.500000e-01  unimportant
   98  48    0.500000e-01  unimportant 

If the mass fractions do not sum up to one, this subroutine will rescale them. In order to make this routine work, also parameter_class::read_initial_composition has to be set to "yes".

See also
parameter_class::seed_file, parameter_class::seed_format
Author
Moritz Reichert
Remarks
This subroutine replaced the previous "read_seed_urs" subroutine of M. Ugliano.

Edited:

  • 01.06.22 - MR
Parameters
[out]iniabinitial abundances

Definition at line 307 of file readini.f90.

Here is the call graph for this function:

◆ read_track_nuclei()

subroutine, private readini::read_track_nuclei
private

Reads in nuclei that will be tracked and where the abundance will be stored in a separate file.

If parameter_class::track_nuclei_every is larger than 0, a file will be created, containing the information of the abundances of the nuclei. This is based on the input parameter parameter_class::track_nuclei_file that contains a path to a file with the format of the nuclei being the same as the format of the sunet file. An example could be:

  c13
 ti44
 ni56 
Note
Verbose level greater than 2 will create the file debug_track_nuclei.dat to debug this subroutine.
Author
Moritz Reichert

Edited:

  • 09.06.17

Definition at line 180 of file readini.f90.

Here is the call graph for this function:

◆ readini_finalize()

subroutine, public readini::readini_finalize

Readini cleanup subroutine.

At the moment, the subroutine does nothing.

Author
Moritz Reichert

Edited:

  • MR 28.02.2017
  • OK 18.06.2017

Definition at line 1414 of file readini.f90.

◆ readini_init()

subroutine, public readini::readini_init

Initialize the readini module and open files.

Depending on user defined parameters (e.g., parameter_class::custom_snapshots) this subroutine calls several subroutines to read input files. After reading the files, it calls consistency_check to make a couple of checks and raise an error if something is not correct.

Author
: Moritz Reichert

Edited: 28.02.17

Definition at line 45 of file readini.f90.

Here is the call graph for this function:

◆ temp_unit_conversion()

real(r_kind) function, private readini::temp_unit_conversion ( character(len=*), intent(in)  unit)
private

Function to convert from a given temperature unit to GK.

Example

b = temp_unit_conversion("K")

b will be \( 10^{-9} \) as it is the conversion from K to GK.

Note
New temperature units can be defined here.
Author
Moritz Reichert
Parameters
[in]unitInput unit
Returns
Conversion factor to obtain GK

Definition at line 597 of file readini.f90.

Here is the call graph for this function:

◆ time_unit_conversion()

real(r_kind) function, private readini::time_unit_conversion ( character(len=*), intent(in)  unit)
private

Function to convert from a given time unit to seconds.

Example

b = time_unit_conversion("ms")

b will be \( 10^{-3} \) as it is the conversion from ms to s.

Note
New time units can be defined here.
Author
Moritz Reichert
Parameters
[in]unitInput unit
Returns
Conversion factor to obtain seconds

Definition at line 557 of file readini.f90.

Here is the call graph for this function: