Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1# 

2# LSST Data Management System 

3# Copyright 2008-2018 LSST/AURA. 

4# 

5# This product includes software developed by the 

6# LSST Project (http://www.lsst.org/). 

7# 

8# This program is free software: you can redistribute it and/or modify 

9# it under the terms of the GNU General Public License as published by 

10# the Free Software Foundation, either version 3 of the License, or 

11# (at your option) any later version. 

12# 

13# This program is distributed in the hope that it will be useful, 

14# but WITHOUT ANY WARRANTY; without even the implied warranty of 

15# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 

16# GNU General Public License for more details. 

17# 

18# You should have received a copy of the LSST License Statement and 

19# the GNU General Public License along with this program. If not, 

20# see <http://www.lsstcorp.org/LegalNotices/>. 

21# 

22 

23from astropy import units 

24from astropy.units import cds 

25import numpy as np 

26 

27import lsst.geom 

28from lsst.afw.coord.weather import Weather 

29 

30__all__ = ["refraction", "differentialRefraction"] 

31 

32# The differential refractive index is the difference of the refractive index from 1., 

33# multiplied by 1E8 to simplify the notation and equations. 

34deltaRefractScale = 1.0E8 

35 

36 

37def refraction(wavelength, elevation, observatory, weather=None): 

38 """Calculate overall refraction under atmospheric and observing conditions. 

39 

40 Parameters 

41 ---------- 

42 wavelength : `float` 

43 wavelength is in nm (valid for 230.2 < wavelength < 2058.6) 

44 elevation : `lsst.geom.Angle` 

45 Elevation of the observation, as an Angle. 

46 observatory : `lsst.afw.coord.Observatory` 

47 Class containing the longitude, latitude, 

48 and altitude of the observatory. 

49 weather : `lsst.afw.coord.Weather`, optional 

50 Class containing the measured temperature, pressure, and humidity 

51 at the observatory during an observation 

52 If omitted, typical conditions for the observatory's elevation will be calculated. 

53 

54 Returns 

55 ------- 

56 refraction : `lsst.geom.Angle` 

57 The angular refraction for light of the given wavelength, 

58 under the given observing conditions. 

59 

60 Notes 

61 ----- 

62 The calculation is taken from [1]_. 

63 

64 References 

65 ---------- 

66 .. [1] R. C. Stone, "An Accurate Method for Computing Atmospheric 

67 Refraction," Publications of the Astronomical Society of the Pacific, 

68 vol. 108, p. 1051, 1996. 

69 """ 

70 if wavelength < 230.2: 

71 raise ValueError("Refraction calculation is valid for wavelengths between 230.2 and 2058.6 nm.") 

72 if wavelength > 2058.6: 

73 raise ValueError("Refraction calculation is valid for wavelengths between 230.2 and 2058.6 nm.") 

74 latitude = observatory.getLatitude() 

75 altitude = observatory.getElevation() 

76 if weather is None: 

77 weather = defaultWeather(altitude*units.meter) 

78 reducedN = deltaN(wavelength, weather)/deltaRefractScale 

79 temperature = extractTemperature(weather, useKelvin=True) 

80 atmosScaleheightRatio = 4.5908E-6*temperature.to_value(units.Kelvin) 

81 

82 # Account for oblate Earth 

83 # This replicates equation 10 of Stone 1996 

84 relativeGravity = (1. + 0.005302*np.sin(latitude.asRadians())**2. 

85 - 0.00000583*np.sin(2.*latitude.asRadians())**2. - 0.000000315*altitude) 

86 

87 # Calculate the tangent of the zenith angle. 

88 tanZ = np.tan(np.pi/2. - elevation.asRadians()) 

89 atmosTerm1 = reducedN*relativeGravity*(1. - atmosScaleheightRatio) 

90 atmosTerm2 = reducedN*relativeGravity*(atmosScaleheightRatio - reducedN/2.) 

91 result = float(atmosTerm1*tanZ + atmosTerm2*tanZ**3.)*lsst.geom.radians 

92 return result 

93 

94 

95def differentialRefraction(wavelength, wavelengthRef, elevation, observatory, weather=None): 

96 """Calculate the differential refraction between two wavelengths. 

97 

98 Parameters 

99 ---------- 

100 wavelength : `float` 

101 wavelength is in nm (valid for 230.2 < wavelength < 2058.6) 

102 wavelengthRef : `float` 

103 Reference wavelength, typically the effective wavelength of a filter. 

104 elevation : `lsst.geom.Angle` 

105 Elevation of the observation, as an Angle. 

106 observatory : `lsst.afw.coord.Observatory` 

107 Class containing the longitude, latitude, 

108 and altitude of the observatory. 

109 weather : `lsst.afw.coord.Weather`, optional 

110 Class containing the measured temperature, pressure, and humidity 

111 at the observatory during an observation 

112 If omitted, typical conditions for the observatory's elevation will be calculated. 

113 

114 Returns 

115 ------- 

116 differentialRefraction : `lsst.geom.Angle` 

117 The refraction at `wavelength` minus the refraction at `wavelengthRef`. 

118 """ 

119 refractionStart = refraction(wavelength, elevation, observatory, weather=weather) 

120 refractionEnd = refraction(wavelengthRef, elevation, observatory, weather=weather) 

121 return refractionStart - refractionEnd 

122 

123 

124def deltaN(wavelength, weather): 

125 """Calculate the differential refractive index of air. 

126 

127 Parameters 

128 ---------- 

129 wavelength : `float` 

130 wavelength is in nanometers 

131 weather : `lsst.afw.coord.Weather` 

132 Class containing the measured temperature, pressure, and humidity 

133 at the observatory during an observation 

134 

135 Returns 

136 ------- 

137 deltaN : `float` 

138 The difference of the refractive index of air from 1., 

139 calculated as (n_air - 1)*10^8 

140 

141 Notes 

142 ----- 

143 The differential refractive index is the difference of 

144 the refractive index from 1., multiplied by 1E8 to simplify 

145 the notation and equations. Calculated as (n_air - 1)*10^8 

146 

147 This replicates equation 14 of [1]_ 

148 

149 References 

150 ---------- 

151 .. [1] R. C. Stone, "An Accurate Method for Computing Atmospheric 

152 Refraction," Publications of the Astronomical Society of the Pacific, 

153 vol. 108, p. 1051, 1996. 

154 """ 

155 waveNum = 1E3/wavelength # want wave number in units 1/micron 

156 dryAirTerm = 2371.34 + (683939.7/(130. - waveNum**2.)) + (4547.3/(38.9 - waveNum**2.)) 

157 wetAirTerm = 6487.31 + 58.058*waveNum**2. - 0.71150*waveNum**4. + 0.08851*waveNum**6. 

158 return (dryAirTerm*densityFactorDry(weather) + wetAirTerm*densityFactorWater(weather)) 

159 

160 

161def densityFactorDry(weather): 

162 """Calculate dry air pressure term to refractive index calculation. 

163 

164 Parameters 

165 ---------- 

166 weather : `lsst.afw.coord.Weather` 

167 Class containing the measured temperature, pressure, and humidity 

168 at the observatory during an observation 

169 

170 Returns 

171 ------- 

172 densityFactor : `float` 

173 Returns the relative density of dry air 

174 at the given pressure and temperature. 

175 

176 Notes 

177 ----- 

178 This replicates equation 15 of [1]_ 

179 

180 References 

181 ---------- 

182 .. [1] R. C. Stone, "An Accurate Method for Computing Atmospheric 

183 Refraction," Publications of the Astronomical Society of the Pacific, 

184 vol. 108, p. 1051, 1996. 

185 """ 

186 temperature = extractTemperature(weather, useKelvin=True) 

187 waterVaporPressure = humidityToPressure(weather) 

188 airPressure = weather.getAirPressure()*units.pascal 

189 dryPressure = airPressure - waterVaporPressure 

190 eqn = dryPressure.to_value(cds.mbar)*(57.90E-8 - 9.3250E-4/temperature.to_value(units.Kelvin) 

191 + 0.25844/temperature.to_value(units.Kelvin)**2.) 

192 densityFactor = (1. + eqn)*dryPressure.to_value(cds.mbar)/temperature.to_value(units.Kelvin) 

193 return densityFactor 

194 

195 

196def densityFactorWater(weather): 

197 """Calculate water vapor pressure term to refractive index calculation. 

198 

199 Parameters 

200 ---------- 

201 weather : `lsst.afw.coord.Weather` 

202 Class containing the measured temperature, pressure, and humidity 

203 at the observatory during an observation 

204 

205 Returns 

206 ------- 

207 densityFactor : `float` 

208 Returns the relative density of water vapor 

209 at the given pressure and temperature. 

210 

211 Notes 

212 ----- 

213 This replicates equation 16 of [1]_ 

214 

215 References 

216 ---------- 

217 .. [1] R. C. Stone, "An Accurate Method for Computing Atmospheric 

218 Refraction," Publications of the Astronomical Society of the Pacific, 

219 vol. 108, p. 1051, 1996. 

220 """ 

221 temperature = extractTemperature(weather, useKelvin=True) 

222 waterVaporPressure = humidityToPressure(weather) 

223 densityEqn1 = (-2.37321E-3 + 2.23366/temperature.to_value(units.Kelvin) 

224 - 710.792/temperature.to_value(units.Kelvin)**2. 

225 + 7.75141E-4/temperature.to_value(units.Kelvin)**3.) 

226 densityEqn2 = waterVaporPressure.to_value(cds.mbar)*(1. + 3.7E-4*waterVaporPressure.to_value(cds.mbar)) 

227 relativeDensity = waterVaporPressure.to_value(cds.mbar)/temperature.to_value(units.Kelvin) 

228 densityFactor = (1 + densityEqn2*densityEqn1)*relativeDensity 

229 

230 return densityFactor 

231 

232 

233def humidityToPressure(weather): 

234 """Convert humidity and temperature to water vapor pressure. 

235 

236 Parameters 

237 ---------- 

238 weather : `lsst.afw.coord.Weather` 

239 Class containing the measured temperature, pressure, and humidity 

240 at the observatory during an observation 

241 

242 Returns 

243 ------- 

244 pressure : `astropy.units.Quantity` 

245 The water vapor pressure in Pascals 

246 calculated from the given humidity and temperature. 

247 

248 Notes 

249 ----- 

250 This replicates equations 18 & 20 of [1]_ 

251 

252 References 

253 ---------- 

254 .. [1] R. C. Stone, "An Accurate Method for Computing Atmospheric 

255 Refraction," Publications of the Astronomical Society of the Pacific, 

256 vol. 108, p. 1051, 1996. 

257 """ 

258 humidity = weather.getHumidity() 

259 x = np.log(humidity/100.0) 

260 temperature = extractTemperature(weather) 

261 temperatureEqn1 = (temperature + 238.3*units.Celsius)*x + 17.2694*temperature 

262 temperatureEqn2 = (temperature + 238.3*units.Celsius)*(17.2694 - x) - 17.2694*temperature 

263 dewPoint = 238.3*temperatureEqn1/temperatureEqn2 

264 waterVaporPressure = (4.50874 + 0.341724*dewPoint + 0.0106778*dewPoint**2 + 0.184889E-3*dewPoint**3 

265 + 0.238294E-5*dewPoint**4 + 0.203447E-7*dewPoint**5)*133.32239*units.pascal 

266 

267 return waterVaporPressure 

268 

269 

270def extractTemperature(weather, useKelvin=False): 

271 """Thin wrapper to return the measured temperature from an observation. 

272 

273 Parameters 

274 ---------- 

275 weather : `lsst.afw.coord.Weather` 

276 Class containing the measured temperature, pressure, and humidity 

277 at the observatory during an observation 

278 useKelvin : bool, optional 

279 Set to True to return the temperature in Kelvin instead of Celsius 

280 This is needed because Astropy can't easily convert 

281 between Kelvin and Celsius. 

282 

283 Returns 

284 ------- 

285 temperature : `astropy.units.Quantity` 

286 The temperature in Celsius, unless `useKelvin` is set. 

287 """ 

288 temperature = weather.getAirTemperature()*units.Celsius 

289 if useKelvin: 

290 temperature = temperature.to(units.Kelvin, equivalencies=units.temperature()) 

291 return temperature 

292 

293 

294def defaultWeather(altitude): 

295 """Set default local weather conditions if they are missing. 

296 

297 Parameters 

298 ---------- 

299 weather : `lsst.afw.coord.Weather` 

300 Class containing the measured temperature, pressure, and humidity 

301 at the observatory during an observation 

302 altitude : `astropy.units.Quantity` 

303 The altitude of the observatory, in meters. 

304 

305 Returns 

306 ------- 

307 default : `lsst.afw.coord.Weather` 

308 Updated Weather class with any `nan` values replaced by defaults. 

309 """ 

310 if isinstance(altitude, units.quantity.Quantity): 

311 altitude2 = altitude 

312 else: 

313 altitude2 = altitude*units.meter 

314 p0 = 101325.*units.pascal # sea level air pressure 

315 g = 9.80665*units.meter/units.second**2 # typical gravitational acceleration at sea level 

316 R0 = 8.31447*units.Joule/(units.mol*units.Kelvin) # gas constant 

317 T0 = 19.*units.Celsius # Typical sea-level temperature 

318 lapseRate = -6.5*units.Celsius/units.km # Typical rate of change of temperature with altitude 

319 M = 0.0289644*units.kg/units.mol # molar mass of dry air 

320 

321 temperature = T0 + lapseRate*altitude2 

322 temperatureK = temperature.to(units.Kelvin, equivalencies=units.temperature()) 

323 pressure = p0*np.exp(-(g*M*altitude2)/(R0*temperatureK)) 

324 humidity = 40. # Typical humidity at many observatory sites. 

325 weather = Weather((temperature/units.Celsius).value, (pressure/units.pascal).value, humidity) 

326 return weather