]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
fdc1aac234b2d0a7004f6405957a3e9e40896df3
[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   if (evtData.fEventID == 10)
211     CALLGRIND_START_INSTRUMENTATION;
212
213   if(!IsDataEvent())return 0;
214
215   HLTDebug( "NofBlocks %i", evtData.fBlockCnt );
216   // Process an event
217   AliHLTUInt32_t totalSize = 0, offset = 0;
218
219   //implement a usage of the following
220   //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
221   //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
222   //   void *triggerData = trigData.fData;
223   //HLTDebug( "Trigger data received. Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
224
225   // Loop over all input blocks in the event
226   AliHLTComponentDataType expectedDataType = (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
227   for ( UInt_t iBlock = 0; iBlock < evtData.fBlockCnt; iBlock++ )
228     {      
229       const AliHLTComponentBlockData &block = blocks[iBlock];
230       // lets not use the internal TRD data types here : AliHLTTRDDefinitions::fgkDDLRawDataType
231       // which is depreciated - we use HLT global defs instead
232       //      if ( block.fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD) )
233       AliHLTComponentDataType inputDataType = block.fDataType;
234       if ( inputDataType != expectedDataType)
235         {
236           HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) Wrong received datatype: %s - required datatype: %s; Skipping",
237                     iBlock, evtData.fBlockCnt,
238                     evtData.fEventID, evtData.fEventID, 
239                     DataType2Text(inputDataType).c_str(), 
240                     DataType2Text(expectedDataType).c_str());
241           if(block.fDataType == kAliHLTDataTypeEOR)
242             CALLGRIND_STOP_INSTRUMENTATION;
243           continue;
244         }
245       else 
246         {
247           HLTDebug("We get the right data type: Block # %i/%i; Event 0x%08LX (%Lu) Received datatype: %s; Block Size: %i",
248                    iBlock, evtData.fBlockCnt,
249                    evtData.fEventID, evtData.fEventID, 
250                    DataType2Text(inputDataType).c_str(),
251                    block.fSize);
252         }
253       
254 #ifndef NDEBUG
255       unsigned long constBase;
256       double inputMultiplier;
257       GetOutputDataSize(constBase,inputMultiplier);
258       if(size<(constBase+block.fSize*inputMultiplier)){
259         HLTWarning("Memory Block given might be too small: %i < %i; Event %Lu", size, constBase+block.fSize*inputMultiplier, evtData.fEventID);
260       }
261 #endif
262
263       // fMemReader->Reset();
264       fMemReader->SetMemory((UChar_t*) block.fPtr, block.fSize);
265
266       AliHLTUInt32_t spec = block.fSpecification;
267       
268       Int_t id = AliHLTTRDUtils::GetSM(spec) + 1024;
269
270       fMemReader->SetEquipmentID(id);
271       
272       fClusterizer->SetMemBlock(outputPtr+offset);
273       Bool_t bclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
274       if(bclustered)
275         {
276           HLTDebug("Clustered successfully");
277         }
278       else
279         {
280           HLTError("Clustering ERROR");
281           return -1;
282         }
283
284       AliHLTUInt32_t addedSize;
285       if(fReconstructor->IsProcessingTracklets()){
286         addedSize = fClusterizer->GetAddedTrSize();
287         totalSize += fClusterizer->GetTrMemBlockSize();  //if IsProcessingTracklets() is enabled we always reserve a data block of size GetTrMemBlockSize() for the tracklets
288         if (addedSize > 0){
289           // Using low-level interface 
290           // with interface classes
291           if ( totalSize > size )
292             {
293               HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
294                        totalSize, size );
295               return EMSGSIZE;
296             }
297
298           // Fill block 
299           AliHLTComponentBlockData bd;
300           FillBlockData( bd );
301           bd.fOffset = offset;
302           bd.fSize = addedSize;
303           bd.fSpecification = block.fSpecification;
304           bd.fDataType = AliHLTTRDDefinitions::fgkMCMtrackletDataType;
305           outputBlocks.push_back( bd );
306           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);
307         }
308         offset = totalSize;
309       }
310
311       addedSize = fClusterizer->GetAddedClSize();
312       if (addedSize > 0){
313         
314         Int_t* nTimeBins = (Int_t*)(outputPtr+offset+fClusterizer->GetAddedClSize());
315         *nTimeBins = fClusterizer->GetNTimeBins();
316         addedSize += sizeof(*nTimeBins);
317
318         totalSize += addedSize;
319         if ( totalSize > size )
320           {
321             HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
322                      totalSize, size );
323             return EMSGSIZE;
324           }
325
326         // Fill block 
327         AliHLTComponentBlockData bd;
328         FillBlockData( bd );
329         bd.fOffset = offset;
330         bd.fSize = addedSize;
331         bd.fSpecification = block.fSpecification;
332         bd.fDataType = AliHLTTRDDefinitions::fgkClusterDataType;
333         outputBlocks.push_back( bd );
334         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);
335         offset = totalSize;
336       }
337       else{
338         HLTDebug("Array of clusters is empty!");
339       }
340     }
341   fReconstructor->SetClusters(0x0);
342
343   size = totalSize;
344   HLTDebug("Event is done. size written to the output is %i", size);
345   return 0;
346 }
347
348 void AliHLTTRDClusterizerComponent::PrintObject( TClonesArray* inClustersArray)
349 {
350   AliTRDcluster* cluster=0x0;
351   
352   for (Int_t i=0; i < inClustersArray->GetEntriesFast(); i++){
353     cluster = dynamic_cast<AliTRDcluster*>(inClustersArray->At(i));
354     HLTDebug("cluster[%i]",i);
355     HLTDebug("  PadCol = %i; PadRow = %i; PadTime = %i", cluster->GetPadCol(), cluster->GetPadRow(), cluster->GetPadTime());
356     HLTDebug("  Detector = %i, Amplitude = %f, Center = %f", cluster->GetDetector(), cluster->GetQ(), cluster->GetCenter());
357     HLTDebug("  LocalTimeBin =  %i; NPads = %i; maskedPosition: %s, status: %s", cluster->GetLocalTimeBin(), cluster->GetNPads(),cluster->GetPadMaskedPosition(),cluster->GetPadMaskedPosition());
358   }
359   
360 }
361
362 int AliHLTTRDClusterizerComponent::Configure(const char* arguments){
363   int iResult=0;
364   if (!arguments) return iResult;
365   
366   TString allArgs=arguments;
367   TString argument;
368   int bMissingParam=0;
369
370   TObjArray* pTokens=allArgs.Tokenize(" ");
371   if (pTokens) {
372     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
373       argument=((TObjString*)pTokens->At(i))->GetString();
374       if (argument.IsNull()) continue;
375       
376       if (argument.CompareTo("output_percentage")==0) {
377         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
378         HLTInfo("Setting output percentage to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
379         fOutputPercentage=((TObjString*)pTokens->At(i))->GetString().Atoi();
380         continue;
381       } 
382       else if (argument.CompareTo("-geometry")==0) {
383         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
384         HLTInfo("Setting geometry to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
385         fgeometryFileName=((TObjString*)pTokens->At(i))->GetString();
386         continue;
387       } 
388       else if (argument.CompareTo("-lowflux")==0) {
389         fRecoParamType = 0;
390         HLTInfo("Low flux reconstruction selected");
391         continue;
392       }
393       else if (argument.CompareTo("-highflux")==0) {
394         fRecoParamType = 1;
395         HLTInfo("High flux reconstruction selected");
396         continue;
397       }
398       else if (argument.CompareTo("-cosmics")==0) {
399         fRecoParamType = 2;
400         HLTInfo("Cosmics reconstruction selected");
401         continue;
402       }
403       else if (argument.CompareTo("-simulation")==0) {
404         fRecoDataType = 0;
405         HLTInfo("Awaiting simulated data");
406         continue;
407       }
408       else if (argument.CompareTo("-experiment")==0) {
409         fRecoDataType = 1;
410         HLTInfo("Awaiting real data");
411         continue;
412       }
413       else if (argument.CompareTo("-processTracklets")==0) {
414         fProcessTracklets = kTRUE;
415         HLTInfo("Writing L1 tracklets to output");
416         continue;
417       }
418       else if (argument.CompareTo("-noZS")==0) {
419         fOutputPercentage = 10;
420         HLTInfo("Awaiting non zero surpressed data");
421         continue;
422       }
423       else if (argument.CompareTo("-HLTflag")==0) {
424         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
425         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
426         if (toCompareTo.CompareTo("yes")==0){
427           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
428           fHLTflag=kTRUE;
429         }
430         else if (toCompareTo.CompareTo("no")==0){
431           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
432           fHLTflag=kFALSE;
433         }
434         else {
435           HLTError("unknown argument for HLTflag: %s", toCompareTo.Data());
436           iResult=-EINVAL;
437           break;
438         }
439         continue;
440       }
441       else if (argument.CompareTo("-faststreamer")==0) {
442         fHLTstreamer = kTRUE;
443         HLTInfo("Useing fast raw streamer");
444         continue;
445       }
446       else if (argument.CompareTo("-nofaststreamer")==0) {
447         fHLTstreamer = kFALSE;
448         HLTInfo("Don't use fast raw streamer");
449         continue;
450       }
451       else if (argument.CompareTo("-tailcancellation")==0) {
452         fTC = kTRUE;
453         HLTInfo("Useing tailcancellation");
454         continue;
455       }
456       else if (argument.CompareTo("-rawver")==0) {
457         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
458         HLTInfo("Raw data version is: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
459         fRawDataVersion=((TObjString*)pTokens->At(i))->GetString().Atoi();
460         continue;
461       } 
462       else if (argument.CompareTo("-highLevelOutput")==0) {
463         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
464         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
465         if (toCompareTo.CompareTo("yes")==0){
466           HLTWarning("Setting highLevelOutput to: %s", toCompareTo.Data());
467           fHighLevelOutput=kTRUE;
468         }
469         else if (toCompareTo.CompareTo("no")==0){
470           HLTInfo("Setting highLevelOutput to: %s", toCompareTo.Data());
471           fHighLevelOutput=kFALSE;
472         }
473         else {
474           HLTError("unknown argument for highLevelOutput: %s", toCompareTo.Data());
475           iResult=-EINVAL;
476           break;
477         }
478         continue;
479       }
480       else if (argument.CompareTo("-emulateHLToutput")==0) {
481         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
482         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
483         if (toCompareTo.CompareTo("yes")==0){
484           HLTWarning("Setting emulateHLToutput to: %s", toCompareTo.Data());
485           fEmulateHLTClusters=kTRUE;
486         }
487         else if (toCompareTo.CompareTo("no")==0){
488           HLTInfo("Setting emulateHLToutput to: %s", toCompareTo.Data());
489           fEmulateHLTClusters=kFALSE;
490         }
491         else {
492           HLTError("unknown argument for emulateHLToutput: %s", toCompareTo.Data());
493           iResult=-EINVAL;
494           break;
495         }
496         continue;
497       }
498       else if (argument.CompareTo("-yPosMethod")==0) {
499         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
500         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
501         if (toCompareTo.CompareTo("COG")==0){
502           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
503           fyPosMethod=0;
504         }
505         else if (toCompareTo.CompareTo("LUT")==0){
506           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
507           fyPosMethod=1;
508         }
509         else if (toCompareTo.CompareTo("Gauss")==0){
510           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
511           fyPosMethod=2;
512         }
513         else {
514           HLTError("unknown argument for yPosMethod: %s", toCompareTo.Data());
515           iResult=-EINVAL;
516           break;
517         }
518         continue;
519       }
520       
521       else {
522         HLTError("unknown argument: %s", argument.Data());
523         iResult=-EINVAL;
524         break;
525       }
526     }
527     delete pTokens;
528   }
529   if (bMissingParam) {
530     HLTError("missing parameter for argument %s", argument.Data());
531     iResult=-EINVAL;
532   }
533   if(iResult>=0){
534     iResult=SetParams();
535   }
536   return iResult;
537 }
538
539 int AliHLTTRDClusterizerComponent::SetParams()
540 {
541   Int_t iResult=0;
542   if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
543     HLTError("DefaultStorage is not set in CDBManager");
544     return -EINVAL;
545   }
546   if(AliCDBManager::Instance()->GetRun()<0){
547     HLTError("Run Number is not set in CDBManager");
548     return -EINVAL;
549   }
550   HLTInfo("CDB default storage: %s; RunNo: %i", (AliCDBManager::Instance()->GetDefaultStorage()->GetBaseFolder()).Data(), AliCDBManager::Instance()->GetRun());
551
552   if(!AliGeomManager::GetGeometry()){
553     if(fgeometryFileName.CompareTo("")==0 || !TFile::Open(fgeometryFileName.Data())){
554       HLTInfo("Loading standard geometry file");
555       AliGeomManager::LoadGeometry();
556     }else{
557       HLTWarning("Loading NON-standard geometry file");
558       AliGeomManager::LoadGeometry(fgeometryFileName.Data());
559     }
560     if(!AliGeomManager::GetGeometry()){
561       HLTError("Could not load geometry");
562       return -EINVAL;
563     }
564     HLTInfo("Applying Alignment from CDB object");
565     AliGeomManager::ApplyAlignObjsFromCDB("TRD");
566   }
567   else{
568     HLTInfo("Geometry Already Loaded!");
569   }
570
571   if(fReconstructor->GetRecoParam()){
572     fRecoParam = new AliTRDrecoParam(*fReconstructor->GetRecoParam());
573     HLTInfo("RecoParam already set!");
574   }else{
575     if(fRecoParamType == 0){
576       HLTDebug("Low flux params init.");
577       fRecoParam = AliTRDrecoParam::GetLowFluxParam();
578     }
579     if(fRecoParamType == 1){
580       HLTDebug("High flux params init.");
581       fRecoParam = AliTRDrecoParam::GetHighFluxParam();
582     }
583     if(fRecoParamType == 2){
584       HLTDebug("Cosmic Test params init.");
585       fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
586     }
587   }
588
589   if (!fRecoParam)
590     {
591       HLTError("No reco params initialized. Sniffing big trouble!");
592       return -EINVAL;
593     }
594
595   if(fTC){fRecoParam->SetTailCancelation(kTRUE); HLTDebug("Enableing Tail Cancelation"); }
596   else{fRecoParam->SetTailCancelation(kFALSE); HLTDebug("Disableing Tail Cancelation"); }
597
598   switch(fyPosMethod){
599   case 0: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kFALSE); break;
600   case 1: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kTRUE); break;
601   case 2: fRecoParam->SetGAUS(kTRUE); fRecoParam->SetLUT(kFALSE); break;
602   }
603
604   fRecoParam->SetStreamLevel(AliTRDrecoParam::kClusterizer, 0);
605   fReconstructor->SetRecoParam(fRecoParam);
606
607   TString recoOptions="!cw";
608   if(fHLTflag)
609     recoOptions += ",hlt";
610   if(fProcessTracklets) recoOptions += ",tp";
611   else  recoOptions += ",!tp";
612
613   HLTDebug("Reconstructor options are: %s",recoOptions.Data());
614   fReconstructor->SetOption(recoOptions.Data());
615
616   if (fRecoDataType < 0 || fRecoDataType > 1)
617     {
618       HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
619       fRecoDataType = 0;
620     }
621
622   if (fRecoDataType == 0)
623     {
624       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDsimStream);
625       HLTDebug("Data type expected is SIMULATION!");
626     }
627
628   if (fRecoDataType == 1)
629     {
630       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDrealStream);
631       HLTDebug("Data type expected is EXPERIMENT!");
632     }
633
634 #ifndef HAVE_NOT_ALITRD_RAWSTREAM_r39608
635   if(fHLTstreamer){
636     AliTRDrawStreamBase::SetRawStreamVersion("default");
637     HLTDebug("fast rawstreamer used");
638   }else{
639     AliTRDrawStreamBase::SetRawStreamVersion("FAST");
640     HLTDebug("old rawstreamer used");
641   }
642 #else
643   if(fHLTstreamer){
644       AliTRDrawStreamBase::SetRawStreamVersion("FAST");
645       HLTDebug("fast rawstreamer used");  
646     }
647 #endif
648
649   if(!fClusterizer){
650     fClusterizer = new AliHLTTRDClusterizer("TRDCclusterizer", "TRDCclusterizer");  
651     HLTDebug("TRDClusterizer at 0x%x", fClusterizer);
652   }
653
654   fClusterizer->SetRawVersion(fRawDataVersion);
655
656   return iResult;
657 }
658
659 int AliHLTTRDClusterizerComponent::Reconfigure(const char* cdbEntry, const char* chainId)
660 {
661   // see header file for class documentation
662
663   int iResult=0;
664   const char* path="HLT/ConfigTRD/ClusterizerComponent";
665   const char* defaultNotify="";
666   if (cdbEntry) {
667     path=cdbEntry;
668     defaultNotify=" (default)";
669   }
670   if (path) {
671     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
672     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
673     if (pEntry) {
674       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
675       if (pString) {
676         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
677         iResult=Configure(pString->GetString().Data());
678       } else {
679         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
680       }
681     } else {
682       HLTError("cannot fetch object \"%s\" from CDB", path);
683     }
684   }
685
686   return iResult;
687 }
688
689 void AliHLTTRDClusterizerComponent::GetOCDBObjectDescription(TMap* const targetMap){
690   // Get a list of OCDB object description needed for the particular component
691   if (!targetMap) return;
692   targetMap->Add(new TObjString("HLT/ConfigTRD/ClusterizerComponent"), new TObjString("component arguments"));
693   targetMap->Add(new TObjString("TRD/Calib/ChamberGainFactor"), new TObjString("gain factor of chambers"));
694   targetMap->Add(new TObjString("TRD/Calib/ChamberT0"), new TObjString("T0 of chambers"));
695   targetMap->Add(new TObjString("TRD/Calib/ChamberVdrift"), new TObjString("drift velocity of chambers"));
696   targetMap->Add(new TObjString("TRD/Calib/DetNoise"), new TObjString("noise of chambers"));
697   targetMap->Add(new TObjString("TRD/Calib/LocalGainFactor"), new TObjString("per pad gain factor"));
698   targetMap->Add(new TObjString("TRD/Calib/LocalT0"), new TObjString("per pad T0"));
699   targetMap->Add(new TObjString("TRD/Calib/LocalVdrift"), new TObjString("per pad drift velocity"));
700   targetMap->Add(new TObjString("TRD/Calib/PadNoise"), new TObjString("per pad noise"));
701   targetMap->Add(new TObjString("TRD/Calib/PadStatus"), new TObjString("pad status"));
702   targetMap->Add(new TObjString("TRD/Calib/PRFWidth"), new TObjString("pad response function"));
703 }