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