source: trunk/Mars/mgui/MHexagon.cc@ 9937

Last change on this file since 9937 was 9385, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 10.6 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MHexagon.cc,v 1.32 2009-03-04 18:45:26 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: Software Development, 2000-2009
24!
25!
26\* ======================================================================== */
27
28//////////////////////////////////////////////////////////////////////////////
29//
30// MHexagon
31//
32// Class Version 2:
33// - added fPhi
34//
35//////////////////////////////////////////////////////////////////////////////
36
37#include "MHexagon.h"
38
39#include <fstream>
40#include <iostream>
41
42#include <TVector2.h> // TVector2
43#include <TVirtualPad.h> // gPad
44#include <TOrdCollection.h> // TOrdCollection
45
46#include "MLog.h"
47#include "MLogManip.h"
48
49#include "MMath.h"
50
51ClassImp(MHexagon);
52
53using namespace std;
54
55const Double_t MHexagon::fgCos60 = 0.5; // TMath::Cos(60/TMath::RadToDeg());
56const Double_t MHexagon::fgSin60 = sqrt(3.)/2; // TMath::Sin(60/TMath::RadToDeg());
57
58const Double_t MHexagon::fgDx[6] = { fgCos60, 0., -fgCos60, -fgCos60, 0., fgCos60 };
59const Double_t MHexagon::fgDy[6] = { fgSin60/3, fgSin60*2/3, fgSin60/3, -fgSin60/3, -fgSin60*2/3, -fgSin60/3 };
60
61// ------------------------------------------------------------------------
62//
63// default constructor for MHexagon
64//
65MHexagon::MHexagon()
66{
67}
68
69// ------------------------------------------------------------------------
70//
71// normal constructor for MHexagon
72//
73MHexagon::MHexagon(Float_t x, Float_t y, Float_t d, Float_t phi)
74: TAttLine(1, 1, 1), TAttFill(0, 1001), fX(x), fY(y), fD(d), fPhi(phi)
75{
76}
77
78// ------------------------------------------------------------------------
79//
80// copy constructor for MHexagon
81//
82MHexagon::MHexagon(const MHexagon &hexagon) : TObject(hexagon), TAttLine(hexagon), TAttFill(hexagon)
83{
84 fX = hexagon.fX;
85 fY = hexagon.fY;
86 fD = hexagon.fD;
87 fPhi = hexagon.fPhi;
88}
89
90// ------------------------------------------------------------------------
91//
92// copy this hexagon to hexagon
93//
94void MHexagon::Copy(TObject &obj)
95#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
96const
97#endif
98{
99 MHexagon &hex = (MHexagon&) obj;
100
101 TObject::Copy(obj);
102 TAttLine::Copy(hex);
103 TAttFill::Copy(hex);
104
105 hex.fX = fX;
106 hex.fY = fY;
107 hex.fD = fD;
108 hex.fPhi = fPhi;
109}
110
111// ------------------------------------------------------------------------
112//
113// compute the distance of a point (px,py) to the Hexagon
114// this functions needed for graphical primitives, that
115// means without this function you are not able to interact
116// with the graphical primitive with the mouse!!!
117//
118// All calcutations are running in pixel coordinates
119//
120Int_t MHexagon::DistancetoPrimitive(Int_t px, Int_t py, Float_t conv)
121{
122 //FIXME: Rotation phi missing!
123
124 //
125 // compute the distance of the Point to the center of the Hexagon
126 //
127 const Int_t pxhex = gPad->XtoAbsPixel(fX*conv);
128 const Int_t pyhex = gPad->YtoAbsPixel(fY*conv);
129
130 //const Double_t x = TMath::Abs(px-pxhex);
131 //const Double_t y = TMath::Abs(py-pyhex);
132
133 TVector2 v(TMath::Abs(px-pxhex), TMath::Abs(py-pyhex));
134 // FIXME: fPhi or -fPhi?
135 v = v.Rotate(-fPhi); // FIXME: Replace with a precalculates sin/cos vector
136
137 const Double_t x = TMath::Abs(v.X());
138 const Double_t y = TMath::Abs(v.Y());
139
140 const Double_t disthex = TMath::Sqrt(x*x + y*y);
141
142 if (disthex==0)
143 return 0;
144
145 const Double_t cosa = x / disthex;
146 const Double_t sina = y / disthex;
147
148 //
149 // comput the distance of hexagon center to pixel border
150 //
151 const Double_t dx = fD/2 * cosa;
152 const Double_t dy = fD/2 * sina;
153
154 const Int_t pxborder = gPad->XtoAbsPixel((fX + dx)*conv);
155 const Int_t pyborder = gPad->YtoAbsPixel((fY + dy)*conv);
156
157 const Double_t bx = pxhex-pxborder;
158 const Double_t by = pyhex-pyborder;
159
160 const Double_t distborder = TMath::Sqrt(bx*bx + by*by);
161
162 //
163 // compute the distance from the border of Pixel
164 // here in the first implementation is just circle inside
165 //
166 return distborder < disthex ? (int)((disthex-distborder)/conv+1) : 0;
167}
168/*
169// ------------------------------------------------------------------------
170//
171// compute the distance of a point (px,py) to the Hexagon center in world
172// coordinates. Return -1 if inside.
173//
174Float_t MHexagon::DistanceToPrimitive(Float_t px, Float_t py) const
175{
176 //FIXME: Rotation phi missing!
177
178 //
179 // compute the distance of the Point to the center of the Hexagon
180 //
181 //const Double_t dx = px-fX;
182 //const Double_t dy = py-fY;
183
184 TVector2 v(px-fX, py-fY);
185 // FIXME: fPhi or -fPhi?
186 v = v.Rotate(-fPhi); // FIXME: Replace with a precalculates sin/cos vector
187
188 const Double_t dx = v.X();
189 const Double_t dy = v.Y();
190
191 const Double_t disthex = TMath::Sqrt(dx*dx + dy*dy);
192
193 //
194 // Now check if point is outside of hexagon; just check x coordinate
195 // in three coordinate systems: the default one, in which two sides of
196 // the hexagon are paralel to the y axis (see camera displays) and two
197 // more, rotated with respect to that one by +- 60 degrees.
198 //
199
200 if (TMath::Abs(dx) > fD/2.)
201 return disthex;
202
203 const Double_t dx2 = dx*fgCos60 + dy*fgSin60;
204
205 if (TMath::Abs(dx2) > fD/2.)
206 return disthex;
207
208 const Double_t dx3 = dx*fgCos60 - dy*fgSin60;
209
210 if (TMath::Abs(dx3) > fD/2.)
211 return disthex;
212
213 return -1;
214}
215*/
216/*
217Float_t MHexagon::DistanceToPrimitive(Float_t px, Float_t py)
218{
219 //
220 // compute the distance of the Point to the center of the Hexagon
221 //
222 const Double_t dx = px-fX;
223 const Double_t dy = py-fY;
224
225 const Double_t disthex = TMath::Sqrt(dx*dx + dy*dy);
226
227 //
228 // compute the distance from the border of Pixel
229 // here in the first implementation is just circle outside
230 //
231 return fD*0.5772 < disthex ? disthex-fD*0.5772 : -1;
232 //
233 // FIXME: this approximate method results in some photons being
234 // assigned to two or even three different pixels.
235 //
236}
237*/
238
239// ------------------------------------------------------------------------
240//
241// Draw this ellipse with new coordinate
242//
243void MHexagon::DrawHexagon(Float_t x, Float_t y, Float_t d, Float_t phi)
244{
245 MHexagon *newhexagon = new MHexagon(x, y, d, phi);
246
247 TAttLine::Copy(*newhexagon);
248 TAttFill::Copy(*newhexagon);
249
250 newhexagon->SetBit(kCanDelete);
251 newhexagon->AppendPad();
252}
253
254/*
255// ------------------------------------------------------------------------
256//
257// This is the first test of implementing a clickable interface
258// for one pixel
259//
260void MHexagon::ExecuteEvent(Int_t event, Int_t px, Int_t py)
261{
262 switch (event)
263 {
264 case kButton1Down:
265 cout << endl << "kButton1Down" << endl;
266 SetFillColor(2);
267 gPad->Modified();
268 break;
269
270 case kButton1Double:
271 SetFillColor(0);
272 gPad->Modified();
273 break;
274 // case kMouseEnter:
275 // printf ("\n Mouse inside object \n" ) ;
276 // break;
277 }
278}
279*/
280
281// ------------------------------------------------------------------------
282//
283// list this hexagon with its attributes
284//
285void MHexagon::ls(const Option_t *) const
286{
287 TROOT::IndentLevel();
288 Print();
289}
290
291// ------------------------------------------------------------------------
292//
293// paint this hexagon with its attribute
294//
295void MHexagon::Paint(Option_t *)
296{
297 PaintHexagon(fX, fY, fD, fPhi);
298}
299
300// ------------------------------------------------------------------------
301//
302// draw this hexagon with the coordinates
303//
304void MHexagon::PaintHexagon(Float_t inX, Float_t inY, Float_t inD, Float_t phi)
305{
306 //
307 // calculate the positions of the pixel corners
308 //
309 Double_t x[7], y[7];
310 for (Int_t i=0; i<7; i++)
311 {
312 TVector2 v(fgDx[i%6], fgDy[i%6]);
313
314 v = v.Rotate(phi); // FIXME: Replace with a precalculates sin/cos vector
315
316 x[i] = inX + v.X()*inD;
317 y[i] = inY + v.Y()*inD;
318 }
319
320 TAttLine::Modify(); // Change line attributes only if neccessary
321 TAttFill::Modify(); // Change fill attributes only if neccessary
322
323 //
324 // paint the pixel
325 //
326 if (GetFillColor())
327 gPad->PaintFillArea(6, x, y);
328
329 if (GetLineStyle())
330 gPad->PaintPolyLine(7, x, y);
331}
332
333// ------------------------------------------------------------------------
334//
335// print/dump this hexagon with its attributes
336//
337void MHexagon::Print(Option_t *) const
338{
339 cout << "MHexagon: ";
340 cout << "x=" << fX << "mm y=" << fY << "mm r=" << fD << "mm phi=" << TMath::RadToDeg() << "deg" << endl;
341
342 cout << " Line:";
343 cout << " Color=" << GetLineColor() << ",";
344 cout << " Style=" << GetLineStyle() << ",";
345 cout << " Width=" << GetLineWidth() << endl;
346 cout << " Fill:";
347 cout << " Color=" << GetFillColor() << ",";
348 cout << " Style=" << GetFillStyle() << endl;
349}
350
351// ------------------------------------------------------------------------
352//
353// Save primitive as a C++ statement(s) on output stream out
354//
355void MHexagon::SavePrimitive(ostream &out, Option_t *)
356{
357#if ROOT_VERSION_CODE >= ROOT_VERSION(5,12,00)
358 if (gROOT->ClassSaved(Class()))
359 out << " ";
360 else
361 out << " MHexagon *";
362
363 out << "hexagon = new MHexagon(" << fX << "," << fY << "," << fD << "," << fPhi << ");" << endl;
364
365 SaveFillAttributes(out, "hexagon");
366 SaveLineAttributes(out, "hexagon");
367
368 out << " hexagon->Draw();" << endl;
369#endif
370}
371
372void MHexagon::SavePrimitive(ofstream &out, Option_t *)
373{
374#if ROOT_VERSION_CODE < ROOT_VERSION(5,12,00)
375 if (gROOT->ClassSaved(Class()))
376 out << " ";
377 else
378 out << " MHexagon *";
379
380 out << "hexagon = new MHexagon(" << fX << "," << fY << "," << fD << "," << fPhi << ");" << endl;
381
382 SaveFillAttributes(out, "hexagon");
383 SaveLineAttributes(out, "hexagon");
384
385 out << " hexagon->Draw();" << endl;
386#else
387 MHexagon::SavePrimitive(static_cast<ostream&>(out), "");
388#endif
389}
Note: See TracBrowser for help on using the repository browser.