]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/Alieve/JetPlaneGL.cxx
Add some class docs.
[u/mrichter/AliRoot.git] / EVE / Alieve / JetPlaneGL.cxx
CommitLineData
4673ff03 1// $Header$
2
3#include "JetPlaneGL.h"
4#include <Alieve/JetPlane.h>
5#include <Reve/GLUtilNS.h>
6
7#include <TGLRnrCtx.h>
8#include <TGLSelectRecord.h>
9#include <TGLIncludes.h>
10
11#include "TGLUtil.h"
12#include "TGLAxis.h"
13
14#include <TColor.h>
15#include <TStyle.h>
16#include <TROOT.h>
17
18using namespace Reve;
19using namespace Alieve;
20
21//______________________________________________________________________
22// JetPlaneGL
23//
24
25ClassImp(JetPlaneGL)
26
27JetPlaneGL::JetPlaneGL() : TGLObject(), fM(0)
28{
29 fDLCache = kFALSE; // Disable display list -- axis pain.
30}
31
32JetPlaneGL::~JetPlaneGL()
33{}
34
35/**************************************************************************/
36
37Bool_t JetPlaneGL::SetModel(TObject* obj, const Option_t* /*opt*/)
38{
39 if(SetModelCheckClass(obj, Alieve::JetPlane::Class())) {
40 fM = dynamic_cast<JetPlane*>(obj);
41 return kTRUE;
42 }
43 return kFALSE;
44}
45
46void JetPlaneGL::SetBBox()
47{
48 // !! This ok if master sub-classed from TAttBBox
49 SetAxisAlignedBBox(((JetPlane*)fExternalObj)->AssertBBox());
50}
51
52/**************************************************************************/
53
54void JetPlaneGL::DirectDraw(TGLRnrCtx & /*rnrCtx*/) const
55{
56
57 Float_t minEta = (fM->fMinEta)*(fM->fEtaScale);
58 Float_t maxEta = (fM->fMaxEta)*(fM->fEtaScale);
59 Float_t minPhi = (fM->fMinPhi)*(fM->fPhiScale);
60 Float_t maxPhi = (fM->fMaxPhi)*(fM->fPhiScale);
61 Float_t phiCoord, etaCoord, dPhi, dEta;
62 Double_t eta, phi, E, x, y;
63
64 // Show frame for Eta-Phi coordinates
65
66 glBegin(GL_LINE_LOOP);
67 glVertex3f( minEta, minPhi, 0);
68 glVertex3f( maxEta, minPhi, 0);
69 glVertex3f( maxEta, maxPhi, 0);
70 glVertex3f( minEta, maxPhi, 0);
71 glEnd();
72
73 // Show grid in Eta-Phi coordinates
74
75 dPhi = (maxPhi-minPhi)/(fM->fNPhiDiv - 1);
76 dEta = (maxEta-minEta)/(fM->fNEtaDiv - 1);
77
78 for (int count=1; count < fM->fNPhiDiv-1; count++)
79 {
80 phiCoord = minPhi + count*dPhi;
81 glBegin(GL_LINES);
82 glVertex3f( minEta, phiCoord, 0);
83 glVertex3f( maxEta, phiCoord, 0);
84 glEnd();
85 }
86
87 for (int count=1; count < fM->fNEtaDiv-1; count++)
88 {
89 etaCoord = minEta + count*dEta;
90 glBegin(GL_LINES);
91 glVertex3f( etaCoord, minPhi, 0);
92 glVertex3f( etaCoord, maxPhi, 0);
93 glEnd();
94 }
95
96 // Show axis tick marks and labels
97
98 {
99 GLUtilNS::GL_Capability_Switch lights_off(GL_LIGHTING, false);
100
101 TGLAxis ap;
102 ap.SetLineColor(fM->fGridColor);
103 ap.SetTextColor(fM->fGridColor);
104 TGLVector3 start, end;
105
106 start.Set(minEta, minPhi, 0);
107 end.Set(maxEta, minPhi, 0);
108 ap.PaintGLAxis(start.CArr(), end.CArr(), fM->fMinEta, fM->fMaxEta, 205);
109
110 start.Set(maxEta, minPhi, 0);
111 end.Set(maxEta, maxPhi, 0);
112 ap.PaintGLAxis(start.CArr(), end.CArr(), fM->fMinPhi, fM->fMaxPhi, 205);
113
114 }
115
116 // Finding the maximum energy
117
118 std::vector<AliAODTrack>::iterator k = fM->fTracks.begin();
119 std::vector<AliAODJet>::iterator j = fM->fJets.begin();
120
121 Double_t eJetMax = 0., eTrackMax = 0., eMax;
122
123 while (j != fM->fJets.end())
124 {
125 if (j->E() > eJetMax) eJetMax = j->E();
126 ++j;
127 }
128
129 while (k != fM->fTracks.end())
130 {
131 if (k->E() > eTrackMax) eTrackMax = k->E();
132 ++k;
133 }
134
135 eMax = eJetMax > eTrackMax ? eJetMax : eTrackMax;
136
137 // Rendering Jets and Tracks in the Eta-Phi plane
138
139 Int_t nCol = gStyle->GetNumberOfColors();
140 Float_t col[4] = { 0., 0., 0., 0.75};
141
142 glBlendFunc(GL_SRC_ALPHA,GL_ONE);
143 TGLUtil::SetDrawQuality(6);
144
145
146 glPushName(0);
147 if (fM->fRnrJets)
148 {
149 glEnable(GL_BLEND); // Turn Blending On
150 glDisable(GL_DEPTH_TEST); // Turn Depth Testing Off
151 UInt_t jetid = 0;
152
153
154 j = fM->fJets.begin();
155
156 while (j != fM->fJets.end())
157 {
158 eta = j->Eta();
159 phi = j->Phi();
160 E = j->E();
161
162 x = eta*(fM->fEtaScale);
163 y = phi*(fM->fPhiScale);
164
165 Int_t colBin = TMath::Min((Int_t) ((nCol-2)*E*
166 TMath::Power(10.,fM->fEnergyColorScale)/(eMax)),nCol-2);
167 Int_t colIdx = gStyle->GetColorPalette(colBin);
168 TColor* c = gROOT->GetColor(colIdx);
169
170 if(c)
171 {
172 col[0] = c->GetRed();
173 col[1] = c->GetGreen();
174 col[2] = c->GetBlue();
175 }
176
177 glLoadName(jetid);
178 TGLUtil::DrawLine( TGLVertex3(x,y,0.),
179 TGLVector3(0.,0.,TMath::Log(E + 1.)*fM->fEnergyScale),
180 TGLUtil::kLineHeadArrow, 25.0, col);
181
182 ++j; ++jetid;
183 }
184 }
185
186 col[3] = 1.0;
187
188 glPushName(1);
189 if(fM->fRnrTracks)
190 {
191 glDisable(GL_BLEND); // Turn Blending Off
192 glEnable(GL_DEPTH_TEST); // Turn Depth Testing On
193 UInt_t trackid = 0;
194
195
196 k = fM->fTracks.begin();
197
198 while (k != fM->fTracks.end())
199 {
200 eta = k->Eta();
201 phi = k->Phi();
202 E = k->E();
203
204 if (E < 0.)
205 {
206 // printf(" WARNING: Particle with negative energy has been found.\n");
207 // printf(" PARTICLE NOT DISPLAYED. TrackID: %i\n", trackid);
208 ++k; ++trackid;
209 continue;
210 }
211
212 x = eta*(fM->fEtaScale);
213 y = phi*(fM->fPhiScale);
214
215 Int_t colBin = TMath::Min((Int_t) ((nCol-2)*E*
216 TMath::Power(10.,fM->fEnergyColorScale)/(eMax)),nCol-2);
217 Int_t colIdx = gStyle->GetColorPalette(colBin);
218 TColor* c = gROOT->GetColor(colIdx);
219
220 if(c)
221 {
222 col[0] = c->GetRed();
223 col[1] = c->GetGreen();
224 col[2] = c->GetBlue();
225 }
226
227 glLoadName(trackid);
228 TGLUtil::DrawLine( TGLVertex3(x,y,0.),
229 TGLVector3(0.,0.,TMath::Log(E + 1.)*fM->fEnergyScale),
230 TGLUtil::kLineHeadArrow, 5.0, col);
231
232 ++k; ++trackid;
233 }
234 }
235
236 glPopName();
237 TGLUtil::ResetDrawQuality();
238
239
240
241}
242
243/**************************************************************************/
244
245void JetPlaneGL::ProcessSelection(TGLRnrCtx & /*rnrCtx*/, TGLSelectRecord & rec)
246{
247 // printf("beep %u\n", rec.GetN());
248 // rec.Print();
249 static Int_t jet1State;
250 static Int_t jet2State;
251 static Int_t track1State;
252 static Int_t track2State;
253
254 if (fM->fOneSelection)
255 {
256 printf("\n");
257
258 if (rec.GetN() == 2)
259 {
260 AliAODJet v = fM->fJets[rec.GetItem(1)];
261 printf("Jet 4-momentum: %f, %f, %f, %f \n", v.Px(),v.Py(),v.Pz(),v.Pt() );
262 printf("Eta-Phi values: %f, %f\n", v.Eta(), v.Phi());
263 }
264
265 if (rec.GetN() == 3)
266 {
267 AliAODTrack v = fM->fTracks[rec.GetItem(2)];
268 printf("Track 4-momentum: %f, %f, %f, %f \n", v.Px(),v.Py(),v.Pz(),v.Pt() );
269 printf("Eta-Phi values: %f, %f\n", v.Eta(), v.Phi());
270 }
271 }
272
273 if (fM->fTwoSelection)
274 {
275
276 if ( fM->fSelectionFlag == 1)
277 {
278 if (rec.GetN() == 2)
279 {
280 fM->SetJet1(&(fM->fJets[rec.GetItem(1)]));
281 jet1State = 1;
282 track1State = 0;
283 }
284
285 if (rec.GetN() == 3)
286 {
287 fM->SetTrack1(&(fM->fTracks[rec.GetItem(1)]));
288 jet1State = 0;
289 track1State = 1;
290 }
291
292 fM->SetSelectionFlag(2);
293
294 return;
295 }
296
297 if ( fM->fSelectionFlag == 2)
298 {
299 printf("\n");
300 if (rec.GetN() == 2)
301 {
302 fM->SetJet2(&(fM->fJets[rec.GetItem(1)]));
303 jet2State = 1;
304 track2State = 0;
305 }
306
307 if (rec.GetN() == 3)
308 {
309 fM->SetTrack2(&(fM->fTracks[rec.GetItem(1)]));
310 jet2State = 0;
311 track2State = 1;
312 }
313
314 printf("Jet: %i, Track: %i \n", jet1State, track1State);
315 printf("Jet: %i, Track: %i \n\n", jet2State, track2State);
316
317 if(jet1State && jet2State)
318 {
319 Double_t Eta1, Eta2, Phi1, Phi2, d;
320
321 Eta1 = (fM->GetJet1()).Eta();
322 Eta2 = (fM->GetJet2()).Eta();
323 Phi1 = (fM->GetJet1()).Phi();
324 Phi2 = (fM->GetJet2()).Phi();
325
326 d = TMath::Sqrt(TMath::Power(Eta2-Eta1,2) + TMath::Power(Phi2-Phi1,2));
327
328 printf("Eta-Phi: %f, %f\n", Eta1, Phi1 );
329 printf("Eta-Phi: %f, %f\n", Eta2, Phi2 );
330 printf("Eta-Phi Distance: %f\n", d);
331 }
332
333 fM->SetSelectionFlag(1);
334 }
335
336 }
337
338}
339
340
341
342