pyTOPSScrape.parse package

Submodules

pyTOPSScrape.parse.abundance module

Author: Thomas M. Boudreaux

Created: May 2021

Last Modified: September 2022

Module responsible for the parsing and handeling of chemical composition files in the form of

#STD [Fe/H] [alpha/Fe] [C/Fe] [N/Fe] [O/Fe] [r/Fe] [s/Fe] C/O X Y,Z
F -1.13 0.32 -0.43 -0.28 0.31 -1.13 -1.13 0.10 0.7584 0.2400,1.599E-03
#H He Li Be B C N O F Ne
12.00 10.898 -0.08 0.25 1.57 6.87 6.42 7.87 3.43 7.12
#Na Mg Al Si P S Cl Ar K Ca
5.11 6.86 5.21 6.65 4.28 6.31 -1.13 5.59 3.90 5.21
#Sc Ti V Cr Mn Fe Co Ni Cu Zn
2.02 3.82 2.80 4.51 4.30 6.37 3.86 5.09 3.06 2.30
#Ga Ge As Se Br Kr Rb Sr Y Zr
0.78 1.39 0.04 1.08 0.28 0.99 0.26 0.61 1.08 1.45
#Nb Mo Tc Ru Rh Pd Ag Cd In Sn
-0.80 -0.38 -99.00 -0.51 -1.35 -0.69 -1.32 -0.55 -1.46 -0.22
#Sb Te I Xe Cs Ba La Ce Pr Nd
-1.25 -0.08 -0.71 -0.02 -1.18 1.05 -0.03 0.45 -1.54 0.29
#Pm Sm Eu Gd Tb Dy Ho Er Tm Yb
-99.00 -1.30 -0.61 -1.19 -1.96 -1.16 -1.78 -1.34 -2.16 -1.42
#Lu Hf Ta W Re Os Ir Pt Au Hg
-2.16 -1.41 -2.38 -1.41 -2.00 -0.86 -0.88 -0.64 -1.34 -1.09
#Tl Pb Bi Po At Rn Fr Ra Ac Th
-1.36 -0.51 -1.61 -99.00 -99.00 -99.00 -99.00 -99.00 -99.00 -2.20
#Pa U
-99.00 -2.80

Where each number is a(i) for the ith element and lines starting with # are comments.

pyTOPSScrape.parse.abundance.a_to_mfrac(a, amass, X)

Convert \(a(i)\) for the \(i^{th}\) element to a mass fraction using the expression

\[a(i) = \log(1.008) + \log(F_{i}) - \left[\log(X) + \log(m_{i})\right] + 12\]

Or, equivilenetly, to go from \(a(i)\) to mass fraction

\[F_{i} = \left[\frac{X m_{i}}{1.008}\right]\times 10^{a(i)-12}\]

Where \(F_{i}\) is the mass fraction of the \(i^{th}\) element, \(X\) is the Hydrogen mass fraction, and \(m_{i}\) is the ith element mass in hydrogen masses.

Parameters
  • a (float) – \(a(i)\) for the \(i^{th}\) element. For example for He chem might be 10.93. For Hydrogen it would definititionally be 12.

  • amass (float) – Mass of \(i^{th}\) element given in atomic mass units.

  • X (float) – Hydrogen mass fraction

Returns

mf – Mass fraction of \(i^{th}\) element.

Return type

float

pyTOPSScrape.parse.abundance.get_atomic_masses()

Return a dict of atomic masses from Hydrogen all the way to plutonium

Returns

amasses – Dicionary of atomic masses in atomic mass units indexed by elemental symbol.

Return type

dict of floats

pyTOPSScrape.parse.abundance.get_base_composition(aTablePath: str) Tuple[list, float, float, float]

For some abundance path return the “base” composition, this is mainly to be used for headers.

Parameters

aTablePath (str) – Path to the abundance table in the form as described in the parseChemFile module documentation

Returns

  • list – list of the composition in the form [(‘Element’,massFrac,numberFrac),…]

  • float – Hydrogen mass fraction

  • float – Helium mass fraction

  • float – Metal mass fraction

pyTOPSScrape.parse.abundance.mfrac_to_a(mfrac, amass, X, Y)

Convert mass fracition of a given element to a for that element at a given hydrogen mass fraction using the equation

\[a(i) = \log(1.008) + \log(F_{i}) - \left[\log(X) + \log(m_{i}) \right]\]

Where \(F_{i}\) is the mass fraction for the \(i^{th}\) element and \(m_{i}\) is the mass fraction for the \(i^{th}\) element.

Parameters
  • mfrac (float) – Mass fraction of the ith element.

  • amass (float) – Mass of the ith element in atomic mass units.

  • X (float) – Hydrogen mass fraction

  • Y (float) – Helium mass fraction, will be used as reference if X = 0

Returns

a – a for the ith element

Return type

float

pyTOPSScrape.parse.abundance.open_and_parse(path)

Open and parse the contents of a chemical composition file

Parameters

path (str) – Path to open file

Returns

parsed

Dictionary with two indexes.

  • Abundance Ratio

    Includes the indexes:

    • STD (str)

    • [Fe/H] (float)

    • [alpha/Fe] (float)

    • [C/Fe] (float)

    • [N/Fe] (float)

    • [O/Fe] (float)

    • [r/Fe] (float)

    • [s/Fe] (float)

    • C/O (float)

    • X (float)

    • Y (float)

    • Z (float)

  • RelativeAbundance

    Includes an index for each chemical symbol given in the file format definition provided in the module documentation. These are all floats.

Return type

dict

pyTOPSScrape.parse.abundance.open_chm_file(path)

Open a chemical composition file (format defined in the module documentation). Split the contents by line then remove all lines which start with #. Finally split each line by both whitespace and commas.

Parameters

path (str) – Path to file to open

Returns

contents – List of list of strings. The outter index selects the row, the inner index selectes the column within the row.

Return type

list

pyTOPSScrape.parse.abundance.parse(contents: list) dict

Parse chem file in the format described in the module documentation.

The abuundance ratios and abundances on the first row are added to a dict under the key [‘AbundanceRatio’] and sub indexed by the comments above each entry (Note that these are not read; rather, they are assumed to be the same in every file). The subsequent values (on all other rows) are added to the same dict under the key [‘RelativeAbundance’] and sub indexed by their chemical symbols.

Parameters

contents (list) – List of list of strings. The outter index selects the row, the inner index selected the column in the row and at each coordinate is a string which can be cast as a float. The one exception is that string at 0,0 is a charectar.

Returns

extracted

Dictionary with two indexes.

  • Abundance Ratio

    Includes the indexes:

    • STD (str)

    • [Fe/H] (float)

    • [alpha/Fe] (float)

    • [C/Fe] (float)

    • [N/Fe] (float)

    • [O/Fe] (float)

    • [r/Fe] (float)

    • [s/Fe] (float)

    • C/O (float)

    • X (float)

    • Y (float)

    • Z (float)

  • RelativeAbundance

    Includes an index for each chemical symbol given in the file format from the module documentation. These are all floats.

Return type

dict

pyTOPSScrape.parse.abundance.parse_abundance_map(path: str) numpy.ndarray

Parse Hydrogen, Helium, and metal mass fraction out of a csv where each row is one composition, the first column is X, second is Y, and the third is Z. Comments may be included in the file if the first non white space character on the line is a hash.

Parameters

path (str) – Path to the abundance map. This should be an ascii file where each row contains X, Y, and Z (comma delimited with no white space). Each row will define one set of tables to be queried with the idea being that the entire file describes the entire set of tables to be queried.

Returns

pContents – numpy array of all the compositions of length n where n is the number of rows whos first non white space character was not a hash. For a DSEP n=126. Along the second axis the first column is X, the second is Y, and the third is Z.

Return type

np.ndarray(shape=(n,3))

pyTOPSScrape.parse.opal module

pyTOPSScrape.parse.opal.load_opal(path: str) Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]

Load both the opacity table as well as the LogR and LogT for each table. LogT and LogR are not actually taken from the table; rather, they are assumed to be the correct LogT and LogR and are simply constructed from the general LogR and LogT constuctors already present in pysep.

Parameters

path (str) – path to opacity table to load.

Returns

  • LogT (np.ndarray(dtype=’float64’)) – Numpy array of shape (70,) with each required LogT value in it.

  • LogR (np.ndarray(dtype=’float64’)) – Numpy array of shape (19,) wich each required LogR value in it.

  • p (np.ndarray(dtype=’float64’)) – array of shape (126, 70, 19) where the first axis is the composition axis for the 126 compositions which dsep expects, the second is the temperature axis, and the third is the R axis.

pyTOPSScrape.parse.opal.parse_OPAL_opacity_table(path: str) numpy.ndarray

Parse the 126 tables out of a properly formated opacity table which dsep can understand. This idetifies all lines starting with Table # after the summary section and uses those to index where the tables begin. Given that dsep opacity tables are not square and that numpy can only handel rectangular data all rows are padded to the length of the longest row with np.nan. Therefore, nan(s) should be interprited as locations where the opacity table was undefined.

Parameters

path (str) – path to opacity table

Returns

p – array of shape (126, 70, 19) where the first axis is the composition axis for the 126 compositions which dsep expects, the second is the temperature axis, and the third is the R axis.

Return type

np.ndarray

pyTOPSScrape.parse.tops module

Author: Thomas M. Boudreaux

Created: September 2021

Last Modified: September 2022

Functions for parsing raw output from TOPS webform

Examples

Say you have some raw output from the TOPS webform saved to a file called “rawOutputTest.dat”, then you may parse that with

>>> rho, LogT, RMO = load_tops("rawOutputTest.dat")
pyTOPSScrape.parse.tops.extract_composition_path(path: str) Tuple[float, float, float]

Given the name of a TOPS return file (named in the format OP:n_X_Y_Z.dat) extract X, Y, and Z

Parameters

path (string) – path to TOPS return file

Returns

  • X (float) – Hydrogen mass fraction

  • Y (float) – Helium mass fraction

  • Z (Metal mass fraction)

pyTOPSScrape.parse.tops.load_tops(TOPSTable: str, n: int = 100) Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]

Given the path to a file queried from the TOPS webform put it into a computer usable form of 3 arrays. One array of mass density, one of LogT and one of log Rossland Mean Opacity

Parameters
  • TOPSTable (string) – Path to file queried from TOPS webform

  • n (int, default=100) – The size of density grid used in TOPS query form.

Returns

  • rho (np.ndarray(shape=n)) – Array of mass densities (in cgs) parsed from TOPS table.

  • LogT (np.ndarray(shape=m)) – Array of temperatures (in Kelvin) parsed from TOPS table.

  • OPALTableInit (np.ndarray(shape=(m,n))) – Array of Rossland Mean Opacities parsed from TOPS table.

Module contents