]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/EVE/AliHLTEveHLT.cxx
Upated calo class to use clustes from ESD
[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 "AliEveHLTEventManager.h"
25 #include "AliHLTGlobalTriggerDecision.h"
26 #include "TEveManager.h"
27 #include "TEvePointSet.h"
28 #include "TEveTrack.h"
29 #include "TCanvas.h"
30 #include "AliESDEvent.h"
31 #include "TEveTrackPropagator.h"
32 #include "AliEveTrack.h"
33 #include "TEveVSDStructs.h"
34 #include "TString.h"
35 #include "TPCLib/tracking-ca/AliHLTTPCCATrackParam.h"
36 #include "TPCLib/tracking-ca/AliHLTTPCCATrackConvertor.h"
37 #include "AliEveMagField.h"
38 #include "TH1F.h"
39 #include "TH2F.h"
40 #include "TThread.h"
41
42 ClassImp(AliHLTEveHLT)
43
44 AliHLTEveHLT::AliHLTEveHLT() : 
45   AliHLTEveBase("TPC tracks"), 
46   fTrueField(kFALSE),
47   fUseIpOnFailedITS(kFALSE),
48   fUseRkStepper(kFALSE),
49   fTrackList(NULL),
50   fOldTrackList(NULL),
51   fPointSetVertex(NULL),
52   fTrCanvas(NULL),
53   fHistPt(NULL), 
54   fHistP(NULL), 
55   fHistEta(NULL),
56   fHistTheta(NULL),
57   fHistPhi(NULL),
58   fHistnClusters(NULL),
59   fHistMult(NULL)
60 {
61   // Constructor.
62   //CreateHistograms();
63 }
64 ///___________________________________________________________________
65 AliHLTEveHLT::~AliHLTEveHLT()
66 {
67   //Destructor, not implemented
68   if(fTrackList)
69     delete fTrackList;
70   fTrackList = NULL;
71 }
72 ///____________________________________________________________________
73 void AliHLTEveHLT::ProcessEsdEvent( AliESDEvent * esd ) {
74   //See header file for documentation
75   if(!fTrackList) CreateTrackList();
76   if(!fPointSetVertex) CreateVertexPointSet();
77   ProcessEsdEvent(esd, fTrackList);
78     
79 }
80 ///_____________________________________________________________________
81 void AliHLTEveHLT::ProcessBlock(AliHLTHOMERBlockDesc * block) {
82   //See header file for documentation
83   if ( ! block->GetDataType().CompareTo("ALIESDV0") ) {
84     if(!fTrackList) CreateTrackList();
85     if(!fPointSetVertex) CreateVertexPointSet();
86     ProcessEsdBlock(block, fTrackList);
87   } 
88   
89   else if ( ! block->GetDataType().CompareTo("ROOTTOBJ") ) {
90     //processROOTTOBJ( block, gHLTText );
91   } 
92
93   else if ( ! block->GetDataType().CompareTo("HLTRDLST") ) {
94     cout << "ignoring hlt rdlst"<<endl;
95     //processHLTRDLST( block );
96   } 
97
98   else if ( ! block->GetDataType().CompareTo("GLOBTRIG") ) {
99     ProcessGlobalTrigger( block );
100   } 
101
102   else if ( !block->GetDataType().CompareTo("ROOTHIST") ) {      
103     if( !fCanvas ) { 
104       fCanvas = CreateCanvas("Primary Vertex", "Primary Vertex");
105       fCanvas->Divide(3, 2);
106     }
107     ProcessHistograms( block , fCanvas);
108   }  
109 }
110
111 ///____________________________________________________________________________
112 void AliHLTEveHLT::UpdateElements() {
113   //See header file for documentation
114   //if(fCanvas) fCanvas->Update();
115   //DrawHistograms();
116   if(fTrackList) fTrackList->ElementChanged();
117   if(fPointSetVertex) fPointSetVertex->ResetBBox();
118
119 }
120 ///_________________________________________________________________________________
121 void AliHLTEveHLT::ResetElements(){
122   //See header file for documentation
123   
124   cout << "destroy"<<endl;
125   if(fTrackList) {
126     fTrackList->Destroy();
127     fTrackList = NULL;
128     // RemoveElement(fTrackList);
129     // TThread * destructor = new TThread(DestroyGarbage, (void*) this);
130     // destructor->Run();
131     // fTrackList = NULL;
132   }
133   
134   if(fPointSetVertex) fPointSetVertex->Reset();
135   cout<< "reset done"<<endl;
136   fHistoCount = 0;
137
138 }
139
140 ///_____________________________________________________________________________________
141 void * AliHLTEveHLT::DestroyGarbage(void * arg) {
142   AliHLTEveHLT * hlt = reinterpret_cast<AliHLTEveHLT*>(arg);
143   if(hlt) hlt->DestroyOldTrackList();
144   return (void*)0;
145 }
146 ///_____________________________________________________________________________________
147 void AliHLTEveHLT::DestroyOldTrackList() {
148   cout << "Destroying the old tracklist's elements"<<endl;
149   fOldTrackList->DestroyElements();
150   cout << "Destroying the old tracklist itself"<<endl;
151   fOldTrackList->Destroy();
152 }
153
154 ///_____________________________________________________________________________________
155 void AliHLTEveHLT::ProcessHistograms(AliHLTHOMERBlockDesc * block, TCanvas * canvas) {
156   //See header file for documentation
157   if ( ! block->GetClassName().CompareTo("TH1F")) {
158     TH1F* histo = reinterpret_cast<TH1F*>(block->GetTObject());
159     if( histo ){
160       TString name(histo->GetName());
161       if( !name.CompareTo("primVertexZ") ){
162         canvas->cd(2);
163         histo->Draw();
164       }else if( !name.CompareTo("primVertexX") ){
165         canvas->cd(3);
166         histo->Draw();
167       }else if( !name.CompareTo("primVertexY") ){
168         canvas->cd(4);
169         histo->Draw();
170       }
171     }
172   }  else if ( ! block->GetClassName().CompareTo("TH2F")) {
173     TH2F *hista = reinterpret_cast<TH2F*>(block->GetTObject());
174     if (hista ){
175        TString name(hista->GetName());
176        if( !name.CompareTo("primVertexXY")) {      
177          canvas->cd(1);
178          hista->Draw();
179        }
180     }
181   }
182   canvas->cd();
183 }
184
185 ///________________________________________________________________________________________
186 void AliHLTEveHLT::CreateTrackList() {
187   //See header file for documentation
188   fTrackList = new TEveTrackList("ESD Tracks");
189   fTrackList->SetMainColor(6);
190   AddElement(fTrackList);
191 }
192
193 ///_________________________________________________________________________________________
194 void AliHLTEveHLT::CreateVertexPointSet() {
195   //See header file for documentation
196   fPointSetVertex = new TEvePointSet("Primary Vertex");
197   fPointSetVertex->SetMainColor(6);
198   fPointSetVertex->SetMarkerStyle((Style_t)kFullStar);
199
200   AddElement(fPointSetVertex);
201 }
202
203 ///________________________________________________________________________
204 void AliHLTEveHLT::ProcessGlobalTrigger( AliHLTHOMERBlockDesc * block ) {
205   //See header file for documentation
206   AliHLTGlobalTriggerDecision * decision = dynamic_cast<AliHLTGlobalTriggerDecision*>(block->GetTObject());
207   decision->Print();
208
209 }
210
211
212 ///______________________________________________________________________
213 void AliHLTEveHLT::ProcessEsdBlock( AliHLTHOMERBlockDesc * block, TEveTrackList * cont ) {
214   //See header file for documentation
215
216   AliESDEvent* esd = (AliESDEvent *) (block->GetTObject());
217   if (!esd) return;
218   
219   ProcessEsdEvent(esd, cont);
220 }
221 ///___________________________________________________________________________________
222 void AliHLTEveHLT::ProcessEsdEvent(AliESDEvent * esd, TEveTrackList * cont) {
223
224   esd->GetStdContent();
225
226   cout << "ProcessESDEvent() :"<< esd->GetEventNumberInFile()<< "  " << esd->GetNumberOfCaloClusters() << " tracks : " << esd->GetNumberOfTracks() << endl;
227
228   //fEventManager->SetRunNumber(esd->GetRunNumber());
229
230   Double_t vertex[3];
231   const AliESDVertex * esdVertex = esd->GetPrimaryVertex();
232   
233   if(esdVertex) {
234     esdVertex->GetXYZ(vertex);
235     fPointSetVertex->SetNextPoint(vertex[0], vertex[1], vertex[2]);
236   }
237   
238   SetUpTrackPropagator(cont->GetPropagator(),-0.1*esd->GetMagneticField(), 520);
239
240   for (Int_t iter = 0; iter < esd->GetNumberOfTracks(); ++iter) {
241     AliEveTrack* track = dynamic_cast<AliEveTrack*>(MakeEsdTrack(esd->GetTrack(iter), cont));
242     cont->AddElement(track);
243    
244     // fHistPt->Fill(esd->GetTrack(iter)->Pt());   // KK
245     // fHistP->Fill(esd->GetTrack(iter)->P()*esd->GetTrack(iter)->Charge());
246     // fHistEta->Fill(esd->GetTrack(iter)->Eta());
247     // fHistTheta->Fill(esd->GetTrack(iter)->Theta()*TMath::RadToDeg());
248     // fHistPhi->Fill(esd->GetTrack(iter)->Phi()*TMath::RadToDeg());
249     // if(esd->GetTrack(iter)->GetStatus()&AliESDtrack::kTPCin || (esd->GetTrack(iter)->GetStatus()&AliESDtrack::kTPCin && esd->GetTrack(iter)->GetStatus()&AliESDtrack::kITSin)){
250     //    fHistnClusters->Fill(esd->GetTrack(iter)->GetTPCNcls());  
251     //}
252   }
253   
254 //fHistMult->Fill(esd->GetNumberOfTracks()); // KK
255   
256   
257   cont->SetTitle(Form("N=%d", esd->GetNumberOfTracks()) );
258   cont->MakeTracks();
259   
260 }
261
262
263
264 void AliHLTEveHLT::DrawHistograms(){
265   //See header file for documentation
266
267   if (!fTrCanvas) {
268     fTrCanvas = CreateCanvas("TPC Tr QA", "TPC Track QA");
269     fTrCanvas->Divide(4, 2);
270   }
271
272   Int_t icd = 1;
273   fTrCanvas->cd(icd++);
274   fHistPt->Draw();
275   fTrCanvas->cd(icd++);
276   fHistP->Draw();
277   fTrCanvas->cd(icd++);
278   fHistEta->Draw();
279   fTrCanvas->cd(icd++);
280   fHistTheta->Draw();
281   fTrCanvas->cd(icd++);
282   fHistPhi->Draw();
283   fTrCanvas->cd(icd++);
284   fHistnClusters->Draw();
285   fTrCanvas->cd(icd++);
286   fHistMult->Draw();
287   fTrCanvas->cd();
288
289   fTrCanvas->Update();
290
291 }
292
293 AliEveTrack* AliHLTEveHLT::MakeEsdTrack (AliESDtrack *at, TEveTrackList* cont) {
294   //See header file for documentation
295
296
297   const double kCLight = 0.000299792458;
298   double bz = - kCLight*10.*( cont->GetPropagator()->GetMagField(0,0,0).fZ);
299
300   Bool_t innerTaken = kFALSE;
301   if ( ! at->IsOn(AliESDtrack::kITSrefit) && fUseIpOnFailedITS)
302   {
303     //tp = at->GetInnerParam();
304     innerTaken = kTRUE;
305   }
306
307   // Add inner/outer track parameters as path-marks.
308
309   Double_t     pbuf[3], vbuf[3];
310
311   AliExternalTrackParam trackParam = *at;
312
313   // take parameters constrained to vertex (if they are)
314
315   if( at->GetConstrainedParam() ){
316     trackParam = *at->GetConstrainedParam();
317   }
318   else if( at->GetInnerParam() ){
319     trackParam = *(at->GetInnerParam());
320   }
321   if( at->GetStatus()&AliESDtrack::kTRDin ){
322     // transport to TRD in
323     trackParam = *at;
324     trackParam.PropagateTo( 290.45, -10.*( cont->GetPropagator()->GetMagField(0,0,0).fZ) );
325   }
326
327   TEveRecTrack rt;
328   {
329     rt.fLabel  = at->GetLabel();
330     rt.fIndex  = (Int_t) at->GetID();
331     rt.fStatus = (Int_t) at->GetStatus();
332     rt.fSign   = (Int_t) trackParam.GetSign();  
333     trackParam.GetXYZ(vbuf);
334     trackParam.GetPxPyPz(pbuf);    
335     rt.fV.Set(vbuf);
336     rt.fP.Set(pbuf);
337     Double_t ep = at->GetP(), mc = at->GetMass();
338     rt.fBeta = ep/TMath::Sqrt(ep*ep + mc*mc);
339   }
340
341   AliEveTrack* track = new AliEveTrack(&rt, cont->GetPropagator());
342   track->SetAttLineAttMarker(cont);
343   track->SetName(Form("AliEveTrack %d", at->GetID()));
344   track->SetElementTitle(CreateTrackTitle(at));
345   track->SetSourceObject(at);
346
347
348   // Set reference points along the trajectory
349   // and the last point
350
351   { 
352     TEvePathMark startPoint(TEvePathMark::kReference);
353     trackParam.GetXYZ(vbuf);
354     trackParam.GetPxPyPz(pbuf);    
355     startPoint.fV.Set(vbuf);
356     startPoint.fP.Set(pbuf);
357     rt.fV.Set(vbuf);
358     rt.fP.Set(pbuf);
359     Double_t ep = at->GetP(), mc = at->GetMass();
360     rt.fBeta = ep/TMath::Sqrt(ep*ep + mc*mc);
361
362     track->AddPathMark( startPoint );    
363   }
364
365
366   if( at->GetTPCPoints(2)>80 ){
367   
368     //
369     // use AliHLTTPCCATrackParam propagator 
370     // since AliExternalTrackParam:PropagateTo()
371     // has an offset at big distances
372     //
373     
374     AliHLTTPCCATrackParam t;
375     AliHLTTPCCATrackConvertor::SetExtParam( t, trackParam );
376     
377     Double_t x0 = trackParam.GetX();
378     Double_t dx = at->GetTPCPoints(2) - x0;
379     
380     //
381     // set a reference at the half of trajectory for better drawing
382     //
383     
384     for( double dxx=dx/2; TMath::Abs(dxx)>=1.; dxx*=.9 ){
385       if( !t.TransportToX(x0+dxx, bz, .999 ) ) continue;
386       AliHLTTPCCATrackConvertor::GetExtParam( t, trackParam, trackParam.GetAlpha() ); 
387       trackParam.GetXYZ(vbuf);
388       trackParam.GetPxPyPz(pbuf);
389       TEvePathMark midPoint(TEvePathMark::kReference);
390       midPoint.fV.Set(vbuf);
391       midPoint.fP.Set(pbuf);    
392       track->AddPathMark( midPoint );
393       break;
394     }
395     
396     //
397     // Set a reference at the end of the trajectory
398     // and a "decay point", to let the event display know where the track ends
399     //
400     
401     for( ; TMath::Abs(dx)>=1.; dx*=.9 ){
402       if( !t.TransportToX(x0+dx, bz, .999 ) ) continue;
403       AliHLTTPCCATrackConvertor::GetExtParam( t, trackParam, trackParam.GetAlpha() ); 
404       trackParam.GetXYZ(vbuf);
405       trackParam.GetPxPyPz(pbuf);
406       TEvePathMark endPoint(TEvePathMark::kReference);
407       TEvePathMark decPoint(TEvePathMark::kDecay);
408       endPoint.fV.Set(vbuf);
409       endPoint.fP.Set(pbuf);
410       decPoint.fV.Set(vbuf);
411       decPoint.fP.Set(pbuf);
412       track->AddPathMark( endPoint );
413       track->AddPathMark( decPoint );
414       break;
415     }  
416   }
417
418   if (at->IsOn(AliESDtrack::kTPCrefit))
419   {
420     if ( ! innerTaken)
421     {
422       AddTrackParamToTrack(track, at->GetInnerParam());
423     }
424     AddTrackParamToTrack(track, at->GetOuterParam());
425   }
426   return track;
427 }
428
429 void AliHLTEveHLT::SetUpTrackPropagator(TEveTrackPropagator* trkProp, Float_t magF, Float_t maxR) {
430   //See header file for documentation
431
432   if (fTrueField) {
433     trkProp->SetMagFieldObj(new AliEveMagField);
434   
435   } else {
436     trkProp->SetMagField(magF);
437   }
438  
439   if (fUseRkStepper) {
440     trkProp->SetStepper(TEveTrackPropagator::kRungeKutta);
441   }
442
443   trkProp->SetMaxR(maxR);
444 }
445
446
447 void AliHLTEveHLT::AddTrackParamToTrack(AliEveTrack* track, const AliExternalTrackParam* tp) {
448   //See header file for documentation
449
450   if (tp == 0)
451     return;
452
453   Double_t pbuf[3], vbuf[3];
454   tp->GetXYZ(vbuf);
455   tp->GetPxPyPz(pbuf);
456
457   TEvePathMark pm(TEvePathMark::kReference);
458   pm.fV.Set(vbuf);
459   pm.fP.Set(pbuf);
460   track->AddPathMark(pm);
461 }
462
463
464
465 TString AliHLTEveHLT::CreateTrackTitle(AliESDtrack* t) {
466   // Add additional track parameters as a path-mark to track.
467
468   TString s;
469
470   Int_t label = t->GetLabel(), index = t->GetID();
471   TString idx(index == kMinInt ? "<undef>" : Form("%d", index));
472   TString lbl(label == kMinInt ? "<undef>" : Form("%d", label));
473
474   Double_t p[3], v[3];
475   t->GetXYZ(v);
476   t->GetPxPyPz(p);
477   Double_t pt    = t->Pt();
478   Double_t ptsig = TMath::Sqrt(t->GetSigma1Pt2());
479   Double_t ptsq  = pt*pt;
480   Double_t ptm   = pt / (1.0 + pt*ptsig);
481   Double_t ptM   = pt / (1.0 - pt*ptsig);
482
483   s = Form("Index=%s, Label=%s\nChg=%d, Pdg=%d\n"
484            "pT = %.3f + %.3f - %.3f [%.3f]\n"
485            "P  = (%.3f, %.3f, %.3f)\n"
486            "V  = (%.3f, %.3f, %.3f)\n",
487            idx.Data(), lbl.Data(), t->Charge(), 0,
488            pt, ptM - pt, pt - ptm, ptsig*ptsq,
489            p[0], p[1], p[2],
490            v[0], v[1], v[2]);
491
492   Int_t   o;
493   s += "Det (in,out,refit,pid):\n";
494   o  = AliESDtrack::kITSin;
495   s += Form("ITS (%d,%d,%d,%d)  ",  t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
496   o  = AliESDtrack::kTPCin;
497   s += Form("TPC(%d,%d,%d,%d)\n",   t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
498   o  = AliESDtrack::kTRDin;
499   s += Form("TRD(%d,%d,%d,%d) ",    t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
500   o  = AliESDtrack::kTOFin;
501   s += Form("TOF(%d,%d,%d,%d)\n",   t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
502   o  = AliESDtrack::kHMPIDout;
503   s += Form("HMPID(out=%d,pid=%d)\n", t->IsOn(o), t->IsOn(o<<1));
504   s += Form("ESD pid=%d", t->IsOn(AliESDtrack::kESDpid));
505
506   if (t->IsOn(AliESDtrack::kESDpid))
507   {
508     Double_t pid[5];
509     t->GetESDpid(pid);
510     s += Form("\n[%.2f %.2f %.2f %.2f %.2f]", pid[0], pid[1], pid[2], pid[3], pid[4]);
511   }
512
513   return s;
514 }
515
516