]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
Changes done by Konstatin Antipin in order to speed up TRD code
[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
41 #include "AliGeomManager.h"
42 #include "AliTRDReconstructor.h"
43 #include "AliCDBManager.h"
44 #include "AliTRDclusterizerHLT.h"
45 #include "AliTRDrecoParam.h"
46 #include "AliTRDrawStreamBase.h"
47
48 #include "AliRawReaderMemory.h"
49
50 #include <cstdlib>
51 #include <cerrno>
52 #include <string>
53
54 // this is a global object used for automatic component registration, do not use this
55 AliHLTTRDClusterizerComponent gAliHLTTRDClusterizerComponent;
56
57 ClassImp(AliHLTTRDClusterizerComponent);
58    
59 AliHLTTRDClusterizerComponent::AliHLTTRDClusterizerComponent()
60   : AliHLTProcessor()
61   , fOutputPercentage(100) // By default we copy to the output exactly what we got as input  
62   , fStrorageDBpath("local://$ALICE_ROOT")
63   , fClusterizer(NULL)
64   , fRecoParam(NULL)
65   , fCDB(NULL)
66   , fMemReader(NULL)
67   , fGeometryFileName("")
68   , fGeometryFile(NULL)
69   , fReconstructor(NULL)
70 {
71   // Default constructor
72
73   fGeometryFileName = getenv("ALICE_ROOT");
74   fGeometryFileName += "/HLT/TRD/geometry.root";
75 }
76
77 AliHLTTRDClusterizerComponent::~AliHLTTRDClusterizerComponent()
78 {
79   // Destructor
80   // Work is Done in DoDeInit()
81 }
82 const char* AliHLTTRDClusterizerComponent::GetComponentID()
83 {
84   // Return the component ID const char *
85   return "TRDClusterizer"; // The ID of this component
86 }
87
88 void AliHLTTRDClusterizerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
89 {
90   // Get the list of input data
91   list.clear(); // We do not have any requirements for our input data type(s).
92   list.push_back( AliHLTTRDDefinitions::fgkDDLRawDataType );
93 }
94
95 AliHLTComponent_DataType AliHLTTRDClusterizerComponent::GetOutputDataType()
96 {
97   // Get the output data type
98   return AliHLTTRDDefinitions::fgkClusterDataType;
99 }
100
101 void AliHLTTRDClusterizerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
102 {
103   // Get the output data size
104   constBase = 0;
105   inputMultiplier = ((double)fOutputPercentage)/100.0;
106 }
107
108 AliHLTComponent* AliHLTTRDClusterizerComponent::Spawn()
109 {
110   // Spawn function, return new instance of this class
111   return new AliHLTTRDClusterizerComponent;
112 };
113
114 int AliHLTTRDClusterizerComponent::DoInit( int argc, const char** argv )
115 {
116   // perform initialization. We check whether our relative output size is specified in the arguments.
117   fOutputPercentage = 100;
118   Int_t iRawDataVersion = 2;
119   int i = 0;
120   char* cpErr;
121   Bool_t bWriteClusters = kTRUE;
122   
123
124   Int_t iRecoParamType = -1; // default will be the low flux
125
126   // the data type will become obsolete as soon as the formats are established
127   Int_t iRecoDataType = -1; // default will be simulation
128   
129   while ( i < argc )
130     {
131       HLTDebug("argv[%d] == %s", i, argv[i] );
132       if ( !strcmp( argv[i], "output_percentage" ) )
133         {
134           if ( i+1>=argc )
135             {
136               HLTError("Missing output_percentage parameter");
137               return ENOTSUP;
138             }
139           HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
140           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
141           if ( *cpErr )
142             {
143               HLTError("Cannot convert output_percentage parameter '%s'", argv[i+1] );
144               return EINVAL;
145             }
146           HLTInfo("Output percentage set to %lu %%", fOutputPercentage );
147           i += 2;
148         }
149       else if ( strcmp( argv[i], "-cdb" ) == 0)
150         {
151           if ( i+1 >= argc )
152             {
153               HLTError("Missing -cdb argument");
154               return ENOTSUP;         
155             }
156           fStrorageDBpath = argv[i+1];
157           HLTInfo("DB storage is %s", fStrorageDBpath.c_str() );          
158           i += 2;
159           
160         }      
161
162       else if ( strcmp( argv[i], "-lowflux" ) == 0)
163         {
164           iRecoParamType = 0;     
165           HLTDebug("Low flux reco selected.");
166           i++;
167           
168         }      
169
170       else if ( strcmp( argv[i], "-highflux" ) == 0)
171         {
172           iRecoParamType = 1;     
173           HLTDebug("High flux reco selected.");
174           i++;
175           
176         }      
177
178       else if ( strcmp( argv[i], "-cosmics" ) == 0)
179         {
180           iRecoParamType = 2;     
181           HLTDebug("Cosmic test reco selected.");
182           i++;
183           
184         }      
185
186       // raw data type - sim or experiment
187       else if ( strcmp( argv[i], "-simulation" ) == 0)
188         {
189           iRecoDataType = 0;
190           i++;
191           
192         }
193
194       else if ( strcmp( argv[i], "-experiment" ) == 0)
195         {
196           iRecoDataType = 1;
197           i++;
198           
199         }
200
201       else if ( strcmp( argv[i], "-rawver" ) == 0)
202         {
203           if ( i+1 >= argc )
204             {
205               HLTError("Missing -rawver argument");
206               return ENOTSUP;         
207             }
208           iRawDataVersion = atoi( argv[i+1] );
209           HLTInfo("Raw data version is %d", iRawDataVersion );    
210           i += 2;
211           
212         }      
213
214       else if ( strcmp( argv[i], "-geometry" ) == 0)
215         {
216           if ( i+1 >= argc )
217             {
218               HLTError("Missing -geometry argument");
219               return ENOTSUP;         
220             }
221           fGeometryFileName = argv[i+1];
222           HLTInfo("GeomFile storage is %s", fGeometryFileName.c_str() );          
223           i += 2;
224           
225         }      
226       else{
227         HLTError("Unknown option '%s'", argv[i] );
228         return EINVAL;
229       }
230       
231     }
232
233   // THE "REAL" INIT COMES HERE
234
235   if (iRecoParamType < 0 || iRecoParamType > 2)
236     {
237       HLTWarning("No reco param selected. Use -lowflux or -highflux flag. Defaulting to low flux.");
238       iRecoParamType = 0;
239     }
240
241   if (iRecoParamType == 0)
242     {
243       fRecoParam = AliTRDrecoParam::GetLowFluxParam();
244       HLTDebug("Low flux params init.");
245     }
246
247   if (iRecoParamType == 1)
248     {
249       fRecoParam = AliTRDrecoParam::GetHighFluxParam();
250       HLTDebug("High flux params init.");
251     }
252
253   if (iRecoParamType == 2)
254     {
255       fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
256       HLTDebug("Cosmic Test params init.");
257     }
258
259   if (fRecoParam == 0)
260     {
261       HLTError("No reco params initialized. Sniffing big trouble!");
262       return -1;
263     }
264
265   fReconstructor = new AliTRDReconstructor();
266   fReconstructor->SetRecoParam(fRecoParam);
267   //fReconstructor->SetStreamLevel(0, AliTRDReconstructor::kClusterizer); // default value
268   fReconstructor->SetOption("!cw,sl_cf_0");
269
270   // init the raw data type to be used...
271   // the switch here will become obsolete as soon as the data structures is fixed 
272   // both: in sim and reality
273   if (iRecoDataType < 0 || iRecoDataType > 1)
274     {
275       HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
276       iRecoDataType = 0;
277     }
278
279   if (iRecoDataType == 0)
280     {
281       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDsimStream);
282       HLTDebug("Data type expected is SIMULATION!");
283     }
284
285   if (iRecoDataType == 1)
286     {
287       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDrealStream);
288       HLTDebug("Data type expected is EXPERIMENT!");
289     }
290
291   // the DATA BASE STUFF
292   fCDB = AliCDBManager::Instance();
293   if (!fCDB)
294     {
295       HLTError("Could not get CDB instance", "fCDB 0x%x", fCDB);
296     }
297   else
298     {
299       fCDB->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
300       fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
301       HLTDebug("CDB instance; fCDB 0x%x", fCDB);
302     }
303
304   fGeometryFile = TFile::Open(fGeometryFileName.c_str());
305   if (fGeometryFile)
306     {
307       AliGeomManager::LoadGeometry(fGeometryFileName.c_str());
308     }
309   else
310     {
311       HLTError("Unable to open file. FATAL!");
312       return -1;
313     }
314
315   fMemReader = new AliRawReaderMemory;
316
317   fClusterizer = new AliTRDclusterizerHLT("TRDCclusterizer", "TRDCclusterizer");
318   fClusterizer->SetReconstructor(fReconstructor);
319
320   fClusterizer->SetRawVersion(iRawDataVersion);
321   fClusterizer->InitClusterTree();
322   return 0;
323 }
324
325 int AliHLTTRDClusterizerComponent::DoDeinit()
326 {
327   // Deinitialization of the component
328   delete fMemReader;
329   fMemReader = 0;
330   delete fClusterizer;
331   fClusterizer = 0;
332   
333   fReconstructor->SetClusters(0x0);
334   delete fReconstructor;
335   fReconstructor = 0x0;
336   return 0;
337
338   if (fGeometryFile)
339     {
340       fGeometryFile->Close();
341       delete fGeometryFile;
342       fGeometryFile = 0;
343     }
344
345   if (fCDB)
346     {
347       HLTDebug("destroy fCDB");
348       fCDB->Destroy();
349       fCDB = 0;
350     }
351
352   if (fRecoParam)
353     {
354       HLTDebug("Deleting fRecoParam");
355       delete fRecoParam;
356       fRecoParam = 0;
357     }
358 }
359
360 int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponentEventData& evtData, 
361                                             const AliHLTComponentBlockData* blocks, 
362                                             AliHLTComponent_TriggerData& /*trigData*/, 
363                                             AliHLTUInt8_t* /*outputPtr*/, 
364                                             AliHLTUInt32_t& size, 
365                                             vector<AliHLTComponent_BlockData>& /*outputBlocks*/ )
366 {
367   // Process an event
368   HLTDebug( "NofBlocks %lu", evtData.fBlockCnt );
369   // Process an event
370
371   //implement a usage of the following
372   //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
373   //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
374   //   void *triggerData = trigData.fData;
375   //HLTDebug( "Trigger data received. Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
376
377   // Loop over all input blocks in the event
378   AliHLTComponentDataType expectedDataType = (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
379   for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
380     {
381       // lets not use the internal TRD data types here : AliHLTTRDDefinitions::fgkDDLRawDataType
382       // which is depreciated - we use HLT global defs instead
383       //      if ( blocks[i].fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD) )
384       AliHLTComponentDataType inputDataType = blocks[i].fDataType;
385       if ( inputDataType != expectedDataType)
386         {
387           HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s; Skipping",
388                     i, evtData.fBlockCnt,
389                     evtData.fEventID, evtData.fEventID, 
390                     DataType2Text(inputDataType).c_str(), 
391                     DataType2Text(expectedDataType).c_str());
392           continue;
393         }
394       else 
395         {
396           HLTDebug("We get the right data type!");
397         }
398       
399       //      fMemReader->Reset();
400       fMemReader->SetMemory((UChar_t*) blocks[i].fPtr, blocks[i].fSize);
401
402       AliHLTUInt32_t spec = blocks[i].fSpecification;
403       
404       Int_t id = 1024;
405       
406       for ( Int_t ii = 0; ii < 18 ; ii++ ) {
407         if ( spec & 0x00000001 ) {
408           id += ii;
409           break;
410         }
411         spec = spec >> 1 ;
412       }
413       
414       fMemReader->SetEquipmentID( id ); 
415
416       fClusterizer->ResetTree();  
417       Bool_t iclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
418       if (iclustered == kTRUE)
419         {
420           HLTDebug( "Clustered successfully");
421         }
422       else
423         {
424           HLTError("Clustering ERROR");
425           return -1;
426         }
427
428       // put the tree into output blocks
429       //fcTree->Print();
430       TClonesArray *clusterArray = fClusterizer->RecPoints();
431       fClusterizer->SetClustersOwner(kFALSE);
432         
433       if (clusterArray){
434         HLTDebug("clusterArray: Entries - %i; Size - %i",clusterArray->GetEntriesFast(),sizeof(clusterArray));
435         PushBack(clusterArray, AliHLTTRDDefinitions::fgkClusterDataType, blocks[i].fSpecification);
436       }
437       else 
438         HLTWarning("Array of clusters is empty!");
439         
440       if (clusterArray){
441         clusterArray->Delete();
442         delete clusterArray;
443       }
444       
445
446     }//  for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
447
448   
449   
450   
451   size=0; // this function did not write data to the buffer directly
452   return 0;
453 }