]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSClusterFinderComponent.cxx
Added test for OCDB objects in DoInit and updated documentation (Gaute)
[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
73   // see header file for class documentation
74   // or
75   // refer to README to build package
76   // or
77   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
78
79   switch(fModeSwitch){
80   case kClusterFinderSPD:
81     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSPD;
82     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD;
83     break;
84   case kClusterFinderSDD:        
85     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSDD;
86     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD;
87     break;
88   case kClusterFinderSSD:
89     fInputDataType  = kAliHLTDataTypeDDLRaw | kAliHLTDataOriginITSSSD;
90     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD;
91     break;
92   case kClusterFinderDigits:
93     fInputDataType  = kAliHLTDataTypeAliTreeD|kAliHLTDataOriginITS;
94     fOutputDataType = kAliHLTDataTypeClusters|kAliHLTDataOriginITS;
95     break;
96   default:
97     HLTFatal("unknown cluster finder");
98   }
99 }
100
101 AliHLTITSClusterFinderComponent::~AliHLTITSClusterFinderComponent() 
102 {
103   // see header file for class documentation
104   delete fRawReader;
105   delete fDettype;
106   delete fgeomInit;  
107   delete fSPD;
108   delete fSSD;
109 }
110
111 // Public functions to implement AliHLTComponent's interface.
112 // These functions are required for the registration process
113
114 const char* AliHLTITSClusterFinderComponent::GetComponentID()
115 {
116   // see header file for class documentation
117   switch(fModeSwitch){
118   case kClusterFinderSPD:
119     return "ITSClusterFinderSPD";
120     break;
121   case kClusterFinderSDD:        
122     return "ITSClusterFinderSDD";
123     break;
124   case kClusterFinderSSD:
125     return "ITSClusterFinderSSD";
126     break;
127   case kClusterFinderDigits:
128     return "ITSClusterFinderDigits";
129     break;
130   }
131   return "";
132 }
133
134 void AliHLTITSClusterFinderComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) 
135 {
136   // see header file for class documentation
137   list.clear(); 
138   list.push_back( fInputDataType );
139 }
140
141 AliHLTComponentDataType AliHLTITSClusterFinderComponent::GetOutputDataType() 
142 {
143   // see header file for class documentation
144   return fOutputDataType;
145 }
146
147 void AliHLTITSClusterFinderComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
148   // see header file for class documentation
149   constBase = 0;
150   inputMultiplier = 100;
151 }
152
153 AliHLTComponent* AliHLTITSClusterFinderComponent::Spawn() {
154   // see header file for class documentation
155   return new AliHLTITSClusterFinderComponent(fModeSwitch);
156 }
157         
158 Int_t AliHLTITSClusterFinderComponent::DoInit( int argc, const char** argv ) {
159   // see header file for class documentation
160   /*
161   fStatTime = 0;
162   fStatTimeAll = 0;
163   fStatTimeC = 0;
164   fStatTimeAllC = 0;
165   fStatNEv = 0;
166   */
167   
168   Int_t runNo = GetRunNo();
169   AliCDBStorage* store = AliCDBManager::Instance()->GetDefaultStorage();
170   if (!store) {
171     return NULL;
172   }
173
174   bool cdbOK = true;
175   //OCDB for SPD
176   if(store->GetLatestVersion("ITS/Calib/SPDNoisy", runNo)<0){
177     HLTError("SPDNoisy is not found in SPD/Calib");
178     cdbOK = false;
179   }
180   if(store->GetLatestVersion("ITS/Calib/SPDDead", runNo)<0){
181     HLTError("SPDDead is not found in SPD/Calib");
182     cdbOK = false;
183   }
184   if(store->GetLatestVersion("TRIGGER/SPD/PITConditions", runNo)<0){
185     HLTError("PITConditions is not found in TRIGGER/SPD");
186     cdbOK = false;
187   }
188   
189   //OCDB for SDD
190   if(store->GetLatestVersion("ITS/Calib/CalibSDD", runNo)<0){
191     HLTError("CalibSDD is not found in ITS/Calib");
192     cdbOK = false;
193   }
194   if(store->GetLatestVersion("ITS/Calib/RespSDD", runNo)<0){
195     HLTError("RespSDD is not found in ITS/Calib");
196     cdbOK = false;
197   }
198   if(store->GetLatestVersion("ITS/Calib/DriftSpeedSDD", runNo)<0){
199     HLTError("DriftSpeedSDD is not found in ITS/Calib");
200     cdbOK = false;
201   }
202   if(store->GetLatestVersion("ITS/Calib/DDLMapSDD", runNo)<0){
203     HLTError("DDLMapSDD is not found in ITS/Calib");
204     cdbOK = false;
205   }
206   if(store->GetLatestVersion("ITS/Calib/MapsTimeSDD", runNo)<0){
207     HLTError("MapsTimeSDD is not found in ITS/Calib");
208     cdbOK = false;
209   }
210
211   //OCDB for SSD
212   if(store->GetLatestVersion("ITS/Calib/NoiseSSD", runNo)<0){
213     HLTError("NoiseSSD is not found in ITS/Calib");
214     cdbOK = false;
215   }
216   if(store->GetLatestVersion("ITS/Calib/GainSSD", runNo)<0){
217     HLTError("GainSSD is not found in ITS/Calib");
218     cdbOK = false;
219   }
220   if(store->GetLatestVersion("ITS/Calib/BadChannelsSSD", runNo)<0){
221     HLTError("BadChannelsSSD is not found in ITS/Calib");
222     cdbOK = false;
223   }
224   
225   //General reconstruction
226   if(store->GetLatestVersion("GRP/CTP/Scalers", runNo)<0){
227     HLTError("Scalers is not found in GRP/CTP/");
228     cdbOK = false;
229   }
230   if(!cdbOK){return NULL;}
231
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");
237   }
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");
243   }
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");
249   }
250   else if(fModeSwitch==kClusterFinderDigits) {
251     //tR = new TTree();
252   }
253   else{
254      HLTFatal("No mode set for clusterfindercomponent");
255   }
256
257   //Removed the warning for loading default RecoParam in HLT
258   AliITSRecoParam *par = AliITSRecoParam::GetLowFluxParam();
259   AliITSReconstructor *rec = new AliITSReconstructor();
260   rec->SetRecoParam(par);
261   
262   AliCDBManager* man = AliCDBManager::Instance();
263   if (!man->IsDefaultStorageSet()){
264     HLTError("Default CDB storage has not been set !");
265     return -ENOENT;
266   }
267
268   if(AliGeomManager::GetGeometry()==NULL){
269     AliGeomManager::LoadGeometry();
270   }
271
272  //fgeomInit = new AliITSInitGeometry(kvSPD02,2);
273   fgeomInit = new AliITSInitGeometry(kvPPRasymmFMD,2);
274   //fgeomInit->InitAliITSgeom(fgeom);
275   fgeom = fgeomInit->CreateAliITSgeom();
276  
277   fNModules = fgeom->GetIndexMax();
278
279   fClusters = new TClonesArray*[fNModules]; 
280   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
281     fClusters[iModule] = NULL;
282   } 
283
284   //set dettype
285   fDettype = new AliITSDetTypeRec();
286   fDettype->SetITSgeom(fgeom); 
287   fDettype->SetDefaults();
288   fDettype->SetDefaultClusterFindersV2(kTRUE); 
289
290   if ( fRawReader )
291     return -EINPROGRESS;
292
293   fRawReader = new AliRawReaderMemory();
294   fSPD = new AliHLTITSClusterFinderSPD( fDettype );
295   fSSD = new AliHLTITSClusterFinderSSD( fDettype, fRawReader );
296
297   TString arguments = "";
298   for ( int i = 0; i < argc; i++ ) {
299     if ( !arguments.IsNull() ) arguments += " ";
300     arguments += argv[i];
301   }
302
303   tD = NULL;
304   tR = NULL;
305
306   return Configure( arguments.Data() );
307 }
308
309 Int_t AliHLTITSClusterFinderComponent::DoDeinit() {
310   // see header file for class documentation
311
312   if ( fRawReader )
313     delete fRawReader;
314   fRawReader = NULL;
315
316   if ( fDettype )
317     delete fDettype;
318   fDettype = NULL;
319
320   if ( fgeomInit )
321     delete fgeomInit;
322   fgeomInit = NULL;
323   
324   delete fSPD;
325   fSPD = 0;
326
327   delete fSSD;
328   fSSD = 0;
329
330   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
331     if(fClusters[iModule] != NULL){
332       fClusters[iModule]->Delete();
333       delete fClusters[iModule];
334     }
335     fClusters[iModule] = NULL;
336   } 
337   
338   fUseOfflineFinder = 0;
339
340   return 0;
341 }
342
343 // #include "TStopwatch.h"
344
345 int AliHLTITSClusterFinderComponent::DoEvent
346 (
347  const AliHLTComponentEventData& evtData,
348  const AliHLTComponentBlockData* blocks,
349  AliHLTComponentTriggerData& /*trigData*/,
350  AliHLTUInt8_t* outputPtr,
351  AliHLTUInt32_t& size,
352   vector<AliHLTComponentBlockData>& outputBlocks )
353 {
354  // see header file for class documentation
355   
356   AliHLTUInt32_t maxBufferSize = size;
357   size = 0; // output size
358
359   if (!IsDataEvent()) return 0;
360   
361   if ( evtData.fBlockCnt<=0 )
362     {
363       HLTDebug("no blocks in event" );
364       return 0;
365     }
366   
367   // TStopwatch timer;
368   Int_t ret = 0;
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 ) );
373       if(!tD){
374         HLTFatal("No Digit Tree found");
375         return -1;
376       }
377       tR = new TTree();
378       fDettype->SetTreeAddressD(tD);
379       fDettype->MakeBranch(tR,"R");
380       fDettype->SetTreeAddressR(tR);
381       Option_t *opt="All";
382       fDettype->DigitsToRecPoints(tD,tR,0,opt,kTRUE);
383       
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);
391         }
392       }
393       
394       if(tR){
395         tR->Delete();
396       }
397       UInt_t nClusters=fclusters.size();
398       
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);
402         ret = -ENOSPC;      
403         break;          
404       }
405       if( nClusters>0 ){
406
407         RecPointToSpacePoint(outputPtr,size);
408         
409         AliHLTComponentBlockData bd;
410         FillBlockData( bd );
411         bd.fOffset = size;
412         bd.fSize = bufferSize;
413         bd.fSpecification = 0x00000000;
414         bd.fDataType = GetOutputDataType();
415         outputBlocks.push_back( bd );
416         size += bufferSize;
417         
418         fclusters.clear();      
419       }
420     }
421   }
422   else{
423       
424     // -- Loop over blocks
425     for( const AliHLTComponentBlockData* iter = GetFirstInputBlock(fInputDataType); iter != NULL; iter = GetNextInputBlock() ) {
426       
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());
432       
433       // -- Check for the correct data type
434       if ( iter->fDataType != (fInputDataType) )  
435         continue;
436       
437       // -- Get equipment ID out of specification
438       AliHLTUInt32_t spec = iter->fSpecification;
439       
440       Int_t id = fId;
441       for ( Int_t ii = 0; ii < fNddl ; ii++ ) {   //number of ddl's
442         if ( spec & 0x00000001 ) {
443           id += ii;
444           break;
445         }
446         spec = spec >> 1 ;
447       }
448       
449       // -- Set equipment ID to the raw reader
450       
451       if(!fRawReader->AddBuffer((UChar_t*) iter->fPtr, iter->fSize, id)){
452         HLTWarning("Could not add buffer");
453       }
454       // TStopwatch timer1;
455       
456       if(fModeSwitch==kClusterFinderSPD && !fUseOfflineFinder){ fSPD->RawdataToClusters( fRawReader, fclusters ); }
457       else if(fModeSwitch==kClusterFinderSSD && !fUseOfflineFinder){ fSSD->RawdataToClusters( fclusters ); }
458       else{
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);
467             }
468             fClusters[i]->Delete();
469             delete fClusters[i];
470           }
471           fClusters[i] = NULL;
472         }     
473       }
474       
475       // timer1.Stop();
476       // fStatTime+=timer1.RealTime();
477       // fStatTimeC+=timer1.CpuTime();
478       
479       fRawReader->ClearBuffers();    
480          
481       UInt_t nClusters=fclusters.size();
482       
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);
486         ret = -ENOSPC;      
487         break;          
488       }
489       if( nClusters>0 ){
490
491         RecPointToSpacePoint(outputPtr,size);
492
493         AliHLTComponentBlockData bd;
494         FillBlockData( bd );
495         bd.fOffset = size;
496         bd.fSize = bufferSize;
497         bd.fSpecification = iter->fSpecification;
498         bd.fDataType = GetOutputDataType();
499         outputBlocks.push_back( bd );
500         size += bufferSize;
501         
502         fclusters.clear();      
503       }
504       
505     } // input blocks
506   }
507   /*
508   timer.Stop();
509   
510   fStatTimeAll+=timer.RealTime();
511   fStatTimeAllC+=timer.CpuTime();
512   fStatNEv++;
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;
516   */
517
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(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);      
608     clustIdx++;
609   }
610 }