LTE

class sbpy.activity.LTE[source]

Bases: object

LTE Methods for calculating production_rate

Methods Summary

cdensity_Bockelee(integrated_flux, mol_data)

Basic equation relating column density with observed integrated flux without the need for an initial column density to be given.

from_Drahus(integrated_flux, mol_data, ephemobj)

Returns production rate based on Drahus 2012 model referenced.

Methods Documentation

cdensity_Bockelee(integrated_flux, mol_data)[source]

Basic equation relating column density with observed integrated flux without the need for an initial column density to be given. This is found in equation 10 in https://ui.adsabs.harvard.edu/abs/2004come.book..391B and is derived from data from JPLSpec, feel free to use your own column density to calculate production rate or use this function with your own molecular data as long as you are aware of the needed data.

Parameters:
integrated_fluxQuantity

Integrated flux of emission line.

mol_datasbpy.data.phys
sbpy.data.phys object that contains AT LEAST the following data:
Transition frequency in MHz
Einstein Coefficient (1/s)

This function will calculate the column density from Bockelee-Morvan et al. 2004 and append it to the phys object as ‘Column Density’ or any of its alternative field names. The values above can either be given by the user or obtained from the functions einstein_coeff and beta_factor Keywords that can be used for these values are found under fieldnames documentation. We recommend the use of the JPL Molecular Spectral Catalog and the use of from_jplspec to obtain these values in order to maintain consistency. Yet, if you wish to use your own molecular data, it is possible. 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.

Returns:
Column Densityastropy.units.Quantity

Column density from Bockelee-Morvan et al. 2004 as astropy Quantity (1/m^2)

from_Drahus(integrated_flux, mol_data, ephemobj, vgas=<Quantity 1. km / s>, aper=<Quantity 25. m>, b=1.2)[source]

Returns production rate based on Drahus 2012 model referenced. Does not include photodissociation, good for first guesses for more computationally intensive methods or for the Haser model under sbpy.activity.gas.productionrate.from_Haser

Parameters:
integrated_fluxQuantity

Line integral derived from spectral data in Kelvins * km/s

mol_datasbpy.data.phys
sbpy.data.phys object that contains the following data:
Transition frequency in MHz
Temperature in Kelvins
Partition function at designated temperature (unitless)
Upper state degeneracy (unitless)
Upper level energy in Joules
Degrees of freedom (unitless)
Einstein Coefficient (1/s)

These fields can be given by the user directly or calculated using from_jplspec, einstein_coeff, Keywords that can be used for these values are found under fieldnames documentation. We recommend the use of the JPL Molecular Spectral Catalog and the use of from_jplspec to obtain these values in order to maintain consistency. Yet, if you wish to use your own molecular data, it is possible. 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.

ephemobjsbpy.data.ephem

sbpy.data.ephem object holding ephemeride information including distance from comet to Sun [‘r’] and from comet to observer [‘delta’]

vgasQuantity

Gas velocity approximation in km / s. Default is 1 km / s

aperQuantity

Telescope aperture in meters. Default is 25 m

bint

Dimensionless factor intrinsic to every antenna. Typical value, and the default for this model, is 1.22. See references for more information on this parameter.

Returns:
qQuantity

Production rate, not including photodissociation

References

Drahus et al. September 2012. The Sources of HCN and CH3OH and the Rotational Temperature in Comet 103P/Hartley 2 from Time-resolved Millimeter Spectroscopy. The Astrophysical Journal, Volume 756, Issue 1.

Examples

>>> import astropy.units as u
>>> from astropy.time import Time
>>> from sbpy.data import Ephem, Phys
>>> from sbpy.activity import LTE, einstein_coeff, intensity_conversion
>>> temp_estimate = 47. * u.K
>>> target = '103P'
>>> vgas = 0.8 * u.km / u.s
>>> aper = 30 * u.m
>>> b = 1.13
>>> mol_tag = 27001
>>> transition_freq = (265.886434 * u.GHz).to('MHz')
>>> integrated_flux = 1.22 * u.K * u.km / u.s
>>> time = Time('2010-11-3 00:48:06', format='iso')
>>> ephemobj = Ephem.from_horizons(
...     target, epochs=time, closest_apparition=True,
...     id_type='designation') 
>>> mol_data = Phys.from_jplspec(temp_estimate, transition_freq,
...                              mol_tag) 
>>> intl = intensity_conversion(mol_data) 
>>> mol_data.apply([intl.value] * intl.unit,
...                name='intl') 
>>> au = einstein_coeff(mol_data) 
>>> mol_data.apply([au.value] * au.unit,
...                name='eincoeff') 
>>> lte = LTE()
>>> q = lte.from_Drahus(integrated_flux, mol_data,
...                     ephemobj, vgas, aper, b=b) 
>>> q  
<MaskedQuantity 1.09899965e+25 1 / s>