]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
- coverity fixes for unused but set variables in various places
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDClusterizerComponent.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors:                                                       *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /** @file   AliHLTTRDClusterizerComponent.cxx
20     @author Theodor Rascanu
21     @date   
22     @brief  A TRDClusterizer processing component for the HLT. 
23 */
24
25 // see header file for class documentation                                   //
26 // or                                                                        //
27 // refer to README to build package                                          //
28 // or                                                                        //
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt                          //
30
31 #include "TTree.h"
32 #include "TFile.h"
33 #include "TBranch.h"
34
35 #include "AliHLTTRDClusterizerComponent.h"
36 #include "AliHLTTRDDefinitions.h"
37 #include "AliHLTTRDClusterizer.h"
38 #include "AliHLTTRDUtils.h"
39
40 #include "AliGeomManager.h"
41 #include "AliTRDReconstructor.h"
42 #include "AliCDBManager.h"
43 #include "AliCDBStorage.h"
44 #include "AliCDBEntry.h"
45 #include "AliTRDrecoParam.h"
46 #include "AliTRDrawStream.h"
47 #include "AliTRDcluster.h"
48
49 #include "AliRawReaderMemory.h"
50
51 #ifdef HAVE_VALGRIND_CALLGRIND_H
52 #include <valgrind/callgrind.h>
53 #else
54 #define CALLGRIND_START_INSTRUMENTATION (void)0
55 #define CALLGRIND_STOP_INSTRUMENTATION (void)0
56 #endif
57
58 #include <cstdlib>
59 #include <cerrno>
60 #include <string>
61
62 using namespace std;
63
64 ClassImp(AliHLTTRDClusterizerComponent)
65    
66 AliHLTTRDClusterizerComponent::AliHLTTRDClusterizerComponent()
67 : AliHLTProcessor(),
68   fOutputPercentage(100),
69   fOutputConst(0),
70   fClusterizer(NULL),
71   fRecoParam(NULL),
72   fMemReader(NULL),
73   fReconstructor(NULL),
74   fRecoParamType(-1),
75   fRecoDataType(-1),
76   fRawDataVersion(2),
77   fyPosMethod(1),
78   fgeometryFileName(""),
79   fProcessTracklets(kFALSE),
80   fHLTstreamer(kTRUE),
81   fTC(kFALSE),
82   fHLTflag(kTRUE),
83   fHighLevelOutput(kFALSE),
84   fEmulateHLTClusters(kFALSE)
85 {
86   // Default constructor
87
88 }
89
90 AliHLTTRDClusterizerComponent::~AliHLTTRDClusterizerComponent()
91 {
92   // Destructor
93   // Work is Done in DoDeInit()
94 }
95
96
97 const char* AliHLTTRDClusterizerComponent::GetComponentID()
98 {
99   // Return the component ID const char *
100   return "TRDClusterizer"; // The ID of this component
101 }
102
103 void AliHLTTRDClusterizerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
104 {
105   // Get the list of input data
106   list.clear(); // We do not have any requirements for our input data type(s).
107   list.push_back(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
108 }
109
110 AliHLTComponentDataType AliHLTTRDClusterizerComponent::GetOutputDataType()
111 {
112   // Get the output data type
113   return kAliHLTMultipleDataType;
114 }
115
116 int AliHLTTRDClusterizerComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
117 {
118   // Get the output data type
119   tgtList.clear();
120   tgtList.push_back(AliHLTTRDDefinitions::fgkClusterDataType);
121   tgtList.push_back(AliHLTTRDDefinitions::fgkMCMtrackletDataType);
122   return tgtList.size();
123 }
124
125
126 void AliHLTTRDClusterizerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
127 {
128   // Get the output data size
129   constBase = fOutputConst;
130   inputMultiplier = ((double)fOutputPercentage)*4/100.0;
131 }
132
133 AliHLTComponent* AliHLTTRDClusterizerComponent::Spawn()
134 {
135   // Spawn function, return new instance of this class
136   return new AliHLTTRDClusterizerComponent;
137 };
138
139 int AliHLTTRDClusterizerComponent::DoInit( int argc, const char** argv )
140 {
141   // perform initialization. We check whether our relative output size is specified in the arguments.
142   int iResult=0;
143   
144   fReconstructor = new AliTRDReconstructor();
145   HLTDebug("TRDReconstructor at 0x%x", fReconstructor);
146
147   TString configuration="";
148   TString argument="";
149   for (int i=0; i<argc && iResult>=0; i++) {
150     argument=argv[i];
151     if (!configuration.IsNull()) configuration+=" ";
152     configuration+=argument;
153   }
154
155   if (!configuration.IsNull()) {
156     iResult=Configure(configuration.Data());
157   } else {
158     iResult=Reconfigure(NULL, NULL);
159   }
160
161   if(!fClusterizer){
162     HLTFatal("Clusterizer was not initialized!");
163     return -1;
164   }
165
166   if(iResult<0) return iResult;
167
168   fMemReader = new AliRawReaderMemory;
169   fClusterizer->SetReconstructor(fReconstructor);
170   fClusterizer->SetUseLabels(kFALSE);
171
172   if(fReconstructor->IsProcessingTracklets())
173     fOutputConst = fClusterizer->GetTrMemBlockSize();
174
175   return iResult;
176 }
177
178 int AliHLTTRDClusterizerComponent::DoDeinit()
179 {
180   // Deinitialization of the component
181   delete fMemReader;
182   fMemReader = 0;
183   delete fClusterizer;
184   fClusterizer = 0;
185   
186   //fReconstructor->SetClusters(0x0);
187   delete fReconstructor;
188   fReconstructor = 0x0;
189
190   if (fRecoParam)
191     {
192       HLTDebug("Deleting fRecoParam");
193       delete fRecoParam;
194       fRecoParam = 0;
195     }
196
197   return 0;
198 }
199
200 int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponentEventData& evtData, 
201                                             const AliHLTComponentBlockData* blocks, 
202                                             AliHLTComponent_TriggerData& /*trigData*/, 
203                                             AliHLTUInt8_t* outputPtr, 
204                                             AliHLTUInt32_t& size, 
205                                             vector<AliHLTComponent_BlockData>& outputBlocks )
206 {
207   // Process an event
208
209 #ifdef HAVE_VALGRIND_CALLGRIND_H
210   if (evtData.fEventID == 10)
211     CALLGRIND_START_INSTRUMENTATION;
212
213   if(GetFirstInputBlock(kAliHLTDataTypeEOR))
214     CALLGRIND_STOP_INSTRUMENTATION;
215 #endif
216
217   if(!IsDataEvent())return 0;
218
219   HLTDebug( "NofBlocks %i", evtData.fBlockCnt );
220   // Process an event
221   AliHLTUInt32_t totalSize = 0, offset = 0;
222
223   //implement a usage of the following
224   //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
225   //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
226   //   void *triggerData = trigData.fData;
227   //HLTDebug( "Trigger data received. Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
228
229   // Loop over all input blocks in the event
230   AliHLTComponentDataType expectedDataType = (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
231   for ( UInt_t iBlock = 0; iBlock < evtData.fBlockCnt; iBlock++ )
232     {      
233       const AliHLTComponentBlockData &block = blocks[iBlock];
234       // lets not use the internal TRD data types here : AliHLTTRDDefinitions::fgkDDLRawDataType
235       // which is depreciated - we use HLT global defs instead
236       //      if ( block.fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD) )
237       AliHLTComponentDataType inputDataType = block.fDataType;
238       if ( inputDataType != expectedDataType)
239         {
240           HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) Wrong received datatype: %s - required datatype: %s; Skipping",
241                     iBlock, evtData.fBlockCnt,
242                     evtData.fEventID, evtData.fEventID, 
243                     DataType2Text(inputDataType).c_str(), 
244                     DataType2Text(expectedDataType).c_str());
245           continue;
246         }
247       else 
248         {
249           HLTDebug("We get the right data type: Block # %i/%i; Event 0x%08LX (%Lu) Received datatype: %s; Block Size: %i",
250                    iBlock, evtData.fBlockCnt,
251                    evtData.fEventID, evtData.fEventID, 
252                    DataType2Text(inputDataType).c_str(),
253                    block.fSize);
254         }
255       
256 #ifndef NDEBUG
257       unsigned long constBase;
258       double inputMultiplier;
259       GetOutputDataSize(constBase,inputMultiplier);
260       if(size<(constBase+block.fSize*inputMultiplier)){
261         HLTWarning("Memory Block given might be too small: %i < %i; Event %Lu", size, constBase+block.fSize*inputMultiplier, evtData.fEventID);
262       }
263 #endif
264
265       // fMemReader->Reset();
266       fMemReader->SetMemory((UChar_t*) block.fPtr, block.fSize);
267
268       AliHLTUInt32_t spec = block.fSpecification;
269       
270       Int_t id = AliHLTTRDUtils::GetSM(spec) + 1024;
271
272       fMemReader->SetEquipmentID(id);
273       
274       fClusterizer->SetMemBlock(outputPtr+offset);
275       Bool_t bclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
276       if(bclustered)
277         {
278           HLTDebug("Clustered successfully");
279         }
280       else
281         {
282           HLTError("Clustering ERROR");
283           return -1;
284         }
285
286       AliHLTUInt32_t addedSize;
287       if(fReconstructor->IsProcessingTracklets()){
288         addedSize = fClusterizer->GetAddedTrSize();
289         totalSize += fClusterizer->GetTrMemBlockSize();  //if IsProcessingTracklets() is enabled we always reserve a data block of size GetTrMemBlockSize() for the tracklets
290         if (addedSize > 0){
291           // Using low-level interface 
292           // with interface classes
293           if ( totalSize > size )
294             {
295               HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
296                        totalSize, size );
297               return EMSGSIZE;
298             }
299
300           // Fill block 
301           AliHLTComponentBlockData bd;
302           FillBlockData( bd );
303           bd.fOffset = offset;
304           bd.fSize = addedSize;
305           bd.fSpecification = block.fSpecification;
306           bd.fDataType = AliHLTTRDDefinitions::fgkMCMtrackletDataType;
307           outputBlocks.push_back( bd );
308           HLTDebug( "BD ptr 0x%x, offset %i, size %i, dataType %s, spec 0x%x ", bd.fPtr, bd.fOffset, bd.fSize, DataType2Text(bd.fDataType).c_str(), spec);
309         }
310         offset = totalSize;
311       }
312
313       addedSize = fClusterizer->GetAddedClSize();
314       if (addedSize > 0){
315         
316         Int_t* nTimeBins = (Int_t*)(outputPtr+offset+fClusterizer->GetAddedClSize());
317         *nTimeBins = fClusterizer->GetNTimeBins();
318         addedSize += sizeof(*nTimeBins);
319
320         totalSize += addedSize;
321         if ( totalSize > size )
322           {
323             HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
324                      totalSize, size );
325             return EMSGSIZE;
326           }
327
328         // Fill block 
329         AliHLTComponentBlockData bd;
330         FillBlockData( bd );
331         bd.fOffset = offset;
332         bd.fSize = addedSize;
333         bd.fSpecification = block.fSpecification;
334         bd.fDataType = AliHLTTRDDefinitions::fgkClusterDataType;
335         outputBlocks.push_back( bd );
336         HLTDebug( "BD ptr 0x%x, offset %i, size %i, dataType %s, spec 0x%x ", bd.fPtr, bd.fOffset, bd.fSize, DataType2Text(bd.fDataType).c_str(), spec);
337         offset = totalSize;
338       }
339       else{
340         HLTDebug("Array of clusters is empty!");
341       }
342     }
343   //fReconstructor->SetClusters(0x0);
344
345   size = totalSize;
346   HLTDebug("Event is done. size written to the output is %i", size);
347   return 0;
348 }
349
350 void AliHLTTRDClusterizerComponent::PrintObject(
351 #ifdef __DEBUG
352                 TClonesArray* inClustersArray
353 #else
354                 TClonesArray*
355 #endif
356         )
357 {
358 #ifdef __DEBUG
359   AliTRDcluster* cluster=0x0;
360   
361   for (Int_t i=0; i < inClustersArray->GetEntriesFast(); i++){
362     cluster = dynamic_cast<AliTRDcluster*>(inClustersArray->At(i));
363     HLTDebug("cluster[%i]",i);
364     HLTDebug("  PadCol = %i; PadRow = %i; PadTime = %i", cluster->GetPadCol(), cluster->GetPadRow(), cluster->GetPadTime());
365     HLTDebug("  Detector = %i, Amplitude = %f, Center = %f", cluster->GetDetector(), cluster->GetQ(), cluster->GetCenter());
366     HLTDebug("  LocalTimeBin =  %i; NPads = %i; maskedPosition: %s, status: %s", cluster->GetLocalTimeBin(), cluster->GetNPads(),cluster->GetPadMaskedPosition(),cluster->GetPadMaskedPosition());
367   }
368 #endif
369 }
370
371 int AliHLTTRDClusterizerComponent::Configure(const char* arguments){
372   int iResult=0;
373   if (!arguments) return iResult;
374   
375   TString allArgs=arguments;
376   TString argument;
377   int bMissingParam=0;
378
379   TObjArray* pTokens=allArgs.Tokenize(" ");
380   if (pTokens) {
381     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
382       argument=((TObjString*)pTokens->At(i))->GetString();
383       if (argument.IsNull()) continue;
384       
385       if (argument.CompareTo("output_percentage")==0) {
386         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
387         HLTInfo("Setting output percentage to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
388         fOutputPercentage=((TObjString*)pTokens->At(i))->GetString().Atoi();
389         continue;
390       } 
391       else if (argument.CompareTo("-geometry")==0) {
392         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
393         HLTInfo("Setting geometry to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
394         fgeometryFileName=((TObjString*)pTokens->At(i))->GetString();
395         continue;
396       } 
397       else if (argument.CompareTo("-lowflux")==0) {
398         fRecoParamType = 0;
399         HLTInfo("Low flux reconstruction selected");
400         continue;
401       }
402       else if (argument.CompareTo("-highflux")==0) {
403         fRecoParamType = 1;
404         HLTInfo("High flux reconstruction selected");
405         continue;
406       }
407       else if (argument.CompareTo("-cosmics")==0) {
408         fRecoParamType = 2;
409         HLTInfo("Cosmics reconstruction selected");
410         continue;
411       }
412       else if (argument.CompareTo("-simulation")==0) {
413         fRecoDataType = 0;
414         HLTInfo("Awaiting simulated data");
415         continue;
416       }
417       else if (argument.CompareTo("-experiment")==0) {
418         fRecoDataType = 1;
419         HLTInfo("Awaiting real data");
420         continue;
421       }
422       else if (argument.CompareTo("-processTracklets")==0) {
423         fProcessTracklets = kTRUE;
424         HLTInfo("Writing L1 tracklets to output");
425         continue;
426       }
427       else if (argument.CompareTo("-noZS")==0) {
428         fOutputPercentage = 10;
429         HLTInfo("Awaiting non zero surpressed data");
430         continue;
431       }
432       else if (argument.CompareTo("-HLTflag")==0) {
433         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
434         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
435         if (toCompareTo.CompareTo("yes")==0){
436           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
437           fHLTflag=kTRUE;
438         }
439         else if (toCompareTo.CompareTo("no")==0){
440           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
441           fHLTflag=kFALSE;
442         }
443         else {
444           HLTError("unknown argument for HLTflag: %s", toCompareTo.Data());
445           iResult=-EINVAL;
446           break;
447         }
448         continue;
449       }
450       else if (argument.CompareTo("-faststreamer")==0) {
451         fHLTstreamer = kTRUE;
452         HLTInfo("Useing fast raw streamer");
453         continue;
454       }
455       else if (argument.CompareTo("-nofaststreamer")==0) {
456         fHLTstreamer = kFALSE;
457         HLTInfo("Don't use fast raw streamer");
458         continue;
459       }
460       else if (argument.CompareTo("-tailcancellation")==0) {
461         fTC = kTRUE;
462         HLTInfo("Useing tailcancellation");
463         continue;
464       }
465       else if (argument.CompareTo("-rawver")==0) {
466         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
467         HLTInfo("Raw data version is: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
468         fRawDataVersion=((TObjString*)pTokens->At(i))->GetString().Atoi();
469         continue;
470       } 
471       else if (argument.CompareTo("-highLevelOutput")==0) {
472         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
473         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
474         if (toCompareTo.CompareTo("yes")==0){
475           HLTWarning("Setting highLevelOutput to: %s", toCompareTo.Data());
476           fHighLevelOutput=kTRUE;
477         }
478         else if (toCompareTo.CompareTo("no")==0){
479           HLTInfo("Setting highLevelOutput to: %s", toCompareTo.Data());
480           fHighLevelOutput=kFALSE;
481         }
482         else {
483           HLTError("unknown argument for highLevelOutput: %s", toCompareTo.Data());
484           iResult=-EINVAL;
485           break;
486         }
487         continue;
488       }
489       else if (argument.CompareTo("-emulateHLToutput")==0) {
490         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
491         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
492         if (toCompareTo.CompareTo("yes")==0){
493           HLTWarning("Setting emulateHLToutput to: %s", toCompareTo.Data());
494           fEmulateHLTClusters=kTRUE;
495         }
496         else if (toCompareTo.CompareTo("no")==0){
497           HLTInfo("Setting emulateHLToutput to: %s", toCompareTo.Data());
498           fEmulateHLTClusters=kFALSE;
499         }
500         else {
501           HLTError("unknown argument for emulateHLToutput: %s", toCompareTo.Data());
502           iResult=-EINVAL;
503           break;
504         }
505         continue;
506       }
507       else if (argument.CompareTo("-yPosMethod")==0) {
508         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
509         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
510         if (toCompareTo.CompareTo("COG")==0){
511           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
512           fyPosMethod=0;
513         }
514         else if (toCompareTo.CompareTo("LUT")==0){
515           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
516           fyPosMethod=1;
517         }
518         else if (toCompareTo.CompareTo("Gauss")==0){
519           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
520           fyPosMethod=2;
521         }
522         else {
523           HLTError("unknown argument for yPosMethod: %s", toCompareTo.Data());
524           iResult=-EINVAL;
525           break;
526         }
527         continue;
528       }
529       
530       else {
531         HLTError("unknown argument: %s", argument.Data());
532         iResult=-EINVAL;
533         break;
534       }
535     }
536     delete pTokens;
537   }
538   if (bMissingParam) {
539     HLTError("missing parameter for argument %s", argument.Data());
540     iResult=-EINVAL;
541   }
542   if(iResult>=0){
543     iResult=SetParams();
544   }
545   return iResult;
546 }
547
548 int AliHLTTRDClusterizerComponent::SetParams()
549 {
550   Int_t iResult=0;
551   if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
552     HLTError("DefaultStorage is not set in CDBManager");
553     return -EINVAL;
554   }
555   if(AliCDBManager::Instance()->GetRun()<0){
556     HLTError("Run Number is not set in CDBManager");
557     return -EINVAL;
558   }
559   HLTInfo("CDB default storage: %s; RunNo: %i", (AliCDBManager::Instance()->GetDefaultStorage()->GetBaseFolder()).Data(), AliCDBManager::Instance()->GetRun());
560
561   if(!AliGeomManager::GetGeometry()){
562     if(fgeometryFileName.CompareTo("")==0 || !TFile::Open(fgeometryFileName.Data())){
563       HLTInfo("Loading standard geometry file");
564       AliGeomManager::LoadGeometry();
565     }else{
566       HLTWarning("Loading NON-standard geometry file");
567       AliGeomManager::LoadGeometry(fgeometryFileName.Data());
568     }
569     if(!AliGeomManager::GetGeometry()){
570       HLTError("Could not load geometry");
571       return -EINVAL;
572     }
573     HLTInfo("Applying Alignment from CDB object");
574     AliGeomManager::ApplyAlignObjsFromCDB("TRD");
575   }
576   else{
577     HLTInfo("Geometry Already Loaded!");
578   }
579
580   if(fReconstructor->GetRecoParam()){
581     fRecoParam = new AliTRDrecoParam(*fReconstructor->GetRecoParam());
582     HLTInfo("RecoParam already set!");
583   }else{
584     if(fRecoParamType == 0){
585 #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
586       if(fHLTflag){
587         HLTInfo("Low flux HLT params init.");
588         fRecoParam = AliTRDrecoParam::GetLowFluxHLTParam();
589       }else
590 #endif
591         {
592           HLTInfo("Low flux params init.");
593           fRecoParam = AliTRDrecoParam::GetLowFluxParam();
594         }
595     }
596     if(fRecoParamType == 1){
597 #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
598       if(fHLTflag){
599         HLTInfo("High flux HLT params init.");
600         fRecoParam = AliTRDrecoParam::GetHighFluxHLTParam();
601       }else
602 #endif
603         {
604           HLTInfo("High flux params init.");
605           fRecoParam = AliTRDrecoParam::GetHighFluxParam();
606         }
607     }
608     if(fRecoParamType == 2){
609       HLTInfo("Cosmic Test params init.");
610       fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
611     }
612   }
613
614   if (!fRecoParam)
615     {
616       HLTError("No reco params initialized. Sniffing big trouble!");
617       return -EINVAL;
618     }
619
620   if(fTC){fRecoParam->SetTailCancelation(kTRUE); HLTDebug("Enableing Tail Cancelation"); }
621   else{fRecoParam->SetTailCancelation(kFALSE); HLTDebug("Disableing Tail Cancelation"); }
622
623   switch(fyPosMethod){
624   case 0: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kFALSE); break;
625   case 1: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kTRUE); break;
626   case 2: fRecoParam->SetGAUS(kTRUE); fRecoParam->SetLUT(kFALSE); break;
627   }
628
629   fRecoParam->SetStreamLevel(AliTRDrecoParam::kClusterizer, 0);
630   fReconstructor->SetRecoParam(fRecoParam);
631
632   if(!fClusterizer){
633     fClusterizer = new AliHLTTRDClusterizer("TRDCclusterizer", "TRDCclusterizer");  
634     HLTDebug("TRDClusterizer at 0x%x", fClusterizer);
635   }
636
637   TString recoOptions="!cw";
638   if(fHLTflag){
639     recoOptions += ",hlt";
640     
641     // we transfer clusters that do no contain the XYZ coodrinates (AliHLTTRDCluster),
642     // thus this coordinate transformation ist useless
643 #ifndef HAVE_NOT_ALITRD_CLUSTERIZER_r42837
644     fClusterizer->SetSkipTransform();
645 #endif
646   }
647   if(fProcessTracklets) recoOptions += ",tp";
648   else  recoOptions += ",!tp";
649
650   HLTDebug("Reconstructor options are: %s",recoOptions.Data());
651   fReconstructor->SetOption(recoOptions.Data());
652
653   if (fRecoDataType < 0 || fRecoDataType > 1)
654     {
655       HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
656       fRecoDataType = 0;
657     }
658
659   if (fRecoDataType == 0)
660     {
661       HLTDebug("Data type expected is SIMULATION!");
662     }
663
664   if (fRecoDataType == 1)
665     {
666       HLTDebug("Data type expected is EXPERIMENT!");
667     }
668
669   fClusterizer->SetRawVersion(fRawDataVersion);
670
671   return iResult;
672 }
673
674 int AliHLTTRDClusterizerComponent::Reconfigure(const char* cdbEntry, const char* chainId)
675 {
676   // see header file for class documentation
677
678   int iResult=0;
679   const char* path="HLT/ConfigTRD/ClusterizerComponent";
680   const char* defaultNotify="";
681   if (cdbEntry) {
682     path=cdbEntry;
683     defaultNotify=" (default)";
684   }
685   if (path) {
686     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
687     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
688     if (pEntry) {
689       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
690       if (pString) {
691         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
692         iResult=Configure(pString->GetString().Data());
693       } else {
694         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
695       }
696     } else {
697       HLTError("cannot fetch object \"%s\" from CDB", path);
698     }
699   }
700
701   return iResult;
702 }
703
704 void AliHLTTRDClusterizerComponent::GetOCDBObjectDescription(TMap* const targetMap){
705   // Get a list of OCDB object description needed for the particular component
706   if (!targetMap) return;
707   targetMap->Add(new TObjString("HLT/ConfigTRD/ClusterizerComponent"), new TObjString("component arguments"));
708   targetMap->Add(new TObjString("TRD/Calib/ChamberGainFactor"), new TObjString("gain factor of chambers"));
709   targetMap->Add(new TObjString("TRD/Calib/ChamberT0"), new TObjString("T0 of chambers"));
710   targetMap->Add(new TObjString("TRD/Calib/ChamberVdrift"), new TObjString("drift velocity of chambers"));
711   targetMap->Add(new TObjString("TRD/Calib/DetNoise"), new TObjString("noise of chambers"));
712   targetMap->Add(new TObjString("TRD/Calib/LocalGainFactor"), new TObjString("per pad gain factor"));
713   targetMap->Add(new TObjString("TRD/Calib/LocalT0"), new TObjString("per pad T0"));
714   targetMap->Add(new TObjString("TRD/Calib/LocalVdrift"), new TObjString("per pad drift velocity"));
715   targetMap->Add(new TObjString("TRD/Calib/PadNoise"), new TObjString("per pad noise"));
716   targetMap->Add(new TObjString("TRD/Calib/PadStatus"), new TObjString("pad status"));
717   targetMap->Add(new TObjString("TRD/Calib/PRFWidth"), new TObjString("pad response function"));
718 }