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