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 Ovrebekk <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 Ovrebekk <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"
45 #include "AliRunLoader.h"
46 #include "AliLoader.h"
52 #include "TObjString.h"
55 /** ROOT macro for the implementation of ROOT specific class methods */
56 ClassImp(AliHLTITSClusterFinderComponent);
58 AliHLTITSClusterFinderComponent::AliHLTITSClusterFinderComponent(int mode)
61 fInputDataType(kAliHLTVoidDataType),
62 fOutputDataType(kAliHLTVoidDataType),
82 fBenchmark(GetComponentID()),
84 fInputMultiplierDigits(20),
87 // see header file for class documentation
89 // refer to README to build package
91 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
94 case kClusterFinderSPD:
95 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;
96 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
98 case kClusterFinderSDD:
99 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;
100 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
102 case kClusterFinderSSD:
103 fInputDataType = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;
104 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
106 case kClusterFinderDigits:
107 fInputDataType = kAliHLTDataTypeAliTreeD|kAliHLTDataOriginITS;
108 fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITS;
111 HLTFatal("unknown cluster finder");
115 AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent()
117 // see header file for class documentation
125 // Public functions to implement AliHLTComponent's interface.
126 // These functions are required for the registration process
128 const char* AliHLTITSClusterFinderComponent::GetComponentID()
130 // see header file for class documentation
132 case kClusterFinderSPD:
133 return "ITSClusterFinderSPD";
135 case kClusterFinderSDD:
136 return "ITSClusterFinderSDD";
138 case kClusterFinderSSD:
139 return "ITSClusterFinderSSD";
141 case kClusterFinderDigits:
142 return "ITSClusterFinderDigits";
148 void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
150 // see header file for class documentation
152 list.push_back( fInputDataType );
155 AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType()
157 // see header file for class documentation
158 return fOutputDataType;
161 void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
162 // see header file for class documentation
163 constBase = fOutputSizeOffset;
165 case kClusterFinderDigits:
166 inputMultiplier = fInputMultiplierDigits;
168 case kClusterFinderSPD:
169 case kClusterFinderSDD:
170 case kClusterFinderSSD:
172 inputMultiplier = 20;
176 AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
177 // see header file for class documentation
178 return new AliHLTITSClusterFinderComponent(fModeSwitch);
181 Int_t AliHLTITSClusterFinderComponent::DoInit( int argc, const char** argv ) {
182 // see header file for class documentation
184 fBenchmark.SetTimer(0,"total");
185 fBenchmark.SetTimer(1,"reco");
195 Int_t runNo = GetRunNo();
196 AliCDBStorage* store = AliCDBManager::Instance()->GetDefaultStorage();
203 if(store->GetLatestVersion("ITS/Calib/SPDNoisy", runNo)<0){
204 HLTError("SPDNoisy is not found in SPD/Calib");
207 if(store->GetLatestVersion("ITS/Calib/SPDDead", runNo)<0){
208 HLTError("SPDDead is not found in SPD/Calib");
211 if(store->GetLatestVersion("TRIGGER/SPD/PITConditions", runNo)<0){
212 HLTError("PITConditions is not found in TRIGGER/SPD");
217 if(store->GetLatestVersion("ITS/Calib/CalibSDD", runNo)<0){
218 HLTError("CalibSDD is not found in ITS/Calib");
221 if(store->GetLatestVersion("ITS/Calib/RespSDD", runNo)<0){
222 HLTError("RespSDD is not found in ITS/Calib");
225 if(store->GetLatestVersion("ITS/Calib/DriftSpeedSDD", runNo)<0){
226 HLTError("DriftSpeedSDD is not found in ITS/Calib");
229 if(store->GetLatestVersion("ITS/Calib/DDLMapSDD", runNo)<0){
230 HLTError("DDLMapSDD is not found in ITS/Calib");
233 if(store->GetLatestVersion("ITS/Calib/MapsTimeSDD", runNo)<0){
234 HLTError("MapsTimeSDD is not found in ITS/Calib");
239 if(store->GetLatestVersion("ITS/Calib/NoiseSSD", runNo)<0){
240 HLTError("NoiseSSD is not found in ITS/Calib");
243 if(store->GetLatestVersion("ITS/Calib/GainSSD", runNo)<0){
244 HLTError("GainSSD is not found in ITS/Calib");
247 if(store->GetLatestVersion("ITS/Calib/BadChannelsSSD", runNo)<0){
248 HLTError("BadChannelsSSD is not found in ITS/Calib");
252 //General reconstruction
253 if(store->GetLatestVersion("GRP/CTP/Scalers", runNo)<0){
254 HLTError("Scalers is not found in GRP/CTP/");
257 if(!cdbOK){return -EIO;}
259 if(fModeSwitch==kClusterFinderSPD) {
260 HLTDebug("using ClusterFinder for SPD");
261 //fNModules=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1) + AliITSgeomTGeo::GetNDetectors(2)*AliITSgeomTGeo::GetNLadders(2);
262 fId=AliHLTDAQ::DdlIDOffset("ITSSPD");
263 fNddl=AliHLTDAQ::NumberOfDdls("ITSSPD");
265 else if(fModeSwitch==kClusterFinderSDD) {
266 HLTDebug("using ClusterFinder for SDD");
267 //fNModules=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3) + AliITSgeomTGeo::GetNDetectors(4)*AliITSgeomTGeo::GetNLadders(4);
268 fId=AliHLTDAQ::DdlIDOffset("ITSSDD");
269 fNddl=AliHLTDAQ::NumberOfDdls("ITSSDD");
271 else if(fModeSwitch==kClusterFinderSSD) {
272 HLTDebug("using ClusterFinder for SSD");
273 //fNModules=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5) + AliITSgeomTGeo::GetNDetectors(6)*AliITSgeomTGeo::GetNLadders(6);
274 fId=AliHLTDAQ::DdlIDOffset("ITSSSD");
275 fNddl=AliHLTDAQ::NumberOfDdls("ITSSSD");
277 else if(fModeSwitch==kClusterFinderDigits) {
281 HLTFatal("No mode set for clusterfindercomponent");
284 //Removed the warning for loading default RecoParam in HLT
285 AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
286 AliITSReconstructor *rec = new AliITSReconstructor();
287 rec->SetRecoParam(par);
289 AliCDBManager* man = AliCDBManager::Instance();
290 if (!man->IsDefaultStorageSet()){
291 HLTError("Default CDB storage has not been set !");
295 if(AliGeomManager::GetGeometry()==NULL){
296 AliGeomManager::LoadGeometry();
299 fgeomInit = new AliITSInitGeometry();
300 //fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
301 //fgeomInit->InitAliITSgeom(fgeom);
302 fgeom = fgeomInit->CreateAliITSgeom();
304 fNModules = fgeom->GetIndexMax();
306 for(Int_t i=0;i<6;i++)modperlay[i]=AliITSgeomTGeo::GetNDetectors(1+i)*AliITSgeomTGeo::GetNLadders(1+i);
307 fSPDNModules=modperlay[0]+modperlay[1];
308 fSDDNModules=modperlay[2]+modperlay[3];
309 fSSDNModules=modperlay[4]+modperlay[5];
311 if(fModeSwitch==kClusterFinderSPD) {
313 fLastModule=fSPDNModules;
315 else if(fModeSwitch==kClusterFinderSDD) {
316 fFirstModule=fSPDNModules;
317 fLastModule=fFirstModule + fSDDNModules;
319 else if(fModeSwitch==kClusterFinderSSD) {
320 fFirstModule=fSPDNModules + fSDDNModules;
321 fLastModule=fFirstModule + fSSDNModules;
325 fDettype = new AliITSDetTypeRec();
326 fDettype->SetITSgeom(fgeom);
327 fDettype->SetDefaults();
328 fDettype->SetDefaultClusterFindersV2(kTRUE);
333 fRawReader = new AliRawReaderMemory();
334 fSPD = new AliHLTITSClusterFinderSPD( fDettype );
335 fSSD = new AliHLTITSClusterFinderSSD( fDettype, fRawReader );
337 TString arguments = "";
338 for ( int i = 0; i < argc; i++ ) {
339 if ( !arguments.IsNull() ) arguments += " ";
340 arguments += argv[i];
346 return Configure( arguments.Data() );
349 Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
350 // see header file for class documentation
370 fUseOfflineFinder = 0;
373 fpLoader->UnloadDigits();
380 int AliHLTITSClusterFinderComponent::DoEvent
382 const AliHLTComponentEventData& evtData,
383 const AliHLTComponentBlockData* /*blocks*/,
384 AliHLTComponentTriggerData& /*trigData*/,
385 AliHLTUInt8_t* outputPtr,
386 AliHLTUInt32_t& size,
387 vector<AliHLTComponentBlockData>& outputBlocks )
389 // see header file for class documentation
391 AliHLTUInt32_t maxBufferSize = size;
392 size = 0; // output size
394 if (!IsDataEvent()) return 0;
396 if ( evtData.fBlockCnt<=0 )
398 HLTDebug("no blocks in event" );
402 AliHLTUInt32_t totalInputSize=0;
403 fBenchmark.StartNewEvent();
405 for( const AliHLTComponentBlockData *i= GetFirstInputBlock(fInputDataType); i!=NULL; i=GetNextInputBlock() ){
406 totalInputSize+=i->fSize;
407 fBenchmark.AddInput(i->fSize);
412 if(fModeSwitch==kClusterFinderDigits) {
414 for ( const TObject *iter = GetFirstInputObject(fInputDataType); iter != NULL; iter = GetNextInputObject() ) {
415 tD = dynamic_cast<TTree*>(const_cast<TObject*>( iter ) );
417 HLTFatal("No Digit Tree found");
420 // 2010-04-17 very crude workaround: TTree objects are difficult to send
421 // The actual case: Running ITS and TPC reconstruction fails at the second event
422 // to read the ITS digits from the TreeD
424 // Reason: reading fails in TBranch::GetBasket, there a new basket is opened from
425 // a TFile object. The function TBranch::GetFile returns the file object from
426 // an internal fDirectory (TDirectory) object. This file is at the second event
427 // set to the TPC.Digits.root. The internal mismatch creates a seg fault
429 // Investigation: TBranch::Streamer uses a crude assignment after creating the
431 // fDirectory = gDirectory;
432 // gDirectory is obviously not set correctly. Setting the directory to a TFile
433 // object for the ITS digits helps to fix the internal mess. Tried also to set
434 // the Directory for the TreeD to NULL (This has only effect if ones sets it
435 // to something not NULL first, and then to NULL). But then no content, i.e.
436 // ITS clusters could be retrieved.
438 // Conclusion: TTree objects are hardly to be sent via TMessage, there are direct
439 // links to the file required anyhow.
441 // 2011-01-28 hotfix reloaded: accessing the files like that fails if there are
442 // multiple digit files because of a large number of events. New ugly fix is to
443 // use the global runloader instance to get hold on the digits tree.
445 AliRunLoader* pRunLoader=AliRunLoader::Instance();
447 HLTError("failed to get global runloader instance");
450 // get the specific loader for the module
452 const char* loaderType="ITSLoader";
453 fpLoader=pRunLoader->GetLoader(loaderType);
455 HLTError("can not get loader \"%s\" from runloader", loaderType);
458 // prepare the loader
459 fpLoader->LoadDigits("read");
461 pRunLoader->GetEvent(GetEventCount());
463 tD=fpLoader->TreeD();
466 fDettype->SetTreeAddressD(tD);
467 fDettype->MakeBranch(tR,"R");
468 fDettype->SetTreeAddressR(tR);
471 fDettype->DigitsToRecPoints(tD,tR,0,opt,0);
473 TClonesArray * fRecPoints = NULL;
474 tR->SetBranchAddress("ITSRecPoints",&fRecPoints);
475 for(Int_t treeEntry=0;treeEntry<tR->GetEntries();treeEntry++){
476 tR->GetEntry(treeEntry);
477 fnClusters += fRecPoints->GetEntries();
480 UInt_t bufferSize = fnClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
481 if( size + bufferSize > maxBufferSize ){
482 //HLTWarning( "Output buffer size exceed (buffer size %d, required size %d)", maxBufferSize, size+bufferSize);
483 if (totalInputSize>0) {
484 fInputMultiplierDigits=(float)(size+bufferSize)/totalInputSize;
485 fInputMultiplierDigits+=1.;
487 fOutputSizeOffset=totalInputSize;
488 fInputMultiplierDigits=1.;
495 RecPointToSpacePoint(outputPtr,size);
497 AliHLTComponentBlockData bd;
500 bd.fSize = bufferSize;
501 bd.fSpecification = 0x00000000;
502 bd.fDataType = GetOutputDataType();
503 outputBlocks.push_back( bd );
505 fBenchmark.AddOutput(bd.fSize);
511 AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
513 // -- Loop over blocks
514 for( const AliHLTComponentBlockData* iter = GetFirstInputBlock(fInputDataType); iter != NULL; iter = GetNextInputBlock() ) {
516 if(fUseOfflineFinder){
517 if(fModeSwitch==kClusterFinderSPD){rpc->ResetSPD();}
518 if(fModeSwitch==kClusterFinderSSD){rpc->ResetSSD();}
520 if(fModeSwitch==kClusterFinderSDD){rpc->ResetSDD();}
522 // -- Debug output of datatype --
523 HLTDebug("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
524 evtData.fEventID, evtData.fEventID,
525 DataType2Text(iter->fDataType).c_str(),
526 DataType2Text(fInputDataType).c_str());
528 // -- Check for the correct data type
529 if ( iter->fDataType != (fInputDataType) )
532 // -- Get equipment ID out of specification
533 AliHLTUInt32_t spec = iter->fSpecification;
536 for ( Int_t ii = 0; ii < fNddl ; ii++ ) { //number of ddl's
537 if ( spec & 0x00000001 ) {
544 // -- Set equipment ID to the raw reader
547 HLTWarning("The fRawReader pointer is NULL");
551 if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
552 HLTWarning("Could not add buffer");
559 if(fModeSwitch==kClusterFinderSPD && !fUseOfflineFinder){ fSPD->RawdataToClusters( fRawReader, fclusters ); }
560 else if(fModeSwitch==kClusterFinderSSD && !fUseOfflineFinder){ fSSD->RawdataToClusters( fclusters ); }
562 if(fModeSwitch==kClusterFinderSPD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,"SPD");}
563 if(fModeSwitch==kClusterFinderSSD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,"SSD");}
564 if(fModeSwitch==kClusterFinderSDD) {fDettype->DigitsToRecPoints(fRawReader,"SDD");}
565 TClonesArray* clusters = NULL;
566 for(int i=fFirstModule;i<fLastModule;i++){
567 clusters = rpc->UncheckedGetClusters(i);
568 if(clusters != NULL){
569 fnClusters += clusters->GetEntriesFast();
575 if(fModeSwitch==kClusterfinderSPD){
576 fnClusters = rpc->GetNClustersInLayerFast(1) + rpc->GetNClustersInLayerFast(2);
578 if(fModeSwitch==kClusterfinderSDD){
579 fnClusters = rpc->GetNClustersInLayerFast(3) + rpc->GetNClustersInLayerFast(4);
581 if(fModeSwitch==kClusterfinderSSD){
582 fnClusters = rpc->GetNClustersInLayerFast(5) + rpc->GetNClustersInLayerFast(6);
588 fRawReader->ClearBuffers();
590 UInt_t nClusters=fclusters.size();
591 if(nClusters>0){fnClusters = nClusters;}
593 UInt_t bufferSize = fnClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
594 if( size + bufferSize > maxBufferSize ){
595 HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
601 RecPointToSpacePoint(outputPtr,size);
603 AliHLTComponentBlockData bd;
606 bd.fSize = bufferSize;
607 bd.fSpecification = iter->fSpecification;
608 bd.fDataType = GetOutputDataType();
609 outputBlocks.push_back( bd );
611 fBenchmark.AddOutput(bd.fSize);
612 if(nClusters>0){fclusters.clear();}
619 HLTInfo(fBenchmark.GetStatistics());
623 int AliHLTITSClusterFinderComponent::Configure(const char* arguments)
628 if (!arguments) return iResult;
630 TString allArgs=arguments;
633 TObjArray* pTokens=allArgs.Tokenize(" ");
636 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
637 argument=((TObjString*)pTokens->At(i))->GetString();
638 if (argument.IsNull()) continue;
639 if (argument.CompareTo("-use-offline-finder")==0) {
640 fUseOfflineFinder = 1;
641 HLTInfo("Off-line ClusterFinder will be used");
645 else if (argument.CompareTo("")==0) {
651 HLTError("unknown argument %s", argument.Data());
662 int AliHLTITSClusterFinderComponent::Reconfigure(const char* cdbEntry, const char* chainId)
664 // see header file for class documentation
670 case kClusterFinderSPD:
671 path = "HLT/ConfigITS/ITSClusterFinderSPD";
673 case kClusterFinderSDD:
674 path = "HLT/ConfigITS/ITSClusterFinderSDD";
676 case kClusterFinderSSD:
677 path = "HLT/ConfigITS/ITSClusterFinderSSD";
679 case kClusterFinderDigits:
683 HLTFatal("unknown cluster finder");
686 const char* defaultNotify="";
689 defaultNotify=" (default)";
692 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
693 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
695 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
697 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
698 iResult=Configure(pString->GetString().Data());
700 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
703 HLTError("can not fetch object \"%s\" from CDB", path);
709 void AliHLTITSClusterFinderComponent::GetOCDBObjectDescription( TMap* const targetMap)
711 // Get a list of OCDB object description.
712 if (!targetMap) return;
714 targetMap->Add(new TObjString("ITS/Calib/SPDNoisy"),new TObjString("Calibration object for SPD" ));
715 targetMap->Add(new TObjString("ITS/Calib/SPDDead"),new TObjString("Calibration object for SPD" ));
716 targetMap->Add(new TObjString("TRIGGER/SPD/PITConditions"),new TObjString("Calibration object for SPD" ));
718 targetMap->Add(new TObjString("ITS/Calib/CalibSDD"),new TObjString("Calibration object for SDD" ));
719 targetMap->Add(new TObjString("ITS/Calib/RespSDD"),new TObjString("Calibration object for SDD" ));
720 targetMap->Add(new TObjString("ITS/Calib/DriftSpeedSDD"),new TObjString("Calibration object for SDD" ));
721 targetMap->Add(new TObjString("ITS/Calib/DDLMapSDD"),new TObjString("Calibration object for SDD" ));
722 targetMap->Add(new TObjString("ITS/Calib/MapsTimeSDD"),new TObjString("Calibration object for SDD" ));
724 targetMap->Add(new TObjString("ITS/Calib/NoiseSSD"),new TObjString("Calibration object for SSD" ));
725 targetMap->Add(new TObjString("ITS/Calib/GainSSD"),new TObjString("Calibration object for SSD" ));
726 targetMap->Add(new TObjString("ITS/Calib/BadChannelsSSD"),new TObjString("Calibration object for SSD" ));
727 //General reconstruction
728 targetMap->Add(new TObjString("GRP/CTP/Scalers"),new TObjString("General reconstruction object" ));
732 void AliHLTITSClusterFinderComponent::RecPointToSpacePoint(AliHLTUInt8_t* outputPtr,AliHLTUInt32_t& size){
733 AliHLTITSClusterData *outputClusters = reinterpret_cast<AliHLTITSClusterData*>(outputPtr + size);
734 outputClusters->fSpacePointCnt=fnClusters;
736 if(fModeSwitch==kClusterFinderDigits) {
737 TClonesArray * fRecPoints = NULL;
738 tR->SetBranchAddress("ITSRecPoints",&fRecPoints);
739 for(Int_t treeEntry=0;treeEntry<tR->GetEntries();treeEntry++){
740 tR->GetEntry(treeEntry);
741 fnClusters += fRecPoints->GetEntries();
742 for(Int_t tCloneEntry=0;tCloneEntry<fRecPoints->GetEntries();tCloneEntry++){
743 AliITSRecPoint *recpoint=(AliITSRecPoint*)fRecPoints->At(tCloneEntry);
744 RecpointToOutput(outputClusters,recpoint,clustIdx);
748 else if(fclusters.size()>0){
749 for(UInt_t i=0;i<fclusters.size();i++){
750 AliITSRecPoint *recpoint = (AliITSRecPoint*) &(fclusters[i]);
751 RecpointToOutput(outputClusters,recpoint,clustIdx);
755 AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
756 TClonesArray* clusters = NULL;
757 for(Int_t i=fFirstModule;i<fLastModule;i++){
758 clusters = rpc->UncheckedGetClusters(i);
759 for(Int_t j=0;j<clusters->GetEntriesFast();j++){
760 AliITSRecPoint *recpoint = (AliITSRecPoint*) clusters->At(j);
761 RecpointToOutput(outputClusters,recpoint,clustIdx);
767 void AliHLTITSClusterFinderComponent::RecpointToOutput(AliHLTITSClusterData *outputClusters, AliITSRecPoint *recpoint, int &clustIdx){
768 outputClusters->fSpacePoints[clustIdx].fY=recpoint->GetY();
769 outputClusters->fSpacePoints[clustIdx].fZ=recpoint->GetZ();
770 outputClusters->fSpacePoints[clustIdx].fSigmaY2=recpoint->GetSigmaY2();
771 outputClusters->fSpacePoints[clustIdx].fSigmaZ2=recpoint->GetSigmaZ2();
772 outputClusters->fSpacePoints[clustIdx].fSigmaYZ=recpoint->GetSigmaYZ();
773 outputClusters->fSpacePoints[clustIdx].fQ=recpoint->GetQ();
774 outputClusters->fSpacePoints[clustIdx].fNy=recpoint->GetNy();
775 outputClusters->fSpacePoints[clustIdx].fNz=recpoint->GetNz();
776 outputClusters->fSpacePoints[clustIdx].fLayer=recpoint->GetLayer();
777 outputClusters->fSpacePoints[clustIdx].fIndex=recpoint->GetDetectorIndex() | recpoint->GetPindex() | recpoint->GetNindex();
778 outputClusters->fSpacePoints[clustIdx].fTracks[0]=recpoint->GetLabel(0);
779 outputClusters->fSpacePoints[clustIdx].fTracks[1]=recpoint->GetLabel(1);
780 outputClusters->fSpacePoints[clustIdx].fTracks[2]=recpoint->GetLabel(2);