![]() |
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... | |
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
|
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 \)).
Definition at line 1361 of file readini.f90.
|
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:
Created: 17.12.20
Definition at line 882 of file readini.f90.
|
private |
Function to convert from a given density unit to g/ccm.
b will be \( 10^{-6} \) as it is the conversion from g/m^3 to g/ccm.
[in] | unit | Input unit |
Definition at line 658 of file readini.f90.
|
private |
Function to convert from a given distance unit to km.
b will be \( 1 \) as it is the conversion from km to km.
[in] | unit | Input unit |
Definition at line 698 of file readini.f90.
|
private |
Function to convert from a given neutrino energy unit to neutrino temperatures in MeV.
b will be \( 3.15^{-1} \) as it is the conversion from neutrino energies to neutrino temperatures (assuming average energies).
[in] | unit | Input unit |
Definition at line 773 of file readini.f90.
|
private |
Function to convert from a given neutrino luminosity unit to erg/s.
b will be \( 1 \) as it is the conversion from erg/s to erg/s.
[in] | unit | Input unit |
Definition at line 808 of file readini.f90.
|
private |
Function to convert from a given neutrino temperature unit to MeV.
b will be \( 1 \) as it is the conversion from Mev to Mev.
Definition at line 737 of file readini.f90.
|
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.
Definition at line 89 of file readini.f90.
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".
Edited:
[out] | iniab | initial abundances |
Definition at line 326 of file readini.f90.
|
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
Edited:
Definition at line 187 of file readini.f90.
subroutine, public readini::readini_finalize |
Readini cleanup subroutine.
At the moment, the subroutine does nothing.
Edited:
Definition at line 1433 of file readini.f90.
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.
Edited: 28.02.17
Definition at line 45 of file readini.f90.
|
private |
Function to convert from a given temperature unit to GK.
b will be \( 10^{-9} \) as it is the conversion from K to GK.
[in] | unit | Input unit |
Definition at line 616 of file readini.f90.
|
private |
Function to convert from a given time unit to seconds.
b will be \( 10^{-3} \) as it is the conversion from ms to s.
[in] | unit | Input unit |
Definition at line 576 of file readini.f90.