]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDTrackerV1Component.cxx
HLT TRD bugfix: correct calculation of the total charge of the cluster (Theo)
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDTrackerV1Component.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   AliHLTTRDTrackerV1Component.cxx
20     @author Theodor Rascanu
21     @date   
22     @brief  A TRDTrackerV1 processing component for the HLT.
23 */
24
25 #if __GNUC__ >= 3
26 using namespace std;
27 #endif
28
29 #include "AliHLTTRDTrackerV1Component.h"
30 #include "AliHLTTRDDefinitions.h"
31 #include "AliHLTTRDTrack.h"
32 #include "AliHLTTRDUtils.h"
33 #include "AliHLTTRDCluster.h"
34
35 #include "TFile.h"
36 #include "TChain.h"
37
38 #include "AliGeomManager.h"
39 #include "AliCDBManager.h"
40 #include "AliCDBStorage.h"
41 #include "AliCDBEntry.h"
42 #include "AliESDEvent.h"
43 #include "AliESDfriend.h"
44
45 #include "AliTRDcalibDB.h"
46 #include "AliTRDReconstructor.h"
47 #include "AliTRDtrackerV1.h"
48 #include "AliTRDrecoParam.h"
49
50 #include <cstdlib>
51 #include <cerrno>
52 #include <string>
53
54 ClassImp(AliHLTTRDTrackerV1Component)
55
56 void AliHLTTRDTrackerV1Component::AliHLTTRDESDEvent::CreateStdContent()
57 {
58   TClonesArray* tracksArray = new TClonesArray("AliESDtrack",0);
59   tracksArray->SetName(AliESDEvent::fgkESDListName[AliESDEvent::kTracks]);
60   AddObject(tracksArray);
61   GetStdContent();
62 }
63
64 void AliHLTTRDTrackerV1Component::AliHLTTRDESDEvent::Streamer(TBuffer &/*R__b*/)
65 {
66   AliFatal("class is for internal us only and not for streaming");
67 }
68
69 AliHLTTRDTrackerV1Component::AliHLTTRDTrackerV1Component():
70   AliHLTProcessor(),
71   fOutputPercentage(100), // By default we copy to the output exactly what we got as input 
72   fTracker(NULL),
73   fRecoParam(NULL),
74   fReconstructor(NULL),
75   fESD(NULL),
76   fClusterArray(NULL),
77   fRecoParamType(-1),
78   fNtimeBins(-1),
79   fPIDmethod(1),
80   fgeometryFileName(""),
81   fHLTflag(kTRUE),
82   fOutputV1Tracks(kTRUE),
83   fHighLevelOutput(kFALSE),
84   fEmulateHLTTracks(kFALSE),
85   fImproveTracklets(kFALSE)
86 {
87   // Default constructor
88
89 }
90
91 AliHLTTRDTrackerV1Component::~AliHLTTRDTrackerV1Component()
92 {
93   // Destructor
94 }
95
96 const char* AliHLTTRDTrackerV1Component::GetComponentID()
97 {
98   // Return the component ID const char *
99   return "TRDTrackerV1"; // The ID of this component
100 }
101
102 void AliHLTTRDTrackerV1Component::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
103 {
104   // Get the list of input data  
105   list.clear(); // We do not have any requirements for our input data type(s).
106   list.push_back(AliHLTTRDDefinitions::fgkClusterDataType);
107 }
108
109 AliHLTComponentDataType AliHLTTRDTrackerV1Component::GetOutputDataType()
110 {
111   // Get the output data type
112   return kAliHLTMultipleDataType;
113 }
114
115 int AliHLTTRDTrackerV1Component::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
116 {
117   // Get the output data types
118   tgtList.clear();
119   tgtList.push_back(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
120   tgtList.push_back(AliHLTTRDDefinitions::fgkTracksDataType);
121   return tgtList.size();
122 }
123
124 void AliHLTTRDTrackerV1Component::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
125 {
126   // Get the output data size
127   constBase = 0;
128   inputMultiplier = fOutputV1Tracks ? 2*((double)fOutputPercentage)/100.0 : 0.5*((double)fOutputPercentage)/100.0;
129   if(sizeof(AliHLTTRDClustersArray::cluster_type) == sizeof(AliHLTTRDCluster)) inputMultiplier *= 28.0/8;
130 }
131
132 // Spawn function, return new instance of this class
133 AliHLTComponent* AliHLTTRDTrackerV1Component::Spawn()
134 {
135   return new AliHLTTRDTrackerV1Component;
136 };
137
138
139 int AliHLTTRDTrackerV1Component::DoInit( int argc, const char** argv )
140 {
141   // perform initialization. We check whether our relative output size is specified in the arguments.
142   int iResult=0;
143
144   fReconstructor = new AliTRDReconstructor();
145   HLTDebug("TRDReconstructor at 0x%x", fReconstructor);
146   fESD = new AliHLTTRDESDEvent();
147   fESD->CreateStdContent();
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(iResult<0) return iResult;
164
165   fTracker = new AliTRDtrackerV1();
166   HLTDebug("TRDTracker at 0x%x", fTracker);
167   fTracker->SetReconstructor(fReconstructor);
168
169   fClusterArray = new TClonesArray("AliTRDcluster"); // would be nice to allocate memory for all clusters here.
170
171   return iResult;
172 }
173
174 int AliHLTTRDTrackerV1Component::DoDeinit()
175 {
176   // Deinitialization of the component
177
178   fTracker->SetClustersOwner(kFALSE);
179   delete fTracker;
180   fTracker = NULL;
181
182   fClusterArray->Delete();
183   delete fClusterArray;
184   fClusterArray = NULL;
185   
186   // We need to set clusters in Reconstructor to null to prevent from 
187   // double deleting, since we delete TClonesArray by ourself.
188   fReconstructor->SetClusters(0x0);
189   delete fReconstructor;
190   fReconstructor = NULL;
191   delete fESD;
192   fESD = NULL;
193   
194   AliTRDcalibDB::Terminate();
195
196   return 0;
197 }
198
199 int AliHLTTRDTrackerV1Component::DoEvent( const AliHLTComponentEventData& evtData, 
200                                           const AliHLTComponentBlockData* blocks, 
201                                           AliHLTComponent_TriggerData& /*trigData*/, 
202                                           AliHLTUInt8_t* outputPtr, 
203                                           AliHLTUInt32_t& size, 
204                                           vector<AliHLTComponent_BlockData>& outputBlocks )
205 {
206   // Process an event
207
208   HLTDebug("NofBlocks %i", evtData.fBlockCnt );
209   
210   AliHLTUInt32_t totalSize = 0, offset = 0;
211
212   AliHLTComponentDataType expectedDataType = AliHLTTRDDefinitions::fgkClusterDataType;
213   for ( unsigned long iBlock = 0; iBlock < evtData.fBlockCnt; iBlock++ ) 
214     {
215       const AliHLTComponentBlockData &block = blocks[iBlock];
216       AliHLTComponentDataType inputDataType = block.fDataType;
217
218       if(inputDataType != expectedDataType)
219         {
220           HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) Wrong received datatype: %s - Skipping",
221                     iBlock, evtData.fBlockCnt-1,
222                     evtData.fEventID, evtData.fEventID, 
223                     DataType2Text(inputDataType).c_str());
224           continue;
225         }
226       else {
227         HLTDebug("We get the right data type: Block # %i/%i; Event 0x%08LX (%Lu) Received datatype: %s; Block Size: %i",
228                  iBlock, evtData.fBlockCnt-1,
229                  evtData.fEventID, evtData.fEventID, 
230                  DataType2Text(inputDataType).c_str(),
231                  block.fSize);
232       }
233
234 #ifndef NDEBUG
235       unsigned long constBase;
236       double inputMultiplier;
237       GetOutputDataSize(constBase,inputMultiplier);
238       if(size<(constBase+block.fSize*inputMultiplier)){
239         HLTWarning("Memory Block given might be too small: %i < %f; Event %Lu", size, constBase+block.fSize*inputMultiplier, evtData.fEventID);
240       }
241 #endif
242
243       fESD->Reset();
244       //fESD->SetMagneticField(GetBz());
245
246       AliHLTTRDUtils::ReadClusters(fClusterArray, block.fPtr, block.fSize, &fNtimeBins);
247       HLTDebug("Reading number of time bins from input block. Setting number of timebins to %d", fNtimeBins);
248       AliTRDtrackerV1::SetNTimeBins(fNtimeBins);
249
250       HLTDebug("TClonesArray of clusters: nbEntries = %i", fClusterArray->GetEntriesFast());
251       fTracker->LoadClusters(fClusterArray);
252
253       fTracker->Clusters2Tracks(fESD);
254
255       Int_t nTracks = fESD->GetNumberOfTracks();
256       HLTInfo("Number of tracks  == %d ==", nTracks);  
257
258       TClonesArray* trdTracks;
259       trdTracks = fTracker->GetListOfTracks();
260       
261       if(fHighLevelOutput){
262         if(fEmulateHLTTracks && trdTracks){
263           // TClonesArray* oldArr = trdTracks;
264           trdTracks = new TClonesArray(*trdTracks);
265           AliHLTTRDUtils::EmulateHLTTracks(trdTracks);
266           // if(oldArr->At(0)){
267           //   HLTInfo("Old Track:");
268           //   ((AliTRDtrackV1*)oldArr->At(0))->Print("a");
269           //   HLTInfo("\nNew Track:");
270           //   ((AliTRDtrackV1*)trdTracks->At(0))->Print("a");
271           // }
272         }
273
274         TObjString strg;
275         strg.String() += fNtimeBins;
276         if(trdTracks)
277           PushBack(trdTracks, AliHLTTRDDefinitions::fgkHiLvlTracksDataType, block.fSpecification);
278         else{
279           TClonesArray temp("AliTRDtrackV1");
280           PushBack(&temp, AliHLTTRDDefinitions::fgkHiLvlTracksDataType, block.fSpecification);
281         }
282         PushBack(&strg, AliHLTTRDDefinitions::fgkHiLvlTracksDataType, block.fSpecification);
283
284         if(fEmulateHLTTracks && trdTracks){
285           trdTracks->Delete();
286           delete trdTracks;
287         }
288       }
289       else if(nTracks>0){
290         HLTDebug("We have an output ESDEvent: 0x%x with %i tracks", fESD, nTracks);
291         AliHLTUInt32_t addedSize = AliHLTTRDUtils::AddESDToOutput(fESD, outputPtr+offset);
292         totalSize += addedSize;
293           
294         // Fill block 
295         AliHLTComponentBlockData bd;
296         FillBlockData( bd );
297         //bd.fPtr = outputPtr;
298         bd.fOffset = offset;
299         bd.fSize = addedSize;
300         bd.fSpecification = block.fSpecification;
301         bd.fDataType = kAliHLTDataTypeTrack | kAliHLTDataOriginTRD;
302         outputBlocks.push_back( bd );
303         HLTDebug("BD ptr 0x%x, offset %i, size %i, datav1Type %s, spec 0x%x ", bd.fPtr, bd.fOffset, bd.fSize, DataType2Text(bd.fDataType).c_str(), bd.fSpecification);
304         offset = totalSize;
305
306         if (fOutputV1Tracks && trdTracks){
307           HLTDebug("We have an output array: pointer to trdTracks = 0x%x, nbEntries = %i", trdTracks, trdTracks->GetEntriesFast());
308           
309           addedSize = AliHLTTRDUtils::AddTracksToOutput(trdTracks, outputPtr+offset, fNtimeBins);
310           totalSize += addedSize;
311           
312           // Fill block 
313           FillBlockData( bd );
314           //bd.fPtr = outputPtr;
315           bd.fOffset = offset;
316           bd.fSize = addedSize;
317           bd.fSpecification = block.fSpecification;
318           bd.fDataType = AliHLTTRDDefinitions::fgkTracksDataType;
319           outputBlocks.push_back( bd );
320           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(), bd.fSpecification);
321           offset = totalSize;
322         }
323       }
324
325       HLTDebug("totalSize: %i", totalSize);
326       
327 //       if ( totalSize > allocSize )
328 //      {
329 //        HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
330 //        totalSize, size );
331 //        return EMSGSIZE;
332 //      }
333
334       //here we are deleting clusters (but not the TClonesArray itself)
335       fTracker->UnloadClusters();
336       AliTRDReconstructor::SetClusters(0x0);
337       fClusterArray->Delete();
338       
339     }
340       
341   size = totalSize;
342   HLTDebug("Event is done. size written to the output is %i", size);
343   return 0;
344 }
345
346 int AliHLTTRDTrackerV1Component::Configure(const char* arguments){
347   int iResult=0;
348   if (!arguments) return iResult;
349   
350   TString allArgs=arguments;
351   TString argument;
352   int bMissingParam=0;
353
354   TObjArray* pTokens=allArgs.Tokenize(" ");
355   if (pTokens) {
356     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
357       argument=((TObjString*)pTokens->At(i))->GetString();
358       if (argument.IsNull()) continue;
359       
360       if (argument.CompareTo("output_percentage")==0) {
361         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
362         HLTInfo("Setting output percentage to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
363         fOutputPercentage=((TObjString*)pTokens->At(i))->GetString().Atoi();
364         continue;
365       } 
366       else if (argument.CompareTo("-solenoidBz")==0) {
367         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
368         HLTWarning("argument -solenoidBz is deprecated, magnetic field set up globally (%f)", GetBz());
369         continue;
370       } 
371       else if (argument.CompareTo("-geometry")==0) {
372         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
373         HLTInfo("Setting geometry to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
374         fgeometryFileName=((TObjString*)pTokens->At(i))->GetString();
375         continue;
376       } 
377       else if (argument.CompareTo("-lowflux")==0) {
378         fRecoParamType = 0;
379         HLTInfo("Low flux reconstruction selected");
380         continue;
381       }
382       else if (argument.CompareTo("-highflux")==0) {
383         fRecoParamType = 1;
384         HLTInfo("High flux reconstruction selected");
385         continue;
386       }
387       else if (argument.CompareTo("-cosmics")==0) {
388         fRecoParamType = 2;
389         HLTInfo("Cosmics reconstruction selected");
390         continue;
391       }
392       else if (argument.CompareTo("-HLTflag")==0) {
393         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
394         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
395         if (toCompareTo.CompareTo("yes")==0){
396           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
397           fHLTflag=kTRUE;
398         }
399         else if (toCompareTo.CompareTo("no")==0){
400           HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
401           fHLTflag=kFALSE;
402         }
403         else {
404           HLTError("unknown argument for HLTflag: %s", toCompareTo.Data());
405           iResult=-EINVAL;
406           break;
407         }
408         continue;
409       }
410       else if (argument.CompareTo("-outputV1Tracks")==0) {
411         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
412         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
413         if (toCompareTo.CompareTo("yes")==0){
414           HLTInfo("Setting OutputV1Tracks to: %s", toCompareTo.Data());
415           fOutputV1Tracks=kTRUE;
416         }
417         else if (toCompareTo.CompareTo("no")==0){
418           HLTInfo("Setting OutputV1Tracks to: %s", toCompareTo.Data());
419           fOutputV1Tracks=kFALSE;
420         }
421         else {
422           HLTError("unknown argument for OutputV1Tracks: %s", toCompareTo.Data());
423           iResult=-EINVAL;
424           break;
425         }
426         continue;
427       }
428       else if (argument.CompareTo("-highLevelOutput")==0) {
429         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
430         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
431         if (toCompareTo.CompareTo("yes")==0){
432           HLTWarning("Setting highLevelOutput to: %s", toCompareTo.Data());
433           fHighLevelOutput=kTRUE;
434         }
435         else if (toCompareTo.CompareTo("no")==0){
436           HLTInfo("Setting highLevelOutput to: %s", toCompareTo.Data());
437           fHighLevelOutput=kFALSE;
438         }
439         else {
440           HLTError("unknown argument for highLevelOutput: %s", toCompareTo.Data());
441           iResult=-EINVAL;
442           break;
443         }
444         continue;
445       }
446       else if (argument.CompareTo("-emulateHLToutput")==0) {
447         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
448         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
449         if (toCompareTo.CompareTo("yes")==0){
450           HLTWarning("Setting emulateHLToutput to: %s", toCompareTo.Data());
451           fEmulateHLTTracks=kTRUE;
452         }
453         else if (toCompareTo.CompareTo("no")==0){
454           HLTInfo("Setting emulateHLToutput to: %s", toCompareTo.Data());
455           fEmulateHLTTracks=kFALSE;
456         }
457         else {
458           HLTError("unknown argument for emulateHLToutput: %s", toCompareTo.Data());
459           iResult=-EINVAL;
460           break;
461         }
462         continue;
463       }
464       else if (argument.CompareTo("-PIDmethod")==0) {
465         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
466         TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
467         if (toCompareTo.CompareTo("LH")==0){
468           HLTInfo("Setting PID method to: %s", toCompareTo.Data());
469           fPIDmethod=0;
470         }
471         else if (toCompareTo.CompareTo("NN")==0){
472           HLTInfo("Setting PID method to: %s", toCompareTo.Data());
473           fPIDmethod=1;
474         }
475         else if (toCompareTo.CompareTo("TM")==0){
476           HLTInfo("Setting PID method to: %s", toCompareTo.Data());
477           fPIDmethod=2;
478         }
479         else {
480           HLTError("unknown argument for PID method: %s", toCompareTo.Data());
481           iResult=-EINVAL;
482           break;
483         }
484         continue;
485       } 
486       
487       else {
488         HLTError("unknown argument: %s", argument.Data());
489         iResult=-EINVAL;
490         break;
491       }
492     }
493     delete pTokens;
494   }
495   if (bMissingParam) {
496     HLTError("missing parameter for argument %s", argument.Data());
497     iResult=-EINVAL;
498   }
499   if(iResult>=0){
500     iResult=SetParams();
501   }
502   return iResult;
503 }
504
505 int AliHLTTRDTrackerV1Component::SetParams()
506 {
507   Int_t iResult=0;
508   if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
509     HLTError("DefaultStorage is not set in CDBManager");
510     return -EINVAL;
511   }
512   if(AliCDBManager::Instance()->GetRun()<0){
513     HLTError("Run Number is not set in CDBManager");
514     return -EINVAL;
515   }
516   HLTInfo("CDB default storage: %s; RunNo: %i", (AliCDBManager::Instance()->GetDefaultStorage()->GetBaseFolder()).Data(), AliCDBManager::Instance()->GetRun());
517
518   if(!AliGeomManager::GetGeometry()){
519     if(fgeometryFileName.CompareTo("")==0 || !TFile::Open(fgeometryFileName.Data())){
520       HLTInfo("Loading standard geometry file");
521       AliGeomManager::LoadGeometry();
522     }else{
523       HLTWarning("Loading NON-standard geometry file");
524       AliGeomManager::LoadGeometry(fgeometryFileName.Data());
525     }
526     if(!AliGeomManager::GetGeometry()){
527       HLTError("Could not load geometry");
528       return -EINVAL;
529     }
530     HLTInfo("Applying Alignment from CDB object");
531     AliGeomManager::ApplyAlignObjsFromCDB("TRD");
532   }
533   else{
534     HLTInfo("Geometry Already Loaded!");
535   }
536   
537   if(fReconstructor->GetRecoParam()){
538     fRecoParam = new AliTRDrecoParam(*fReconstructor->GetRecoParam());
539     HLTInfo("RecoParam already set!");
540   }else{
541     if(fRecoParamType == 0){
542 #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
543       if(fHLTflag){
544         HLTInfo("Low flux HLT params init.");
545         fRecoParam = AliTRDrecoParam::GetLowFluxHLTParam();
546       }else
547 #endif
548         {
549           HLTInfo("Low flux params init.");
550           fRecoParam = AliTRDrecoParam::GetLowFluxParam();
551         }
552     }
553     if(fRecoParamType == 1){
554 #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
555       if(fHLTflag){
556         HLTInfo("High flux HLT params init.");
557         fRecoParam = AliTRDrecoParam::GetHighFluxHLTParam();
558       }else
559 #endif
560         {
561           HLTInfo("High flux params init.");
562           fRecoParam = AliTRDrecoParam::GetHighFluxParam();
563         }
564     }
565     if(fRecoParamType == 2){
566       HLTInfo("Cosmic Test params init.");
567       fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
568     }
569   }
570
571   if(!fRecoParam)
572     {
573       HLTError("No reco params initialized. Sniffing big trouble!");
574       return -EINVAL;
575     }
576
577   switch(fPIDmethod){
578   case 0: fRecoParam->SetPIDNeuralNetwork(kFALSE); break;
579   case 1: fRecoParam->SetPIDNeuralNetwork(kTRUE); break;
580   case 2: fRecoParam->SetPIDNeuralNetwork(kFALSE); break;
581   }
582
583   fRecoParam->SetImproveTracklets(fImproveTracklets);
584
585   fRecoParam->SetStreamLevel(AliTRDrecoParam::kTracker, 0);
586   fReconstructor->SetRecoParam(fRecoParam);
587
588   TString recoOptions="sa,!cw";
589   
590   if(fHLTflag)
591     recoOptions += ",hlt";
592
593   HLTDebug("Reconstructor options are: %s",recoOptions.Data());
594   fReconstructor->SetOption(recoOptions.Data());
595
596   return iResult;
597 }
598
599 int AliHLTTRDTrackerV1Component::Reconfigure(const char* cdbEntry, const char* chainId)
600 {
601   // see header file for class documentation
602
603   int iResult=0;
604   const char* path="HLT/ConfigTRD/TrackerV1Component";
605   const char* defaultNotify="";
606   if (cdbEntry) {
607     path=cdbEntry;
608     defaultNotify=" (default)";
609   }
610   if (path) {
611     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
612     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
613     if (pEntry) {
614       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
615       if (pString) {
616         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
617         iResult=Configure(pString->GetString().Data());
618       } else {
619         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
620       }
621     } else {
622       HLTError("cannot fetch object \"%s\" from CDB", path);
623     }
624   }
625
626   return iResult;
627
628 }
629
630 int AliHLTTRDTrackerV1Component::ReadPreprocessorValues(const char* modules)
631 {
632   // see header file for class documentation
633   
634   int iResult = 0;
635   TString str(modules);
636   if(str.Contains("HLT") || str.Contains("TRD") || str.Contains("GRP")){
637   
638   }  
639   return iResult;
640 }
641
642 void AliHLTTRDTrackerV1Component::GetOCDBObjectDescription(TMap* const targetMap){
643   // Get a list of OCDB object description needed for the particular component
644   if (!targetMap) return;
645   targetMap->Add(new TObjString("HLT/ConfigTRD/TrackerV1Component"), new TObjString("component arguments"));
646   targetMap->Add(new TObjString("TRD/Calib/ChamberGainFactor"), new TObjString("gain factor of chambers"));
647   targetMap->Add(new TObjString("TRD/Calib/ChamberT0"), new TObjString("T0 of chambers"));
648   targetMap->Add(new TObjString("TRD/Calib/ChamberVdrift"), new TObjString("drift velocity of chambers"));
649   targetMap->Add(new TObjString("TRD/Calib/DetNoise"), new TObjString("noise of chambers"));
650   targetMap->Add(new TObjString("TRD/Calib/LocalGainFactor"), new TObjString("per pad gain factor"));
651   targetMap->Add(new TObjString("TRD/Calib/LocalT0"), new TObjString("per pad T0"));
652   targetMap->Add(new TObjString("TRD/Calib/LocalVdrift"), new TObjString("per pad drift velocity"));
653   targetMap->Add(new TObjString("TRD/Calib/PadNoise"), new TObjString("per pad noise"));
654   targetMap->Add(new TObjString("TRD/Calib/PadStatus"), new TObjString("pad status"));
655   targetMap->Add(new TObjString("TRD/Calib/PRFWidth"), new TObjString("pad response function"));
656   targetMap->Add(new TObjString("TRD/Calib/ChamberStatus"), new TObjString("status of chambers"));
657   targetMap->Add(new TObjString("TRD/Calib/PIDLQ"), new TObjString("likelyhood PID"));
658   targetMap->Add(new TObjString("TRD/Calib/PIDNN"), new TObjString("neuronal network PID"));
659   targetMap->Add(new TObjString("TRD/Calib/PIDThresholds"), new TObjString("threshold for PID"));
660 }