]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSClusterFinderComponent.cxx
bug fix https://savannah.cern.ch/bugs/?79511
[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 Ovrebekk <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 Ovrebekk <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 "AliCDBStorage.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 #include "AliHLTITSClusterFinderSPD.h"
42 #include "AliHLTITSClusterFinderSSD.h"
43 #include "TMap.h"
44 #include "AliITSRecPointContainer.h"
45 #include "AliRunLoader.h"
46 #include "AliLoader.h"
47
48 #include <cstdlib>
49 #include <cerrno>
50 #include "TFile.h"
51 #include "TString.h"
52 #include "TObjString.h"
53 #include <sys/time.h>
54
55 /** ROOT macro for the implementation of ROOT specific class methods */
56 ClassImp(AliHLTITSClusterFinderComponent);
57
58 AliHLTITSClusterFinderComponent::AliHLTITSClusterFinderComponent(int mode)
59   :
60   fModeSwitch(mode),
61   fInputDataType(kAliHLTVoidDataType),
62   fOutputDataType(kAliHLTVoidDataType),
63   fUseOfflineFinder(0),
64   fNModules(0),
65   fId(0),
66   fNddl(0),
67   fRawReader(NULL),
68   fDettype(NULL),
69   fgeom(NULL),
70   fgeomInit(NULL),
71   fSPD(NULL),
72   fSSD(NULL),
73   tD(NULL),
74   tR(NULL),
75   fSPDNModules(0),
76   fSDDNModules(0),
77   fSSDNModules(0),
78   fFirstModule(0),
79   fLastModule(0),
80   fnClusters(0),
81   fclusters(),
82   fBenchmark(GetComponentID()),
83   fOutputSizeOffset(0),
84   fInputMultiplierDigits(20),
85   fpLoader(NULL)
86
87   // see header file for class documentation
88   // or
89   // refer to README to build package
90   // or
91   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
92
93   switch(fModeSwitch){
94   case kClusterFinderSPD:
95     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;
96     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
97     break;
98   case kClusterFinderSDD:        
99     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;
100     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
101     break;
102   case kClusterFinderSSD:
103     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;
104     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
105     break;
106   case kClusterFinderDigits:
107     fInputDataType  = kAliHLTDataTypeAliTreeD|kAliHLTDataOriginITS;
108     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITS;
109     break;
110   default:
111     HLTFatal("unknown cluster finder");
112   }
113 }
114
115 AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent() 
116 {
117   // see header file for class documentation
118   delete fRawReader;
119   delete fDettype;
120   delete fgeomInit;  
121   delete fSPD;
122   delete fSSD;
123 }
124
125 // Public functions to implement AliHLTComponent's interface.
126 // These functions are required for the registration process
127
128 const char* AliHLTITSClusterFinderComponent::GetComponentID()
129 {
130   // see header file for class documentation
131   switch(fModeSwitch){
132   case kClusterFinderSPD:
133     return "ITSClusterFinderSPD";
134     break;
135   case kClusterFinderSDD:        
136     return "ITSClusterFinderSDD";
137     break;
138   case kClusterFinderSSD:
139     return "ITSClusterFinderSSD";
140     break;
141   case kClusterFinderDigits:
142     return "ITSClusterFinderDigits";
143     break;
144   }
145   return "";
146 }
147
148 void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) 
149 {
150   // see header file for class documentation
151   list.clear(); 
152   list.push_back( fInputDataType );
153 }
154
155 AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType() 
156 {
157   // see header file for class documentation
158   return fOutputDataType;
159 }
160
161 void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
162   // see header file for class documentation
163   constBase = fOutputSizeOffset;
164   switch(fModeSwitch){
165   case kClusterFinderDigits:
166     inputMultiplier = fInputMultiplierDigits;
167     break;
168   case kClusterFinderSPD:
169   case kClusterFinderSDD:        
170   case kClusterFinderSSD:
171   default:
172     inputMultiplier = 20;
173   }
174 }
175
176 AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
177   // see header file for class documentation
178   return new AliHLTITSClusterFinderComponent(fModeSwitch);
179 }
180         
181 Int_t AliHLTITSClusterFinderComponent::DoInit( int argc, const char** argv ) {
182   // see header file for class documentation
183   fBenchmark.Reset();
184   fBenchmark.SetTimer(0,"total");
185   fBenchmark.SetTimer(1,"reco");
186
187   /*
188   fStatTime = 0;
189   fStatTimeAll = 0;
190   fStatTimeC = 0;
191   fStatTimeAllC = 0;
192   fStatNEv = 0;
193   */
194   
195   Int_t runNo = GetRunNo();
196   AliCDBStorage* store = AliCDBManager::Instance()->GetDefaultStorage();
197   if (!store) {
198     return -ENOENT;
199   }
200
201   bool cdbOK = true;
202   //OCDB for SPD
203   if(store->GetLatestVersion("ITS/Calib/SPDNoisy", runNo)<0){
204     HLTError("SPDNoisy is not found in SPD/Calib");
205     cdbOK = false;
206   }
207   if(store->GetLatestVersion("ITS/Calib/SPDDead", runNo)<0){
208     HLTError("SPDDead is not found in SPD/Calib");
209     cdbOK = false;
210   }
211   if(store->GetLatestVersion("TRIGGER/SPD/PITConditions", runNo)<0){
212     HLTError("PITConditions is not found in TRIGGER/SPD");
213     cdbOK = false;
214   }
215   
216   //OCDB for SDD
217   if(store->GetLatestVersion("ITS/Calib/CalibSDD", runNo)<0){
218     HLTError("CalibSDD is not found in ITS/Calib");
219     cdbOK = false;
220   }
221   if(store->GetLatestVersion("ITS/Calib/RespSDD", runNo)<0){
222     HLTError("RespSDD is not found in ITS/Calib");
223     cdbOK = false;
224   }
225   if(store->GetLatestVersion("ITS/Calib/DriftSpeedSDD", runNo)<0){
226     HLTError("DriftSpeedSDD is not found in ITS/Calib");
227     cdbOK = false;
228   }
229   if(store->GetLatestVersion("ITS/Calib/DDLMapSDD", runNo)<0){
230     HLTError("DDLMapSDD is not found in ITS/Calib");
231     cdbOK = false;
232   }
233   if(store->GetLatestVersion("ITS/Calib/MapsTimeSDD", runNo)<0){
234     HLTError("MapsTimeSDD is not found in ITS/Calib");
235     cdbOK = false;
236   }
237
238   //OCDB for SSD
239   if(store->GetLatestVersion("ITS/Calib/NoiseSSD", runNo)<0){
240     HLTError("NoiseSSD is not found in ITS/Calib");
241     cdbOK = false;
242   }
243   if(store->GetLatestVersion("ITS/Calib/GainSSD", runNo)<0){
244     HLTError("GainSSD is not found in ITS/Calib");
245     cdbOK = false;
246   }
247   if(store->GetLatestVersion("ITS/Calib/BadChannelsSSD", runNo)<0){
248     HLTError("BadChannelsSSD is not found in ITS/Calib");
249     cdbOK = false;
250   }
251   
252   //General reconstruction
253   if(store->GetLatestVersion("GRP/CTP/Scalers", runNo)<0){
254     HLTError("Scalers is not found in GRP/CTP/");
255     cdbOK = false;
256   }
257   if(!cdbOK){return -EIO;}
258
259   if(fModeSwitch==kClusterFinderSPD) {
260     HLTDebug("using ClusterFinder for SPD");
261     //fNModules=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1) + AliITSgeomTGeo::GetNDetectors(2)*AliITSgeomTGeo::GetNLadders(2);
262     fId=AliHLTDAQ::DdlIDOffset("ITSSPD");
263     fNddl=AliHLTDAQ::NumberOfDdls("ITSSPD");
264   }
265   else if(fModeSwitch==kClusterFinderSDD) {
266     HLTDebug("using ClusterFinder for SDD");
267     //fNModules=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3) + AliITSgeomTGeo::GetNDetectors(4)*AliITSgeomTGeo::GetNLadders(4);
268     fId=AliHLTDAQ::DdlIDOffset("ITSSDD");
269     fNddl=AliHLTDAQ::NumberOfDdls("ITSSDD");
270   }
271   else if(fModeSwitch==kClusterFinderSSD) {
272     HLTDebug("using ClusterFinder for SSD");
273     //fNModules=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5) + AliITSgeomTGeo::GetNDetectors(6)*AliITSgeomTGeo::GetNLadders(6);
274     fId=AliHLTDAQ::DdlIDOffset("ITSSSD");
275     fNddl=AliHLTDAQ::NumberOfDdls("ITSSSD");
276   }
277   else if(fModeSwitch==kClusterFinderDigits) {
278     tR = new TTree();
279   }
280   else{
281      HLTFatal("No mode set for clusterfindercomponent");
282   }
283
284   //Removed the warning for loading default RecoParam in HLT
285   AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
286   AliITSReconstructor *rec = new AliITSReconstructor();
287   rec->SetRecoParam(par);
288   
289   AliCDBManager* man = AliCDBManager::Instance();
290   if (!man->IsDefaultStorageSet()){
291     HLTError("Default CDB storage has not been set !");
292     return -ENOENT;
293   }
294
295   if(AliGeomManager::GetGeometry()==NULL){
296     AliGeomManager::LoadGeometry();
297   }
298
299   fgeomInit = new AliITSInitGeometry();
300   //fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
301   //fgeomInit->InitAliITSgeom(fgeom);
302   fgeom = fgeomInit->CreateAliITSgeom();
303  
304   fNModules = fgeom->GetIndexMax();
305   Int_t modperlay[6];
306   for(Int_t i=0;i<6;i++)modperlay[i]=AliITSgeomTGeo::GetNDetectors(1+i)*AliITSgeomTGeo::GetNLadders(1+i);
307   fSPDNModules=modperlay[0]+modperlay[1];
308   fSDDNModules=modperlay[2]+modperlay[3];
309   fSSDNModules=modperlay[4]+modperlay[5];
310   
311   if(fModeSwitch==kClusterFinderSPD) {
312     fFirstModule=0;
313     fLastModule=fSPDNModules;
314   }
315   else if(fModeSwitch==kClusterFinderSDD) {
316      fFirstModule=fSPDNModules;
317      fLastModule=fFirstModule + fSDDNModules;
318   }
319   else if(fModeSwitch==kClusterFinderSSD) {
320     fFirstModule=fSPDNModules + fSDDNModules;
321     fLastModule=fFirstModule + fSSDNModules;
322   }
323  
324   //set dettype
325   fDettype = new AliITSDetTypeRec();
326   fDettype->SetITSgeom(fgeom); 
327   fDettype->SetDefaults();
328   fDettype->SetDefaultClusterFindersV2(kTRUE); 
329
330   if ( fRawReader )
331     return -EINPROGRESS;
332
333   fRawReader = new AliRawReaderMemory();
334   fSPD = new AliHLTITSClusterFinderSPD( fDettype );
335   fSSD = new AliHLTITSClusterFinderSSD( fDettype, fRawReader );
336
337   TString arguments = "";
338   for ( int i = 0; i < argc; i++ ) {
339     if ( !arguments.IsNull() ) arguments += " ";
340     arguments += argv[i];
341   }
342
343   tD = NULL;
344   //tR = NULL;
345
346   return Configure( arguments.Data() );
347 }
348
349 Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
350   // see header file for class documentation
351
352   if ( fRawReader )
353     delete fRawReader;
354   fRawReader = NULL;
355
356   if ( fDettype )
357     delete fDettype;
358   fDettype = NULL;
359
360   if ( fgeomInit )
361     delete fgeomInit;
362   fgeomInit = NULL;
363   
364   delete fSPD;
365   fSPD = 0;
366
367   delete fSSD;
368   fSSD = 0;
369
370   fUseOfflineFinder = 0;
371
372   if (fpLoader) {
373     fpLoader->UnloadDigits();
374   }
375   fpLoader=NULL;
376
377   return 0;
378 }
379
380 int AliHLTITSClusterFinderComponent::DoEvent
381 (
382  const AliHLTComponentEventData& evtData,
383  const AliHLTComponentBlockData* /*blocks*/,
384  AliHLTComponentTriggerData& /*trigData*/,
385  AliHLTUInt8_t* outputPtr,
386  AliHLTUInt32_t& size,
387   vector<AliHLTComponentBlockData>& outputBlocks )
388 {
389  // see header file for class documentation
390   
391   AliHLTUInt32_t maxBufferSize = size;
392   size = 0; // output size
393
394   if (!IsDataEvent()) return 0;
395   
396   if ( evtData.fBlockCnt<=0 )
397     {
398       HLTDebug("no blocks in event" );
399       return 0;
400     }
401
402   AliHLTUInt32_t totalInputSize=0;
403   fBenchmark.StartNewEvent();
404   fBenchmark.Start(0);
405   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(fInputDataType); i!=NULL; i=GetNextInputBlock() ){
406     totalInputSize+=i->fSize;
407     fBenchmark.AddInput(i->fSize);
408   }
409
410   Int_t ret = 0;
411
412   if(fModeSwitch==kClusterFinderDigits) {
413
414     for ( const TObject *iter = GetFirstInputObject(fInputDataType); iter != NULL; iter = GetNextInputObject() ) {  
415       tD = dynamic_cast<TTree*>(const_cast<TObject*>( iter ) );
416       if(!tD){
417         HLTFatal("No Digit Tree found");
418         return -1;
419       }
420       // 2010-04-17 very crude workaround: TTree objects are difficult to send
421       // The actual case: Running ITS and TPC reconstruction fails at the second event
422       // to read the ITS digits from the TreeD
423       //
424       // Reason: reading fails in TBranch::GetBasket, there a new basket is opened from
425       // a TFile object. The function TBranch::GetFile returns the file object from
426       // an internal fDirectory (TDirectory) object. This file is at the second event
427       // set to the TPC.Digits.root. The internal mismatch creates a seg fault
428       //
429       // Investigation: TBranch::Streamer uses a crude assignment after creating the
430       // TBranch object
431       //    fDirectory = gDirectory;
432       // gDirectory is obviously not set correctly. Setting the directory to a TFile
433       // object for the ITS digits helps to fix the internal mess. Tried also to set
434       // the Directory for the TreeD to NULL (This has only effect if ones sets it 
435       // to something not NULL first, and then to NULL). But then no content, i.e.
436       // ITS clusters could be retrieved.
437       //
438       // Conclusion: TTree objects are hardly to be sent via TMessage, there are direct
439       // links to the file required anyhow.
440       //
441       // 2011-01-28 hotfix reloaded: accessing the files like that fails if there are
442       // multiple digit files because of a large number of events. New ugly fix is to
443       // use the global runloader instance to get hold on the digits tree.
444       fnClusters = 0;
445       AliRunLoader* pRunLoader=AliRunLoader::Instance();
446       if (!pRunLoader) {
447         HLTError("failed to get global runloader instance");
448         return -ENOSYS;
449       }
450       // get the specific loader for the module
451       if (!fpLoader) {
452         const char* loaderType="ITSLoader";
453         fpLoader=pRunLoader->GetLoader(loaderType);
454         if (!fpLoader) {
455           HLTError("can not get loader \"%s\" from runloader", loaderType);
456           return -ENOSYS;
457         }
458         // prepare the loader
459         fpLoader->LoadDigits("read");
460       }
461       pRunLoader->GetEvent(GetEventCount());
462
463       tD=fpLoader->TreeD();
464       tR->Reset();
465       tR->SetDirectory(0);
466       fDettype->SetTreeAddressD(tD);
467       fDettype->MakeBranch(tR,"R");
468       fDettype->SetTreeAddressR(tR);
469       Option_t *opt="All";
470       fBenchmark.Start(1);
471       fDettype->DigitsToRecPoints(tD,tR,0,opt,0);
472       fBenchmark.Stop(1);
473       TClonesArray * fRecPoints = NULL;
474       tR->SetBranchAddress("ITSRecPoints",&fRecPoints);
475       for(Int_t treeEntry=0;treeEntry<tR->GetEntries();treeEntry++){
476         tR->GetEntry(treeEntry);
477         fnClusters += fRecPoints->GetEntries();
478       }
479
480       UInt_t bufferSize = fnClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
481       if( size + bufferSize > maxBufferSize ){
482         //HLTWarning( "Output buffer size exceed (buffer size %d, required size %d)", maxBufferSize, size+bufferSize);
483         if (totalInputSize>0) {
484           fInputMultiplierDigits=(float)(size+bufferSize)/totalInputSize;
485           fInputMultiplierDigits+=1.;
486         } else {
487           fOutputSizeOffset=totalInputSize;
488           fInputMultiplierDigits=1.;
489         }
490         ret = -ENOSPC;      
491         break;          
492       }
493       if( fnClusters>0 ){
494         fBenchmark.Start(1);
495         RecPointToSpacePoint(outputPtr,size);
496         fBenchmark.Stop(1);
497         AliHLTComponentBlockData bd;
498         FillBlockData( bd );
499         bd.fOffset = size;
500         bd.fSize = bufferSize;
501         bd.fSpecification = 0x00000000;
502         bd.fDataType = GetOutputDataType();
503         outputBlocks.push_back( bd );
504         size += bufferSize;
505         fBenchmark.AddOutput(bd.fSize);
506       }
507     }
508   }
509   else{
510
511     AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
512
513     // -- Loop over blocks
514     for( const AliHLTComponentBlockData* iter = GetFirstInputBlock(fInputDataType); iter != NULL; iter = GetNextInputBlock() ) {
515
516       if(fUseOfflineFinder){
517         if(fModeSwitch==kClusterFinderSPD){rpc->ResetSPD();}
518         if(fModeSwitch==kClusterFinderSSD){rpc->ResetSSD();}
519       }
520       if(fModeSwitch==kClusterFinderSDD){rpc->ResetSDD();}
521       
522       // -- Debug output of datatype --
523       HLTDebug("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
524                evtData.fEventID, evtData.fEventID, 
525                DataType2Text(iter->fDataType).c_str(), 
526                DataType2Text(fInputDataType).c_str());
527       
528       // -- Check for the correct data type
529       if ( iter->fDataType != (fInputDataType) )  
530         continue;
531       
532       // -- Get equipment ID out of specification
533       AliHLTUInt32_t spec = iter->fSpecification;
534       
535       Int_t id = fId;
536       for ( Int_t ii = 0; ii < fNddl ; ii++ ) {   //number of ddl's
537         if ( spec & 0x00000001 ) {
538           id += ii;
539           break;
540         }
541         spec = spec >> 1 ;
542       }
543       
544       // -- Set equipment ID to the raw reader
545       
546       if(!fRawReader){
547         HLTWarning("The fRawReader pointer is NULL");
548         continue;
549       }
550
551       if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
552         HLTWarning("Could not add buffer");
553       }
554
555       fBenchmark.Start(1);
556
557       fnClusters = 0;
558
559       if(fModeSwitch==kClusterFinderSPD && !fUseOfflineFinder){ fSPD->RawdataToClusters( fRawReader, fclusters ); }
560       else if(fModeSwitch==kClusterFinderSSD && !fUseOfflineFinder){ fSSD->RawdataToClusters( fclusters ); }
561       else{
562         if(fModeSwitch==kClusterFinderSPD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,"SPD");}
563         if(fModeSwitch==kClusterFinderSSD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,"SSD");}
564         if(fModeSwitch==kClusterFinderSDD) {fDettype->DigitsToRecPoints(fRawReader,"SDD");}
565         TClonesArray* clusters = NULL;
566         for(int i=fFirstModule;i<fLastModule;i++){
567           clusters = rpc->UncheckedGetClusters(i);
568           if(clusters != NULL){
569             fnClusters += clusters->GetEntriesFast();
570           }
571         }     
572       }
573       
574       /* Do not work
575          if(fModeSwitch==kClusterfinderSPD){
576          fnClusters = rpc->GetNClustersInLayerFast(1) + rpc->GetNClustersInLayerFast(2);
577          }
578          if(fModeSwitch==kClusterfinderSDD){
579          fnClusters = rpc->GetNClustersInLayerFast(3) + rpc->GetNClustersInLayerFast(4);
580          }
581          if(fModeSwitch==kClusterfinderSSD){
582          fnClusters = rpc->GetNClustersInLayerFast(5) + rpc->GetNClustersInLayerFast(6);
583          }
584       */
585
586       fBenchmark.Stop(1);
587       
588       fRawReader->ClearBuffers();    
589       
590       UInt_t nClusters=fclusters.size();
591       if(nClusters>0){fnClusters = nClusters;}
592
593       UInt_t bufferSize = fnClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
594       if( size + bufferSize > maxBufferSize ){
595         HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
596         ret = -ENOSPC;      
597         break;          
598       }
599       if( fnClusters>0 ){
600
601         RecPointToSpacePoint(outputPtr,size);
602
603         AliHLTComponentBlockData bd;
604         FillBlockData( bd );
605         bd.fOffset = size;
606         bd.fSize = bufferSize;
607         bd.fSpecification = iter->fSpecification;
608         bd.fDataType = GetOutputDataType();
609         outputBlocks.push_back( bd );
610         size += bufferSize;
611         fBenchmark.AddOutput(bd.fSize);
612         if(nClusters>0){fclusters.clear();}     
613       }
614       
615     } // input blocks
616   }
617
618   fBenchmark.Stop(0);
619   HLTInfo(fBenchmark.GetStatistics());
620   return ret;
621 }
622
623 int AliHLTITSClusterFinderComponent::Configure(const char* arguments)
624 {
625   
626   int iResult=0;
627   
628   if (!arguments) return iResult;
629   
630   TString allArgs=arguments;
631   TString argument;
632   
633   TObjArray* pTokens=allArgs.Tokenize(" ");
634   
635   if (pTokens) {
636     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
637       argument=((TObjString*)pTokens->At(i))->GetString();
638       if (argument.IsNull()) continue;      
639       if (argument.CompareTo("-use-offline-finder")==0) {
640         fUseOfflineFinder = 1;
641         HLTInfo("Off-line ClusterFinder will be used");
642         continue;
643       }
644       /*
645       else if (argument.CompareTo("")==0) {
646         HLTInfo("");
647         continue;
648       }
649       */
650       else {
651         HLTError("unknown argument %s", argument.Data());
652         iResult=-EINVAL;
653         break;
654       }
655     }
656     delete pTokens;
657   }
658   
659   return iResult;
660 }
661
662 int AliHLTITSClusterFinderComponent::Reconfigure(const char* cdbEntry, const char* chainId)
663 {
664   // see header file for class documentation
665   int iResult=0;
666   
667   const char* path="";
668
669   switch(fModeSwitch){
670   case kClusterFinderSPD:
671     path = "HLT/ConfigITS/ITSClusterFinderSPD";
672     break;
673   case kClusterFinderSDD:        
674     path = "HLT/ConfigITS/ITSClusterFinderSDD";
675     break;
676   case kClusterFinderSSD:
677     path = "HLT/ConfigITS/ITSClusterFinderSSD";
678     break;
679   case kClusterFinderDigits:
680     path = "";
681     break;
682   default:
683     HLTFatal("unknown cluster finder");
684   }
685
686   const char* defaultNotify="";
687   if (cdbEntry) {
688     path=cdbEntry;
689     defaultNotify=" (default)";
690   }
691   if (path) {
692     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
693     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
694     if (pEntry) {
695       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
696       if (pString) {
697         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
698         iResult=Configure(pString->GetString().Data());
699       } else {
700         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
701       }
702     } else {
703       HLTError("can not fetch object \"%s\" from CDB", path);
704     }
705   }
706   
707   return iResult;
708 }
709 void AliHLTITSClusterFinderComponent::GetOCDBObjectDescription( TMap* const targetMap)
710 {
711   // Get a list of OCDB object description.
712   if (!targetMap) return;
713   //SPD
714   targetMap->Add(new TObjString("ITS/Calib/SPDNoisy"),new TObjString("Calibration object for SPD" ));
715   targetMap->Add(new TObjString("ITS/Calib/SPDDead"),new TObjString("Calibration object for SPD" ));
716   targetMap->Add(new TObjString("TRIGGER/SPD/PITConditions"),new TObjString("Calibration object for SPD" ));
717   //SDD
718   targetMap->Add(new TObjString("ITS/Calib/CalibSDD"),new TObjString("Calibration object for SDD" ));
719   targetMap->Add(new TObjString("ITS/Calib/RespSDD"),new TObjString("Calibration object for SDD" ));
720   targetMap->Add(new TObjString("ITS/Calib/DriftSpeedSDD"),new TObjString("Calibration object for SDD" ));
721   targetMap->Add(new TObjString("ITS/Calib/DDLMapSDD"),new TObjString("Calibration object for SDD" ));
722   targetMap->Add(new TObjString("ITS/Calib/MapsTimeSDD"),new TObjString("Calibration object for SDD" ));
723   //SSD
724   targetMap->Add(new TObjString("ITS/Calib/NoiseSSD"),new TObjString("Calibration object for SSD" ));
725   targetMap->Add(new TObjString("ITS/Calib/GainSSD"),new TObjString("Calibration object for SSD" ));
726   targetMap->Add(new TObjString("ITS/Calib/BadChannelsSSD"),new TObjString("Calibration object for SSD" ));
727   //General reconstruction
728   targetMap->Add(new TObjString("GRP/CTP/Scalers"),new TObjString("General reconstruction object" ));
729 }
730
731
732 void AliHLTITSClusterFinderComponent::RecPointToSpacePoint(AliHLTUInt8_t* outputPtr,AliHLTUInt32_t& size){
733   AliHLTITSClusterData *outputClusters = reinterpret_cast<AliHLTITSClusterData*>(outputPtr + size);
734   outputClusters->fSpacePointCnt=fnClusters;    
735   int clustIdx=0;
736   if(fModeSwitch==kClusterFinderDigits) {
737     TClonesArray * fRecPoints = NULL;
738     tR->SetBranchAddress("ITSRecPoints",&fRecPoints);
739     for(Int_t treeEntry=0;treeEntry<tR->GetEntries();treeEntry++){
740       tR->GetEntry(treeEntry);
741       fnClusters += fRecPoints->GetEntries();
742       for(Int_t tCloneEntry=0;tCloneEntry<fRecPoints->GetEntries();tCloneEntry++){
743         AliITSRecPoint *recpoint=(AliITSRecPoint*)fRecPoints->At(tCloneEntry);
744         RecpointToOutput(outputClusters,recpoint,clustIdx);
745       }
746     } 
747   }
748   else if(fclusters.size()>0){
749     for(UInt_t i=0;i<fclusters.size();i++){
750       AliITSRecPoint *recpoint = (AliITSRecPoint*) &(fclusters[i]);
751       RecpointToOutput(outputClusters,recpoint,clustIdx);
752     }
753   }
754   else{
755     AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
756     TClonesArray* clusters = NULL;
757     for(Int_t i=fFirstModule;i<fLastModule;i++){
758       clusters = rpc->UncheckedGetClusters(i);
759       for(Int_t j=0;j<clusters->GetEntriesFast();j++){
760         AliITSRecPoint *recpoint = (AliITSRecPoint*) clusters->At(j);
761         RecpointToOutput(outputClusters,recpoint,clustIdx);
762       }
763     } 
764   }
765 }
766
767 void AliHLTITSClusterFinderComponent::RecpointToOutput(AliHLTITSClusterData *outputClusters, AliITSRecPoint *recpoint, int &clustIdx){
768   outputClusters->fSpacePoints[clustIdx].fY=recpoint->GetY();
769   outputClusters->fSpacePoints[clustIdx].fZ=recpoint->GetZ();
770   outputClusters->fSpacePoints[clustIdx].fSigmaY2=recpoint->GetSigmaY2();
771   outputClusters->fSpacePoints[clustIdx].fSigmaZ2=recpoint->GetSigmaZ2();
772   outputClusters->fSpacePoints[clustIdx].fSigmaYZ=recpoint->GetSigmaYZ();
773   outputClusters->fSpacePoints[clustIdx].fQ=recpoint->GetQ();
774   outputClusters->fSpacePoints[clustIdx].fNy=recpoint->GetNy();
775   outputClusters->fSpacePoints[clustIdx].fNz=recpoint->GetNz();
776   outputClusters->fSpacePoints[clustIdx].fLayer=recpoint->GetLayer();
777   outputClusters->fSpacePoints[clustIdx].fIndex=recpoint->GetDetectorIndex() | recpoint->GetPindex() | recpoint->GetNindex();
778   outputClusters->fSpacePoints[clustIdx].fTracks[0]=recpoint->GetLabel(0);
779   outputClusters->fSpacePoints[clustIdx].fTracks[1]=recpoint->GetLabel(1);
780   outputClusters->fSpacePoints[clustIdx].fTracks[2]=recpoint->GetLabel(2);      
781   clustIdx++;
782 }