source: trunk/MagicSoft/Mars/mgui/MHexagon.cc@ 8199

Last change on this file since 8199 was 8178, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 14.1 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MHexagon.cc,v 1.30 2006-10-30 12:46:13 tbretz Exp $
3! --------------------------------------------------------------------------
4!
5! *
6! * This file is part of MARS, the MAGIC Analysis and Reconstruction
7! * Software. It is distributed to you in the hope that it can be a useful
8! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
9! * It is distributed WITHOUT ANY WARRANTY.
10! *
11! * Permission to use, copy, modify and distribute this software and its
12! * documentation for any purpose is hereby granted without fee,
13! * provided that the above copyright notice appear in all copies and
14! * that both that copyright notice and this permission notice appear
15! * in supporting documentation. It is provided "as is" without express
16! * or implied warranty.
17! *
18!
19!
20! Author(s): Thomas Bretz 12/2000 <mailto:tbretz@astro.uni-wuerzburg.de>
21! Author(s): Harald Kornmayer 1/2001
22!
23! Copyright: MAGIC Software Development, 2000-2006
24!
25!
26\* ======================================================================== */
27
28//////////////////////////////////////////////////////////////////////////////
29//
30// MHexagon
31//
32//////////////////////////////////////////////////////////////////////////////
33
34#include "MHexagon.h"
35
36#include <fstream>
37#include <iostream>
38
39#include <TVector2.h> // TVector2
40#include <TVirtualPad.h> // gPad
41#include <TOrdCollection.h> // TOrdCollection
42
43#include "MMath.h"
44#include "MGeomPix.h" // GetX
45
46ClassImp(MHexagon);
47
48using namespace std;
49
50const Double_t MHexagon::fgCos60 = 0.5; // TMath::Cos(60/TMath::RadToDeg());
51const Double_t MHexagon::fgSin60 = sqrt(3.)/2; // TMath::Sin(60/TMath::RadToDeg());
52
53const Double_t MHexagon::fgDx[6] = { fgCos60, 0., -fgCos60, -fgCos60, 0., fgCos60 };
54const Double_t MHexagon::fgDy[6] = { fgSin60/3, fgSin60*2/3, fgSin60/3, -fgSin60/3, -fgSin60*2/3, -fgSin60/3 };
55
56// ------------------------------------------------------------------------
57//
58// default constructor for MHexagon
59//
60MHexagon::MHexagon()
61{
62}
63
64// ------------------------------------------------------------------------
65//
66// normal constructor for MHexagon
67//
68MHexagon::MHexagon(Float_t x, Float_t y, Float_t d)
69: TAttLine(1, 1, 1), TAttFill(0, 1001), fX(x), fY(y), fD(d)
70{
71}
72
73// ------------------------------------------------------------------------
74//
75// normal constructor for MHexagon
76//
77MHexagon::MHexagon(const MGeomPix &pix)
78: TAttLine(1, 1, 1), TAttFill(0, 1001)
79{
80 fX = pix.GetX();
81 fY = pix.GetY();
82 fD = pix.GetD();
83}
84
85// ------------------------------------------------------------------------
86//
87// copy constructor for MHexagon
88//
89MHexagon::MHexagon(const MHexagon &hexagon) : TObject(hexagon), TAttLine(hexagon), TAttFill(hexagon)
90{
91 fX = hexagon.fX;
92 fY = hexagon.fY;
93 fD = hexagon.fD;
94}
95
96// ------------------------------------------------------------------------
97//
98// copy this hexagon to hexagon
99//
100void MHexagon::Copy(TObject &obj)
101#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
102const
103#endif
104{
105 MHexagon &hex = (MHexagon&) obj;
106
107 TObject::Copy(obj);
108 TAttLine::Copy(hex);
109 TAttFill::Copy(hex);
110
111 hex.fX = fX;
112 hex.fY = fY;
113 hex.fD = fD;
114}
115
116// ------------------------------------------------------------------------
117//
118// compute the distance of a point (px,py) to the Hexagon
119// this functions needed for graphical primitives, that
120// means without this function you are not able to interact
121// with the graphical primitive with the mouse!!!
122//
123// All calcutations are running in pixel coordinates
124//
125Int_t MHexagon::DistancetoPrimitive(Int_t px, Int_t py, Float_t conv)
126{
127 //
128 // compute the distance of the Point to the center of the Hexagon
129 //
130 const Int_t pxhex = gPad->XtoAbsPixel(fX*conv);
131 const Int_t pyhex = gPad->YtoAbsPixel(fY*conv);
132
133 const Double_t x = TMath::Abs(px-pxhex);
134 const Double_t y = TMath::Abs(py-pyhex);
135
136 const Double_t disthex = TMath::Sqrt(x*x + y*y);
137
138 if (disthex==0)
139 return 0;
140
141 const Double_t cosa = x / disthex;
142 const Double_t sina = y / disthex;
143
144 //
145 // comput the distance of hexagon center to pixel border
146 //
147 const Double_t dx = fD/2 * cosa;
148 const Double_t dy = fD/2 * sina;
149
150 const Int_t pxborder = gPad->XtoAbsPixel((fX + dx)*conv);
151 const Int_t pyborder = gPad->YtoAbsPixel((fY + dy)*conv);
152
153 const Double_t bx = pxhex-pxborder;
154 const Double_t by = pyhex-pyborder;
155
156 const Double_t distborder = TMath::Sqrt(bx*bx + by*by);
157
158 //
159 // compute the distance from the border of Pixel
160 // here in the first implementation is just circle inside
161 //
162 return distborder < disthex ? (int)((disthex-distborder)/conv+1) : 0;
163}
164
165// ------------------------------------------------------------------------
166//
167// compute the distance of a point (px,py) to the Hexagon center in world
168// coordinates. Return -1 if inside.
169//
170Float_t MHexagon::DistanceToPrimitive(Float_t px, Float_t py) const
171{
172 //
173 // compute the distance of the Point to the center of the Hexagon
174 //
175 const Double_t dx = px-fX;
176 const Double_t dy = py-fY;
177 const Double_t disthex = TMath::Sqrt(dx*dx + dy*dy);
178
179 //
180 // Now check if point is outside of hexagon; just check x coordinate
181 // in three coordinate systems: the default one, in which two sides of
182 // the hexagon are paralel to the y axis (see camera displays) and two
183 // more, rotated with respect to that one by +- 60 degrees.
184 //
185
186 if (TMath::Abs(dx) > fD/2.)
187 return disthex;
188
189 const Double_t dx2 = dx*fgCos60 + dy*fgSin60;
190
191 if (TMath::Abs(dx2) > fD/2.)
192 return disthex;
193
194 const Double_t dx3 = dx*fgCos60 - dy*fgSin60;
195
196 if (TMath::Abs(dx3) > fD/2.)
197 return disthex;
198
199 return -1;
200}
201
202/*
203Float_t MHexagon::DistanceToPrimitive(Float_t px, Float_t py)
204{
205 //
206 // compute the distance of the Point to the center of the Hexagon
207 //
208 const Double_t dx = px-fX;
209 const Double_t dy = py-fY;
210
211 const Double_t disthex = TMath::Sqrt(dx*dx + dy*dy);
212
213 //
214 // compute the distance from the border of Pixel
215 // here in the first implementation is just circle outside
216 //
217 return fD*0.5772 < disthex ? disthex-fD*0.5772 : -1;
218 //
219 // FIXME: this approximate method results in some photons being
220 // assigned to two or even three different pixels.
221 //
222}
223*/
224
225// ------------------------------------------------------------------------
226//
227// Draw this ellipse with new coordinate
228//
229void MHexagon::DrawHexagon(Float_t x, Float_t y, Float_t d)
230{
231 MHexagon *newhexagon = new MHexagon(x, y, d);
232
233 TAttLine::Copy(*newhexagon);
234 TAttFill::Copy(*newhexagon);
235
236 newhexagon->SetBit(kCanDelete);
237 newhexagon->AppendPad();
238}
239
240/*
241// ------------------------------------------------------------------------
242//
243// This is the first test of implementing a clickable interface
244// for one pixel
245//
246void MHexagon::ExecuteEvent(Int_t event, Int_t px, Int_t py)
247{
248 switch (event)
249 {
250 case kButton1Down:
251 cout << endl << "kButton1Down" << endl;
252 SetFillColor(2);
253 gPad->Modified();
254 break;
255
256 case kButton1Double:
257 SetFillColor(0);
258 gPad->Modified();
259 break;
260 // case kMouseEnter:
261 // printf ("\n Mouse inside object \n" ) ;
262 // break;
263 }
264}
265*/
266
267// ------------------------------------------------------------------------
268//
269// list this hexagon with its attributes
270//
271void MHexagon::ls(const Option_t *) const
272{
273 TROOT::IndentLevel();
274 Print();
275}
276
277// ------------------------------------------------------------------------
278//
279// paint this hexagon with its attribute
280//
281void MHexagon::Paint(Option_t *)
282{
283 PaintHexagon(fX, fY, fD);
284}
285
286// ------------------------------------------------------------------------
287//
288// draw this hexagon with the coordinates
289//
290void MHexagon::PaintHexagon(Float_t inX, Float_t inY, Float_t inD)
291{
292 //
293 // calculate the positions of the pixel corners
294 //
295 Double_t x[7], y[7];
296 for (Int_t i=0; i<7; i++)
297 {
298 x[i] = inX + fgDx[i%6]*inD;
299 y[i] = inY + fgDy[i%6]*inD;
300 }
301
302 TAttLine::Modify(); // Change line attributes only if neccessary
303 TAttFill::Modify(); // Change fill attributes only if neccessary
304
305 //
306 // paint the pixel
307 //
308 if (GetFillColor())
309 gPad->PaintFillArea(6, x, y);
310
311 if (GetLineStyle())
312 gPad->PaintPolyLine(7, x, y);
313}
314
315// ------------------------------------------------------------------------
316//
317// print/dump this hexagon with its attributes
318//
319void MHexagon::Print(Option_t *) const
320{
321 cout << "MHexagon: ";
322 cout << "x=" << fX << "mm y=" << fY << "mm r=" << fD << "mm" << endl;
323
324 cout << " Line:";
325 cout << " Color=" << GetLineColor() << ",";
326 cout << " Style=" << GetLineStyle() << ",";
327 cout << " Width=" << GetLineWidth() << endl;
328 cout << " Fill:";
329 cout << " Color=" << GetFillColor() << ",";
330 cout << " Style=" << GetFillStyle() << endl;
331}
332
333// ------------------------------------------------------------------------
334//
335// Save primitive as a C++ statement(s) on output stream out
336//
337void MHexagon::SavePrimitive(ostream &out, Option_t *)
338{
339#if ROOT_VERSION_CODE >= ROOT_VERSION(5,12,00)
340 if (gROOT->ClassSaved(Class()))
341 out << " ";
342 else
343 out << " MHexagon *";
344
345 out << "hexagon = new MHexagon(" << fX << "," << fY << "," << fD << ");" << endl;
346
347 SaveFillAttributes(out, "hexagon");
348 SaveLineAttributes(out, "hexagon");
349
350 out << " hexagon->Draw();" << endl;
351#endif
352}
353
354void MHexagon::SavePrimitive(ofstream &out, Option_t *)
355{
356#if ROOT_VERSION_CODE < ROOT_VERSION(5,12,00)
357 if (gROOT->ClassSaved(Class()))
358 out << " ";
359 else
360 out << " MHexagon *";
361
362 out << "hexagon = new MHexagon(" << fX << "," << fY << "," << fD << ");" << endl;
363
364 SaveFillAttributes(out, "hexagon");
365 SaveLineAttributes(out, "hexagon");
366
367 out << " hexagon->Draw();" << endl;
368#else
369 MHexagon::SavePrimitive(static_cast<ostream&>(out), "");
370#endif
371}
372
373// ------------------------------------------------------------------------
374//
375// Small helper class to allow fast sorting of TVector2 by angle
376//
377class HVector2 : public TVector2
378{
379 Double_t fAngle;
380public:
381 HVector2() : TVector2() { }
382 HVector2(const TVector2 &v) : TVector2(v) { }
383 HVector2(Double_t x, Double_t y) : TVector2 (x, y) { }
384 void CalcAngle() { fAngle = Phi(); }
385 Bool_t IsSortable() const { return kTRUE; }
386 Int_t Compare(const TObject *obj) const { return ((HVector2*)obj)->fAngle>fAngle ? 1 : -1; }
387 Double_t GetAngle() const { return fAngle; }
388};
389
390// ------------------------------------------------------------------------
391//
392// Calculate the edge points of the intersection area of two hexagons.
393// The points are added as TVector2 to the TOrdCollection.
394// The user is responsible of delete the objects.
395//
396void MHexagon::GetIntersectionBorder(TOrdCollection &col, const MHexagon &hex) const
397{
398 Bool_t inside0[6], inside1[6];
399
400 HVector2 v0[6];
401 HVector2 v1[6];
402
403 Int_t cnt0=0;
404 Int_t cnt1=0;
405
406 // Calculate teh edges of each hexagon and whether this edge
407 // is inside the other hexgon or not
408 for (int i=0; i<6; i++)
409 {
410 const Double_t x0 = fX+fgDx[i]*fD;
411 const Double_t y0 = fY+fgDy[i]*fD;
412
413 const Double_t x1 = hex.fX+fgDx[i]*hex.fD;
414 const Double_t y1 = hex.fY+fgDy[i]*hex.fD;
415
416 v0[i].Set(x0, y0);
417 v1[i].Set(x1, y1);
418
419 inside0[i] = hex.DistanceToPrimitive(x0, y0) < 0;
420 inside1[i] = DistanceToPrimitive(x1, y1) < 0;
421
422 if (inside0[i])
423 {
424 col.Add(new HVector2(v0[i]));
425 cnt0++;
426 }
427 if (inside1[i])
428 {
429 col.Add(new HVector2(v1[i]));
430 cnt1++;
431 }
432 }
433
434 if (cnt0==0 || cnt1==0)
435 return;
436
437 // No calculate which vorder lines intersect
438 Bool_t iscross0[6], iscross1[6];
439 for (int i=0; i<6; i++)
440 {
441 iscross0[i] = (inside0[i] && !inside0[(i+1)%6]) || (!inside0[i] && inside0[(i+1)%6]);
442 iscross1[i] = (inside1[i] && !inside1[(i+1)%6]) || (!inside1[i] && inside1[(i+1)%6]);
443 }
444
445 // Calculate the border points of our intersection area
446 for (int i=0; i<6; i++)
447 {
448 // Skip non intersecting lines
449 if (!iscross0[i])
450 continue;
451
452 for (int j=0; j<6; j++)
453 {
454 // Skip non intersecting lines
455 if (!iscross1[j])
456 continue;
457
458 const TVector2 p = MMath::GetIntersectionPoint(v0[i], v0[(i+1)%6], v1[j], v1[(j+1)%6]);
459 if (hex.DistanceToPrimitive(p.X(), p.Y())<1e-9)
460 col.Add(new HVector2(p));
461 }
462 }
463}
464
465// ------------------------------------------------------------------------
466//
467// Calculate the overlapping area of the two hexagons.
468//
469Double_t MHexagon::CalcOverlapArea(const MHexagon &cam) const
470{
471 TOrdCollection col;
472 col.SetOwner();
473
474 GetIntersectionBorder(col, cam);
475
476 // Check if there is an intersection to proceed with
477 const Int_t n = col.GetEntries();
478 if (n==0)
479 return 0;
480
481 // Calculate the center of gravity
482 TVector2 cog;
483
484 TIter Next(&col);
485 HVector2 *v=0;
486 while ((v=(HVector2*)Next()))
487 cog += *v;
488 cog /= n;
489
490 // Shift the figure to its center-og-gravity and
491 // calculate the angle of the connection line between the
492 // border points of our intersesction area and its cog
493 Next.Reset();
494 while ((v=(HVector2*)Next()))
495 {
496 *v -= cog;
497 v->CalcAngle();
498 }
499
500 // Sort these points by this angle
501 col.Sort();
502
503 // Now sum up the area of all the triangles between two
504 // following points and the cog.
505 Double_t A = 0;
506 for (int i=0; i<n; i++)
507 {
508 // Vectors from cog to two nearby border points
509 const HVector2 *v1 = (HVector2*)col.At(i);
510 const HVector2 *v2 = (HVector2*)col.At((i+1)%n);
511
512 // Angle between both vectors
513 const Double_t a = fmod(v1->GetAngle()-v2->GetAngle()+TMath::TwoPi(), TMath::TwoPi());
514
515 // Length of both vectors
516 const Double_t d1 = v1->Mod();
517 const Double_t d2 = v2->Mod();
518
519 A += d1*d2/2*sin(a);
520 }
521 return A;
522}
Note: See TracBrowser for help on using the repository browser.