fiasco currently only supports version 8 of the CHIANTI database.

Source code for fiasco.io.sources.continuum_sources

"""
Source classes for continuum files
"""
import astropy.units as u
import h5py
import numpy as np
import pathlib
import plasmapy

from fiasco.io.generic import GenericParser

__all__ = ['GffguParser', 'GffintParser', 'KlgfbParser', 'VernerParser', 'ItohParser',
           'HSeqParser', 'HeSeqParser']


[docs] class GffguParser(GenericParser): """ Free-free Gaunt factor as a function of scaled frequency and energy """ filetype = 'gffgu' dtypes = [float, float, float] units = [u.dimensionless_unscaled, u.dimensionless_unscaled, u.dimensionless_unscaled] headings = ['u', 'gamma_squared', 'gaunt_factor'] descriptions = ['scaled frequency', 'scaled temperature', 'free-free Gaunt factor'] def __init__(self, filename, **kwargs): super().__init__(filename, **kwargs) self.full_path = pathlib.Path(kwargs.get('full_path', self.ascii_dbase_root / 'continuum' / filename)) self.body_index = 5
[docs] def preprocessor(self, table, line, index): if index >= self.body_index and '--' not in line: super().preprocessor(table, line, index)
[docs] def to_hdf5(self, hf, df, **kwargs): grp_name = '/'.join(['continuum', self.filetype]) if grp_name not in hf: grp = hf.create_group(grp_name) grp.attrs['chianti_version'] = df.meta['chianti_version'] grp.attrs['footer'] = df.meta['footer'] else: grp = hf[grp_name] for name in df.colnames: col = df[name] if type(col) == u.Quantity: data = col.value else: data = col.data if '<U' in data.dtype.str: numchar = data.dtype.str[2:] data = data.astype(f'|S{numchar}') if name in grp: ds = grp[name] else: if data.dtype == np.dtype('O'): ragged_dtype = h5py.special_dtype(vlen=np.dtype('float64')) ds = grp.create_dataset(name, data=data, dtype=ragged_dtype) else: ds = grp.create_dataset(name, data=data, dtype=data.dtype) if col.unit is None: ds.attrs['unit'] = 'SKIP' else: ds.attrs['unit'] = col.unit.to_string() ds.attrs['description'] = df.meta['descriptions'][name]
[docs] class GffintParser(GffguParser): """ Total free-free Gaunt factor as a function of scaled temperature. """ filetype = 'gffint' dtypes = [float, float, float, float, float] units = [u.dimensionless_unscaled, u.dimensionless_unscaled, u.dimensionless_unscaled, u.dimensionless_unscaled, u.dimensionless_unscaled] headings = ['log_gamma_squared', 'gaunt_factor', 's1', 's2', 's3'] descriptions = ['log scaled temperature', 'total free-free Gaunt factor', 'spline coefficient', 'spline coefficient', 'spline coefficient'] def __init__(self, filename, **kwargs): super().__init__(filename, **kwargs) self.body_index = 4
[docs] class KlgfbParser(GenericParser): """ Free-bound gaunt factor as a function of photon energy for several different energy levels. """ filetype = 'klgfb' dtypes = [int, int, float, float] units = [None, None, u.dimensionless_unscaled, u.dimensionless_unscaled] headings = ['n', 'l', 'log_pe', 'log_gaunt_factor'] descriptions = ['principal quantum number', 'orbital angular momentum number', 'log photon energy divided by ionization potential', 'log free-bound Gaunt factor'] def __init__(self, filename, **kwargs): super().__init__(filename, **kwargs) self.full_path = pathlib.Path(kwargs.get('full_path', self.ascii_dbase_root / 'continuum' / filename))
[docs] def preprocessor(self, table, line, index): if index == 0: pass elif index == 1: self._photon_energy = np.array(line.strip().split(), dtype=float) else: line = line.strip().split() gf = np.array(line[2:], dtype=float) table.append(line[:2] + [self._photon_energy, gf])
[docs] def to_hdf5(self, hf, df, **kwargs): grp_name = '/'.join(['continuum', self.filetype]) if grp_name not in hf: grp = hf.create_group(grp_name) grp.attrs['chianti_version'] = df.meta['chianti_version'] grp.attrs['footer'] = df.meta['footer'] else: grp = hf[grp_name] for name in df.colnames: col = df[name] if type(col) == u.Quantity: data = col.value else: data = col.data if '<U' in data.dtype.str: numchar = data.dtype.str[2:] data = data.astype(f'|S{numchar}') if name in grp: ds = grp[name] else: if data.dtype == np.dtype('O'): ragged_dtype = h5py.special_dtype(vlen=np.dtype('float64')) ds = grp.create_dataset(name, data=data, dtype=ragged_dtype) else: ds = grp.create_dataset(name, data=data, dtype=data.dtype) if col.unit is None: ds.attrs['unit'] = 'SKIP' else: ds.attrs['unit'] = col.unit.to_string() ds.attrs['description'] = df.meta['descriptions'][name]
[docs] class VernerParser(GenericParser): """ Fit parameters for calculating partial photoionization cross-sections using the method of :cite:t:`verner_analytic_1995`. """ filetype = 'verner_short' dtypes = [int, int, int, int, float, float, float, float, float, float] units = [None, None, None, None, u.eV, u.eV, u.megabarn, u.dimensionless_unscaled, u.dimensionless_unscaled, u.dimensionless_unscaled] headings = ['Z', 'n_electrons', 'n', 'l', 'E_thresh', 'E_0_fit', 'sigma_0', 'y_a_fit', 'P_fit', 'y_w_fit'] descriptions = ['atomic number', 'number of electrons', 'principal quantum number', 'orbital angular momentum number', 'threshold energy below which cross-section is 0', 'E_0 fit parameter', 'nominal value of cross-section', 'y_a fit parameter', 'P fit parameter', 'y_w fit parameter'] def __init__(self, filename, **kwargs): super().__init__(filename, **kwargs) self.full_path = pathlib.Path(kwargs.get('full_path', self.ascii_dbase_root / 'continuum' / filename))
[docs] def to_hdf5(self, hf, df, **kwargs): for row in df: el = plasmapy.particles.atomic_symbol(int(row['Z'])).lower() stage = row['Z'] - row['n_electrons'] + 1 grp_name = f'{el}/{el}_{stage}/continuum/{self.filetype}' if grp_name not in hf: grp = hf.create_group(grp_name) grp.attrs['chianti_version'] = df.meta['chianti_version'] grp.attrs['footer'] = df.meta['footer'] else: grp = hf[grp_name] for col in row.colnames: if col == 'Z' or col == 'n_electrons': continue ds = grp.create_dataset(col, data=row[col]) ds.attrs['description'] = df.meta['descriptions'][col] if not hasattr(row[col], 'unit'): ds.attrs['unit'] = 'SKIP' else: ds.attrs['unit'] = row[col].unit.to_string()
[docs] class ItohParser(GenericParser): """ Fit parameters for calculating relativistic free-free Gaunt factor using the method of :cite:t:`itoh_relativistic_2000`. """ filetype = 'itoh' dtypes = [int, float] units = [None, u.dimensionless_unscaled] headings = ['Z', 'a'] descriptions = ['atomic number', 'fit coefficient'] def __init__(self, filename, **kwargs): super().__init__(filename, **kwargs) self.full_path = pathlib.Path(kwargs.get('full_path', self.ascii_dbase_root / 'continuum' / filename))
[docs] def preprocessor(self, table, line, index): a_matrix = np.array(line.strip().split()).reshape((11, 11)) table.append([index+1] + [a_matrix])
[docs] def to_hdf5(self, hf, df, **kwargs): for row in df: el = plasmapy.particles.atomic_symbol(int(row['Z'])).lower() grp_name = f'{el}/continuum/{self.filetype}' if grp_name not in hf: grp = hf.create_group(grp_name) grp.attrs['chianti_version'] = df.meta['chianti_version'] grp.attrs['footer'] = df.meta['footer'] else: grp = hf[grp_name] for col in row.colnames: if col == 'Z': continue ds = grp.create_dataset(col, data=row[col]) ds.attrs['description'] = df.meta['descriptions'][col] if not hasattr(row[col], 'unit'): ds.attrs['unit'] = 'SKIP' else: ds.attrs['unit'] = row[col].unit.to_string()
[docs] class HSeqParser(GenericParser): r""" Parameters for calculating two-photon continuum for hydrogen-like ions Notes ----- * The parameter :math: `\psi_{\text{norm}}` (called :math: `A_{\text{sum}}` in CHIANTI) is a normalization factor of the integral of the spectral distribution function :math:`\psi(y)` from 0 to 1, such that :math: `\frac{1}{\psi_{\text{norm}}} \int_{0}^{1} \psi(y) dy = 2`. This normalization is only used For hydrogenic ions. """ filetype = 'hseq_2photon' dtypes = [int, float, int, float, float, float] units = [None, u.dimensionless_unscaled, None, 1/u.s, u.dimensionless_unscaled, u.dimensionless_unscaled] headings = ['Z', 'y', 'Z_0', 'A', 'psi_norm', 'psi'] descriptions = ['atomic number', 'fraction of energy carried by one of the two photons', 'nominal atomic number', 'radiative decay rate', 'normalization of the integral of psi from 0 to 1', 'spectral distribution function'] def __init__(self, filename, **kwargs): super().__init__(filename, **kwargs) self.full_path = pathlib.Path(kwargs.get('full_path', self.ascii_dbase_root / 'continuum' / filename))
[docs] def preprocessor(self, table, line, index): if index == 0: self._y0 = np.array(line.strip().split(), dtype=float) elif index == 1: self._z0 = np.array(line.strip().split(), dtype=float) else: line = line.strip().split() table.append([line[0]] + [self._y0] + [self._z0] + line[1:3] + [np.array(line[3:])])
[docs] def to_hdf5(self, hf, df, **kwargs): for row in df: el = plasmapy.particles.atomic_symbol(int(row['Z'])).lower() grp_name = f'{el}/continuum/{self.filetype}' if grp_name not in hf: grp = hf.create_group(grp_name) grp.attrs['chianti_version'] = df.meta['chianti_version'] grp.attrs['footer'] = df.meta['footer'] else: grp = hf[grp_name] for col in row.colnames: if col == 'Z': continue ds = grp.create_dataset(col, data=row[col]) ds.attrs['description'] = df.meta['descriptions'][col] if not hasattr(row[col], 'unit'): ds.attrs['unit'] = 'SKIP' else: ds.attrs['unit'] = row[col].unit.to_string()
[docs] class HeSeqParser(GenericParser): """ Parameters for calculating two-photon continuum for helium-like ions. """ filetype = 'heseq_2photon' dtypes = [int, float, float, float] units = [None, u.dimensionless_unscaled, 1/u.s, u.dimensionless_unscaled] headings = ['Z', 'y', 'A', 'psi'] descriptions = ['atomic number', 'fraction of energy carried by one of the two photons', 'radiative decay rate', 'spectral distribution function'] def __init__(self, filename, **kwargs): super().__init__(filename, **kwargs) self.full_path = pathlib.Path(kwargs.get('full_path', self.ascii_dbase_root / 'continuum' / filename))
[docs] def preprocessor(self, table, line, index): if index == 0: self._y0 = np.array(line.strip().split(), dtype=float) else: line = line.strip().split() table.append([line[0]] + [self._y0] + line[1:2] + [np.array(line[2:])])
[docs] def to_hdf5(self, hf, df, **kwargs): for row in df: el = plasmapy.particles.atomic_symbol(int(row['Z'])).lower() grp_name = f'{el}/continuum/{self.filetype}' if grp_name not in hf: grp = hf.create_group(grp_name) grp.attrs['chianti_version'] = df.meta['chianti_version'] grp.attrs['footer'] = df.meta['footer'] else: grp = hf[grp_name] for col in row.colnames: if col == 'Z': continue ds = grp.create_dataset(col, data=row[col]) ds.attrs['description'] = df.meta['descriptions'][col] if not hasattr(row[col], 'unit'): ds.attrs['unit'] = 'SKIP' else: ds.attrs['unit'] = row[col].unit.to_string()