casa  $Rev:20696$
 All Classes Namespaces Files Functions Variables
opacity.py
Go to the documentation of this file.
00001 __all__ = ["model", "skydip"]
00002 import os
00003 import math
00004 from asap.scantable import scantable
00005 from asap.asapmath import merge
00006 from asap.asapfitter import fitter
00007 from asap.selector import selector
00008 from asap._asap import atmosphere
00009 from asap import rcParams
00010 
00011 class model(object):
00012     def _to_pascals(self, val):
00013         if val > 2000:
00014             return val
00015         return val*100
00016 
00017     def __init__(self, temperature=288, pressure=101325., humidity=0.5,
00018                  elevation=700.):
00019         """
00020         This class implements opacity/atmospheric brightness temperature model
00021         equivalent to the model available in MIRIAD. The actual math is a
00022         convertion of the Fortran code written by Bob Sault for MIRIAD.
00023         It implements a simple model of the atmosphere and Liebe's model (1985)
00024         of the complex refractive index of air.
00025 
00026         The model of the atmosphere is one with an exponential fall-off in
00027         the water vapour content (scale height of 1540 m) and a temperature
00028         lapse rate of 6.5 mK/m. Otherwise the atmosphere obeys the ideal gas
00029         equation and hydrostatic equilibrium.
00030 
00031         Note, the model includes atmospheric lines up to 800 GHz, but was not
00032         rigorously tested above 100 GHz and for instruments located at
00033         a significant elevation. For high-elevation sites it may be necessary to
00034         adjust scale height and lapse rate.
00035 
00036         Parameters:
00037             temperature:    air temperature at the observatory (K)
00038             pressure:       air pressure at the sea level if the observatory
00039                             elevation is set to non-zero value (note, by
00040                             default is set to 700m) or at the observatory
00041                             ground level if the elevation is set to 0. (The
00042                             value is in Pascals or hPa, default 101325 Pa
00043             humidity:       air humidity at the observatory (fractional),
00044                             default is 0.5
00045             elevation:     observatory elevation about sea level (in meters)
00046         """
00047         self._atm = atmosphere(temperature, self._to_pascals(pressure),
00048                                humidity)
00049         self.set_observatory_elevation(elevation)
00050 
00051     def get_opacities(self, freq, elevation=None):
00052         """Get the opacity value(s) for the given frequency(ies).
00053         If no elevation is given the opacities for the zenith are returned.
00054         If an elevation is specified refraction is also taken into account.
00055         Parameters:
00056             freq:       a frequency value in Hz, or a list of frequency values.
00057                         One opacity value per frequency is returned as a scalar
00058                         or list.
00059             elevation:  the elevation at which to compute the opacity. If `None`
00060                         is given (default) the zenith opacity is returned.
00061 
00062 
00063         """
00064         func = None
00065         if isinstance(freq, (list, tuple)):
00066             if elevation is None:
00067                 return self._atm.zenith_opacities(freq)
00068             else:
00069                 elevation *= math.pi/180.
00070                 return self._atm.opacities(freq, elevation)
00071         else:
00072             if elevation is None:
00073                 return self._atm.zenith_opacity(freq)
00074             else:
00075                 elevation *= math.pi/180.
00076                 return self._atm.opacity(freq, elevation)
00077 
00078     def set_weather(self, temperature, pressure, humidity):
00079         """Update the model using the given environmental parameters.
00080         Parameters:
00081             temperature:    air temperature at the observatory (K)
00082             pressure:       air pressure at the sea level if the observatory
00083                             elevation is set to non-zero value (note, by
00084                             default is set to 700m) or at the observatory
00085                             ground level if the elevation is set to 0. (The
00086                             value is in Pascals or hPa, default 101325 Pa
00087             humidity:       air humidity at the observatory (fractional),
00088                             default is 0.5
00089         """
00090         pressure = self._to_pascals(pressure)
00091         self._atm.set_weather(temperature, pressure, humidity)
00092 
00093     def set_observatory_elevation(self, elevation):
00094         """Update the model using the given the observatory elevation
00095         Parameters:
00096             elevation:  the elevation at which to compute the opacity. If `None`
00097                         is given (default) the zenith opacity is returned.
00098         """
00099         self._atm.set_observatory_elevation(elevation)
00100 
00101 
00102 def _import_data(data):
00103     if not isinstance(data, (list,tuple)):
00104         if isinstance(data, scantable):
00105             return data
00106         elif isinstance(data, str):
00107             return scantable(data)
00108     tables = []
00109     for d in data:
00110         if isinstance(d, scantable):
00111             tables.append(d)
00112         elif isinstance(d, str):
00113             if os.path.exists(d):
00114                 tables.append(scantable(d))
00115             else:
00116                 raise IOError("Data file doesn't exists")
00117         else:
00118             raise TypeError("data is not a scantable or valid file")
00119     return merge(tables)
00120 
00121 def skydip(data, averagepol=True, tsky=300., plot=False,
00122            temperature=288, pressure=101325., humidity=0.5):
00123     """Determine the opacity from a set of 'skydip' obervations.
00124     This can be any set of observations over a range of elevations,
00125     but will ususally be a dedicated (set of) scan(s).
00126     Return a list of 'n' opacities for 'n' IFs. In case of averagepol
00127     being 'False' a list of 'n*m' elements where 'm' is the number of
00128     polarisations, e.g.
00129     nIF = 3, nPol = 2 => [if0pol0, if0pol1, if1pol0, if1pol1, if2pol0, if2pol1]
00130 
00131     The opacity is determined by fitting a first order polynomial to:
00132 
00133 
00134         Tsys(airmass) = p0 + airmass*p1
00135 
00136     where
00137 
00138         airmass = 1/sin(elevation)
00139 
00140         tau =  p1/Tsky
00141 
00142     Parameters:
00143         data:       a list of file names or scantables or a single
00144                     file name or scantable.
00145         averagepol: Return the average of the opacities for the polarisations
00146                     This might be useful to set to 'False' if one polarisation
00147                     is corrupted (Mopra). If set to 'False', an opacity value
00148                     per polarisation is returned.
00149         tsky:       The sky temperature (default 300.0K). This might
00150                     be read from the data in the future.
00151         plot:       Plot each fit (airmass vs. Tsys). Default is 'False'
00152     """
00153     # quiten output
00154     verbose = rcParams["verbose"]
00155     rcParams["verbose"] = False
00156     try:
00157         if plot:
00158             from matplotlib import pylab
00159         scan = _import_data(data)
00160         f = fitter()
00161         f.set_function(poly=1)
00162         sel = selector()
00163         basesel = scan.get_selection()
00164         inos = scan.getifnos()
00165         pnos = scan.getpolnos()
00166         opacities = []
00167         om = model(temperature, pressure, humidity)
00168         for ino in inos:
00169             sel.set_ifs(ino)
00170             opacity = []
00171             fits = []
00172             airms = []
00173             tsyss = []
00174             if plot:
00175                 pylab.cla()
00176                 pylab.ioff()
00177                 pylab.clf()
00178                 pylab.xlabel("Airmass")
00179                 pylab.ylabel(r"$T_{sys}$")
00180             for pno in pnos:
00181                 sel.set_polarisations(pno)
00182                 scan.set_selection(basesel+sel)
00183                 freq = scan.get_coordinate(0).get_reference_value()/1e9
00184                 freqstr = "%0.4f GHz" % freq
00185                 tsys = scan.get_tsys()
00186                 elev = scan.get_elevation()
00187                 airmass = [ 1./math.sin(i) for i in elev ]
00188                 airms.append(airmass)
00189                 tsyss.append(tsys)
00190                 f.set_data(airmass, tsys)
00191                 f.fit()
00192                 fits.append(f.get_fit())
00193                 params = f.get_parameters()["params"]
00194                 opacity.append(params[1]/tsky)
00195             if averagepol:
00196                 opacities.append(sum(opacity)/len(opacity))
00197             else:
00198                 opacities += opacity
00199             if plot:
00200                 colors = ['b','g','k']
00201                 n = len(airms)
00202                 for i in range(n):
00203                     pylab.plot(airms[i], tsyss[i], 'o', color=colors[i])
00204                     pylab.plot(airms[i], fits[i], '-', color=colors[i])
00205                     pylab.figtext(0.7,0.3-(i/30.0),
00206                                       r"$\tau_{fit}=%0.2f$" % opacity[i],
00207                                       color=colors[i])
00208                 if averagepol:
00209                     pylab.figtext(0.7,0.3-(n/30.0),
00210                                       r"$\tau_{avg}=%0.2f$" % opacities[-1],
00211                                       color='r')
00212                     n +=1
00213                 pylab.figtext(0.7,0.3-(n/30.0),
00214                               r"$\tau_{model}=%0.2f$" % om.get_opacities(freq*1e9),
00215                               color='grey')
00216 
00217                 pylab.title("IF%d : %s" % (ino, freqstr))
00218 
00219                 pylab.ion()
00220                 pylab.draw()
00221                 raw_input("Hit <return> for next fit...")
00222             sel.reset()
00223 
00224         scan.set_selection(basesel)
00225         if plot:
00226             pylab.close()
00227         return opacities
00228     finally:
00229         rcParams["verbose"] = verbose