source: branches/Mars_MC/mastro/MAstroCamera.cc@ 17782

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