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