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