]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterizerComponent.cxx
112e95118d2c2dba3ff549b69acaa6e45f0289ff
[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 #if __GNUC__ >= 3
25 using namespace std;
26 #endif
27
28 #include "AliHLTTRDClusterizerComponent.h"
29 #include "AliHLTTRDDefinitions.h"
30
31 #include "AliCDBManager.h"
32 #include "AliTRDclusterizerV1HLT.h"
33 #include "AliRawReaderMemory.h"
34
35 #include "TTree.h"
36 #include "TBranch.h"
37
38 #include <cstdlib>
39 #include <cerrno>
40 #include <string>
41
42 // this is a global object used for automatic component registration, do not use this
43 AliHLTTRDClusterizerComponent gAliHLTTRDClusterizerComponent;
44
45 ClassImp(AliHLTTRDClusterizerComponent)
46     
47   AliHLTTRDClusterizerComponent::AliHLTTRDClusterizerComponent()
48 {
49   fOutputPercentage = 100; // By default we copy to the output exactly what we got as input
50   
51 }
52
53 AliHLTTRDClusterizerComponent::~AliHLTTRDClusterizerComponent()
54 {
55 }
56
57 const char* AliHLTTRDClusterizerComponent::GetComponentID()
58 {
59   return "TRDClusterizer"; // The ID of this component
60 }
61
62 void AliHLTTRDClusterizerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
63 {
64   list.clear(); // We do not have any requirements for our input data type(s).
65   list.push_back( AliHLTTRDDefinitions::fgkDDLRawDataType );
66 }
67
68 AliHLTComponent_DataType AliHLTTRDClusterizerComponent::GetOutputDataType()
69 {
70   return AliHLTTRDDefinitions::fgkClusterDataType;
71 }
72
73 void AliHLTTRDClusterizerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
74 {
75   constBase = 0;
76   inputMultiplier = ((double)fOutputPercentage)/100.0;
77 }
78
79
80
81 // Spawn function, return new instance of this class
82 AliHLTComponent* AliHLTTRDClusterizerComponent::Spawn()
83 {
84   return new AliHLTTRDClusterizerComponent;
85 };
86
87 int AliHLTTRDClusterizerComponent::DoInit( int argc, const char** argv )
88 {
89   // perform initialization. We check whether our relative output size is specified in the arguments.
90   fOutputPercentage = 100;
91   Int_t fRawDataVersion = 2;
92   int i = 0;
93   char* cpErr;
94   while ( i < argc )
95     {
96       Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoInit", "Arguments", "argv[%d] == %s", i, argv[i] );
97       if ( !strcmp( argv[i], "output_percentage" ) )
98         {
99           if ( i+1>=argc )
100             {
101               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Missing Argument", "Missing output_percentage parameter");
102               return ENOTSUP;
103             }
104           Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoInit", "Arguments", "argv[%d+1] == %s", i, argv[i+1] );
105           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
106           if ( *cpErr )
107             {
108               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Wrong Argument", "Cannot convert output_percentage parameter '%s'", argv[i+1] );
109               return EINVAL;
110             }
111           Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoInit", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
112           i += 2;
113           continue;
114         }
115
116       if ( strcmp( argv[i], "-cdb" ) == 0)
117         {
118           if ( i+1 >= argc )
119             {
120               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Missing Argument", "Missing -cdb argument");
121               return ENOTSUP;         
122             }
123           fStrorageDBpath = argv[i+1];
124           Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoInit", "DB storage set", "DB storage is %s", fStrorageDBpath.c_str() );   
125           i += 2;
126           continue;
127         }      
128
129       if ( strcmp( argv[i], "-rawver" ) == 0)
130         {
131           if ( i+1 >= argc )
132             {
133               Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Missing Argument", "Missing -rawver argument");
134               return ENOTSUP;         
135             }
136           fRawDataVersion = atoi( argv[i+1] );
137           Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoInit", "Raw Data", "Version is %d", fRawDataVersion );    
138           i += 2;
139           continue;
140         }      
141
142       Logging(kHLTLogError, "HLT::TRDClusterizer::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
143       return EINVAL;
144     }
145
146   cdb = AliCDBManager::Instance();
147   if (!cdb)
148     {
149       Logging(kHLTLogError, "HLT::TRDCalibration::DoInit", "Could not get CDB instance", "cdb 0x%x", cdb);
150     }
151   else
152     {
153       cdb->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
154       cdb->SetDefaultStorage(fStrorageDBpath.c_str());
155       Logging(kHLTLogDebug, "HLT::TRDCalibration::DoInit", "CDB instance", "cdb 0x%x", cdb);
156     }
157
158   rmem = new AliRawReaderMemory;
159   clusterizer = new AliTRDclusterizerV1HLT("TRDCclusterizer", "TRDCclusterizer");
160   clusterizer->SetRawDataVersion(fRawDataVersion);
161   clusterizer->InitClusterTree();
162   return 0;
163 }
164
165 int AliHLTTRDClusterizerComponent::DoDeinit()
166 {
167   delete rmem;
168   rmem = 0;
169   delete clusterizer;
170   clusterizer = 0;
171   return 0;
172
173   if (cdb)
174     {
175       Logging( kHLTLogDebug, "HLT::TRDCalibration::DoDeinit", "destroy", "cdb");
176       cdb->Destroy();
177       cdb = 0;
178     }
179 }
180
181 int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponent_EventData& evtData, const AliHLTComponent_BlockData* blocks, 
182                                             AliHLTComponent_TriggerData& trigData, AliHLTUInt8_t* outputPtr, 
183                                             AliHLTUInt32_t& size, vector<AliHLTComponent_BlockData>& outputBlocks )
184 {
185   Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
186   Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "BLOCKS", "NofBlocks %lu", evtData.fBlockCnt );
187   // Process an event
188   unsigned long totalSize = 0;
189   AliHLTUInt32_t fDblock_Specification = 0;
190
191   // Loop over all input blocks in the event
192   for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
193     {
194       char tmp1[14], tmp2[14];
195       DataType2Text( blocks[i].fDataType, tmp1 );
196       DataType2Text( AliHLTTRDDefinitions::fgkDDLRawDataType, tmp2 );      
197       Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoEvent", "Event received", 
198                "Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
199                evtData.fEventID, evtData.fEventID, tmp1, tmp2 );
200
201       if ( blocks[i].fDataType != AliHLTTRDDefinitions::fgkDDLRawDataType ) 
202         {
203           Logging (kHLTLogError, "HLT::TRDClusterizer::DoEvent", "COMPARE FAILED", "type=%d is type=%d",
204                    blocks[i].fDataType, AliHLTTRDDefinitions::fgkDDLRawDataType);
205           continue;
206         }
207       fDblock_Specification = blocks[i].fSpecification;
208       unsigned long blockSize = blocks[i].fSize;
209       totalSize += blockSize;
210     }
211
212   void *memBufIn = calloc(totalSize, 1);
213   AliHLTUInt8_t *pBuf = (AliHLTUInt8_t *)memBufIn;
214   if (memBufIn == NULL)
215     {
216       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "MEMORY", "Unable to allocate %lu bytes", totalSize);
217       return -1;
218     }
219
220   // Make the memory continuous
221   unsigned long copied = 0;
222   for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
223     {
224       if ( blocks[i].fDataType != AliHLTTRDDefinitions::fgkDDLRawDataType ) 
225         continue;
226
227       void *pos = (void*)(pBuf + copied);
228       void *copyret = memcpy(pos, blocks[i].fPtr, blocks[i].fSize);
229       if (copyret < 0)
230         {
231           Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "MEMORY", "Unable to copy %lu bytes", blocks[i].fSize);
232           return -1;
233         }
234       copied += blocks[i].fSize;
235     }
236
237   Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "COPY STATS", "total=%lu copied=%lu", totalSize, copied);
238
239   rmem->Reset();
240   rmem->SetMemory((UChar_t*)memBufIn, totalSize);
241   //rmem->Reset();
242   Bool_t ihead = rmem->ReadHeader();
243   if (ihead == kTRUE)
244     {
245       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "HEADER", "Header read successfully");
246     }
247   else
248     {
249       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "HEADER", "Header read ERROR");
250       //return -1; -- not FATAL
251     }
252
253   clusterizer->ResetTree();
254   Bool_t ireadD = clusterizer->ReadDigits(rmem);
255   if (ireadD == kTRUE)
256     {
257       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "DIGITS", "Digits read successfully");
258     }
259   else
260     {
261       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "DIGITS", "Digits read ERROR");
262       return -1;
263     }
264
265   Bool_t iclustered = clusterizer->MakeClusters();
266   if (iclustered == kTRUE)
267     {
268       Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "CLUSTERS", "Clustered successfully");
269     }
270   else
271     {
272       Logging( kHLTLogError, "HLT::TRDClusterizer::DoEvent", "CLUSTERS", "Clustering ERROR");
273       return -1;
274     }
275
276   free(memBufIn);
277
278   //UInt_t memBufOutSize = 0;
279   //   void *memBufOut = clusterizer->WriteClustersToMemory(memBufOut, memBufOutSize);
280   Int_t iNclusters = clusterizer->GetNclusters();
281   Logging( kHLTLogInfo, "HLT::TRDClusterizer::DoEvent", "COUNT", "N of Clusters = %d", iNclusters);
282
283   // put the tree into output blocks of TObjArrays
284   TTree *fcTree = clusterizer->GetClusterTree();
285   TList *lt = (TList*)fcTree->GetListOfBranches();
286   TIter it(lt);
287   it.Reset();
288   TBranch *tb = 0;
289   while ((tb = (TBranch*)it.Next()) != 0)
290     {
291       TObjArray *clusters = 0;
292       tb->SetAddress(&clusters);
293       for (Int_t icb = 0; icb < tb->GetEntries(); icb++)
294         {
295           tb->GetEntry(icb);
296           PushBack(clusters, AliHLTTRDDefinitions::fgkClusterDataType, fDblock_Specification);
297         }
298     }
299
300   return 0;
301 }