]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
Using common HLT track data format for TRD tracks (Theodor)
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDClusterizerComponent.cxx
1 // $Id$
2
3 /**************************************************************************
4  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Authors: Matthias Richter <Matthias.Richter@ift.uib.no>                *
7  *          Timm Steinbeck <timm@kip.uni-heidelberg.de>                   *
8  *          for The ALICE Off-line 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 Timm Steinbeck, Matthias Richter
21     @date   
22     @brief  A TRDClusterizer processing component for the HLT. */
23
24 // see header file for class documentation                                   //
25 // or                                                                        //
26 // refer to README to build package                                          //
27 // or                                                                        //
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt                          //
29
30 #if __GNUC__ >= 3
31 using namespace std;
32 #endif
33
34 #include "TTree.h"
35 #include "TFile.h"
36 #include "TBranch.h"
37
38 #include "AliHLTTRDClusterizerComponent.h"
39 #include "AliHLTTRDDefinitions.h"
40 #include "AliHLTTRDCluster.h"
41
42 #include "AliGeomManager.h"
43 #include "AliTRDReconstructor.h"
44 #include "AliCDBManager.h"
45 #include "AliHLTTRDClusterizer.h"
46 #include "AliTRDrecoParam.h"
47 #include "AliTRDrawStreamBase.h"
48 #include "AliTRDcluster.h"
49
50 #include "AliRawReaderMemory.h"
51
52 #ifdef HAVE_VALGRIND_CALLGRIND_H
53 #include <valgrind/callgrind.h>
54 #else
55 #define CALLGRIND_START_INSTRUMENTATION do { } while (0)
56 #define CALLGRIND_STOP_INSTRUMENTATION do { } while (0)
57 #endif
58
59 #include <cstdlib>
60 #include <cerrno>
61 #include <string>
62
63 ClassImp(AliHLTTRDClusterizerComponent)
64    
65 AliHLTTRDClusterizerComponent::AliHLTTRDClusterizerComponent():
66   AliHLTProcessor(),
67   fOutputPercentage(100), // By default we copy to the output exactly what we got as input
68   fStrorageDBpath("local://$ALICE_ROOT/OCDB"),
69   fClusterizer(NULL),
70   fRecoParam(NULL),
71   fCDB(NULL),
72   fMemReader(NULL),
73   fReconstructor(NULL),
74   fGeometryFileName("")
75 {
76   // Default constructor
77
78   fGeometryFileName = getenv("ALICE_ROOT");
79   fGeometryFileName += "/HLT/TRD/geometry.root";
80 }
81
82 AliHLTTRDClusterizerComponent::~AliHLTTRDClusterizerComponent()
83 {
84   // Destructor
85   // Work is Done in DoDeInit()
86 }
87
88
89 const char* AliHLTTRDClusterizerComponent::GetComponentID()
90 {
91   // Return the component ID const char *
92   return "TRDClusterizer"; // The ID of this component
93 }
94
95 void AliHLTTRDClusterizerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
96 {
97   // Get the list of input data
98   list.clear(); // We do not have any requirements for our input data type(s).
99   list.push_back( AliHLTTRDDefinitions::fgkDDLRawDataType );
100 }
101
102 AliHLTComponent_DataType AliHLTTRDClusterizerComponent::GetOutputDataType()
103 {
104   // Get the output data type
105   return AliHLTTRDDefinitions::fgkClusterDataType;
106 }
107
108 void AliHLTTRDClusterizerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
109 {
110   // Get the output data size
111   constBase = 0;
112   inputMultiplier = ((double)fOutputPercentage)/100.0;
113 }
114
115 AliHLTComponent* AliHLTTRDClusterizerComponent::Spawn()
116 {
117   // Spawn function, return new instance of this class
118   return new AliHLTTRDClusterizerComponent;
119 };
120
121 int AliHLTTRDClusterizerComponent::DoInit( int argc, const char** argv )
122 {
123   // perform initialization. We check whether our relative output size is specified in the arguments.
124   fOutputPercentage = 100;
125   Int_t iRawDataVersion = 2;
126   int i = 0;
127   char* cpErr;
128
129   Int_t iRecoParamType = -1; // default will be the low flux
130
131   // the data type will become obsolete as soon as the formats are established
132   Int_t iRecoDataType = -1; // default will be simulation
133   Int_t iyPosMethod = 1;     // 0=COG 1=LUT 2=Gauss 
134   
135   while ( i < argc )
136     {
137       HLTDebug("argv[%d] == %s", i, argv[i] );
138       if ( !strcmp( argv[i], "output_percentage" ) )
139         {
140           if ( i+1>=argc )
141             {
142               HLTError("Missing output_percentage parameter");
143               return ENOTSUP;
144             }
145           HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
146           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
147           if ( *cpErr )
148             {
149               HLTError("Cannot convert output_percentage parameter '%s'", argv[i+1] );
150               return EINVAL;
151             }
152           HLTInfo("Output percentage set to %lu %%", fOutputPercentage );
153           i += 2;
154         }
155       else if ( strcmp( argv[i], "-cdb" ) == 0)
156         {
157           if ( i+1 >= argc )
158             {
159               HLTError("Missing -cdb argument");
160               return ENOTSUP;         
161             }
162           fStrorageDBpath = argv[i+1];
163           HLTInfo("DB storage is %s", fStrorageDBpath.c_str() );          
164           i += 2;
165         }
166       else if ( strcmp( argv[i], "-lowflux" ) == 0)
167         {
168           iRecoParamType = 0;     
169           HLTDebug("Low flux reco selected.");
170           i++;
171         }
172       else if ( strcmp( argv[i], "-highflux" ) == 0)
173         {
174           iRecoParamType = 1;     
175           HLTDebug("High flux reco selected.");
176           i++;
177         }
178       else if ( strcmp( argv[i], "-cosmics" ) == 0)
179         {
180           iRecoParamType = 2;     
181           HLTDebug("Cosmic test reco selected.");
182           i++;
183         }
184       // raw data type - sim or experiment
185       else if ( strcmp( argv[i], "-simulation" ) == 0)
186         {
187           iRecoDataType = 0;
188           i++;
189         }
190       else if ( strcmp( argv[i], "-experiment" ) == 0)
191         {
192           iRecoDataType = 1;
193           i++;
194         }
195       else if ( strcmp( argv[i], "-rawver" ) == 0)
196         {
197           if ( i+1 >= argc )
198             {
199               HLTError("Missing -rawver argument");
200               return ENOTSUP;         
201             }
202           iRawDataVersion = atoi( argv[i+1] );
203           HLTInfo("Raw data version is %d", iRawDataVersion );    
204           i += 2;
205           
206         }      
207
208       else if ( strcmp( argv[i], "-geometry" ) == 0)
209         {
210           if ( i+1 >= argc )
211             {
212               HLTError("Missing -geometry argument");
213               return ENOTSUP;         
214             }
215           fGeometryFileName = argv[i+1];
216           HLTInfo("GeomFile storage is %s", fGeometryFileName.c_str() );          
217           i += 2;
218         }
219
220       else if ( strcmp( argv[i], "-yPosMethod" ) == 0)
221         {
222           if ( i+1 >= argc )
223             {
224               HLTError("Missing -yPosMethod argument");
225               return ENOTSUP;         
226             }
227           if( strcmp(argv[i], "COG") )
228             iyPosMethod=0;
229           else if( strcmp(argv[i], "LUT") )
230             iyPosMethod=1;
231           else if( strcmp(argv[i], "Gauss") )
232             iyPosMethod=2;
233           else {
234             HLTError("Unknown -yPosMethod argument");
235             return ENOTSUP;           
236           }
237           i += 2;
238         }
239       
240       else{
241         HLTError("Unknown option '%s'", argv[i] );
242         return EINVAL;
243       }
244       
245     }
246
247   // THE "REAL" INIT COMES HERE
248
249   if (iRecoParamType < 0 || iRecoParamType > 2)
250     {
251       HLTWarning("No reco param selected. Use -lowflux or -highflux flag. Defaulting to low flux.");
252       iRecoParamType = 0;
253     }
254
255   if (iRecoParamType == 0)
256     {
257       fRecoParam = AliTRDrecoParam::GetLowFluxParam();
258       HLTDebug("Low flux params init.");
259     }
260
261   if (iRecoParamType == 1)
262     {
263       fRecoParam = AliTRDrecoParam::GetHighFluxParam();
264       HLTDebug("High flux params init.");
265     }
266
267   if (iRecoParamType == 2)
268     {
269       fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
270       HLTDebug("Cosmic Test params init.");
271     }
272
273   if (fRecoParam == 0)
274     {
275       HLTError("No reco params initialized. Sniffing big trouble!");
276       return -1;
277     }
278
279   fReconstructor = new AliTRDReconstructor();
280   fReconstructor->SetRecoParam(fRecoParam);
281   fReconstructor->SetStreamLevel(0, AliTRDReconstructor::kClusterizer); // default value
282   HLTInfo("Not writing clusters. I.e. output is a TClonesArray of clusters");
283   TString recoOptions="hlt,!cw,sl_cf_0";
284   switch(iRecoDataType){
285   case 0: recoOptions += ",tc"; break;
286   case 1: recoOptions += ",!tc"; break;
287   }
288   switch(iyPosMethod){
289   case 0: recoOptions += ",!gs,!lut"; break;
290   case 1: recoOptions += ",!gs,lut"; break;
291   case 2: recoOptions += ",gs,!lut"; break;
292   }
293   
294   fReconstructor->SetOption(recoOptions.Data());
295   
296   // init the raw data type to be used...
297   // the switch here will become obsolete as soon as the data structures is fixed 
298   // both: in sim and reality
299   if (iRecoDataType < 0 || iRecoDataType > 1)
300     {
301       HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
302       iRecoDataType = 0;
303     }
304
305   if (iRecoDataType == 0)
306     {
307       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDsimStream);
308       HLTDebug("Data type expected is SIMULATION!");
309     }
310
311   if (iRecoDataType == 1)
312     {
313       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDrealStream);
314       HLTDebug("Data type expected is EXPERIMENT!");
315     }
316
317   // the DATA BASE STUFF
318   fCDB = AliCDBManager::Instance();
319   if (!fCDB)
320     {
321       HLTError("Could not get CDB instance", "fCDB 0x%x", fCDB);
322     }
323   else
324     {
325       fCDB->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
326       fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
327       HLTDebug("CDB instance; fCDB 0x%x", fCDB);
328     }
329
330   if((AliGeomManager::GetGeometry()) == NULL){
331     
332     if ( TFile::Open(fGeometryFileName.c_str())) {
333       AliGeomManager::LoadGeometry(fGeometryFileName.c_str());
334     }
335     else {
336       HLTError("Cannot load geometry from file %s",fGeometryFileName.c_str());
337       return EINVAL;
338     }
339   }
340   else
341     HLTInfo("Geometry Already Loaded");
342   
343   fMemReader = new AliRawReaderMemory;
344
345   fClusterizer = new AliHLTTRDClusterizer("TRDCclusterizer", "TRDCclusterizer");
346   fClusterizer->SetReconstructor(fReconstructor);
347   fClusterizer->SetUseLabels(kFALSE);
348   fClusterizer->SetRawVersion(iRawDataVersion);
349   return 0;
350 }
351
352 int AliHLTTRDClusterizerComponent::DoDeinit()
353 {
354   // Deinitialization of the component
355   delete fMemReader;
356   fMemReader = 0;
357   delete fClusterizer;
358   fClusterizer = 0;
359   
360   fReconstructor->SetClusters(0x0);
361   delete fReconstructor;
362   fReconstructor = 0x0;
363   return 0;
364
365
366   if (fCDB)
367     {
368       HLTDebug("destroy fCDB");
369       fCDB->Destroy();
370       fCDB = 0;
371     }
372
373   if (fRecoParam)
374     {
375       HLTDebug("Deleting fRecoParam");
376       delete fRecoParam;
377       fRecoParam = 0;
378     }
379 }
380
381 int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponentEventData& evtData, 
382                                             const AliHLTComponentBlockData* blocks, 
383                                             AliHLTComponent_TriggerData& /*trigData*/, 
384                                             AliHLTUInt8_t* outputPtr, 
385                                             AliHLTUInt32_t& size, 
386                                             vector<AliHLTComponent_BlockData>& outputBlocks )
387 {
388   // Process an event
389   HLTDebug( "NofBlocks %lu", evtData.fBlockCnt );
390   // Process an event
391   AliHLTUInt32_t totalSize = 0, offset = 0;
392
393   //implement a usage of the following
394   //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
395   //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
396   //   void *triggerData = trigData.fData;
397   //HLTDebug( "Trigger data received. Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
398
399   // Loop over all input blocks in the event
400   AliHLTComponentDataType expectedDataType = (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
401   for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
402     {
403       if (evtData.fEventID == 1)
404         CALLGRIND_START_INSTRUMENTATION;
405       
406       const AliHLTComponentBlockData &block = blocks[i];
407       // lets not use the internal TRD data types here : AliHLTTRDDefinitions::fgkDDLRawDataType
408       // which is depreciated - we use HLT global defs instead
409       //      if ( block.fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD) )
410       AliHLTComponentDataType inputDataType = block.fDataType;
411       if ( inputDataType != expectedDataType)
412         {
413           HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s; Skipping",
414                     i, evtData.fBlockCnt,
415                     evtData.fEventID, evtData.fEventID, 
416                     DataType2Text(inputDataType).c_str(), 
417                     DataType2Text(expectedDataType).c_str());
418           continue;
419         }
420       else 
421         {
422         HLTDebug("We get the right data type: Block # %i/%i; Event 0x%08LX (%Lu) Received datatype: %s",
423                     i, evtData.fBlockCnt,
424                     evtData.fEventID, evtData.fEventID, 
425                     DataType2Text(inputDataType).c_str());
426         }
427       
428       //      fMemReader->Reset();
429       fMemReader->SetMemory((UChar_t*) block.fPtr, block.fSize);
430
431       AliHLTUInt32_t spec = block.fSpecification;
432       
433       Int_t id = 1024;
434       
435       for ( Int_t ii = 0; ii < 18 ; ii++ ) {
436         if ( spec & 0x00000001 ) {
437           id += ii;
438           break;
439         }
440         spec = spec >> 1 ;
441       }
442
443       fMemReader->SetEquipmentID( id ); 
444       
445       fClusterizer->SetMemBlock((AliHLTUInt8_t*)(outputPtr+offset));
446       Bool_t iclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
447       if (iclustered == kTRUE)
448         {
449           HLTDebug( "Clustered successfully");
450         }
451       else
452         {
453           HLTError("Clustering ERROR");
454           return -1;
455         }
456
457       // put the tree into output
458       //fcTree->Print();
459       
460       AliHLTUInt32_t addedSize = fClusterizer->GetAddedSize();
461         if (addedSize > 0){
462           // Using low-level interface 
463           // with interface classes
464           totalSize += addedSize;
465           if ( totalSize > size )
466             {
467               HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
468                        totalSize, size );
469               return EMSGSIZE;
470             }
471                 
472           // Fill block 
473           AliHLTComponentBlockData bd;
474           FillBlockData( bd );
475           bd.fOffset = offset;
476           bd.fSize = addedSize;
477           //bd.fSpecification = spec;
478           bd.fSpecification = gkAliEventTypeData;
479           bd.fDataType = AliHLTTRDDefinitions::fgkClusterDataType;
480           outputBlocks.push_back( bd );
481           HLTDebug( "BD fPtr 0x%x, fOffset %i, size %i, dataType %s, spec 0x%x ", bd.fPtr, bd.fOffset, bd.fSize, DataType2Text(bd.fDataType).c_str(), spec);
482           offset = totalSize;
483               
484         }
485         else 
486           HLTWarning("Array of clusters is empty!");
487     }
488   fReconstructor->SetClusters(0x0);
489
490   size = totalSize;
491   HLTDebug("Event is done. size written to the output is %i", size);
492   return 0;
493 }
494
495 void AliHLTTRDClusterizerComponent::PrintObject( TClonesArray* inClustersArray)
496 {
497   AliTRDcluster* cluster=0x0;
498   
499   for (Int_t i=0; i < inClustersArray->GetEntriesFast(); i++){
500     cluster = dynamic_cast<AliTRDcluster*>(inClustersArray->At(i));
501     HLTDebug("cluster[%i]",i);
502     HLTDebug("  PadCol = %i; PadRow = %i; PadTime = %i", cluster->GetPadCol(), cluster->GetPadRow(), cluster->GetPadTime());
503     HLTDebug("  Detector = %i, Amplitude = %f, Center = %f", cluster->GetDetector(), cluster->GetQ(), cluster->GetCenter());
504     HLTDebug("  LocalTimeBin =  %i; NPads = %i; maskedPosition: %s, status: %s", cluster->GetLocalTimeBin(), cluster->GetNPads(),cluster->GetPadMaskedPosition(),cluster->GetPadMaskedPosition());
505   }
506   
507 }
508