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