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

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