]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/ITS/AliHLTITSClusterFinderComponent.cxx
Related to the report #60402: protection in case of timestamp=0 with warnings added
[u/mrichter/AliRoot.git] / HLT / ITS / AliHLTITSClusterFinderComponent.cxx
CommitLineData
6b7742a2 1// $Id$
2//**************************************************************************
3//* This file is property of and copyright by the ALICE HLT Project *
4//* ALICE Experiment at CERN, All rights reserved. *
5//* *
6//* Primary Authors: Gaute Øvrebekk <st05886@alf.uib.no> *
7//* for The ALICE HLT Project. *
8//* *
9//* Permission to use, copy, modify and distribute this software and its *
10//* documentation strictly for non-commercial purposes is hereby granted *
11//* without fee, provided that the above copyright notice appears in all *
12//* copies and that both the copyright notice and this permission notice *
13//* appear in the supporting documentation. The authors make no claims *
14//* about the suitability of this software for any purpose. It is *
15//* provided "as is" without express or implied warranty. *
16//**************************************************************************
17
18/** @file AliHLTITSClusterFinderComponent.cxx
19 @author Gaute Øvrebekk <st05886@alf.uib.no>
20 @date
21 @brief Component to run offline clusterfinders
22*/
23
24#if __GNUC__>= 3
25using namespace std;
26#endif
27
28#include "AliHLTITSClusterFinderComponent.h"
29
30#include "AliCDBEntry.h"
31#include "AliCDBManager.h"
32#include "AliHLTDataTypes.h"
33#include "AliITSgeomTGeo.h"
34#include "AliITSRecPoint.h"
35#include "AliHLTITSSpacePointData.h"
36#include "AliHLTITSClusterDataFormat.h"
37#include <AliHLTDAQ.h>
38#include "AliGeomManager.h"
060f29ad 39#include "AliITSRecoParam.h"
40#include "AliITSReconstructor.h"
3f61e9ce 41#include "AliHLTITSClusterFinderSPD.h"
a2a2a7ce 42#include "AliHLTITSClusterFinderSSD.h"
6b7742a2 43
44#include <cstdlib>
45#include <cerrno>
46#include "TString.h"
47#include "TObjString.h"
48#include <sys/time.h>
49
50/** ROOT macro for the implementation of ROOT specific class methods */
51ClassImp(AliHLTITSClusterFinderComponent);
52
53AliHLTITSClusterFinderComponent::AliHLTITSClusterFinderComponent(int mode)
54 :
55 fModeSwitch(mode),
a2a2a7ce 56 fInputDataType(kAliHLTVoidDataType),
57 fOutputDataType(kAliHLTVoidDataType),
58 fUseOfflineFinder(0),
6b7742a2 59 fNModules(0),
60 fId(0),
61 fNddl(0),
d293cef8 62 fClusters(NULL),
6b7742a2 63 fRawReader(NULL),
64 fDettype(NULL),
65 fgeom(NULL),
3f61e9ce 66 fgeomInit(NULL),
a2a2a7ce 67 fSPD(NULL),
68 fSSD(NULL)
6b7742a2 69{
70 // see header file for class documentation
71 // or
72 // refer to README to build package
73 // or
74 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
a2a2a7ce 75
76 switch(fModeSwitch){
77 case kClusterFinderSPD:
78 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;
79 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
80 break;
81 case kClusterFinderSDD:
82 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;
83 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
84 break;
85 case kClusterFinderSSD:
86 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;
87 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
88 break;
89 default:
6b7742a2 90 HLTFatal("unknown cluster finder");
91 }
92}
93
3f61e9ce 94AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent()
95{
6b7742a2 96 // see header file for class documentation
3f61e9ce 97 delete fRawReader;
98 delete fDettype;
99 delete fgeomInit;
100 delete fSPD;
a2a2a7ce 101 delete fSSD;
6b7742a2 102}
103
104// Public functions to implement AliHLTComponent's interface.
105// These functions are required for the registration process
106
107const char* AliHLTITSClusterFinderComponent::GetComponentID()
108{
109 // see header file for class documentation
110 switch(fModeSwitch){
111 case kClusterFinderSPD:
112 return "ITSClusterFinderSPD";
113 break;
114 case kClusterFinderSDD:
115 return "ITSClusterFinderSDD";
116 break;
117 case kClusterFinderSSD:
118 return "ITSClusterFinderSSD";
119 break;
120 }
121 return "";
122}
123
a2a2a7ce 124void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
125{
6b7742a2 126 // see header file for class documentation
127 list.clear();
a2a2a7ce 128 list.push_back( fInputDataType );
6b7742a2 129}
130
a2a2a7ce 131AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType()
132{
6b7742a2 133 // see header file for class documentation
a2a2a7ce 134 return fOutputDataType;
6b7742a2 135}
136
137void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
138 // see header file for class documentation
139 constBase = 0;
140 inputMultiplier = 100;
141}
142
143AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
144 // see header file for class documentation
145 return new AliHLTITSClusterFinderComponent(fModeSwitch);
146}
147
1ea949fd 148Int_t AliHLTITSClusterFinderComponent::DoInit( int argc, const char** argv ) {
6b7742a2 149 // see header file for class documentation
a2a2a7ce 150 /*
3f61e9ce 151 fStatTime = 0;
152 fStatTimeAll = 0;
153 fStatTimeC = 0;
154 fStatTimeAllC = 0;
155 fStatNEv = 0;
a2a2a7ce 156 */
6b7742a2 157 if(fModeSwitch==kClusterFinderSPD) {
158 HLTDebug("using ClusterFinder for SPD");
d293cef8 159 //fNModules=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1) + AliITSgeomTGeo::GetNDetectors(2)*AliITSgeomTGeo::GetNLadders(2);
6b7742a2 160 fId=AliHLTDAQ::DdlIDOffset("ITSSPD");
161 fNddl=AliHLTDAQ::NumberOfDdls("ITSSPD");
162 }
163 else if(fModeSwitch==kClusterFinderSDD) {
164 HLTDebug("using ClusterFinder for SDD");
d293cef8 165 //fNModules=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3) + AliITSgeomTGeo::GetNDetectors(4)*AliITSgeomTGeo::GetNLadders(4);
6b7742a2 166 fId=AliHLTDAQ::DdlIDOffset("ITSSDD");
167 fNddl=AliHLTDAQ::NumberOfDdls("ITSSDD");
168 }
169 else if(fModeSwitch==kClusterFinderSSD) {
170 HLTDebug("using ClusterFinder for SSD");
d293cef8 171 //fNModules=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5) + AliITSgeomTGeo::GetNDetectors(6)*AliITSgeomTGeo::GetNLadders(6);
6b7742a2 172 fId=AliHLTDAQ::DdlIDOffset("ITSSSD");
173 fNddl=AliHLTDAQ::NumberOfDdls("ITSSSD");
174 }
175 else{
176 HLTFatal("No mode set for clusterfindercomponent");
177 }
178
060f29ad 179 //Removed the warning for loading default RecoParam in HLT
180 AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
181 AliITSReconstructor *rec = new AliITSReconstructor();
182 rec->SetRecoParam(par);
183
6b7742a2 184 AliCDBManager* man = AliCDBManager::Instance();
185 if (!man->IsDefaultStorageSet()){
186 HLTError("Default CDB storage has not been set !");
187 return -ENOENT;
188 }
189
190 if(AliGeomManager::GetGeometry()==NULL){
191 AliGeomManager::LoadGeometry();
192 }
3f61e9ce 193
194 //fgeomInit = new AliITSInitGeometry(kvSPD02,2);
6b7742a2 195 fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
196 //fgeomInit->InitAliITSgeom(fgeom);
197 fgeom = fgeomInit->CreateAliITSgeom();
3f61e9ce 198
d293cef8 199 fNModules = fgeom->GetIndexMax();
200
201 fClusters = new TClonesArray*[fNModules];
202 for (Int_t iModule = 0; iModule < fNModules; iModule++) {
203 fClusters[iModule] = NULL;
204 }
205
6b7742a2 206 //set dettype
207 fDettype = new AliITSDetTypeRec();
3f61e9ce 208 fDettype->SetITSgeom(fgeom);
6b7742a2 209 fDettype->SetDefaults();
86e746ed 210 fDettype->SetDefaultClusterFindersV2(kTRUE);
211
6b7742a2 212 if ( fRawReader )
213 return -EINPROGRESS;
214
215 fRawReader = new AliRawReaderMemory();
3f61e9ce 216 fSPD = new AliHLTITSClusterFinderSPD( fDettype );
a2a2a7ce 217 fSSD = new AliHLTITSClusterFinderSSD( fDettype, fRawReader );
6b7742a2 218
1ea949fd 219 TString arguments = "";
220 for ( int i = 0; i < argc; i++ ) {
221 if ( !arguments.IsNull() ) arguments += " ";
222 arguments += argv[i];
223 }
224
225 return Configure( arguments.Data() );
6b7742a2 226}
227
228Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
229 // see header file for class documentation
230
231 if ( fRawReader )
232 delete fRawReader;
233 fRawReader = NULL;
234
6b7742a2 235 if ( fDettype )
236 delete fDettype;
237 fDettype = NULL;
238
239 if ( fgeomInit )
240 delete fgeomInit;
241 fgeomInit = NULL;
3f61e9ce 242
243 delete fSPD;
244 fSPD = 0;
6b7742a2 245
a2a2a7ce 246 delete fSSD;
247 fSSD = 0;
248
d293cef8 249 for (Int_t iModule = 0; iModule < fNModules; iModule++) {
250 if(fClusters[iModule] != NULL){
251 fClusters[iModule]->Delete();
252 delete fClusters[iModule];
253 }
254 fClusters[iModule] = NULL;
255 }
256
a2a2a7ce 257 fUseOfflineFinder = 0;
258
6b7742a2 259 return 0;
260}
261
a2a2a7ce 262// #include "TStopwatch.h"
3f61e9ce 263
a2a2a7ce 264int AliHLTITSClusterFinderComponent::DoEvent
265(
266 const AliHLTComponentEventData& evtData,
267 const AliHLTComponentBlockData* blocks,
268 AliHLTComponentTriggerData& /*trigData*/,
269 AliHLTUInt8_t* outputPtr,
270 AliHLTUInt32_t& size,
271 vector<AliHLTComponentBlockData>& outputBlocks )
272{
273 // see header file for class documentation
274
275 AliHLTUInt32_t maxBufferSize = size;
276 size = 0; // output size
6b7742a2 277
a2a2a7ce 278 if (!IsDataEvent()) return 0;
279
6b7742a2 280 if ( evtData.fBlockCnt<=0 )
a2a2a7ce 281 {
282 HLTDebug("no blocks in event" );
283 return 0;
284 }
6b7742a2 285
a2a2a7ce 286 // TStopwatch timer;
287
288 Int_t ret = 0;
3f61e9ce 289
6b7742a2 290 // -- Loop over blocks
a2a2a7ce 291 for( const AliHLTComponentBlockData* iter = GetFirstInputBlock(fInputDataType); iter != NULL; iter = GetNextInputBlock() ) {
6b7742a2 292
293 // -- Debug output of datatype --
294 HLTDebug("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
295 evtData.fEventID, evtData.fEventID,
296 DataType2Text(iter->fDataType).c_str(),
df7307a2 297 DataType2Text(fInputDataType).c_str());
6b7742a2 298
299 // -- Check for the correct data type
a2a2a7ce 300 if ( iter->fDataType != (fInputDataType) )
6b7742a2 301 continue;
302
303 // -- Get equipment ID out of specification
304 AliHLTUInt32_t spec = iter->fSpecification;
305
306 Int_t id = fId;
307 for ( Int_t ii = 0; ii < fNddl ; ii++ ) { //number of ddl's
308 if ( spec & 0x00000001 ) {
309 id += ii;
310 break;
311 }
312 spec = spec >> 1 ;
313 }
314
315 // -- Set equipment ID to the raw reader
316
317 if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
318 HLTWarning("Could not add buffer");
319 }
a2a2a7ce 320 // TStopwatch timer1;
321
3f61e9ce 322 std::vector<AliITSRecPoint> vclusters;
6b7742a2 323
a2a2a7ce 324 if(fModeSwitch==kClusterFinderSPD && !fUseOfflineFinder){ fSPD->RawdataToClusters( fRawReader, vclusters ); }
325 else if(fModeSwitch==kClusterFinderSSD && !fUseOfflineFinder){ fSSD->RawdataToClusters( vclusters ); }
3f61e9ce 326 else{
a2a2a7ce 327 if(fModeSwitch==kClusterFinderSPD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SPD");}
328 if(fModeSwitch==kClusterFinderSSD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SSD");}
3f61e9ce 329 if(fModeSwitch==kClusterFinderSDD) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SDD");}
3f61e9ce 330 for(int i=0;i<fNModules;i++){
331 if(fClusters[i] != NULL){
332 for(int j=0;j<fClusters[i]->GetEntriesFast();j++){
333 AliITSRecPoint *recpoint = (AliITSRecPoint*) (fClusters[i]->At(j));
334 vclusters.push_back(*recpoint);
335 }
336 fClusters[i]->Delete();
337 delete fClusters[i];
338 }
339 fClusters[i] = NULL;
340 }
6b7742a2 341 }
342
a2a2a7ce 343 // timer1.Stop();
344 // fStatTime+=timer1.RealTime();
345 // fStatTimeC+=timer1.CpuTime();
3f61e9ce 346
a2a2a7ce 347 fRawReader->ClearBuffers();
348
3f61e9ce 349 UInt_t nClusters=vclusters.size();
350
6b7742a2 351 UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
a2a2a7ce 352 if( size + bufferSize > maxBufferSize ){
353 HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
354 ret = -ENOSPC;
355 break;
356 }
357 //cout<<"event "<<fStatNEv<<", nclu="<<nClusters<<":"<<endl;
358 if( nClusters>0 ){
359 AliHLTITSClusterData *outputClusters = reinterpret_cast<AliHLTITSClusterData*>(outputPtr + size);
360 outputClusters->fSpacePointCnt=nClusters;
361 int clustIdx=0;
362 for(int i=0;i<vclusters.size();i++){
363 AliITSRecPoint *recpoint = (AliITSRecPoint*) &(vclusters[i]);
364 //cout<<recpoint->GetDetectorIndex()<<" "<<recpoint->GetY()<<" "<<recpoint->GetZ()<<endl;
365 outputClusters->fSpacePoints[clustIdx].fY=recpoint->GetY();
366 outputClusters->fSpacePoints[clustIdx].fZ=recpoint->GetZ();
367 outputClusters->fSpacePoints[clustIdx].fSigmaY2=recpoint->GetSigmaY2();
368 outputClusters->fSpacePoints[clustIdx].fSigmaZ2=recpoint->GetSigmaZ2();
369 outputClusters->fSpacePoints[clustIdx].fSigmaYZ=recpoint->GetSigmaYZ();
370 outputClusters->fSpacePoints[clustIdx].fQ=recpoint->GetQ();
371 outputClusters->fSpacePoints[clustIdx].fNy=recpoint->GetNy();
372 outputClusters->fSpacePoints[clustIdx].fNz=recpoint->GetNz();
373 outputClusters->fSpacePoints[clustIdx].fLayer=recpoint->GetLayer();
374 outputClusters->fSpacePoints[clustIdx].fIndex=recpoint->GetDetectorIndex() | recpoint->GetPindex() | recpoint->GetNindex();
375 outputClusters->fSpacePoints[clustIdx].fTracks[0]=recpoint->GetLabel(0);
376 outputClusters->fSpacePoints[clustIdx].fTracks[1]=recpoint->GetLabel(1);
377 outputClusters->fSpacePoints[clustIdx].fTracks[2]=recpoint->GetLabel(2);
378 clustIdx++;
379 }
380 AliHLTComponentBlockData bd;
381 FillBlockData( bd );
382 bd.fOffset = size;
383 bd.fSize = bufferSize;
384 bd.fSpecification = iter->fSpecification;
385 bd.fDataType = GetOutputDataType();
386 outputBlocks.push_back( bd );
387 size += bufferSize;
6b7742a2 388 }
3f61e9ce 389
a2a2a7ce 390 } // input blocks
3f61e9ce 391
a2a2a7ce 392 /*
3f61e9ce 393 timer.Stop();
a2a2a7ce 394
3f61e9ce 395 fStatTimeAll+=timer.RealTime();
396 fStatTimeAllC+=timer.CpuTime();
397 fStatNEv++;
a2a2a7ce 398 if( fStatNEv%1000==0 && fStatTimeAll>0.0 && fStatTime>0.0 && fStatTimeAllC>0.0 && fStatTimeC>0.0)
399 cout<<fStatTimeAll/fStatNEv*1.e3<<" "<<fStatTime/fStatNEv*1.e3<<" "
400 <<fStatTimeAllC/fStatNEv*1.e3<<" "<<fStatTimeC/fStatNEv*1.e3<<" ms"<<endl;
401 */
3f61e9ce 402
a2a2a7ce 403 return ret;
6b7742a2 404}
405
406int AliHLTITSClusterFinderComponent::Configure(const char* arguments)
407{
408
409 int iResult=0;
410
411 if (!arguments) return iResult;
412
413 TString allArgs=arguments;
414 TString argument;
415
416 TObjArray* pTokens=allArgs.Tokenize(" ");
417
418 if (pTokens) {
419 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
420 argument=((TObjString*)pTokens->At(i))->GetString();
a2a2a7ce 421 if (argument.IsNull()) continue;
422 if (argument.CompareTo("-use-offline-finder")==0) {
423 fUseOfflineFinder = 1;
424 HLTInfo("Off-line ClusterFinder will be used");
6b7742a2 425 continue;
426 }
a2a2a7ce 427 /*
6b7742a2 428 else if (argument.CompareTo("")==0) {
429 HLTInfo("");
430 continue;
431 }
432 */
433 else {
434 HLTError("unknown argument %s", argument.Data());
435 iResult=-EINVAL;
436 break;
437 }
438 }
439 delete pTokens;
440 }
441
442 return iResult;
443}
444
445int AliHLTITSClusterFinderComponent::Reconfigure(const char* cdbEntry, const char* chainId)
446{
447 // see header file for class documentation
448 int iResult=0;
449
450 const char* path="HLT/ConfigITS/ClusterFinderComponent";
451 const char* defaultNotify="";
452 if (cdbEntry) {
453 path=cdbEntry;
454 defaultNotify=" (default)";
455 }
456 if (path) {
457 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
458 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
459 if (pEntry) {
460 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
461 if (pString) {
462 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
463 iResult=Configure(pString->GetString().Data());
464 } else {
465 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
466 }
467 } else {
468 HLTError("can not fetch object \"%s\" from CDB", path);
469 }
470 }
471
472 return iResult;
473}
474