]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSClusterFinderComponent.cxx
dynamic adjustment of the output buffer size estimator
[u/mrichter/AliRoot.git] / HLT / ITS / AliHLTITSClusterFinderComponent.cxx
1 // $Id$
2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project        * 
4 //* ALICE Experiment at CERN, All rights reserved.                         *
5 //*                                                                        *
6 //* Primary Authors: Gaute Øvrebekk <st05886@alf.uib.no>                   *
7 //*                  for The ALICE HLT Project.                            *
8 //*                                                                        *
9 //* Permission to use, copy, modify and distribute this software and its   *
10 //* documentation strictly for non-commercial purposes is hereby granted   *
11 //* without fee, provided that the above copyright notice appears in all   *
12 //* copies and that both the copyright notice and this permission notice   *
13 //* appear in the supporting documentation. The authors make no claims     *
14 //* about the suitability of this software for any purpose. It is          *
15 //* provided "as is" without express or implied warranty.                  *
16 //**************************************************************************
17
18 /** @file   AliHLTITSClusterFinderComponent.cxx
19     @author Gaute Øvrebekk <st05886@alf.uib.no>
20     @date   
21     @brief  Component to run offline clusterfinders
22 */
23
24 #if __GNUC__>= 3
25 using namespace std;
26 #endif
27
28 #include "AliHLTITSClusterFinderComponent.h" 
29
30 #include "AliCDBEntry.h"
31 #include "AliCDBManager.h"
32 #include "AliHLTDataTypes.h"
33 #include "AliITSgeomTGeo.h"
34 #include "AliITSRecPoint.h"
35 #include "AliHLTITSSpacePointData.h"
36 #include "AliHLTITSClusterDataFormat.h"
37 #include <AliHLTDAQ.h>
38 #include "AliGeomManager.h"
39 #include "AliITSRecoParam.h"
40 #include "AliITSReconstructor.h"
41
42 #include <cstdlib>
43 #include <cerrno>
44 #include "TString.h"
45 #include "TObjString.h"
46 #include <sys/time.h>
47
48 /** ROOT macro for the implementation of ROOT specific class methods */
49 ClassImp(AliHLTITSClusterFinderComponent);
50
51 AliHLTITSClusterFinderComponent::AliHLTITSClusterFinderComponent(int mode)
52   :
53   fModeSwitch(mode),
54   fNModules(0),
55   fId(0),
56   fNddl(0),
57   fClusters(NULL),
58   fRawReader(NULL),
59   fDettype(NULL),
60   fgeom(NULL),
61   fgeomInit(NULL)
62
63   // see header file for class documentation
64   // or
65   // refer to README to build package
66   // or
67   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
68   if (fModeSwitch!=kClusterFinderSPD &&
69       fModeSwitch!=kClusterFinderSDD &&
70       fModeSwitch!=kClusterFinderSSD) {
71     HLTFatal("unknown cluster finder");
72   }
73 }
74
75 AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent() {
76   // see header file for class documentation
77 }
78
79 // Public functions to implement AliHLTComponent's interface.
80 // These functions are required for the registration process
81
82 const char* AliHLTITSClusterFinderComponent::GetComponentID()
83 {
84   // see header file for class documentation
85   switch(fModeSwitch){
86   case kClusterFinderSPD:
87     return "ITSClusterFinderSPD";
88     break;
89   case kClusterFinderSDD:        
90     return "ITSClusterFinderSDD";
91     break;
92   case kClusterFinderSSD:
93     return "ITSClusterFinderSSD";
94     break;
95   }
96   return "";
97 }
98
99 void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
100   // see header file for class documentation
101   list.clear(); 
102   switch(fModeSwitch){
103   case kClusterFinderSPD:
104     list.push_back( kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD );
105     break;
106   case kClusterFinderSDD:        
107     list.push_back( kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD );
108     break;
109   case kClusterFinderSSD:
110     list.push_back( kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD );
111     break;
112   }
113 }
114
115 AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType() {
116   // see header file for class documentation
117   switch(fModeSwitch){
118   case kClusterFinderSPD:
119     return kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
120     break;
121   case kClusterFinderSDD:        
122     return kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
123     break;
124   case kClusterFinderSSD:
125     return kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
126     break;
127   }
128   return kAliHLTVoidDataType;
129 }
130
131 void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
132   // see header file for class documentation
133   constBase = 0;
134   inputMultiplier = 100;
135 }
136
137 AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
138   // see header file for class documentation
139   return new AliHLTITSClusterFinderComponent(fModeSwitch);
140 }
141         
142 Int_t AliHLTITSClusterFinderComponent::DoInit( int /*argc*/, const char** /*argv*/ ) {
143   // see header file for class documentation
144   
145   if(fModeSwitch==kClusterFinderSPD) {
146     HLTDebug("using ClusterFinder for SPD");
147     //fNModules=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1) + AliITSgeomTGeo::GetNDetectors(2)*AliITSgeomTGeo::GetNLadders(2);
148     fId=AliHLTDAQ::DdlIDOffset("ITSSPD");
149     fNddl=AliHLTDAQ::NumberOfDdls("ITSSPD");
150   }
151   else if(fModeSwitch==kClusterFinderSDD) {
152     HLTDebug("using ClusterFinder for SDD");
153     //fNModules=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3) + AliITSgeomTGeo::GetNDetectors(4)*AliITSgeomTGeo::GetNLadders(4);
154     fId=AliHLTDAQ::DdlIDOffset("ITSSDD");
155     fNddl=AliHLTDAQ::NumberOfDdls("ITSSDD");
156   }
157   else if(fModeSwitch==kClusterFinderSSD) {
158     HLTDebug("using ClusterFinder for SSD");
159     //fNModules=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5) + AliITSgeomTGeo::GetNDetectors(6)*AliITSgeomTGeo::GetNLadders(6);
160     fId=AliHLTDAQ::DdlIDOffset("ITSSSD");
161     fNddl=AliHLTDAQ::NumberOfDdls("ITSSSD");
162   }
163   else{
164      HLTFatal("No mode set for clusterfindercomponent");
165   }
166
167   //Removed the warning for loading default RecoParam in HLT
168   AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
169   AliITSReconstructor *rec = new AliITSReconstructor();
170   rec->SetRecoParam(par);
171   
172   AliCDBManager* man = AliCDBManager::Instance();
173   if (!man->IsDefaultStorageSet()){
174     HLTError("Default CDB storage has not been set !");
175     return -ENOENT;
176   }
177
178   if(AliGeomManager::GetGeometry()==NULL){
179     AliGeomManager::LoadGeometry();
180   }
181   //fgeomInit = new AliITSInitGeometry(kvSPD02,2);
182   fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
183   //fgeomInit->InitAliITSgeom(fgeom);
184   fgeom = fgeomInit->CreateAliITSgeom();
185   
186   fNModules = fgeom->GetIndexMax();
187
188   fClusters = new TClonesArray*[fNModules]; 
189   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
190     fClusters[iModule] = NULL;
191   } 
192
193   //set dettype
194   fDettype = new AliITSDetTypeRec();
195   fDettype->SetITSgeom(fgeom);
196   fDettype->SetDefaultClusterFindersV2(kTRUE);
197   fDettype->SetDefaults();
198   
199   if ( fRawReader )
200     return -EINPROGRESS;
201
202   fRawReader = new AliRawReaderMemory();
203
204   return 0;
205 }
206
207 Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
208   // see header file for class documentation
209
210   if ( fRawReader )
211     delete fRawReader;
212   fRawReader = NULL;
213
214   if ( fDettype )
215     delete fDettype;
216   fDettype = NULL;
217
218   if ( fgeomInit )
219     delete fgeomInit;
220   fgeomInit = NULL;
221
222   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
223     if(fClusters[iModule] != NULL){
224       fClusters[iModule]->Delete();
225       delete fClusters[iModule];
226     }
227     fClusters[iModule] = NULL;
228   } 
229   
230   return 0;
231 }
232
233 Int_t AliHLTITSClusterFinderComponent::DoEvent( const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& /*trigData*/)
234 {  // see header file for class documentation
235
236   // -- Iterator over Data Blocks --
237   const AliHLTComponentBlockData* iter = NULL;
238   
239   if (!IsDataEvent()) return 0;
240
241   if ( evtData.fBlockCnt<=0 )
242       {
243         HLTDebug("no blocks in event" );
244         return 0;
245       }
246   AliHLTComponentDataType datatype;
247   if(fModeSwitch==kClusterFinderSPD) { datatype = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;}
248   else if(fModeSwitch==kClusterFinderSDD) {datatype = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;}
249   else if(fModeSwitch==kClusterFinderSSD) {datatype = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;}
250   
251   // -- Loop over blocks
252   for ( iter = GetFirstInputBlock(datatype); iter != NULL; iter = GetNextInputBlock() ) {
253   
254     // -- Debug output of datatype --
255     HLTDebug("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
256                evtData.fEventID, evtData.fEventID, 
257                DataType2Text(iter->fDataType).c_str(), 
258                DataType2Text(datatype).c_str());
259     
260     // -- Check for the correct data type
261     if ( iter->fDataType != (datatype) )  
262       continue;
263     
264     // -- Get equipment ID out of specification
265     AliHLTUInt32_t spec = iter->fSpecification;
266   
267     Int_t id = fId;
268     for ( Int_t ii = 0; ii < fNddl ; ii++ ) {   //number of ddl's
269       if ( spec & 0x00000001 ) {
270         id += ii;
271         break;
272       }
273       spec = spec >> 1 ;
274     }
275     
276     // -- Set equipment ID to the raw reader
277
278     if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
279       HLTWarning("Could not add buffer");
280     }
281     
282     if(fModeSwitch==kClusterFinderSPD) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SPD");}
283     else if(fModeSwitch==kClusterFinderSDD) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SDD");}
284     else if(fModeSwitch==kClusterFinderSSD) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SSD");}
285     
286     UInt_t nClusters=0;
287     for(int i=0;i<fNModules;i++){
288       if(fClusters[i] != NULL){
289         nClusters += fClusters[i]->GetEntries(); 
290       }
291     }
292     
293     UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
294     AliHLTUInt8_t *buffer = new AliHLTUInt8_t[bufferSize];
295     AliHLTITSClusterData *outputClusters = reinterpret_cast<AliHLTITSClusterData*>(buffer);
296     outputClusters->fSpacePointCnt=nClusters;
297     
298     int clustIdx=0;
299     for(int i=0;i<fNModules;i++){
300       if(fClusters[i] != NULL){
301         for(int j=0;j<fClusters[i]->GetEntries();j++){
302           AliITSRecPoint *recpoint = (AliITSRecPoint*) fClusters[i]->At(j);
303           outputClusters->fSpacePoints[clustIdx].fY=recpoint->GetY();
304           outputClusters->fSpacePoints[clustIdx].fZ=recpoint->GetZ();
305           outputClusters->fSpacePoints[clustIdx].fSigmaY2=recpoint->GetSigmaY2();
306           outputClusters->fSpacePoints[clustIdx].fSigmaZ2=recpoint->GetSigmaZ2();
307           outputClusters->fSpacePoints[clustIdx].fSigmaYZ=recpoint->GetSigmaYZ();
308           outputClusters->fSpacePoints[clustIdx].fQ=recpoint->GetQ();
309           outputClusters->fSpacePoints[clustIdx].fNy=recpoint->GetNy();
310           outputClusters->fSpacePoints[clustIdx].fNz=recpoint->GetNz();
311           outputClusters->fSpacePoints[clustIdx].fLayer=recpoint->GetLayer();
312           outputClusters->fSpacePoints[clustIdx].fIndex=recpoint->GetDetectorIndex() | recpoint->GetPindex() | recpoint->GetNindex();
313           outputClusters->fSpacePoints[clustIdx].fTracks[0]=recpoint->GetLabel(0);
314           outputClusters->fSpacePoints[clustIdx].fTracks[1]=recpoint->GetLabel(1);
315           outputClusters->fSpacePoints[clustIdx].fTracks[2]=recpoint->GetLabel(2);
316           
317           clustIdx++;
318         }
319       }
320     }
321     
322     if(fModeSwitch==kClusterFinderSPD) {PushBack(buffer,bufferSize,kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD,iter->fSpecification);}
323     else if(fModeSwitch==kClusterFinderSDD) {PushBack(buffer,bufferSize,kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD,iter->fSpecification);}
324     else if(fModeSwitch==kClusterFinderSSD) {PushBack(buffer,bufferSize,kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD,iter->fSpecification);}
325     
326     delete buffer; 
327     fRawReader->ClearBuffers();
328     
329     for (Int_t iModule = 0; iModule < fNModules; iModule++) {
330       if(fClusters[iModule] != NULL){
331         fClusters[iModule]->Delete();
332         delete fClusters[iModule];
333       }
334       fClusters[iModule] = NULL;
335     }  
336         
337   } //  for ( ndx = 0; ndx < evtData.fBlockCnt; ndx++ ) {    
338   
339   return 0;
340 }
341
342 int AliHLTITSClusterFinderComponent::Configure(const char* arguments)
343 {
344   
345   int iResult=0;
346   
347   if (!arguments) return iResult;
348   
349   TString allArgs=arguments;
350   TString argument;
351   
352   TObjArray* pTokens=allArgs.Tokenize(" ");
353   
354   if (pTokens) {
355     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
356       argument=((TObjString*)pTokens->At(i))->GetString();
357       if (argument.IsNull()) continue;
358       /*
359       if (argument.CompareTo("")==0) {
360         HLTInfo("");
361         continue;
362       }
363
364       else if (argument.CompareTo("")==0) {
365         HLTInfo("");
366         continue;
367       }
368       */
369       else {
370         HLTError("unknown argument %s", argument.Data());
371         iResult=-EINVAL;
372         break;
373       }
374     }
375     delete pTokens;
376   }
377   
378   return iResult;
379 }
380
381 int AliHLTITSClusterFinderComponent::Reconfigure(const char* cdbEntry, const char* chainId)
382 {
383   // see header file for class documentation
384   int iResult=0;
385   
386   const char* path="HLT/ConfigITS/ClusterFinderComponent";
387   const char* defaultNotify="";
388   if (cdbEntry) {
389     path=cdbEntry;
390     defaultNotify=" (default)";
391   }
392   if (path) {
393     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
394     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
395     if (pEntry) {
396       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
397       if (pString) {
398         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
399         iResult=Configure(pString->GetString().Data());
400       } else {
401         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
402       }
403     } else {
404       HLTError("can not fetch object \"%s\" from CDB", path);
405     }
406   }
407   
408   return iResult;
409 }
410