source: trunk/MagicSoft/Mars/mastro/MAstroCatalog.cc@ 3705

Last change on this file since 3705 was 3704, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 31.2 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 03/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2002-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MAstroCatalog
28//
29// THIS IMPLEMENTATION IS PRELIMINARY AND WILL BE MERGED WITH
30// SOME PARTS OF THE DRIVE SOFTWARE SOON!
31//
32// To display a starfiield you must have a supported catalog, then do:
33//
34// MTime time;
35// // Time for which to get the picture
36// time.Set(2004, 2, 28, 20, 14, 7);
37// // Current observatory
38// MObservatory magic1;
39// // Right Ascension [h] and declination [deg] of source
40// // Currently 'perfect' pointing is assumed
41// const Double_t ra = MAstro::Hms2Rad(5, 34, 31.9);
42// const Double_t dec = MAstro::Dms2Rad(22, 0, 52.0);
43// MAstroCatalog stars;
44// // Magnitude up to which the stars are loaded from the catalog
45// stars.SetLimMag(6);
46// // Radius of FOV around the source position to load the stars
47// stars.SetRadiusFOV(3);
48// // Source position
49// stars.SetRaDec(ra, dec);
50// // Catalog to load (here: Bright Star Catalog V5)
51// stars.ReadBSC("bsc5.dat");
52// // Obersavatory and time to also get local coordinate information
53// stars.SetObservatory(magic1);
54// stars.SetTime(time);
55// // Enable interactive GUI
56// stars.SetGuiActive();
57// //Clone the catalog due to the validity range of the instance
58// TObject *o = stars.Clone();
59// o->SetBit(kCanDelete);
60// o->Draw();
61//
62// If no time and/or Obervatory location is given no local coordinate
63// information is displayed.
64//
65// The conversion from sky coordinates to local coordinates is done using
66// MAstroSky2Local which does a simple rotation of the coordinate system.
67// This is inaccurate in the order of 30arcsec due to ignorance of all
68// astrometrical corrections (nutation, precission, abberation, ...)
69//
70// GUI: If the gui is interactive you can use the cursor keys to change
71// the position you are looking at and with plus/minus you
72// can (un)zoom the FOV (Field Of View)
73//
74// ToDo:
75// - replace MVetcor3 by a more convinient class. Maybe use TExMap, too.
76//
77//////////////////////////////////////////////////////////////////////////////
78#include "MAstroCatalog.h"
79
80#include <errno.h>
81#include <fstream>
82#include <stdlib.h>
83#include <limits.h> // INT_MAX (Suse 7.3/gcc 2.95)
84
85#include <KeySymbols.h>
86
87#include <TPad.h> // TPad::GetMaxPickDistance
88#include <TLine.h>
89#include <TMarker.h>
90#include <TCanvas.h>
91#include <TArrayI.h>
92#include <TGToolTip.h>
93#include <TRotation.h>
94#include <TPaveText.h>
95
96#include "MLog.h"
97#include "MLogManip.h"
98
99#include "MTime.h"
100#include "MAstro.h"
101#include "MAstroSky2Local.h"
102#include "MObservatory.h"
103
104#undef DEBUG
105//#define DEBUG
106
107#ifdef DEBUG
108#include <TStopwatch.h>
109#endif
110
111ClassImp(MVector3);
112ClassImp(MAstroCatalog);
113
114using namespace std;
115
116// --------------------------------------------------------------------------
117//
118// Default Constructor. Set Default values:
119// fLimMag = 99
120// fRadiusFOV = 90
121//
122MAstroCatalog::MAstroCatalog() : fLimMag(99), fRadiusFOV(90), fToolTip(0), fObservatory(0), fTime(0)
123{
124 fList.SetOwner();
125 fToolTip = new TGToolTip(0, "", 0);
126}
127
128// --------------------------------------------------------------------------
129//
130// Destructor. Delete fTime, fObservatory. Disconnect signal. delete tooltip.
131// Delete Map with gui primitives
132//
133MAstroCatalog::~MAstroCatalog()
134{
135 // First disconnect the EventInfo...
136 // FIXME: There must be an easier way!
137 TIter Next(gROOT->GetListOfCanvases());
138 TCanvas *c;
139 while ((c=(TCanvas*)Next()))
140 c->Disconnect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", this,
141 "EventInfo(Int_t,Int_t,Int_t,TObject*)");
142
143 // Now delete the data members
144 if (fTime)
145 delete fTime;
146 if (fObservatory)
147 delete fObservatory;
148
149 fToolTip->Hide();
150 delete fToolTip;
151
152 DeleteMap();
153}
154
155// --------------------------------------------------------------------------
156//
157// Snippet to for reading ctalog files.
158//
159TString MAstroCatalog::FindToken(TString &line, Char_t tok)
160{
161 Ssiz_t token = line.First(tok);
162 if (token<0)
163 {
164 const TString copy(line);
165 line = "";
166 return copy;
167 }
168
169 const TString res = line(0, token);
170 line.Remove(0, token+1);
171 return res;
172}
173
174// --------------------------------------------------------------------------
175//
176// return int correspoding to TSubString
177//
178Int_t MAstroCatalog::atoi(const TSubString &sub)
179{
180 return atoi(TString(sub));
181}
182
183// --------------------------------------------------------------------------
184//
185// return float correspoding to TSubString
186//
187Float_t MAstroCatalog::atof(const TSubString &sub)
188{
189 return atof(TString(sub));
190}
191
192// --------------------------------------------------------------------------
193//
194// return int correspoding to TString
195//
196Int_t MAstroCatalog::atoi(const TString &s)
197{
198 return std::atoi(s);
199}
200
201// --------------------------------------------------------------------------
202//
203// return float correspoding to TString
204//
205Float_t MAstroCatalog::atof(const TString &s)
206{
207 return std::atof(s);
208}
209
210// --------------------------------------------------------------------------
211//
212// Read from a xephem catalog, set bit kHasChahanged.
213// Already read data is not deleted. To delete the stored data call
214// Delete().
215//
216Int_t MAstroCatalog::ReadXephem(TString catalog)
217{
218 SetBit(kHasChanged);
219
220 gLog << inf << "Reading Xephem catalog: " << catalog << endl;
221
222 ifstream fin(catalog);
223 if (!fin)
224 {
225 gLog << err << "Cannot open file " << catalog << ": ";
226 gLog << strerror(errno) << endl;
227 return 0;
228 }
229
230 Int_t add =0;
231 Int_t cnt =0;
232 Int_t line=0;
233
234 Double_t maxmag=0;
235
236 while (1)
237 {
238 TString row;
239 row.ReadLine(fin);
240 if (!fin)
241 break;
242
243 line++;
244
245 if (row[0]=='#')
246 continue;
247
248 TString line(row);
249
250 TString name = FindToken(line);
251 TString dummy = FindToken(line);
252 TString r = FindToken(line);
253 TString d = FindToken(line);
254 TString m = FindToken(line);
255 TString epoch = FindToken(line);
256
257 if (name.IsNull() || r.IsNull() || d.IsNull() || m.IsNull() || epoch.IsNull())
258 {
259 gLog << warn << "Invalid Entry Line #" << line << ": " << row << endl;
260 continue;
261 }
262
263 cnt++;
264
265 const Double_t mag = atof(m);
266
267 maxmag = TMath::Max(maxmag, mag);
268
269 if (mag>fLimMag)
270 continue;
271
272 if (epoch!="2000")
273 {
274 gLog << warn << "Epoch != 2000... skipped." << endl;
275 continue;
276 }
277
278 Double_t ra0, dec0;
279 MAstro::Coordinate2Angle(r, ra0);
280 MAstro::Coordinate2Angle(d, dec0);
281
282 ra0 *= TMath::Pi()/12;
283 dec0 *= TMath::Pi()/180;
284
285 MVector3 *star=new MVector3;
286 star->SetRaDec(ra0, dec0, mag);
287 star->SetName(name);
288 if (star->Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
289 {
290 delete star;
291 continue;
292 }
293
294 fList.Add(star);
295 add++;
296 }
297 gLog << inf << "Read " << add << " out of " << cnt << " (Total max mag=" << maxmag << ")" << endl;
298
299 return add;
300}
301
302// --------------------------------------------------------------------------
303//
304// Read from a NGC2000 catalog. set bit kHasChanged
305// Already read data is not deleted. To delete the stored data call
306// Delete().
307//
308Int_t MAstroCatalog::ReadNGC2000(TString catalog)
309{
310 SetBit(kHasChanged);
311
312 gLog << inf << "Reading NGC2000 catalog: " << catalog << endl;
313
314 ifstream fin(catalog);
315 if (!fin)
316 {
317 gLog << err << "Cannot open file " << catalog << ": ";
318 gLog << strerror(errno) << endl;
319 return 0;
320 }
321
322 Int_t add=0;
323 Int_t cnt=0;
324 Int_t n =0;
325
326 Double_t maxmag=0;
327
328 while (1)
329 {
330 TString row;
331 row.ReadLine(fin);
332 if (!fin)
333 break;
334
335 cnt++;
336
337 const Int_t rah = atoi(row(13, 2));
338 const Float_t ram = atof(row(16, 4));
339 const Char_t decs = row(22);
340 const Int_t decd = atoi(row(23, 2));
341 const Int_t decm = atoi(row(26, 2));
342 const TString m = row(43, 4);
343
344 if (m.Strip().IsNull())
345 continue;
346
347 n++;
348
349 const Double_t mag = atof(m);
350
351 maxmag = TMath::Max(maxmag, mag);
352
353 if (mag>fLimMag)
354 continue;
355
356 const Double_t ra = MAstro::Hms2Rad(rah, (int)ram, fmod(ram, 1)*60);
357 const Double_t dec = MAstro::Dms2Rad(decd, decm, 0, decs);
358
359 MVector3 *star=new MVector3;
360 star->SetRaDec(ra, dec, mag);
361 star->SetName(row(0, 8));
362 if (star->Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
363 {
364 delete star;
365 continue;
366 }
367
368 fList.Add(star);
369 add++;
370 }
371
372 gLog << inf << "Read " << add << " out of " << n << " (Total max mag=" << maxmag << ")" << endl;
373
374 return add;
375}
376
377// --------------------------------------------------------------------------
378//
379// Read from a Bright Star Catalog catalog. set bit kHasChanged
380// Already read data is not deleted. To delete the stored data call
381// Delete().
382//
383Int_t MAstroCatalog::ReadBSC(TString catalog)
384{
385 SetBit(kHasChanged);
386
387 gLog << inf << "Reading Bright Star Catalog (BSC5) catalog: " << catalog << endl;
388
389 ifstream fin(catalog);
390 if (!fin)
391 {
392 gLog << err << "Cannot open file " << catalog << ": ";
393 gLog << strerror(errno) << endl;
394 return 0;
395 }
396
397 Int_t add=0;
398 Int_t cnt=0;
399 Int_t n =0;
400
401 Double_t maxmag=0;
402
403 while (1)
404 {
405 TString row;
406 row.ReadLine(fin);
407 if (!fin)
408 break;
409
410 cnt++;
411
412 const Int_t rah = atoi(row(75, 2));
413 const Int_t ram = atoi(row(77, 2));
414 const Float_t ras = atof(row(79, 4));
415 const Char_t decsgn = row(83);
416 const Int_t decd = atoi(row(84, 2));
417 const Int_t decm = atoi(row(86, 2));
418 const Int_t decs = atoi(row(88, 2));
419 const TString m = row(102, 5);
420
421 if (m.Strip().IsNull())
422 continue;
423
424 n++;
425
426 const Double_t mag = atof(m.Data());
427
428 maxmag = TMath::Max(maxmag, mag);
429
430 if (mag>fLimMag)
431 continue;
432
433 const Double_t ra = MAstro::Hms2Rad(rah, ram, ras);
434 const Double_t dec = MAstro::Dms2Rad(decd, decm, decs, decsgn);
435
436 MVector3 *star=new MVector3;
437 star->SetRaDec(ra, dec, mag);
438 star->SetName(row(4,9));
439 if (star->Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
440 {
441 delete star;
442 continue;
443 }
444
445 fList.Add(star);
446 add++;
447 }
448
449 gLog << inf << "Read " << add << " out of " << n << " (Total max mag=" << maxmag << ")" << endl;
450
451 return add;
452}
453
454// --------------------------------------------------------------------------
455//
456// Set Range of pad. If something has changed create and draw new primitives.
457// Paint all gui primitives.
458//
459void MAstroCatalog::Paint(Option_t *o)
460{
461 SetRangePad(o);
462
463 if (TestBit(kHasChanged))
464 DrawPrimitives(o);
465
466 PaintMap();
467}
468
469// --------------------------------------------------------------------------
470//
471// Draw a black marker at the position of the star. Create a corresponding
472// tooltip with the coordinates.
473// x, y: Pad Coordinates to draw star
474// v: Sky position (Ra/Dec) of the star
475// transparent: Draw marker or tooltip only
476// txt: additional tooltip text
477//
478void MAstroCatalog::DrawStar(Double_t x, Double_t y, const TVector3 &v, Bool_t transparent, const char *txt)
479{
480 const Double_t ra = v.Phi()*TMath::RadToDeg()/15;
481 const Double_t dec = (TMath::Pi()/2-v.Theta())*TMath::RadToDeg();
482
483 TString str = v.GetName();
484 str += Form(": Ra=%.2fh", ra);
485 str += Form(" Dec=%.1fd", dec);
486 str += Form(" Mag=%.1f", -2.5*log10(v.Mag()));
487 if (txt)
488 str += Form(" (%s)", txt);
489
490 // draw star on the camera display
491 TMarker *tip=new TMarker(x, y, transparent ? kDot : kFullDotMedium);;
492 tip->SetMarkerColor(kBlack);
493 AddMap(tip, new TString(str));
494}
495
496// --------------------------------------------------------------------------
497//
498// Set pad as modified.
499//
500void MAstroCatalog::Update(Bool_t upd)
501{
502 if (gPad && TestBit(kMustCleanup))
503 {
504 SetBit(kHasChanged);
505 gPad->Modified();
506 if (upd)
507 gPad->Update();
508 }
509}
510
511// --------------------------------------------------------------------------
512//
513// Set the observation time. Necessary to use local coordinate
514// system. The MTime object is cloned.
515//
516void MAstroCatalog::SetTime(const MTime &time)
517{
518 if (fTime)
519 delete fTime;
520 fTime=(MTime*)time.Clone();
521}
522
523// --------------------------------------------------------------------------
524//
525// Set the observatory location. Necessary to use local coordinate
526// system. The MObservatory object is cloned.
527//
528void MAstroCatalog::SetObservatory(const MObservatory &obs)
529{
530 if (fObservatory)
531 delete fObservatory;
532 fObservatory=(MObservatory*)obs.Clone();
533}
534
535// --------------------------------------------------------------------------
536//
537// Convert the vector to pad coordinates. After conversion
538// the x- coordinate of the vector must be the x coordinate
539// of the pad - the same for y. If the coordinate is inside
540// the current draw area return kTRUE, otherwise kFALSE.
541// If it is an invalid coordinate return kERROR
542//
543Int_t MAstroCatalog::ConvertToPad(const TVector3 &w0, TVector2 &v) const
544{
545 TVector3 w(w0);
546
547 // Stretch such, that the Z-component is alwas the same. Now
548 // X and Y contains the intersection point between the star-light
549 // and the plain of a virtual plain screen (ccd...)
550 if (TestBit(kPlainScreen))
551 w *= 1./w(2);
552
553 w *= TMath::RadToDeg(); // FIXME: *conversion factor?
554 v.Set(TestBit(kMirrorX) ? -w(0) : w(0),
555 TestBit(kMirrorY) ? -w(1) : w(1));
556
557 if (w(2)<0)
558 return kERROR;
559
560 return v.X()>gPad->GetX1() && v.Y()>gPad->GetY1() &&
561 v.X()<gPad->GetX2() && v.Y()<gPad->GetY2();
562}
563
564// --------------------------------------------------------------------------
565//
566// Convert theta/phi coordinates of v by TRotation into new coordinate
567// system and convert the coordinated to pad by ConvertToPad.
568// The result is retunred in v.
569//
570Int_t MAstroCatalog::Convert(const TRotation &rot, TVector2 &v) const
571{
572 MVector3 w;
573 w.SetMagThetaPhi(1, v.Y(), v.X());
574 w *= rot;
575
576 return ConvertToPad(w, v);
577}
578
579// --------------------------------------------------------------------------
580//
581// Draw a line from v to v+(dx,dy) using Convert/ConvertToPad to get the
582// corresponding pad coordinates.
583//
584Bool_t MAstroCatalog::DrawLine(const TVector2 &v, Double_t dx, Double_t dy, const TRotation &rot, Int_t type)
585{
586 const TVector2 add(dx*TMath::DegToRad(), dy*TMath::DegToRad());
587
588 TVector2 v0 = v;
589 TVector2 v1 = v+add;
590
591 const Int_t rc0 = Convert(rot, v0);
592 const Int_t rc1 = Convert(rot, v1);
593
594 // Both are kFALSE or both are kERROR
595 if ((rc0|rc1)==kFALSE || (rc0&rc1)==kERROR)
596 return kFALSE;
597
598 TLine *line = new TLine(v0.X(), v0.Y(), v1.X(), v1.Y());
599 line->SetLineStyle(kDashDotted); //kDashed, kDotted, kDashDotted
600 line->SetLineColor(kWhite+type*2);
601 AddMap(line);
602
603 const TVector2 deg = v*TMath::RadToDeg();
604 TString txt = type==1 ?
605 Form("Ra=%.2fh Dec=%.1fd", fmod(deg.X()/15+48, 24), fmod(90-deg.Y()+270,180)-90) :
606 Form("Zd=%.1fd Az=%.1fd", fmod(deg.Y()+270,180)-90, fmod(deg.X()+720, 360));
607
608 TMarker *tip=new TMarker(v0.X(), v0.Y(), kDot);
609 tip->SetMarkerColor(kWhite+type*2);
610 AddMap(tip, new TString(txt));
611
612 return kTRUE;
613}
614
615// --------------------------------------------------------------------------
616//
617// Use "local" draw option to align the display to the local
618// coordinate system instead of the sky coordinate system.
619// dx, dy are arrays storing recuresively all touched points
620// stepx, stepy are the step-size of the current grid.
621//
622void MAstroCatalog::Draw(const TVector2 &v0, const TRotation &rot, TArrayI &dx, TArrayI &dy, Int_t stepx, Int_t stepy, Int_t type)
623{
624 // Calculate the end point
625 const TVector2 v1 = v0 + TVector2(dx[0]*TMath::DegToRad(), dy[0]*TMath::DegToRad());
626
627 // Check whether the point has already been touched.
628 Int_t idx[] = {1, 1, 1, 1};
629
630 Int_t dirs[4][2] = { {0, stepy}, {stepx, 0}, {0, -stepy}, {-stepx, 0} };
631
632 // Check for ambiguities.
633 for (int i=0; i<dx.GetSize(); i++)
634 {
635 for (int j=0; j<4; j++)
636 {
637 const Bool_t rcx0 = (dx[i]+720)%360==(dx[0]+dirs[j][0]+720)%360;
638 const Bool_t rcy0 = (dy[i]+360)%180==(dy[0]+dirs[j][1]+360)%180;
639 if (rcx0&&rcy0)
640 idx[j] = 0;
641 }
642 }
643
644 // Enhance size of array by 1, copy current
645 // position as last entry
646 dx.Set(dx.GetSize()+1);
647 dy.Set(dy.GetSize()+1);
648
649 dx[dx.GetSize()-1] = dx[0];
650 dy[dy.GetSize()-1] = dy[0];
651
652 // Store current positon
653 const Int_t d[2] = { dx[0], dy[0] };
654
655 for (int i=0; i<4; i++)
656 if (idx[i])
657 {
658 // Calculate new position
659 dx[0] = d[0]+dirs[i][0];
660 dy[0] = d[1]+dirs[i][1];
661
662 // Draw corresponding line and iterate through grid
663 if (DrawLine(v1, dirs[i][0], dirs[i][1], rot, type))
664 Draw(v0, rot, dx, dy, stepx, stepy, type);
665
666 dx[0]=d[0];
667 dy[0]=d[1];
668 }
669}
670
671// --------------------------------------------------------------------------
672//
673// Draw a grid recursively around the point v0 (either Ra/Dec or Zd/Az)
674// The points in the grid are converted by a TRotation and CovertToPad
675// to pad coordinates. The type arguemnts is neccessary to create the
676// correct tooltip (Ra/Dec, Zd/Az) at the grid-points.
677// From the pointing position the step-size of teh gris is caluclated.
678//
679void MAstroCatalog::DrawGrid(const TVector3 &v0, const TRotation &rot, Int_t type)
680{
681 TArrayI dx(1);
682 TArrayI dy(1);
683
684 // align to 1deg boundary
685 TVector2 v(v0.Phi()*TMath::RadToDeg(), v0.Theta()*TMath::RadToDeg());
686 v.Set((Float_t)TMath::Nint(v.X()), (Float_t)TMath::Nint(v.Y()));
687
688 // calculate stepsizes based on visible FOV
689 Int_t stepx = 1;
690
691 if (v.Y()<fRadiusFOV || v.Y()>180-fRadiusFOV)
692 stepx=36;
693 else
694 {
695 // This is a rough estimate how many degrees are visible
696 const Float_t m = log(fRadiusFOV/180.)/log(90./(fRadiusFOV+1)+1);
697 const Float_t t = log(180.)-m*log(fRadiusFOV);
698 const Float_t f = m*log(90-fabs(90-v.Y()))+t;
699 const Int_t nx = (Int_t)(exp(f)+0.5);
700 stepx = nx<4 ? 1 : nx/4;
701 if (stepx>36)
702 stepx=36;
703 }
704
705 const Int_t ny = (Int_t)(fRadiusFOV+1);
706 Int_t stepy = ny<4 ? 1 : ny/4;
707
708 // align stepsizes to be devisor or 180 and 90
709 while (180%stepx)
710 stepx++;
711 while (90%stepy)
712 stepy++;
713
714 // align to step-size boundary (search for the nearest one)
715 Int_t dv = 1;
716 while ((int)(v.X())%stepx)
717 {
718 v.Set(v.X()+dv, v.Y());
719 dv = -TMath::Sign(TMath::Abs(dv)+1, dv);
720 }
721
722 dv = 1;
723 while ((int)(v.Y())%stepy)
724 {
725 v.Set(v.X(), v.Y()+dv);
726 dv = -TMath::Sign(TMath::Abs(dv)+1, dv);
727 }
728
729 // draw...
730 v *= TMath::DegToRad();
731
732 Draw(v, rot, dx, dy, stepx, stepy, type);
733}
734
735// --------------------------------------------------------------------------
736//
737// Get a rotation matrix which aligns the pointing position
738// to the center of the x,y plain
739//
740TRotation MAstroCatalog::AlignCoordinates(const TVector3 &v) const
741{
742 TRotation trans;
743 trans.RotateZ(-v.Phi());
744 trans.RotateY(-v.Theta());
745 trans.RotateZ(-TMath::Pi()/2);
746 return trans;
747}
748
749// --------------------------------------------------------------------------
750//
751// Return the rotation matrix which converts either sky or
752// local coordinates to coordinates which pole is the current
753// pointing direction.
754//
755TRotation MAstroCatalog::GetGrid(Bool_t local)
756{
757 const Bool_t enable = fTime && fObservatory;
758
759 // If sky coordinate view is requested get rotation matrix and
760 // draw corresponding sky-grid and if possible local grid
761 if (!local)
762 {
763 const TRotation trans(AlignCoordinates(fRaDec));
764
765 DrawGrid(fRaDec, trans, 1);
766
767 if (enable)
768 {
769 const MAstroSky2Local rot(*fTime, *fObservatory);
770 DrawGrid(rot*fRaDec, trans*rot.Inverse(), 2);
771 }
772
773 // Return the correct rotation matrix
774 return trans;
775 }
776
777 // If local coordinate view is requested get rotation matrix and
778 // draw corresponding sky-grid and if possible local grid
779 if (local && enable)
780 {
781 const MAstroSky2Local rot(*fTime, *fObservatory);
782
783 const TRotation trans(AlignCoordinates(rot*fRaDec));
784
785 DrawGrid(fRaDec, trans*rot, 1);
786 DrawGrid(rot*fRaDec, trans, 2);
787
788 // Return the correct rotation matrix
789 return trans*rot;
790 }
791
792 return TRotation();
793}
794
795// --------------------------------------------------------------------------
796//
797// Create the title for the pad.
798//
799TString MAstroCatalog::GetPadTitle() const
800{
801 const Double_t ra = fRaDec.Phi()*TMath::RadToDeg();
802 const Double_t dec = (TMath::Pi()/2-fRaDec.Theta())*TMath::RadToDeg();
803
804 TString txt;
805 txt += Form("\\alpha=%.2fh ", fmod(ra/15+48, 24));
806 txt += Form("\\delta=%.1f\\circ ", fmod(dec+270,180)-90);
807 txt += Form("/ FOV=%.1f\\circ", fRadiusFOV);
808
809 if (!fTime || !fObservatory)
810 return txt;
811
812 const MAstroSky2Local rot(*fTime, *fObservatory);
813 const TVector3 loc = rot*fRaDec;
814
815 const Double_t rho = rot.RotationAngle(fRaDec.Phi(), TMath::Pi()/2-fRaDec.Theta());
816
817 const Double_t zd = TMath::RadToDeg()*loc.Theta();
818 const Double_t az = TMath::RadToDeg()*loc.Phi();
819
820 txt.Prepend("#splitline{");
821 txt += Form(" \\theta=%.1f\\circ ", fmod(zd+270,180)-90);
822 txt += Form("\\phi=%.1f\\circ ", fmod(az+720, 360));
823 txt += Form(" / \\rho=%.1f\\circ", rho*TMath::RadToDeg());
824 txt += "}{<";
825 txt += fTime->GetSqlDateTime();
826 txt += ">}";
827 return txt;
828}
829
830// --------------------------------------------------------------------------
831//
832// To overlay the catalog make sure, that in any case you are using
833// the 'same' option.
834//
835// If you want to overlay this on top of any picture which is created
836// by derotation of the camera plain you have to use the 'mirror' option
837// the compensate the mirroring of the image in the camera plain.
838//
839// If you have already compensated this by x=-x and y=-y when creating
840// the histogram you can simply overlay the catalog.
841//
842// To overlay the catalog on a 2D histogram the histogram must have
843// units of degrees (which are plain, like you directly convert the
844// camera units by multiplication to degrees)
845//
846// To be 100% exact you must use the option 'plain' which assumes a plain
847// screen. This is not necessary for the MAGIC-camera because the
848// difference between both is less than 1e-3.
849//
850// You should always be aware of the fact, that the shown stars and the
851// displayed grid is the ideal case, like a reflection on a virtual
852// perfectly aligned central mirror. In reality the star-positions are
853// smeared to the edge of the camera the more the distance to the center
854// is, such that the center of gravity of the light distribution might
855// be more far away from the center than the display shows.
856//
857//
858void MAstroCatalog::AddPrimitives(TString o)
859{
860 const Bool_t same = o.Contains("same", TString::kIgnoreCase);
861 const Bool_t local = o.Contains("local", TString::kIgnoreCase);
862 const Bool_t mirx = o.Contains("mirrorx", TString::kIgnoreCase);
863 const Bool_t miry = o.Contains("mirrory", TString::kIgnoreCase);
864 const Bool_t mirror = o.Contains("mirror", TString::kIgnoreCase) && !mirx && !miry;
865
866 // X is vice versa, because ra is defined anti-clockwise
867 mirx || mirror ? ResetBit(kMirrorX) : SetBit(kMirrorX);
868 miry || mirror ? SetBit(kMirrorY) : ResetBit(kMirrorY);
869
870 const TRotation rot(GetGrid(local));
871
872 TIter Next(&fList);
873 TVector3 *v=0;
874 while ((v=(TVector3*)Next()))
875 {
876 // FIXME: Check Magnitude!
877 TVector2 s(v->Phi(), v->Theta());
878 if (Convert(rot, s)==kTRUE)
879 DrawStar(s.X(), s.Y(), *v, kFALSE);
880 }
881
882 if (!same)
883 {
884 TPaveText *pv = new TPaveText(0.01, 0.90, 0.63, 0.99, "brNDC");
885 pv->AddText(GetPadTitle());
886 AddMap(pv);
887 }
888
889 TMarker *mk=new TMarker(0, 0, kMultiply);
890 mk->SetMarkerColor(kBlack);
891 mk->SetMarkerSize(1.5);
892 AddMap(mk);
893}
894
895// --------------------------------------------------------------------------
896//
897// Do nothing if 'same' option given.
898// Otherwise set pad-range such that x- and y- coordinates have the same
899// step-size
900//
901void MAstroCatalog::SetRangePad(Option_t *o)
902{
903 if (TString(o).Contains("same", TString::kIgnoreCase))
904 return;
905
906 const Double_t edge = fRadiusFOV/TMath::Sqrt(2.);
907 gPad->Range(-edge, -edge, edge, edge);
908
909 const Float_t w = gPad->GetWw();
910 const Float_t h = gPad->GetWh();
911
912 if (w<h)
913 gPad->Range(-edge, -edge*h/w, edge, edge*h/w);
914 else
915 gPad->Range(-edge*w/h, -edge, edge*w/h, edge);
916}
917
918// --------------------------------------------------------------------------
919//
920// First delete all gui elements.
921// Set the correct range of the pad.
922// Create all gui primitives
923// If "this" is not in the current pad add it to the current pad.
924// Reset vit kHasChanged
925//
926void MAstroCatalog::DrawPrimitives(Option_t *o)
927{
928 DeleteMap();
929
930 SetRangePad(o);
931
932#ifdef DEBUG
933 TStopwatch clk;
934 clk.Start();
935#endif DEBUG
936 AddPrimitives(o);
937#ifdef DEBUG
938 clk.Stop();
939 clk.Print();
940#endif DEBUG
941
942 // Append to a possible second pad
943 if (!gPad->GetListOfPrimitives()->FindObject(this))
944 AppendPad(o);
945
946 ResetBit(kHasChanged);
947}
948
949// --------------------------------------------------------------------------
950//
951// Call Paint() of all gui elements
952//
953void MAstroCatalog::PaintMap()
954{
955 Long_t key, val;
956 TExMapIter map(&fMapG);
957 while (map.Next(key, val))
958 ((TObject*)key)->Paint();
959}
960
961// --------------------------------------------------------------------------
962//
963// Append "this" to current pad
964// set bit kHasChanged to recreate all gui elements
965// Connect signal
966//
967void MAstroCatalog::Draw(Option_t *o)
968{
969 // Append to first pad
970 AppendPad(o);
971
972 // If contents have not previously changed make sure that
973 // all primitives are recreated.
974 SetBit(kHasChanged);
975
976 // Connect all TCanvas::ProcessedEvent to this->EventInfo
977 // This means, that after TCanvas has processed an event
978 // EventInfo of this class is called, see TCanvas::HandleInput
979 gPad->GetCanvas()->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)",
980 "MAstroCatalog", this,
981 "EventInfo(Int_t,Int_t,Int_t,TObject*)");
982}
983
984// --------------------------------------------------------------------------
985//
986// This function was connected to all created canvases. It is used
987// to redirect GetObjectInfo into our own status bar.
988//
989// The 'connection' is done in AddTab
990//
991void MAstroCatalog::EventInfo(Int_t event, Int_t px, Int_t py, TObject *selected)
992{
993 TCanvas *c = (TCanvas*)gTQSender;
994
995 gPad = c ? c->GetSelectedPad() : NULL;
996 if (!gPad)
997 return;
998
999 // Try to find a corresponding object with kCannotPick set and
1000 // an available TString (for a tool tip)
1001 TString *str=0;
1002 if (!selected || selected==this)
1003 {
1004 Long_t key, val;
1005 TExMapIter map(&fMapG);
1006 while (map.Next(key, val))
1007 {
1008 if (!val)
1009 continue;
1010
1011 TObject *o=(TObject*)key;
1012 if (o->DistancetoPrimitive(px, py)>TPad::GetMaxPickDistance())
1013 continue;
1014
1015 selected = o;
1016 str = (TString*)val;
1017 break;
1018 }
1019 }
1020
1021 if (!selected)
1022 return;
1023
1024 // Handle some gui events
1025 switch (event)
1026 {
1027 case kMouseMotion:
1028 if (!fToolTip->IsMapped() && str)
1029 ShowToolTip(px, py, *str);
1030 break;
1031
1032 case kMouseLeave:
1033 if (fToolTip->IsMapped())
1034 fToolTip->Hide();
1035 break;
1036
1037 case kKeyPress:
1038 ExecuteEvent(kKeyPress, px, py);
1039 break;
1040 }
1041}
1042
1043// --------------------------------------------------------------------------
1044//
1045// Handle keyboard events.
1046//
1047void MAstroCatalog::ExecuteEventKbd(Int_t keycode, Int_t keysym)
1048{
1049 Double_t dra =0;
1050 Double_t ddec=0;
1051
1052 switch (keysym)
1053 {
1054 case kKey_Left:
1055 dra = -TMath::DegToRad();
1056 break;
1057 case kKey_Right:
1058 dra = +TMath::DegToRad();
1059 break;
1060 case kKey_Up:
1061 ddec = +TMath::DegToRad();
1062 break;
1063 case kKey_Down:
1064 ddec = -TMath::DegToRad();
1065 break;
1066 case kKey_Plus:
1067 SetRadiusFOV(fRadiusFOV+1);
1068 break;
1069 case kKey_Minus:
1070 SetRadiusFOV(fRadiusFOV-1);
1071 break;
1072
1073 default:
1074 return;
1075 }
1076
1077 const Double_t r = fRaDec.Phi();
1078 const Double_t d = TMath::Pi()/2-fRaDec.Theta();
1079
1080 SetRaDec(r+dra, d+ddec);
1081
1082 gPad->Update();
1083}
1084
1085// ------------------------------------------------------------------------
1086//
1087// Execute a gui event on the camera
1088//
1089void MAstroCatalog::ExecuteEvent(Int_t event, Int_t mp1, Int_t mp2)
1090{
1091 if (!TestBit(kGuiActive))
1092 return;
1093
1094 if (event==kKeyPress)
1095 ExecuteEventKbd(mp1, mp2);
1096}
1097
1098// --------------------------------------------------------------------------
1099//
1100// Calculate distance to primitive by checking all gui elements
1101//
1102Int_t MAstroCatalog::DistancetoPrimitive(Int_t px, Int_t py)
1103{
1104 Int_t min = INT_MAX;
1105
1106 Long_t key, val;
1107 TExMapIter map(&fMapG);
1108 while (map.Next(key, val))
1109 {
1110 TObject *o=(TObject*)key;
1111
1112 const Int_t d = o->DistancetoPrimitive(px, py);
1113
1114 if (d<TPad::GetMaxPickDistance())
1115 return 0;
1116
1117 if (d<min)
1118 min=d;
1119 }
1120
1121 return min;
1122}
1123
1124// --------------------------------------------------------------------------
1125void MAstroCatalog::ShowToolTip(Int_t px, Int_t py, const char *txt)
1126{
1127 Int_t x=0;
1128 Int_t y=0;
1129
1130 const Window_t id1 = gVirtualX->GetWindowID(gPad->GetCanvasID());
1131 const Window_t id2 = fToolTip->GetParent()->GetId();
1132
1133 Window_t id3;
1134 gVirtualX->TranslateCoordinates(id1, id2, px, py, x, y, id3);
1135
1136 // Show tool tip
1137 fToolTip->SetText(txt);
1138 fToolTip->Show(x+4, y+4);
1139}
1140
1141/*
1142void MAstroCatalog::RecursiveRemove(TObject *obj)
1143{
1144 ULong_t hash;
1145 Long_t key, val;
1146
1147 TExMapIter map(&fMapG);
1148 while (map.Next(hash, key, val))
1149 {
1150 if (key != (Long_t)obj)
1151 continue;
1152
1153 fMapG.Remove(hash, key);
1154 delete (TObject*)(key);
1155 if (val)
1156 delete (TString*)(val);
1157 break;
1158 }
1159}
1160*/
Note: See TracBrowser for help on using the repository browser.