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