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