]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
implementation of SOR and EOR events (not yet enabled), bugfix in DataBuffer: data...
[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 "AliCDBManager.h"
42 #include "AliTRDclusterizerHLT.h"
43 #include "AliRawReaderMemory.h"
44
45 #include <cstdlib>
46 #include <cerrno>
47 #include <string>
48
49 // this is a global object used for automatic component registration, do not use this
50 AliHLTTRDClusterizerComponent gAliHLTTRDClusterizerComponent;
51
52 ClassImp(AliHLTTRDClusterizerComponent);
53    
54 AliHLTTRDClusterizerComponent::AliHLTTRDClusterizerComponent()
55   : AliHLTProcessor()
56   , fOutputPercentage(100) // By default we copy to the output exactly what we got as input  
57   , fStrorageDBpath("local://$ALICE_ROOT")
58   , fClusterizer(NULL)
59   , fCDB(NULL)
60   , fMemReader(NULL)
61   , fGeometryFileName("")
62   , fGeometryFile(NULL)
63   , fGeoManager(NULL)
64 {
65   // Default constructor
66
67   fGeometryFileName = getenv("ALICE_ROOT");
68   fGeometryFileName += "/HLT/TRD/geometry.root";
69 }
70
71 AliHLTTRDClusterizerComponent::~AliHLTTRDClusterizerComponent()
72 {
73   // Destructor
74   ;
75 }
76
77 const char* AliHLTTRDClusterizerComponent::GetComponentID()
78 {
79   // Return the component ID const char *
80   return "TRDClusterizer"; // The ID of this component
81 }
82
83 void AliHLTTRDClusterizerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
84 {
85   // Get the list of input data
86   list.clear(); // We do not have any requirements for our input data type(s).
87   list.push_back( AliHLTTRDDefinitions::fgkDDLRawDataType );
88 }
89
90 AliHLTComponent_DataType AliHLTTRDClusterizerComponent::GetOutputDataType()
91 {
92   // Get the output data type
93   return AliHLTTRDDefinitions::fgkClusterDataType;
94 }
95
96 void AliHLTTRDClusterizerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
97 {
98   // Get the output data size
99   constBase = 0;
100   inputMultiplier = ((double)fOutputPercentage)/100.0;
101 }
102
103 AliHLTComponent* AliHLTTRDClusterizerComponent::Spawn()
104 {
105   // Spawn function, return new instance of this class
106   return new AliHLTTRDClusterizerComponent;
107 };
108
109 int AliHLTTRDClusterizerComponent::DoInit( int argc, const char** argv )
110 {
111   // perform initialization. We check whether our relative output size is specified in the arguments.
112   fOutputPercentage = 100;
113   Int_t fRawDataVersion = 2;
114   int i = 0;
115   char* cpErr;
116   while ( i < argc )
117     {
118       Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoInit", "Arguments", "argv[%d] == %s", i, argv[i] );
119       if ( !strcmp( argv[i], "output_percentage" ) )
120         {
121           if ( i+1>=argc )
122             {
123               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Missing Argument", "Missing output_percentage parameter");
124               return ENOTSUP;
125             }
126           Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoInit", "Arguments", "argv[%d+1] == %s", i, argv[i+1] );
127           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
128           if ( *cpErr )
129             {
130               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Wrong Argument", "Cannot convert output_percentage parameter '%s'", argv[i+1] );
131               return EINVAL;
132             }
133           Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoInit", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
134           i += 2;
135           continue;
136         }
137
138       if ( strcmp( argv[i], "-cdb" ) == 0)
139         {
140           if ( i+1 >= argc )
141             {
142               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Missing Argument", "Missing -cdb argument");
143               return ENOTSUP;         
144             }
145           fStrorageDBpath = argv[i+1];
146           Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoInit", "DB storage set", "DB storage is %s", fStrorageDBpath.c_str() );   
147           i += 2;
148           continue;
149         }      
150
151       if ( strcmp( argv[i], "-rawver" ) == 0)
152         {
153           if ( i+1 >= argc )
154             {
155               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Missing Argument", "Missing -rawver argument");
156               return ENOTSUP;         
157             }
158           fRawDataVersion = atoi( argv[i+1] );
159           Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoInit", "Raw Data", "Version is %d", fRawDataVersion );    
160           i += 2;
161           continue;
162         }      
163
164       if ( strcmp( argv[i], "-geometry" ) == 0)
165         {
166           if ( i+1 >= argc )
167             {
168               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing -geometry argument");
169               return ENOTSUP;         
170             }
171           fGeometryFileName = argv[i+1];
172           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "GeomFile storage set", "GeomFile storage is %s", 
173                    fGeometryFileName.c_str() );   
174           i += 2;
175           continue;
176         }      
177
178       Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
179       return EINVAL;
180     }
181
182   fCDB = AliCDBManager::Instance();
183   if (!fCDB)
184     {
185       Logging(kHLTLogError, "HLT::TRDCalibration::DoInit", "Could not get CDB instance", "fCDB 0x%x", fCDB);
186     }
187   else
188     {
189       fCDB->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
190       fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
191       Logging(kHLTLogDebug, "HLT::TRDCalibration::DoInit", "CDB instance", "fCDB 0x%x", fCDB);
192     }
193
194   fGeometryFile = TFile::Open(fGeometryFileName.c_str());
195   if (fGeometryFile)
196     {
197       fGeoManager = (TGeoManager *)fGeometryFile->Get("Geometry");
198     }
199   else
200     {
201       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "fGeometryFile", "Unable to open file. FATAL!");
202       return -1;
203     }
204
205   fMemReader = new AliRawReaderMemory;
206
207   fClusterizer = new AliTRDclusterizerHLT("TRDCclusterizer", "TRDCclusterizer");
208   fClusterizer->SetRawVersion(fRawDataVersion);
209   fClusterizer->InitClusterTree();
210   return 0;
211 }
212
213 int AliHLTTRDClusterizerComponent::DoDeinit()
214 {
215   // Deinitialization of the component
216   delete fMemReader;
217   fMemReader = 0;
218   delete fClusterizer;
219   fClusterizer = 0;
220   return 0;
221
222   if (fGeometryFile)
223     {
224       fGeometryFile->Close();
225       delete fGeometryFile;
226       fGeometryFile = 0;
227     }
228
229   if (fCDB)
230     {
231       Logging( kHLTLogDebug, "HLT::TRDCalibration::DoDeinit", "destroy", "fCDB");
232       fCDB->Destroy();
233       fCDB = 0;
234     }
235 }
236
237 int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponent_EventData& evtData, const AliHLTComponent_BlockData* blocks, 
238                                             AliHLTComponent_TriggerData& trigData, AliHLTUInt8_t* outputPtr, 
239                                             AliHLTUInt32_t& size, vector<AliHLTComponent_BlockData>& outputBlocks )
240 {
241   // Process an event
242 //   Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
243   Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoEvent", "BLOCKS", "NofBlocks %lu", evtData.fBlockCnt );
244   // Process an event
245   unsigned long totalSize = 0;
246   AliHLTUInt32_t dBlockSpecification = 0;
247
248   //implement a usage of the following
249 //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
250 //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
251 //   void *triggerData = trigData.fData;
252   Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoEvent", "Trigger data received", 
253            "Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
254   Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoEvent", "Output status", 
255            "Output pointer at 0x%x Output vector blocks at 0x%x", outputPtr, &outputBlocks);
256
257   // Loop over all input blocks in the event
258   for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
259     {
260       char tmp1[14], tmp2[14];
261       DataType2Text( blocks[i].fDataType, tmp1 );
262       DataType2Text( AliHLTTRDDefinitions::fgkDDLRawDataType, tmp2 );      
263       Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoEvent", "Event received", 
264                "Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
265                evtData.fEventID, evtData.fEventID, tmp1, tmp2 );
266
267       if ( blocks[i].fDataType != AliHLTTRDDefinitions::fgkDDLRawDataType ) 
268         {
269           Logging (kHLTLogError, "HLT::TRDClusterizer::DoEvent", "COMPARE FAILED", "type=%d is type=%d",
270                    blocks[i].fDataType, AliHLTTRDDefinitions::fgkDDLRawDataType);
271           continue;
272         }
273       dBlockSpecification = blocks[i].fSpecification;
274 //       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "CHECKSPEC", "fDblock_Spec %d %d", i, dBlockSpecification);
275       unsigned long blockSize = blocks[i].fSize;
276       totalSize += blockSize;
277     }
278
279   void *memBufIn = calloc(totalSize, 1);
280   AliHLTUInt8_t *pBuf = (AliHLTUInt8_t *)memBufIn;
281   if (memBufIn == NULL)
282     {
283       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "MEMORY", "Unable to allocate %lu bytes", totalSize);
284       return -1;
285     }
286
287   // Make the memory continuous
288   unsigned long copied = 0;
289   for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
290     {
291       if ( blocks[i].fDataType != AliHLTTRDDefinitions::fgkDDLRawDataType ) 
292         continue;
293
294       void *pos = (void*)(pBuf + copied);
295       void *copyret = memcpy(pos, blocks[i].fPtr, blocks[i].fSize);
296       if (copyret < 0)
297         {
298           Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "MEMORY", "Unable to copy %lu bytes", blocks[i].fSize);
299           return -1;
300         }
301       copied += blocks[i].fSize;
302     }
303
304 //   Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "COPY STATS", "total=%lu copied=%lu", totalSize, copied);
305
306   fMemReader->Reset();
307   fMemReader->SetMemory((UChar_t*)memBufIn, totalSize);
308   //fMemReader->SelectEquipment(0, 1024, 1041);
309   fMemReader->SetEquipmentID(1024);
310   //fMemReader->Reset();
311 //   Bool_t ihead = fMemReader->ReadHeader();
312 //   if (ihead == kTRUE)
313 //     {
314 //       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "HEADER", "Header read successfully");
315 //     }
316 //   else
317 //     {
318 //       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "HEADER", "Header read ERROR");
319 //       //return -1; -- not FATAL
320 //     }
321
322   fClusterizer->ResetTree();  
323   Bool_t iclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
324   if (iclustered == kTRUE)
325     {
326       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "CLUSTERS", "Clustered successfully");
327     }
328   else
329     {
330       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "CLUSTERS", "Clustering ERROR");
331       return -1;
332     }
333
334 //   AliRawReaderMemory reader;
335 //   reader.Reset();
336 //   reader.SetMemory((UChar_t*)memBufIn, totalSize);
337 //   //reader->Reset();
338 //   Bool_t ihead = reader.ReadHeader();
339 // //   if (ihead == kTRUE)
340 // //     {
341 // //       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "HEADER", "Header read successfully");
342 // //     }
343 // //   else
344 // //     {
345 // //       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "HEADER", "Header read ERROR");
346 // //       //return -1; -- not FATAL
347 // //     }
348
349 //   fClusterizer->ResetTree();
350   
351 //   Bool_t iclustered = fClusterizer->Raw2ClustersChamber(&reader);
352 //   if (iclustered == kTRUE)
353 //     {
354 //       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "CLUSTERS", "Clustered successfully");
355 //     }
356 //   else
357 //     {
358 //       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "CLUSTERS", "Clustering ERROR");
359 //       return -1;
360 //     }
361   
362   free(memBufIn);
363   
364   // put the tree into output blocks of TObjArrays
365   TTree *fcTree = fClusterizer->GetClusterTree();
366   
367   PushBack(fcTree, AliHLTTRDDefinitions::fgkClusterDataType, dBlockSpecification);
368   
369   Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "DONE", "Output size %d", size);
370   return 0;
371 }