]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCFileHandler.cxx
- use_newio switch removed from libAliHLTTPC, always on; not removed in
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCFileHandler.cxx
1 // @(#) $Id$
2 // Original: AliHLTFileHandler.cxx,v 1.49 2005/06/23 17:46:55 hristov 
3
4 // Author: Uli Frankenfeld <mailto:franken@fi.uib.no>, Anders Vestbo <mailto:vestbo$fi.uib.no>, C. Loizides <mailto:loizides@ikf.uni-frankfurt.de>
5 //*-- Copyright &copy ALICE HLT Group 
6
7 #include <TClonesArray.h>
8 #include <TSystem.h>
9 #include <TMath.h>
10
11 #include <AliRunLoader.h>
12 #include <TObject.h>
13 #include <TFile.h>
14 #include <TTree.h>
15 #include <AliTPCParamSR.h>
16 #include <AliTPCDigitsArray.h>
17 #include <AliTPCClustersArray.h>
18 #include <AliTPCcluster.h>
19 #include <AliTPCClustersRow.h>
20 #include <AliSimDigits.h>
21
22 #include "AliHLTTPCLogging.h"
23 #include "AliHLTTPCTransform.h"
24 #include "AliHLTTPCDigitData.h"
25 //#include "AliHLTTPCTrackSegmentData.h"
26 #include "AliHLTTPCSpacePointData.h"
27 //#include "AliHLTTPCTrackArray.h"
28 #include "AliHLTTPCFileHandler.h"
29
30 #if __GNUC__ >= 3
31 using namespace std;
32 #endif
33
34 /** \class AliHLTTPCFileHandler
35 <pre>
36 //_____________________________________________________________
37 // AliHLTTPCFileHandler
38 //
39 // The HLT ROOT <-> binary files handling class
40 //
41 // This class provides the interface between AliROOT files,
42 // and HLT binary files. It should be used for converting 
43 // TPC data stored in AliROOT format (outputfile from a simulation),
44 // into the data format currently used by in the HLT framework. 
45 // This enables the possibility to always use the same data format, 
46 // whether you are using a binary file as an input, or a AliROOT file.
47 //
48 // For example on how to create binary files from a AliROOT simulation,
49 // see example macro exa/Binary.C.
50 //
51 // For reading a AliROOT file into HLT format in memory, do the following:
52 //
53 // AliHLTTPCFileHandler file;
54 // file.Init(slice,patch);
55 // file.SetAliInput("galice.root");
56 // AliHLTTPCDigitRowData *dataPt = (AliHLTTPCDigitRowData*)file.AliDigits2Memory(nrows,eventnr);
57 // 
58 // All the data are then stored in memory and accessible via the pointer dataPt.
59 // Accesing the data is then identical to the example 1) showed in AliHLTTPCMemHandler class.
60 //
61 // For converting the data back, and writing it to a new AliROOT file do:
62 //
63 // AliHLTTPCFileHandler file;
64 // file.Init(slice,patch);
65 // file.SetAliInput("galice.root");
66 // file.Init(slice,patch,NumberOfRowsInPatch);
67 // file.AliDigits2RootFile(dataPt,"new_galice.root");
68 // file.CloseAliInput();
69 </pre>
70 */
71
72 ClassImp(AliHLTTPCFileHandler)
73
74 AliHLTTPCFileHandler::AliHLTTPCFileHandler(Bool_t b)
75   :
76   fInAli(NULL),
77   fUseRunLoader(kFALSE),
78   fParam(NULL),
79   fMC(NULL),
80   fDigits(NULL),
81   fDigitsTree(NULL),
82   fIndexCreated(kFALSE),
83   fUseStaticIndex(b)
84 {
85   //Default constructor
86
87   for(Int_t i=0;i<AliHLTTPCTransform::GetNSlice();i++)
88     for(Int_t j=0;j<AliHLTTPCTransform::GetNRows();j++) 
89       fIndex[i][j]=-1;
90
91   if(fUseStaticIndex&&!fgStaticIndexCreated) CleanStaticIndex();
92 }
93
94 AliHLTTPCFileHandler::AliHLTTPCFileHandler(const AliHLTTPCFileHandler& ref)
95   :
96   fInAli(NULL),
97   fUseRunLoader(kFALSE),
98   fParam(NULL),
99   fMC(NULL),
100   fDigits(NULL),
101   fDigitsTree(NULL),
102   fIndexCreated(kFALSE),
103   fUseStaticIndex(ref.fUseStaticIndex)
104 {
105   HLTFatal("copy constructor untested");
106 }
107
108 AliHLTTPCFileHandler& AliHLTTPCFileHandler::operator=(const AliHLTTPCFileHandler&)
109
110   HLTFatal("assignment operator untested");
111   return *this;
112 }
113
114 AliHLTTPCFileHandler::~AliHLTTPCFileHandler()
115 {
116   //Destructor
117   if(fMC) CloseMCOutput();
118   FreeDigitsTree();
119   if(fInAli) CloseAliInput();
120 }
121
122 // of course on start up the index is not created
123 Bool_t AliHLTTPCFileHandler::fgStaticIndexCreated=kFALSE;
124 Int_t  AliHLTTPCFileHandler::fgStaticIndex[36][159]; 
125
126 void AliHLTTPCFileHandler::CleanStaticIndex() 
127
128   // use this static call to clean static index after
129   // running over one event
130   for(Int_t i=0;i<AliHLTTPCTransform::GetNSlice();i++){
131     for(Int_t j=0;j<AliHLTTPCTransform::GetNRows();j++)
132       fgStaticIndex[i][j]=-1;
133   }
134   fgStaticIndexCreated=kFALSE;
135 }
136
137 Int_t AliHLTTPCFileHandler::SaveStaticIndex(Char_t *prefix,Int_t event) 
138
139   // use this static call to store static index after
140   if(!fgStaticIndexCreated) return -1;
141
142   Char_t fname[1024];
143   if(prefix)
144     sprintf(fname,"%s-%d.txt",prefix,event);
145   else
146     sprintf(fname,"TPC.Digits.staticindex-%d.txt",event);
147
148   ofstream file(fname,ios::trunc);
149   if(!file.good()) return -1;
150
151   for(Int_t i=0;i<AliHLTTPCTransform::GetNSlice();i++){
152     for(Int_t j=0;j<AliHLTTPCTransform::GetNRows();j++)
153       file << fgStaticIndex[i][j] << " ";
154     file << endl;
155   }
156   file.close();
157   return 0;
158 }
159
160 Int_t AliHLTTPCFileHandler::LoadStaticIndex(Char_t *prefix,Int_t event) 
161
162   // use this static call to store static index after
163   if(fgStaticIndexCreated){
164       LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::LoadStaticIndex","Inxed")
165         <<"Static index already created, will overwrite"<<ENDLOG;
166       CleanStaticIndex();
167   }
168
169   Char_t fname[1024];
170   if(prefix)
171     sprintf(fname,"%s-%d.txt",prefix,event);
172   else
173     sprintf(fname,"TPC.Digits.staticindex-%d.txt",event);
174
175   ifstream file(fname);
176   if(!file.good()) return -1;
177
178   for(Int_t i=0;i<AliHLTTPCTransform::GetNSlice();i++){
179     for(Int_t j=0;j<AliHLTTPCTransform::GetNRows();j++)
180       file >> fgStaticIndex[i][j];
181   }
182   file.close();
183
184   fgStaticIndexCreated=kTRUE;
185   return 0;
186 }
187
188 void AliHLTTPCFileHandler::FreeDigitsTree()
189
190   //free digits tree
191   if(!fDigitsTree)
192     {
193       LOG(AliHLTTPCLog::kInformational,"AliHLTTPCFileHandler::FreeDigitsTree()","Pointer")
194         <<"Cannot free digitstree, it is not present"<<ENDLOG;
195       return;
196     }
197   delete fDigits;
198   fDigits=0;
199   fDigitsTree=0;
200
201   for(Int_t i=0;i<AliHLTTPCTransform::GetNSlice();i++){
202     for(Int_t j=0;j<AliHLTTPCTransform::GetNRows();j++)
203       fIndex[i][j]=-1;
204   }
205   fIndexCreated=kFALSE;
206 }
207
208 Bool_t AliHLTTPCFileHandler::SetMCOutput(Char_t *name)
209
210   //set mc input
211   fMC = fopen(name,"w");
212   if(!fMC){
213     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::SetMCOutput","File Open")
214       <<"Pointer to File = 0x0 "<<ENDLOG;
215     return kFALSE;
216   }
217   return kTRUE;
218 }
219
220 Bool_t AliHLTTPCFileHandler::SetMCOutput(FILE *file)
221
222   //set mc output
223   fMC = file;
224   if(!fMC){
225     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::SetMCOutput","File Open")
226       <<"Pointer to File = 0x0 "<<ENDLOG;
227     return kFALSE;
228   }
229   return kTRUE;
230 }
231
232 void AliHLTTPCFileHandler::CloseMCOutput()
233
234   //close mc output
235   if(!fMC){
236     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::CloseMCOutPut","File Close")
237       <<"Nothing to Close"<<ENDLOG;
238     return;
239   }
240   fclose(fMC);
241   fMC =0;
242 }
243
244 Bool_t AliHLTTPCFileHandler::SetAliInput()
245
246   //set ali input
247   fInAli->CdGAFile();
248   fParam = (AliTPCParam*)gFile->Get("75x40_100x60_150x60");
249   if(!fParam){
250     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::SetAliInput","File")
251       <<"No TPC parameters found in \""<<gFile->GetName()
252       <<"\", creating standard parameters "
253       <<"which might not be what you want!"<<ENDLOG;
254     fParam = new AliTPCParamSR;
255   }
256   if(!fParam){ 
257     LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::SetAliInput","File Open")
258       <<"No AliTPCParam "<<AliHLTTPCTransform::GetParamName()<<" in File "<<gFile->GetName()<<ENDLOG;
259     return kFALSE;
260   }
261
262   return kTRUE;
263 }
264
265 Bool_t AliHLTTPCFileHandler::SetAliInput(Char_t *name)
266
267   //Open the AliROOT file with name.
268   fInAli= AliRunLoader::Open(name);
269   if(!fInAli){
270     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::SetAliInput","File Open")
271     <<"Pointer to fInAli = 0x0 "<<ENDLOG;
272     return kFALSE;
273   }
274   return SetAliInput();
275 }
276
277 Bool_t AliHLTTPCFileHandler::SetAliInput(AliRunLoader *runLoader)
278
279   //set ali input as runloader
280   fInAli=runLoader;
281   fUseRunLoader = kTRUE;
282   if(!fInAli){
283     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::SetAliInput","File Open")
284     <<"Pointer to AliRunLoader = 0x0 "<<ENDLOG;
285     return kFALSE;
286   }
287   return SetAliInput();
288 }
289
290 void AliHLTTPCFileHandler::CloseAliInput()
291
292   //close ali input
293   if(fUseRunLoader) return;
294   if(!fInAli){
295     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::CloseAliInput","RunLoader")
296       <<"Nothing to Close"<<ENDLOG;
297     return;
298   }
299
300   delete fInAli;
301   fInAli = 0;
302 }
303
304 Bool_t AliHLTTPCFileHandler::IsDigit(Int_t event)
305 {
306   //Check if there is a TPC digit tree in the current file.
307   //Return kTRUE if tree was found, and kFALSE if not found.
308   
309   if(!fInAli){
310     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::IsDigit","File")
311     <<"Pointer to fInAli = 0x0 "<<ENDLOG;
312     return kTRUE;  //maybe you are using binary input which is Digits!
313   }
314   AliLoader* tpcLoader = fInAli->GetLoader("TPCLoader");
315   if(!tpcLoader){
316     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandlerNewIO::IsDigit","File")
317     <<"Pointer to AliLoader for TPC = 0x0 "<<ENDLOG;
318     return kFALSE;
319   }
320   fInAli->GetEvent(event);
321   tpcLoader->LoadDigits();
322   TTree *t=tpcLoader->TreeD();
323   if(t){
324     LOG(AliHLTTPCLog::kInformational,"AliHLTTPCFileHandlerNewIO::IsDigit","File Type")
325     <<"Found Digit Tree -> Use Fast Cluster Finder"<<ENDLOG;
326     return kTRUE;
327   }
328   else{
329     LOG(AliHLTTPCLog::kInformational,"AliHLTTPCFileHandlerNewIO::IsDigit","File Type")
330     <<"No Digit Tree -> Use Cluster Tree"<<ENDLOG;
331     return kFALSE;
332   }
333 }
334
335 ///////////////////////////////////////// Digit IO  
336 Bool_t AliHLTTPCFileHandler::AliDigits2Binary(Int_t event,Bool_t altro)
337 {
338   //save alidigits as binary
339   Bool_t out = kTRUE;
340   UInt_t nrow;
341   AliHLTTPCDigitRowData* data = 0;
342   if(altro)
343     data = AliAltroDigits2Memory(nrow,event);
344   else
345     data = AliDigits2Memory(nrow,event);
346   out = Memory2Binary(nrow,data);
347   Free();
348   return out;
349 }
350
351 Bool_t AliHLTTPCFileHandler::AliDigits2CompBinary(Int_t event,Bool_t altro)
352 {
353   //Convert AliROOT TPC data, into HLT data format.
354   //event specifies the event you want in the aliroot file.
355   
356   Bool_t out = kTRUE;
357   UInt_t ndigits=0;
358   AliHLTTPCDigitRowData *digits=0;
359   if(altro)
360     digits = AliAltroDigits2Memory(ndigits,event);
361   else
362     digits = AliDigits2Memory(ndigits,event);
363   out = Memory2CompBinary(ndigits,digits);
364   Free();
365   return out;
366 }
367
368 Bool_t AliHLTTPCFileHandler::CreateIndex()
369 {
370   //create the access index or copy from static index
371   fIndexCreated=kFALSE;
372
373   if(!fgStaticIndexCreated || !fUseStaticIndex) { //we have to create index 
374     LOG(AliHLTTPCLog::kInformational,"AliHLTTPCFileHandler::CreateIndex","Index")
375       <<"Starting to create index, this can take a while."<<ENDLOG;
376
377     for(Int_t n=0; n<fDigitsTree->GetEntries(); n++) {
378       Int_t sector, row;
379       Int_t lslice,lrow;
380       fDigitsTree->GetEvent(n);
381       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
382       if(!AliHLTTPCTransform::Sector2Slice(lslice,lrow,sector,row)){
383         LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::CreateIndex","Slice/Row")
384           <<AliHLTTPCLog::kDec<<"Index could not be created. Wrong values "
385           <<sector<<" "<<row<<ENDLOG;
386         return kFALSE;
387       }
388       if(fIndex[lslice][lrow]==-1) {
389         fIndex[lslice][lrow]=n;
390       }
391     }
392     if(fUseStaticIndex) { // create static index
393       for(Int_t i=0;i<AliHLTTPCTransform::GetNSlice();i++){
394         for(Int_t j=0;j<AliHLTTPCTransform::GetNRows();j++)
395           fgStaticIndex[i][j]=fIndex[i][j];
396       }
397       fgStaticIndexCreated=kTRUE; //remember that index has been created
398     }
399
400   LOG(AliHLTTPCLog::kInformational,"AliHLTTPCFileHandler::CreateIndex","Index")
401     <<"Index successfully created."<<ENDLOG;
402
403   } else if(fUseStaticIndex) { //simply copy static index
404     for(Int_t i=0;i<AliHLTTPCTransform::GetNSlice();i++){
405       for(Int_t j=0;j<AliHLTTPCTransform::GetNRows();j++)
406         fIndex[i][j]=fgStaticIndex[i][j];
407     }
408
409   LOG(AliHLTTPCLog::kInformational,"AliHLTTPCFileHandler::CreateIndex","Index")
410     <<"Index successfully taken from static copy."<<ENDLOG;
411   }
412   fIndexCreated=kTRUE;
413   return kTRUE;
414 }
415
416 AliHLTTPCDigitRowData * AliHLTTPCFileHandler::AliDigits2Memory(UInt_t & nrow,Int_t event)
417 {
418   //Read data from AliROOT file into memory, and store it in the HLT data format.
419   //Returns a pointer to the data.
420
421   AliHLTTPCDigitRowData *data = 0;
422   nrow=0;
423   
424   if(!fInAli){
425     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliDigits2Memory","File")
426     <<"No Input avalible: Pointer to fInAli == NULL"<<ENDLOG;
427     return 0; 
428   }
429
430   if(!fDigitsTree)
431     if(!GetDigitsTree(event)) return 0;
432
433   UShort_t dig;
434   Int_t time,pad,sector,row;
435   Int_t lslice,lrow;
436   Int_t nrows=0;
437   Int_t ndigitcount=0;
438   Int_t entries = (Int_t)fDigitsTree->GetEntries();
439   if(entries==0) {
440     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliDigits2Memory","ndigits")
441       <<"No TPC digits (entries==0)!"<<ENDLOG;
442     nrow = (UInt_t)(fRowMax-fRowMin+1);
443     Int_t size = nrow*sizeof(AliHLTTPCDigitRowData);
444     data=(AliHLTTPCDigitRowData*) Allocate(size);
445     AliHLTTPCDigitRowData *tempPt = data;
446     for(Int_t r=fRowMin;r<=fRowMax;r++){
447       tempPt->fRow = r;
448       tempPt->fNDigit = 0;
449       tempPt++;
450     }
451     return data;
452   }
453
454   Int_t * ndigits = new Int_t[fRowMax+1];
455   Float_t xyz[3];
456
457   for(Int_t r=fRowMin;r<=fRowMax;r++){
458     Int_t n=fIndex[fSlice][r];
459     if(n!=-1){ //data on that row
460       fDigitsTree->GetEvent(n);
461       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
462       AliHLTTPCTransform::Sector2Slice(lslice,lrow,sector,row);
463
464       if(lrow!=r){
465         LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliDigits2Memory","Row")
466           <<AliHLTTPCLog::kDec<<"Rows in slice " << fSlice << " dont match "<<lrow<<" "<<r<<ENDLOG;
467         continue;
468       }
469
470       ndigits[lrow] = 0;
471       fDigits->First();
472       do {
473         time=fDigits->CurrentRow();
474         pad=fDigits->CurrentColumn();
475         dig = fDigits->GetDigit(time,pad);
476         if(dig <= fParam->GetZeroSup()) continue;
477         if(dig >= AliHLTTPCTransform::GetADCSat())
478           dig = AliHLTTPCTransform::GetADCSat();
479       
480         AliHLTTPCTransform::Raw2Local(xyz,sector,row,pad,time);
481         //      if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
482         //        continue; // why 230???
483
484         ndigits[lrow]++; //for this row only
485         ndigitcount++;   //total number of digits to be published
486
487       } while (fDigits->Next());
488       //cout << lrow << " " << ndigits[lrow] << " - " << ndigitcount << endl;
489     }
490     nrows++;
491   }
492
493   Int_t size = sizeof(AliHLTTPCDigitData)*ndigitcount
494     + nrows*sizeof(AliHLTTPCDigitRowData);
495
496   LOG(AliHLTTPCLog::kDebug,"AliHLTTPCFileHandler::AliDigits2Memory","Digits")
497     <<AliHLTTPCLog::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
498   
499   data=(AliHLTTPCDigitRowData*) Allocate(size);
500   nrow = (UInt_t)nrows;
501   AliHLTTPCDigitRowData *tempPt = data;
502
503   for(Int_t r=fRowMin;r<=fRowMax;r++){
504     Int_t n=fIndex[fSlice][r];
505     tempPt->fRow = r;
506     tempPt->fNDigit = 0;
507
508     if(n!=-1){//data on that row
509       fDigitsTree->GetEvent(n);
510       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
511       AliHLTTPCTransform::Sector2Slice(lslice,lrow,sector,row);
512       if(lrow!=r){
513         LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliDigits2Memory","Row")
514           <<AliHLTTPCLog::kDec<<"Rows on slice " << fSlice << " dont match "<<lrow<<" "<<r<<ENDLOG;
515         continue;
516       }
517
518       tempPt->fNDigit = ndigits[lrow];
519
520       Int_t localcount=0;
521       fDigits->First();
522       do {
523         time=fDigits->CurrentRow();
524         pad=fDigits->CurrentColumn();
525         dig = fDigits->GetDigit(time,pad);
526         if (dig <= fParam->GetZeroSup()) continue;
527         if(dig >= AliHLTTPCTransform::GetADCSat())
528           dig = AliHLTTPCTransform::GetADCSat();
529
530         //Exclude data outside cone:
531         AliHLTTPCTransform::Raw2Local(xyz,sector,row,pad,time);
532         //      if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
533         //        continue; // why 230???
534
535         if(localcount >= ndigits[lrow])
536           LOG(AliHLTTPCLog::kFatal,"AliHLTTPCFileHandler::AliDigits2Binary","Memory")
537             <<AliHLTTPCLog::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
538             <<ndigits[lrow]<<ENDLOG;
539         
540         tempPt->fDigitData[localcount].fCharge=dig;
541         tempPt->fDigitData[localcount].fPad=pad;
542         tempPt->fDigitData[localcount].fTime=time;
543 #ifdef do_mc
544         tempPt->fDigitData[localcount].fTrackID[0] = fDigits->GetTrackID(time,pad,0);
545         tempPt->fDigitData[localcount].fTrackID[1] = fDigits->GetTrackID(time,pad,1);
546         tempPt->fDigitData[localcount].fTrackID[2] = fDigits->GetTrackID(time,pad,2);
547 #endif
548         localcount++;
549       } while (fDigits->Next());
550     }
551     Byte_t *tmp = (Byte_t*)tempPt;
552     Int_t size = sizeof(AliHLTTPCDigitRowData)
553                                       + ndigits[lrow]*sizeof(AliHLTTPCDigitData);
554     tmp += size;
555     tempPt = (AliHLTTPCDigitRowData*)tmp;
556   }
557   delete [] ndigits;
558   return data;
559 }
560
561 AliHLTTPCDigitRowData * AliHLTTPCFileHandler::AliAltroDigits2Memory(UInt_t & nrow,Int_t event,Bool_t eventmerge)
562 {
563   //Read data from AliROOT file into memory, and store it in the HLT data format.
564   //Returns a pointer to the data.
565   //This functions filter out single timebins, which is noise. The timebins which
566   //are removed are timebins which have the 4 zero neighbours; 
567   //(pad-1,time),(pad+1,time),(pad,time-1),(pad,time+1).
568   
569   AliHLTTPCDigitRowData *data = 0;
570   nrow=0;
571   
572   if(!fInAli){
573     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliAltroDigits2Memory","File")
574     <<"No Input avalible: Pointer to TFile == NULL"<<ENDLOG;
575     return 0; 
576   }
577   if(eventmerge == kTRUE && event >= 1024)
578     {
579       LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliAltroDigits2Memory","TrackIDs")
580         <<"Too many events if you want to merge!"<<ENDLOG;
581       return 0;
582     }
583   delete fDigits;
584   fDigits=0;
585   /* Dont understand why we have to do 
586      reload the tree, but otherwise the code crashes */
587   fDigitsTree=0;
588   if(!GetDigitsTree(event)) return 0;
589
590   UShort_t dig;
591   Int_t time,pad,sector,row;
592   Int_t nrows=0;
593   Int_t ndigitcount=0;
594   Int_t entries = (Int_t)fDigitsTree->GetEntries();
595   if(entries==0) {
596     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliAltroDigits2Memory","ndigits")
597       <<"No TPC digits (entries==0)!"<<ENDLOG;
598     nrow = (UInt_t)(fRowMax-fRowMin+1);
599     Int_t size = nrow*sizeof(AliHLTTPCDigitRowData);
600     data=(AliHLTTPCDigitRowData*) Allocate(size);
601     AliHLTTPCDigitRowData *tempPt = data;
602     for(Int_t r=fRowMin;r<=fRowMax;r++){
603       tempPt->fRow = r;
604       tempPt->fNDigit = 0;
605       tempPt++;
606     }
607     return data;
608   }
609   Int_t * ndigits = new Int_t[fRowMax+1];
610   Int_t lslice,lrow;
611   Int_t zerosupval=AliHLTTPCTransform::GetZeroSup();
612   Float_t xyz[3];
613
614   for(Int_t r=fRowMin;r<=fRowMax;r++){
615     Int_t n=fIndex[fSlice][r];
616
617     ndigits[r] = 0;
618
619     if(n!=-1){//data on that row
620       fDigitsTree->GetEvent(n);
621       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
622       AliHLTTPCTransform::Sector2Slice(lslice,lrow,sector,row);
623       //cout << lslice << " " << fSlice << " " << lrow << " " << r << " " << sector << " " << row << endl;
624       if((lslice!=fSlice)||(lrow!=r)){
625         LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliAltroDigits2Memory","Row")
626           <<AliHLTTPCLog::kDec<<"Rows on slice " << fSlice << " dont match "<<lrow<<" "<<r<<ENDLOG;
627         continue;
628       }
629
630       fDigits->ExpandBuffer();
631       fDigits->ExpandTrackBuffer();
632       for(Int_t i=0; i<fDigits->GetNCols(); i++){
633         for(Int_t j=0; j<fDigits->GetNRows(); j++){
634           pad=i;
635           time=j;
636           dig = fDigits->GetDigitFast(time,pad);
637           if(dig <= zerosupval) continue;
638           if(dig >= AliHLTTPCTransform::GetADCSat())
639             dig = AliHLTTPCTransform::GetADCSat();
640
641           //Check for single timebins, and remove them because they are noise for sure.
642           if(i>0 && i<fDigits->GetNCols()-1 && j>0 && j<fDigits->GetNRows()-1)
643             if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
644                fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
645                fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
646                fDigits->GetDigitFast(time,pad+1)<=zerosupval)
647               continue;
648               
649           //Boundaries:
650           if(i==0) //pad==0
651             {
652               if(j < fDigits->GetNRows()-1 && j > 0) 
653                 {
654                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
655                      fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
656                      fDigits->GetDigitFast(time,pad+1)<=zerosupval)
657                     continue;
658                 }
659               else if(j > 0)
660                 {
661                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
662                      fDigits->GetDigitFast(time,pad+1)<=zerosupval)
663                     continue;
664                 }
665             }
666           if(j==0)
667             {
668               if(i < fDigits->GetNCols()-1 && i > 0)
669                 {
670                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
671                      fDigits->GetDigitFast(time,pad+1)<=zerosupval &&
672                      fDigits->GetDigitFast(time+1,pad)<=zerosupval)
673                     continue;
674                 }
675               else if(i > 0)
676                 {
677                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
678                      fDigits->GetDigitFast(time+1,pad)<=zerosupval)
679                     continue;
680                 }
681             }
682
683           if(i==fDigits->GetNCols()-1)
684             {
685               if(j>0 && j<fDigits->GetNRows()-1)
686                 {
687                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
688                      fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
689                      fDigits->GetDigitFast(time,pad-1)<=zerosupval)
690                     continue;
691                 }
692               else if(j==0 && j<fDigits->GetNRows()-1)
693                 {
694                   if(fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
695                      fDigits->GetDigitFast(time,pad-1)<=zerosupval)
696                     continue;
697                 }
698               else 
699                 {
700                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
701                      fDigits->GetDigitFast(time,pad-1)<=zerosupval)
702                     continue;
703                 }
704             }
705         
706           if(j==fDigits->GetNRows()-1)
707             {
708               if(i>0 && i<fDigits->GetNCols()-1)
709                 {
710                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
711                      fDigits->GetDigitFast(time,pad+1)<=zerosupval &&
712                      fDigits->GetDigitFast(time-1,pad)<=zerosupval)
713                     continue;
714                 }
715               else if(i==0 && fDigits->GetNCols()-1)
716                 {
717                   if(fDigits->GetDigitFast(time,pad+1)<=zerosupval &&
718                      fDigits->GetDigitFast(time-1,pad)<=zerosupval)
719                     continue;
720                 }
721               else 
722                 {
723                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
724                      fDigits->GetDigitFast(time-1,pad)<=zerosupval)
725                     continue;
726                 }
727             }
728
729           AliHLTTPCTransform::Raw2Local(xyz,sector,row,pad,time);
730           //      if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
731           //      continue; 
732               
733           ndigits[lrow]++; //for this row only
734           ndigitcount++;   //total number of digits to be published
735         }
736       }
737     }
738     nrows++;
739   }
740   
741   Int_t size = sizeof(AliHLTTPCDigitData)*ndigitcount
742     + nrows*sizeof(AliHLTTPCDigitRowData);
743
744   LOG(AliHLTTPCLog::kDebug,"AliHLTTPCFileHandler::AliAltroDigits2Memory","Digits")
745     <<AliHLTTPCLog::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
746   
747   data=(AliHLTTPCDigitRowData*) Allocate(size);
748   nrow = (UInt_t)nrows;
749   AliHLTTPCDigitRowData *tempPt = data;
750  
751   for(Int_t r=fRowMin;r<=fRowMax;r++){
752     Int_t n=fIndex[fSlice][r];
753     tempPt->fRow = r;
754     tempPt->fNDigit = 0;
755     if(n!=-1){ //no data on that row
756       fDigitsTree->GetEvent(n);
757       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
758       AliHLTTPCTransform::Sector2Slice(lslice,lrow,sector,row);
759
760       if(lrow!=r){
761         LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliAltroDigits2Memory","Row")
762           <<AliHLTTPCLog::kDec<<"Rows on slice " << fSlice << " dont match "<<lrow<<" "<<r<<ENDLOG;
763         continue;
764       }
765
766       tempPt->fNDigit = ndigits[lrow];
767
768       Int_t localcount=0;
769       fDigits->ExpandBuffer();
770       fDigits->ExpandTrackBuffer();
771       for(Int_t i=0; i<fDigits->GetNCols(); i++){
772         for(Int_t j=0; j<fDigits->GetNRows(); j++){
773           pad=i;
774           time=j;
775           dig = fDigits->GetDigitFast(time,pad);
776           if(dig <= zerosupval) continue;
777           if(dig >= AliHLTTPCTransform::GetADCSat())
778             dig = AliHLTTPCTransform::GetADCSat();
779               
780           //Check for single timebins, and remove them because they are noise for sure.
781           if(i>0 && i<fDigits->GetNCols()-1 && j>0 && j<fDigits->GetNRows()-1)
782             if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
783                fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
784                fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
785                fDigits->GetDigitFast(time,pad+1)<=zerosupval)
786               continue;
787           
788           //Boundaries:
789           if(i==0) //pad ==0
790             {
791               if(j < fDigits->GetNRows()-1 && j > 0) 
792                 {
793                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
794                      fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
795                      fDigits->GetDigitFast(time,pad+1)<=zerosupval)
796                     continue;
797                 }
798               else if(j > 0)
799                 {
800                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
801                      fDigits->GetDigitFast(time,pad+1)<=zerosupval)
802                     continue;
803                 }
804             }
805           if(j==0)
806             {
807               if(i < fDigits->GetNCols()-1 && i > 0)
808                 {
809                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
810                      fDigits->GetDigitFast(time,pad+1)<=zerosupval &&
811                      fDigits->GetDigitFast(time+1,pad)<=zerosupval)
812                     continue;
813                 }
814               else if(i > 0)
815                 {
816                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
817                      fDigits->GetDigitFast(time+1,pad)<=zerosupval)
818                     continue;
819                 }
820             }
821         
822           if(i == fDigits->GetNCols()-1)
823             {
824               if(j>0 && j<fDigits->GetNRows()-1)
825                 {
826                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
827                      fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
828                      fDigits->GetDigitFast(time,pad-1)<=zerosupval)
829                     continue;
830                 }
831               else if(j==0 && j<fDigits->GetNRows()-1)
832                 {
833                   if(fDigits->GetDigitFast(time+1,pad)<=zerosupval &&
834                      fDigits->GetDigitFast(time,pad-1)<=zerosupval)
835                     continue;
836                 }
837               else 
838                 {
839                   if(fDigits->GetDigitFast(time-1,pad)<=zerosupval &&
840                      fDigits->GetDigitFast(time,pad-1)<=zerosupval)
841                     continue;
842                 }
843             }
844           if(j==fDigits->GetNRows()-1)
845             {
846               if(i>0 && i<fDigits->GetNCols()-1)
847                 {
848                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
849                      fDigits->GetDigitFast(time,pad+1)<=zerosupval &&
850                      fDigits->GetDigitFast(time-1,pad)<=zerosupval)
851                     continue;
852                 }
853               else if(i==0 && fDigits->GetNCols()-1)
854                 {
855                   if(fDigits->GetDigitFast(time,pad+1)<=zerosupval &&
856                      fDigits->GetDigitFast(time-1,pad)<=zerosupval)
857                     continue;
858                 }
859               else 
860                 {
861                   if(fDigits->GetDigitFast(time,pad-1)<=zerosupval &&
862                      fDigits->GetDigitFast(time-1,pad)<=zerosupval)
863                     continue;
864                 }
865             }
866         
867           AliHLTTPCTransform::Raw2Local(xyz,sector,row,pad,time);
868           //      if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
869           //        continue;
870           
871           if(localcount >= ndigits[lrow])
872             LOG(AliHLTTPCLog::kFatal,"AliHLTTPCFileHandler::AliAltroDigits2Binary","Memory")
873               <<AliHLTTPCLog::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
874               <<ndigits[lrow]<<ENDLOG;
875         
876           tempPt->fDigitData[localcount].fCharge=dig;
877           tempPt->fDigitData[localcount].fPad=pad;
878           tempPt->fDigitData[localcount].fTime=time;
879 #ifdef do_mc
880           tempPt->fDigitData[localcount].fTrackID[0] = (fDigits->GetTrackIDFast(time,pad,0)-2);
881           tempPt->fDigitData[localcount].fTrackID[1] = (fDigits->GetTrackIDFast(time,pad,1)-2);
882           tempPt->fDigitData[localcount].fTrackID[2] = (fDigits->GetTrackIDFast(time,pad,2)-2);
883           if(eventmerge == kTRUE) //careful track mc info will be touched
884             {//Event are going to be merged, so event number is stored in the upper 10 bits.
885               tempPt->fDigitData[localcount].fTrackID[0] += 128; //leave some room
886               tempPt->fDigitData[localcount].fTrackID[1] += 128; //for neg. numbers
887               tempPt->fDigitData[localcount].fTrackID[2] += 128;
888               tempPt->fDigitData[localcount].fTrackID[0] += ((event&0x3ff)<<22);
889               tempPt->fDigitData[localcount].fTrackID[1] += ((event&0x3ff)<<22);
890               tempPt->fDigitData[localcount].fTrackID[2] += ((event&0x3ff)<<22);
891             }
892 #endif
893           localcount++;
894         }
895       }
896     }
897     Byte_t *tmp = (Byte_t*)tempPt;
898     Int_t size = sizeof(AliHLTTPCDigitRowData)
899       + ndigits[r]*sizeof(AliHLTTPCDigitData);
900     tmp += size;
901     tempPt = (AliHLTTPCDigitRowData*)tmp;
902   }
903   delete [] ndigits;
904   return data;
905 }
906  
907 Bool_t AliHLTTPCFileHandler::GetDigitsTree(Int_t event)
908 {
909   //Connects to the TPC digit tree in the AliROOT file.
910   AliLoader* tpcLoader = fInAli->GetLoader("TPCLoader");
911   if(!tpcLoader){
912     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::GetDigitsTree","File")
913     <<"Pointer to AliLoader for TPC = 0x0 "<<ENDLOG;
914     return kFALSE;
915   }
916   fInAli->GetEvent(event);
917   tpcLoader->LoadDigits();
918   fDigitsTree = tpcLoader->TreeD();
919   if(!fDigitsTree) 
920     {
921       LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::GetDigitsTree","Digits Tree")
922         <<AliHLTTPCLog::kHex<<"Error getting digitstree "<<(void*)fDigitsTree<<ENDLOG;
923       return kFALSE;
924     }
925   fDigitsTree->GetBranch("Segment")->SetAddress(&fDigits);
926
927   if(!fIndexCreated) return CreateIndex();
928   else return kTRUE;
929 }
930
931 void AliHLTTPCFileHandler::AliDigits2RootFile(AliHLTTPCDigitRowData *rowPt,Char_t *newDigitsfile)
932 {
933   //Write the data stored in rowPt, into a new AliROOT file.
934   //The data is stored in the AliROOT format 
935   //This is specially a nice thing if you have modified data, and wants to run it  
936   //through the offline reconstruction chain.
937   //The arguments is a pointer to the data, and the name of the new AliROOT file.
938   //Remember to pass the original AliROOT file (the one that contains the original
939   //simulated data) to this object, in order to retrieve the MC id's of the digits.
940
941   if(!fInAli)
942     {
943       LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliDigits2RootFile","File")
944         <<"No rootfile "<<ENDLOG;
945       return;
946     }
947   if(!fParam)
948     {
949       LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliDigits2RootFile","File")
950         <<"No parameter object. Run on rootfile "<<ENDLOG;
951       return;
952     }
953
954   //Get the original digitstree:
955   AliLoader* tpcLoader = fInAli->GetLoader("TPCLoader");
956   if(!tpcLoader){
957     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliDigits2RootFile","File")
958     <<"Pointer to AliLoader for TPC = 0x0 "<<ENDLOG;
959     return;
960   }
961   tpcLoader->LoadDigits();
962   TTree *t=tpcLoader->TreeD();
963
964   AliTPCDigitsArray *oldArray = new AliTPCDigitsArray();
965   oldArray->Setup(fParam);
966   oldArray->SetClass("AliSimDigits");
967
968   Bool_t ok = oldArray->ConnectTree(t);
969   if(!ok)
970     {
971       LOG(AliHLTTPCLog::kError,"AliHLTTPCFileHandler::AliDigits2RootFile","File")
972         << "No digits tree object" << ENDLOG;
973       return;
974     }
975
976   tpcLoader->SetDigitsFileName(newDigitsfile);
977   tpcLoader->MakeDigitsContainer();
978     
979   //setup a new one, or connect it to the existing one:
980   AliTPCDigitsArray *arr = new AliTPCDigitsArray(); 
981   arr->SetClass("AliSimDigits");
982   arr->Setup(fParam);
983   arr->MakeTree(tpcLoader->TreeD());
984
985   Int_t digcounter=0,trackID[3];
986
987   for(Int_t i=fRowMin; i<=fRowMax; i++)
988     {
989       
990       if((Int_t)rowPt->fRow != i) 
991         LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliDigits2RootFile","Data")
992           <<"Mismatching row numbering "<<(Int_t)rowPt->fRow<<" "<<i<<ENDLOG;
993             
994       Int_t sector,row;
995       AliHLTTPCTransform::Slice2Sector(fSlice,i,sector,row);
996       
997       AliSimDigits *oldDig = (AliSimDigits*)oldArray->LoadRow(sector,row);
998       AliSimDigits * dig = (AliSimDigits*)arr->CreateRow(sector,row);
999       oldDig->ExpandBuffer();
1000       oldDig->ExpandTrackBuffer();
1001       dig->ExpandBuffer();
1002       dig->ExpandTrackBuffer();
1003       
1004       if(!oldDig)
1005         LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliDigits2RootFile","Data")
1006           <<"No padrow " << sector << " " << row <<ENDLOG;
1007
1008       AliHLTTPCDigitData *digPt = rowPt->fDigitData;
1009       digcounter=0;
1010       for(UInt_t j=0; j<rowPt->fNDigit; j++)
1011         {
1012           Short_t charge = (Short_t)digPt[j].fCharge;
1013           Int_t pad = (Int_t)digPt[j].fPad;
1014           Int_t time = (Int_t)digPt[j].fTime;
1015           
1016           if(charge == 0) //Only write the digits that has not been removed
1017             {
1018               LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliDigits2RootFile","Data")
1019                 <<"Zero charge" <<ENDLOG;
1020               continue;
1021             }
1022
1023           digcounter++;
1024           
1025           //Tricks to get and set the correct track id's. 
1026           for(Int_t t=0; t<3; t++)
1027             {
1028               Int_t label = oldDig->GetTrackIDFast(time,pad,t);
1029               if(label > 1)
1030                 trackID[t] = label - 2;
1031               else if(label==0)
1032                 trackID[t] = -2;
1033               else
1034                 trackID[t] = -1;
1035             }
1036           
1037           dig->SetDigitFast(charge,time,pad);
1038           
1039           for(Int_t t=0; t<3; t++)
1040             ((AliSimDigits*)dig)->SetTrackIDFast(trackID[t],time,pad,t);
1041           
1042         }
1043       //cout<<"Wrote "<<digcounter<<" on row "<<i<<endl;
1044       UpdateRowPointer(rowPt);
1045       arr->StoreRow(sector,row);
1046       arr->ClearRow(sector,row);  
1047       oldArray->ClearRow(sector,row);
1048     }
1049
1050   char treeName[100];
1051   sprintf(treeName,"TreeD_%s_0",fParam->GetTitle());
1052   
1053   arr->GetTree()->SetName(treeName);
1054   arr->GetTree()->AutoSave();
1055   tpcLoader->WriteDigits("OVERWRITE");
1056   delete arr;
1057   delete oldArray;
1058 }
1059
1060 ///////////////////////////////////////// Point IO  
1061 Bool_t AliHLTTPCFileHandler::AliPoints2Binary(Int_t eventn)
1062 {
1063   //points to binary
1064   Bool_t out = kTRUE;
1065   UInt_t npoint;
1066   AliHLTTPCSpacePointData *data = AliPoints2Memory(npoint,eventn);
1067   out = Memory2Binary(npoint,data);
1068   Free();
1069   return out;
1070 }
1071
1072 AliHLTTPCSpacePointData * AliHLTTPCFileHandler::AliPoints2Memory(UInt_t & npoint,Int_t eventn)
1073 {
1074   //points to memory
1075   AliHLTTPCSpacePointData *data = 0;
1076   npoint=0;
1077   if(!fInAli){
1078     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliPoints2Memory","File")
1079     <<"No Input avalible: no object fInAli"<<ENDLOG;
1080     return 0;
1081   }
1082
1083   TDirectory *savedir = gDirectory;
1084   AliLoader* tpcLoader = fInAli->GetLoader("TPCLoader");
1085   if(!tpcLoader){
1086     LOG(AliHLTTPCLog::kWarning,"AliHLTTPCFileHandler::AliPoints2Memory","File")
1087     <<"Pointer to AliLoader for TPC = 0x0 "<<ENDLOG;
1088     return 0;
1089   }
1090   fInAli->GetEvent(eventn);
1091   tpcLoader->LoadRecPoints();
1092
1093   AliTPCClustersArray carray;
1094   carray.Setup(fParam);
1095   carray.SetClusterType("AliTPCcluster");
1096   Bool_t clusterok = carray.ConnectTree(tpcLoader->TreeR());
1097
1098   if(!clusterok) return 0;
1099
1100   AliTPCClustersRow ** clusterrow = 
1101                new AliTPCClustersRow*[ (int)carray.GetTree()->GetEntries()];
1102   Int_t *rows = new int[ (int)carray.GetTree()->GetEntries()];
1103   Int_t *sects = new int[  (int)carray.GetTree()->GetEntries()];
1104   Int_t sum=0;
1105
1106   Int_t lslice,lrow;
1107   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
1108     AliSegmentID *s = carray.LoadEntry(i);
1109     Int_t sector,row;
1110     fParam->AdjustSectorRow(s->GetID(),sector,row);
1111     rows[i] = row;
1112     sects[i] = sector;
1113     clusterrow[i] = 0;
1114     AliHLTTPCTransform::Sector2Slice(lslice,lrow,sector,row);
1115     if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
1116     clusterrow[i] = carray.GetRow(sector,row);
1117     if(clusterrow[i])
1118       sum+=clusterrow[i]->GetArray()->GetEntriesFast();
1119   }
1120   UInt_t size = sum*sizeof(AliHLTTPCSpacePointData);
1121
1122   LOG(AliHLTTPCLog::kDebug,"AliHLTTPCFileHandler::AliPoints2Memory","File")
1123   <<AliHLTTPCLog::kDec<<"Found "<<sum<<" SpacePoints"<<ENDLOG;
1124
1125   data = (AliHLTTPCSpacePointData *) Allocate(size);
1126   npoint = sum;
1127   UInt_t n=0; 
1128   Int_t pat=fPatch;
1129   if(fPatch==-1)
1130     pat=0;
1131   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
1132     if(!clusterrow[i]) continue;
1133     Int_t row = rows[i];
1134     Int_t sector = sects[i];
1135     AliHLTTPCTransform::Sector2Slice(lslice,lrow,sector,row);
1136     Int_t entriesInRow = clusterrow[i]->GetArray()->GetEntriesFast();
1137     for(Int_t j = 0;j<entriesInRow;j++){
1138       AliTPCcluster *c = (AliTPCcluster*)(*clusterrow[i])[j];
1139       data[n].fZ = c->GetZ();
1140       data[n].fY = c->GetY();
1141       data[n].fX = fParam->GetPadRowRadii(sector,row);
1142       data[n].fCharge = (UInt_t)c->GetQ();
1143       data[n].fID = n+((fSlice&0x7f)<<25)+((pat&0x7)<<22);//uli
1144       data[n].fPadRow = lrow;
1145       data[n].fSigmaY2 = c->GetSigmaY2();
1146       data[n].fSigmaZ2 = c->GetSigmaZ2();
1147 #ifdef do_mc
1148       data[n].fTrackID[0] = c->GetLabel(0);
1149       data[n].fTrackID[1] = c->GetLabel(1);
1150       data[n].fTrackID[2] = c->GetLabel(2);
1151 #endif
1152       if(fMC) fprintf(fMC,"%d %d\n",data[n].fID,c->GetLabel(0));
1153       n++;
1154     }
1155   }
1156   for(Int_t i=0;i<carray.GetTree()->GetEntries();i++){
1157     Int_t row = rows[i];
1158     Int_t sector = sects[i];
1159     if(carray.GetRow(sector,row))
1160       carray.ClearRow(sector,row);
1161   }
1162
1163   delete [] clusterrow;
1164   delete [] rows;
1165   delete [] sects;
1166   savedir->cd();   
1167
1168   return data;
1169 }
1170