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