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