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