]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCTrackHistoComponent.cxx
registration of AliHLTTPCTrackHistoComponent and AliHLTTPCTrackDumpComponent from...
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCTrackHistoComponent.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Gaute Ovrebekk <ovrebekk@ift.uib.no>                  *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /** @file   AliHLTTPCTrackHistoComponent.cxx
20     @author Gaute Ovrebekk, Matthias Richter
21     @date   
22     @brief  The TPC conformal mapping tracker component.
23 */
24
25
26 #if __GNUC__>= 3
27 using namespace std;
28 #endif
29
30 #include "AliHLTTPCTrackHistoComponent.h"
31 #include "AliHLTTPCTransform.h"
32 #include "AliHLTTPCClusterDataFormat.h"
33 #include "AliHLTTPCTrackletDataFormat.h"
34 #include "AliHLTTPCMemHandler.h"
35 #include "AliHLTTPCDefinitions.h"
36 #include <TFile.h>
37 #include <TString.h>
38 #include "TH1F.h"
39 #include "TObjString.h"
40 #include "TObjArray.h"
41
42 //#include "AliHLTTPC.h"
43 //#include <stdlib.h>
44 //#include <cerrno>
45
46 AliHLTTPCTrackHistoComponent gAliHLTTPCTrackHistoComponent;
47
48 /** ROOT macro for the implementation of ROOT specific class methods */
49 ClassImp(AliHLTTPCTrackHistoComponent)
50
51 AliHLTTPCTrackHistoComponent::AliHLTTPCTrackHistoComponent()
52 :
53 fHistoNClustersOnTracks(NULL),                                  
54   fHistoChargeAllClusters(NULL),                                                                         
55   fHistoChargeUsedClusters(NULL),                                                                        
56   fHistoPT(NULL),                                                                                  
57   fHistoResidual(NULL),                                                                            
58   fHistoTgl(NULL),                                                                                 
59   fHistoNClusters(NULL),
60   fHistoNUsedClusters(NULL),
61   fHistoNTracks(NULL),
62   fHistoQMaxAllClusters(NULL),
63   fHistoQMaxUsedClusters(NULL),
64   fPlotAll(kFALSE),                                                
65   fPlotNClustersOnTracks(kFALSE),                                                
66   fPlotChargeClusters(kFALSE),                                                                                     
67   fPlotChargeUsedClusters(kFALSE),                                                                                 
68   fPlotPT(kFALSE),                                                                                                 
69   fPlotResidual(kFALSE),                                                                                           
70   fPlotTgl(kFALSE),
71   fPlotNClusters(kFALSE),  
72   fPlotNUsedClusters(kFALSE), 
73   fPlotNTracks(kFALSE),
74   fPlotQMaxClusters(kFALSE),
75   fPlotQMaxUsedClusters(kFALSE),
76   fResetPlots(kFALSE),
77   fClusters(),
78   fTracks()
79 {
80   
81   // see header file for class documentation
82   // or
83   // refer to README to build package
84   // or
85   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
86
87 }
88
89 AliHLTTPCTrackHistoComponent::~AliHLTTPCTrackHistoComponent()
90 {
91   // see header file for class documentation
92 }
93
94 // Public functions to implement AliHLTComponent's interface.
95 // These functions are required for the registration process
96
97 const char* AliHLTTPCTrackHistoComponent::GetComponentID()
98 {
99   // see header file for class documentation
100   
101   return "TPCTrackHisto";
102 }
103
104 void AliHLTTPCTrackHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
105 {
106   // see header file for class documentation
107   list.clear();
108   list.push_back( AliHLTTPCDefinitions::fgkClustersDataType );
109   list.push_back( AliHLTTPCDefinitions::fgkTrackSegmentsDataType );
110   list.push_back( AliHLTTPCDefinitions::fgkTracksDataType );
111 }
112
113 AliHLTComponentDataType AliHLTTPCTrackHistoComponent::GetOutputDataType()
114 {
115   // see header file for class documentation
116   return kAliHLTDataTypeHistogram;
117
118 }
119
120 void AliHLTTPCTrackHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
121 {
122   // see header file for class documentation
123   // XXX TODO: Find more realistic values.
124   constBase = 0;
125   inputMultiplier = 1;
126 }
127
128 AliHLTComponent* AliHLTTPCTrackHistoComponent::Spawn()
129 {
130   // see header file for class documentation
131   return new AliHLTTPCTrackHistoComponent;
132 }
133
134 int AliHLTTPCTrackHistoComponent::DoInit( int argc, const char** argv )
135 {
136   fHistoNClustersOnTracks = new TH1F("fHistoNClustersOnTracks","Number of Clusters on Tracks",160,0,160);                                   
137   fHistoChargeAllClusters = new TH1F("fHistoChargeAllClusters","Charge of All Clusters",4000,0,4000);                                                  
138   fHistoChargeUsedClusters = new TH1F("fHistoChargeUsedClusters","Charge of Clusters used on Tracks",4000,0,4000);                                 
139   fHistoPT = new TH1F("fHistoPT","pT of Tracks",100,0,10);                                                                    
140   fHistoResidual = new TH1F("fHistoResidual","Residuals",360,0,360);    //change. Testing                                                               
141   fHistoTgl = new TH1F("fHistoTgl","Tgl of Tracks",900,0,90);  
142   fHistoNClusters = new TH1F("fHistoNClusters","Total number of Clusters in Event",3000,0,3000);  
143   fHistoNUsedClusters = new TH1F("fHistoNUsedClusters","Number of Used Cluster in event",3000,0,3000);  
144   fHistoNTracks = new TH1F("fHistoNTracks","Number of Tracks in Event",10,0,10);  
145   fHistoQMaxAllClusters = new TH1F("fHistoQMaxAllClusters","Charge of All Clusters",4000,0,4000);
146   fHistoQMaxUsedClusters = new TH1F("fHistoQMaxUsedClusters","Charge of Clusters used on Tracks",4000,0,4000);
147   fPlotAll=kFALSE;                                                
148   fPlotNClustersOnTracks=kFALSE;                                                
149   fPlotChargeClusters=kFALSE;                                                                                     
150   fPlotChargeUsedClusters=kFALSE;                                                                                 
151   fPlotPT=kFALSE;                                                                                                 
152   fPlotResidual=kFALSE;                                                                                           
153   fPlotTgl=kFALSE;                 
154   fPlotNClusters=kFALSE;    
155   fPlotNUsedClusters=kFALSE;
156   fPlotNTracks=kFALSE;        
157   fPlotQMaxClusters=kFALSE;
158   fPlotQMaxUsedClusters=kFALSE;
159   fResetPlots=kFALSE;
160
161   int iResult=0;
162   TString configuration="";
163   TString argument="";
164   for (int i=0; i<argc && iResult>=0; i++) {
165     argument=argv[i];
166     if (!configuration.IsNull()) configuration+=" ";
167     configuration+=argument;
168   }
169   
170   if (!configuration.IsNull()) {
171     iResult=Configure(configuration.Data());
172   }  
173   return iResult; 
174 }
175   
176 int AliHLTTPCTrackHistoComponent::DoDeinit()
177 {
178   // see header file for class documentation
179   
180   delete fHistoNClustersOnTracks;                                 
181   delete fHistoChargeAllClusters;                                                                         
182   delete fHistoChargeUsedClusters;                                                                        
183   delete fHistoPT;                                                                                
184   delete fHistoResidual;                                                                            
185   delete fHistoTgl;        
186   delete fHistoNClusters;
187   delete fHistoNUsedClusters;
188   delete fHistoNTracks;
189   delete fHistoQMaxAllClusters;
190   delete fHistoQMaxUsedClusters;
191   
192   return 0;
193 }
194
195 int AliHLTTPCTrackHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
196 {
197   const AliHLTComponentBlockData* iter = NULL;
198
199   if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) )
200     return 0;
201
202   Int_t TotalTrack = 0;
203
204   //Reading Merged Tracks
205   for ( iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkTracksDataType); iter != NULL; iter = GetNextInputBlock() ) {
206     if(iter->fDataType!=AliHLTTPCDefinitions::fgkTracksDataType){continue;}
207     ReadTracks(iter,TotalTrack);  
208   }
209
210   //Reading Tracks form slice
211   for ( iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkTrackSegmentsDataType); iter != NULL; iter = GetNextInputBlock() ) {
212     if(iter->fDataType!=AliHLTTPCDefinitions::fgkTrackSegmentsDataType){continue;}
213     ReadTracks(iter,TotalTrack);
214   }
215
216   int TotalSpacePoint = 0;
217   int nClustersUsed=0;
218  
219   for ( iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); iter != NULL; iter = GetNextInputBlock() ) {
220     
221     if(iter->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType){continue;}
222
223     //AliHLTUInt8_t slice = AliHLTTPCDefinitions::GetMinSliceNr( *iter );
224     //AliHLTUInt8_t patch = AliHLTTPCDefinitions::GetMinPatchNr( *iter );
225
226     //HLTDebug ( "Input Data - TPC cluster - Slice/Patch: %d/%d.", slice, patch );
227     const AliHLTTPCClusterData* clusterData = (const AliHLTTPCClusterData*) iter->fPtr;
228     Int_t nSpacepoint = (Int_t) clusterData->fSpacePointCnt;    
229     TotalSpacePoint += nSpacepoint;
230     //HLTInfo("TrackHisto found %d Spacepoints in slice %d patch %d", nSpacepoint, slice, patch);
231     AliHLTTPCSpacePointData *clusters = (AliHLTTPCSpacePointData*) clusterData->fSpacePoints;
232
233     for(int i=0;i<nSpacepoint;i++){
234       UInt_t idCluster = clusters[i].fID;
235       Int_t sliceCl = (idCluster>>25) & 0x7f;
236       Int_t patchCl = (idCluster>>22) & 0x7;
237       UInt_t pos = idCluster&0x3fffff;
238       if(fPlotChargeClusters || fPlotAll){fHistoChargeAllClusters->Fill(clusters[i].fCharge);}
239       if(fPlotQMaxClusters || fPlotAll){fHistoQMaxAllClusters->Fill(clusters[i].fMaxQ);}
240       for(UInt_t id=0;id<fTrackClusterID[sliceCl][patchCl].size();id++){
241         if(fTrackClusterID[sliceCl][patchCl][id]==pos){
242           clusters[i].fUsed=kTRUE;
243           nClustersUsed++;
244           if(fPlotChargeUsedClusters || fPlotAll){fHistoChargeUsedClusters->Fill(clusters[i].fCharge);}
245           if(fPlotQMaxUsedClusters || fPlotAll){fHistoQMaxUsedClusters->Fill(clusters[i].fMaxQ);}
246         }
247       } 
248       fClusters.push_back(clusters[i]);
249     }
250   } 
251   
252   fHistoNClusters->Fill(TotalSpacePoint);
253   if(TotalTrack>0){fHistoNUsedClusters->Fill(nClustersUsed);}
254   fHistoNTracks->Fill(TotalTrack);
255
256   HLTInfo("TrackHisto found %d Spacepoints",TotalSpacePoint);
257   HLTInfo("TrackHisto found %d Tracks",TotalTrack);
258   
259   PushHisto();
260
261   fClusters.clear();
262   fTracks.clear();
263   for(UInt_t i=0;i<36;i++){
264     for(UInt_t j=0;j<6;j++){ 
265       fTrackClusterID[i][j].clear();
266     }
267   }
268
269   return 0;
270 }
271  
272  int AliHLTTPCTrackHistoComponent::Configure(const char* arguments)
273  {
274    
275    int iResult=0;
276    if (!arguments) return iResult;
277    
278    TString allArgs=arguments;
279    TString argument;
280    
281    TObjArray* pTokens=allArgs.Tokenize(" ");
282    if (pTokens) {
283      for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
284        argument=((TObjString*)pTokens->At(i))->GetString();
285        if (argument.IsNull()) continue;
286
287        if (argument.CompareTo("-plot-All")==0) {
288         HLTInfo("Ploting All Histograms for Tracks");
289         fPlotAll = kTRUE;
290         fPlotNClustersOnTracks=kTRUE;                                                
291         fPlotChargeClusters=kTRUE;                                                                                     
292         fPlotChargeUsedClusters=kTRUE;                                                                                 
293         fPlotPT=kTRUE;                                                                                                 
294         fPlotResidual=kTRUE;                                                                                           
295         fPlotTgl=kTRUE;
296         fPlotNClusters=kTRUE;
297         fPlotNUsedClusters=kTRUE;
298         fPlotNTracks=kTRUE;
299         fPlotQMaxClusters=kTRUE;
300         fPlotQMaxUsedClusters=kTRUE;
301             
302         continue;
303        }
304        else if (argument.CompareTo("-plot-nClusters")==0) {     
305          HLTInfo("Ploting Number of clusters Used on Tracks");
306          fPlotNClustersOnTracks = kTRUE;
307          continue;
308        }
309        else if (argument.CompareTo("-plot-ChargeClusters")==0) {
310          HLTInfo("Ploting Charge of All Clusters");
311          fPlotChargeClusters = kTRUE;
312          continue;
313        }
314        else if (argument.CompareTo("-plot-ChargeUsedClusters")==0) {
315          HLTInfo("Ploting Charge of Clusters Used on Tracks");
316          fPlotChargeUsedClusters = kTRUE; 
317          continue;
318        }
319        else if (argument.CompareTo("-plot-pT")==0) {
320          HLTInfo("Ploting pT of Tracks");
321          fPlotPT=kTRUE;
322          continue;
323        }
324        else if (argument.CompareTo("-plot-Residuals")==0) {
325          HLTInfo("Ploting Residuals");
326          fPlotResidual=kTRUE; 
327          continue;
328        }
329        else if (argument.CompareTo("-plot-Tgl")==0) {
330          HLTInfo("Ploting Tgl of Tracks");
331          fPlotTgl=kTRUE;  
332          continue;
333        }
334        else if (argument.CompareTo("-plot-NClusters")==0) {
335          HLTInfo("Ploting Number of Clusters");
336          fPlotNClusters=kTRUE;  
337          continue;
338        }
339        else if (argument.CompareTo("-plot-NUsedClusters")==0) {
340          HLTInfo("Ploting Number Of Used Clusters");
341          fPlotNUsedClusters=kTRUE;  
342          continue;
343        }
344        else if (argument.CompareTo("-plot-NTracks")==0) {
345          HLTInfo("Ploting Number Of Tracks");
346          fPlotNTracks=kTRUE;  
347          continue;
348        }
349        else if (argument.CompareTo("-plot-QMaxAll")==0) {
350          HLTInfo("Ploting QMax for All Clusters");
351          fPlotQMaxClusters=kTRUE;  
352          continue;
353        }
354        else if (argument.CompareTo("-plot-QMaxUsed")==0) {
355          HLTInfo("Ploting QMax for Used Clusters");
356          fPlotQMaxUsedClusters=kTRUE;  
357          continue;
358        }
359        else if (argument.CompareTo("-reset-plots")==0) {
360          HLTInfo("Reseting plots");
361          fResetPlots=kTRUE;  
362          continue;
363        }
364        else {
365          HLTError("unknown argument %s", argument.Data());
366          iResult=-EINVAL;
367          break;
368        }
369      }
370      delete pTokens;
371    }
372    
373    return iResult;
374  }
375  
376 void AliHLTTPCTrackHistoComponent::ReadTracks(const AliHLTComponentBlockData* iter,Int_t &tt){
377
378   //AliHLTUInt8_t slice = AliHLTTPCDefinitions::GetMinSliceNr( *iter );
379   //AliHLTUInt8_t patch = AliHLTTPCDefinitions::GetMinPatchNr( *iter );
380   
381   //HLTDebug ( "Input Data - TPC cluster - Slice/Patch: %d/%d.", slice, patch );
382   const AliHLTTPCTrackletData* trackData = (const AliHLTTPCTrackletData*) iter->fPtr;
383   AliHLTUInt32_t nTracks = trackData->fTrackletCnt;
384   tt += nTracks;
385   //HLTInfo("TrackHisto found %d Tracks in slice %d patch %d", nTracks, slice, patch);
386   AliHLTTPCTrackSegmentData *tracks = (AliHLTTPCTrackSegmentData*) trackData->fTracklets;
387   
388   for(AliHLTUInt32_t i=0;i<nTracks;i++){
389     fTracks.push_back(tracks[i]);
390     UInt_t nHits = tracks->fNPoints;
391     if(fPlotNClustersOnTracks || fPlotAll){fHistoNClustersOnTracks->Fill(nHits);}
392     if(fPlotPT || fPlotAll){fHistoPT->Fill(tracks[i].fPt);}
393     if(fPlotTgl || fPlotAll){fHistoTgl->Fill(tracks[i].fTgl);}
394     const UInt_t *hitnum = tracks->fPointIDs;
395     for(UInt_t h=0; h<nHits; h++){
396       UInt_t idTrack = hitnum[h];
397       Int_t sliceTrack = (idTrack>>25) & 0x7f;
398       Int_t patchTrack = (idTrack>>22) & 0x7;
399       UInt_t pos = idTrack&0x3fffff;
400       fTrackClusterID[sliceTrack][patchTrack].push_back(pos);
401     }
402     UChar_t *tmpP = (UChar_t*)tracks;
403     tmpP += sizeof(AliHLTTPCTrackSegmentData)+tracks->fNPoints*sizeof(UInt_t);
404     tracks = (AliHLTTPCTrackSegmentData*)tmpP;
405   } 
406 }
407
408 void AliHLTTPCTrackHistoComponent::PushHisto(){
409
410   if(fPlotNClustersOnTracks || fPlotAll){
411     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
412     PushBack( (TObject*) fHistoNClustersOnTracks,kAliHLTDataTypeHistogram, fSpecification);   
413   }
414   if(fPlotChargeClusters || fPlotAll){
415     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
416     PushBack( (TObject*) fHistoChargeAllClusters,kAliHLTDataTypeHistogram, fSpecification);   
417   }
418   if(fPlotChargeUsedClusters || fPlotAll){
419     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
420     PushBack( (TObject*) fHistoChargeUsedClusters,kAliHLTDataTypeHistogram, fSpecification);   
421   }
422   if(fPlotPT || fPlotAll){
423     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
424     PushBack( (TObject*) fHistoPT,kAliHLTDataTypeHistogram, fSpecification);   
425   }
426   if(fPlotResidual || fPlotAll){
427     for(unsigned int i=0;i<fTracks.size();i++){
428       fHistoResidual->Fill(fTracks[i].fPsi); //Not rigth. Change here and x in Histo. Just for test.
429     }
430     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
431     PushBack( (TObject*) fHistoResidual,kAliHLTDataTypeHistogram, fSpecification);   
432   }
433   if(fPlotTgl || fPlotAll){
434     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
435     PushBack( (TObject*) fHistoTgl,kAliHLTDataTypeHistogram, fSpecification);   
436   }
437   if(fPlotNClusters || fPlotAll){
438     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
439     PushBack( (TObject*) fHistoNClusters,kAliHLTDataTypeHistogram, fSpecification);
440   }
441   if(fPlotNUsedClusters || fPlotAll){
442     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
443     PushBack( (TObject*) fHistoNUsedClusters,kAliHLTDataTypeHistogram, fSpecification);
444   }
445   if(fPlotNTracks || fPlotAll){
446     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
447     PushBack( (TObject*) fHistoNTracks,kAliHLTDataTypeHistogram, fSpecification);
448   }
449   if(fPlotQMaxUsedClusters || fPlotAll){
450     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
451     PushBack( (TObject*) fHistoQMaxAllClusters,kAliHLTDataTypeHistogram, fSpecification);
452   }
453   if(fPlotQMaxUsedClusters || fPlotAll){
454     AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5);
455     PushBack( (TObject*) fHistoQMaxUsedClusters,kAliHLTDataTypeHistogram, fSpecification);
456   }
457   
458   if(fResetPlots){
459     fHistoNClustersOnTracks->Reset();                                 
460     fHistoChargeAllClusters->Reset();                                                                         
461     fHistoChargeUsedClusters->Reset();                                                                        
462     fHistoPT->Reset();                                                                                
463     fHistoResidual->Reset();                                                                            
464     fHistoTgl->Reset();        
465     fHistoNClusters->Reset();     
466     fHistoNUsedClusters->Reset();
467     fHistoNTracks->Reset();
468     fHistoQMaxAllClusters->Reset();
469     fHistoQMaxUsedClusters->Reset();
470   }
471 }