]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
fixing warnings
[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(
349 #ifdef __DEBUG
350                 TClonesArray* inClustersArray
351 #else
352                 TClonesArray*
353 #endif
354         )
355 {
356 #ifdef __DEBUG
357   AliTRDcluster* cluster=0x0;
358   
359   for (Int_t i=0; i < inClustersArray->GetEntriesFast(); i++){
360     cluster = dynamic_cast<AliTRDcluster*>(inClustersArray->At(i));
361     HLTDebug("cluster[%i]",i);
362     HLTDebug("  PadCol = %i; PadRow = %i; PadTime = %i", cluster->GetPadCol(), cluster->GetPadRow(), cluster->GetPadTime());
363     HLTDebug("  Detector = %i, Amplitude = %f, Center = %f", cluster->GetDetector(), cluster->GetQ(), cluster->GetCenter());
364     HLTDebug("  LocalTimeBin =  %i; NPads = %i; maskedPosition: %s, status: %s", cluster->GetLocalTimeBin(), cluster->GetNPads(),cluster->GetPadMaskedPosition(),cluster->GetPadMaskedPosition());
365   }
366 #endif
367 }
368
369 int AliHLTTRDClusterizerComponent::Configure(const char* arguments){
370   int iResult=0;
371   if (!arguments) return iResult;
372   
373   TString allArgs=arguments;
374   TString argument;
375   int bMissingParam=0;
376
377   TObjArray* pTokens=allArgs.Tokenize(" ");
378   if (pTokens) {
379     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
380       argument=((TObjString*)pTokens->At(i))->GetString();
381       if (argument.IsNull()) continue;
382       
383       if (argument.CompareTo("output_percentage")==0) {
384         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
385         HLTInfo("Setting output percentage to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
386         fOutputPercentage=((TObjString*)pTokens->At(i))->GetString().Atoi();
387         continue;
388       } 
389       else if (argument.CompareTo("-geometry")==0) {
390         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
391         HLTInfo("Setting geometry to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
392         fgeometryFileName=((TObjString*)pTokens->At(i))->GetString();
393         continue;
394       } 
395       else if (argument.CompareTo("-lowflux")==0) {
396         fRecoParamType = 0;
397         HLTInfo("Low flux reconstruction selected");
398         continue;
399       }
400       else if (argument.CompareTo("-highflux")==0) {
401         fRecoParamType = 1;
402         HLTInfo("High flux reconstruction selected");
403         continue;
404       }
405       else if (argument.CompareTo("-cosmics")==0) {
406         fRecoParamType = 2;
407         HLTInfo("Cosmics reconstruction selected");
408         continue;
409       }
410       else if (argument.CompareTo("-simulation")==0) {
411         fRecoDataType = 0;
412         HLTInfo("Awaiting simulated data");
413         continue;
414       }
415       else if (argument.CompareTo("-experiment")==0) {
416         fRecoDataType = 1;
417         HLTInfo("Awaiting real data");
418         continue;
419       }
420       else if (argument.CompareTo("-processTracklets")==0) {
421         fProcessTracklets = kTRUE;
422         HLTInfo("Writing L1 tracklets to output");
423         continue;
424       }
425       else if (argument.CompareTo("-noZS")==0) {
426         fOutputPercentage = 10;
427         HLTInfo("Awaiting non zero surpressed data");
428         continue;
429       }
430       else if (argument.CompareTo("-HLTflag")==0) {
431         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
432         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
433         if (toCompareTo.CompareTo("yes")==0){
434           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
435           fHLTflag=kTRUE;
436         }
437         else if (toCompareTo.CompareTo("no")==0){
438           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
439           fHLTflag=kFALSE;
440         }
441         else {
442           HLTError("unknown argument for HLTflag: %s", toCompareTo.Data());
443           iResult=-EINVAL;
444           break;
445         }
446         continue;
447       }
448       else if (argument.CompareTo("-faststreamer")==0) {
449         fHLTstreamer = kTRUE;
450         HLTInfo("Useing fast raw streamer");
451         continue;
452       }
453       else if (argument.CompareTo("-nofaststreamer")==0) {
454         fHLTstreamer = kFALSE;
455         HLTInfo("Don't use fast raw streamer");
456         continue;
457       }
458       else if (argument.CompareTo("-tailcancellation")==0) {
459         fTC = kTRUE;
460         HLTInfo("Useing tailcancellation");
461         continue;
462       }
463       else if (argument.CompareTo("-rawver")==0) {
464         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
465         HLTInfo("Raw data version is: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
466         fRawDataVersion=((TObjString*)pTokens->At(i))->GetString().Atoi();
467         continue;
468       } 
469       else if (argument.CompareTo("-highLevelOutput")==0) {
470         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
471         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
472         if (toCompareTo.CompareTo("yes")==0){
473           HLTWarning("Setting highLevelOutput to: %s", toCompareTo.Data());
474           fHighLevelOutput=kTRUE;
475         }
476         else if (toCompareTo.CompareTo("no")==0){
477           HLTInfo("Setting highLevelOutput to: %s", toCompareTo.Data());
478           fHighLevelOutput=kFALSE;
479         }
480         else {
481           HLTError("unknown argument for highLevelOutput: %s", toCompareTo.Data());
482           iResult=-EINVAL;
483           break;
484         }
485         continue;
486       }
487       else if (argument.CompareTo("-emulateHLToutput")==0) {
488         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
489         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
490         if (toCompareTo.CompareTo("yes")==0){
491           HLTWarning("Setting emulateHLToutput to: %s", toCompareTo.Data());
492           fEmulateHLTClusters=kTRUE;
493         }
494         else if (toCompareTo.CompareTo("no")==0){
495           HLTInfo("Setting emulateHLToutput to: %s", toCompareTo.Data());
496           fEmulateHLTClusters=kFALSE;
497         }
498         else {
499           HLTError("unknown argument for emulateHLToutput: %s", toCompareTo.Data());
500           iResult=-EINVAL;
501           break;
502         }
503         continue;
504       }
505       else if (argument.CompareTo("-yPosMethod")==0) {
506         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
507         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
508         if (toCompareTo.CompareTo("COG")==0){
509           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
510           fyPosMethod=0;
511         }
512         else if (toCompareTo.CompareTo("LUT")==0){
513           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
514           fyPosMethod=1;
515         }
516         else if (toCompareTo.CompareTo("Gauss")==0){
517           HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
518           fyPosMethod=2;
519         }
520         else {
521           HLTError("unknown argument for yPosMethod: %s", toCompareTo.Data());
522           iResult=-EINVAL;
523           break;
524         }
525         continue;
526       }
527       
528       else {
529         HLTError("unknown argument: %s", argument.Data());
530         iResult=-EINVAL;
531         break;
532       }
533     }
534     delete pTokens;
535   }
536   if (bMissingParam) {
537     HLTError("missing parameter for argument %s", argument.Data());
538     iResult=-EINVAL;
539   }
540   if(iResult>=0){
541     iResult=SetParams();
542   }
543   return iResult;
544 }
545
546 int AliHLTTRDClusterizerComponent::SetParams()
547 {
548   Int_t iResult=0;
549   if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
550     HLTError("DefaultStorage is not set in CDBManager");
551     return -EINVAL;
552   }
553   if(AliCDBManager::Instance()->GetRun()<0){
554     HLTError("Run Number is not set in CDBManager");
555     return -EINVAL;
556   }
557   HLTInfo("CDB default storage: %s; RunNo: %i", (AliCDBManager::Instance()->GetDefaultStorage()->GetBaseFolder()).Data(), AliCDBManager::Instance()->GetRun());
558
559   if(!AliGeomManager::GetGeometry()){
560     if(fgeometryFileName.CompareTo("")==0 || !TFile::Open(fgeometryFileName.Data())){
561       HLTInfo("Loading standard geometry file");
562       AliGeomManager::LoadGeometry();
563     }else{
564       HLTWarning("Loading NON-standard geometry file");
565       AliGeomManager::LoadGeometry(fgeometryFileName.Data());
566     }
567     if(!AliGeomManager::GetGeometry()){
568       HLTError("Could not load geometry");
569       return -EINVAL;
570     }
571     HLTInfo("Applying Alignment from CDB object");
572     AliGeomManager::ApplyAlignObjsFromCDB("TRD");
573   }
574   else{
575     HLTInfo("Geometry Already Loaded!");
576   }
577
578   if(fReconstructor->GetRecoParam()){
579     fRecoParam = new AliTRDrecoParam(*fReconstructor->GetRecoParam());
580     HLTInfo("RecoParam already set!");
581   }else{
582     if(fRecoParamType == 0){
583 #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
584       if(fHLTflag){
585         HLTInfo("Low flux HLT params init.");
586         fRecoParam = AliTRDrecoParam::GetLowFluxHLTParam();
587       }else
588 #endif
589         {
590           HLTInfo("Low flux params init.");
591           fRecoParam = AliTRDrecoParam::GetLowFluxParam();
592         }
593     }
594     if(fRecoParamType == 1){
595 #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
596       if(fHLTflag){
597         HLTInfo("High flux HLT params init.");
598         fRecoParam = AliTRDrecoParam::GetHighFluxHLTParam();
599       }else
600 #endif
601         {
602           HLTInfo("High flux params init.");
603           fRecoParam = AliTRDrecoParam::GetHighFluxParam();
604         }
605     }
606     if(fRecoParamType == 2){
607       HLTInfo("Cosmic Test params init.");
608       fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
609     }
610   }
611
612   if (!fRecoParam)
613     {
614       HLTError("No reco params initialized. Sniffing big trouble!");
615       return -EINVAL;
616     }
617
618   if(fTC){fRecoParam->SetTailCancelation(kTRUE); HLTDebug("Enableing Tail Cancelation"); }
619   else{fRecoParam->SetTailCancelation(kFALSE); HLTDebug("Disableing Tail Cancelation"); }
620
621   switch(fyPosMethod){
622   case 0: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kFALSE); break;
623   case 1: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kTRUE); break;
624   case 2: fRecoParam->SetGAUS(kTRUE); fRecoParam->SetLUT(kFALSE); break;
625   }
626
627   fRecoParam->SetStreamLevel(AliTRDrecoParam::kClusterizer, 0);
628   fReconstructor->SetRecoParam(fRecoParam);
629
630   TString recoOptions="!cw";
631   if(fHLTflag)
632     recoOptions += ",hlt";
633   if(fProcessTracklets) recoOptions += ",tp";
634   else  recoOptions += ",!tp";
635
636   HLTDebug("Reconstructor options are: %s",recoOptions.Data());
637   fReconstructor->SetOption(recoOptions.Data());
638
639   if (fRecoDataType < 0 || fRecoDataType > 1)
640     {
641       HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
642       fRecoDataType = 0;
643     }
644
645   if (fRecoDataType == 0)
646     {
647       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDsimStream);
648       HLTDebug("Data type expected is SIMULATION!");
649     }
650
651   if (fRecoDataType == 1)
652     {
653       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDrealStream);
654       HLTDebug("Data type expected is EXPERIMENT!");
655     }
656
657 #ifndef HAVE_NOT_ALITRD_RAWSTREAM_r39608
658   if(fHLTstreamer){
659     AliTRDrawStreamBase::SetRawStreamVersion("default");
660     HLTDebug("fast rawstreamer used");
661   }else{
662     AliTRDrawStreamBase::SetRawStreamVersion("FAST");
663     HLTDebug("old rawstreamer used");
664   }
665 #else
666   if(fHLTstreamer){
667       AliTRDrawStreamBase::SetRawStreamVersion("FAST");
668       HLTDebug("fast rawstreamer used");  
669     }
670 #endif
671
672   if(!fClusterizer){
673     fClusterizer = new AliHLTTRDClusterizer("TRDCclusterizer", "TRDCclusterizer");  
674     HLTDebug("TRDClusterizer at 0x%x", fClusterizer);
675   }
676
677   fClusterizer->SetRawVersion(fRawDataVersion);
678
679   return iResult;
680 }
681
682 int AliHLTTRDClusterizerComponent::Reconfigure(const char* cdbEntry, const char* chainId)
683 {
684   // see header file for class documentation
685
686   int iResult=0;
687   const char* path="HLT/ConfigTRD/ClusterizerComponent";
688   const char* defaultNotify="";
689   if (cdbEntry) {
690     path=cdbEntry;
691     defaultNotify=" (default)";
692   }
693   if (path) {
694     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
695     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
696     if (pEntry) {
697       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
698       if (pString) {
699         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
700         iResult=Configure(pString->GetString().Data());
701       } else {
702         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
703       }
704     } else {
705       HLTError("cannot fetch object \"%s\" from CDB", path);
706     }
707   }
708
709   return iResult;
710 }
711
712 void AliHLTTRDClusterizerComponent::GetOCDBObjectDescription(TMap* const targetMap){
713   // Get a list of OCDB object description needed for the particular component
714   if (!targetMap) return;
715   targetMap->Add(new TObjString("HLT/ConfigTRD/ClusterizerComponent"), new TObjString("component arguments"));
716   targetMap->Add(new TObjString("TRD/Calib/ChamberGainFactor"), new TObjString("gain factor of chambers"));
717   targetMap->Add(new TObjString("TRD/Calib/ChamberT0"), new TObjString("T0 of chambers"));
718   targetMap->Add(new TObjString("TRD/Calib/ChamberVdrift"), new TObjString("drift velocity of chambers"));
719   targetMap->Add(new TObjString("TRD/Calib/DetNoise"), new TObjString("noise of chambers"));
720   targetMap->Add(new TObjString("TRD/Calib/LocalGainFactor"), new TObjString("per pad gain factor"));
721   targetMap->Add(new TObjString("TRD/Calib/LocalT0"), new TObjString("per pad T0"));
722   targetMap->Add(new TObjString("TRD/Calib/LocalVdrift"), new TObjString("per pad drift velocity"));
723   targetMap->Add(new TObjString("TRD/Calib/PadNoise"), new TObjString("per pad noise"));
724   targetMap->Add(new TObjString("TRD/Calib/PadStatus"), new TObjString("pad status"));
725   targetMap->Add(new TObjString("TRD/Calib/PRFWidth"), new TObjString("pad response function"));
726 }