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"
44 #include "AliITSRecPointContainer.h"
50 #include "TObjString.h"
53 /** ROOT macro for the implementation of ROOT specific class methods */
54 ClassImp(AliHLTITSClusterFinderComponent);
56 AliHLTITSClusterFinderComponent::AliHLTITSClusterFinderComponent(int mode)
59 fInputDataType(kAliHLTVoidDataType),
60 fOutputDataType(kAliHLTVoidDataType),
77 fBenchmark(GetComponentID())
79 // see header file for class documentation
81 // refer to README to build package
83 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
86 case kClusterFinderSPD:
87 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;
88 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
90 case kClusterFinderSDD:
91 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;
92 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
94 case kClusterFinderSSD:
95 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;
96 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
98 case kClusterFinderDigits:
99 fInputDataType = kAliHLTDataTypeAliTreeD|kAliHLTDataOriginITS;
100 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITS;
103 HLTFatal("unknown cluster finder");
107 AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent()
109 // see header file for class documentation
117 // Public functions to implement AliHLTComponent's interface.
118 // These functions are required for the registration process
120 const char* AliHLTITSClusterFinderComponent::GetComponentID()
122 // see header file for class documentation
124 case kClusterFinderSPD:
125 return "ITSClusterFinderSPD";
127 case kClusterFinderSDD:
128 return "ITSClusterFinderSDD";
130 case kClusterFinderSSD:
131 return "ITSClusterFinderSSD";
133 case kClusterFinderDigits:
134 return "ITSClusterFinderDigits";
140 void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
142 // see header file for class documentation
144 list.push_back( fInputDataType );
147 AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType()
149 // see header file for class documentation
150 return fOutputDataType;
153 void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
154 // see header file for class documentation
156 inputMultiplier = 100;
159 AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
160 // see header file for class documentation
161 return new AliHLTITSClusterFinderComponent(fModeSwitch);
164 Int_t AliHLTITSClusterFinderComponent::DoInit( int argc, const char** argv ) {
165 // see header file for class documentation
167 fBenchmark.SetTimer(0,"total");
168 fBenchmark.SetTimer(1,"reco");
178 Int_t runNo = GetRunNo();
179 AliCDBStorage* store = AliCDBManager::Instance()->GetDefaultStorage();
186 if(store->GetLatestVersion("ITS/Calib/SPDNoisy", runNo)<0){
187 HLTError("SPDNoisy is not found in SPD/Calib");
190 if(store->GetLatestVersion("ITS/Calib/SPDDead", runNo)<0){
191 HLTError("SPDDead is not found in SPD/Calib");
194 if(store->GetLatestVersion("TRIGGER/SPD/PITConditions", runNo)<0){
195 HLTError("PITConditions is not found in TRIGGER/SPD");
200 if(store->GetLatestVersion("ITS/Calib/CalibSDD", runNo)<0){
201 HLTError("CalibSDD is not found in ITS/Calib");
204 if(store->GetLatestVersion("ITS/Calib/RespSDD", runNo)<0){
205 HLTError("RespSDD is not found in ITS/Calib");
208 if(store->GetLatestVersion("ITS/Calib/DriftSpeedSDD", runNo)<0){
209 HLTError("DriftSpeedSDD is not found in ITS/Calib");
212 if(store->GetLatestVersion("ITS/Calib/DDLMapSDD", runNo)<0){
213 HLTError("DDLMapSDD is not found in ITS/Calib");
216 if(store->GetLatestVersion("ITS/Calib/MapsTimeSDD", runNo)<0){
217 HLTError("MapsTimeSDD is not found in ITS/Calib");
222 if(store->GetLatestVersion("ITS/Calib/NoiseSSD", runNo)<0){
223 HLTError("NoiseSSD is not found in ITS/Calib");
226 if(store->GetLatestVersion("ITS/Calib/GainSSD", runNo)<0){
227 HLTError("GainSSD is not found in ITS/Calib");
230 if(store->GetLatestVersion("ITS/Calib/BadChannelsSSD", runNo)<0){
231 HLTError("BadChannelsSSD is not found in ITS/Calib");
235 //General reconstruction
236 if(store->GetLatestVersion("GRP/CTP/Scalers", runNo)<0){
237 HLTError("Scalers is not found in GRP/CTP/");
240 if(!cdbOK){return NULL;}
242 if(fModeSwitch==kClusterFinderSPD) {
243 HLTDebug("using ClusterFinder for SPD");
244 //fNModules=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1) + AliITSgeomTGeo::GetNDetectors(2)*AliITSgeomTGeo::GetNLadders(2);
245 fId=AliHLTDAQ::DdlIDOffset("ITSSPD");
246 fNddl=AliHLTDAQ::NumberOfDdls("ITSSPD");
248 else if(fModeSwitch==kClusterFinderSDD) {
249 HLTDebug("using ClusterFinder for SDD");
250 //fNModules=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3) + AliITSgeomTGeo::GetNDetectors(4)*AliITSgeomTGeo::GetNLadders(4);
251 fId=AliHLTDAQ::DdlIDOffset("ITSSDD");
252 fNddl=AliHLTDAQ::NumberOfDdls("ITSSDD");
254 else if(fModeSwitch==kClusterFinderSSD) {
255 HLTDebug("using ClusterFinder for SSD");
256 //fNModules=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5) + AliITSgeomTGeo::GetNDetectors(6)*AliITSgeomTGeo::GetNLadders(6);
257 fId=AliHLTDAQ::DdlIDOffset("ITSSSD");
258 fNddl=AliHLTDAQ::NumberOfDdls("ITSSSD");
260 else if(fModeSwitch==kClusterFinderDigits) {
264 HLTFatal("No mode set for clusterfindercomponent");
267 //Removed the warning for loading default RecoParam in HLT
268 AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
269 AliITSReconstructor *rec = new AliITSReconstructor();
270 rec->SetRecoParam(par);
272 AliCDBManager* man = AliCDBManager::Instance();
273 if (!man->IsDefaultStorageSet()){
274 HLTError("Default CDB storage has not been set !");
278 if(AliGeomManager::GetGeometry()==NULL){
279 AliGeomManager::LoadGeometry();
282 fgeomInit = new AliITSInitGeometry();
283 //fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
284 //fgeomInit->InitAliITSgeom(fgeom);
285 fgeom = fgeomInit->CreateAliITSgeom();
287 fNModules = fgeom->GetIndexMax();
289 for(Int_t i=0;i<6;i++)modperlay[i]=AliITSgeomTGeo::GetNDetectors(1+i)*AliITSgeomTGeo::GetNLadders(1+i);
290 fSPDNModules=modperlay[0]+modperlay[1];
291 fSDDNModules=modperlay[2]+modperlay[3];
292 fSSDNModules=modperlay[4]+modperlay[5];
294 if(fModeSwitch==kClusterFinderSPD) {
296 fLastModule=fSPDNModules;
298 else if(fModeSwitch==kClusterFinderSDD) {
299 fFirstModule=fSPDNModules;
300 fLastModule=fFirstModule + fSDDNModules;
302 else if(fModeSwitch==kClusterFinderSSD) {
303 fFirstModule=fSPDNModules + fSDDNModules;
304 fLastModule=fFirstModule + fSSDNModules;
308 fDettype = new AliITSDetTypeRec();
309 fDettype->SetITSgeom(fgeom);
310 fDettype->SetDefaults();
311 fDettype->SetDefaultClusterFindersV2(kTRUE);
316 fRawReader = new AliRawReaderMemory();
317 fSPD = new AliHLTITSClusterFinderSPD( fDettype );
318 fSSD = new AliHLTITSClusterFinderSSD( fDettype, fRawReader );
320 TString arguments = "";
321 for ( int i = 0; i < argc; i++ ) {
322 if ( !arguments.IsNull() ) arguments += " ";
323 arguments += argv[i];
329 return Configure( arguments.Data() );
332 Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
333 // see header file for class documentation
353 fUseOfflineFinder = 0;
358 int AliHLTITSClusterFinderComponent::DoEvent
360 const AliHLTComponentEventData& evtData,
361 const AliHLTComponentBlockData* /*blocks*/,
362 AliHLTComponentTriggerData& /*trigData*/,
363 AliHLTUInt8_t* outputPtr,
364 AliHLTUInt32_t& size,
365 vector<AliHLTComponentBlockData>& outputBlocks )
367 // see header file for class documentation
369 AliHLTUInt32_t maxBufferSize = size;
370 size = 0; // output size
372 if (!IsDataEvent()) return 0;
374 if ( evtData.fBlockCnt<=0 )
376 HLTDebug("no blocks in event" );
380 fBenchmark.StartNewEvent();
382 for( const AliHLTComponentBlockData *i= GetFirstInputBlock(fInputDataType); i!=NULL; i=GetNextInputBlock() ){
383 fBenchmark.AddInput(i->fSize);
388 if(fModeSwitch==kClusterFinderDigits) {
390 for ( const TObject *iter = GetFirstInputObject(fInputDataType); iter != NULL; iter = GetNextInputObject() ) {
391 tD = dynamic_cast<TTree*>(const_cast<TObject*>( iter ) );
393 HLTFatal("No Digit Tree found");
396 // 2010-04-17 very crude workaround: TTree objects are difficult to send
397 // The actual case: Running ITS and TPC reconstruction fails at the second event
398 // to read the ITS digits from the TreeD
400 // Reason: reading fails in TBranch::GetBasket, there a new basket is opened from
401 // a TFile object. The function TBranch::GetFile returns the file object from
402 // an internal fDirectory (TDirectory) object. This file is at the second event
403 // set to the TPC.Digits.root. The internal mismatch creates a seg fault
405 // Investigation: TBranch::Streamer uses a crude assignment after creating the
407 // fDirectory = gDirectory;
408 // gDirectory is obviously not set correctly. Setting the directory to a TFile
409 // object for the ITS digits helps to fix the internal mess. Tried also to set
410 // the Directory for the TreeD to NULL (This has only effect if ones sets it
411 // to something not NULL first, and then to NULL). But then no content, i.e.
412 // ITS clusters could be retrieved.
414 // Conclusion: TTree objects are hardly to be sent via TMessage, there are direct
415 // links to the file required anyhow.
416 TFile* dummy=new TFile("ITS.Digits.root");
417 tD->SetDirectory(dummy);
420 fDettype->SetTreeAddressD(tD);
421 fDettype->MakeBranch(tR,"R");
422 fDettype->SetTreeAddressR(tR);
425 fDettype->DigitsToRecPoints(tD,tR,0,opt,0);
427 TClonesArray * fRecPoints = NULL;
428 tR->SetBranchAddress("ITSRecPoints",&fRecPoints);
429 for(Int_t treeEntry=0;treeEntry<tR->GetEntries();treeEntry++){
430 tR->GetEntry(treeEntry);
431 for(Int_t tCloneEntry=0;tCloneEntry<fRecPoints->GetEntries();tCloneEntry++){
432 AliITSRecPoint *recpoint=(AliITSRecPoint*)fRecPoints->At(tCloneEntry);
433 fclusters.push_back(*recpoint);
443 UInt_t nClusters=fclusters.size();
445 UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
446 if( size + bufferSize > maxBufferSize ){
447 HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
453 RecPointToSpacePoint(outputPtr,size);
455 AliHLTComponentBlockData bd;
458 bd.fSize = bufferSize;
459 bd.fSpecification = 0x00000000;
460 bd.fDataType = GetOutputDataType();
461 outputBlocks.push_back( bd );
463 fBenchmark.AddOutput(bd.fSize);
470 AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
472 if(fUseOfflineFinder){
473 if(fModeSwitch==kClusterFinderSPD){rpc->ResetSPD();}
474 if(fModeSwitch==kClusterFinderSSD){rpc->ResetSSD();}
476 if(fModeSwitch==kClusterFinderSDD){rpc->ResetSDD();}
478 // -- Loop over blocks
479 for( const AliHLTComponentBlockData* iter = GetFirstInputBlock(fInputDataType); iter != NULL; iter = GetNextInputBlock() ) {
481 // -- Debug output of datatype --
482 HLTDebug("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
483 evtData.fEventID, evtData.fEventID,
484 DataType2Text(iter->fDataType).c_str(),
485 DataType2Text(fInputDataType).c_str());
487 // -- Check for the correct data type
488 if ( iter->fDataType != (fInputDataType) )
491 // -- Get equipment ID out of specification
492 AliHLTUInt32_t spec = iter->fSpecification;
495 for ( Int_t ii = 0; ii < fNddl ; ii++ ) { //number of ddl's
496 if ( spec & 0x00000001 ) {
503 // -- Set equipment ID to the raw reader
506 HLTWarning("The fRawReader pointer is NULL");
510 if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
511 HLTWarning("Could not add buffer");
516 if(fModeSwitch==kClusterFinderSPD && !fUseOfflineFinder){ fSPD->RawdataToClusters( fRawReader, fclusters ); }
517 else if(fModeSwitch==kClusterFinderSSD && !fUseOfflineFinder){ fSSD->RawdataToClusters( fclusters ); }
519 if(fModeSwitch==kClusterFinderSPD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,"SPD");}
520 if(fModeSwitch==kClusterFinderSSD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,"SSD");}
521 if(fModeSwitch==kClusterFinderSDD) {fDettype->DigitsToRecPoints(fRawReader,"SDD");}
522 //AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
523 TClonesArray* clusters = NULL;
524 for(int i=fFirstModule;i<fLastModule;i++){
525 clusters = rpc->UncheckedGetClusters(i);
526 if(clusters != NULL){
527 for(int j=0;j<clusters->GetEntriesFast();j++){
528 AliITSRecPoint *recpoint = (AliITSRecPoint*) (clusters->At(j));
529 fclusters.push_back(*recpoint);
535 if(fUseOfflineFinder){
536 if(fModeSwitch==kClusterFinderSPD){rpc->ResetSPD();}
537 if(fModeSwitch==kClusterFinderSSD){rpc->ResetSSD();}
539 if(fModeSwitch==kClusterFinderSDD){rpc->ResetSDD();}
543 fRawReader->ClearBuffers();
545 UInt_t nClusters=fclusters.size();
547 UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
548 if( size + bufferSize > maxBufferSize ){
549 HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
555 RecPointToSpacePoint(outputPtr,size);
557 AliHLTComponentBlockData bd;
560 bd.fSize = bufferSize;
561 bd.fSpecification = iter->fSpecification;
562 bd.fDataType = GetOutputDataType();
563 outputBlocks.push_back( bd );
565 fBenchmark.AddOutput(bd.fSize);
573 HLTInfo(fBenchmark.GetStatistics());
577 int AliHLTITSClusterFinderComponent::Configure(const char* arguments)
582 if (!arguments) return iResult;
584 TString allArgs=arguments;
587 TObjArray* pTokens=allArgs.Tokenize(" ");
590 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
591 argument=((TObjString*)pTokens->At(i))->GetString();
592 if (argument.IsNull()) continue;
593 if (argument.CompareTo("-use-offline-finder")==0) {
594 fUseOfflineFinder = 1;
595 HLTInfo("Off-line ClusterFinder will be used");
599 else if (argument.CompareTo("")==0) {
605 HLTError("unknown argument %s", argument.Data());
616 int AliHLTITSClusterFinderComponent::Reconfigure(const char* cdbEntry, const char* chainId)
618 // see header file for class documentation
621 const char* path="HLT/ConfigITS/ClusterFinderComponent";
622 const char* defaultNotify="";
625 defaultNotify=" (default)";
628 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
629 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
631 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
633 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
634 iResult=Configure(pString->GetString().Data());
636 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
639 HLTError("can not fetch object \"%s\" from CDB", path);
645 void AliHLTITSClusterFinderComponent::GetOCDBObjectDescription( TMap* const targetMap)
647 // Get a list of OCDB object description.
648 if (!targetMap) return;
650 targetMap->Add(new TObjString("ITS/Calib/SPDNoisy"),new TObjString("Calibration object for SPD" ));
651 targetMap->Add(new TObjString("ITS/Calib/SPDDead"),new TObjString("Calibration object for SPD" ));
652 targetMap->Add(new TObjString("TRIGGER/SPD/PITConditions"),new TObjString("Calibration object for SPD" ));
654 targetMap->Add(new TObjString("ITS/Calib/CalibSDD"),new TObjString("Calibration object for SDD" ));
655 targetMap->Add(new TObjString("ITS/Calib/RespSDD"),new TObjString("Calibration object for SDD" ));
656 targetMap->Add(new TObjString("ITS/Calib/DriftSpeedSDD"),new TObjString("Calibration object for SDD" ));
657 targetMap->Add(new TObjString("ITS/Calib/DDLMapSDD"),new TObjString("Calibration object for SDD" ));
658 targetMap->Add(new TObjString("ITS/Calib/MapsTimeSDD"),new TObjString("Calibration object for SDD" ));
660 targetMap->Add(new TObjString("ITS/Calib/NoiseSSD"),new TObjString("Calibration object for SSD" ));
661 targetMap->Add(new TObjString("ITS/Calib/GainSSD"),new TObjString("Calibration object for SSD" ));
662 targetMap->Add(new TObjString("ITS/Calib/BadChannelsSSD"),new TObjString("Calibration object for SSD" ));
663 //General reconstruction
664 targetMap->Add(new TObjString("GRP/CTP/Scalers"),new TObjString("General reconstruction object" ));
668 void AliHLTITSClusterFinderComponent::RecPointToSpacePoint(AliHLTUInt8_t* outputPtr,AliHLTUInt32_t& size){
669 AliHLTITSClusterData *outputClusters = reinterpret_cast<AliHLTITSClusterData*>(outputPtr + size);
670 outputClusters->fSpacePointCnt=fclusters.size();
672 for(UInt_t i=0;i<fclusters.size();i++){
673 AliITSRecPoint *recpoint = (AliITSRecPoint*) &(fclusters[i]);
674 outputClusters->fSpacePoints[clustIdx].fY=recpoint->GetY();
675 outputClusters->fSpacePoints[clustIdx].fZ=recpoint->GetZ();
676 outputClusters->fSpacePoints[clustIdx].fSigmaY2=recpoint->GetSigmaY2();
677 outputClusters->fSpacePoints[clustIdx].fSigmaZ2=recpoint->GetSigmaZ2();
678 outputClusters->fSpacePoints[clustIdx].fSigmaYZ=recpoint->GetSigmaYZ();
679 outputClusters->fSpacePoints[clustIdx].fQ=recpoint->GetQ();
680 outputClusters->fSpacePoints[clustIdx].fNy=recpoint->GetNy();
681 outputClusters->fSpacePoints[clustIdx].fNz=recpoint->GetNz();
682 outputClusters->fSpacePoints[clustIdx].fLayer=recpoint->GetLayer();
683 outputClusters->fSpacePoints[clustIdx].fIndex=recpoint->GetDetectorIndex() | recpoint->GetPindex() | recpoint->GetNindex();
684 outputClusters->fSpacePoints[clustIdx].fTracks[0]=recpoint->GetLabel(0);
685 outputClusters->fSpacePoints[clustIdx].fTracks[1]=recpoint->GetLabel(1);
686 outputClusters->fSpacePoints[clustIdx].fTracks[2]=recpoint->GetLabel(2);