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