IonCollection#
- class fiasco.IonCollection(*args)[source]#
Bases:
objectContainer for holding multiple
Ioninstances.This container is most useful when needing to group many ions together in order to perform some aggregate calculation like a radiative loss curve or a composite spectrum made up of ions from many different elements.
- Parameters:
*ions (
fiasco.Ionorfiasco.IonCollection) – Entries can be either ions or collections of ion.
Attributes Summary
Methods Summary
bound_bound_radiative_loss(density, **kwargs)Calculate the radiative loss rate from bound-bound emission (line emission) integrated over wavelength.
free_bound(wavelength, **kwargs)Compute combined free-bound continuum emission.
Calculate the radiative loss rate from free-bound emission (collisional recombination) integrated over wavelength.
free_free(wavelength)Compute combined free-free continuum emission (bremsstrahlung).
free_free_radiative_loss([use_itoh])Calculate the radiative loss rate from free-free emission (bremsstrahlung) integrated over wavelength.
radiative_loss(density[, use_itoh])Calculate the total wavelength-integrated radiative loss rate including the bound-bound, free-bound, and free-free emission contributions
spectrum(density, emission_measure[, ...])Calculate spectrum for multiple ions
two_photon(wavelength, electron_density, ...)Compute the two-photon continuum emission.
Attributes Documentation
- temperature#
Methods Documentation
- bound_bound_radiative_loss(density, **kwargs)[source]#
Calculate the radiative loss rate from bound-bound emission (line emission) integrated over wavelength.
- free_bound(wavelength: Unit('Angstrom'), **kwargs)[source]#
Compute combined free-bound continuum emission.
Note
Both abundance and ionization fraction are included here.
The combined free-bound continuum is given by,
\[P_{fb}(\lambda,T) = \sum_{X,k}\mathrm{Ab}(X)f(X_{k+1})C_{fb, X_k}(\lambda,T)\]where \(\mathrm{Ab}(X)\) is the abundance of element \(X\), \(f(X_{k+1})\) is the ionization equilibrium of the recombining ion \(X_{k+1}\), and \(C_{fb, X_k}(\lambda,T)\) is the free-bound emission of the recombined ion \(X_k\) as computed by
fiasco.Ion.free_bound. The sum is taken over all ions in the collection.- Parameters:
wavelength (
Quantity)
See also
- free_bound_radiative_loss()[source]#
Calculate the radiative loss rate from free-bound emission (collisional recombination) integrated over wavelength.
- Returns:
rad_loss (
Quantity) – The bolometric free-bound radiative loss rate per unit emission measure
- free_free(wavelength: Unit('Angstrom'))[source]#
Compute combined free-free continuum emission (bremsstrahlung).
Note
Both abundance and ionization fraction are included here
The combined free-free continuum is given by,
\[P_{ff}(\lambda,T) = \sum_{X,k}\mathrm{Ab}(X)f(X_{k})C_{ff, X_k}(\lambda,T)\]where \(\mathrm{Ab}(X)\) is the abundance of element \(X\), \(f(X_{k})\) is the ionization equilibrium of the ion, and \(C_{ff, X_k}(\lambda,T)\) is the free-free emission of the ion as computed by
fiasco.Ion.free_free. The sum is taken over all ions in the collection.- Parameters:
wavelength (
Quantity)
See also
- free_free_radiative_loss(use_itoh=False)[source]#
Calculate the radiative loss rate from free-free emission (bremsstrahlung) integrated over wavelength.
- radiative_loss(density: Unit('1 / cm3'), use_itoh=False, **kwargs)[source]#
Calculate the total wavelength-integrated radiative loss rate including the bound-bound, free-bound, and free-free emission contributions
Note
The calculation does not include two-photon continuum emission, which is also neglected in the CHIANTI IDL routines.
- spectrum(density: Unit('1 / cm3'), emission_measure: Unit('1 / cm5'), wavelength_range=None, bin_width=None, kernel=None, **kwargs)[source]#
Calculate spectrum for multiple ions
Warning
This function is still experimental and may be removed or significantly refactored in future releases.
- Parameters:
density (
Quantity) – Electron number densityemission_measure (
Quantity) – Column emission measurewavelength_range (
Quantity, optional) – Tuple of bounds on which transitions to include. Default includes allbin_width (
Quantity, optional) – Wavelength resolution to bin intensity values. Default to 1/100 of rangekernel (
Model1DKernel, optional) – Convolution kernel for computing spectrum. Default is gaussian kernel with thermal width
- Returns:
See also
fiasco.Ion.spectrumCompute spectrum for a single ion
- two_photon(wavelength: Unit('Angstrom'), electron_density: Unit('1 / cm3'), **kwargs)[source]#
Compute the two-photon continuum emission.
Note
Both abundance and ionization equilibrium are included here.
The combined two-photon continuum is given by
\[P_{2p}(\lambda,T,n_{e}) = \sum_{X,k}\mathrm{Ab}(X)f(X_{k})C_{2p, X_k}(\lambda,T,n_{e})\]where \(\mathrm{Ab}(X)\) is the abundance of element \(X\), \(f(X_{k})\) is the ionization equilibrium of the emitting ion \(X_{k}\), and \(C_{fb, X_k}(\lambda,T)\) is the two-photon emission of the ion \(X_k\) as computed by
fiasco.Ion.two_photon. The sum is taken over all ions in the collection.See also