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