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

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

from builtins import zip 

from builtins import object 

import os 

import numpy 

from astropy.io import fits 

 

from lsst.sims.utils.CodeUtilities import sims_clean_up 

from lsst.sims.utils import _galacticFromEquatorial 

from functools import reduce 

 

__all__ = ["EBVmap", "EBVbase"] 

 

 

def interp1D(z1, z2, offset): 

""" 1D interpolation on a grid""" 

 

zPrime = (z2-z1)*offset + z1 

 

return zPrime 

 

 

class EBVmap(object): 

'''Class for describing a map of EBV 

 

Images are read in from a fits file and assume a ZEA projection 

''' 

 

def __del__(self): 

self.hdulist.close() 

 

def readMapFits(self, fileName): 

""" read a fits file containing the ebv data""" 

 

self._file_name = fileName 

 

self.hdulist = fits.open(fileName) 

self.header = self.hdulist[0].header 

self.data = self.hdulist[0].data 

self.nr = self.data.shape[0] 

self.nc = self.data.shape[1] 

 

# read WCS information 

self.cd11 = self.header['CD1_1'] 

self.cd22 = self.header['CD2_2'] 

self.cd12 = 0. 

self.cd21 = 0. 

 

self.crpix1 = self.header['CRPIX1'] 

self.crval1 = self.header['CRVAL1'] 

 

self.crpix2 = self.header['CRPIX2'] 

self.crval2 = self.header['CRVAL2'] 

 

# read projection information 

self.nsgp = self.header['LAM_NSGP'] 

self.scale = self.header['LAM_SCAL'] 

self.lonpole = self.header['LONPOLE'] 

 

def xyFromSky(self, gLon, gLat): 

""" convert long, lat angles to pixel x y 

 

input angles are in radians but the conversion assumes radians 

 

@param [in] gLon galactic longitude in radians 

 

@param [in] gLat galactic latitude in radians 

 

@param [out] x is the x pixel coordinate 

 

@param [out] y is the y pixel coordinate 

 

""" 

 

rad2deg = 180./numpy.pi 

 

# use the SFD approach to define xy pixel positions 

# ROTATION - Equn (4) - degenerate case 

if (self.crval2 > 89.9999): 

theta = gLat*rad2deg 

phi = gLon*rad2deg + 180.0 + self.lonpole - self.crval1 

elif (self.crval2 < -89.9999): 

theta = -gLat*rad2deg 

phi = self.lonpole + self.crval1 - gLon*rad2deg 

else: 

# Assume it's an NGP projection ... 

theta = gLat*rad2deg 

phi = gLon*rad2deg + 180.0 + self.lonpole - self.crval1 

 

# Put phi in the range [0,360) degrees 

phi = phi - 360.0 * numpy.floor(phi/360.0) 

 

# FORWARD MAP PROJECTION - Equn (26) 

Rtheta = 2.0 * rad2deg * numpy.sin((0.5 / rad2deg) * (90.0 - theta)) 

 

# Equns (10), (11) 

xr = Rtheta * numpy.sin(phi / rad2deg) 

yr = - Rtheta * numpy.cos(phi / rad2deg) 

 

# SCALE FROM PHYSICAL UNITS - Equn (3) after inverting the matrix 

denom = self.cd11 * self.cd22 - self.cd12 * self.cd21 

x = (self.cd22 * xr - self.cd12 * yr) / denom + (self.crpix1 - 1.0) 

y = (self.cd11 * yr - self.cd21 * xr) / denom + (self.crpix2 - 1.0) 

 

return x, y 

 

def generateEbv(self, glon, glat, interpolate = False): 

""" 

Calculate EBV with option for interpolation 

 

@param [in] glon galactic longitude in radians 

 

@param [in] galactic latitude in radians 

 

@param [out] ebvVal the scalar value of EBV extinction 

 

""" 

 

# calculate pixel values 

x, y = self.xyFromSky(glon, glat) 

 

ix = (x + 0.5).astype(int) 

iy = (y + 0.5).astype(int) 

 

unity = numpy.ones(len(ix), dtype=int) 

 

if (interpolate): 

 

# find the indices of the pixels bounding the point of interest 

ixLow = numpy.minimum(ix, (self.nc - 2)*unity) 

ixHigh = ixLow + 1 

dx = x - ixLow 

 

iyLow = numpy.minimum(iy, (self.nr - 2)*unity) 

iyHigh = iyLow + 1 

dy = y - iyLow 

 

# interpolate the EBV value at the point of interest by interpolating 

# first in x and then in y 

x1 = numpy.array([self.data[ii][jj] for (ii, jj) in zip(iyLow, ixLow)]) 

x2 = numpy.array([self.data[ii][jj] for (ii, jj) in zip(iyLow, ixHigh)]) 

xLow = interp1D(x1, x2, dx) 

 

x1 = numpy.array([self.data[ii][jj] for (ii, jj) in zip(iyHigh, ixLow)]) 

x2 = numpy.array([self.data[ii][jj] for (ii, jj) in zip(iyHigh, ixHigh)]) 

xHigh = interp1D(x1, x2, dx) 

 

ebvVal = interp1D(xLow, xHigh, dy) 

 

else: 

ebvVal = numpy.array([self.data[ii][jj] for (ii, jj) in zip(iy, ix)]) 

 

return ebvVal 

 

def xyIntFromSky(self, gLong, gLat): 

x, y = self.xyFromSky(gLong, gLat) 

ix = int(x + 0.5) 

iy = int(y + 0.5) 

 

return ix, iy 

 

 

class EBVbase(object): 

""" 

This class will give users access to calculateEbv oustide of the framework of a catalog. 

 

To find the value of EBV at a point on the sky, create an instance of this object, and 

then call calculateEbv passing the coordinates of interest as kwargs 

 

e.g. 

 

ebvObject = EBVbase() 

ebvValue = ebvObject.calculateEbv(galacticCoordinates = myGalacticCoordinates) 

 

or 

 

ebvValue = ebvObject.calculateEbv(equatorialCoordinates = myEquatorialCoordinates) 

 

where myGalacticCoordinates is a 2-d numpy array where the first row is galactic longitude 

and the second row is galactic latitude. 

 

myEquatorialCoordinates is a 2-d numpy array where the first row is RA and the second row 

is dec 

 

All coordinates are in radians 

 

You can also specify dust maps in the northern and southern galactic hemispheres, but 

there are default values that the code will automatically load (see the class variables 

below). 

 

The information regarding where the dust maps are located is stored in 

member variables ebvDataDir, ebvMapNorthName, ebvMapSouthName 

 

The actual dust maps (when loaded) are stored in ebvMapNorth and ebvMapSouth 

""" 

 

# these variables will tell the mixin where to get the dust maps 

ebvDataDir = os.environ.get("SIMS_MAPS_DIR") 

ebvMapNorthName = "DustMaps/SFD_dust_4096_ngp.fits" 

ebvMapSouthName = "DustMaps/SFD_dust_4096_sgp.fits" 

ebvMapNorth = None 

ebvMapSouth = None 

 

# A dict to hold every open instance of an EBVmap. 

# Since this is being declared outside of the constructor, 

# it will be a class member, which means that, every time 

# an EBVmap is added to the cache, all EBVBase instances will 

# know about it. 

_ebv_map_cache = {} 

 

# the set_xxxx routines below will allow the user to point elsewhere for the dust maps 

def set_ebvMapNorth(self, word): 

""" 

This allows the user to pick a new northern SFD map file 

""" 

self.ebvMapNorthName = word 

 

def set_ebvMapSouth(self, word): 

""" 

This allows the user to pick a new southern SFD map file 

""" 

self.ebvMapSouthName = word 

 

# these routines will load the dust maps for the galactic north and south hemispheres 

def _load_ebv_map(self, file_name): 

""" 

Load the EBV map specified by file_name. If that map has already been loaded, 

just return the map stored in self._ebv_map_cache. If it must be loaded, store 

it in the cache. 

""" 

if file_name in self._ebv_map_cache: 

return self._ebv_map_cache[file_name] 

 

ebv_map = EBVmap() 

ebv_map.readMapFits(file_name) 

self._ebv_map_cache[file_name] = ebv_map 

return ebv_map 

 

def load_ebvMapNorth(self): 

""" 

This will load the northern SFD map 

""" 

file_name = os.path.join(self.ebvDataDir, self.ebvMapNorthName) 

self.ebvMapNorth = self._load_ebv_map(file_name) 

return None 

 

def load_ebvMapSouth(self): 

""" 

This will load the southern SFD map 

""" 

file_name = os.path.join(self.ebvDataDir, self.ebvMapSouthName) 

self.ebvMapSouth = self._load_ebv_map(file_name) 

return None 

 

def calculateEbv(self, galacticCoordinates=None, equatorialCoordinates=None, northMap=None, southMap=None, 

interp=False): 

""" 

For an array of Gal long, lat calculate E(B-V) 

 

 

@param [in] galacticCoordinates is a numpy.array; the first row is galactic longitude, 

the second row is galactic latitude in radians 

 

@param [in] equatorialCoordinates is a numpy.array; the first row is RA, the second row is Dec in 

radians 

 

@param [in] northMap the northern dust map 

 

@param [in] southMap the southern dust map 

 

@param [in] interp is a boolean determining whether or not to interpolate the EBV value 

 

@param [out] ebv is a list of EBV values for all of the gLon, gLat pairs 

 

""" 

 

# raise an error if the coordinates are specified in both systems 

if galacticCoordinates is not None: 

if equatorialCoordinates is not None: 

raise RuntimeError("Specified both galacticCoordinates and " 

"equatorialCoordinates in calculateEbv") 

 

# convert (ra,dec) into gLon, gLat 

if galacticCoordinates is None: 

 

# raise an error if you did not specify ra or dec 

if equatorialCoordinates is None: 

raise RuntimeError("Must specify coordinates in calculateEbv") 

 

galacticCoordinates = numpy.array(_galacticFromEquatorial(equatorialCoordinates[0, :], 

equatorialCoordinates[1, :])) 

 

if northMap is None: 

if self.ebvMapNorth is None: 

self.load_ebvMapNorth() 

 

northMap = self.ebvMapNorth 

 

if southMap is None: 

if self.ebvMapSouth is None: 

self.load_ebvMapSouth() 

 

southMap = self.ebvMapSouth 

 

ebv = None 

 

if galacticCoordinates.shape[1] > 0: 

 

ebv = numpy.zeros(len(galacticCoordinates[0, :])) 

 

# identify (by index) which points are in the galactic northern hemisphere 

# and which points are in the galactic southern hemisphere 

# taken from 

# http://stackoverflow.com/questions/4578590/python-equivalent-of-filter-getting-two-output-lists-i-e-partition-of-a-list 

inorth, isouth = reduce(lambda x, y: x[not y[1] > 0.0].append(y[0]) or x, 

enumerate(galacticCoordinates[1, :]), ([], [])) 

 

nSet = galacticCoordinates[:, inorth] 

sSet = galacticCoordinates[:, isouth] 

 

ebvNorth = northMap.generateEbv(nSet[0, :], nSet[1, :], interpolate=interp) 

ebvSouth = southMap.generateEbv(sSet[0, :], sSet[1, :], interpolate=interp) 

 

for (i, ee) in zip(inorth, ebvNorth): 

ebv[i] = ee 

 

for (i, ee) in zip(isouth, ebvSouth): 

ebv[i] = ee 

 

return ebv 

 

 

sims_clean_up.targets.append(EBVbase._ebv_map_cache)