Logo coherent WaveBurst  
Library Reference Guide
Logo
skymap.cc
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 Sergey Klimenko, Gabriele Vedovato
3 #
4 # This program is free software: you can redistribute it and/or modify
5 # it under the terms of the GNU General Public License as published by
6 # the Free Software Foundation, either version 3 of the License, or
7 # (at your option) any later version.
8 #
9 # This program is distributed in the hope that it will be useful,
10 # but WITHOUT ANY WARRANTY; without even the implied warranty of
11 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 # GNU General Public License for more details.
13 #
14 # You should have received a copy of the GNU General Public License
15 # along with this program. If not, see <https://www.gnu.org/licenses/>.
16 */
17 
18 
19 // ++++++++++++++++++++++++++++++++++++++++++++++
20 // S. Klimenko, University of Florida
21 // WAT sky map class
22 // ++++++++++++++++++++++++++++++++++++++++++++++
23 
24 #define SKYMAP_CC
25 #include "skymap.hh"
26 #include "wavearray.hh"
27 #include "TObjArray.h"
28 #include "TObjString.h"
29 #include "TSystem.h"
30 #include "time.hh"
31 
32 using namespace std;
33 
34 #ifdef _USE_HEALPIX
35 unsigned int healpix_repcount (tsize npix)
36  {
37  if (npix<1024) return 1;
38  if ((npix%1024)==0) return 1024;
39  return isqrt (npix/12);
40  }
41 #endif
42 
43 ClassImp(skymap)
44 
45 // constructors
46 
48 
49  deg2rad = PI/180.;
50  rad2deg = 180./PI;
51  healpix = NULL;
52  healpix_order = 0;
53 }
54 
55 skymap::skymap(double sms,double t1,double t2,double p1,double p2)
56 {
57  int i,j;
58  int ntheta, nphi;
59  double c,s,p,d;
60  double a = PI/180.;
61 
62  if(t1<0 || t1>180. || t2<0 || t2>180. || t2<t1) {
63  cout<<"skymap(): invalid theta parameters"<<endl;
64  return;
65  }
66 
67  if(p1<0 || p1>360. || p2<0 || p2>360. || p2<p1) {
68  cout<<"skymap(): invalid phi parameters"<<endl;
69  return;
70  }
71 
72  this->sms = sms; // step on phi and theta
73  theta_1 = t1; // theta range
74  theta_2 = t2; // set theta range
75  phi_1 = p1; // phi range
76  phi_2 = p2; // phi range
77 
78  ntheta = 2*size_t((t2-t1)/sms/2.)+1; // divisions on theta
79  d = ntheta>1 ? (t2-t1)/(ntheta-1) : sms;
80 
81  this->index.reserve(ntheta); // reserve index array
82  this->value.reserve(ntheta); // reserve value array
83 
84  for(i=-ntheta/2; i<=ntheta/2; i++) {
85  c = cos(((t2+t1)/2.+i*d)*a);
86  s = sin(((t2+t1)/2.+i*d)*a);
87  c = s>0. ? (cos(sms*a)-c*c)/s/s : 2.;
88  p = c>0&&c<1. ? acos(c)/a : 361.; // step on phi
89  nphi = 2*size_t((p2-p1)/p/2.)+1; // divisions on phi
90 
91  std::vector<double> v;
92 
93  v.reserve(nphi);
94  for(j=0; j<nphi; j++) v.push_back(0.);
95  index.push_back(nphi); // initialize index array
96  value.push_back(v); // initialize value array
97  }
98  mIndex = mTheta = mPhi = -1;
99  gps = -1.;
100 
101  deg2rad = PI/180.;
102  rad2deg = 180./PI;
103  healpix = NULL;
104  healpix_order = 0;
105 
106 }
107 
109 {
110 #ifdef _USE_HEALPIX
111  int i,j;
112  int ntheta, nphi;
113 
114  sms = 0.0; // step on phi and theta
115  theta_1 = 0.0; // theta range
116  theta_2 = 180.0; // set theta range
117  phi_1 = 0.0; // phi range
118  phi_2 = 360.0; // phi range
119 
120  healpix = new Healpix_Map<double>(); //Healpix MAP
121  read_Healpix_map_from_fits(ifile,*healpix,1,2);
122  this->healpix_order = healpix->Order();
123 
124  ntheta = 1; // divisions on theta
125  for(i=0;i<=healpix->Order();i++) ntheta=2*ntheta+1;
126 
127  this->index.reserve(ntheta); // reserve index array
128  this->value.reserve(ntheta); // reserve value array
129  for(i=1; i<=ntheta; i++) {
130 
131  int ring=i; int startpix; int ringpix; double costheta; double sintheta; bool shifted=false;
132  healpix->get_ring_info(ring, startpix, ringpix, costheta, sintheta, shifted);
133  nphi=ringpix; // divisions on phi
134 
135  std::vector<double> v;
136 
137  v.reserve(nphi);
138  for(j=0; j<nphi; j++) v.push_back(0.);
139  index.push_back(nphi); // initialize index array
140  value.push_back(v); // initialize value array
141  }
142  mIndex = mTheta = mPhi = -1;
143  gps = -1.;
144 
145  // fill skymap
146  int L = healpix->Npix();
147  for(int l=0;l<L;l++) this->set(l,(*healpix)[l]);
148 
149  deg2rad = PI/180.;
150  rad2deg = 180./PI;
151 #else
152  cout << "skymap::skymap(char* ifile) : not available - healpix not enabled" << endl;
153  exit(1);
154 #endif
155 }
156 
157 skymap::skymap(int healpix_order)
158 {
159 #ifdef _USE_HEALPIX
160  int i,j;
161  int ntheta, nphi;
162 
163  sms = 0.0; // step on phi and theta
164  theta_1 = 0.0; // theta range
165  theta_2 = 180.0; // set theta range
166  phi_1 = 0.0; // phi range
167  phi_2 = 360.0; // phi range
168 
169  ntheta = 1; // divisions on theta
170  for(i=0;i<=healpix_order;i++) ntheta=2*ntheta+1;
171 
172  healpix = new Healpix_Map<double>(healpix_order,RING); //Healpix MAP
173  int L = healpix->Npix();
174  this->healpix_order = healpix->Order();
175 
176 // TEST
177 /*
178  char pixDir[8][3] = {"SW", "W", "NW", "N", "NE", "E", "SE", "S"};
179  int ring=0; int startpix; int ringpix; double costheta; double sintheta; bool shifted=false;
180  for(int i=1;i<=ntheta;i++) {
181  ring=i;
182  healpix->get_ring_info(ring, startpix, ringpix, costheta, sintheta, shifted);
183  cout << "ring : " << ring << " startpix : " << startpix << " ringpix : " << ringpix << endl;
184  fix_arr< int, 8 > result;
185  for(int j=startpix;j<startpix+ringpix;j++) {
186  healpix->neighbors (j, result);
187  for(int k=0;k<8;k++) cout << "pixel : " << j << " pixDir : " << pixDir[k] << " " << result[k] << endl;
188  cout << endl;
189  }
190  }
191 */
192 // END TEST
193 
194  this->index.reserve(ntheta); // reserve index array
195  this->value.reserve(ntheta); // reserve value array
196  for(i=1; i<=ntheta; i++) {
197 
198  int ring=i; int startpix; int ringpix; double costheta; double sintheta; bool shifted=false;
199  healpix->get_ring_info(ring, startpix, ringpix, costheta, sintheta, shifted);
200  nphi=ringpix; // divisions on phi
201 
202  std::vector<double> v;
203 
204  v.reserve(nphi);
205  for(j=0; j<nphi; j++) v.push_back(0.);
206  index.push_back(nphi); // initialize index array
207  value.push_back(v); // initialize value array
208  }
209  mIndex = mTheta = mPhi = -1;
210  gps = -1.;
211 
212  deg2rad = PI/180.;
213  rad2deg = 180./PI;
214 #else
215  cout << "skymap::skymap(int healpix_order) : not available - healpix not enabled" << endl;
216  exit(1);
217 #endif
218 }
219 
221 {
222  healpix = NULL;
223  healpix_order = 0;
224 
225  TFile *rfile = TFile::Open(ifile);
226  if (rfile==NULL) {
227  cout << "skymap::skymap - Error : No " << ifile.Data() << " file !!!" << endl;
228  exit(1);
229  }
230  // read skymap object
231  if(rfile->Get(name)!=NULL) {
232  *this = *(skymap*)rfile->Get(name);
233  } else {
234  cout << "skymap::skymap - Error : skymap is not contained in root file " << ifile.Data() << endl;
235  exit(1);
236  }
237  rfile->Close();
238 }
239 
240 skymap::skymap(const skymap& value)
241 {
242  healpix = NULL;
243  healpix_order = 0;
244 
245  *this = value;
246 }
247 
248 // destructor
249 
251  if(healpix!=NULL) delete healpix;
252 }
253 
254 //: operators
255 
256 skymap& skymap::operator=(const skymap& a)
257 {
258  if(healpix!=NULL) {delete healpix;healpix=NULL;}
259 
260  if(a.healpix!=NULL) {
261 #ifdef _USE_HEALPIX
262  healpix = new Healpix_Map<double>(a.healpix->Order(),a.healpix->Scheme());
263  this->healpix_order = healpix->Order();
264 #endif
265  } else {
266  healpix = NULL;
267  this->healpix_order = 0;
268  }
269 
270  sms = a.sms; // step on phi and theta
271  theta_1 = a.theta_1; // theta range
272  theta_2 = a.theta_2; // set theta range
273  phi_1 = a.phi_1; // phi range
274  phi_2 = a.phi_2; // phi range
275  value = a.value;
276  index = a.index;
277  mTheta = a.mTheta;
278  mPhi = a.mPhi;
279  mIndex = a.mIndex;
280  gps = a.gps;
281 
282  deg2rad = PI/180.;
283  rad2deg = 180./PI;
284 
285  return *this;
286 }
287 
288 skymap& skymap::operator+=(const skymap& a)
289 {
290  size_t i,j;
291  size_t n = value.size();
292  size_t m;
293 
294  if( a.theta_1 != theta_1 ||
295  a.theta_2 != theta_2 ||
296  a.phi_1 != phi_1 ||
297  a.phi_2 != phi_2 ||
298  a.sms != sms ||
299  a.value.size() != n )
300  {
301  cout<<"skymap::operator+ - incompatible skymaps"<<endl;
302  }
303 
304 
305  for(i=0; i<n; i++) {
306  m = value[i].size();
307  if(m != a.value[i].size()) {
308  cout<<"skymap::operator+ - incompatible skymaps"<<endl;
309  break;
310  }
311  for(j=0; j<m; j++) { value[i][j] += a.value[i][j]; }
312  }
313 
314  return *this;
315 }
316 
317 skymap& skymap::operator-=(const skymap& a)
318 {
319  size_t i,j;
320  size_t n = value.size();
321  size_t m;
322 
323  if( a.theta_1 != theta_1 ||
324  a.theta_2 != theta_2 ||
325  a.phi_1 != phi_1 ||
326  a.phi_2 != phi_2 ||
327  a.sms != sms ||
328  a.value.size() != n )
329  {
330  cout<<"skymap::operator- - incompatible skymaps"<<endl;
331  }
332 
333 
334  for(i=0; i<n; i++) {
335  m = value[i].size();
336  if(m != a.value[i].size()) {
337  cout<<"skymap::operator- - incompatible skymaps"<<endl;
338  break;
339  }
340  for(j=0; j<m; j++) { value[i][j] -= a.value[i][j]; }
341  }
342 
343  return *this;
344 }
345 
346 skymap& skymap::operator*=(const skymap& a)
347 {
348  size_t i,j,m;
349  size_t n = value.size();
350 
351  if( a.theta_1 != theta_1 ||
352  a.theta_2 != theta_2 ||
353  a.phi_1 != phi_1 ||
354  a.phi_2 != phi_2 ||
355  a.sms != sms ||
356  a.value.size() != n )
357  {
358  cout<<"skymap::operator* - incompatible skymaps"<<endl;
359  }
360 
361 
362  for(i=0; i<n; i++) {
363  m = value[i].size();
364  if(m != a.value[i].size()) {
365  cout<<"skymap::operator* - incompatible skymaps"<<endl;
366  break;
367  }
368  for(j=0; j<m; j++) { value[i][j] *= a.value[i][j]; }
369  }
370  return *this;
371 }
372 
373 skymap& skymap::operator/=(const skymap& a)
374 {
375  size_t i,j,m;
376  size_t n = value.size();
377 
378  if( a.theta_1 != theta_1 ||
379  a.theta_2 != theta_2 ||
380  a.phi_1 != phi_1 ||
381  a.phi_2 != phi_2 ||
382  a.sms != sms ||
383  a.value.size() != n )
384  {
385  cout<<"skymap::operator/ - incompatible skymaps"<<endl;
386  }
387 
388 
389  for(i=0; i<n; i++) {
390  m = value[i].size();
391  if(m != a.value[i].size()) {
392  cout<<"skymap::operator/ - incompatible skymaps"<<endl;
393  break;
394  }
395  for(j=0; j<m; j++) { value[i][j] /= a.value[i][j]!=0. ? a.value[i][j] : 1.; }
396  }
397  return *this;
398 }
399 
400 skymap& skymap::operator=(const double a)
401 {
402  size_t i,j;
403  size_t n = value.size();
404  size_t m;
405 
406  for(i=0; i<n; i++) {
407  m = value[i].size();
408  for(j=0; j<m; j++) { this->value[i][j] = a; }
409  }
410  return *this;
411 }
412 
413 skymap& skymap::operator*=(const double a)
414 {
415  size_t i,j;
416  size_t n = value.size();
417  size_t m;
418 
419  for(i=0; i<n; i++) {
420  m = value[i].size();
421  for(j=0; j<m; j++) { this->value[i][j] *= a; }
422  }
423  return *this;
424 }
425 
426 skymap& skymap::operator+=(const double a)
427 {
428  size_t i,j,m;
429  size_t n = value.size();
430 
431  for(i=0; i<n; i++) {
432  m = value[i].size();
433  for(j=0; j<m; j++) { this->value[i][j] += a; }
434  }
435  return *this;
436 }
437 
438 double skymap::max()
439 {
440  size_t i,j,m;
441  size_t k = 0;
442  size_t n = value.size();
443  double a = -1.e100;
444  double x;
445 
446  mTheta = mPhi = mIndex = -1;
447 
448  for(i=0; i<n; i++) {
449  m = value[i].size();
450  for(j=0; j<m; j++) {
451  x = this->value[i][j];
452  if(x>a) { a=x; mTheta=i; mPhi=j; mIndex=k; }
453  k++;
454  }
455  }
456  return a;
457 }
458 
459 double skymap::min()
460 {
461  size_t i,j,m;
462  size_t k = 0;
463  size_t n = value.size();
464  double a = 1.e100;
465  double x;
466 
467  mTheta = mPhi = mIndex = -1;
468 
469  for(i=0; i<n; i++) {
470  m = value[i].size();
471  for(j=0; j<m; j++) {
472  x = this->value[i][j];
473  if(x<a) { a=x; mTheta=i; mPhi=j; mIndex=k; }
474  k++;
475  }
476  }
477  return a;
478 }
479 
480 double skymap::mean()
481 {
482  size_t i,j,m;
483  size_t n = value.size();
484  size_t k = 0;
485  double a = 0.;
486 
487  for(i=0; i<n; i++) {
488  m = value[i].size();
489  k += m;
490  for(j=0; j<m; j++) {
491  a += this->value[i][j];
492  }
493  }
494  return a/double(k);
495 }
496 
497 double skymap::fraction(double t)
498 {
499  size_t i,j,m;
500  size_t n = value.size();
501  size_t k = 0;
502  double a = 0.;
503 
504  for(i=0; i<n; i++) {
505  m = value[i].size();
506  k += m;
507  for(j=0; j<m; j++) {
508  if(this->value[i][j]>t) a += 1.;
509  }
510  }
511  return a/double(k);
512 }
513 
514 double skymap::norm(double a)
515 {
516  size_t i,j,m;
517  size_t n = value.size();
518  double s = 0;
519  double x;
520 
521  (*this) += 1.-this->max();
522  if(a==0.) return a;
523 
524  for(i=0; i<n; i++) {
525  m = value[i].size();
526  for(j=0; j<m; j++) {
527  x = exp((this->value[i][j]-1.)*fabs(a));
528  s += x;
529  if(a>0.) this->value[i][j] = x;
530  }
531  }
532  if(a>0.) (*this) *= 1./s;
533  return s;
534 }
535 
537 {
538  size_t i,j;
539  size_t n = value.size();
540  size_t m = this->size();
541  size_t l = 0;
542 
543  if(index.size() != m) index.resize(m);
544  index = 0;
545 
546  for(i=0; i<n; i++) {
547  m = value[i].size();
548  for(j=0; j<m; j++) {
549  if(i&1) index.data[l] = 1;
550  if(j&1 && k==4) index.data[l] = 1;
551  l++;
552  }
553  }
554  return;
555 }
556 
557 #ifdef _USE_HEALPIX
558 //: dump skymap into fits file
559 void skymap::Dump2fits(const char* file, double gps_obs, const char* configur,
560  const char* TTYPE1, const char* TUNIT1, char coordsys)
561 {
562 
563  if(healpix!=NULL) {
564  int L = healpix->Npix();
565  for(int l=0;l<L;l++) (*healpix)[l]=this->get(l);
566  try {
567 
568  TString fName = file;
569  if(fName.EndsWith(".gz")) fName.Resize(fName.Sizeof()-4);
570 
571  if(!fName.EndsWith(".fits")) {
572  cout << "skymap::Dump2fits Error : wrong file extension" << endl;
573  cout << "fits file must ends with '.fits' or '.fits.gz'" << endl;
574  exit(1);
575  }
576 
577  // if already exist the delete file
578  Long_t id,size,flags,mt;
579  int estat = gSystem->GetPathInfo(fName.Data(),&id,&size,&flags,&mt);
580  if (estat==0) {
581  char cmd[1024];
582  sprintf(cmd,"rm %s",fName.Data());
583  cout << cmd << endl;
584  int err=gSystem->Exec(cmd);
585  if(err) {
586  cout << "skymap::Dump2fits Error : failed to remove file" << endl;
587  exit(1);
588  }
589  }
590 
591  // write_Healpix_map_to_fits(file,*healpix,PLANCK_FLOAT32);
592  fitshandle out;
593  out.create(fName.Data());
594  PDT datatype = PLANCK_FLOAT32;
595 
596  // prepare_Healpix_fitsmap (out, *healpix, datatype, colname);
597  arr<string> colname(1);
598  colname[0] = (TString(TTYPE1)!="") ? TTYPE1 : "unknown ";
599  string tunit1 = (TString(TUNIT1)!="") ? TUNIT1 : "unknown ";
600  vector<fitscolumn> cols;
601  int repcount = healpix_repcount (healpix->Npix());
602  for (tsize m=0; m<colname.size(); ++m)
603  cols.push_back (fitscolumn (colname[m],tunit1,repcount, datatype));
604  out.insert_bintab(cols);
605  out.set_key ("PIXTYPE",string("HEALPIX"),"HEALPIX pixelisation");
606  string ordering = (healpix->Scheme()==RING) ? "RING" : "NESTED";
607  out.set_key ("ORDERING",ordering, "Pixel ordering scheme, either RING or NESTED");
608  out.set_key ("NSIDE",healpix->Nside(),"Resolution parameter for HEALPIX");
609  out.set_key ("FIRSTPIX",0,"First pixel # (0 based)");
610  out.set_key ("LASTPIX",healpix->Npix()-1,"Last pixel # (0 based)");
611  out.set_key ("INDXSCHM",string("IMPLICIT"), "Indexing: IMPLICIT or EXPLICIT");
612  if(coordsys=='c' || coordsys=='C')
613  out.set_key ("COORDSYS",string("C "),"Pixelisation coordinate system");
614  out.set_key ("CREATOR",string("CWB "),"Program that created this file");
615  if(TString(configur)!="")
616  out.set_key ("CONFIGUR",string(configur),"software configuration used to process the data");
617 
618  // convert gps_obs to (YYYY-MM-DDThh:mm:ss.xx UT)
619  if(gps_obs>0) {
620  wat::Time date(gps_obs);
621  TString sdate = date.GetDateString();
622  sdate.Resize(19);
623  sdate.ReplaceAll(" ","T");
624  char date_obs[32];
625  sprintf(date_obs,"%s.%d",sdate.Data(),int(100*(date.GetNSec()/1000000000.)));
626  out.set_key ("DATE-OBS",string(date_obs),"UTC date of the observation");
627  char mjd_obs[32];
628  sprintf(mjd_obs,"%.9f",date.GetModJulianDate());
629  out.set_key ("MJD-OBS",string(mjd_obs),"modified Julian date of the observation");
630  }
631 
632  // UTC date of file creation
633  wat::Time date("now");
634  TString sdate = date.GetDateString();
635  sdate.Resize(19);
636  sdate.ReplaceAll(" ","T");
637  char date_obs[32];
638  sprintf(date_obs,"%s.%d",sdate.Data(),int(100*(date.GetNSec()/1000000000.)));
639  out.set_key ("DATE",string(date_obs),"UTC date of file creation");
640 
641  out.write_column(1,healpix->Map());
642 
643  out.close();
644 
645  // if file ends with .gz then the file zipped
646  if(TString(file).EndsWith(".gz")) {
647  char cmd[1024];
648  sprintf(cmd,"gzip -f %s",fName.Data());
649  cout << cmd << endl;
650  int err=gSystem->Exec(cmd);
651  if(err) {
652  cout << "skymap::Dump2fits Error : gzip error" << endl;
653  exit(1);
654  }
655  }
656 
657  } catch(...) {}
658  return;
659  } else {
660  cout << "skymap::Dump2fits Error : healpix not initialized" << endl;
661  exit(1);
662  }
663 }
664 #endif
665 
666 //: dump skymap into root file
667 void skymap::DumpObject(char* file)
668 {
669  TFile* rfile = new TFile(file, "RECREATE");
670  this->Write("skymap"); // write this object
671  rfile->Close();
672 }
673 
674 //: dump skymap into binary file
675 void skymap::DumpBinary(char* file, int mode)
676 {
677  size_t i,j,m;
678  size_t n = value.size();
679  size_t k = 0;
680 
681  for(i=0; i<n; i++) k += value[i].size();
682  wavearray<float> x(k);
683 
684  k = 0;
685  for(i=0; i<n; i++) {
686  m = value[i].size();
687  for(j=0; j<m; j++) {
688  x.data[k++] = this->value[i][j];
689  }
690  }
691 
692  x.DumpBinary(file,mode);
693  return;
694 }
695 
696 //: get skymap value at index l (access as a linear array)
697 // and fill in mTheta, mPhi and mIndex fields
698 //!param: sky index
699 double skymap::get(size_t l) {
700  size_t i,m;
701  size_t n = value.size();
702  size_t k = 0;
703 
704  mIndex = l;
705  for(i=0; i<n; i++) {
706  m = value[i].size();
707  k += m;
708  if(k <= l) continue;
709  mTheta = i;
710  mPhi = m-int(k-l);
711  return this->value[mTheta][mPhi];
712  }
713  mTheta = mPhi = mIndex = -1;
714  return 0.;
715 }
716 
717 //: get sky index at theta,phi
718 //!param: theta
719 //!param: phi
720 size_t skymap::getSkyIndex(double th, double ph) {
721  size_t k;
722 
723  if(healpix!=NULL) {
724 #ifdef _USE_HEALPIX
725  pointing P(th*deg2rad,ph*deg2rad);
726  k = healpix->ang2pix(P);
727  this->get(k);
728  return mIndex;
729 #endif
730  }
731 
732  double g;
733  g = (value.size()-1)/(theta_2-theta_1);
734  if(g>0.) {
735  if(th<theta_1) th = theta_1;
736  if(th>theta_2) th = theta_2;
737  mTheta = int(g*(th-theta_1)+0.5);
738  }
739  else { mTheta = 0; }
740 
741  g = value[mTheta].size()/(phi_2-phi_1);
742  if(ph< phi_1) ph = phi_2 - 0.5/g;
743  if(ph>=phi_2) ph = phi_1 + 0.5/g;
744  mPhi = int(g*(ph-phi_1));
745 
746  k = mPhi;
747  for(int i=0; i<mTheta; i++) k += value[i].size();
748  mIndex = k;
749  return k;
750 }
751 
752 #ifdef _USE_HEALPIX
753 //:
754 wavearray<int> skymap::neighbors(int index)
755 {
756 /*! Returns the neighboring pixels of pixel index in neighbors.
757  On exit, neighbors contains (in this order)
758  the pixel numbers of the SW, W, NW, N, NE, E, SE and S neighbor
759  of the pixel index. If a neighbor does not exist (this can only be the case
760  for the W, N, E and S neighbors), its entry is set to -1.
761 */
762 
763  wavearray<int> neighbors(8);
764 
765  if(healpix!=NULL) {
766  fix_arr< int, 8 > result;
767  healpix->neighbors (index, result);
768  for(int k=0;k<8;k++) neighbors[k]=result[k];
769  }
770 
771  return neighbors;
772 }
773 #endif
774 
775 #ifdef _USE_HEALPIX
776 void skymap::median(double radius) {
777 //
778 // Applies a median filter with the desired radius
779 //
780 // Input: radius - the radius (in degrees) of the disc
781 
782  if(healpix==NULL) {
783  cout << "skymap::median Error : healpix not initialized" << endl;
784  exit(1);
785  }
786 
787  radius*=deg2rad;
788 
789  // fill healpix map with skymap values
790  int L = healpix->Npix();
791  for(int l=0;l<L;l++) (*healpix)[l]=this->get(l);
792 
793  // This method works in both RING and NEST schemes, but is
794  // considerably faster in the RING scheme.
795  if (healpix->Scheme()==NEST) healpix->swap_scheme();
796 
797  vector<int> list;
798  vector<float> list2;
799  for (int l=0; l<L; ++l) {
800  // Returns the numbers of all pixels that lie at least partially within
801  // radius of healpix->pix2ang(l) in list. It may also return a few pixels
802  // which do not lie in the disk at all
803  healpix->query_disc(healpix->pix2ang(l),radius,list);
804  list2.resize(list.size());
805  for (tsize i=0; i<list.size(); ++i) list2[i] = (*healpix)[list[i]];
806  double x=::median(list2.begin(),list2.end());
807  // fill skymap with healpix map values
808  this->set(l,x);
809  }
810 
811  return;
812 }
813 #endif
814 
815 #ifdef _USE_HEALPIX
816 void skymap::smoothing(double fwhm, int nlmax, int num_iter) {
817 //
818 // Applies a convolution with a Gaussian beam with an FWHM of fwhm_arcmin arcmin to alm.
819 //
820 // Input: fwhm - Gaussian beam with an FWHM (degrees)
821 // note : If fwhm<0, a deconvolution with -fwhm is performed
822 // nlmax - l max of spherical harmonic coefficients.
823 // num_iter - the number of iterations (default=0 is identical to map2alm())
824 
825  if (fwhm<0)
826  cout << "NOTE: negative FWHM supplied, doing a deconvolution..." << endl;
827 
828  if(healpix==NULL) {
829  cout << "skymap::smoothing Error : healpix not initialized" << endl;
830  exit(1);
831  }
832 
833  int L = healpix->Npix();
834 
835  // fill healpix map with skymap values
836  for(int l=0;l<L;l++) (*healpix)[l]=this->get(l);
837  double avg=healpix->average();
838  healpix->Add(double(-avg));
839 
840  // get alm coefficients
841  wat::Alm alm = getAlm(nlmax, num_iter);
842  // applies gaussian smoothing to alm
843  alm.smoothWithGauss(fwhm);
844  // converts a set of alm to a HEALPix map.
845  setAlm(alm);
846 
847  // fill skymap with healpix map values
848  healpix->Add(double(avg));
849  for(int l=0;l<L;l++) this->set(l,(*healpix)[l]);
850 
851  return;
852 }
853 #endif
854 
855 #ifdef _USE_HEALPIX
856 void skymap::rotate(double psi, double theta, double phi, int nlmax, int num_iter) {
857 //
858 // Rotates alm through the Euler angles psi, theta and phi.
859 // The Euler angle convention is right handed, rotations are active.
860 // - psi is the first rotation about the z-axis (vertical)
861 // - then theta about the ORIGINAL (unrotated) y-axis
862 // - then phi about the ORIGINAL (unrotated) z-axis (vertical)
863 // relates Alm */
864 
865  if(healpix==NULL) {
866  cout << "skymap::smoothing Error : healpix not initialized" << endl;
867  exit(1);
868  }
869 
870  // get alm coefficients
871  wat::Alm alm = getAlm(nlmax, num_iter);
872  // applies gaussian smoothing to alm
873  alm.rotate(psi, theta, phi);
874  // converts a set of alm to a HEALPix map.
875  setAlm(alm);
876 
877  return;
878 }
879 #endif
880 
881 #ifdef _USE_HEALPIX
882 void skymap::setAlm(wat::Alm alm) {
883 //
884 // Converts a set of a_lm to a HEALPix map.
885 //
886 // Input : alm - spherical harmonic coefficients.
887 
888  if(healpix==NULL) {
889  cout << "skymap::setAlm Error : healpix not initialized" << endl;
890  exit(1);
891  }
892 
893  // comvert wat::alm to alm
894  Alm<xcomplex<double> > _alm(alm.Lmax(),alm.Mmax());
895  for(int l=0;l<=_alm.Lmax();l++) {
896  int limit = TMath::Min(l,_alm.Mmax());
897  for (int m=0; m<=limit; m++)
898  _alm(l,m)=complex<double>(alm(l,m).real(),alm(l,m).imag());
899  }
900  // converts a set of alm to a HEALPix map.
901  alm2map(_alm,*healpix);
902 
903  // fill skymap with healpix map values
904  int L = healpix->Npix();
905  for(int l=0;l<L;l++) this->set(l,(*healpix)[l]);
906 
907  return;
908 }
909 #endif
910 
911 #ifdef _USE_HEALPIX
912 wat::Alm
913 skymap::getAlm(int nlmax, int num_iter) {
914 //
915 // Converts a Healpix map to a set of a_lms, using an iterative scheme
916 // which is more accurate than plain map2alm().
917 //
918 // Input: nlmax - l max of spherical harmonic coefficients.
919 // num_iter - the number of iterations (default=0 is identical to map2alm())
920 
921 
922  if(healpix==NULL) {
923  cout << "skymap::getAlm Error : healpix not initialized" << endl;
924  exit(1);
925  }
926 
927  if(nlmax<0) nlmax=0;
928 
929  // fill healpix map with skymap values
930  int L = healpix->Npix();
931  for(int l=0;l<L;l++) (*healpix)[l]=this->get(l);
932 
933  // create weight array and set to 1
934  // weight array containing the weights for the individual rings of the map
935  arr<double> weight(2*healpix->Nside());
936  for(int i=0;i<2*healpix->Nside();i++) weight[i]=1.0;
937 
938  // get alm coefficients
939  Alm<xcomplex<double> > alm(nlmax,nlmax);
940  if (healpix->Scheme()==NEST) healpix->swap_scheme();
941  // Converts a Healpix map to a set of a_lms, using an iterative scheme
942  map2alm_iter(*healpix,alm,num_iter,weight);
943 
944  // comverts alm to wat::alm
945  wat::Alm _alm = alm;
946 
947  return _alm;
948 }
949 #endif
950 
951 #ifdef _USE_HEALPIX
952 void
953 skymap::resample(int order) {
954 //
955 // Resample the Healpix map to the healpix order
956 // the resampling is done using the spherical harmonic coefficients
957 //
958 // Input: order - healpix order of the resampled skymap
959 
960  if(healpix==NULL) {
961  cout << "skymap::resample Error : healpix not initialized" << endl;
962  exit(1);
963  }
964 
965  if(order == getOrder()) return; // nothing to do
966 
967  // get alm coefficients
968  int nlmax = 2*healpix->Nside();
969  wat::Alm alm = getAlm(nlmax);
970 
971  // skymap with new healpix order
972  *this = skymap(order);
973 
974  int _nlmax = 2*healpix->Nside();
975  wat::Alm _alm(_nlmax,_nlmax);
976 
977  // fill new _alm with original alm
978  int min_nlmax = TMath::Min(nlmax,_nlmax);
979  for(int l=0;l<=min_nlmax;l++) {
980  int limit = TMath::Min(l,alm.Mmax());
981  for (int m=0; m<=limit; m++) _alm(l,m)=alm(l,m);
982  }
983 
984  // set new _alm coefficients
985  setAlm(_alm);
986 
987  return;
988 }
989 #endif
990 
991 #ifdef _USE_HEALPIX
992 int
993 skymap::getRings() {
994 //
995 // return the number of division in theta
996 
997  if(healpix==NULL) {
998  cout << "skymap::getRings Error : healpix not initialized" << endl;
999  exit(1);
1000  }
1001 
1002  int nrings = 1; // divisions on theta
1003  for(int i=0;i<=healpix->Order();i++) nrings=2*nrings+1;
1004 
1005  return nrings;
1006 }
1007 #endif
1008 
1009 #ifdef _USE_HEALPIX
1010 int
1011 skymap::getRingPixels(int ring) {
1012 //
1013 // return the number of pixels in the ring
1014 
1015  if(healpix==NULL) {
1016  cout << "skymap::getRingPixels Error : healpix not initialized" << endl;
1017  exit(1);
1018  }
1019 
1020  if(ring<1 || ring>getRings()) {
1021  cout << "skymap::getRingPixels Error : ring index " << ring
1022  << " not allowed [1:" << getRings() << "]" << endl;
1023  exit(1);
1024  }
1025 
1026  int startpix; int ringpix; double costheta; double sintheta; bool shifted=false;
1027  healpix->get_ring_info(ring, startpix, ringpix, costheta, sintheta, shifted);
1028  return ringpix;
1029 }
1030 #endif
1031 
1032 #ifdef _USE_HEALPIX
1033 int
1034 skymap::getStartRingPixel(int ring) {
1035 //
1036 // return the start pixel index in the ring
1037 
1038  if(healpix==NULL) {
1039  cout << "skymap::getStartRingPixel Error : healpix not initialized" << endl;
1040  exit(1);
1041  }
1042 
1043  if(ring<1 || ring>getRings()) {
1044  cout << "skymap::getRingPixels Error : ring index " << ring
1045  << " not allowed [1:" << getRings() << "]" << endl;
1046  exit(1);
1047  }
1048 
1049  int startpix; int ringpix; double costheta; double sintheta; bool shifted=false;
1050  healpix->get_ring_info(ring, startpix, ringpix, costheta, sintheta, shifted);
1051  return startpix;
1052 }
1053 #endif
1054 
1055 #ifdef _USE_HEALPIX
1056 int
1057 skymap::getEulerCharacteristic(double threshold) {
1058 //
1059 // return the get the Euler characteristic for pixels value > threshold
1060 
1061  if(healpix==NULL) {
1062  cout << "skymap::getEulerCharacteristic Error : healpix not initialized" << endl;
1063  exit(1);
1064  }
1065 
1066 
1067  int nV=0; // number of vertices (= number of pixels)
1068  int nE=0; // number of edges
1069  int nF=0; // number of faces
1070 
1071  int F1=0;
1072  int F2=0;
1073  int F3=0;
1074  int F4=0;
1075 
1077 
1078  double a;
1079  for(int l=0;l<size();l++) {
1080  a=get(l);
1081  if(a>threshold){
1082  nV++;
1083  index = neighbors(l);
1084  a=get(index[0]);
1085  if(a>threshold) {nE++;F1++;F4++;}
1086  a=get(index[1]);
1087  if(a>threshold) {F1++;}
1088  a=get(index[2]);
1089  if(a>threshold) {nE++;F1++;F2++;}
1090  a=get(index[3]);
1091  if(a>threshold) {F2++;}
1092  a=get(index[4]);
1093  if(a>threshold) {nE++;F2++;F3++;}
1094  a=get(index[5]);
1095  if(a>threshold) {F3++;}
1096  a=get(index[6]);
1097  if(a>threshold) {nE++;F3++;F4++;}
1098  a=get(index[7]);
1099  if(a>threshold) {F4++;}
1100  nF+=F1/3+F2/3+F3/3+F4/3;
1101  F1=0; F2=0; F3=0; F4=0;
1102  }
1103  }
1104  nE/=2;nF/=4;
1105 
1106  return nV-nE+nF;
1107 }
1108 #endif
1109 
1111 
1112  TString name = fname;
1113  name.ReplaceAll(" ","");
1114  TObjArray* token = TString(fname).Tokenize(TString("."));
1115  TObjString* ext_tok = (TObjString*)token->At(token->GetEntries()-1);
1116  TString ext = ext_tok->GetString();
1117  if(ext=="dat") {
1118  //: dump skymap into binary file
1119  DumpBinary(fname,0);
1120  } else
1121 #ifdef _USE_HEALPIX
1122  if(ext=="fits") {
1123  //: dump skymap into fits file
1124  Dump2fits(fname);
1125  } else
1126 #endif
1127  if(ext=="root") {
1128  //: dump skymap object into root file
1129  DumpObject(fname);
1130  } else {
1131  cout << "skymap operator (>>) : file type " << ext.Data() << " not supported" << endl;
1132  }
1133  return fname;
1134 }
1135 
1136 //______________________________________________________________________________
1137 void skymap::Streamer(TBuffer &R__b)
1138 {
1139  // Stream an object of class skymap.
1140 
1141  UInt_t R__s, R__c;
1142  if (R__b.IsReading()) {
1143  Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
1144  TNamed::Streamer(R__b);
1145  {
1146  vector<vectorD> &R__stl = value;
1147  R__stl.clear();
1148  TClass *R__tcl1 = TBuffer::GetClass(typeid(vector<double,allocator<double> >));
1149  if (R__tcl1==0) {
1150  Error("value streamer","Missing the TClass object for vector<double,allocator<double> >!");
1151  return;
1152  }
1153  int R__i, R__n;
1154  R__b >> R__n;
1155  R__stl.reserve(R__n);
1156  for (R__i = 0; R__i < R__n; R__i++) {
1157  vector<double,allocator<double> > R__t;
1158  R__b.StreamObject(&R__t,R__tcl1);
1159  R__stl.push_back(R__t);
1160  }
1161  }
1162  {
1163  vector<int> &R__stl = index;
1164  R__stl.clear();
1165  int R__i, R__n;
1166  R__b >> R__n;
1167  R__stl.reserve(R__n);
1168  for (R__i = 0; R__i < R__n; R__i++) {
1169  int R__t;
1170  R__b >> R__t;
1171  R__stl.push_back(R__t);
1172  }
1173  }
1174  if(R__v > 2) R__b >> sms;
1175  R__b >> theta_1;
1176  R__b >> theta_2;
1177  R__b >> phi_1;
1178  R__b >> phi_2;
1179  R__b >> gps;
1180  R__b >> mTheta;
1181  R__b >> mPhi;
1182  R__b >> mIndex;
1183  if(R__v > 1) {
1184  R__b >> healpix_order;
1185  if(healpix_order > 0) {
1186 #ifdef _USE_HEALPIX
1187  healpix = new Healpix_Map<double>(healpix_order,RING);
1188 #endif
1189  } else {
1190  healpix = NULL;
1191  healpix_order = 0;
1192  }
1193  }
1194  R__b.CheckByteCount(R__s, R__c, skymap::IsA());
1195  } else {
1196  R__c = R__b.WriteVersion(skymap::IsA(), kTRUE);
1197  TNamed::Streamer(R__b);
1198  {
1199  vector<vectorD> &R__stl = value;
1200  int R__n=(&R__stl) ? int(R__stl.size()) : 0;
1201  R__b << R__n;
1202  if(R__n) {
1203  TClass *R__tcl1 = TBuffer::GetClass(typeid(vector<double,allocator<double> >));
1204  if (R__tcl1==0) {
1205  Error("value streamer","Missing the TClass object for vector<double,allocator<double> >!");
1206  return;
1207  }
1208  vector<vectorD>::iterator R__k;
1209  for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1210  R__b.StreamObject((vector<double,allocator<double> >*)&(*R__k),R__tcl1);
1211  }
1212  }
1213  }
1214  {
1215  vector<int> &R__stl = index;
1216  int R__n=(&R__stl) ? int(R__stl.size()) : 0;
1217  R__b << R__n;
1218  if(R__n) {
1219  vector<int>::iterator R__k;
1220  for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1221  R__b << (*R__k);
1222  }
1223  }
1224  }
1225  R__b << sms;
1226  R__b << theta_1;
1227  R__b << theta_2;
1228  R__b << phi_1;
1229  R__b << phi_2;
1230  R__b << gps;
1231  R__b << mTheta;
1232  R__b << mPhi;
1233  R__b << mIndex;
1234  R__b << healpix_order;
1235  R__b.SetByteCount(R__c, kTRUE);
1236  }
1237 }
1238 
wavearray< double > t(hp.size())
void smoothWithGauss(double fwhm)
Definition: alm.hh:207
#define F1
Definition: Regression_H1.C:13
int Mmax() const
Returns the maximum m.
Definition: alm.hh:72
double norm(double=0.)
Definition: skymap.cc:514
static double g(double e)
Definition: GNGen.cc:116
TString GetDateString()
Definition: time.cc:461
TH1 * t1
par [0] value
int mPhi
Definition: skymap.hh:327
wavearray< double > a(hp.size())
int mTheta
Definition: skymap.hh:326
WSeries< float > v[nIFO]
Definition: cwb_net.C:80
skymap & operator-=(const skymap &)
Definition: skymap.cc:317
par [0] name
int n
Definition: cwb_net.C:28
TString("c")
double gps
Definition: skymap.hh:325
ofstream out
Definition: cwb_merge.C:214
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\ layers : "<< nLAYERS<< "\ dF(hz) : "<< dF<< "\ dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1) *itime+ifreq;double time=itime *dT;double freq=(ifreq >0) ? ifreq *dF :dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
INT_4U GetNSec()
Definition: time.cc:654
float theta
double median(std::vector< double > &vec)
Definition: wavegraph.cc:485
TH2F * ph
skymap & operator*=(const skymap &)
Definition: skymap.cc:346
virtual ~skymap()
Definition: skymap.cc:250
double min()
Definition: skymap.cc:459
Long_t flags
STL namespace.
std::vector< int > index
Definition: skymap.hh:319
double theta_1
Definition: skymap.hh:321
Long_t size
int m
Definition: cwb_net.C:28
double max()
Definition: skymap.cc:438
void downsample(wavearray< short > &, size_t=4)
Definition: skymap.cc:536
int j
Definition: cwb_net.C:28
i drho i
TList * list
double theta_2
Definition: skymap.hh:322
double mean()
Definition: skymap.cc:480
skymap()
Definition: skymap.cc:47
#define PI
Definition: watfun.hh:32
size_t mode
watplot p2(const_cast< char *>("TFMap2"))
size_t getSkyIndex(double th, double ph)
param: theta param: phi
Definition: skymap.cc:720
virtual size_t size() const
Definition: wavearray.hh:145
pTF [i] DumpBinary(file)
tlive_fix Write()
float phi
std::vector< vectorD > value
Definition: skymap.hh:318
float psi
int Lmax() const
Returns the maximum l.
Definition: alm.hh:70
#define F2
Definition: Regression_H1.C:14
double deg2rad
double fraction(double)
Definition: skymap.cc:497
skymap & operator+=(const skymap &)
Definition: skymap.cc:288
i() int(T_cor *100))
void DumpObject(char *)
Definition: skymap.cc:667
watplot p1(const_cast< char *>("TFMap1"))
int l
double sms
Definition: skymap.hh:320
char fname[1024]
char * operator>>(char *fname)
Definition: skymap.cc:1110
int mIndex
Definition: skymap.hh:328
int k
TObjArray * token
Definition: skymap.hh:63
void DumpBinary(char *, int)
Definition: skymap.cc:675
TFile * ifile
int npix
virtual void DumpBinary(const char *, int=0)
Definition: wavearray.cc:353
Definition: alm.hh:193
s s
Definition: cwb_net.C:155
double rad2deg
double gps
void rotate(double psi, double theta, double phi)
Definition: alm.hh:231
wavearray< int > index
double fabs(const Complex &x)
Definition: numpy.cc:55
skymap & operator=(const skymap &)
Definition: skymap.cc:256
double GetModJulianDate()
Definition: time.cc:591
char cmd[1024]
int estat
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
Long_t mt
double phi_2
Definition: skymap.hh:324
DataType_t * data
Definition: wavearray.hh:319
Long_t id
double phi_1
Definition: skymap.hh:323
double get(size_t i)
param: sky index
Definition: skymap.cc:699
bool * healpix
Definition: skymap.hh:337
char fName[256]
virtual void resize(unsigned int)
Definition: wavearray.cc:463
skymap & operator/=(const skymap &)
Definition: skymap.cc:373
exit(0)
size_t healpix