]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
correcting compilation warning and making change of AliShuttleInterface (rev 29388)
[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   , fReconstructor(NULL)
68   , fGeometryFileName("")
69   , fGeometryFile(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 = kFALSE;
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 if ( strcmp( argv[i], "-writeClusters" ) == 0)
227         {
228           bWriteClusters = kTRUE;
229           i++;
230         }
231       else{
232         HLTError("Unknown option '%s'", argv[i] );
233         return EINVAL;
234       }
235       
236     }
237
238   // THE "REAL" INIT COMES HERE
239
240   if (iRecoParamType < 0 || iRecoParamType > 2)
241     {
242       HLTWarning("No reco param selected. Use -lowflux or -highflux flag. Defaulting to low flux.");
243       iRecoParamType = 0;
244     }
245
246   if (iRecoParamType == 0)
247     {
248       fRecoParam = AliTRDrecoParam::GetLowFluxParam();
249       HLTDebug("Low flux params init.");
250     }
251
252   if (iRecoParamType == 1)
253     {
254       fRecoParam = AliTRDrecoParam::GetHighFluxParam();
255       HLTDebug("High flux params init.");
256     }
257
258   if (iRecoParamType == 2)
259     {
260       fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
261       HLTDebug("Cosmic Test params init.");
262     }
263
264   if (fRecoParam == 0)
265     {
266       HLTError("No reco params initialized. Sniffing big trouble!");
267       return -1;
268     }
269
270   fReconstructor = new AliTRDReconstructor();
271   fReconstructor->SetRecoParam(fRecoParam);
272   fReconstructor->SetStreamLevel(0, AliTRDReconstructor::kClusterizer); // default value
273   if (bWriteClusters)
274     {
275     HLTInfo("Writing clusters. I.e. output is a TTree with clusters");
276   fReconstructor->SetOption("cw,sl_cf_0");
277     }
278   else
279     {
280       HLTInfo("Not witing clusters. I.e. output is a TClonesArray of clusters");
281       fReconstructor->SetOption("!cw,sl_cf_0");
282     }
283   
284   
285   // init the raw data type to be used...
286   // the switch here will become obsolete as soon as the data structures is fixed 
287   // both: in sim and reality
288   if (iRecoDataType < 0 || iRecoDataType > 1)
289     {
290       HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
291       iRecoDataType = 0;
292     }
293
294   if (iRecoDataType == 0)
295     {
296       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDsimStream);
297       HLTDebug("Data type expected is SIMULATION!");
298     }
299
300   if (iRecoDataType == 1)
301     {
302       AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDrealStream);
303       HLTDebug("Data type expected is EXPERIMENT!");
304     }
305
306   // the DATA BASE STUFF
307   fCDB = AliCDBManager::Instance();
308   if (!fCDB)
309     {
310       HLTError("Could not get CDB instance", "fCDB 0x%x", fCDB);
311     }
312   else
313     {
314       fCDB->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
315       fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
316       HLTDebug("CDB instance; fCDB 0x%x", fCDB);
317     }
318
319   fGeometryFile = TFile::Open(fGeometryFileName.c_str());
320   if (fGeometryFile)
321     {
322       AliGeomManager::LoadGeometry(fGeometryFileName.c_str());
323     }
324   else
325     {
326       HLTError("Unable to open file. FATAL!");
327       return -1;
328     }
329
330   fMemReader = new AliRawReaderMemory;
331
332   fClusterizer = new AliTRDclusterizerHLT("TRDCclusterizer", "TRDCclusterizer");
333   fClusterizer->SetReconstructor(fReconstructor);
334
335   fClusterizer->SetRawVersion(iRawDataVersion);
336   fClusterizer->InitClusterTree();
337   return 0;
338 }
339
340 int AliHLTTRDClusterizerComponent::DoDeinit()
341 {
342   // Deinitialization of the component
343   delete fMemReader;
344   fMemReader = 0;
345   delete fClusterizer;
346   fClusterizer = 0;
347   
348   fReconstructor->SetClusters(0x0);
349   delete fReconstructor;
350   fReconstructor = 0x0;
351   return 0;
352
353   if (fGeometryFile)
354     {
355       fGeometryFile->Close();
356       delete fGeometryFile;
357       fGeometryFile = 0;
358     }
359
360   if (fCDB)
361     {
362       HLTDebug("destroy fCDB");
363       fCDB->Destroy();
364       fCDB = 0;
365     }
366
367   if (fRecoParam)
368     {
369       HLTDebug("Deleting fRecoParam");
370       delete fRecoParam;
371       fRecoParam = 0;
372     }
373 }
374
375 int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponentEventData& evtData, 
376                                             const AliHLTComponentBlockData* blocks, 
377                                             AliHLTComponent_TriggerData& /*trigData*/, 
378                                             AliHLTUInt8_t* /*outputPtr*/, 
379                                             AliHLTUInt32_t& size, 
380                                             vector<AliHLTComponent_BlockData>& /*outputBlocks*/ )
381 {
382   // Process an event
383   HLTDebug( "NofBlocks %lu", evtData.fBlockCnt );
384   // Process an event
385   Bool_t bWriteClusters = fReconstructor->IsWritingClusters();
386
387   //implement a usage of the following
388   //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
389   //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
390   //   void *triggerData = trigData.fData;
391   //HLTDebug( "Trigger data received. Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
392
393   // Loop over all input blocks in the event
394   AliHLTComponentDataType expectedDataType = (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
395   for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
396     {
397       // lets not use the internal TRD data types here : AliHLTTRDDefinitions::fgkDDLRawDataType
398       // which is depreciated - we use HLT global defs instead
399       //      if ( blocks[i].fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD) )
400       AliHLTComponentDataType inputDataType = blocks[i].fDataType;
401       if ( inputDataType != expectedDataType)
402         {
403           HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s; Skipping",
404                     i, evtData.fBlockCnt,
405                     evtData.fEventID, evtData.fEventID, 
406                     DataType2Text(inputDataType).c_str(), 
407                     DataType2Text(expectedDataType).c_str());
408           continue;
409         }
410       else 
411         {
412           HLTDebug("We get the right data type!");
413         }
414       
415       //      fMemReader->Reset();
416       fMemReader->SetMemory((UChar_t*) blocks[i].fPtr, blocks[i].fSize);
417
418       AliHLTUInt32_t spec = blocks[i].fSpecification;
419       
420       Int_t id = 1024;
421       
422       for ( Int_t ii = 0; ii < 18 ; ii++ ) {
423         if ( spec & 0x00000001 ) {
424           id += ii;
425           break;
426         }
427         spec = spec >> 1 ;
428       }
429       
430       fMemReader->SetEquipmentID( id ); 
431
432       fClusterizer->ResetTree();  
433       Bool_t iclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
434       if (iclustered == kTRUE)
435         {
436           HLTDebug( "Clustered successfully");
437         }
438       else
439         {
440           HLTError("Clustering ERROR");
441           return -1;
442         }
443
444       // put the tree into output blocks
445       //fcTree->Print();
446       if (bWriteClusters)
447         {
448           TTree *fcTree = fClusterizer->GetClusterTree();
449           if (fcTree)
450             {
451               HLTDebug("fcTree: Entries - %i; Size - %i",fcTree->GetEntriesFast(),sizeof(fcTree));
452               PushBack(fcTree, AliHLTTRDDefinitions::fgkClusterDataType, blocks[i].fSpecification);
453             }
454           fClusterizer->RecPoints()->Delete();
455           
456         }
457       else
458         {
459           TClonesArray *clustersArray = fClusterizer->RecPoints();
460           fClusterizer->SetClustersOwner(kFALSE);
461
462           if (clustersArray){
463             clustersArray->BypassStreamer(kFALSE);
464             HLTDebug("clustersArray: Entries - %i; Size - %i",clustersArray->GetEntriesFast(),sizeof(clustersArray));
465             PushBack(clustersArray, AliHLTTRDDefinitions::fgkClusterDataType, blocks[i].fSpecification);
466             clustersArray->Delete();
467             delete clustersArray;
468           }
469           else 
470             HLTWarning("Array of clusters is empty!");
471         }
472       
473         
474       
475
476     }//  for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
477
478   
479   
480   
481   size=0; // this function did not write data to the buffer directly
482   return 0;
483 }