]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveJetPlaneGL.cxx
Temporary fix to avoid xrootd thrashing
[u/mrichter/AliRoot.git] / EVE / EveDet / 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 <EveDet/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 ClassImp(AliEveJetPlaneGL)
30
31 AliEveJetPlaneGL::AliEveJetPlaneGL() : TGLObject(), fM(0)
32 {
33   fDLCache = kFALSE; // Disable display list -- axis pain.
34 }
35
36 AliEveJetPlaneGL::~AliEveJetPlaneGL()
37 {}
38
39 /******************************************************************************/
40
41 Bool_t AliEveJetPlaneGL::SetModel(TObject* obj, const Option_t* /*opt*/)
42 {
43   if(SetModelCheckClass(obj, AliEveJetPlane::Class())) {
44     fM = dynamic_cast<AliEveJetPlane*>(obj);
45     return kTRUE;
46   }
47   return kFALSE;
48 }
49
50 void AliEveJetPlaneGL::SetBBox()
51 {
52   // !! This ok if master sub-classed from TAttBBox
53   SetAxisAlignedBBox(((AliEveJetPlane*)fExternalObj)->AssertBBox());
54 }
55
56 /******************************************************************************/
57
58 void AliEveJetPlaneGL::DirectDraw(TGLRnrCtx & /*rnrCtx*/) const
59 {
60
61   Float_t minEta = (fM->fMinEta)*(fM->fEtaScale);
62   Float_t maxEta = (fM->fMaxEta)*(fM->fEtaScale);
63   Float_t minPhi = (fM->fMinPhi)*(fM->fPhiScale);
64   Float_t maxPhi = (fM->fMaxPhi)*(fM->fPhiScale);
65   Float_t phiCoord, etaCoord, dPhi, dEta;
66   Double_t eta, phi, E, x, y;
67
68   // Show frame for Eta-Phi coordinates
69
70   glBegin(GL_LINE_LOOP);
71   glVertex3f( minEta, minPhi, 0);
72   glVertex3f( maxEta, minPhi, 0);
73   glVertex3f( maxEta, maxPhi, 0);
74   glVertex3f( minEta, maxPhi, 0);
75   glEnd();
76
77   // Show grid in Eta-Phi coordinates
78
79   dPhi = (maxPhi-minPhi)/(fM->fNPhiDiv - 1);
80   dEta = (maxEta-minEta)/(fM->fNEtaDiv - 1);
81
82   for (int count=1; count < fM->fNPhiDiv-1; count++)
83   {
84     phiCoord = minPhi + count*dPhi;
85     glBegin(GL_LINES);
86     glVertex3f( minEta, phiCoord, 0);
87     glVertex3f( maxEta, phiCoord, 0);
88     glEnd();
89   }
90
91   for (int count=1; count < fM->fNEtaDiv-1; count++)
92   {
93     etaCoord = minEta + count*dEta;
94     glBegin(GL_LINES);
95     glVertex3f( etaCoord, minPhi, 0);
96     glVertex3f( etaCoord, maxPhi, 0);
97     glEnd();
98   }
99
100   // Show axis tick marks and labels
101
102   {
103     TGLCapabilitySwitch lights_off(GL_LIGHTING, false);
104
105     TGLAxis ap;
106     ap.SetLineColor(fM->fGridColor);
107     ap.SetTextColor(fM->fGridColor);
108     TGLVector3 start, end;
109
110     start.Set(minEta, minPhi, 0);
111     end.Set(maxEta, minPhi, 0);
112     ap.PaintGLAxis(start.CArr(), end.CArr(), fM->fMinEta, fM->fMaxEta, 205);
113
114     start.Set(maxEta, minPhi, 0);
115     end.Set(maxEta, maxPhi, 0);
116     ap.PaintGLAxis(start.CArr(), end.CArr(), fM->fMinPhi, fM->fMaxPhi, 205);
117
118   }
119
120   // Finding the maximum energy
121
122   std::vector<AliAODTrack>::iterator k = fM->fTracks.begin();
123   std::vector<AliAODJet>::iterator   j = fM->fJets.begin();
124
125   Double_t eJetMax = 0., eTrackMax = 0., eMax;
126
127   while (j != fM->fJets.end())
128   {
129     if (j->E() > eJetMax) eJetMax = j->E();
130     ++j;
131   }
132
133   while (k != fM->fTracks.end())
134   {
135     if (k->E() > eTrackMax) eTrackMax = k->E();
136     ++k;
137   }
138
139   eMax = eJetMax > eTrackMax ? eJetMax : eTrackMax;
140
141   // Rendering Jets and Tracks in the Eta-Phi plane
142
143   Int_t     nCol = gStyle->GetNumberOfColors();
144   Float_t col[4] = { 0., 0., 0., 0.75};
145
146   glBlendFunc(GL_SRC_ALPHA,GL_ONE);
147   TGLUtil::SetDrawQuality(6);
148
149
150   glPushName(0);
151   if (fM->fRnrJets)
152   {
153     glEnable(GL_BLEND);            // Turn Blending On
154     glDisable(GL_DEPTH_TEST); // Turn Depth Testing Off
155     UInt_t jetid = 0;
156
157
158     j = fM->fJets.begin();
159
160     while (j != fM->fJets.end())
161     {
162       eta = j->Eta();
163       phi = j->Phi();
164       E   = j->E();
165
166       x = eta*(fM->fEtaScale);
167       y = phi*(fM->fPhiScale);
168
169       Int_t colBin = TMath::Min((Int_t) ((nCol-2)*E*
170                                          TMath::Power(10.,fM->fEnergyColorScale)/(eMax)),nCol-2);
171       Int_t colIdx = gStyle->GetColorPalette(colBin);
172       TColor* c    = gROOT->GetColor(colIdx);
173
174       if(c)
175       {
176         col[0] = c->GetRed();
177         col[1] = c->GetGreen();
178         col[2] = c->GetBlue();
179       }
180
181       glLoadName(jetid);
182       TGLUtil::DrawLine( TGLVertex3(x,y,0.),
183                          TGLVector3(0.,0.,TMath::Log(E + 1.)*fM->fEnergyScale),
184                          TGLUtil::kLineHeadArrow, 25.0, col);
185
186       ++j; ++jetid;
187     }
188   }
189
190   col[3] = 1.0;
191
192   glPushName(1);
193   if(fM->fRnrTracks)
194   {
195     glDisable(GL_BLEND);                   // Turn Blending Off
196     glEnable(GL_DEPTH_TEST);     // Turn Depth Testing On
197     UInt_t trackid = 0;
198
199
200     k = fM->fTracks.begin();
201
202     while (k != fM->fTracks.end())
203     {
204       eta = k->Eta();
205       phi = k->Phi();
206       E   = k->E();
207
208       if (E < 0.)
209       {
210         //                      printf(" WARNING: Particle with negative energy has been found.\n");
211         //                      printf(" PARTICLE NOT DISPLAYED. TrackID: %i\n", trackid);
212         ++k; ++trackid;
213         continue;
214       }
215
216       x = eta*(fM->fEtaScale);
217       y = phi*(fM->fPhiScale);
218
219       Int_t colBin = TMath::Min((Int_t) ((nCol-2)*E*
220                                          TMath::Power(10.,fM->fEnergyColorScale)/(eMax)),nCol-2);
221       Int_t colIdx = gStyle->GetColorPalette(colBin);
222       TColor* c    = gROOT->GetColor(colIdx);
223
224       if(c)
225       {
226         col[0] = c->GetRed();
227         col[1] = c->GetGreen();
228         col[2] = c->GetBlue();
229       }
230
231       glLoadName(trackid);
232       TGLUtil::DrawLine( TGLVertex3(x,y,0.),
233                          TGLVector3(0.,0.,TMath::Log(E + 1.)*fM->fEnergyScale),
234                          TGLUtil::kLineHeadArrow, 5.0, col);
235
236       ++k; ++trackid;
237     }
238   }
239
240   glPopName();
241   TGLUtil::ResetDrawQuality();
242
243
244
245 }
246
247 /******************************************************************************/
248
249 void AliEveJetPlaneGL::ProcessSelection(TGLRnrCtx & /*rnrCtx*/, TGLSelectRecord & rec)
250 {
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)];
272       printf("TEveTrack 4-momentum: %f, %f, %f, %f \n", v.Px(),v.Py(),v.Pz(),v.Pt() );
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
318       printf("Jet: %i, TEveTrack: %i \n", jet1State, track1State);
319       printf("Jet: %i, TEveTrack: %i \n\n", jet2State, track2State);
320
321       if(jet1State && jet2State)
322       {
323         Double_t Eta1, Eta2, Phi1, Phi2, d;
324
325         Eta1 = (fM->GetJet1()).Eta();
326         Eta2 = (fM->GetJet2()).Eta();
327         Phi1 = (fM->GetJet1()).Phi();
328         Phi2 = (fM->GetJet2()).Phi();
329
330         d = TMath::Sqrt(TMath::Power(Eta2-Eta1,2) + TMath::Power(Phi2-Phi1,2));
331
332         printf("Eta-Phi: %f, %f\n", Eta1, Phi1 );
333         printf("Eta-Phi: %f, %f\n", Eta2, Phi2 );
334         printf("Eta-Phi Distance: %f\n", d);
335       }
336
337       fM->SetSelectionFlag(1);
338     }
339
340   }
341
342 }
343
344
345
346