]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCTrackHistoComponent.cxx
coverity warnings 15388 10083 10082 fixed
[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, Kalliopi Kanaki
21     @date   
22     @brief  A histogram component with track and associated cluster properties 
23 */
24
25 #if __GNUC__>= 3
26 using namespace std;
27 #endif
28
29 #include "AliHLTTPCTrackHistoComponent.h"
30 #include "AliHLTTPCTransform.h"
31 #include "AliHLTTPCClusterDataFormat.h"
32 #include "AliHLTTPCTrackletDataFormat.h"
33 #include "AliHLTTPCMemHandler.h"
34 #include "AliHLTTPCDefinitions.h"
35 #include "AliHLTGlobalBarrelTrack.h"
36 #include "AliHLTExternalTrackParam.h"
37 #include "AliHLTDataTypes.h"
38
39 #include <TFile.h>
40 #include <TString.h>
41 #include "TNtuple.h"
42 #include "TH1F.h"
43 #include "TProfile.h"
44 #include "TObjString.h"
45 #include "TObjArray.h"
46
47
48 #include "AliHLTTPCTrack.h"
49
50
51 /** ROOT macro for the implementation of ROOT specific class methods */
52 ClassImp(AliHLTTPCTrackHistoComponent)
53
54 AliHLTTPCTrackHistoComponent::AliHLTTPCTrackHistoComponent()
55     :
56   //, fReset(0)
57     fNEvents(0)
58   , fNtotTracks(0)
59   , fEvtMod(20)
60   , fBufferSize(5000)
61   , fdEdx(kFALSE)
62   , fMeanMultiplicity(NULL)
63   , fMultiplicity(NULL)
64   //, fdNdEta(NULL)
65   , fDeDxVsP(NULL)
66   , fClusters(NULL)
67   , fTracks(NULL)
68   //, fNClusterVsXY(NULL)
69   //, fChargeVsXY(NULL)
70   //, fClustersArray(NULL)
71   //, fNSpacePoints(NULL)
72 {
73   // see header file for class documentation
74   // or
75   // refer to README to build package
76   // or
77   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
78 }
79
80 const char* AliHLTTPCTrackHistoComponent::fgkOCDBEntry="HLT/ConfigTPC/TPCTrackHisto";
81
82 AliHLTTPCTrackHistoComponent::~AliHLTTPCTrackHistoComponent(){
83 // see header file for class documentation
84 }
85
86 // Public functions to implement AliHLTComponent's interface.
87 // These functions are required for the registration process
88
89 const char* AliHLTTPCTrackHistoComponent::GetComponentID(){
90 // see header file for class documentation
91   
92   return "TPCTrackHisto";
93 }
94
95 void AliHLTTPCTrackHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list){
96 // see header file for class documentation
97   
98   list.clear();
99   list.push_back(AliHLTTPCDefinitions::fgkClustersDataType|kAliHLTDataOriginTPC);
100   list.push_back(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
101   //list.push_back(kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC);
102 }
103
104 AliHLTComponentDataType AliHLTTPCTrackHistoComponent::GetOutputDataType(){
105 // see header file for class documentation
106   return kAliHLTMultipleDataType;
107 }
108
109 int AliHLTTPCTrackHistoComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList){
110 // see header file for class documentation
111   tgtList.clear();
112   tgtList.push_back(kAliHLTDataTypeTNtuple|kAliHLTDataOriginTPC);
113   tgtList.push_back(kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
114   return tgtList.size();
115 }
116
117 void AliHLTTPCTrackHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ){
118 // see header file for class documentation
119   
120   constBase = 5000;
121   inputMultiplier = 1;// XXX TODO: Find more realistic value
122 }
123
124 AliHLTComponent* AliHLTTPCTrackHistoComponent::Spawn(){
125 // see header file for class documentation
126   return new AliHLTTPCTrackHistoComponent;
127 }
128
129 int AliHLTTPCTrackHistoComponent::DoInit( int argc, const char** argv ){
130 // see header file for class documentation
131  
132   fClusters = new TNtuple("fClusters", "fClusters", "charge:qmax:residualY:residualZ"); 
133   fTracks   = new TNtuple("fTracks",  "fTracks",  "pt:eta:psi:nclusters"); 
134  
135   fClusters->SetCircular(fBufferSize);
136   fTracks->SetCircular(fBufferSize);
137  
138   fMultiplicity     = new TH1F("fMultiplicity",     "Track multiplicity per event",     1000,           0, 1000);
139   fMeanMultiplicity = new TH1F("fMeanMultiplicity", "Mean track multiplicity vs. #evt", 10000/fEvtMod, 0, 10000);
140   fDeDxVsP          = new TProfile("fDeDxVsP",      "E deposition per unit length vs. p",100, 0, 100);
141   fDeDxVsP->SetXTitle("p (GeV/c)");
142   //fdNdEta = new TH1F("fdNdEta", "dN/d#eta",100,-10,10);
143   
144   // first configure the default
145   int iResult=0;
146   if (iResult>=0) iResult=ConfigureFromCDBTObjString(fgkOCDBEntry);
147
148   // configure from the command line parameters if specified
149   if (iResult>=0 && argc>0)  iResult=ConfigureFromArgumentString(argc, argv);
150   
151   return iResult;
152 }
153   
154 int AliHLTTPCTrackHistoComponent::DoDeinit(){
155 // see header file for class documentation
156   
157   delete fClusters;
158   delete fTracks;
159   
160   delete fMultiplicity;
161   delete fMeanMultiplicity;
162   delete fDeDxVsP;
163   //delete fdNdEta;
164   
165   return 0;
166 }
167
168 int AliHLTTPCTrackHistoComponent::Reconfigure(const char* cdbEntry, const char* /*chainId*/){
169 // see header file for class documentation
170
171   // configure from the specified antry or the default one
172   const char* entry=cdbEntry;
173   if (!entry || entry[0]==0) {
174      entry=fgkOCDBEntry;
175   }
176
177   return ConfigureFromCDBTObjString(entry);
178 }
179
180
181 int AliHLTTPCTrackHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/){
182 // see header file for class documentation
183
184   if(GetFirstInputBlock(kAliHLTDataTypeSOR) || GetFirstInputBlock(kAliHLTDataTypeEOR)) return 0;  
185
186   fNEvents++;
187
188   const AliHLTComponentBlockData *iter = NULL;
189   
190   for(int i=0; i<36; i++){
191      for(int j=0; j<6; j++){
192          fClustersArray[i][j] = NULL;
193          fNSpacePoints[i][j]  = 0;     
194      }  
195   }
196  
197  
198   //----------------- loop over cluster blocks ---------------------//
199   
200   Int_t totalSpacePoints = 0;
201   
202   for(iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); iter != NULL; iter = GetNextInputBlock()){
203             
204       if(iter->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType) continue;
205
206       AliHLTUInt8_t minSlice     = AliHLTTPCDefinitions::GetMinSliceNr(*iter);
207       AliHLTUInt8_t minPartition = AliHLTTPCDefinitions::GetMinPatchNr(*iter);
208       //HLTDebug("Input Data - TPC cluster - slice/partition: %d/%d.", minSlice, minPartition);
209
210       const AliHLTTPCClusterData* clusterData = (const AliHLTTPCClusterData*)iter->fPtr;
211       Int_t nSpacepoint = (Int_t)clusterData->fSpacePointCnt;    
212       totalSpacePoints += nSpacepoint;
213       HLTDebug("TrackHisto component found %d spacepoints in slice %d partition %d", nSpacepoint, minSlice, minPartition);
214       
215       AliHLTTPCSpacePointData *clusters = (AliHLTTPCSpacePointData*)clusterData->fSpacePoints;
216       
217       if(fClustersArray[minSlice][minPartition] != NULL){
218          delete(fClustersArray[minSlice][minPartition]); // ???????
219          fClustersArray[minSlice][minPartition] = NULL;
220       }      
221
222       // fill the array with AliHLTTPCSpacePointData pointers
223       // it will be used in the track loop to access information
224       // for the used clusters only
225       fClustersArray[minSlice][minPartition] = clusters;
226       fNSpacePoints[minSlice][minPartition]  = nSpacepoint;
227       
228       if(nSpacepoint==0) fClustersArray[minSlice][minPartition] = NULL;
229
230   } // end of loop over cluster data blocks
231   
232   HLTDebug("TrackHisto found %d spacepoints",totalSpacePoints);
233   
234   
235   
236   
237   //----------------- loop over merged tracks -------------------//
238
239   Int_t totalTracks = 0;
240   
241   for(iter = GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC); iter != NULL; iter = GetNextInputBlock()){    
242       if(iter->fDataType != (kAliHLTDataTypeTrack|kAliHLTDataOriginTPC)) continue; 
243       ReadTracks(iter,totalTracks);
244   }
245   
246   HLTDebug("TrackHisto found %d tracks", totalTracks);  
247
248   fMultiplicity->Fill(totalTracks);
249   
250   fNtotTracks += totalTracks;
251   
252   if(fNEvents%fEvtMod==0){    
253      fMeanMultiplicity->Fill(fNEvents, Float_t(fNtotTracks)/(fEvtMod));
254      //HLTInfo("-------------- Event number: %d, total tracks accummulated %d", fNEvents, fNtotTracks);
255      fNtotTracks = 0;
256   }
257
258   PushHisto();
259   
260   return 0;
261 }
262  
263 void AliHLTTPCTrackHistoComponent::ReadTracks(const AliHLTComponentBlockData* iter,Int_t &tt){
264 // see header file for class documentation
265
266   //AliHLTUInt8_t slice = AliHLTTPCDefinitions::GetMinSliceNr(*iter);
267   
268   Int_t usedSpacePoints = 0;
269   
270   vector<AliHLTGlobalBarrelTrack> tracksVector;
271   AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(iter->fPtr), iter->fSize, tracksVector);
272   
273   tt = tracksVector.size();
274   
275   for(vector<AliHLTGlobalBarrelTrack>::iterator element=tracksVector.begin();  element!=tracksVector.end(); element++){
276        
277        Double_t trackLength = 0.;
278        if(fdEdx==kTRUE) trackLength = element->GetPathLengthTo( element->GetLastPointX(), 5.0);     
279            
280        UInt_t nHits = element->GetNumberOfPoints();
281        fTracks->Fill( 1./element->OneOverPt(), element->GetSnp(), element->GetTgl(), nHits );  
282        //fdNdEta->Fill(element->GetSnp());
283  
284        Double_t totCharge = 0;
285        const UInt_t *hitnum = element->GetPoints();
286        for(UInt_t i=0; i<element->GetNumberOfPoints(); i++){
287            
288            UInt_t idTrack   = hitnum[i];
289            Int_t sliceTrack = AliHLTTPCSpacePointData::GetSlice(idTrack);
290            Int_t patchTrack = AliHLTTPCSpacePointData::GetPatch(idTrack);
291            UInt_t pos       = AliHLTTPCSpacePointData::GetNumber(idTrack);
292  
293            if( !fClustersArray[sliceTrack][patchTrack] ) continue;             
294            
295            if(sliceTrack<0 || sliceTrack>36 || patchTrack<0 || patchTrack>5 ){
296               HLTError("Corrupted TPC cluster Id: slice %d, patch %d, cluster %d", sliceTrack, patchTrack, idTrack);
297               continue;
298            }
299
300            if(fNSpacePoints[sliceTrack][patchTrack]<=pos ){
301               HLTError("Space point array out of boundaries!");
302               continue;
303            }
304            
305            totCharge += (fClustersArray[sliceTrack][patchTrack])[pos].fCharge; 
306            
307            //Float_t xyz[3]; xyz[0] = xyz[1] = xyz[2] = 0.;
308         
309            //xyz[0] = (fClustersArray[sliceTrack][patchTrack])[pos].fX;
310            //xyz[1] = (fClustersArray[sliceTrack][patchTrack])[pos].fY;
311            //xyz[2] = (fClustersArray[sliceTrack][patchTrack])[pos].fZ;
312         
313            //AliHLTTPCTransform::Local2Global(xyz,slice); 
314            
315            //Double_t p[2]   = { xyz[1], xyz[2] };
316            //Double_t cov[3] = { (fClustersArray[sliceTrack][patchTrack])[pos].fSigmaY2, 0., (fClustersArray[sliceTrack][patchTrack])[pos].fSigmaZ2};  
317            //Double_t *res = element->GetResiduals(p,cov,kFALSE); 
318            
319            //HLTInfo("resy: %f, resz: %f", res[0], res[1]);
320            
321            //if(!res)  res[0] = res[1] = -1000.;
322            //else      fClusters->Fill( (fClustersArray[sliceTrack][patchTrack])[pos].fCharge, (fClustersArray[sliceTrack][patchTrack])[pos].fQMax, res[0], res[1]);
323           
324            fClusters->Fill( (fClustersArray[sliceTrack][patchTrack])[pos].fCharge, (fClustersArray[sliceTrack][patchTrack])[pos].fQMax, -1000., -1000.);
325             
326            usedSpacePoints++;      
327        }        
328   if(fdEdx==kTRUE && trackLength > 0) fDeDxVsP->Fill(element->OneOverPt()*TMath::Sqrt(1.+element->GetTgl()*element->GetTgl()), totCharge/trackLength);       
329   }
330 }
331
332 void AliHLTTPCTrackHistoComponent::PushHisto(){
333 // see header file for class documentation
334     
335     PushBack( (TObject*)fTracks,           kAliHLTDataTypeTNtuple  |kAliHLTDataOriginTPC, 0x0);   
336     PushBack( (TObject*)fClusters,         kAliHLTDataTypeTNtuple  |kAliHLTDataOriginTPC, 0x0);
337     PushBack( (TObject*)fMultiplicity,     kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0);
338     PushBack( (TObject*)fMeanMultiplicity, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0); 
339     PushBack( (TObject*)fDeDxVsP,          kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0);
340     //PushBack( (TObject*)fdNdEta,           kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0);
341 }
342
343 int AliHLTTPCTrackHistoComponent::ScanConfigurationArgument(int argc, const char** argv){
344 // see header file for class documentation
345  
346   if (argc<=0) return 0;
347   int i=0;
348   TString argument=argv[i];
349
350   // -event-modulo
351   if (argument.CompareTo("-event-modulo")==0) {
352     if (++i>=argc) return -EPROTO;
353     argument=argv[i];
354     fEvtMod=argument.Atoi();
355     return 2;
356   }    
357
358   // -buffer-size
359   if (argument.CompareTo("-buffer-size")==0) {
360     if (++i>=argc) return -EPROTO;
361     argument=argv[i];
362     fBufferSize=argument.Atoi();
363     return 2;
364   }    
365   
366   // -dEdx
367   if (argument.CompareTo("-dEdx")==0) {
368     fdEdx=kTRUE;
369     return 1;
370   }    
371   return -EINVAL;
372 }
373
374 void AliHLTTPCTrackHistoComponent::GetOCDBObjectDescription( TMap* const targetMap){
375 // Get a list of OCDB object description needed for the particular component
376   if (!targetMap) return;
377   targetMap->Add(new TObjString("HLT/ConfigTPC/TPCTrackHisto"), new TObjString("component arguments for setting the size of the filled ntuples and the event modulo for the mean multiplicity distribution"));
378 }