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