source: trunk/MagicSoft/Mars/mastro/MAstroCamera.cc@ 4983

Last change on this file since 4983 was 4715, checked in by rwagner, 20 years ago
*** empty log message ***
File size: 16.5 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 expressed
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Robert Wagner, 7/2004 <mailto:rwagner@mppmu.mpg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MAstroCamera
29// ============
30//
31// A tools displaying stars from a catalog in the camera display.
32// PRELIMINARY!!
33//
34//
35// Usage:
36// ------
37// For a usage example see macros/starfield.C
38//
39// To be able to reflect the star-light you need the geometry of the
40// mirror and the distance of the plain screen.
41//
42// You can get the mirror geometry from a MC file and the distance of
43// the screen from MGeomCam.
44//
45//
46// Algorithm:
47// ----------
48// The calculation of the position of the reflection in the camera is
49// done by:
50// - Rotation of the star-field such that the camera is looking into
51// the pointing direction
52// - Calculation of the reflected star-light vector by calling
53// MGeomMirror::GetReflection (which returns the point at which
54// the vector hits the camera plain)
55// - Depending on the draw-option you get each reflected point, the
56// reflection on a virtual ideal mirror or the reflection on each
57// individual mirror
58//
59//
60// GUI:
61// ----
62// * You can use the the cursor keys to change the pointing position
63// and plus/minus to change the time by a quarter of an hour.
64//
65// ToDo:
66// -----
67// * possibility to overwrite distance from mirror to screen
68// * read the mirror geometry directly from the MC input file
69//
70/////////////////////////////////////////////////////////////////////////////
71#include "MAstroCamera.h"
72
73#include <errno.h> // strerror
74#include <fstream> // ifstream
75
76#include <KeySymbols.h> // kKey_*
77
78#include <TH2.h> // TH2D
79#include <TMarker.h> // TMarker
80#include <TVirtualPad.h> // gPad
81
82#include "MLog.h"
83#include "MLogManip.h"
84
85#include "MGeomCam.h"
86#include "MGeomMirror.h"
87
88#include "MTime.h"
89#include "MAstroSky2Local.h"
90#include "MObservatory.h"
91
92#include "../mhist/MHCamera.h" // FIXME: This dependancy is very bad!
93 // HOW TO GET RID OF IT? Move MHCamera to mgeom?
94
95#include "MStarPos.h"
96
97ClassImp(MAstroCamera);
98
99using namespace std;
100
101// --------------------------------------------------------------------------
102//
103// Create a virtual MGeomMirror which is in the center of the coordinate
104// system and has a normal vector in z-direction.
105//
106MAstroCamera::MAstroCamera() : fGeom(0), fMirrors(0)
107{
108 fMirror0 = new MGeomMirror;
109 fMirror0->SetMirrorContent(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1);
110}
111
112// --------------------------------------------------------------------------
113//
114// Delete fGeom, fMirrors and the virtual 0-Mirror fMirror0
115//
116MAstroCamera::~MAstroCamera()
117{
118 if (fGeom)
119 delete fGeom;
120 if (fMirrors)
121 delete fMirrors;
122
123 delete fMirror0;
124
125}
126
127// --------------------------------------------------------------------------
128//
129// Set a list of mirrors. The Mirrors must be of type MGeomMirror and
130// stored in a TClonesArray
131//
132void MAstroCamera::SetMirrors(TClonesArray &arr)
133{
134 if (arr.GetClass()!=MGeomMirror::Class())
135 {
136 cout << "ERROR - TClonesArray doesn't contain objects of type MGeomMirror... ignored." << endl;
137 return;
138 }
139
140 const Int_t n = arr.GetSize();
141
142 if (!fMirrors)
143 fMirrors = new TClonesArray(MGeomMirror::Class(), n);
144
145 fMirrors->ExpandCreate(n);
146
147 for (int i=0; i<n; i++)
148 memcpy((*fMirrors)[i], arr[i], sizeof(MGeomMirror));
149
150}
151
152
153// --------------------------------------------------------------------------
154//
155// Read the mirror geometry from a MC .def file. The following
156// structure is expected:
157//
158// #* TYPE=1 (MAGIC)
159// #* i f sx sy x y z thetan phin
160// #*
161// #* i : number of the mirror
162// #* f : focal distance of that mirror
163// #* sx : curvilinear coordinate of mirror's center in X[cm]
164// #* sy : curvilinear coordinate of mirror's center in X[cm]
165// #* x : x coordinate of the center of the mirror [cm]
166// #* y : y coordinate of the center of the mirror [cm]
167// #* z : z coordinate of the center of the mirror [cm]
168// #* thetan : polar theta angle of the direction where the mirror points to
169// #* phin : polar phi angle of the direction where the mirror points to
170// #* xn : xn coordinate of the normal vector in the center (normalized)
171// #* yn : yn coordinate of the normal vector in the center (normalized)
172// #* zn : zn coordinate of the normal vector in the center (normalized)
173// #
174// define_mirrors
175// 1 1700.9200 25.0002 75.0061 25.0000 75.0000 0.9207 0.02328894 1.24904577 -0.00736394 -0.02209183 0.99972882
176// 2 ...
177//
178void MAstroCamera::SetMirrors(const char *fname)
179{
180 ifstream fin(fname);
181 if (!fin)
182 {
183 gLog << err << "Cannot open file " << fname << ": ";
184 gLog << strerror(errno) << endl;
185 return;
186 }
187
188 TString line;
189 while (1)
190 {
191 line.ReadLine(fin);
192 if (!fin)
193 return;
194
195 line = line.Strip(TString::kBoth);
196
197 if (line.BeginsWith("n_mirrors"))
198 {
199 Int_t n;
200 sscanf(line.Data(), "%*s %d", &n);
201
202 if (!fMirrors)
203 fMirrors = new TClonesArray(MGeomMirror::Class(), n);
204
205 fMirrors->ExpandCreate(n);
206 continue;
207 }
208
209
210 Int_t id;
211 Float_t f, sx, sy, x, y, z, thetan, phin, xn, yn, zn;
212
213 const Int_t n = sscanf(line.Data(), "%d %f %f %f %f %f %f %f %f %f %f %f",
214 &id, &f, &sx, &sy, &x, &y, &z, &thetan,
215 &phin, &xn, &yn, &zn);
216 if (n!=12)
217 continue;
218
219 new ((*fMirrors)[id-1]) MGeomMirror;
220 ((MGeomMirror*)(*fMirrors)[id-1])->SetMirrorContent(id, f, sx, sy, x, y, z, thetan, phin, xn, yn, zn);
221 }
222}
223
224// --------------------------------------------------------------------------
225//
226// Set the camera geometry. The MGeomCam object is cloned.
227//
228void MAstroCamera::SetGeom(const MGeomCam &cam)
229{
230 if (fGeom)
231 delete fGeom;
232
233 fGeom = (MGeomCam*)cam.Clone();
234}
235
236// --------------------------------------------------------------------------
237//
238// Convert To Pad coordinates (see MAstroCatalog)
239//
240Int_t MAstroCamera::ConvertToPad(const TVector3 &w, TVector2 &v) const
241{
242 /*
243 --- Use this to plot the 'mean grid' instead of the grid of a
244 theoretical central mirror ---
245
246 TVector3 spot;
247 const Int_t num = fConfig->GetNumMirror();
248 for (int i=0; i<num; i++)
249 spot += fConfig->GetMirror(i).GetReflection(w, fGeom->GetCameraDist())*1000;
250 spot *= 1./num;
251 */
252
253 const TVector3 spot = fMirror0->GetReflection(w, fGeom->GetCameraDist())*1000;
254 v.Set(spot(0), spot(1));
255
256 const Float_t max = fGeom->GetMaxRadius()*0.70;
257 return v.X()>-max && v.Y()>-max && v.X()<max && v.Y()<max;
258}
259
260// --------------------------------------------------------------------------
261//
262// Find an object with a given name in the list of primitives of this pad.
263//
264TObject *FindObjectInPad(const char *name, TVirtualPad *pad)
265{
266 if (!pad)
267 pad = gPad;
268
269 if (!pad)
270 return NULL;
271
272 TObject *o;
273
274 TIter Next(pad->GetListOfPrimitives());
275 while ((o=Next()))
276 {
277 if (o->InheritsFrom(gROOT->GetClass(name)))
278 return o;
279
280 if (o->InheritsFrom("TPad"))
281 if ((o = FindObjectInPad(name, (TVirtualPad*)o)))
282 return o;
283 }
284 return NULL;
285}
286
287// --------------------------------------------------------------------------
288//
289// Options:
290//
291// '*' Draw the mean of the reflections on all mirrors
292// '.' Draw a dot for the reflection on each individual mirror
293// 'h' To create a TH2D of the star-light which is displayed
294// 'c' Use the underlaying MHCamera as histogram
295// '0' Draw the reflection on a virtual perfect mirror
296// '=' Draw '0' or '*' propotional to the star magnitude
297//
298// If the Pad contains an object MHCamera of type MHCamera it is used.
299// Otherwise a new object is created.
300//
301void MAstroCamera::AddPrimitives(TString o)
302{
303 if (!fTime || !fObservatory || !fMirrors)
304 {
305 cout << "Missing data..." << endl;
306 return;
307 }
308
309 if (o.IsNull())
310 o = "*.";
311
312 const Bool_t hasnull = o.Contains("0", TString::kIgnoreCase);
313 const Bool_t hashist = o.Contains("h", TString::kIgnoreCase);
314 const Bool_t hasmean = o.Contains("*", TString::kIgnoreCase);
315 const Bool_t hasdot = o.Contains(".", TString::kIgnoreCase);
316 const Bool_t usecam = o.Contains("c", TString::kIgnoreCase);
317 const Bool_t resize = o.Contains("=", TString::kIgnoreCase);
318
319 // Get camera
320 MHCamera *camera=(MHCamera*)FindObjectInPad("MHCamera", gPad);
321 if (camera)
322 {
323 if (!camera->GetGeometry() || camera->GetGeometry()->IsA()!=fGeom->IsA())
324 camera->SetGeometry(*fGeom);
325 }
326 else
327 {
328 camera = new MHCamera(*fGeom);
329 camera->SetName("MHCamera");
330 camera->SetStats(0);
331 camera->SetInvDeepBlueSeaPalette();
332 camera->SetBit(kCanDelete);
333 camera->Draw();
334 }
335
336 camera->SetTitle(GetPadTitle());
337
338 gPad->cd(1);
339
340 if (!usecam)
341 {
342 if (camera->GetEntries()==0)
343 camera->SetBit(MHCamera::kNoLegend);
344 }
345 else
346 {
347 camera->Reset();
348 camera->SetYTitle("arb.cur");
349 }
350
351 TH2 *h=0;
352 if (hashist)
353 {
354 TH2F hist("","", 90, -650, 650, 90, -650, 650);
355 hist.SetMinimum(0);
356 h = (TH2*)hist.DrawCopy("samecont1");
357 }
358
359 const TRotation rot(GetGrid(kTRUE));
360
361 MVector3 *radec;
362 TIter Next(&fList);
363
364 while ((radec=(MVector3*)Next()))
365 {
366 const Double_t mag = radec->Magnitude();
367 if (mag>GetLimMag())
368 continue;
369
370 TVector3 star(*radec);
371
372 // Rotate Star into telescope system
373 star *= rot;
374
375 TVector3 mean;
376
377 Int_t num = 0;
378
379 MGeomMirror *mirror = 0;
380 TIter NextM(fMirrors);
381 while ((mirror=(MGeomMirror*)NextM()))
382 {
383 const TVector3 spot = mirror->GetReflection(star, fGeom->GetCameraDist())*1000;
384
385 // calculate mean of all 'stars' hitting the camera plane
386 // by taking the sum of the intersection points between
387 // the light vector and the camera plane
388 mean += spot;
389
390 if (hasdot)
391 {
392 TMarker *m=new TMarker(spot(0), spot(1), 1);
393 TVector3 pos=mirror->GetMirrorCenter();
394 Float_t r=pos.Mag();
395
396 m->SetMarkerColor((Int_t)(6*r)+51);
397 m->SetMarkerStyle(kDot);
398 fMapG.Add(m);
399 }
400 if (h)
401 h->Fill(spot(0), spot(1), pow(10, -mag/2.5));
402
403 if (usecam)
404 camera->Fill(spot(0), spot(1), pow(10, -mag/2.5));
405
406 num++;
407 }
408
409 // transform meters into millimeters (camera display works with mm)
410 mean *= 1./num;
411 DrawStar(mean(0), mean(1), *radec, hasmean?kBlack:-1, Form("x=%.1fmm y=%.1fmm", mean(0), mean(1)), resize);
412 if (hasnull)
413 {
414 TVector3 star(*radec);
415 star *= rot;
416 const TVector3 spot = fMirror0->GetReflection(star, fGeom->GetCameraDist())*1000;
417 DrawStar(spot(0), spot(1), *radec, hasmean?kBlack:-1, Form("x=%.1fmm y=%.1fmm", mean(0), mean(1)), resize);
418 // This can be used to get the abberation...
419 //cout << TMath::Hypot(spot(0), spot(1)) << " " << TMath::Hypot(mean(0)-spot(0), mean(1)-spot(1)) << endl;
420 }
421 }
422}
423
424
425// --------------------------------------------------------------------------
426//
427// Fills a TList with all stars found under current presets
428// (Coordinates, Time, FOV). The list is emptied before the filling is done.
429// Currently, the mean spot when averaging over all mirrors is used.
430//
431void MAstroCamera::FillStarList(TList* list)
432{
433 if (!fTime || !fObservatory || !fMirrors) {
434 cout << "Missing data..." << endl;
435 return;
436 }
437
438 if (!list) {
439 cout << "Missing storage list for stars found..." << endl;
440 return;
441 }
442
443 list->Delete(); // dump old list... (otherwise the same stars would pile up)
444
445 const TRotation rot(GetGrid(kTRUE));
446 MVector3 *radec;
447 TIter Next(&fList);
448
449 while ((radec=(MVector3*)Next())) {
450 const Double_t mag = radec->Magnitude();
451 if (mag>GetLimMag())
452 continue;
453 TVector3 star(*radec);
454 // Rotate Star into telescope system
455 star *= rot;
456 TVector3 mean;
457 Int_t num = 0;
458 MGeomMirror *mirror = 0;
459 TIter NextM(fMirrors);
460 while ((mirror=(MGeomMirror*)NextM())) {
461 const TVector3 spot = mirror->GetReflection(star, fGeom->GetCameraDist())*1000;
462 mean += spot;
463 num++;
464 }
465 mean *= 1./num;
466 MStarPos *starpos = new MStarPos;
467 starpos->SetExpValues(mag,mean(0),mean(1));
468
469 // Do _not_ remove this, name needs to be
470 // transferred to the starpos.
471 TString name = radec->GetName();
472 if (name.Length()==0) {
473 starpos->SetName("unknown");
474 } else {
475 starpos->SetName(radec->GetName());
476 }
477
478 const TVector3 spot = fMirror0->GetReflection(star, fGeom->GetCameraDist())*1000;
479 starpos->SetIdealValues(mag,spot(0),spot(1));
480 list->Add(starpos);
481 }
482}
483
484// ------------------------------------------------------------------------
485//
486// Uses fRaDec as a reference point.
487//
488// Return dZd and dAz corresponding to the distance from x,y to fRaDec
489//
490// Before calling this function you should correct for abberation. In
491// case of MAGIC you can do this by:
492// x /= 1.0713;
493// y /= 1.0713;
494//
495// x [mm]: x coordinate in the camera plane (assuming a perfect mirror)
496// y [mm]: y coordinate in the camera plane (assuming a perfect mirror)
497//
498// We assume (0, 0) to be the center of the FOV
499//
500// dzd [deg]: Delta Zd
501// daz [deg]: Delta Az
502//
503void MAstroCamera::GetDiffZdAz(Double_t x, Double_t y, Double_t &dzd, Double_t &daz)
504{
505 // Reflect the corrected pixel on a perfect mirror
506 TVector3 v(x, y, fGeom->GetCameraDist()*1000);
507 TVector3 spot = fMirror0->GetReflection(v);
508
509 // Derotate into the local coordinate system
510 const MAstroSky2Local rot(*fTime, *fObservatory);
511 const TRotation align(AlignCoordinates(rot*fRaDec).Inverse());
512 spot *= align;
513
514 cout << "Zd="<<spot.Theta()*TMath::RadToDeg() << " ";
515 cout << "Az="<<spot.Phi() *TMath::RadToDeg()+360 << endl;
516
517 // Derotatet the center of the camera
518 TVector3 c(0, 0, 1);
519 c *= align;
520
521 dzd = (spot.Theta()-c.Theta())*TMath::RadToDeg();
522 daz = (spot.Phi() -c.Phi()) *TMath::RadToDeg();
523
524 if (daz> 180) daz -= 360;
525 if (daz<-180) daz += 360;
526}
527
528// ------------------------------------------------------------------------
529//
530// Execute a gui event on the camera
531//
532void MAstroCamera::ExecuteEvent(Int_t event, Int_t mp1, Int_t mp2)
533{
534 // if (mp1>0 && mp2>0)
535 // {
536 // // Calculate World coordinates from pixel
537 // Double_t x = gPad->AbsPixeltoX(mp1);
538 // Double_t y = gPad->AbsPixeltoY(mp2);
539 //
540 // // Correct for abberation
541 // x /= 1.0713;
542 // y /= 1.0713;
543 //
544 // Double_t dzd, daz;
545 // GetDiffZdAz(x, y, dzd, daz);
546 //
547 // cout << "dZd="<< dzd << " " << "dAz="<< daz << endl;
548 // }
549 //
550 // For MAGIC the distance of the mean of the light distribution
551 // to the Mirror0 reflection of the star (Abberation) can be
552 // expressed as: dr = 0.0713*r = r/14.03
553 // +-0.0002
554
555 if (event==kKeyPress && fTime)
556 switch (mp2)
557 {
558 case kKey_Plus:
559 fTime->SetMjd(fTime->GetMjd()+0.25/24);
560 Update(kTRUE);
561 return;
562
563 case kKey_Minus:
564 fTime->SetMjd(fTime->GetMjd()-0.25/24);
565 Update(kTRUE);
566 return;
567 }
568
569 MAstroCatalog::ExecuteEvent(event, mp1, mp2);
570}
Note: See TracBrowser for help on using the repository browser.