from_Haser

sbpy.activity.from_Haser(coma, mol_data, aper=<Quantity 25. m>)[source]

Calculate production rate for GasComa

Parameters:
comasbpy.activity.gas.GasComa

Gas coma model for ratio calculation of production rate, the production rate Q that the gas coma model expects should be an educated first guess. A good way to get this guess would be to use the function from_drahus found under sbpy.activity.gas.LTE The molecule name used for the parent argument of the coma model should be the same name or equivalent JPLSpec identifier used to calculate the total number of molecules.

mol_data: `sbpy.data.phys`

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

Total Number of Molecules (See
total_number for a calculation
of this datum if you don’t wish to provide it yourself)

This field can be given by the user directly or calculated using the necessary combinations of the following functions: from_jplspec, einstein_coeff, beta_factor, and total_number. 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.

aperQuantity

Telescope aperture in meters. Default is 25 m

Returns:
QQuantity

production rate

References

Haser 1957, Bulletin de la Societe Royale des Sciences de Liege 43, 740. Newburn and Johnson 1978, Icarus 35, 360-368.

Examples

>>> import astropy.units as u
>>> from astropy.time import Time
>>> from sbpy.data import Ephem, Phys
>>> from sbpy.activity import (Haser, LTE, photo_timescale, einstein_coeff,
...                            from_Haser)
>>> from sbpy.activity import (intensity_conversion, beta_factor,
...                            total_number)
>>> aper = 10 * u.m
>>> mol_tag = 28001
>>> temp_estimate = 25. * u.K
>>> target = 'C/2016 R2'
>>> b = 0.74
>>> vgas = 0.5 * u.km / u.s
>>> transition_freq = (230.53799 * u.GHz).to('MHz')
>>> integrated_flux = 0.26 * u.K * u.km / u.s
>>> time = Time('2017-12-22 05:24:20', format = 'iso')
>>> ephemobj = Ephem.from_horizons(target,
...                                epochs=time) 
>>> 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') 
>>> beta = beta_factor(mol_data, ephemobj) 
>>> mol_data.apply([beta.value] * beta.unit,
...                name='beta') 
>>> lte = LTE()
>>> cdensity = lte.cdensity_Bockelee(integrated_flux,
...                                  mol_data) 
>>> mol_data.apply([cdensity.value] * cdensity.unit,
...                name='cdensity') 
>>> tnum = total_number(mol_data, aper, b) 
>>> mol_data.apply([tnum.value] * tnum.unit,
...                name='total_number') 
>>> Q_estimate = 2.8*10**(28) / u.s
>>> parent = photo_timescale('CO') * vgas
>>> coma = Haser(Q_estimate, vgas, parent)
>>> Q = from_Haser(coma, mol_data, aper=aper) 
>>> Q 
    <Quantity [9.35795579e+27] 1 / s>