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

Last change on this file since 4240 was 4213, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 39.0 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//
30// THIS IMPLEMENTATION IS PRELIMINARY AND WILL BE MERGED WITH
31// SOME PARTS OF THE DRIVE SOFTWARE SOON!
32//
33//
34// Usage:
35// ------
36//
37// To display a starfield you must have a supported catalog, then do:
38//
39// MTime time;
40// // Time for which to get the picture
41// time.Set(2004, 2, 28, 20, 14, 7);
42// // Current observatory
43// MObservatory magic1;
44// // Right Ascension [h] and declination [deg] of source
45// // Currently 'perfect' pointing is assumed
46// const Double_t ra = MAstro::Hms2Rad(5, 34, 31.9);
47// const Double_t dec = MAstro::Dms2Rad(22, 0, 52.0);
48// MAstroCatalog stars;
49// // Magnitude up to which the stars are loaded from the catalog
50// stars.SetLimMag(6);
51// // Radius of FOV around the source position to load the stars
52// stars.SetRadiusFOV(3);
53// // Source position
54// stars.SetRaDec(ra, dec);
55// // Catalog to load (here: Bright Star Catalog V5)
56// stars.ReadBSC("bsc5.dat");
57// // Obersavatory and time to also get local coordinate information
58// stars.SetObservatory(magic1);
59// stars.SetTime(time);
60// // Enable interactive GUI
61// stars.SetGuiActive();
62// //Clone the catalog due to the validity range of the instance
63// TObject *o = stars.Clone();
64// o->SetBit(kCanDelete);
65// o->Draw();
66//
67// If no time and/or Obervatory location is given no local coordinate
68// information is displayed.
69//
70//
71// Coordinate Transformation:
72// -------------------------
73// The conversion from sky coordinates to local coordinates is done using
74// MAstroSky2Local which does a simple rotation of the coordinate system.
75// This is inaccurate in the order of 30arcsec due to ignorance of all
76// astrometrical corrections (nutation, precission, abberation, ...)
77//
78//
79// GUI:
80// ----
81// * If the gui is interactive you can use the cursor keys to change
82// the position you are looking at and with plus/minus you
83// can (un)zoom the FOV (Field Of View)
84// * The displayed values mean the following:
85// + alpha: Right Ascension
86// + delta: Declination
87// + theta: zenith distance / zenith angle
88// + phi: azimuth angle
89// + rho: angle of rotation sky-coordinate system vs local-
90// coordinate system
91// + time of display
92// * Move the mouse on top of the grid points or the stars to get
93// more setailed information.
94// * Enable the event-info in a canvas to see the current
95// ObjectInfo=tooltip-text
96// * call SetNoToolTips to supress the tooltips
97// * the blue lines are the local coordinat system
98// * the red lines are sky coordinate system
99//
100//
101// ToDo:
102// -----
103// - replace MVetcor3 by a more convinient class. Maybe use TExMap, too.
104// - change tooltips to multi-line tools tips as soon as root
105// supports them
106// - a derived class is missing which supports all astrometrical
107// correction (base on slalib and useable in Cosy)
108//
109//////////////////////////////////////////////////////////////////////////////
110#include "MAstroCatalog.h"
111
112#include <errno.h> // strerror
113#include <fstream> // ifstream
114#include <stdlib.h> // ati, atof
115#include <limits.h> // INT_MAX (Suse 7.3/gcc 2.95)
116
117#include <KeySymbols.h> // kKey_*
118
119#include <TLine.h> // TLine
120#include <TMarker.h> // TMarker
121#include <TCanvas.h> // TCanvas
122#include <TArrayI.h> // TArrayI
123#include <TGToolTip.h> // TGToolTip
124#include <TPaveText.h> // TPaveText
125
126#include "MLog.h"
127#include "MLogManip.h"
128
129#include "MTime.h"
130#include "MString.h"
131#include "MAstro.h"
132#include "MAstroSky2Local.h"
133#include "MObservatory.h"
134
135#undef DEBUG
136//#define DEBUG
137
138#ifdef DEBUG
139#include <TStopwatch.h>
140#endif
141
142ClassImp(MVector3);
143ClassImp(MAstroCatalog);
144
145using namespace std;
146
147void MVector3::WriteBinary(ostream &out) const
148{
149 const Double_t t = Theta();
150 const Double_t p = Phi();
151 const Float_t m = Mag();
152 out.write((char*)&t, 8);
153 out.write((char*)&p, 8);
154 out.write((char*)&m, 4);
155 out << fName << endl;
156}
157
158void MVector3::ReadBinary(istream &in)
159{
160 Double_t t, p;
161 Float_t m;
162
163 in.read((char*)&t, 8);
164 in.read((char*)&p, 8);
165 in.read((char*)&m, 4);
166 fName.ReadLine(in);
167
168 SetMagThetaPhi(m, t, p);
169}
170
171// --------------------------------------------------------------------------
172//
173// Default Constructor. Set Default values:
174// fLimMag = 99
175// fRadiusFOV = 90
176//
177MAstroCatalog::MAstroCatalog() : fLimMag(99), fRadiusFOV(90), fToolTip(0), fObservatory(0), fTime(0)
178{
179 fList.SetOwner();
180 fMapG.SetOwner();
181 fToolTip = gROOT->IsBatch() || !gClient ? 0 : new TGToolTip(0, "", 0);
182}
183
184// --------------------------------------------------------------------------
185//
186// Destructor. Delete fTime, fObservatory. Disconnect signal. delete tooltip.
187// Delete Map with gui primitives
188//
189MAstroCatalog::~MAstroCatalog()
190{
191 // First disconnect the EventInfo...
192 // FIXME: There must be an easier way!
193 TIter Next(gROOT->GetListOfCanvases());
194 TCanvas *c;
195 while ((c=(TCanvas*)Next()))
196 c->Disconnect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", this,
197 "EventInfo(Int_t,Int_t,Int_t,TObject*)");
198
199 // Now delete the data members
200 if (fTime)
201 delete fTime;
202 if (fObservatory)
203 delete fObservatory;
204
205 if (fToolTip)
206 {
207 fToolTip->Hide();
208 delete fToolTip;
209 }
210}
211
212// --------------------------------------------------------------------------
213//
214// Snippet to for reading ctalog files.
215//
216TString MAstroCatalog::FindToken(TString &line, Char_t tok)
217{
218 Ssiz_t token = line.First(tok);
219 if (token<0)
220 {
221 const TString copy(line);
222 line = "";
223 return copy;
224 }
225
226 const TString res = line(0, token);
227 line.Remove(0, token+1);
228 return res;
229}
230
231// --------------------------------------------------------------------------
232//
233// return int correspoding to TSubString
234//
235Int_t MAstroCatalog::atoi(const TSubString &sub)
236{
237 return atoi(TString(sub));
238}
239
240// --------------------------------------------------------------------------
241//
242// return float correspoding to TSubString
243//
244Float_t MAstroCatalog::atof(const TSubString &sub)
245{
246 return atof(TString(sub));
247}
248
249// --------------------------------------------------------------------------
250//
251// return int correspoding to TString
252//
253Int_t MAstroCatalog::atoi(const TString &s)
254{
255 return std::atoi(s);
256}
257
258// --------------------------------------------------------------------------
259//
260// return float correspoding to TString
261//
262Float_t MAstroCatalog::atof(const TString &s)
263{
264 return std::atof(s);
265}
266
267// --------------------------------------------------------------------------
268//
269// Read from a xephem catalog, set bit kHasChahanged.
270// Already read data is not deleted. To delete the stored data call
271// Delete().
272//
273Int_t MAstroCatalog::ReadXephem(TString catalog)
274{
275 SetBit(kHasChanged);
276
277 gLog << inf << "Reading Xephem catalog: " << catalog << endl;
278
279 ifstream fin(catalog);
280 if (!fin)
281 {
282 gLog << err << "Cannot open catalog file " << catalog << ": ";
283 gLog << strerror(errno) << endl;
284 return 0;
285 }
286
287 Int_t add =0;
288 Int_t cnt =0;
289 Int_t line=0;
290
291 Double_t maxmag=0;
292
293 while (1)
294 {
295 TString row;
296 row.ReadLine(fin);
297 if (!fin)
298 break;
299
300 line++;
301
302 if (row[0]=='#')
303 continue;
304
305 TString line(row);
306
307 TString name = FindToken(line);
308 TString dummy = FindToken(line);
309 TString r = FindToken(line);
310 TString d = FindToken(line);
311 TString m = FindToken(line);
312 TString epoch = FindToken(line);
313
314 if (name.IsNull() || r.IsNull() || d.IsNull() || m.IsNull() || epoch.IsNull())
315 {
316 gLog << warn << "Invalid Entry Line #" << line << ": " << row << endl;
317 continue;
318 }
319
320 cnt++;
321
322 const Double_t mag = atof(m);
323
324 maxmag = TMath::Max(maxmag, mag);
325
326 if (mag>fLimMag)
327 continue;
328
329 if (epoch!="2000")
330 {
331 gLog << warn << "Epoch != 2000... skipped." << endl;
332 continue;
333 }
334
335 Double_t ra0, dec0;
336 MAstro::Coordinate2Angle(r, ra0);
337 MAstro::Coordinate2Angle(d, dec0);
338
339 ra0 *= TMath::Pi()/12;
340 dec0 *= TMath::Pi()/180;
341
342 MVector3 *star=new MVector3;
343 star->SetRaDec(ra0, dec0, mag);
344 star->SetName(name);
345 if (star->Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
346 {
347 delete star;
348 continue;
349 }
350
351 fList.Add(star);
352 add++;
353 }
354 gLog << inf << "Read " << add << " out of " << cnt << " (Total max mag=" << maxmag << ")" << endl;
355
356 return add;
357}
358
359// --------------------------------------------------------------------------
360//
361// Read from a NGC2000 catalog. set bit kHasChanged
362// Already read data is not deleted. To delete the stored data call
363// Delete().
364//
365Int_t MAstroCatalog::ReadNGC2000(TString catalog)
366{
367 SetBit(kHasChanged);
368
369 gLog << inf << "Reading NGC2000 catalog: " << catalog << endl;
370
371 ifstream fin(catalog);
372 if (!fin)
373 {
374 gLog << err << "Cannot open catalog file " << catalog << ": ";
375 gLog << strerror(errno) << endl;
376 return 0;
377 }
378
379 Int_t add=0;
380 Int_t cnt=0;
381 Int_t n =0;
382
383 Double_t maxmag=0;
384
385 while (1)
386 {
387 TString row;
388 row.ReadLine(fin);
389 if (!fin)
390 break;
391
392 cnt++;
393
394 const Int_t rah = atoi(row(13, 2));
395 const Float_t ram = atof(row(16, 4));
396 const Char_t decs = row(22);
397 const Int_t decd = atoi(row(23, 2));
398 const Int_t decm = atoi(row(26, 2));
399 const TString m = row(43, 4);
400
401 if (m.Strip().IsNull())
402 continue;
403
404 n++;
405
406 const Double_t mag = atof(m);
407
408 maxmag = TMath::Max(maxmag, mag);
409
410 if (mag>fLimMag)
411 continue;
412
413 const Double_t ra = MAstro::Hms2Rad(rah, (int)ram, fmod(ram, 1)*60);
414 const Double_t dec = MAstro::Dms2Rad(decd, decm, 0, decs);
415
416 MVector3 *star=new MVector3;
417 star->SetRaDec(ra, dec, mag);
418 star->SetName(row(0, 8));
419 if (star->Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
420 {
421 delete star;
422 continue;
423 }
424
425 fList.Add(star);
426 add++;
427 }
428
429 gLog << inf << "Read " << add << " out of " << n << " (Total max mag=" << maxmag << ")" << endl;
430
431 return add;
432}
433
434// --------------------------------------------------------------------------
435//
436// Read from a Bright Star Catalog catalog. set bit kHasChanged
437// Already read data is not deleted. To delete the stored data call
438// Delete().
439//
440Int_t MAstroCatalog::ReadBSC(TString catalog)
441{
442 SetBit(kHasChanged);
443
444 gLog << inf << "Reading Bright Star Catalog (BSC5) catalog: " << catalog << endl;
445
446 ifstream fin(catalog);
447 if (!fin)
448 {
449 gLog << err << "Cannot open catalog file " << catalog << ": ";
450 gLog << strerror(errno) << endl;
451 return 0;
452 }
453
454 Int_t add=0;
455 Int_t cnt=0;
456 Int_t n =0;
457
458 Double_t maxmag=0;
459
460 while (1)
461 {
462 TString row;
463 row.ReadLine(fin);
464 if (!fin)
465 break;
466
467 cnt++;
468
469 const Int_t rah = atoi(row(75, 2));
470 const Int_t ram = atoi(row(77, 2));
471 const Float_t ras = atof(row(79, 4));
472 const Char_t decsgn = row(83);
473 const Int_t decd = atoi(row(84, 2));
474 const Int_t decm = atoi(row(86, 2));
475 const Int_t decs = atoi(row(88, 2));
476 const TString m = row(102, 5);
477
478 if (m.Strip().IsNull())
479 continue;
480
481 n++;
482
483 const Double_t mag = atof(m.Data());
484
485 maxmag = TMath::Max(maxmag, mag);
486
487 if (mag>fLimMag)
488 continue;
489
490 const Double_t ra = MAstro::Hms2Rad(rah, ram, ras);
491 const Double_t dec = MAstro::Dms2Rad(decd, decm, decs, decsgn);
492
493 MVector3 *star=new MVector3;
494 star->SetRaDec(ra, dec, mag);
495 star->SetName(row(4,9));
496 if (star->Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
497 {
498 delete star;
499 continue;
500 }
501
502 fList.Add(star);
503 add++;
504 }
505
506 gLog << inf << "Read " << add << " out of " << n << " (Total max mag=" << maxmag << ")" << endl;
507
508 return add;
509}
510
511// --------------------------------------------------------------------------
512//
513// Read from a ascii heasarc ppm catalog. set bit kHasChanged
514// Already read data is not deleted. To delete the stored data call
515// Delete().
516// If the second argument is given all survived stars are written
517// to a file outname. This files will contain an apropriate compressed
518// file format. You can read such files again using ReadCompressed.
519//
520Int_t MAstroCatalog::ReadHeasarcPPM(TString catalog, TString outname)
521{
522 SetBit(kHasChanged);
523
524 gLog << inf << "Reading Heasarc PPM catalog: " << catalog << endl;
525
526 ifstream fin(catalog);
527 if (!fin)
528 {
529 gLog << err << "Cannot open catalog file " << catalog << ": ";
530 gLog << strerror(errno) << endl;
531 return 0;
532 }
533
534 ofstream *fout = outname.IsNull() ? 0 : new ofstream(outname);
535 if (fout && !*fout)
536 {
537 gLog << warn << "Cannot open output file " << outname << ": ";
538 gLog << strerror(errno) << endl;
539 delete fout;
540 fout = 0;
541 }
542
543 Int_t add=0;
544 Int_t cnt=0;
545 Int_t n =0;
546
547 Double_t maxmag=0;
548
549 while (1)
550 {
551 TString row;
552 row.ReadLine(fin);
553 if (!fin)
554 break;
555
556 cnt++;
557
558 if (!row.BeginsWith("PPM "))
559 continue;
560
561 const TString name = row(0, row.First('|'));
562 row = row(row.First('|')+1, row.Length());
563 row = row(row.First('|')+1, row.Length());
564
565 const TString vmag = row(0, row.First('|'));
566
567 n++;
568 const Double_t mag = atof(vmag.Data());
569 maxmag = TMath::Max(maxmag, mag);
570 if (mag>fLimMag)
571 continue;
572
573 row = row(row.First('|')+1, row.Length());
574 row = row(row.First('|')+1, row.Length());
575
576 row = row(row.First('|')+1, row.Length());
577 row = row(row.First('|')+1, row.Length());
578
579 row = row(row.First('|')+1, row.Length());
580 row = row(row.First('|')+1, row.Length());
581
582 const TString ra = row(0, row.First('|'));
583 row = row(row.First('|')+1, row.Length());
584 const TString de = row(0, row.First('|'));
585 row = row(row.First('|')+1, row.Length());
586
587 Char_t sgn;
588 Int_t d, m;
589 Float_t s;
590 if (sscanf(ra.Data(), "%d %d %f", &d, &m, &s)!=3)
591 {
592 // gLog << "Error loading entry in line " << i << endl;
593 continue;
594 }
595 const Double_t ra0 = MAstro::Hms2Rad(d, m, s);
596
597 if (sscanf(de.Data(), "%c%d %d %f", &sgn, &d, &m, &s)!=4)
598 {
599 // gLog << "Error loading entry in line " << i << endl;
600 continue;
601 }
602 const Double_t de0 = MAstro::Dms2Rad(d, m, s, sgn);
603
604 MVector3 *star=new MVector3;
605 star->SetRaDec(ra0, de0, mag);
606 star->SetName(name);
607 if (star->Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
608 {
609 delete star;
610 continue;
611 }
612
613 fList.Add(star);
614 add++;
615
616 if (fout)
617 star->WriteBinary(*fout);
618 }
619
620 gLog << inf << "Read " << add << " out of " << n << " (Total max mag=" << maxmag << ")" << endl;
621
622 return add;
623}
624
625// --------------------------------------------------------------------------
626//
627// Read from a MAstroCatalog compressed catalog. set bit kHasChanged
628// Already read data is not deleted. To delete the stored data call
629// Delete().
630//
631Int_t MAstroCatalog::ReadCompressed(TString catalog)
632{
633 SetBit(kHasChanged);
634
635 gLog << inf << "Reading MAstroCatalog compressed catalog: " << catalog << endl;
636
637 ifstream fin(catalog);
638 if (!fin)
639 {
640 gLog << err << "Cannot open catalog file " << catalog << ": ";
641 gLog << strerror(errno) << endl;
642 return 0;
643 }
644
645 Int_t add=0;
646 Int_t cnt=0;
647 Int_t n =0;
648
649 Double_t maxmag=0;
650
651 MVector3 entry;
652
653 while (1)
654 {
655 cnt++;
656
657 entry.ReadBinary(fin);
658 if (!fin)
659 break;
660
661 n++;
662
663 const Double_t mag = entry.Magnitude();
664 maxmag = TMath::Max(maxmag, mag);
665 if (mag>fLimMag)
666 continue;
667
668 if (entry.Angle(fRaDec)*TMath::RadToDeg()>fRadiusFOV)
669 continue;
670
671 fList.Add(entry.Clone());
672 add++;
673 }
674
675 gLog << inf << "Read " << add << " out of " << n << " (Total max mag=" << maxmag << ")" << endl;
676
677 return add;
678}
679
680// --------------------------------------------------------------------------
681//
682// Set Range of pad. If something has changed create and draw new primitives.
683// Paint all gui primitives.
684//
685void MAstroCatalog::Paint(Option_t *o)
686{
687 SetRangePad(o);
688
689 if (TestBit(kHasChanged))
690 DrawPrimitives(o);
691
692 fMapG.Paint();
693}
694
695// --------------------------------------------------------------------------
696//
697// Set Range of pad if pad available. If something has changed create
698// and draw new primitives. Paint all gui primitives to the Drawable with
699// Id id. This can be used to be able to
700//
701/*
702void MAstroCatalog::PaintImg(Int_t id, Option_t *o)
703{
704 if (gPad)
705 SetRangePad(o);
706
707 if (TestBit(kHasChanged))
708 {
709 if (id>0)
710 gPad=0;
711 DrawPrimitives(o);
712 }
713
714 fMapG.Paint(id, fRadiusFOV);
715}
716*/
717
718// --------------------------------------------------------------------------
719//
720// Set Range of pad. If something has changed create and draw new primitives.
721// Paint all gui primitives.
722//
723// Because in some kind of multi-threaded environments gPad doesn't stay
724// the same in a single thread (because it might be changed in the same
725// thread inside a gui updating timer for example) we have to secure the
726// usage of gPad with a bit. This is also not multi-thread safe against
727// calling this function, but the function should work well in multi-
728// threaded environments. Never call this function from different threads
729// simultaneously.
730//
731void MAstroCatalog::PaintImg(unsigned char *buf, int w, int h, Option_t *o)
732{
733 if (!o)
734 o = "local mirrorx yellow * =";
735
736 if (TestBit(kHasChanged))
737 {
738 SetBit(kDrawingImage);
739 DrawPrimitives(o);
740 ResetBit(kDrawingImage);
741 }
742
743 fMapG.Paint(buf, w, h, fRadiusFOV);
744}
745
746// --------------------------------------------------------------------------
747//
748// Draw a black marker at the position of the star. Create a corresponding
749// tooltip with the coordinates.
750// x, y: Pad Coordinates to draw star
751// v: Sky position (Ra/Dec) of the star
752// col: Color of marker (<0 mean transparent)
753// txt: additional tooltip text
754// resize: means resize the marker according to the magnitude
755//
756void MAstroCatalog::DrawStar(Double_t x, Double_t y, const TVector3 &v, Int_t col, const char *txt, Bool_t resize)
757{
758 const Double_t ra = v.Phi()*TMath::RadToDeg()/15;
759 const Double_t dec = (TMath::Pi()/2-v.Theta())*TMath::RadToDeg();
760
761 const Double_t mag = -2.5*log10(v.Mag());
762
763 MString str;
764 str.Print("%s: Ra=%.2fh Dec=%.1fd Mag=%.1f", v.GetName(), ra, dec, mag);
765 if (txt)
766 {
767 str += " (";
768 str += txt;
769 str += ")";
770 }
771
772 // draw star on the camera display
773 TMarker *tip=new TMarker(x, y, col<0 ? kDot : kFullDotMedium);;
774 tip->SetMarkerColor(col);
775 fMapG.Add(tip, new TString(str));
776
777 if (!resize)
778 return;
779
780 tip->SetMarkerSize((10 - (mag>1 ? mag : 1))/15);
781 tip->SetMarkerStyle(kCircle);
782}
783
784// --------------------------------------------------------------------------
785//
786// Set pad as modified.
787//
788void MAstroCatalog::Update(Bool_t upd)
789{
790 SetBit(kHasChanged);
791 if (gPad && TestBit(kMustCleanup))
792 {
793 gPad->Modified();
794 if (upd)
795 gPad->Update();
796 }
797}
798
799// --------------------------------------------------------------------------
800//
801// Set the observation time. Necessary to use local coordinate
802// system. The MTime object is cloned.
803//
804void MAstroCatalog::SetTime(const MTime &time)
805{
806 if (fTime)
807 delete fTime;
808 fTime=(MTime*)time.Clone();
809}
810
811// --------------------------------------------------------------------------
812//
813// Set the observatory location. Necessary to use local coordinate
814// system. The MObservatory object is cloned.
815//
816void MAstroCatalog::SetObservatory(const MObservatory &obs)
817{
818 if (fObservatory)
819 delete fObservatory;
820 fObservatory=new MObservatory;
821 obs.Copy(*fObservatory);
822}
823
824// --------------------------------------------------------------------------
825//
826// Convert the vector to pad coordinates. After conversion
827// the x- coordinate of the vector must be the x coordinate
828// of the pad - the same for y. If the coordinate is inside
829// the current draw area return kTRUE, otherwise kFALSE.
830// If it is an invalid coordinate return kERROR
831//
832Int_t MAstroCatalog::ConvertToPad(const TVector3 &w0, TVector2 &v) const
833{
834 TVector3 w(w0);
835
836 // Stretch such, that the Z-component is alwas the same. Now
837 // X and Y contains the intersection point between the star-light
838 // and the plain of a virtual plain screen (ccd...)
839 if (TestBit(kPlainScreen))
840 w *= 1./w(2);
841
842 w *= TMath::RadToDeg(); // FIXME: *conversion factor?
843 v.Set(TestBit(kMirrorX) ? -w(0) : w(0),
844 TestBit(kMirrorY) ? -w(1) : w(1));
845
846 if (w(2)<0)
847 return kERROR;
848
849 if (TestBit(kDrawingImage) || !gPad)
850 return v.Mod2()<fRadiusFOV*fRadiusFOV;
851
852 return v.X()>gPad->GetX1() && v.Y()>gPad->GetY1() &&
853 v.X()<gPad->GetX2() && v.Y()<gPad->GetY2();
854}
855
856// --------------------------------------------------------------------------
857//
858// Convert theta/phi coordinates of v by TRotation into new coordinate
859// system and convert the coordinated to pad by ConvertToPad.
860// The result is retunred in v.
861//
862Int_t MAstroCatalog::Convert(const TRotation &rot, TVector2 &v) const
863{
864 MVector3 w;
865 w.SetMagThetaPhi(1, v.Y(), v.X());
866 w *= rot;
867
868 return ConvertToPad(w, v);
869}
870
871// --------------------------------------------------------------------------
872//
873// Draw a line from v to v+(dx,dy) using Convert/ConvertToPad to get the
874// corresponding pad coordinates.
875//
876Bool_t MAstroCatalog::DrawLine(const TVector2 &v, Int_t dx, Int_t dy, const TRotation &rot, Int_t type)
877{
878 const TVector2 add(dx*TMath::DegToRad(), dy*TMath::DegToRad());
879
880 // Define all lines in the same direction
881 const TVector2 va(dy==1?v:v+add);
882 const TVector2 vb(dy==1?v+add:v);
883
884 TVector2 v0(va);
885 TVector2 v1(vb);
886
887 const Int_t rc0 = Convert(rot, v0);
888 const Int_t rc1 = Convert(rot, v1);
889
890 // Both are kFALSE or both are kERROR
891 if ((rc0|rc1)==kFALSE || (rc0&rc1)==kERROR)
892 return kFALSE;
893
894 TLine *line = new TLine(v0.X(), v0.Y(), v1.X(), v1.Y());
895 line->SetLineStyle(kDashDotted); //kDashed, kDotted, kDashDotted
896 line->SetLineColor(kWhite+type*2);
897 fMapG.Add(line);
898
899 if (dx!=0)
900 return kTRUE;
901
902 const TVector2 deg = va*TMath::RadToDeg();
903
904 MString txt;
905 if (type==1)
906 txt.Print("Ra=%.2fh Dec=%.1fd", fmod(deg.X()/15+48, 24), fmod(90-deg.Y()+270,180)-90);
907 else
908 txt.Print("Zd=%.1fd Az=%.1fd", fmod(deg.Y()+270,180)-90, fmod(deg.X()+720, 360));
909
910 TMarker *tip=new TMarker(v0.X(), v0.Y(), kDot);
911 tip->SetMarkerColor(kWhite+type*2);
912 fMapG.Add(tip, new TString(txt));
913
914 return kTRUE;
915}
916
917// --------------------------------------------------------------------------
918//
919// Use "local" draw option to align the display to the local
920// coordinate system instead of the sky coordinate system.
921// dx, dy are arrays storing recuresively all touched points
922// stepx, stepy are the step-size of the current grid.
923//
924void MAstroCatalog::Draw(const TVector2 &v0, const TRotation &rot, TArrayI &dx, TArrayI &dy, Int_t stepx, Int_t stepy, Int_t type)
925{
926 // Calculate the end point
927 const TVector2 v1 = v0 + TVector2(dx[0]*TMath::DegToRad(), dy[0]*TMath::DegToRad());
928
929 // Check whether the point has already been touched.
930 Int_t idx[] = {1, 1, 1, 1};
931
932 Int_t dirs[4][2] = { {0, stepy}, {stepx, 0}, {0, -stepy}, {-stepx, 0} };
933
934 // Check for ambiguities.
935 for (int i=0; i<dx.GetSize(); i++)
936 {
937 for (int j=0; j<4; j++)
938 {
939 const Bool_t rcx0 = (dx[i]+720)%360==(dx[0]+dirs[j][0]+720)%360;
940 const Bool_t rcy0 = (dy[i]+360)%180==(dy[0]+dirs[j][1]+360)%180;
941 if (rcx0&&rcy0)
942 idx[j] = 0;
943 }
944 }
945
946 // Enhance size of array by 1, copy current
947 // position as last entry
948 dx.Set(dx.GetSize()+1);
949 dy.Set(dy.GetSize()+1);
950
951 dx[dx.GetSize()-1] = dx[0];
952 dy[dy.GetSize()-1] = dy[0];
953
954 // Store current positon
955 const Int_t d[2] = { dx[0], dy[0] };
956
957 for (int i=0; i<4; i++)
958 if (idx[i])
959 {
960 // Calculate new position
961 dx[0] = d[0]+dirs[i][0];
962 dy[0] = d[1]+dirs[i][1];
963
964 // Draw corresponding line and iterate through grid
965 if (DrawLine(v1, dirs[i][0], dirs[i][1], rot, type))
966 Draw(v0, rot, dx, dy, stepx, stepy, type);
967
968 dx[0]=d[0];
969 dy[0]=d[1];
970 }
971}
972
973// --------------------------------------------------------------------------
974//
975// Draw a grid recursively around the point v0 (either Ra/Dec or Zd/Az)
976// The points in the grid are converted by a TRotation and CovertToPad
977// to pad coordinates. The type arguemnts is neccessary to create the
978// correct tooltip (Ra/Dec, Zd/Az) at the grid-points.
979// From the pointing position the step-size of teh gris is caluclated.
980//
981void MAstroCatalog::DrawGrid(const TVector3 &v0, const TRotation &rot, Int_t type)
982{
983 TArrayI dx(1);
984 TArrayI dy(1);
985
986 // align to 1deg boundary
987 TVector2 v(v0.Phi()*TMath::RadToDeg(), v0.Theta()*TMath::RadToDeg());
988 v.Set((Float_t)TMath::Nint(v.X()), (Float_t)TMath::Nint(v.Y()));
989
990 // calculate stepsizes based on visible FOV
991 Int_t stepx = 1;
992
993 if (v.Y()<fRadiusFOV || v.Y()>180-fRadiusFOV)
994 stepx=36;
995 else
996 {
997 // This is a rough estimate how many degrees are visible
998 const Float_t m = log(fRadiusFOV/180.)/log(90./(fRadiusFOV+1)+1);
999 const Float_t t = log(180.)-m*log(fRadiusFOV);
1000 const Float_t f = m*log(90-fabs(90-v.Y()))+t;
1001 const Int_t nx = (Int_t)(exp(f)+0.5);
1002 stepx = nx<4 ? 1 : nx/4;
1003 if (stepx>36)
1004 stepx=36;
1005 }
1006
1007 const Int_t ny = (Int_t)(fRadiusFOV+1);
1008 Int_t stepy = ny<4 ? 1 : ny/4;
1009
1010 // align stepsizes to be devisor or 180 and 90
1011 while (180%stepx)
1012 stepx++;
1013 while (90%stepy)
1014 stepy++;
1015
1016 // align to step-size boundary (search for the nearest one)
1017 Int_t dv = 1;
1018 while ((int)(v.X())%stepx)
1019 {
1020 v.Set(v.X()+dv, v.Y());
1021 dv = -TMath::Sign(TMath::Abs(dv)+1, dv);
1022 }
1023
1024 dv = 1;
1025 while ((int)(v.Y())%stepy)
1026 {
1027 v.Set(v.X(), v.Y()+dv);
1028 dv = -TMath::Sign(TMath::Abs(dv)+1, dv);
1029 }
1030
1031 // draw...
1032 v *= TMath::DegToRad();
1033
1034 Draw(v, rot, dx, dy, stepx, stepy, type);
1035}
1036
1037// --------------------------------------------------------------------------
1038//
1039// Get a rotation matrix which aligns the pointing position
1040// to the center of the x,y plain
1041//
1042TRotation MAstroCatalog::AlignCoordinates(const TVector3 &v) const
1043{
1044 TRotation trans;
1045 trans.RotateZ(-v.Phi());
1046 trans.RotateY(-v.Theta());
1047 trans.RotateZ(-TMath::Pi()/2);
1048 return trans;
1049}
1050
1051// --------------------------------------------------------------------------
1052//
1053// Return the rotation matrix which converts either sky or
1054// local coordinates to coordinates which pole is the current
1055// pointing direction.
1056//
1057TRotation MAstroCatalog::GetGrid(Bool_t local)
1058{
1059 const Bool_t enable = fTime && fObservatory;
1060
1061 // If sky coordinate view is requested get rotation matrix and
1062 // draw corresponding sky-grid and if possible local grid
1063 if (!local)
1064 {
1065 const TRotation trans(AlignCoordinates(fRaDec));
1066
1067 DrawGrid(fRaDec, trans, 1);
1068
1069 if (enable)
1070 {
1071 const MAstroSky2Local rot(*fTime, *fObservatory);
1072 DrawGrid(rot*fRaDec, trans*rot.Inverse(), 2);
1073 }
1074
1075 // Return the correct rotation matrix
1076 return trans;
1077 }
1078
1079 // If local coordinate view is requested get rotation matrix and
1080 // draw corresponding sky-grid and if possible local grid
1081 if (local && enable)
1082 {
1083 const MAstroSky2Local rot(*fTime, *fObservatory);
1084
1085 const TRotation trans(AlignCoordinates(rot*fRaDec));
1086
1087 DrawGrid(fRaDec, trans*rot, 1);
1088 DrawGrid(rot*fRaDec, trans, 2);
1089
1090 // Return the correct rotation matrix
1091 return trans*rot;
1092 }
1093
1094 return TRotation();
1095}
1096
1097// --------------------------------------------------------------------------
1098//
1099// Create the title for the pad.
1100//
1101TString MAstroCatalog::GetPadTitle() const
1102{
1103 const Double_t ra = fRaDec.Phi()*TMath::RadToDeg();
1104 const Double_t dec = (TMath::Pi()/2-fRaDec.Theta())*TMath::RadToDeg();
1105
1106 TString txt;
1107 txt += Form("\\alpha=%.2fh ", fmod(ra/15+48, 24));
1108 txt += Form("\\delta=%.1f\\circ ", fmod(dec+270,180)-90);
1109 txt += Form("/ FOV=%.1f\\circ", fRadiusFOV);
1110
1111 if (!fTime || !fObservatory)
1112 return txt;
1113
1114 const MAstroSky2Local rot(*fTime, *fObservatory);
1115 const TVector3 loc = rot*fRaDec;
1116
1117 const Double_t rho = rot.RotationAngle(fRaDec.Phi(), TMath::Pi()/2-fRaDec.Theta());
1118
1119 const Double_t zd = TMath::RadToDeg()*loc.Theta();
1120 const Double_t az = TMath::RadToDeg()*loc.Phi();
1121
1122 txt.Prepend("#splitline{");
1123 txt += Form(" \\theta=%.1f\\circ ", fmod(zd+270,180)-90);
1124 txt += Form("\\phi=%.1f\\circ ", fmod(az+720, 360));
1125 txt += Form(" / \\rho=%.1f\\circ", rho*TMath::RadToDeg());
1126 txt += "}{<";
1127 txt += fTime->GetSqlDateTime();
1128 txt += ">}";
1129 return txt;
1130}
1131
1132// --------------------------------------------------------------------------
1133//
1134// To overlay the catalog make sure, that in any case you are using
1135// the 'same' option.
1136//
1137// If you want to overlay this on top of any picture which is created
1138// by derotation of the camera plain you have to use the 'mirror' option
1139// the compensate the mirroring of the image in the camera plain.
1140//
1141// If you have already compensated this by x=-x and y=-y when creating
1142// the histogram you can simply overlay the catalog.
1143//
1144// To overlay the catalog on a 2D histogram the histogram must have
1145// units of degrees (which are plain, like you directly convert the
1146// camera units by multiplication to degrees)
1147//
1148// To be 100% exact you must use the option 'plain' which assumes a plain
1149// screen. This is not necessary for the MAGIC-camera because the
1150// difference between both is less than 1e-3.
1151//
1152// You should always be aware of the fact, that the shown stars and the
1153// displayed grid is the ideal case, like a reflection on a virtual
1154// perfectly aligned central mirror. In reality the star-positions are
1155// smeared to the edge of the camera the more the distance to the center
1156// is, such that the center of gravity of the light distribution might
1157// be more far away from the center than the display shows.
1158//
1159// If you want the stars to be displayed as circles with a size
1160// showing their magnitude use "*" as an option.
1161//
1162// Use 'white' to display white instead of black stars
1163// Use 'yellow' to display white instead of black stars
1164//
1165//
1166void MAstroCatalog::AddPrimitives(TString o)
1167{
1168 const Bool_t same = o.Contains("same", TString::kIgnoreCase);
1169 const Bool_t local = o.Contains("local", TString::kIgnoreCase);
1170 const Bool_t mirx = o.Contains("mirrorx", TString::kIgnoreCase);
1171 const Bool_t miry = o.Contains("mirrory", TString::kIgnoreCase);
1172 const Bool_t mirror = o.Contains("mirror", TString::kIgnoreCase) && !mirx && !miry;
1173 const Bool_t size = o.Contains("*", TString::kIgnoreCase);
1174 const Bool_t white = o.Contains("white", TString::kIgnoreCase);
1175 const Bool_t yellow = o.Contains("yellow", TString::kIgnoreCase) && !white;
1176
1177 // X is vice versa, because ra is defined anti-clockwise
1178 mirx || mirror ? ResetBit(kMirrorX) : SetBit(kMirrorX);
1179 miry || mirror ? SetBit(kMirrorY) : ResetBit(kMirrorY);
1180
1181 const TRotation rot(GetGrid(local));
1182
1183 TIter Next(&fList);
1184 MVector3 *v=0;
1185 while ((v=(MVector3*)Next()))
1186 {
1187 if (v->Magnitude()>fLimMag)
1188 continue;
1189
1190 TVector2 s(v->Phi(), v->Theta());
1191 if (Convert(rot, s)==kTRUE)
1192 DrawStar(s.X(), s.Y(), *v, yellow?kYellow:(white?kWhite:kBlack), 0, size);
1193 }
1194
1195 if (!same && !TestBit(kDrawingImage) && gPad)
1196 {
1197 TPaveText *pv = new TPaveText(0.01, 0.90, 0.63, 0.99, "brNDC");
1198 pv->AddText(GetPadTitle());
1199 fMapG.Add(pv);
1200 }
1201
1202 TMarker *mk=new TMarker(0, 0, kMultiply);
1203 mk->SetMarkerColor(white||yellow?kWhite:kBlack);
1204 mk->SetMarkerSize(1.5);
1205 fMapG.Add(mk);
1206}
1207
1208// --------------------------------------------------------------------------
1209//
1210// Do nothing if 'same' option given.
1211// Otherwise set pad-range such that x- and y- coordinates have the same
1212// step-size
1213//
1214void MAstroCatalog::SetRangePad(Option_t *o)
1215{
1216 if (TString(o).Contains("same", TString::kIgnoreCase))
1217 return;
1218
1219 const Double_t edge = fRadiusFOV/TMath::Sqrt(2.);
1220 //gPad->Range(-edge, -edge, edge, edge);
1221
1222 const Float_t w = gPad->GetWw();
1223 const Float_t h = gPad->GetWh();
1224
1225 if (w<h)
1226 gPad->Range(-edge, -edge*h/w, edge, edge*h/w);
1227 else
1228 gPad->Range(-edge*w/h, -edge, edge*w/h, edge);
1229}
1230
1231// --------------------------------------------------------------------------
1232//
1233// First delete all gui elements.
1234// Set the correct range of the pad.
1235// Create all gui primitives
1236// If "this" is not in the current pad add it to the current pad.
1237// Reset vit kHasChanged
1238//
1239void MAstroCatalog::DrawPrimitives(Option_t *o)
1240{
1241 fMapG.Delete();
1242
1243 if (!TestBit(kDrawingImage) && gPad)
1244 SetRangePad(o);
1245
1246#ifdef DEBUG
1247 TStopwatch clk;
1248 clk.Start();
1249#endif
1250 AddPrimitives(o);
1251#ifdef DEBUG
1252 clk.Stop();
1253 clk.Print();
1254#endif
1255
1256 // Append to a possible second pad
1257 if (!TestBit(kDrawingImage) && gPad && !gPad->GetListOfPrimitives()->FindObject(this))
1258 AppendPad(o);
1259
1260 ResetBit(kHasChanged);
1261}
1262
1263// --------------------------------------------------------------------------
1264//
1265// Append "this" to current pad
1266// set bit kHasChanged to recreate all gui elements
1267// Connect signal
1268//
1269void MAstroCatalog::Draw(Option_t *o)
1270{
1271 // Append to first pad
1272 AppendPad(o);
1273
1274 // If contents have not previously changed make sure that
1275 // all primitives are recreated.
1276 SetBit(kHasChanged);
1277
1278 // Connect all TCanvas::ProcessedEvent to this->EventInfo
1279 // This means, that after TCanvas has processed an event
1280 // EventInfo of this class is called, see TCanvas::HandleInput
1281 gPad->GetCanvas()->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)",
1282 "MAstroCatalog", this,
1283 "EventInfo(Int_t,Int_t,Int_t,TObject*)");
1284}
1285
1286// --------------------------------------------------------------------------
1287//
1288// This function was connected to all created canvases. It is used
1289// to redirect GetObjectInfo into our own status bar.
1290//
1291// The 'connection' is done in Draw. It seems that 'connected'
1292// functions must be public.
1293//
1294void MAstroCatalog::EventInfo(Int_t event, Int_t px, Int_t py, TObject *selected)
1295{
1296 TCanvas *c = (TCanvas*)gTQSender;
1297
1298 gPad = c ? c->GetSelectedPad() : NULL;
1299 if (!gPad)
1300 return;
1301
1302
1303 // Try to find a corresponding object with kCannotPick set and
1304 // an available TString (for a tool tip)
1305 TString str;
1306 if (!selected || selected==this)
1307 selected = fMapG.PickObject(px, py, str);
1308
1309 if (!selected)
1310 return;
1311
1312 // Handle some gui events
1313 switch (event)
1314 {
1315 case kMouseMotion:
1316 if (fToolTip && !fToolTip->IsMapped() && !str.IsNull())
1317 ShowToolTip(px, py, str);
1318 break;
1319
1320 case kMouseLeave:
1321 if (fToolTip && fToolTip->IsMapped())
1322 fToolTip->Hide();
1323 break;
1324
1325 case kKeyPress:
1326 ExecuteEvent(kKeyPress, px, py);
1327 break;
1328 }
1329}
1330
1331// --------------------------------------------------------------------------
1332//
1333// Handle keyboard events.
1334//
1335void MAstroCatalog::ExecuteEventKbd(Int_t keycode, Int_t keysym)
1336{
1337 Double_t dra =0;
1338 Double_t ddec=0;
1339
1340 switch (keysym)
1341 {
1342 case kKey_Left:
1343 dra = -TMath::DegToRad();
1344 break;
1345 case kKey_Right:
1346 dra = +TMath::DegToRad();
1347 break;
1348 case kKey_Up:
1349 ddec = +TMath::DegToRad();
1350 break;
1351 case kKey_Down:
1352 ddec = -TMath::DegToRad();
1353 break;
1354 case kKey_Plus:
1355 SetRadiusFOV(fRadiusFOV+1);
1356 break;
1357 case kKey_Minus:
1358 SetRadiusFOV(fRadiusFOV-1);
1359 break;
1360
1361 default:
1362 return;
1363 }
1364
1365 const Double_t r = fRaDec.Phi();
1366 const Double_t d = TMath::Pi()/2-fRaDec.Theta();
1367
1368 SetRaDec(r+dra, d+ddec);
1369
1370 gPad->Update();
1371}
1372
1373// ------------------------------------------------------------------------
1374//
1375// Execute a gui event on the camera
1376//
1377void MAstroCatalog::ExecuteEvent(Int_t event, Int_t mp1, Int_t mp2)
1378{
1379 if (!TestBit(kGuiActive))
1380 return;
1381
1382 if (event==kKeyPress)
1383 ExecuteEventKbd(mp1, mp2);
1384}
1385
1386// --------------------------------------------------------------------------
1387//
1388// Displays a tooltip
1389//
1390void MAstroCatalog::ShowToolTip(Int_t px, Int_t py, const char *txt)
1391{
1392 if (TestBit(kNoToolTips))
1393 return;
1394
1395 Int_t x=0;
1396 Int_t y=0;
1397
1398 const Window_t id1 = gVirtualX->GetWindowID(gPad->GetCanvasID());
1399 const Window_t id2 = fToolTip->GetParent()->GetId();
1400
1401 Window_t id3;
1402 gVirtualX->TranslateCoordinates(id1, id2, px, py, x, y, id3);
1403
1404 // Show tool tip
1405 fToolTip->SetText(txt);
1406 fToolTip->Show(x+4, y+4);
1407}
1408
1409// --------------------------------------------------------------------------
1410//
1411// Calculate distance to primitive by checking all gui elements
1412//
1413Int_t MAstroCatalog::DistancetoPrimitive(Int_t px, Int_t py)
1414{
1415 return fMapG.DistancetoPrimitive(px, py);
1416}
1417
1418// ------------------------------------------------------------------------
1419//
1420// Returns string containing info about the object at position (px,py).
1421// Returned string will be re-used (lock in MT environment).
1422//
1423char *MAstroCatalog::GetObjectInfo(Int_t px, Int_t py) const
1424{
1425 return fMapG.GetObjectInfo(px, py);
1426}
Note: See TracBrowser for help on using the repository browser.