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