Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
skymap.hh
Go to the documentation of this file.
1 //**********************************************************
2 // Wavelet Analysis Tool
3 // Sergey Klimenko, University of Florida
4 // sky map data container for network analysis
5 // used with DMT and ROOT
6 //**********************************************************
7 
8 #ifndef SKYMAP_HH
9 #define SKYMAP_HH
10 
11 #include <iostream>
12 #include <vector>
13 #include "complex"
14 #include "wavearray.hh"
15 #include "watfun.hh"
16 #include "TNamed.h"
17 #include "TFile.h"
18 #include "TString.h"
19 #ifdef _USE_HEALPIX
20 #include "alm.hh"
21 #endif
22 
23 #ifndef __CINT__
24 #ifdef _USE_HEALPIX
25 #include "xcomplex.h"
26 #include "paramfile.h"
27 #include "healpix_data_io.h"
28 #include "alm.h"
29 #include "healpix_map.h"
30 #include "healpix_map_fitsio.h"
31 #include "alm_healpix_tools.h"
32 #include "alm_powspec_tools.h"
33 #include "fitshandle.h"
34 #include "levels_facilities.h"
35 #include "lsconstants.h"
36 #endif
37 #endif
38 
39 typedef std::vector<double> vectorD;
40 
41 // GPS seconds of the J2000.0 epoch (2000 JAN 1 12h UTC).
42 #define EPOCH_J2000_0_GPS 630763213
43 
44 
45 class skymap : public TNamed
46 {
47  public:
48 
49  // constructors
50 
51  //: Default constructor
52  skymap();
53 
54  //: skymap constructor
55  //!param - step on phi and theta
56  //!param - theta begin
57  //!param - theta end
58  //!param - phi begin
59  //!param - phi end
60  skymap(double,double=0.,double=180.,double=0.,double=360.);
61 
62  //: skymap constructor
63  //!param - healpix order
64  skymap(int);
65 
66  //: skymap constructor
67  //!param - fits file
68  skymap(char*);
69 
70  //: skymap constructor
71  //!param ifile - root file name
72  //!param name - object name
73  skymap(TString ifile, TString name="skymap");
74 
75  //: Copy constructor
76  //!param: value - object to copy from
77  skymap(const skymap&);
78 
79  //: destructor
80  virtual ~skymap();
81 
82  // operators
83 
84  skymap& operator = (const skymap&);
85  skymap& operator += (const skymap&);
86  skymap& operator -= (const skymap&);
87  skymap& operator *= (const skymap&);
88  skymap& operator /= (const skymap&);
89  skymap& operator = (const double);
90  skymap& operator *= (const double);
91  skymap& operator += (const double);
92  char* operator >> (char* fname);
93 
94  // accessors
95 
96  //: get skymap value at index i
97  //: sets index nTheta and mPhi
98  //!param: sky index
99  double get(size_t i);
100 
101  //: set skymap value at index i,j
102  //!param: sky index
103  //!param: value to set
104  inline void set(size_t i, double a) {
105  if(int(i) != mIndex) this->get(i); // fill in sky index
106  if(mTheta>=0) value[mTheta][mPhi] = a;
107  }
108 
109  //: add skymap value at index i,j
110  //!param: sky index
111  //!param: value to add
112  inline void add(size_t i ,double a) {
113  if(int(i) != mIndex) this->get(i); // fill in sky index
114  if(mTheta>=0) value[mTheta][mPhi] += a;
115  }
116 
117  //: get skymap size
118  inline size_t size() {
119  size_t n = value.size();
120  size_t m = 0;
121  for(size_t i=0; i<n; i++) m += value[i].size();
122  return m;
123  }
124 
125  //: skymap sizes
126  inline size_t size(size_t k) {
127  if(!k) return value.size();
128  else if(k<=value.size()) return value[k-1].size();
129  else return 0;
130  }
131 
132  //: get sky index at theta,phi
133  //!param: theta
134  //!param: phi
135  size_t getSkyIndex(double th, double ph);
136 
137  //: get skymap value at theta,phi
138  //!param: theta
139  //!param: phi
140  inline double get(double th, double ph) {
141  mIndex = getSkyIndex(th,ph);
142  return get(mIndex);
143  }
144 
145  //: get phi value
146  inline double getPhi(size_t i) {
147  if(healpix!=NULL) {
148 #ifdef _USE_HEALPIX
149  pointing P = healpix->pix2ang(i);
150  return rad2deg*P.phi;
151 #else
152  return 0.;
153 #endif
154  } else {
155  if(int(i) != mIndex) this->get(i); // fill in sky index
156  size_t n = this->value[mTheta].size();
157  if(mTheta>=0 && n>1)
158  return phi_1+(mPhi+0.5)*(phi_2-phi_1)/double(n);
159  else return (phi_2-phi_1)/2.;
160  }
161  }
162 
163  //: get phi step
164  inline double getPhiStep(size_t i) {
165  if(int(i) != mIndex) this->get(i); // fill in sky index
166  size_t n = this->value[mTheta].size();
167  if(mTheta>=0 && n)
168  return (phi_2-phi_1)/double(n);
169  else return 0.;
170  }
171 
172  // get RA angle from EFEC phi angle in degrees
173  // Earth angular velocity defines duration of the siderial day
174  // 1 siderial day = 23h 56m 04.0905s
175  // (from http://www.maa.mhn.de/Scholar/times.html)
176  // phiRA function is implemented by Gabriele Vedovato and compared
177  // with LAL - agrees within 0.01 degrees
178 
179  static inline double phiRA(double ph, double gps, bool inverse=false) {
180  double sidereal_time; // sidereal time in sidereal seconds. (magic)
181  double t = (gps-EPOCH_J2000_0_GPS)/86400./36525.;
182 
183  sidereal_time = (-6.2e-6 * t + 0.093104) * t * t + 67310.54841;
184  sidereal_time += 8640184.812866*t + 3155760000.*t;
185 
186  double gmst = 360.*sidereal_time/86400.;
187  if(inverse) gmst=-gmst;
188 
189  double ra = fmod( ph + gmst, 360. );
190  return ra<0. ? ra+360 : ra;
191 
192  }
193 
194  inline double phi2RA(double ph, double gps) { return phiRA(ph,gps,false); }
195  inline double RA2phi(double ph, double gps) { return phiRA(ph,gps,true); }
196 
197  inline double getRA(size_t i) {
198  if(this->gps<=0.) return 0.;
199 // double omega = 7.292115090e-5; // http://hypertextbook.com/facts/2002/JasonAtkins.shtml
200 // double gps2000 = 630720013; // GPS time at 01/01/2000 00:00:00 UTC
201 // double GMST2000 = 24110.54841; // GMST at UT1=0 on Jan 1, 2000
202  return phi2RA(getPhi(i),this->gps);
203  }
204 
205  //: get theta value
206  inline double getTheta(size_t i) {
207  if(healpix!=NULL) {
208 #ifdef _USE_HEALPIX
209  pointing P = healpix->pix2ang(i);
210  return rad2deg*P.theta;
211 #else
212  return 0.;
213 #endif
214  } else {
215  size_t n = this->value.size()-1;
216  if(int(i) != mIndex) this->get(i); // fill in sky index
217  if(mTheta >= 0 && n)
218  return theta_1+mTheta*(theta_2-theta_1)/double(n);
219  else return (theta_2-theta_1)/2.;
220  }
221  }
222 
223  //: get theta step
224  inline double getThetaStep(size_t i) {
225  size_t n = this->value.size()-1;
226  if(int(i) != mIndex) this->get(i); // fill in sky index
227  if(mTheta >= 0 && n)
228  return (theta_2-theta_1)/double(n);
229  else return 0.;
230  }
231 
232  // get declination angle from EFEC theta angle in degrees
233  inline double getDEC(size_t i) { return -(getTheta(i)-90.); }
234 
235  //: find and return maximum value in skymap
236  //: set mTheta and mPhi to be theta and phi index for the maximum value
237  double max();
238 
239  //: find and return minimum value in skymap
240  //: set mTheta and mPhi to be theta and phi index for the minimum value
241  double min();
242 
243  //: find and return mean value over entire skymap
244  double mean();
245 
246  //: find for what fraction of the sky the statistic > threshold t
247  double fraction(double);
248 
249  //: normalize skymap
250  double norm(double=0.);
251 
252  //: downsample sky index array
253  //: input sky index array
254  //: down-sampling option: 2 or 4
255  void downsample(wavearray<short>&, size_t=4);
256 
257  //: dump skymap into binary file
258  //: parameter: file name
259  //: parameter: append mode
260  void DumpBinary(char*, int);
261 
262 #ifdef _USE_HEALPIX
263  //: dump skymap into fits file
264  //: file : output fits file name [ext : .fits or .fits.gz]
265  //: gps_obs : gps time of the observation
266  //: configur : software configuration used to process the data
267  //: TTYPE1 : label for field 1
268  //: TUNIT1 : physical unit of field 1
269  //: coordsys : Pixelisation coordinate system [c/C : select celestial]
270  void Dump2fits(const char* file, double gps_obs=0, const char configur[]="",
271  const char TTYPE1[]="", const char TUNIT1[]="", char coordsys='x'); // *MENU*
272 #endif
273 
274  //: dump skymap object into root file
275  void DumpObject(char*);
276 
277  // return 1 if healpix else 0
278  int getType() {if(healpix==NULL) return 0; else return 1;}
279 
280  // works only with the healpix scheme !!!
281 #ifdef _USE_HEALPIX
282  wavearray<int> neighbors(int index);
283  void median(double radius);
284  void smoothing(double fwhm, int nlmax=256, int num_iter=0);
285  void rotate(double psi, double theta, double phi, int nlmax=256, int num_iter=0);
286  void setAlm(wat::Alm alm);
287  wat::Alm getAlm(int nlmax, int num_iter=0);
288  void resample(int order);
289  int getRings();
290  int getRingPixels(int ring);
291  int getStartRingPixel(int ring);
292  int getEulerCharacteristic(double threshold);
293 #endif
294 
295  // return the healpix order (if=0 healpix skygrid is disabled)
296  int getOrder() {return healpix_order;}
297 
298 // data members
299 
300  std::vector<vectorD> value; // skymap map array
301  std::vector<int> index; // sample index array
302  double sms; // step on phi and theta
303  double theta_1; // theta range begin
304  double theta_2; // theta range end
305  double phi_1; // phi range begin
306  double phi_2; // phi range end
307  double gps; // gps time
308  int mTheta; // theta index
309  int mPhi; // phi index
310  int mIndex; // sky index
311 
312 private:
313  int healpix_order; // healpix order (if=0 healpix is disabled)
314 
315 #ifndef __CINT__
316 #ifdef _USE_HEALPIX
317  Healpix_Map<double>* healpix;//!
318 #else
319  bool* healpix;
320 #endif
321  double deg2rad;
322  double rad2deg;
323 #endif
324 
325  ClassDef(skymap,5)
326 
327 }; // class skymap
328 
329 using namespace std;
330 
331 namespace {
332 
333 template<typename Iterator> typename iterator_traits<Iterator>::value_type
334  median(Iterator first, Iterator last)
335  {
336  Iterator mid = first+(last-first-1)/2;
337  nth_element(first,mid,last);
338  if ((last-first)&1) return *mid;
339  return typename iterator_traits<Iterator>::value_type
340  (0.5*((*mid)+(*min_element(mid+1,last))));
341  }
342 
343 } // unnamed namespace
344 
345 #endif // SKYMAP_HH
346 
wavearray< double > t(hp.size())
double norm(double=0.)
Definition: skymap.cc:496
int mPhi
Definition: skymap.hh:309
#define EPOCH_J2000_0_GPS
Definition: skymap.hh:42
wavearray< double > a(hp.size())
int mTheta
Definition: skymap.hh:308
skymap & operator-=(const skymap &)
Definition: skymap.cc:299
par[0] name
int n
Definition: cwb_net.C:10
TString("c")
double gps
Definition: skymap.hh:307
float theta
double deg2rad
Definition: skymap.hh:321
double median(std::vector< double > &vec)
Definition: wavegraph.cc:467
delete[] radius
TH2F * ph
skymap & operator*=(const skymap &)
Definition: skymap.cc:328
TRandom3 P
Definition: compare_bkg.C:296
virtual ~skymap()
Definition: skymap.cc:232
double min()
Definition: skymap.cc:441
STL namespace.
double getTheta(size_t i)
Definition: skymap.hh:206
std::vector< int > index
Definition: skymap.hh:301
double getThetaStep(size_t i)
Definition: skymap.hh:224
double theta_1
Definition: skymap.hh:303
int m
Definition: cwb_net.C:10
double getPhiStep(size_t i)
Definition: skymap.hh:164
double max()
Definition: skymap.cc:420
void downsample(wavearray< short > &, size_t=4)
Definition: skymap.cc:518
i drho i
double theta_2
Definition: skymap.hh:304
double mean()
Definition: skymap.cc:462
skymap()
Definition: skymap.cc:29
double rad2deg
Definition: skymap.hh:322
size_t getSkyIndex(double th, double ph)
param: theta param: phi
Definition: skymap.cc:702
float phi
double ra
Definition: ConvertGWGC.C:46
std::vector< vectorD > value
Definition: skymap.hh:300
float psi
double fraction(double)
Definition: skymap.cc:479
skymap & operator+=(const skymap &)
Definition: skymap.cc:270
double getRA(size_t i)
Definition: skymap.hh:197
void DumpObject(char *)
Definition: skymap.cc:649
int getOrder()
Definition: skymap.hh:296
double sms
Definition: skymap.hh:302
double getDEC(size_t i)
Definition: skymap.hh:233
char fname[1024]
Class for storing spherical harmonic coefficients.
char * operator>>(char *fname)
Definition: skymap.cc:1092
int mIndex
Definition: skymap.hh:310
int k
double RA2phi(double ph, double gps)
Definition: skymap.hh:195
std::vector< double > vectorD
Definition: skymap.hh:39
Definition: skymap.hh:45
void DumpBinary(char *, int)
Definition: skymap.cc:657
TFile * ifile
double phi2RA(double ph, double gps)
Definition: skymap.hh:194
Definition: alm.hh:175
double getPhi(size_t i)
Definition: skymap.hh:146
int getType()
Definition: skymap.hh:278
string file
Definition: cwb_online.py:385
skymap & operator=(const skymap &)
Definition: skymap.cc:238
int healpix_order
Definition: skymap.hh:313
size_t size(size_t k)
Definition: skymap.hh:126
void set(size_t i, double a)
param: sky index param: value to set
Definition: skymap.hh:104
double phi_2
Definition: skymap.hh:306
double phi_1
Definition: skymap.hh:305
bool * healpix
Definition: skymap.hh:319
size_t size()
Definition: skymap.hh:118
void add(size_t i, double a)
param: sky index param: value to add
Definition: skymap.hh:112
skymap & operator/=(const skymap &)
Definition: skymap.cc:355
static double phiRA(double ph, double gps, bool inverse=false)
Definition: skymap.hh:179