NonLTE

class sbpy.activity.NonLTE[source]

Bases: object

Class for non LTE production rate models.

Methods Summary

from_pyradex(integrated_flux, mol_data[, ...])

Calculate production rate from the Non-LTE iterative code pyradex Presently, only the LAMDA catalog is supported by this function.

Methods Documentation

static from_pyradex(integrated_flux, mol_data, line_width=<Quantity 1. km / s>, escapeProbGeom='lvg', iter=100, collider_density={'H2': 1980.0000000000002})[source]

Calculate production rate from the Non-LTE iterative code pyradex Presently, only the LAMDA catalog is supported by this function. In the future a function will be provided by sbpy to build your own molecular data file from JPLSpec for use in this function. Collider is assumed to be H2O for the default setting since we are only considering comet production rates. See documentation for pyradex installation tips

Parameters:
integrated_fluxQuantity

The integrated flux in K * km/s

mol_datasbpy.data.phys

sbpy.data.phys object that contains AT LEAST the following data:

mol_tag: molecule of interest as string or int JPLSpec identifier
temp: kinetic temperature in gas coma (unit K)
cdensity : cdensity estimate (can be calculated from cdensity_Bockelee) (unit 1/cm^2)
temp_back: (optional) background temperature in K (default is 2.730 K)
lamda_name: (optional) LAMDA molecule identifier to avoid ambiguity. Lamda.molecule_dict provides list

Keywords that can be used for these values are found under fieldnames documentation. Make sure to inform yourself on the values needed for each function, their units, and their interchangeable keywords as part of the Phys data class.

line_widthQuantity

The FWHM line width (really, the single-zone velocity width to scale the column density by: this is most sensibly interpreted as a velocity gradient (dv/length)) in km/s (default is 1.0 km/s)

escapeProbGeomstr

Which escape probability method to use, available choices are ‘sphere’, ‘lvg’, and ‘slab’

iterint

Number of iterations you wish to perform. Default is 100, more iterations will take more time to do but will offer a better range of results to compare against. The range of guesses is built by having the column density guess and subtracting/adding an order of magnitude for the start and end values of the loop, respectively. i.e. a guess of 1e15 will create a range between 1e14 and 1e16

collider_densitydict

Dictionary of colliders and their densities in cm^-3. Allowed entries are any of the following : h2,oh2,ph2,e,He,H,H+ See Pyradex documentation for more information. Default dictionary is {‘H2’ : 900*2.2} where 900 is the collider density of H2 and 2.2 is the value for the square root of the ratio of reduced masses of H2O/H2 as follows: (mu_H2O/mu_H2)**0.5 = ((18**2/18*2)/((18*2)/(18+2)))**0.5 = 2.2 in order to scale the collisional excitation to H2O as the main collisional partner. (Walker, et al. 2014; de val Borro, et al. 2017; & Schoier, et al. 2004)

Returns:
column densityQuantity

column density to use for the Haser model ratio calculation Note: it is normal for pyradex/RADEX to output warnings depending on the setup the user gives it (such as not having found a molecular data file so it searches for it on LAMDA catalog)

References

Haser 1957, Bulletin de la Societe Royale des Sciences de Liege 43, 740.

Walker, et al., On the Validity of Collider-mass Scaling for Molecular Rotational Excitation, APJ, August 2014.

van der Tak, et al., A computer program for fast non-LTE analysis of interstellar line spectra. With diagnostic plots to interpret observed line intensity ratios. A&A, February 12 2013.

de Val Borro, et al., Measuring molecular abundances in comet C/2014 Q2 (Lovejoy) using the APEX telescope. Monthly Notices of the Royal Astronomical Society, October 27 2017.

Schoier, et al., An atomic and molecular database for analysis of submillimetre line observations. A&A, November 4 2004.

Examples

>>> from sbpy.activity import NonLTE
>>> from sbpy.data import Phys
>>> import astropy.units as u
>>> transition_freq = (177.196 * u.GHz).to(u.MHz)
>>> mol_tag =  29002
>>> cdensity_guess = (1.89*10.**(14) / (u.cm * u.cm))
>>> temp_estimate = 20. * u.K
>>> mol_data = Phys.from_jplspec(temp_estimate, transition_freq,
...                              mol_tag)  
>>> mol_data.apply([cdensity_guess.value] * cdensity_guess.unit,
...                 name= 'cdensity')  
>>> mol_data.apply(['HCO+@xpol'],
...                 name='lamda_name')  
>>> nonLTE = NonLTE()
>>> try:  
...     cdensity = nonLTE.from_pyradex(1.234 * u.K * u.km / u.s,
...                                    mol_data, iter=600,
...                                    collider_density={'H2': 900})  
...     print(cdensity)  
... except ImportError:
...     pass
Closest Integrated Flux:[1.24925956] K km / s
Given Integrated Flux: 1.234 K km / s
[1.06363773e+14] 1 / cm2