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