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