]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/EVE/AliHLTEveHLT.cxx
- fill the number of cluster per track only for TPC or TPC+ITS tracks (0 entries...
[u/mrichter/AliRoot.git] / HLT / EVE / AliHLTEveHLT.cxx
1 /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        *
3  * ALICE Experiment at CERN, All rights reserved.                         *
4  *                                                                        *
5  * Primary Authors: Svein Lindal <slindal@fys.uio.no   >                  *
6  *                  for The ALICE HLT Project.                            *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /// @file   AliHLTEvePhos.cxx
18 /// @author Svein Lindal <slindal@fys.uio.no>
19 /// @brief  HLT class for the HLT EVE display
20
21 #include "AliHLTEveHLT.h"
22 #include "AliHLTHOMERBlockDesc.h"
23 #include "AliHLTEveBase.h"
24 #include "AliEveHOMERManager.h"
25 #include "TEveManager.h"
26 #include "TEvePointSet.h"
27 #include "TEveTrack.h"
28 #include "TCanvas.h"
29 #include "AliESDEvent.h"
30 #include "TEveTrackPropagator.h"
31 #include "AliEveTrack.h"
32 #include "TEveVSDStructs.h"
33 #include "TString.h"
34 #include "TPCLib/tracking-ca/AliHLTTPCCATrackParam.h"
35 #include "TPCLib/tracking-ca/AliHLTTPCCATrackConvertor.h"
36 #include "AliEveMagField.h"
37 #include "TH1F.h"
38 #include "TH2F.h"
39
40 ClassImp(AliHLTEveHLT)
41
42 AliHLTEveHLT::AliHLTEveHLT() : 
43   AliHLTEveBase(), 
44   fTrueField(kFALSE),
45   fUseIpOnFailedITS(kFALSE),
46   fUseRkStepper(kFALSE),
47   fTrackList(NULL),
48   fTrCanvas(NULL),
49   fHistPt(NULL), 
50   fHistP(NULL), 
51   fHistEta(NULL),
52   fHistTheta(NULL),
53   fHistPhi(NULL),
54   fHistnClusters(NULL),
55   fHistMult(NULL)
56 {
57   // Constructor.
58   CreateHistograms();
59
60 }
61
62 AliHLTEveHLT::~AliHLTEveHLT()
63 {
64   //Destructor, not implemented
65   if(fTrackList)
66     delete fTrackList;
67   fTrackList = NULL;
68 }
69
70 void AliHLTEveHLT::CreateHistograms(){
71   //See header file for documentation
72   fHistPt        = new TH1F("fHistPt",       "transverse momentum",    100, 0, 10); // KK   
73   fHistP         = new TH1F("fHistP",        "signed momentum",        100,-7,  7);        
74   fHistEta       = new TH1F("fHistEta",      "pseudorapidity",         100,-2,  2);        
75   fHistTheta     = new TH1F("fHistTheta",    "polar angle",            180, 0,180);   
76   fHistPhi       = new TH1F("fHistPhi",      "azimuthal angle",        180, 0,360);   
77   fHistnClusters = new TH1F("fHistnClusters","TPC clusters per track", 160, 0,160);
78   fHistMult      = new TH1F("fHistMult",     "event track multiplicity",50, 0, 50);    
79   
80   fHistPt   ->SetXTitle("p_{t} (GeV/c)");   // KK
81   fHistP    ->SetXTitle("P*charge (GeV/c)");
82   fHistEta  ->SetXTitle("#eta");
83   fHistTheta->SetXTitle("#theta (degrees)");
84   fHistPhi  ->SetXTitle("#phi (degrees)");
85
86 }
87
88 void AliHLTEveHLT::ProcessBlock(AliHLTHOMERBlockDesc * block) {
89   //See header file for documentation
90   if ( ! block->GetDataType().CompareTo("ALIESDV0") ) {
91     if(!fTrackList) CreateTrackList();
92     ProcessEsdBlock(block, fTrackList);
93   } 
94   
95   else if ( ! block->GetDataType().CompareTo("ROOTTOBJ") ) {
96     //processROOTTOBJ( block, gHLTText );
97   } 
98
99   else if ( ! block->GetDataType().CompareTo("HLTRDLST") ) {
100     //processHLTRDLST( block );
101   } 
102
103   else if ( !block->GetDataType().CompareTo("ROOTHIST") ) {      
104     if( !fCanvas ) { 
105       fCanvas = CreateCanvas("Primary Vertex", "Primary Vertex");
106       fCanvas->Divide(3, 2);
107     }
108     ProcessHistograms( block , fCanvas);
109   }
110   
111 }
112
113
114 void AliHLTEveHLT::UpdateElements() {
115   //See header file for documentation
116   if(fCanvas) fCanvas->Update();
117   DrawHistograms();
118   if(fTrackList) fTrackList->ElementChanged();
119 }
120
121 void AliHLTEveHLT::ResetElements(){
122     //See header file for documentation
123   if(fTrackList) fTrackList->DestroyElements();
124   fHistoCount = 0;
125
126 }
127
128 void AliHLTEveHLT::ProcessHistograms(AliHLTHOMERBlockDesc * block, TCanvas * canvas) {
129   //See header file for documentation
130   if ( ! block->GetClassName().CompareTo("TH1F")) {
131     TH1F* histo = reinterpret_cast<TH1F*>(block->GetTObject());
132     if( histo ){
133       TString name(histo->GetName());
134       if( !name.CompareTo("primVertexZ") ){
135         canvas->cd(2);
136         histo->Draw();
137       }else if( !name.CompareTo("primVertexX") ){
138         canvas->cd(3);
139         histo->Draw();
140       }else if( !name.CompareTo("primVertexY") ){
141         canvas->cd(4);
142         histo->Draw();
143       }
144     }
145   }  else if ( ! block->GetClassName().CompareTo("TH2F")) {
146     TH2F *hista = reinterpret_cast<TH2F*>(block->GetTObject());
147     if (hista ){
148        TString name(hista->GetName());
149        if( !name.CompareTo("primVertexXY")) {      
150          canvas->cd(1);
151          hista->Draw();
152        }
153     }
154   }
155   canvas->cd();
156
157
158
159
160 }
161
162 void AliHLTEveHLT::CreateTrackList() {
163   //See header file for documentation
164   fTrackList = new TEveTrackList("ESD Tracks");
165   fTrackList->SetMainColor(6);
166   gEve->AddElement(fTrackList);
167 }
168
169
170 void AliHLTEveHLT::ProcessEsdBlock( AliHLTHOMERBlockDesc * block, TEveTrackList * cont ) {
171   //See header file for documentation
172
173   AliESDEvent* esd = (AliESDEvent *) (block->GetTObject());
174   esd->GetStdContent();
175   
176   SetUpTrackPropagator(cont->GetPropagator(),-0.1*esd->GetMagneticField(), 520);
177
178   for (Int_t iter = 0; iter < esd->GetNumberOfTracks(); ++iter) {
179     AliEveTrack* track = dynamic_cast<AliEveTrack*>(MakeEsdTrack(esd->GetTrack(iter), cont));
180     cont->AddElement(track);
181    
182     fHistPt->Fill(esd->GetTrack(iter)->Pt());   // KK
183     fHistP->Fill(esd->GetTrack(iter)->P()*esd->GetTrack(iter)->Charge());
184     fHistEta->Fill(esd->GetTrack(iter)->Eta());
185     fHistTheta->Fill(esd->GetTrack(iter)->Theta()*TMath::RadToDeg());
186     fHistPhi->Fill(esd->GetTrack(iter)->Phi()*TMath::RadToDeg());
187     fHistnClusters->Fill(esd->GetTrack(iter)->GetTPCNcls());  
188   }
189   
190   fHistMult->Fill(esd->GetNumberOfTracks()); // KK
191   
192   
193   cont->SetTitle(Form("N=%d", esd->GetNumberOfTracks()) );
194   cont->MakeTracks();
195
196 }
197
198
199 void AliHLTEveHLT::DrawHistograms(){
200   //See header file for documentation
201
202   if (!fTrCanvas) {
203     fTrCanvas = CreateCanvas("TPC Tr QA", "TPC Track QA");
204     fTrCanvas->Divide(4, 2);
205   }
206
207   Int_t icd = 1;
208   fTrCanvas->cd(icd++);
209   fHistPt->Draw();
210   fTrCanvas->cd(icd++);
211   fHistP->Draw();
212   fTrCanvas->cd(icd++);
213   fHistEta->Draw();
214   fTrCanvas->cd(icd++);
215   fHistTheta->Draw();
216   fTrCanvas->cd(icd++);
217   fHistPhi->Draw();
218   fTrCanvas->cd(icd++);
219   fHistnClusters->Draw();
220   fTrCanvas->cd(icd++);
221   fHistMult->Draw();
222   fTrCanvas->cd();
223
224   fTrCanvas->Update();
225
226 }
227
228 AliEveTrack* AliHLTEveHLT::MakeEsdTrack (AliESDtrack *at, TEveTrackList* cont) {
229   //See header file for documentation
230
231
232   const double kCLight = 0.000299792458;
233   double bz = - kCLight*10.*( cont->GetPropagator()->GetMagField(0,0,0).fZ);
234
235   Bool_t innerTaken = kFALSE;
236   if ( ! at->IsOn(AliESDtrack::kITSrefit) && fUseIpOnFailedITS)
237   {
238     //tp = at->GetInnerParam();
239     innerTaken = kTRUE;
240   }
241
242   // Add inner/outer track parameters as path-marks.
243
244   Double_t     pbuf[3], vbuf[3];
245
246   AliExternalTrackParam trackParam = *at;
247
248   // take parameters constrained to vertex (if they are)
249
250   if( at->GetConstrainedParam() ){
251     trackParam = *at->GetConstrainedParam();
252   }
253   else if( at->GetInnerParam() ){
254     trackParam = *(at->GetInnerParam());
255   }
256   if( at->GetStatus()&AliESDtrack::kTRDin ){
257     // transport to TRD in
258     trackParam = *at;
259     trackParam.PropagateTo( 290.45, -10.*( cont->GetPropagator()->GetMagField(0,0,0).fZ) );
260   }
261
262   TEveRecTrack rt;
263   {
264     rt.fLabel  = at->GetLabel();
265     rt.fIndex  = (Int_t) at->GetID();
266     rt.fStatus = (Int_t) at->GetStatus();
267     rt.fSign   = (Int_t) trackParam.GetSign();  
268     trackParam.GetXYZ(vbuf);
269     trackParam.GetPxPyPz(pbuf);    
270     rt.fV.Set(vbuf);
271     rt.fP.Set(pbuf);
272     Double_t ep = at->GetP(), mc = at->GetMass();
273     rt.fBeta = ep/TMath::Sqrt(ep*ep + mc*mc);
274   }
275
276   AliEveTrack* track = new AliEveTrack(&rt, cont->GetPropagator());
277   track->SetAttLineAttMarker(cont);
278   track->SetName(Form("AliEveTrack %d", at->GetID()));
279   track->SetElementTitle(CreateTrackTitle(at));
280   track->SetSourceObject(at);
281
282
283   // Set reference points along the trajectory
284   // and the last point
285
286   { 
287     TEvePathMark startPoint(TEvePathMark::kReference);
288     trackParam.GetXYZ(vbuf);
289     trackParam.GetPxPyPz(pbuf);    
290     startPoint.fV.Set(vbuf);
291     startPoint.fP.Set(pbuf);
292     rt.fV.Set(vbuf);
293     rt.fP.Set(pbuf);
294     Double_t ep = at->GetP(), mc = at->GetMass();
295     rt.fBeta = ep/TMath::Sqrt(ep*ep + mc*mc);
296
297     track->AddPathMark( startPoint );    
298   }
299
300
301   if( at->GetTPCPoints(2)>80 ){
302   
303     //
304     // use AliHLTTPCCATrackParam propagator 
305     // since AliExternalTrackParam:PropagateTo()
306     // has an offset at big distances
307     //
308     
309     AliHLTTPCCATrackParam t;
310     AliHLTTPCCATrackConvertor::SetExtParam( t, trackParam );
311     
312     Double_t x0 = trackParam.GetX();
313     Double_t dx = at->GetTPCPoints(2) - x0;
314     
315     //
316     // set a reference at the half of trajectory for better drawing
317     //
318     
319     for( double dxx=dx/2; TMath::Abs(dxx)>=1.; dxx*=.9 ){
320       if( !t.TransportToX(x0+dxx, bz, .99 ) ) continue;
321       AliHLTTPCCATrackConvertor::GetExtParam( t, trackParam, trackParam.GetAlpha() ); 
322       trackParam.GetXYZ(vbuf);
323       trackParam.GetPxPyPz(pbuf);
324       TEvePathMark midPoint(TEvePathMark::kReference);
325       midPoint.fV.Set(vbuf);
326       midPoint.fP.Set(pbuf);    
327       track->AddPathMark( midPoint );
328       break;
329     }
330     
331     //
332     // Set a reference at the end of the trajectory
333     // and a "decay point", to let the event display know where the track ends
334     //
335     
336     for( ; TMath::Abs(dx)>=1.; dx*=.9 ){
337       if( !t.TransportToX(x0+dx, bz, .99 ) ) continue;
338       AliHLTTPCCATrackConvertor::GetExtParam( t, trackParam, trackParam.GetAlpha() ); 
339       trackParam.GetXYZ(vbuf);
340       trackParam.GetPxPyPz(pbuf);
341       TEvePathMark endPoint(TEvePathMark::kReference);
342       TEvePathMark decPoint(TEvePathMark::kDecay);
343       endPoint.fV.Set(vbuf);
344       endPoint.fP.Set(pbuf);
345       decPoint.fV.Set(vbuf);
346       decPoint.fP.Set(pbuf);
347       track->AddPathMark( endPoint );
348       track->AddPathMark( decPoint );
349       break;
350     }  
351   }
352
353   if (at->IsOn(AliESDtrack::kTPCrefit))
354   {
355     if ( ! innerTaken)
356     {
357       AddTrackParamToTrack(track, at->GetInnerParam());
358     }
359     AddTrackParamToTrack(track, at->GetOuterParam());
360   }
361   return track;
362 }
363
364 void AliHLTEveHLT::SetUpTrackPropagator(TEveTrackPropagator* trkProp, Float_t magF, Float_t maxR) {
365   //See header file for documentation
366
367   if (fTrueField) {
368     trkProp->SetMagFieldObj(new AliEveMagField);
369   
370   } else {
371     trkProp->SetMagField(magF);
372   }
373  
374   if (fUseRkStepper) {
375     trkProp->SetStepper(TEveTrackPropagator::kRungeKutta);
376   }
377
378   trkProp->SetMaxR(maxR);
379 }
380
381
382 void AliHLTEveHLT::AddTrackParamToTrack(AliEveTrack* track, const AliExternalTrackParam* tp) {
383   //See header file for documentation
384
385   if (tp == 0)
386     return;
387
388   Double_t pbuf[3], vbuf[3];
389   tp->GetXYZ(vbuf);
390   tp->GetPxPyPz(pbuf);
391
392   TEvePathMark pm(TEvePathMark::kReference);
393   pm.fV.Set(vbuf);
394   pm.fP.Set(pbuf);
395   track->AddPathMark(pm);
396 }
397
398
399
400 TString AliHLTEveHLT::CreateTrackTitle(AliESDtrack* t) {
401   // Add additional track parameters as a path-mark to track.
402
403   TString s;
404
405   Int_t label = t->GetLabel(), index = t->GetID();
406   TString idx(index == kMinInt ? "<undef>" : Form("%d", index));
407   TString lbl(label == kMinInt ? "<undef>" : Form("%d", label));
408
409   Double_t p[3], v[3];
410   t->GetXYZ(v);
411   t->GetPxPyPz(p);
412   Double_t pt    = t->Pt();
413   Double_t ptsig = TMath::Sqrt(t->GetSigma1Pt2());
414   Double_t ptsq  = pt*pt;
415   Double_t ptm   = pt / (1.0 + pt*ptsig);
416   Double_t ptM   = pt / (1.0 - pt*ptsig);
417
418   s = Form("Index=%s, Label=%s\nChg=%d, Pdg=%d\n"
419            "pT = %.3f + %.3f - %.3f [%.3f]\n"
420            "P  = (%.3f, %.3f, %.3f)\n"
421            "V  = (%.3f, %.3f, %.3f)\n",
422            idx.Data(), lbl.Data(), t->Charge(), 0,
423            pt, ptM - pt, pt - ptm, ptsig*ptsq,
424            p[0], p[1], p[2],
425            v[0], v[1], v[2]);
426
427   Int_t   o;
428   s += "Det (in,out,refit,pid):\n";
429   o  = AliESDtrack::kITSin;
430   s += Form("ITS (%d,%d,%d,%d)  ",  t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
431   o  = AliESDtrack::kTPCin;
432   s += Form("TPC(%d,%d,%d,%d)\n",   t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
433   o  = AliESDtrack::kTRDin;
434   s += Form("TRD(%d,%d,%d,%d) ",    t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
435   o  = AliESDtrack::kTOFin;
436   s += Form("TOF(%d,%d,%d,%d)\n",   t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
437   o  = AliESDtrack::kHMPIDout;
438   s += Form("HMPID(out=%d,pid=%d)\n", t->IsOn(o), t->IsOn(o<<1));
439   s += Form("ESD pid=%d", t->IsOn(AliESDtrack::kESDpid));
440
441   if (t->IsOn(AliESDtrack::kESDpid))
442   {
443     Double_t pid[5];
444     t->GetESDpid(pid);
445     s += Form("\n[%.2f %.2f %.2f %.2f %.2f]", pid[0], pid[1], pid[2], pid[3], pid[4]);
446   }
447
448   return s;
449 }
450