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