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