Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
gskymap.cc
Go to the documentation of this file.
1 
2 #include "gskymap.hh"
3 #include "TPaletteAxis.h"
4 
5 
6 ClassImp(gskymap)
7 
8 //______________________________________________________________________________
9 /* Begin_Html
10 <center><h2>gskymap class</h2></center>
11 This class gskymap is derived from the skymap class
12  It is used to produce the skymap plots.
13 
14 <p>
15 
16 <h3><a name="usage">Usage</a></h3>
17 
18 create gskymap with HEALPix order=7
19 <pre>
20  <a class="funcname" href="#gskymap:gskymap">gskymap</a> gSM((int)7);
21 </pre>
22 set gskymap options
23 <pre>
24  gSM.<a class="funcname" href="#gskymap:SetOptions">SetOptions</a>("cartesian","geographic",1);
25 </pre>
26 set title
27 <pre>
28  gSM.<a class="funcname" href="#gskymap:SetTitle">SetTitle</a>("Projection : cartesian - Coordinates : geographic");
29 </pre>
30 set world map
31 <pre>
32  gSM.<a class="funcname" href="#gskymap:SetWorldMap">SetWorldMap</a>();
33 </pre>
34 draw skymap
35 <pre>
36  gSM.<a class="funcname" href="#gskymap:Draw">Draw</a>(0,"col");
37 </pre>
38 
39 <p>
40 <h3><a name="example">Example</a></h3>
41 <p>
42 The macro <a href="./tutorials/gwat/DrawGskymap.C.html">DrawGskymap.C</a> is an example which shown how to use the gskymap class.<br>
43 The pictures below gives the macro output plots.<br>
44 <p>
45 
46 End_Html
47 
48 Begin_Macro
49 DrawGskymap.C("cartesian")
50 End_Macro
51 Begin_Macro
52 DrawGskymap.C("hammer")
53 End_Macro
54 Begin_Macro
55 DrawGskymap.C("parabolic")
56 End_Macro
57 Begin_Macro
58 DrawGskymap.C("sinusoidal")
59 End_Macro */
60 
61 
62 using namespace std;
63 
64 //______________________________________________________________________________
65 void
66 gskymap::SetOptions(TString projection, TString coordinate, double resolution, bool goff) {
67 //
68 // Set graphical options
69 //
70 // Input: projection - set projection type
71 // - "hammer" (default)
72 // - "sinusoidal"
73 // - "parabolic"
74 // - "cartesian" or ""
75 //
76 // coordinate - set the coordinate type
77 // - "geographic" (default)
78 // - "celestial"
79 // - "cwb"
80 //
81 // resolution - factor used to increase the resolution of the skymap (default=1)
82 // phi = 2*360*resolution points, theta = 2*180*resolution
83 //
84 // goff - true -> turn off the graphical output (default = false)
85 //
86 
87  coordinate.ToUpper();
88 
89  if (resolution<=0) {cout << "gskymap::gskymap: Error sky resolution must be >0" << endl;exit(1);}
90 
91  if(gSystem->Getenv("CWB_GWAT")!=NULL) {
92  worldMapPath=TString(gSystem->Getenv("CWB_GWAT"))+"/data";
93  } else {
94  worldMapPath=".";
95  }
96 
97  // generate a unique base name for canvas and hist2d
98  gRandom->SetSeed(0);
99  char sname[64];sprintf(sname,"gskymap_%d",int(gRandom->Rndm(13)*1.e9));
100  name = sname;
101  this->projection = projection;
102  this->coordinate = coordinate;
103  this->resolution = resolution;
104  this->goff = goff;
105 
106  wtopx=35; wtopy=46; ww=817; wh=472;
107 
108  double phMin = -180;
109  double phMax = 180;
110  double thMin = -90;
111  double thMax = 90;
112  if (coordinate.CompareTo("CWB")==0) {phMin=0;phMax=360;thMin=0;thMax=180;}
113  if (coordinate.CompareTo("CELESTIAL")==0) {phMin=0;phMax=360;thMin=-90;thMax=90;}
114 
115  if(h2!=NULL) delete h2;
116  char nameh[64];sprintf(nameh,"h_%s",name.Data());
117  h2 = new TH2D(nameh,"gskymap", int(360*resolution), phMin, phMax, int(180*resolution), thMin, thMax);
118  h2->SetStats(kFALSE);
119  h2->SetTitleFont(12);
120  h2->SetTitle(title);
121  h2->SetFillColor(kWhite);
122 
123  h2->GetXaxis()->SetNdivisions(-604);
124  h2->GetXaxis()->SetLabelFont(42);
125  h2->GetXaxis()->SetLabelOffset(0.014);
126  h2->GetXaxis()->SetTitleOffset(1.2);
127  h2->GetYaxis()->SetTitleOffset(0.8);
128  h2->GetYaxis()->SetNdivisions(-602);
129  h2->GetYaxis()->SetLabelFont(42);
130  h2->GetYaxis()->SetLabelOffset(0.01);
131  h2->GetZaxis()->SetLabelFont(42);
132  h2->GetZaxis()->SetNoExponent(false);
133 
134  h2->GetXaxis()->SetTitleFont(42);
135  if ((coordinate.CompareTo("CWB")==0)||(coordinate.CompareTo("GEOGRAPHIC")==0)) {
136  h2->GetXaxis()->SetTitle("#phi (deg)");
137  } else if (coordinate.CompareTo("CELESTIAL")==0) {
138  h2->GetXaxis()->SetTitle("RA (deg)");
139  } else {
140  {cout << "gskymap::gskymap: Error : coordinate system not declared" << endl;exit(1);}
141  }
142  h2->GetXaxis()->CenterTitle(true);
143 
144  h2->GetYaxis()->SetTitleFont(42);
145  if ((coordinate.CompareTo("CWB")==0)||(coordinate.CompareTo("GEOGRAPHIC")==0)) {
146  h2->GetYaxis()->SetTitle("#theta (deg)");
147  } else if (coordinate.CompareTo("CELESTIAL")==0) {
148  h2->GetYaxis()->SetTitle("DEC (deg)");
149  } else {
150  {cout << "gskymap::gskymap: Error : coordinate system not declared" << endl;exit(1);}
151  }
152  h2->GetYaxis()->CenterTitle(true);
153 
154  h2->GetZaxis()->SetTitleOffset(0.6);
155  h2->GetZaxis()->SetTitleFont(42);
156  h2->GetZaxis()->SetTitle(zAxisTitle);
157  h2->GetZaxis()->CenterTitle(true);
158 
159  if(this->size()) FillData();
160 }
161 
162 //______________________________________________________________________________
163 void
165 //
166 // create canvas object used for skymap plot
167 //
168 
169 // TCanvas* tcanvas = (TCanvas*) gROOT->FindObject(name);
170 // if (tcanvas!=NULL) {cout << "gskymap::gskymap: Error Canvas Name already exist" << endl;exit(1);}
171  if(canvas!=NULL) delete canvas;
172  if(!goff) {
173  char namec[64];sprintf(namec,"c_%s",name.Data());
174  canvas = new TCanvas(namec, "gskymap", wtopx, wtopy, ww, wh);
175  canvas->Clear();
176  canvas->ToggleEventStatus();
177  canvas->SetGridx(false);
178  canvas->SetGridy(false);
179  canvas->SetLogz(isLogz);
180  canvas->SetFillColor(kWhite);
181  canvas->SetLeftMargin(0.08);
182  canvas->SetBottomMargin(0.13);
183  canvas->SetBorderMode(0);
184  //canvas->SetWindowSize(1200,600);
185  //canvas->SetGrayscale();
186  } else {
187  canvas=NULL;
188  }
189 
190  // remove the red box around canvas
191  gStyle->SetFrameBorderMode(0);
192  gROOT->ForceStyle();
193 
194  gStyle->SetTitleH(0.050);
195  gStyle->SetTitleW(0.95);
196  gStyle->SetTitleY(0.98);
197  gStyle->SetTitleFont(12,"D");
198  gStyle->SetTitleColor(kBlue,"D");
199  gStyle->SetTextFont(12);
200  gStyle->SetTitleFillColor(kWhite);
201  gStyle->SetLineColor(kWhite);
202  gStyle->SetNumberContours(256);
203  gStyle->SetMarkerStyle(7);
204  gStyle->SetMarkerSize(2);
205  gStyle->SetCanvasColor(kWhite);
206  gStyle->SetStatBorderSize(1);
207 }
208 
209 //______________________________________________________________________________
211 //
212 // destructor
213 //
214 
215  if(h2!=NULL) delete h2;
216  if(canvas!=NULL) delete canvas;
217 
218  if (wm_size) delete [] wm_lon;
219  if (wm_size) delete [] wm_lat;
220 
222  ClearWorldMap();
223  ClearGridx();
224  ClearGridy();
225  ClearAxisLabel();
226 }
227 
228 //______________________________________________________________________________
229 void
231 //
232 // fill the internal skymap reading from ascii file
233 //
234 // Input: fname - input file name
235 //
236 // format :
237 //
238 // header : sky_res -> area of a pixel (degrees)
239 // theta_1 -> theta min
240 // theta_2 -> theta max
241 // phi_1 -> phi min
242 // phi_2 -> phi max
243 // lenght -> number of entries
244 //
245 // list of data: sky_index value
246 //
247 
248  char iline[1024];
249  TObjArray* tok;
250 
251  ifstream in;
252  in.open(fname,ios::in);
253  if (!in.good()) {cout << "gskymap::FillData: Error Opening File : " << fname << endl;exit(1);}
254 
255  in.getline(iline,1024);
256  in.getline(iline,1024);
257  if(TString(iline).Contains("skymap")==0)
258  {cout << "gskymap::FillData: " << fname << " is not a skymap file" << endl;exit(1);}
259  in.getline(iline,1024);
260  tok = TString(iline).Tokenize(TString(' '));
261  TObjString* tsms = (TObjString*)tok->At(1);
262  if(tsms->GetString().IsAlpha()==1)
263  {cout << "gskymap::FillData: resolution is not well defined " << endl;exit(1);}
264  double sms = tsms->GetString().Atof();
265  in.getline(iline,1024);
266  tok = TString(iline).Tokenize(TString(' '));
267  TObjString* ttheta_1 = (TObjString*)tok->At(1);
268  if(ttheta_1->GetString().IsAlpha()==1)
269  {cout << "gskymap::FillData: theta_1 is not well defined " << endl;exit(1);}
270  double theta_1 = ttheta_1->GetString().Atof();
271  in.getline(iline,1024);
272  tok = TString(iline).Tokenize(TString(' '));
273  TObjString* ttheta_2 = (TObjString*)tok->At(1);
274  if(ttheta_2->GetString().IsAlpha()==1)
275  {cout << "gskymap::FillData: theta_2 is not well defined " << endl;exit(1);}
276  double theta_2 = ttheta_2->GetString().Atof();
277  in.getline(iline,1024);
278  tok = TString(iline).Tokenize(TString(' '));
279  TObjString* tphi_1 = (TObjString*)tok->At(1);
280  if(tphi_1->GetString().IsAlpha()==1)
281  {cout << "gskymap::FillData: phi_1 is not well defined " << endl;exit(1);}
282  double phi_1 = tphi_1->GetString().Atof();
283  in.getline(iline,1024);
284  tok = TString(iline).Tokenize(TString(' '));
285  TObjString* tphi_2 = (TObjString*)tok->At(1);
286  if(tphi_2->GetString().IsAlpha()==1)
287  {cout << "gskymap::FillData: phi_2 is not well defined " << endl;exit(1);}
288  double phi_2 = tphi_2->GetString().Atof();
289  in.getline(iline,1024);
290  tok = TString(iline).Tokenize(TString(' '));
291  TObjString* tlenght = (TObjString*)tok->At(1);
292  if(tlenght->GetString().IsAlpha()==1)
293  {cout << "gskymap::FillData: lenght is not well defined " << endl;exit(1);}
294  int lenght = tlenght->GetString().Atoi();
295 
296  //sms = double(int(1000*sms))/1000.;
297 
298  //cout << sms << " " << theta_1 << " " << theta_2 << " " << phi_1 << " " << phi_2 << endl;
299 
300  int L = skymap::size();
301  if(L!=lenght)
302  {cout << "gskymap::FillData: lenght is not consistent" << endl;exit(1);}
303 
304  while (1) {
305 
306  in.getline(iline,1024);
307  if (!in.good()) break;
308  //cout << ++cnt << " " << iline << endl;
309  tok = TString(iline).Tokenize(TString(' '));
310 
311  if (tok->GetEntries()==1) continue;
312 
313  TObjString* tindex = (TObjString*)tok->At(0);
314  TObjString* tvalue = (TObjString*)tok->At(1);
315 
316  if (tindex->GetString().IsAlpha()==1)
317  {cout << "gskymap::FillData: index " << tindex->GetString().Data() << " is not well defined " << endl;exit(1);}
318  //if (tvalue->GetString().IsAlpha()==1) continue;
319  // {cout << "gskymap::FillData: value " << tvalue->GetString().Data() << " is not well defined " << endl;exit(1);}
320 
321  int index = tindex->GetString().Atoi();
322  double value = tvalue->GetString().Atof();
323 
324  if(index>=0&&index<L)
325  skymap::set(index,value);
326  else
327  {cout << "gskymap::FillData: index " << index << " not allowed - max " << L-1 << endl;}
328 
329  }
330  in.close();
331 
332  FillData();
333 
334  return;
335 }
336 
337 //______________________________________________________________________________
338 void
340 //
341 // Fill 2D histo with the internal skymap data
342 //
343 
344  int size=int(2*180*resolution*2*360*resolution);
345 
346  double* phi = new double[size];
347  double* theta = new double[size];
348  double* binc = new double[size];
349 
350  int k=0;
351  for(int i=0;i<2*180*resolution;i++) {
352  for(int j=0;j<2*360*resolution;j++) {
353  phi[k]=(double)j/(2*resolution);
354  theta[k]=(double)i/(2*resolution);
355  binc[k]=skymap::get(theta[k],phi[k]);
356  k++;
357  }
358  }
359 
360  if (coordinate.CompareTo("GEOGRAPHIC")==0) {
361  for (int i=0;i<k;i++) CwbToGeographic(phi[i],theta[i],phi[i],theta[i]);
362  }
363  if (coordinate.CompareTo("CELESTIAL")==0) {
364  for (int i=0;i<k;i++) CwbToCelestial(phi[i],theta[i],phi[i],theta[i]);
365  }
366 
367  FillData(k, phi, theta, binc);
368 
369  delete [] phi;
370  delete [] theta;
371  delete [] binc;
372 }
373 
374 //______________________________________________________________________________
375 void
376 gskymap::FillData(int size, double* phi, double* theta, double* binc) {
377 //
378 // Fill 2D histo using input arrays
379 //
380 // Input: size - size of the input arrays
381 // phi - arrays of phi values
382 // theta - arrays of theta values
383 // binc - arrays of amplitudes
384 //
385 
386  for (int i=0;i<size;i++) {
387  double x,y;
388  if (coordinate.CompareTo("CWB")==0) {
389  if(phi[i]<0) phi[i]+=360.;
390  if(phi[i]>360) phi[i]-=360.;
391  if(theta[i]<0) theta[i]+=180.;
392  if(theta[i]>180) theta[i]-=180.;
393  }
394  if (coordinate.CompareTo("GEOGRAPHIC")==0) {
395  if(phi[i]<-180) phi[i]+=360.;
396  if(phi[i]>180) phi[i]-=360.;
397  if(theta[i]<-90) theta[i]+=180.;
398  if(theta[i]>90) theta[i]-=180.;
399  phi[i]+=180.;
400  theta[i]+=90.;
401  }
402  if (coordinate.CompareTo("CELESTIAL")==0) {
403  if(phi[i]<0) phi[i]+=360.;
404  if(phi[i]>360) phi[i]-=360.;
405  if(theta[i]<-90) theta[i]+=180.;
406  if(theta[i]>90) theta[i]-=180.;
407  theta[i]+=90.;
408  }
409  if (projection.CompareTo("hammer")==0) {
410  ProjectHammer((phi[i]-180), (theta[i]-90), x, y);
411  x+=180+1;
412  y+=90+1;
413  } else
414  if (projection.CompareTo("sinusoidal")==0) {
415  ProjectSinusoidal((phi[i]-180), (theta[i]-90), x, y);
416  x+=180+1;
417  y+=90+1;
418  } else
419  if (projection.CompareTo("parabolic")==0) {
420  ProjectParabolic((phi[i]-180), (theta[i]-90), x, y);
421  x+=180+1;
422  y+=90+1;
423  } else {
424  x=phi[i];
425  y=theta[i];
426  }
427  h2->SetBinContent(int(x*resolution)+1,int(y*resolution)+1,binc[i]);
428 
429  if (coordinate.CompareTo("GEOGRAPHIC")==0) {
430  CwbToGeographic(phi[i],theta[i],phi[i],theta[i]);
431  //phi[i]-=180.;
432  //theta[i]+=90.;
433  }
434  if (coordinate.CompareTo("CELESTIAL")==0) {
435  CwbToCelestial(phi[i],theta[i],phi[i],theta[i]);
436  }
437  }
438 }
439 
440 //______________________________________________________________________________
441 void
442 gskymap::Draw(int dpaletteId, Option_t* goption) {
443 //
444 // draw skymap
445 //
446 // Input: dpaletteId - color palette used for the skymap
447 // goption - graphical options
448 //
449 
450  if(goff) return;
451 
452  CreateCanvas();
453 
454  TGaxis::SetMaxDigits(3);
455 
456  if(changed) {FillData();changed=false;}
457 
458  if (dpaletteId==DUMMY_PALETTE_ID) {
459  if (paletteId!=0) {
461  } else {
462  gStyle->SetPalette(1,0);
463  }
464  } else {
465  if (dpaletteId!=0) {
466  SetPlotStyle(dpaletteId);
467  } else {
468  gStyle->SetPalette(1,0);
469  }
470  }
471 
472  canvas->cd();
473  canvas->SetLogz(isLogz);
474  h2->Draw(goption);
475  // change palette's width
476  canvas->Update();
477  TPaletteAxis *palette = (TPaletteAxis*)h2->GetListOfFunctions()->FindObject("palette");
478  if(palette) {
479  palette->SetX1NDC(0.91);
480  palette->SetX2NDC(0.933);
481  palette->SetTitleOffset(0.92);
482  palette->GetAxis()->SetTickSize(0.01);
483  canvas->Modified();
484  }
485 
486  Double_t rlonait[360*4],rlatait[360*4];
487  Double_t kfix;
488  Int_t kk,k;
489  double mdistance=8.; // min distance (degrees) between polyline segments (avoid solid line when saved)
490 
491 // draw lines of constant longitude
492  if (isGridy) {
493  ClearGridy();
494  for (kk=0;kk<=360;kk+=15){
495  kfix=kk;
496  int klen=0;
497  for (k=0;k<180*4;k++){
498  double x,y;
499  double phi = kfix;
500  double theta = k/4.;
501  if (projection.CompareTo("hammer")==0) {
502  ProjectHammer((phi-180), (theta-90), x, y);
503  } else
504  if (projection.CompareTo("sinusoidal")==0) {
505  ProjectSinusoidal((phi-180), (theta-90), x, y);
506  } else
507  if (projection.CompareTo("parabolic")==0) {
508  ProjectParabolic((phi-180), (theta-90), x, y);
509  } else {
510  x=phi-180;
511  y=theta-90;
512  }
513  if (coordinate.CompareTo("CWB")==0) {x+=180;y+=90;}
514  //if (coordinate.CompareTo("CELESTIAL")==0) {x+=180;}
515  if (coordinate.CompareTo("CELESTIAL")==0) {x+=180;x=360-x;} //TO BE CHECKED
516  rlonait[klen] = x;
517  rlatait[klen] = y;
518  double distance = klen>0 ? sqrt(pow(x-rlonait[klen-1],2)+pow(y-rlatait[klen-1],2)) : mdistance;
519  if(distance>=mdistance || k==(180*4-1)) klen++;
520  }
521  TPolyLine *pL = new TPolyLine(klen,rlonait,rlatait);
522  pL->SetLineColor(colorGridx);
523  if(kk==0 || kk==360) pL->SetLineStyle(1); else pL->SetLineStyle(3);
524  pL->Draw();
525  gridyL.push_back(pL);
526  }
527  }
528 
529 /// draw lines of constant latitude
530  if (isGridx) {
531  ClearGridx();
532  //for (kk=10;kk<180;kk+=10){
533  for (kk=15;kk<180;kk+=15){
534  kfix=kk;
535  int klen=0;
536  for (k=0;k<360*4;k++){
537  double x,y;
538  double phi = k/4.;
539  double theta = kfix;
540  if (projection.CompareTo("hammer")==0) {
541  ProjectHammer((phi-180), (theta-90), x, y);
542  } else
543  if (projection.CompareTo("sinusoidal")==0) {
544  ProjectSinusoidal((phi-180), (theta-90), x, y);
545  } else
546  if (projection.CompareTo("parabolic")==0) {
547  ProjectParabolic((phi-180), (theta-90), x, y);
548  } else {
549  x=phi-180;
550  y=theta-90;
551  }
552  if (coordinate.CompareTo("CWB")==0) {x+=180;y+=90;}
553  //if (coordinate.CompareTo("CELESTIAL")==0) {x+=180;}
554  if (coordinate.CompareTo("CELESTIAL")==0) {x+=180;x=360-x;} // TO BE CHECKED
555  rlonait[klen] = x;
556  rlatait[klen] = y;
557  double distance = klen>0 ? sqrt(pow(x-rlonait[klen-1],2)+pow(y-rlatait[klen-1],2)) : mdistance;
558  if(distance>=mdistance || k==(360*4-1)) klen++;
559  }
560  TPolyLine *pL = new TPolyLine(klen,rlonait,rlatait);
561  pL->SetLineColor(colorGridy);
562  if(kk==0 || kk==180) pL->SetLineStyle(1); else pL->SetLineStyle(3);
563  pL->Draw();
564  gridxL.push_back(pL);
565  }
566  }
567 
568 /// draw labels of constant latitude
569 /*
570  if ((projection.CompareTo("")!=0)&&(coordinate.CompareTo("GEOGRAPHIC")==0)) {
571  ClearAxisLabel();
572  double phi[100],theta[100];
573  int k=0;
574  //for (kk=0;kk<=360;kk+=40) {phi[k]=kk-180;theta[k]=0;k++;}
575  //for (kk=0;kk<=180;kk+=10) {phi[k]=-180;theta[k]=kk-90;k++;}
576  for (kk=0;kk<=360;kk+=15) {phi[k]=kk-180;theta[k]=0;k++;}
577  for (kk=0;kk<=180;kk+=15) {phi[k]=-180;theta[k]=kk-90;k++;}
578  char label[16];
579  for(kk=0;kk<k;kk++) {
580  if(theta[kk]!=0)
581  sprintf(label,"%d",(int)theta[kk]);
582  else
583  sprintf(label,"%d",(int)phi[kk]+180);
584  double x,y;
585  if (projection.CompareTo("hammer")==0) {
586  ProjectHammer(phi[kk], theta[kk], x, y);
587  } else
588  if (projection.CompareTo("sinusoidal")==0) {
589  ProjectSinusoidal(phi[kk], theta[kk], x, y);
590  } else
591  if (projection.CompareTo("parabolic")==0) {
592  ProjectParabolic(phi[kk], theta[kk], x, y);
593  } else {
594  x=phi[kk];y=theta[kk];
595  }
596  TText *tP = new TText(x,y,label);
597  tP->SetTextSize(0.04);
598  tP->SetTextColor(kBlack);
599  tP->SetTextFont(42);
600  if(theta[kk]>0) tP->SetTextAlign(31);
601  if(theta[kk]==0) tP->SetTextAlign(32);
602  if(theta[kk]<0) tP->SetTextAlign(33);
603  tP->Draw();
604  axisT.push_back(tP);
605  }
606  }
607 */
608 
609 /// draw line of galactic disk
610  if (gpsGalacticDisk>=0.) {
612 
613  double gphi = gpsGalacticDisk>0. ? skymap::RA2phi(0., gpsGalacticDisk) : 0.;
614 
615  int L=4*360;
616  wavearray<float> th(L);
617  wavearray<float> ph(L);
618  for (int l=0;l<L;l++) {
619  double phi,theta;
620  GalacticToEquatorial(l/4.,0.,phi,theta);
621  ph.data[l]=phi+gphi;
622  th.data[l]=theta;
623  ph.data[l]=fmod(ph.data[l],360);
624  }
625 
626  wavearray<double> x(L);
627  wavearray<double> y(L);
628  int P=0;
629  for (int l=0;l<L;l++) {
630  if((fabs(ph.data[l]-ph.data[l==0?L-1:l-1])<180)&&(l<L-1)) {
631  x.data[P]=ph.data[l]; y.data[P]=th.data[l];
632 
633  if (projection.CompareTo("hammer")==0) {
634  ProjectHammer((x.data[P]-180), (-y.data[P]), x.data[P], y.data[P]);
635  x.data[P]=x.data[P]+180; y.data[P]=-y.data[P];
636  } else if (projection.CompareTo("sinusoidal")==0) {
637  ProjectSinusoidal((x.data[P]-180), (-y.data[P]), x.data[P], y.data[P]);
638  x.data[P]=x.data[P]+180; y.data[P]=-y.data[P];
639  } else if (projection.CompareTo("parabolic")==0) {
640  ProjectParabolic((x.data[P]-180), (-y.data[P]), x.data[P], y.data[P]);
641  x.data[P]=x.data[P]+180; y.data[P]=-y.data[P];
642  } else {
643  }
644  if (coordinate.CompareTo("GEOGRAPHIC")==0)
645  {x.data[P]=x.data[P]-180;}
646  if (coordinate.CompareTo("CWB")==0)
647  {y.data[P]=90+y.data[P];}
648  if (coordinate.CompareTo("CELESTIAL")==0)
649  {x.data[P]=360-x.data[P];y.data[P]=90+y.data[P];}
650  P++;
651  } else {
652  TPolyLine *pL = new TPolyLine(P,x.data,y.data);
653  pL->SetLineColor(colorGalacticDisk);
654  pL->Draw();
655  gdL.push_back(pL);
656  P=0;
657  }
658  }
659  }
660 
661  if (drawWorldMap && (coordinate.CompareTo("CELESTIAL")!=0)) {
662  ClearWorldMap();
663  if (wm_size==0) {
664  // geographic coordinates
666  }
667  for (int n=0;n<wm_size;n++) {
668  double x,y;
669  double phi = wm_lon[n];
670  double theta = wm_lat[n];
671  if (coordinate.CompareTo("CWB")==0) {phi+=180;if(phi>180) phi=phi-360;}
672  if (coordinate.CompareTo("CELESTIAL")==0) {phi+=180;if(phi>180) phi=phi-360;phi=360-phi;} //TOBE CHECK
673  if (projection.CompareTo("hammer")==0) {
674  ProjectHammer(phi, theta, x, y);
675  } else
676  if (projection.CompareTo("sinusoidal")==0) {
677  ProjectSinusoidal(phi, theta, x, y);
678  } else
679  if (projection.CompareTo("parabolic")==0) {
680  ProjectParabolic(phi, theta, x, y);
681  } else {
682  x=phi;
683  y=theta;
684  }
685  if (coordinate.CompareTo("CWB")==0) {GeographicToCwb(x+180,y,x,y);}
686  if (coordinate.CompareTo("CELESTIAL")==0) {CelestialToCwb(x+180,y,x,y);} // TO BE CHECKED
687  TMarker *mP = new TMarker(x,y, 1); // WARNING : style=20,size=0.1 has issue when printed
688  mP->SetMarkerColor(kGray+1);
689  mP->Draw();
690  wmM.push_back(mP);
691  }
692  }
693 
694  if (coordinate.CompareTo("CELESTIAL")==0) {
695  ReverseXAxis(h2);
696  }
697 
698  TGaxis::SetMaxDigits();
699 }
700 
701 //______________________________________________________________________________
702 void
703 gskymap::DrawMarker(double ra, double dec, double gps, int marker, Size_t msize, Color_t mcolor) {
704 //
705 // Draw a marker in the ra,dec coordinates (see ROOT TMarker input parameters)
706 //
707 // Input: ra,dec - celestial coordinates
708 // gps - gps time
709 // mcolor - color of marker
710 // msize - size of marker
711 // mstyle - style of marker
712 //
713 
714  double phi,theta;
715  GeographicToCwb(ra,dec,phi,theta);
716  phi = gps>0 ? skymap::RA2phi(phi, gps) : phi;
717  if(coordinate.CompareTo("GEOGRAPHIC")==0) CwbToGeographic(phi,theta,phi,theta);
718  if(coordinate.CompareTo("CELESTIAL")==0) CwbToCelestial(phi,theta,phi,theta);
719  DrawMarker(phi, theta, marker, msize, mcolor);
720 }
721 
722 //______________________________________________________________________________
723 void
724 gskymap::DrawMarker(double phi, double theta, int marker, Size_t msize, Color_t mcolor) {
725 //
726 // Draw a marker in the phi,theta coordinates (see ROOT TMarker input parameters)
727 //
728 // Input: phi,theta - skymap coordinates
729 // mcolor - color of marker
730 // msize - size of marker
731 // mstyle - style of marker
732 //
733 
734  if(goff) return;
735 
736  double x,y;
737  if (projection.CompareTo("hammer")==0) {
738  //ProjectHammer((phi-180), -(theta-90), x, y);
739  ProjectHammer(phi, theta, x, y);
740  } else
741  if (projection.CompareTo("sinusoidal")==0) {
742  //ProjectSinusoidal((phi-180), -(theta-90), x, y);
743  ProjectSinusoidal(phi, theta, x, y);
744  } else
745  if (projection.CompareTo("parabolic")==0) {
746  //ProjectParabolic((phi-180), -(theta-90), x, y);
747  ProjectParabolic(phi, theta, x, y);
748  } else {
749  x=phi;
750  y=theta;
751  }
752 
753  TMarker *mP = new TMarker(x,y, marker);
754  mP->SetMarkerSize(msize);
755  mP->SetMarkerColor(mcolor);
756  mP->Draw();
757 }
758 
759 //______________________________________________________________________________
760 void
761 gskymap::DrawText(double ra, double dec, double gps, TString text, double tsize, Color_t tcolor) {
762 //
763 // Draw a text to ra,dec coordinates
764 //
765 // Input: ra,dec - celestial coordinates
766 // text - text
767 // tsize - size of text
768 // tcolor - color of text
769 //
770 
771  double phi,theta;
772  GeographicToCwb(ra,dec,phi,theta);
773  phi = gps>0 ? skymap::RA2phi(phi, gps) : phi;
774  if(coordinate.CompareTo("GEOGRAPHIC")==0) CwbToGeographic(phi,theta,phi,theta);
775  if(coordinate.CompareTo("CELESTIAL")==0) CwbToCelestial(phi,theta,phi,theta);
776  DrawText(phi, theta, text, tsize, tcolor);
777 }
778 
779 //______________________________________________________________________________
780 void
781 gskymap::DrawText(double phi, double theta, TString text, double tsize, Color_t tcolor) {
782 //
783 // Draw a text to phi,theta coordinates
784 //
785 // Input: phi,theta - skymap coordinates
786 // text - text
787 // tsize - size of text
788 // tcolor - color of text
789 //
790 
791  if(goff) return;
792 
793  double x,y;
794  if (projection.CompareTo("hammer")==0) {
795  //ProjectHammer((phi-180), -(theta-90), x, y);
796  ProjectHammer(phi, theta, x, y);
797  } else
798  if (projection.CompareTo("sinusoidal")==0) {
799  //ProjectSinusoidal((phi-180), -(theta-90), x, y);
800  ProjectSinusoidal(phi, theta, x, y);
801  } else
802  if (projection.CompareTo("parabolic")==0) {
803  //ProjectParabolic((phi-180), -(theta-90), x, y);
804  ProjectParabolic(phi, theta, x, y);
805  } else {
806  x=phi;
807  y=theta;
808  }
809 
810  TText *tP = new TText(x,y,text);
811  tP->SetTextSize(tsize);
812  tP->SetTextColor(tcolor);
813  tP->Draw();
814 }
815 
816 //______________________________________________________________________________
817 void
818 gskymap::ProjectHammer(Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
819 {
820 //
821 // Static function
822 // Convert Right Ascension, Declination to X,Y using an AITOFF projection.
823 // This procedure can be used to create an all-sky map in Galactic
824 // coordinates with an equal-area Aitoff projection. Output map
825 // coordinates are zero longitude centered.
826 // Also called Hammer-Aitoff projection (first presented by Ernst von Hammer in 1892)
827 // source: GMT
828 // code from Ernst-Jan Buis
829 // http://en.wikipedia.org/wiki/Hammer_projection
830 //
831 // INPUT : Geographic Coordinates
832 // LONGITUDE (W:-Pi,E:+Pi)
833 // LATITUDE (N:+Pi/2,S:-Pi/2)
834 //
835 
836  Double_t x, y;
837 
838  Double_t alpha2 = (l/2)*TMath::DegToRad();
839  Double_t delta = b*TMath::DegToRad();
840  Double_t r2 = TMath::Sqrt(2.);
841  Double_t f = 2*r2/TMath::Pi();
842  Double_t cdec = TMath::Cos(delta);
843  Double_t denom = TMath::Sqrt(1. + cdec*TMath::Cos(alpha2));
844  x = cdec*TMath::Sin(alpha2)*2.*r2/denom;
845  y = TMath::Sin(delta)*r2/denom;
846  x *= TMath::RadToDeg()/f;
847  y *= TMath::RadToDeg()/f;
848  // x *= -1.; // for a skymap swap left<->right
849  Al = x;
850  Ab = y;
851 
852  return;
853 }
854 
855 //______________________________________________________________________________
856 void
857 gskymap::ProjectSinusoidal(Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
858 {
859 //
860 // Static function
861 // code from Ernst-Jan Buis
862 //
863 // INPUT : Geographic Coordinates
864 // LONGITUDE (W:-Pi,E:+Pi)
865 // LATITUDE (N:+Pi/2,S:-Pi/2)
866 //
867 
868  Al = l*cos(b*TMath::DegToRad());
869  Ab = b;
870  return;
871 }
872 
873 //______________________________________________________________________________
874 void
875 gskymap::ProjectParabolic(Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
876 {
877 //
878 // Static function
879 // code from Ernst-Jan Buis
880 //
881 // INPUT : Geographic Coordinates
882 // LONGITUDE (W:-Pi,E:+Pi)
883 // LATITUDE (N:+Pi/2,S:-Pi/2)
884 //
885 
886  Al = l*(2.*TMath::Cos(2*b*TMath::DegToRad()/3) - 1);
887  Ab = 180*TMath::Sin(b*TMath::DegToRad()/3);
888  return;
889 }
890 
891 //______________________________________________________________________________
892 void
893 gskymap::SetPlotStyle(int paletteId) {
894 //
895 // set the color palette used for the skymap plot
896 //
897 // Input: paletteId - palette number (1,2,3,4,5)
898 //
899 // -----------------------------------------------------------------
900 // http://ultrahigh.org/2007/08/20/making-pretty-root-color-palettes/
901 // -----------------------------------------------------------------
902 //
903 
904  const Int_t NRGBs = 5;
905  const Int_t NCont = 255;
906 
907  if (fabs(paletteId)==1) {
908  Double_t stops[NRGBs] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
909  Double_t red[NRGBs] = { 0.00, 0.00, 0.87, 1.00, 0.51 };
910  Double_t green[NRGBs] = { 0.00, 0.81, 1.00, 0.20, 0.00 };
911  Double_t blue[NRGBs] = { 0.51, 1.00, 0.12, 0.00, 0.00 };
912  if (paletteId<0) {
913  TColor::CreateGradientColorTable(NRGBs, stops, blue, green, red, NCont);
914  } else {
915  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
916  }
917  } else
918  if (fabs(paletteId)==2) {
919  Double_t stops[NRGBs] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
920  Double_t red[NRGBs] = { 0.00, 0.00, 0.87, 1.00, 0.51 };
921  //Double_t red[NRGBs] = { 0.00, 0.00, 0.00, 1.00, 1.00 };
922  Double_t green[NRGBs] = { 0.00, 0.81, 1.00, 0.20, 0.00 };
923  //Double_t green[NRGBs] = { 0.00, 1.00, 1.00, 1.00, 0.00 };
924  Double_t blue[NRGBs] = { 1.00, 1.00, 0.00, 0.00, 0.00 };
925  if (paletteId<0) {
926  TColor::CreateGradientColorTable(NRGBs, stops, blue, green, red, NCont);
927  } else {
928  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
929  }
930  } else
931  if (fabs(paletteId)==3) {
932  Double_t stops[NRGBs] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
933  Double_t red[NRGBs] = { 0.00, 0.09, 0.18, 0.09, 0.00 };
934  Double_t green[NRGBs] = { 0.01, 0.02, 0.39, 0.68, 0.97 };
935  Double_t blue[NRGBs] = { 0.17, 0.39, 0.62, 0.79, 0.97 };
936  if (paletteId<0) {
937  TColor::CreateGradientColorTable(NRGBs, stops, blue, green, red, NCont);
938  } else {
939  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
940  }
941  } else
942  if (fabs(paletteId)==4) {
943  Double_t stops[NRGBs] = { 0.00, 0.50, 0.75, 0.875, 1.00 };
944  Double_t red[NRGBs] = { 1.00, 1.00, 1.00, 1.00, 1.00 };
945  Double_t green[NRGBs] = { 1.00, 0.75, 0.50, 0.25, 0.00 };
946  Double_t blue[NRGBs] = { 0.00, 0.00, 0.00, 0.00, 0.00 };
947  if (paletteId<0) {
948  TColor::CreateGradientColorTable(NRGBs, stops, blue, green, red, NCont);
949  } else {
950  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
951  }
952  } else
953  if (fabs(paletteId)==5) { // Greyscale palette
954  Double_t stops[NRGBs] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
955  Double_t red[NRGBs] = { 1.00, 0.84, 0.61, 0.34, 0.00 };
956  Double_t green[NRGBs] = { 1.00, 0.84, 0.61, 0.34, 0.00 };
957  Double_t blue[NRGBs] = { 1.00, 0.84, 0.61, 0.34, 0.00 };
958  if (paletteId<0) {
959  TColor::CreateGradientColorTable(NRGBs, stops, blue, green, red, NCont);
960  } else {
961  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
962  }
963  } else
964  if (fabs(paletteId)==57) { // kBird palette (is defined only in ROOT6)
965  Double_t stops[9] = { 0., 0.125, 0.25, 0.375, 0.5, 0.625, 0.75, 0.875, 1.0};
966  Double_t red[9] = { 0.2082, 0.0592, 0.0780, 0.0232, 0.1802, 0.5301, 0.8186, 0.9956, 0.9764};
967  Double_t green[9] = { 0.1664, 0.3599, 0.5041, 0.6419, 0.7178, 0.7492, 0.7328, 0.7862, 0.9832};
968  Double_t blue[9] = { 0.5293, 0.8684, 0.8385, 0.7914, 0.6425, 0.4662, 0.3499, 0.1968, 0.0539};
969  if (paletteId<0) {
970  TColor::CreateGradientColorTable(9, stops, blue, green, red, NCont);
971  } else {
972  TColor::CreateGradientColorTable(9, stops, red, green, blue, NCont);
973  }
974  }
975 
976  gStyle->SetNumberContours(NCont);
977 
978  return;
979 }
980 
981 //______________________________________________________________________________
982 int
983 gskymap::ReadWorlMapCoastLine(double*& wm_lon, double*& wm_lat) {
984 //
985 // Read the World Coastline map
986 //
987 // return in longitudes/latitudes arrays the world cordinates
988 //
989 // -----------------------------------------------------------------
990 // Coastline map (shapefile) download from :
991 // http://www.naturalearthdata.com/downloads/10m-physical-vectors/10m-coastline/
992 // hape2Text transforms ESRI shapefiles into Ascii files
993 // http://www.zonums.com/shp2text.html
994 // -----------------------------------------------------------------
995 //
996 
997  char ifileName[256];
998  sprintf(ifileName,"%s/%s",worldMapPath.Data(),WorlMapCoastLine);
999 
1000  ifstream in;
1001  in.open(ifileName,ios::in);
1002  if (!in.good()) {cout << "gskymap::ReadWorlMapCoastLine: Error Opening File : " << ifileName << endl;exit(1);}
1003 
1004  int cnt=0;
1005  char iline[1024];
1006  int size = WM_ENTRIES;
1007  wm_lon = new double[size];
1008  wm_lat = new double[size];
1009  while (1) {
1010 
1011  in.getline(iline,1024);
1012  if (!in.good()) break;
1013  //cout << ++cnt << " " << iline << endl;
1014  TObjArray* tok = TString(iline).Tokenize(TString(' '));
1015 
1016  if (tok->GetEntries()==1) continue;
1017 
1018  TObjString* tra = (TObjString*)tok->At(0);
1019  TObjString* tdec = (TObjString*)tok->At(1);
1020 
1021  if (tdec->GetString().IsAlpha()==1) continue;
1022  //cout << tra->GetString().Data() << " " << tdec->GetString().Data() << endl;
1023 
1024  double ra = tra->GetString().Atof();
1025  double dec = tdec->GetString().Atof();
1026 
1027  if (fabs(dec)<0.001) continue;
1028 
1029  wm_lon[cnt] = ra;
1030  wm_lat[cnt] = dec;
1031  cnt++;
1032 
1033  }
1034  cout << cnt << endl;
1035 
1036  in.close();
1037 
1038  drawWorldMap = false;
1039 
1040  return size;
1041 }
1042 
1043 //______________________________________________________________________________
1044 void
1046 //
1047 // clear the axis labels
1048 //
1049 
1050  for(int i=0; i<(int)axisT.size(); i++) {if(axisT[i]) delete axisT[i];}
1051  axisT.clear();
1052  return;
1053 }
1054 
1055 //______________________________________________________________________________
1056 void
1058 //
1059 // clear galactic disk
1060 //
1061 
1062  for(int i=0; i<(int)gdL.size(); i++) {if(gdL[i]) delete gdL[i];}
1063  gdL.clear();
1064  return;
1065 }
1066 
1067 //______________________________________________________________________________
1068 void
1070 //
1071 // clear world map
1072 //
1073 
1074  for(int i=0; i<(int)wmM.size(); i++) {if(wmM[i]) delete wmM[i];}
1075  wmM.clear();
1076  return;
1077 }
1078 
1079 //______________________________________________________________________________
1080 void
1082 //
1083 // clear grids along the x axis
1084 //
1085 
1086  for(int i=0; i<(int)gridxL.size(); i++) {if(gridxL[i]) delete gridxL[i];}
1087  gridxL.clear();
1088 }
1089 
1090 //______________________________________________________________________________
1091 void
1093 //
1094 // clear grids along the y axis
1095 //
1096 
1097  for(int i=0; i<(int)gridyL.size(); i++) {if(gridyL[i]) delete gridyL[i];}
1098  gridyL.clear();
1099  return;
1100 }
1101 
1102 //______________________________________________________________________________
1103 void
1105 //
1106 // save plot to file
1107 //
1108 // Input: pname - output file name
1109 //
1110 
1111  if(canvas==NULL) {cout << "gskymap::Print: Error - canvas not initialized " << endl;exit(1);}
1112 
1113  TGaxis::SetMaxDigits(3);
1114 
1115 /*
1116  if(TString(pname).Contains(".png")!=0) { // fix gray background for png plots
1117  TString gname=pname;
1118  gname.ReplaceAll(".png",".gif");
1119  canvas->Print(gname);
1120  char cmd[1024];
1121  sprintf(cmd,"convert %s %s",gname.Data(),pname.Data());
1122  //cout << cmd << endl;
1123  gSystem->Exec(cmd);
1124  sprintf(cmd,"rm %s",gname.Data());
1125  //cout << cmd << endl;
1126  gSystem->Exec(cmd);
1127  } else {
1128  pname.ReplaceAll(".PNG",".png");
1129  canvas->Print(pname);
1130  }
1131 */
1132 
1133  pname.ReplaceAll(".PNG",".png");
1134  canvas->Print(pname);
1135 
1136  TGaxis::SetMaxDigits();
1137 }
1138 
1139 //______________________________________________________________________________
1140 void
1142 //
1143 // dump skymap to ascii file
1144 //
1145 // Input: fname - output file name
1146 //
1147 // format :
1148 //
1149 // header : sky_res -> area of a pixel (degrees)
1150 // theta_1 -> theta min
1151 // theta_2 -> theta max
1152 // phi_1 -> phi min
1153 // phi_2 -> phi max
1154 // lenght -> number of entries
1155 //
1156 // list of data: sky_index value
1157 //
1158 
1159  int L = skymap::size();
1160 
1161  double x;
1162  FILE* fp; // dump file
1163 
1164  if((fp = fopen(fname, "w")) == NULL) {
1165  cout << "netevent::DumpSkyMap() error: cannot open file " << fname <<". \n";
1166  return;
1167  };
1168  fprintf(fp,"wat %s\n",watversion('f'));
1169  fprintf(fp,"class skymap\n");
1170  x = cos(this->theta_1*PI/180.)-cos(this->theta_2*PI/180.);
1171  x*= (this->phi_2-this->phi_1)*180/PI/L;
1172  x = double(int(1000*x))/1000.;
1173  fprintf(fp,"sky_res %f\n",sqrt(fabs(x)));
1174  fprintf(fp,"theta_1 %f\n",this->theta_1);
1175  fprintf(fp,"theta_2 %f\n",this->theta_2);
1176  fprintf(fp,"phi_1 %f\n",this->phi_1);
1177  fprintf(fp,"phi_2 %f\n",this->phi_2);
1178  fprintf(fp,"lenght %d\n",L);
1179 
1180  for(int l=0;l<L;l++) fprintf(fp,"%d %e\n",l,this->get(l));
1181 
1182  fclose(fp);
1183  return;
1184 }
1185 
1186 //______________________________________________________________________________
1187 void
1189 //
1190 // reverse the x axis
1191 //
1192 
1193  char xlabel[16];
1194  for (int i=1;i<=h2->GetNbinsX();i++) {
1195  //cout << i << " " << h2->GetXaxis()->GetBinCenter(i) << endl;
1196  double xvalue = h2->GetXaxis()->GetBinCenter(i)+h2->GetXaxis()->GetBinWidth(i)/2.;
1197  int nbin=int(h2->GetNbinsX()/4.);
1198  if(i%nbin==0) {
1199  sprintf(xlabel,"%d",int(360-xvalue));
1200  //cout << i << " " << xvalue << endl;
1201  h2->GetXaxis()->SetBinLabel(i,xlabel);
1202  h2->GetXaxis()->SetLabelSize(0.06);
1203  h2->GetXaxis()->LabelsOption("h");
1204  h2->GetXaxis()->SetNdivisions(604);
1205  h2->GetXaxis()->SetLabelFont(42);
1206  h2->GetXaxis()->SetLabelOffset(0.010);
1207  h2->GetXaxis()->SetTitleOffset(1.5);
1208  }
1209  }
1210 }
1211 
1212 //______________________________________________________________________________
1213 void gskymap::DumpObject(const char* file, const char* name)
1214 {
1215 //
1216 // dump skymap into root file
1217 //
1218 
1219  TFile* rfile = new TFile(const_cast<char*>(file), "RECREATE");
1220  this->Write(const_cast<char*>(name)); // write this object
1221  rfile->Close();
1222 }
1223 
1224 //______________________________________________________________________________
1225 void gskymap::LoadObject(const char* file, const char* name)
1226 {
1227 //
1228 // load gskymap object from root file
1229 //
1230 
1231  TFile* rfile = new TFile(const_cast<char*>(file));
1232  gskymap* gsm = (gskymap*)rfile->Get(const_cast<char*>(name));
1233  if(gsm==NULL) {
1234  cout << "gskymap::LoadObject : Error - input file don't contains object gskymap" << endl;
1235  exit(1);
1236  }
1237  *this=*gsm;
1238 
1239  gridxL.clear();
1240  gridyL.clear();
1241  gdL.clear();
1242  wmM.clear();
1243  axisT.clear();
1244 
1245  FillData();
1246  rfile->Close();
1247 }
1248 
1249 //______________________________________________________________________________
1251 {
1252 //
1253 // Used to draw skymap from TBrowser
1254 //
1255 
1256  gridxL.clear();
1257  gridyL.clear();
1258  gdL.clear();
1259  wmM.clear();
1260  axisT.clear();
1261 
1262  FillData();
1263  Draw();
1264 }
1265 
void CreateCanvas()
Definition: gskymap.cc:164
#define NRGBs
TText * text
Definition: cbc_plots.C:717
Color_t colorGridy
Definition: gskymap.hh:202
double * wm_lon
Definition: gskymap.hh:227
Double_t green[nRGBs]
tuple f
Definition: cwb_online.py:91
std::vector< TPolyLine * > gdL
Definition: gskymap.hh:223
TString ifileName
Definition: MergeTrees.C:35
fprintf(stdout,"start=%f duration=%f rate=%f\n", x.start(), x.size()/x.rate(), x.rate())
std::vector< TPolyLine * > gridxL
Definition: gskymap.hh:221
double delta
void CwbToCelestial(double ilongitude, double ilatitude, double &olongitude, double &olatitude, double gps=0)
Definition: skycoord.hh:410
Color_t colorGridx
Definition: gskymap.hh:200
void SetPlotStyle(int paletteId=1)
Definition: gskymap.cc:893
double gpsGalacticDisk
Definition: gskymap.hh:213
par[0] name
int n
Definition: cwb_net.C:10
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
char * watversion(char c='s')
Definition: watversion.hh:21
TString("c")
int palette
Definition: DrawGnetwork2.C:17
Int_t ww
Definition: gskymap.hh:219
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<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t 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
float theta
void DrawMarker(double phi, double theta, int marker, Size_t msize=1, Color_t tcolor=1)
Definition: gskymap.cc:724
TString coordinate
Definition: gskymap.hh:211
void DrawText(double phi, double theta, TString text, double tsize=0.04, Color_t tcolor=1)
Definition: gskymap.cc:781
TH2F * ph
TRandom3 P
Definition: compare_bkg.C:296
Double_t blue[nRGBs]
STL namespace.
std::vector< int > index
Definition: skymap.hh:301
double theta_1
Definition: skymap.hh:303
Long_t size
TCanvas * canvas
Definition: gskymap.hh:194
bool isGridx
Definition: gskymap.hh:199
int j
Definition: cwb_net.C:10
i drho i
int wm_size
Definition: gskymap.hh:208
void Draw(int dpaletteId=0, Option_t *option="colfz")
Definition: gskymap.cc:442
virtual ~gskymap()
Definition: gskymap.cc:210
TH2D * h2
`
Definition: gskymap.hh:195
double theta_2
Definition: skymap.hh:304
float distance
Definition: cbc_plots.C:436
#define NCont
#define PI
Definition: watfun.hh:14
#define DUMMY_PALETTE_ID
Definition: gskymap.hh:47
void ClearGridx()
Definition: gskymap.cc:1081
void ClearGridy()
Definition: gskymap.cc:1092
tlive_fix Write()
float phi
Double_t stops[nRGBs]
double ra
Definition: ConvertGWGC.C:46
std::vector< vectorD > value
Definition: skymap.hh:300
void Plot()
Definition: gskymap.cc:1250
static double r2
Definition: geodesics.cc:8
void ClearGalacticDisk()
Definition: gskymap.cc:1057
bool isGridy
Definition: gskymap.hh:201
gSM SetOptions(cwb_antpat_projection, COORDINATES, RESOLUTION/2)
i() int(T_cor *100))
double Pi
Int_t wtopy
Definition: gskymap.hh:219
std::vector< TText * > axisT
Definition: gskymap.hh:225
double sms
Definition: skymap.hh:302
char fname[1024]
double * wm_lat
Definition: gskymap.hh:228
int k
double RA2phi(double ph, double gps)
Definition: skymap.hh:195
#define WM_ENTRIES
Definition: gskymap.hh:46
void ClearAxisLabel()
Definition: gskymap.cc:1045
bool isLogz
Definition: gskymap.hh:203
WSeries< double > ww
Definition: Regression_H1.C:33
void ReverseXAxis(TH2D *h2)
Definition: gskymap.cc:1188
#define WorlMapCoastLine
Definition: gskymap.hh:45
int paletteId
Definition: gskymap.hh:210
double resolution
Definition: gskymap.hh:206
char title[256]
Definition: SSeriesExample.C:1
void CelestialToCwb(double ilongitude, double ilatitude, double &olongitude, double &olatitude, double gps=0)
Definition: skycoord.hh:419
double gps
TString worldMapPath
Definition: gskymap.hh:209
Int_t wtopx
Definition: gskymap.hh:219
ifstream in
void ProjectParabolic(Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
Definition: gskymap.cc:875
double fabs(const Complex &x)
Definition: numpy.cc:37
string file
Definition: cwb_online.py:385
void ProjectHammer(Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
Definition: gskymap.cc:818
void FillData()
Definition: gskymap.cc:339
void LoadObject(const char *file, const char *name="gskymap")
Definition: gskymap.cc:1225
void ClearWorldMap()
Definition: gskymap.cc:1069
void CwbToGeographic(double ilongitude, double ilatitude, double &olongitude, double &olatitude)
Definition: skycoord.hh:396
int cnt
void DumpObject(const char *file, const char *name="gskymap")
Definition: gskymap.cc:1213
int l
Definition: cbc_plots.C:434
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
TObjString * tra
Definition: ConvertGWGC.C:40
void GalacticToEquatorial(double ilongitude, double ilatitude, double &olongitude, double &olatitude)
Definition: skycoord.hh:29
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
void ProjectSinusoidal(Double_t l, Double_t b, Double_t &Al, Double_t &Ab)
Definition: gskymap.cc:857
DataType_t * data
Definition: wavearray.hh:301
void DumpSkyMap(char *fname)
Definition: gskymap.cc:1141
gskymap() gskymap(double sms, double t1=0., double t2=180., double p1=0., double p2=360.) gskymap(char *ifile changed)
Definition: gskymap.hh:65
TObjString * tdec
Definition: ConvertGWGC.C:41
TString name
Definition: gskymap.hh:218
double phi_1
Definition: skymap.hh:305
double get(size_t i)
param: sky index
Definition: skymap.cc:681
void GeographicToCwb(double ilongitude, double ilatitude, double &olongitude, double &olatitude)
Definition: skycoord.hh:403
bool goff
Definition: gskymap.hh:207
fclose(ftrig)
Color_t colorGalacticDisk
Definition: gskymap.hh:214
int ReadWorlMapCoastLine(double *&wm_lon, double *&wm_lat)
Definition: gskymap.cc:983
size_t size()
Definition: skymap.hh:118
void Print(TString pname)
Definition: gskymap.cc:1104
std::vector< TPolyLine * > gridyL
Definition: gskymap.hh:222
wavearray< double > y
Definition: Test10.C:31
std::vector< TMarker * > wmM
Definition: gskymap.hh:224
TString projection
Definition: gskymap.hh:212
Double_t red[nRGBs]
exit(0)
bool drawWorldMap
Definition: gskymap.hh:205
Int_t wh
Definition: gskymap.hh:219