2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project *
4 //* ALICE Experiment at CERN, All rights reserved. *
6 //* Primary Authors: Gaute Øvrebekk <st05886@alf.uib.no> *
7 //* for The ALICE HLT Project. *
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 //**************************************************************************
18 /** @file AliHLTITSClusterFinderComponent.cxx
19 @author Gaute Øvrebekk <st05886@alf.uib.no>
21 @brief Component to run offline clusterfinders
28 #include "AliHLTITSClusterFinderComponent.h"
30 #include "AliCDBEntry.h"
31 #include "AliCDBManager.h"
32 #include "AliCDBStorage.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"
39 #include "AliITSRecoParam.h"
40 #include "AliITSReconstructor.h"
41 #include "AliHLTITSClusterFinderSPD.h"
42 #include "AliHLTITSClusterFinderSSD.h"
47 #include "TObjString.h"
50 /** ROOT macro for the implementation of ROOT specific class methods */
51 ClassImp(AliHLTITSClusterFinderComponent);
53 AliHLTITSClusterFinderComponent::AliHLTITSClusterFinderComponent(int mode)
56 fInputDataType(kAliHLTVoidDataType),
57 fOutputDataType(kAliHLTVoidDataType),
73 // see header file for class documentation
75 // refer to README to build package
77 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
80 case kClusterFinderSPD:
81 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;
82 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
84 case kClusterFinderSDD:
85 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;
86 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
88 case kClusterFinderSSD:
89 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;
90 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
92 case kClusterFinderDigits:
93 fInputDataType = kAliHLTDataTypeAliTreeD|kAliHLTDataOriginITS;
94 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITS;
97 HLTFatal("unknown cluster finder");
101 AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent()
103 // see header file for class documentation
111 // Public functions to implement AliHLTComponent's interface.
112 // These functions are required for the registration process
114 const char* AliHLTITSClusterFinderComponent::GetComponentID()
116 // see header file for class documentation
118 case kClusterFinderSPD:
119 return "ITSClusterFinderSPD";
121 case kClusterFinderSDD:
122 return "ITSClusterFinderSDD";
124 case kClusterFinderSSD:
125 return "ITSClusterFinderSSD";
127 case kClusterFinderDigits:
128 return "ITSClusterFinderDigits";
134 void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
136 // see header file for class documentation
138 list.push_back( fInputDataType );
141 AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType()
143 // see header file for class documentation
144 return fOutputDataType;
147 void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
148 // see header file for class documentation
150 inputMultiplier = 100;
153 AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
154 // see header file for class documentation
155 return new AliHLTITSClusterFinderComponent(fModeSwitch);
158 Int_t AliHLTITSClusterFinderComponent::DoInit( int argc, const char** argv ) {
159 // see header file for class documentation
168 Int_t runNo = GetRunNo();
169 AliCDBStorage* store = AliCDBManager::Instance()->GetDefaultStorage();
176 if(store->GetLatestVersion("ITS/Calib/SPDNoisy", runNo)<0){
177 HLTError("SPDNoisy is not found in SPD/Calib");
180 if(store->GetLatestVersion("ITS/Calib/SPDDead", runNo)<0){
181 HLTError("SPDDead is not found in SPD/Calib");
184 if(store->GetLatestVersion("TRIGGER/SPD/PITConditions", runNo)<0){
185 HLTError("PITConditions is not found in TRIGGER/SPD");
190 if(store->GetLatestVersion("ITS/Calib/CalibSDD", runNo)<0){
191 HLTError("CalibSDD is not found in ITS/Calib");
194 if(store->GetLatestVersion("ITS/Calib/RespSDD", runNo)<0){
195 HLTError("RespSDD is not found in ITS/Calib");
198 if(store->GetLatestVersion("ITS/Calib/DriftSpeedSDD", runNo)<0){
199 HLTError("DriftSpeedSDD is not found in ITS/Calib");
202 if(store->GetLatestVersion("ITS/Calib/DDLMapSDD", runNo)<0){
203 HLTError("DDLMapSDD is not found in ITS/Calib");
206 if(store->GetLatestVersion("ITS/Calib/MapsTimeSDD", runNo)<0){
207 HLTError("MapsTimeSDD is not found in ITS/Calib");
212 if(store->GetLatestVersion("ITS/Calib/NoiseSSD", runNo)<0){
213 HLTError("NoiseSSD is not found in ITS/Calib");
216 if(store->GetLatestVersion("ITS/Calib/GainSSD", runNo)<0){
217 HLTError("GainSSD is not found in ITS/Calib");
220 if(store->GetLatestVersion("ITS/Calib/BadChannelsSSD", runNo)<0){
221 HLTError("BadChannelsSSD is not found in ITS/Calib");
225 //General reconstruction
226 if(store->GetLatestVersion("GRP/CTP/Scalers", runNo)<0){
227 HLTError("Scalers is not found in GRP/CTP/");
230 if(!cdbOK){return NULL;}
232 if(fModeSwitch==kClusterFinderSPD) {
233 HLTDebug("using ClusterFinder for SPD");
234 //fNModules=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1) + AliITSgeomTGeo::GetNDetectors(2)*AliITSgeomTGeo::GetNLadders(2);
235 fId=AliHLTDAQ::DdlIDOffset("ITSSPD");
236 fNddl=AliHLTDAQ::NumberOfDdls("ITSSPD");
238 else if(fModeSwitch==kClusterFinderSDD) {
239 HLTDebug("using ClusterFinder for SDD");
240 //fNModules=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3) + AliITSgeomTGeo::GetNDetectors(4)*AliITSgeomTGeo::GetNLadders(4);
241 fId=AliHLTDAQ::DdlIDOffset("ITSSDD");
242 fNddl=AliHLTDAQ::NumberOfDdls("ITSSDD");
244 else if(fModeSwitch==kClusterFinderSSD) {
245 HLTDebug("using ClusterFinder for SSD");
246 //fNModules=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5) + AliITSgeomTGeo::GetNDetectors(6)*AliITSgeomTGeo::GetNLadders(6);
247 fId=AliHLTDAQ::DdlIDOffset("ITSSSD");
248 fNddl=AliHLTDAQ::NumberOfDdls("ITSSSD");
250 else if(fModeSwitch==kClusterFinderDigits) {
254 HLTFatal("No mode set for clusterfindercomponent");
257 //Removed the warning for loading default RecoParam in HLT
258 AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
259 AliITSReconstructor *rec = new AliITSReconstructor();
260 rec->SetRecoParam(par);
262 AliCDBManager* man = AliCDBManager::Instance();
263 if (!man->IsDefaultStorageSet()){
264 HLTError("Default CDB storage has not been set !");
268 if(AliGeomManager::GetGeometry()==NULL){
269 AliGeomManager::LoadGeometry();
272 //fgeomInit = new AliITSInitGeometry(kvSPD02,2);
273 fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
274 //fgeomInit->InitAliITSgeom(fgeom);
275 fgeom = fgeomInit->CreateAliITSgeom();
277 fNModules = fgeom->GetIndexMax();
279 fClusters = new TClonesArray*[fNModules];
280 for (Int_t iModule = 0; iModule < fNModules; iModule++) {
281 fClusters[iModule] = NULL;
285 fDettype = new AliITSDetTypeRec();
286 fDettype->SetITSgeom(fgeom);
287 fDettype->SetDefaults();
288 fDettype->SetDefaultClusterFindersV2(kTRUE);
293 fRawReader = new AliRawReaderMemory();
294 fSPD = new AliHLTITSClusterFinderSPD( fDettype );
295 fSSD = new AliHLTITSClusterFinderSSD( fDettype, fRawReader );
297 TString arguments = "";
298 for ( int i = 0; i < argc; i++ ) {
299 if ( !arguments.IsNull() ) arguments += " ";
300 arguments += argv[i];
306 return Configure( arguments.Data() );
309 Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
310 // see header file for class documentation
330 for (Int_t iModule = 0; iModule < fNModules; iModule++) {
331 if(fClusters[iModule] != NULL){
332 fClusters[iModule]->Delete();
333 delete fClusters[iModule];
335 fClusters[iModule] = NULL;
338 fUseOfflineFinder = 0;
343 // #include "TStopwatch.h"
345 int AliHLTITSClusterFinderComponent::DoEvent
347 const AliHLTComponentEventData& evtData,
348 const AliHLTComponentBlockData* blocks,
349 AliHLTComponentTriggerData& /*trigData*/,
350 AliHLTUInt8_t* outputPtr,
351 AliHLTUInt32_t& size,
352 vector<AliHLTComponentBlockData>& outputBlocks )
354 // see header file for class documentation
356 AliHLTUInt32_t maxBufferSize = size;
357 size = 0; // output size
359 if (!IsDataEvent()) return 0;
361 if ( evtData.fBlockCnt<=0 )
363 HLTDebug("no blocks in event" );
369 //std::vector<AliITSRecPoint> vclusters;
370 if(fModeSwitch==kClusterFinderDigits) {
371 for ( const TObject *iter = GetFirstInputObject(fInputDataType); iter != NULL; iter = GetNextInputObject() ) {
372 tD = dynamic_cast<TTree*>(const_cast<TObject*>( iter ) );
374 HLTFatal("No Digit Tree found");
378 fDettype->SetTreeAddressD(tD);
379 fDettype->MakeBranch(tR,"R");
380 fDettype->SetTreeAddressR(tR);
382 fDettype->DigitsToRecPoints(tD,tR,0,opt,kTRUE);
384 TClonesArray * fRecPoints;
385 tR->SetBranchAddress("ITSRecPoints",&fRecPoints);
386 for(Int_t treeEntry=0;treeEntry<tR->GetEntries();treeEntry++){
387 tR->GetEntry(treeEntry);
388 for(Int_t tCloneEntry=0;tCloneEntry<fRecPoints->GetEntries();tCloneEntry++){
389 AliITSRecPoint *recpoint=(AliITSRecPoint*)fRecPoints->At(tCloneEntry);
390 fclusters.push_back(*recpoint);
397 UInt_t nClusters=fclusters.size();
399 UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
400 if( size + bufferSize > maxBufferSize ){
401 HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
407 RecPointToSpacePoint(outputPtr,size);
409 AliHLTComponentBlockData bd;
412 bd.fSize = bufferSize;
413 bd.fSpecification = 0x00000000;
414 bd.fDataType = GetOutputDataType();
415 outputBlocks.push_back( bd );
424 // -- Loop over blocks
425 for( const AliHLTComponentBlockData* iter = GetFirstInputBlock(fInputDataType); iter != NULL; iter = GetNextInputBlock() ) {
427 // -- Debug output of datatype --
428 HLTDebug("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
429 evtData.fEventID, evtData.fEventID,
430 DataType2Text(iter->fDataType).c_str(),
431 DataType2Text(fInputDataType).c_str());
433 // -- Check for the correct data type
434 if ( iter->fDataType != (fInputDataType) )
437 // -- Get equipment ID out of specification
438 AliHLTUInt32_t spec = iter->fSpecification;
441 for ( Int_t ii = 0; ii < fNddl ; ii++ ) { //number of ddl's
442 if ( spec & 0x00000001 ) {
449 // -- Set equipment ID to the raw reader
451 if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
452 HLTWarning("Could not add buffer");
454 // TStopwatch timer1;
456 if(fModeSwitch==kClusterFinderSPD && !fUseOfflineFinder){ fSPD->RawdataToClusters( fRawReader, fclusters ); }
457 else if(fModeSwitch==kClusterFinderSSD && !fUseOfflineFinder){ fSSD->RawdataToClusters( fclusters ); }
459 if(fModeSwitch==kClusterFinderSPD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SPD");}
460 if(fModeSwitch==kClusterFinderSSD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SSD");}
461 if(fModeSwitch==kClusterFinderSDD) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SDD");}
462 for(int i=0;i<fNModules;i++){
463 if(fClusters[i] != NULL){
464 for(int j=0;j<fClusters[i]->GetEntriesFast();j++){
465 AliITSRecPoint *recpoint = (AliITSRecPoint*) (fClusters[i]->At(j));
466 fclusters.push_back(*recpoint);
468 fClusters[i]->Delete();
476 // fStatTime+=timer1.RealTime();
477 // fStatTimeC+=timer1.CpuTime();
479 fRawReader->ClearBuffers();
481 UInt_t nClusters=fclusters.size();
483 UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
484 if( size + bufferSize > maxBufferSize ){
485 HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
491 RecPointToSpacePoint(outputPtr,size);
493 AliHLTComponentBlockData bd;
496 bd.fSize = bufferSize;
497 bd.fSpecification = iter->fSpecification;
498 bd.fDataType = GetOutputDataType();
499 outputBlocks.push_back( bd );
510 fStatTimeAll+=timer.RealTime();
511 fStatTimeAllC+=timer.CpuTime();
513 if( fStatNEv%1000==0 && fStatTimeAll>0.0 && fStatTime>0.0 && fStatTimeAllC>0.0 && fStatTimeC>0.0)
514 cout<<fStatTimeAll/fStatNEv*1.e3<<" "<<fStatTime/fStatNEv*1.e3<<" "
515 <<fStatTimeAllC/fStatNEv*1.e3<<" "<<fStatTimeC/fStatNEv*1.e3<<" ms"<<endl;
521 int AliHLTITSClusterFinderComponent::Configure(const char* arguments)
526 if (!arguments) return iResult;
528 TString allArgs=arguments;
531 TObjArray* pTokens=allArgs.Tokenize(" ");
534 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
535 argument=((TObjString*)pTokens->At(i))->GetString();
536 if (argument.IsNull()) continue;
537 if (argument.CompareTo("-use-offline-finder")==0) {
538 fUseOfflineFinder = 1;
539 HLTInfo("Off-line ClusterFinder will be used");
543 else if (argument.CompareTo("")==0) {
549 HLTError("unknown argument %s", argument.Data());
560 int AliHLTITSClusterFinderComponent::Reconfigure(const char* cdbEntry, const char* chainId)
562 // see header file for class documentation
565 const char* path="HLT/ConfigITS/ClusterFinderComponent";
566 const char* defaultNotify="";
569 defaultNotify=" (default)";
572 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
573 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
575 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
577 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
578 iResult=Configure(pString->GetString().Data());
580 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
583 HLTError("can not fetch object \"%s\" from CDB", path);
589 void AliHLTITSClusterFinderComponent::RecPointToSpacePoint(AliHLTUInt8_t* outputPtr,AliHLTUInt32_t& size){
590 AliHLTITSClusterData *outputClusters = reinterpret_cast<AliHLTITSClusterData*>(outputPtr + size);
591 outputClusters->fSpacePointCnt=fclusters.size();
593 for(int i=0;i<fclusters.size();i++){
594 AliITSRecPoint *recpoint = (AliITSRecPoint*) &(fclusters[i]);
595 outputClusters->fSpacePoints[clustIdx].fY=recpoint->GetY();
596 outputClusters->fSpacePoints[clustIdx].fZ=recpoint->GetZ();
597 outputClusters->fSpacePoints[clustIdx].fSigmaY2=recpoint->GetSigmaY2();
598 outputClusters->fSpacePoints[clustIdx].fSigmaZ2=recpoint->GetSigmaZ2();
599 outputClusters->fSpacePoints[clustIdx].fSigmaYZ=recpoint->GetSigmaYZ();
600 outputClusters->fSpacePoints[clustIdx].fQ=recpoint->GetQ();
601 outputClusters->fSpacePoints[clustIdx].fNy=recpoint->GetNy();
602 outputClusters->fSpacePoints[clustIdx].fNz=recpoint->GetNz();
603 outputClusters->fSpacePoints[clustIdx].fLayer=recpoint->GetLayer();
604 outputClusters->fSpacePoints[clustIdx].fIndex=recpoint->GetDetectorIndex() | recpoint->GetPindex() | recpoint->GetNindex();
605 outputClusters->fSpacePoints[clustIdx].fTracks[0]=recpoint->GetLabel(0);
606 outputClusters->fSpacePoints[clustIdx].fTracks[1]=recpoint->GetLabel(1);
607 outputClusters->fSpacePoints[clustIdx].fTracks[2]=recpoint->GetLabel(2);