]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSClusterFinderComponent.cxx
extra benchmarks are added
[u/mrichter/AliRoot.git] / HLT / ITS / AliHLTITSClusterFinderComponent.cxx
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
25 using namespace std;
26 #endif
27
28 #include "AliHLTITSClusterFinderComponent.h" 
29
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"
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 */
51 ClassImp(AliHLTITSClusterFinderComponent);
52
53 AliHLTITSClusterFinderComponent::AliHLTITSClusterFinderComponent(int mode)
54   :
55   fModeSwitch(mode),
56   fInputDataType(kAliHLTVoidDataType),
57   fOutputDataType(kAliHLTVoidDataType),
58   fUseOfflineFinder(0),
59   fNModules(0),
60   fId(0),
61   fNddl(0),
62   fClusters(NULL),
63   fRawReader(NULL),
64   fDettype(NULL),
65   fgeom(NULL),
66   fgeomInit(NULL),
67   fSPD(NULL),
68   fSSD(NULL),
69   tD(NULL),
70   tR(NULL),
71   fclusters(),
72   fBenchmark(GetComponentID())
73
74   // see header file for class documentation
75   // or
76   // refer to README to build package
77   // or
78   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
79
80   switch(fModeSwitch){
81   case kClusterFinderSPD:
82     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;
83     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
84     break;
85   case kClusterFinderSDD:        
86     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;
87     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
88     break;
89   case kClusterFinderSSD:
90     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;
91     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
92     break;
93   case kClusterFinderDigits:
94     fInputDataType  = kAliHLTDataTypeAliTreeD|kAliHLTDataOriginITS;
95     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITS;
96     break;
97   default:
98     HLTFatal("unknown cluster finder");
99   }
100 }
101
102 AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent() 
103 {
104   // see header file for class documentation
105   delete fRawReader;
106   delete fDettype;
107   delete fgeomInit;  
108   delete fSPD;
109   delete fSSD;
110 }
111
112 // Public functions to implement AliHLTComponent's interface.
113 // These functions are required for the registration process
114
115 const char* AliHLTITSClusterFinderComponent::GetComponentID()
116 {
117   // see header file for class documentation
118   switch(fModeSwitch){
119   case kClusterFinderSPD:
120     return "ITSClusterFinderSPD";
121     break;
122   case kClusterFinderSDD:        
123     return "ITSClusterFinderSDD";
124     break;
125   case kClusterFinderSSD:
126     return "ITSClusterFinderSSD";
127     break;
128   case kClusterFinderDigits:
129     return "ITSClusterFinderDigits";
130     break;
131   }
132   return "";
133 }
134
135 void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) 
136 {
137   // see header file for class documentation
138   list.clear(); 
139   list.push_back( fInputDataType );
140 }
141
142 AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType() 
143 {
144   // see header file for class documentation
145   return fOutputDataType;
146 }
147
148 void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
149   // see header file for class documentation
150   constBase = 0;
151   inputMultiplier = 100;
152 }
153
154 AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
155   // see header file for class documentation
156   return new AliHLTITSClusterFinderComponent(fModeSwitch);
157 }
158         
159 Int_t AliHLTITSClusterFinderComponent::DoInit( int argc, const char** argv ) {
160   // see header file for class documentation
161   fBenchmark.Reset();
162   fBenchmark.SetTimer(0,"total");
163   fBenchmark.SetTimer(1,"reco");
164
165   /*
166   fStatTime = 0;
167   fStatTimeAll = 0;
168   fStatTimeC = 0;
169   fStatTimeAllC = 0;
170   fStatNEv = 0;
171   */
172   
173   Int_t runNo = GetRunNo();
174   AliCDBStorage* store = AliCDBManager::Instance()->GetDefaultStorage();
175   if (!store) {
176     return NULL;
177   }
178
179   bool cdbOK = true;
180   //OCDB for SPD
181   if(store->GetLatestVersion("ITS/Calib/SPDNoisy", runNo)<0){
182     HLTError("SPDNoisy is not found in SPD/Calib");
183     cdbOK = false;
184   }
185   if(store->GetLatestVersion("ITS/Calib/SPDDead", runNo)<0){
186     HLTError("SPDDead is not found in SPD/Calib");
187     cdbOK = false;
188   }
189   if(store->GetLatestVersion("TRIGGER/SPD/PITConditions", runNo)<0){
190     HLTError("PITConditions is not found in TRIGGER/SPD");
191     cdbOK = false;
192   }
193   
194   //OCDB for SDD
195   if(store->GetLatestVersion("ITS/Calib/CalibSDD", runNo)<0){
196     HLTError("CalibSDD is not found in ITS/Calib");
197     cdbOK = false;
198   }
199   if(store->GetLatestVersion("ITS/Calib/RespSDD", runNo)<0){
200     HLTError("RespSDD is not found in ITS/Calib");
201     cdbOK = false;
202   }
203   if(store->GetLatestVersion("ITS/Calib/DriftSpeedSDD", runNo)<0){
204     HLTError("DriftSpeedSDD is not found in ITS/Calib");
205     cdbOK = false;
206   }
207   if(store->GetLatestVersion("ITS/Calib/DDLMapSDD", runNo)<0){
208     HLTError("DDLMapSDD is not found in ITS/Calib");
209     cdbOK = false;
210   }
211   if(store->GetLatestVersion("ITS/Calib/MapsTimeSDD", runNo)<0){
212     HLTError("MapsTimeSDD is not found in ITS/Calib");
213     cdbOK = false;
214   }
215
216   //OCDB for SSD
217   if(store->GetLatestVersion("ITS/Calib/NoiseSSD", runNo)<0){
218     HLTError("NoiseSSD is not found in ITS/Calib");
219     cdbOK = false;
220   }
221   if(store->GetLatestVersion("ITS/Calib/GainSSD", runNo)<0){
222     HLTError("GainSSD is not found in ITS/Calib");
223     cdbOK = false;
224   }
225   if(store->GetLatestVersion("ITS/Calib/BadChannelsSSD", runNo)<0){
226     HLTError("BadChannelsSSD is not found in ITS/Calib");
227     cdbOK = false;
228   }
229   
230   //General reconstruction
231   if(store->GetLatestVersion("GRP/CTP/Scalers", runNo)<0){
232     HLTError("Scalers is not found in GRP/CTP/");
233     cdbOK = false;
234   }
235   if(!cdbOK){return NULL;}
236
237   if(fModeSwitch==kClusterFinderSPD) {
238     HLTDebug("using ClusterFinder for SPD");
239     //fNModules=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1) + AliITSgeomTGeo::GetNDetectors(2)*AliITSgeomTGeo::GetNLadders(2);
240     fId=AliHLTDAQ::DdlIDOffset("ITSSPD");
241     fNddl=AliHLTDAQ::NumberOfDdls("ITSSPD");
242   }
243   else if(fModeSwitch==kClusterFinderSDD) {
244     HLTDebug("using ClusterFinder for SDD");
245     //fNModules=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3) + AliITSgeomTGeo::GetNDetectors(4)*AliITSgeomTGeo::GetNLadders(4);
246     fId=AliHLTDAQ::DdlIDOffset("ITSSDD");
247     fNddl=AliHLTDAQ::NumberOfDdls("ITSSDD");
248   }
249   else if(fModeSwitch==kClusterFinderSSD) {
250     HLTDebug("using ClusterFinder for SSD");
251     //fNModules=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5) + AliITSgeomTGeo::GetNDetectors(6)*AliITSgeomTGeo::GetNLadders(6);
252     fId=AliHLTDAQ::DdlIDOffset("ITSSSD");
253     fNddl=AliHLTDAQ::NumberOfDdls("ITSSSD");
254   }
255   else if(fModeSwitch==kClusterFinderDigits) {
256     //tR = new TTree();
257   }
258   else{
259      HLTFatal("No mode set for clusterfindercomponent");
260   }
261
262   //Removed the warning for loading default RecoParam in HLT
263   AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
264   AliITSReconstructor *rec = new AliITSReconstructor();
265   rec->SetRecoParam(par);
266   
267   AliCDBManager* man = AliCDBManager::Instance();
268   if (!man->IsDefaultStorageSet()){
269     HLTError("Default CDB storage has not been set !");
270     return -ENOENT;
271   }
272
273   if(AliGeomManager::GetGeometry()==NULL){
274     AliGeomManager::LoadGeometry();
275   }
276
277  //fgeomInit = new AliITSInitGeometry(kvSPD02,2);
278   fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
279   //fgeomInit->InitAliITSgeom(fgeom);
280   fgeom = fgeomInit->CreateAliITSgeom();
281  
282   fNModules = fgeom->GetIndexMax();
283
284   fClusters = new TClonesArray*[fNModules]; 
285   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
286     fClusters[iModule] = NULL;
287   } 
288
289   //set dettype
290   fDettype = new AliITSDetTypeRec();
291   fDettype->SetITSgeom(fgeom); 
292   fDettype->SetDefaults();
293   fDettype->SetDefaultClusterFindersV2(kTRUE); 
294
295   if ( fRawReader )
296     return -EINPROGRESS;
297
298   fRawReader = new AliRawReaderMemory();
299   fSPD = new AliHLTITSClusterFinderSPD( fDettype );
300   fSSD = new AliHLTITSClusterFinderSSD( fDettype, fRawReader );
301
302   TString arguments = "";
303   for ( int i = 0; i < argc; i++ ) {
304     if ( !arguments.IsNull() ) arguments += " ";
305     arguments += argv[i];
306   }
307
308   tD = NULL;
309   tR = NULL;
310
311   return Configure( arguments.Data() );
312 }
313
314 Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
315   // see header file for class documentation
316
317   if ( fRawReader )
318     delete fRawReader;
319   fRawReader = NULL;
320
321   if ( fDettype )
322     delete fDettype;
323   fDettype = NULL;
324
325   if ( fgeomInit )
326     delete fgeomInit;
327   fgeomInit = NULL;
328   
329   delete fSPD;
330   fSPD = 0;
331
332   delete fSSD;
333   fSSD = 0;
334
335   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
336     if(fClusters[iModule] != NULL){
337       fClusters[iModule]->Delete();
338       delete fClusters[iModule];
339     }
340     fClusters[iModule] = NULL;
341   } 
342   
343   fUseOfflineFinder = 0;
344
345   return 0;
346 }
347
348 int AliHLTITSClusterFinderComponent::DoEvent
349 (
350  const AliHLTComponentEventData& evtData,
351  const AliHLTComponentBlockData* /*blocks*/,
352  AliHLTComponentTriggerData& /*trigData*/,
353  AliHLTUInt8_t* outputPtr,
354  AliHLTUInt32_t& size,
355   vector<AliHLTComponentBlockData>& outputBlocks )
356 {
357  // see header file for class documentation
358   
359   AliHLTUInt32_t maxBufferSize = size;
360   size = 0; // output size
361
362   if (!IsDataEvent()) return 0;
363   
364   if ( evtData.fBlockCnt<=0 )
365     {
366       HLTDebug("no blocks in event" );
367       return 0;
368     }
369
370   fBenchmark.StartNewEvent();
371   fBenchmark.Start(0);
372   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(fInputDataType); i!=NULL; i=GetNextInputBlock() ){
373     fBenchmark.AddInput(i->fSize);
374   }
375
376   Int_t ret = 0;
377
378   if(fModeSwitch==kClusterFinderDigits) {
379
380     for ( const TObject *iter = GetFirstInputObject(fInputDataType); iter != NULL; iter = GetNextInputObject() ) {  
381       tD = dynamic_cast<TTree*>(const_cast<TObject*>( iter ) );
382       if(!tD){
383         HLTFatal("No Digit Tree found");
384         return -1;
385       }
386       tR = new TTree();
387       fDettype->SetTreeAddressD(tD);
388       fDettype->MakeBranch(tR,"R");
389       fDettype->SetTreeAddressR(tR);
390       Option_t *opt="All";
391       fBenchmark.Start(1);
392       fDettype->DigitsToRecPoints(tD,tR,0,opt,kTRUE);
393       fBenchmark.Stop(1);
394       TClonesArray * fRecPoints;
395       tR->SetBranchAddress("ITSRecPoints",&fRecPoints);
396       for(Int_t treeEntry=0;treeEntry<tR->GetEntries();treeEntry++){
397         tR->GetEntry(treeEntry);
398         for(Int_t tCloneEntry=0;tCloneEntry<fRecPoints->GetEntries();tCloneEntry++){
399           AliITSRecPoint *recpoint=(AliITSRecPoint*)fRecPoints->At(tCloneEntry);
400           fclusters.push_back(*recpoint);
401         }
402       }
403       
404       if(tR){
405         tR->Delete();
406       }
407       UInt_t nClusters=fclusters.size();
408       
409       UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
410       if( size + bufferSize > maxBufferSize ){
411         HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
412         ret = -ENOSPC;      
413         break;          
414       }
415       if( nClusters>0 ){
416         fBenchmark.Start(1);
417         RecPointToSpacePoint(outputPtr,size);
418         fBenchmark.Stop(1);
419         AliHLTComponentBlockData bd;
420         FillBlockData( bd );
421         bd.fOffset = size;
422         bd.fSize = bufferSize;
423         bd.fSpecification = 0x00000000;
424         bd.fDataType = GetOutputDataType();
425         outputBlocks.push_back( bd );
426         size += bufferSize;
427         fBenchmark.AddOutput(bd.fSize);
428         fclusters.clear();      
429       }
430     }
431   }
432   else{
433       
434     // -- Loop over blocks
435     for( const AliHLTComponentBlockData* iter = GetFirstInputBlock(fInputDataType); iter != NULL; iter = GetNextInputBlock() ) {
436       
437       // -- Debug output of datatype --
438       HLTDebug("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
439                evtData.fEventID, evtData.fEventID, 
440                DataType2Text(iter->fDataType).c_str(), 
441                DataType2Text(fInputDataType).c_str());
442       
443       // -- Check for the correct data type
444       if ( iter->fDataType != (fInputDataType) )  
445         continue;
446       
447       // -- Get equipment ID out of specification
448       AliHLTUInt32_t spec = iter->fSpecification;
449       
450       Int_t id = fId;
451       for ( Int_t ii = 0; ii < fNddl ; ii++ ) {   //number of ddl's
452         if ( spec & 0x00000001 ) {
453           id += ii;
454           break;
455         }
456         spec = spec >> 1 ;
457       }
458       
459       // -- Set equipment ID to the raw reader
460       
461       if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
462         HLTWarning("Could not add buffer");
463       }
464
465       fBenchmark.Start(1);
466
467       if(fModeSwitch==kClusterFinderSPD && !fUseOfflineFinder){ fSPD->RawdataToClusters( fRawReader, fclusters ); }
468       else if(fModeSwitch==kClusterFinderSSD && !fUseOfflineFinder){ fSSD->RawdataToClusters( fclusters ); }
469       else{
470         if(fModeSwitch==kClusterFinderSPD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SPD");}
471         if(fModeSwitch==kClusterFinderSSD && fUseOfflineFinder) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SSD");}
472         if(fModeSwitch==kClusterFinderSDD) {fDettype->DigitsToRecPoints(fRawReader,fClusters,"SDD");}
473         for(int i=0;i<fNModules;i++){
474           if(fClusters[i] != NULL){
475             for(int j=0;j<fClusters[i]->GetEntriesFast();j++){
476               AliITSRecPoint *recpoint = (AliITSRecPoint*) (fClusters[i]->At(j));
477               fclusters.push_back(*recpoint);
478             }
479             fClusters[i]->Delete();
480             delete fClusters[i];
481           }
482           fClusters[i] = NULL;
483         }     
484       }
485       fBenchmark.Stop(1);
486       
487       fRawReader->ClearBuffers();    
488          
489       UInt_t nClusters=fclusters.size();
490       
491       UInt_t bufferSize = nClusters * sizeof(AliHLTITSSpacePointData) + sizeof(AliHLTITSClusterData);
492       if( size + bufferSize > maxBufferSize ){
493         HLTWarning( "Output buffer size exceed (buffer size %d, current size %d)", maxBufferSize, size+bufferSize);
494         ret = -ENOSPC;      
495         break;          
496       }
497       if( nClusters>0 ){
498
499         RecPointToSpacePoint(outputPtr,size);
500
501         AliHLTComponentBlockData bd;
502         FillBlockData( bd );
503         bd.fOffset = size;
504         bd.fSize = bufferSize;
505         bd.fSpecification = iter->fSpecification;
506         bd.fDataType = GetOutputDataType();
507         outputBlocks.push_back( bd );
508         size += bufferSize;
509         fBenchmark.AddOutput(bd.fSize);
510         fclusters.clear();      
511       }
512       
513     } // input blocks
514   }
515
516   fBenchmark.Stop(0);
517   HLTInfo(fBenchmark.GetStatistics());
518   return ret;
519 }
520
521 int AliHLTITSClusterFinderComponent::Configure(const char* arguments)
522 {
523   
524   int iResult=0;
525   
526   if (!arguments) return iResult;
527   
528   TString allArgs=arguments;
529   TString argument;
530   
531   TObjArray* pTokens=allArgs.Tokenize(" ");
532   
533   if (pTokens) {
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");
540         continue;
541       }
542       /*
543       else if (argument.CompareTo("")==0) {
544         HLTInfo("");
545         continue;
546       }
547       */
548       else {
549         HLTError("unknown argument %s", argument.Data());
550         iResult=-EINVAL;
551         break;
552       }
553     }
554     delete pTokens;
555   }
556   
557   return iResult;
558 }
559
560 int AliHLTITSClusterFinderComponent::Reconfigure(const char* cdbEntry, const char* chainId)
561 {
562   // see header file for class documentation
563   int iResult=0;
564   
565   const char* path="HLT/ConfigITS/ClusterFinderComponent";
566   const char* defaultNotify="";
567   if (cdbEntry) {
568     path=cdbEntry;
569     defaultNotify=" (default)";
570   }
571   if (path) {
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);
574     if (pEntry) {
575       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
576       if (pString) {
577         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
578         iResult=Configure(pString->GetString().Data());
579       } else {
580         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
581       }
582     } else {
583       HLTError("can not fetch object \"%s\" from CDB", path);
584     }
585   }
586   
587   return iResult;
588 }
589 void AliHLTITSClusterFinderComponent::RecPointToSpacePoint(AliHLTUInt8_t* outputPtr,AliHLTUInt32_t& size){
590   AliHLTITSClusterData *outputClusters = reinterpret_cast<AliHLTITSClusterData*>(outputPtr + size);
591   outputClusters->fSpacePointCnt=fclusters.size();    
592   int clustIdx=0;
593   for(UInt_t 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);      
608     clustIdx++;
609   }
610 }