]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/src/AliL3FileHandler.cxx
Possible to read several events from same rootfile
[u/mrichter/AliRoot.git] / HLT / src / AliL3FileHandler.cxx
1 //$Id$
2
3 // Author: Uli Frankenfeld <mailto:franken@fi.uib.no>
4 //*-- Copyright &copy Uli 
5
6 #include <math.h>
7 #include <iostream.h>
8 #include <TObject.h>
9 #include <TFile.h>
10 #include <TTree.h>
11 #include <stdio.h>
12
13 #include "AliL3Transform.h"
14 #include "AliL3Logging.h"
15 #include "AliL3MemHandler.h"
16 #include "AliL3FileHandler.h"
17
18 #include "AliTPCDigitsArray.h"
19 #include "AliTPCClustersArray.h"
20 #include "AliTPCcluster.h"
21 #include "AliTPCClustersRow.h"
22 #include "AliTPCParam.h"
23 #include "AliSimDigits.h"
24
25 #include "AliL3DigitData.h"
26 #include "AliL3TrackSegmentData.h"
27 #include "AliL3SpacePointData.h"
28 #include "AliL3TrackArray.h"
29
30 //_____________________________________________________________
31 // AliL3FileHandler
32 //
33 // The L3 Binary File handler 
34 //
35
36 ClassImp(AliL3FileHandler)
37
38 AliL3FileHandler::AliL3FileHandler(){
39   //Default constructor
40   fInAli = 0;
41   fParam = 0;
42   fTransformer = 0;
43   fMC =0;
44 }
45
46
47 AliL3FileHandler::~AliL3FileHandler(){
48   //Destructor
49   if(fTransformer) delete fTransformer;
50   if(fMC) CloseMCOutput();
51   if(fInAli) CloseAliInput();
52 }
53
54 Bool_t AliL3FileHandler::SetMCOutput(char *name){
55   fMC = fopen(name,"w");
56   if(!fMC){
57     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
58     <<"Pointer to File = 0x0 "<<ENDLOG;
59     return kFALSE;
60   }
61   return kTRUE;
62 }
63
64 Bool_t AliL3FileHandler::SetMCOutput(FILE *file){
65   fMC = file;
66   if(!fMC){
67     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
68     <<"Pointer to File = 0x0 "<<ENDLOG;
69     return kFALSE;
70   }
71   return kTRUE;
72 }
73
74 void AliL3FileHandler::CloseMCOutput(){
75   if(!fMC){
76     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseMCOutPut","File Close")
77     <<"Nothing to Close"<<ENDLOG;
78     return;
79   }
80   fclose(fMC);
81   fMC =0;
82 }
83
84 Bool_t AliL3FileHandler::SetAliInput(){
85   if(!fInAli->IsOpen()){
86     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
87     <<"Ali File "<<fInAli->GetName()<<" does not exist"<<ENDLOG;
88     return kFALSE;
89   }
90   fParam = (AliTPCParam*)fInAli->Get("75x40_100x60");
91   if(!fParam){ 
92     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
93     <<"No AliTPCParam 75x40_100x60 in File "<<fInAli->GetName()<<ENDLOG;
94      return kFALSE;
95   }
96   fTransformer = new AliL3Transform();
97   return kTRUE;
98 }
99
100 Bool_t AliL3FileHandler::SetAliInput(char *name){
101   fInAli= new TFile(name,"READ");
102   if(!fInAli){
103     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
104     <<"Pointer to TFile = 0x0 "<<ENDLOG;
105     return kFALSE;
106   }
107   return SetAliInput();
108 }
109
110 Bool_t AliL3FileHandler::SetAliInput(TFile *file){
111   fInAli=file;
112   if(!fInAli){
113     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
114     <<"Pointer to TFile = 0x0 "<<ENDLOG;
115     return kFALSE;
116   }
117   return SetAliInput();
118 }
119
120 void AliL3FileHandler::CloseAliInput(){
121   if(!fInAli){
122     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseAliInput","File Close")
123     <<"Nothing to Close"<<ENDLOG;
124      return;
125   }
126   if(fInAli->IsOpen()) fInAli->Close();
127   delete fInAli;
128   fInAli = 0;
129     
130 }
131
132 Bool_t AliL3FileHandler::IsDigit(){
133   if(!fInAli){
134     LOG(AliL3Log::kWarning,"AliL3FileHandler::IsDigit","File")
135     <<"Pointer to TFile = 0x0 "<<ENDLOG;
136     return kTRUE;  //may you are use binary input which is Digits!!
137   }
138   TTree *t=(TTree*)fInAli->Get("TreeD_75x40_100x60_0");
139   if(t){
140     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
141     <<"Found Digit Tree -> Use Fast Cluster Finder"<<ENDLOG;
142     return kTRUE;
143   }
144   else{
145     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
146     <<"No Digit Tree -> Use Cluster Tree"<<ENDLOG;
147     return kFALSE;
148   }
149 }
150
151 ///////////////////////////////////////// Digit IO  
152 Bool_t AliL3FileHandler::AliDigits2Binary(Int_t event){
153   Bool_t out = kTRUE;
154   UInt_t nrow;
155   AliL3DigitRowData* data = AliDigits2Memory(nrow,event);
156   out = Memory2Binary(nrow,data);
157   Free();
158   return out;
159 }
160
161
162 Bool_t AliL3FileHandler::AliDigits2CompBinary(Int_t event){
163   Bool_t out = kTRUE;
164   UInt_t ndigits=0;
165   AliL3DigitRowData* digits = AliDigits2Memory(ndigits,event);
166   out = Memory2CompBinary(ndigits,digits);
167   Free();
168   return out;
169 }
170
171
172 AliL3DigitRowData * AliL3FileHandler::AliDigits2Memory(UInt_t & nrow,Int_t event){
173   AliL3DigitRowData *data = 0;
174   nrow=0;
175   if(!fInAli){
176     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
177     <<"No Input avalible: no object TFile"<<ENDLOG;
178     return 0; 
179   }
180   if(!fInAli->IsOpen()){
181     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
182     <<"No Input avalible: TFile not opend"<<ENDLOG;
183     return 0;
184   }
185
186   
187   TDirectory *savedir = gDirectory;
188   fInAli->cd();
189   //TTree *t=(TTree*)fInAli->Get("TreeD_75x40_100x60_0");
190   Char_t dname[100];
191   sprintf(dname,"TreeD_75x40_100x60_%d",event);
192   TTree *t=(TTree*)fInAli->Get(dname);
193   if(!t){
194     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Binary","AliRoot")
195     <<"No Digit Tree inside!"<<ENDLOG;
196     return 0;
197   }
198   AliSimDigits digarr, *dummy=&digarr;
199   t->GetBranch("Segment")->SetAddress(&dummy);
200   UShort_t dig;
201   Int_t time,pad,sector,row;
202   Int_t nrows=0;
203   Int_t ndigitcount=0;
204   Int_t entries = (Int_t)t->GetEntries();
205   Int_t ndigits[entries];
206   Int_t lslice,lrow;
207   for(Int_t n=0; n<t->GetEntries(); n++)
208     {
209       t->GetEvent(n);
210       fParam->AdjustSectorRow(digarr.GetID(),sector,row);
211       fTransformer->Sector2Slice(lslice,lrow,sector,row);
212       if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
213       //      if(fSlice != lslice) continue;
214
215       Float_t xyz[3];
216       ndigits[lrow] = 0;
217       digarr.First();
218       do {
219         time=digarr.CurrentRow();
220         pad=digarr.CurrentColumn();
221         dig = digarr.GetDigit(time,pad);
222         if(dig<=fParam->GetZeroSup()) continue;
223         if(time < fParam->GetMaxTBin()-1 && time > 0)
224           if(digarr.GetDigit(time+1,pad) <= fParam->GetZeroSup()
225              && digarr.GetDigit(time-1,pad) <= fParam->GetZeroSup())
226             continue;
227
228         fTransformer->Raw2Local(xyz,sector,row,pad,time);
229         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
230           continue;
231
232         ndigits[lrow]++; //for this row only
233         ndigitcount++;  //total number of digits to be published
234
235       } while (digarr.Next());
236
237       nrows++;
238     }
239   Int_t size = sizeof(AliL3DigitData)*ndigitcount
240     + nrows*sizeof(AliL3DigitRowData);
241
242   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2Memory","Digits")
243   <<AliL3Log::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
244
245   data=(AliL3DigitRowData*) Allocate(size);
246   nrow = (UInt_t)nrows;
247   AliL3DigitRowData *tempPt = data;
248   for(Int_t n=0; n<t->GetEntries(); n++)
249     {
250       t->GetEvent(n);
251
252       Float_t xyz[3];
253       fParam->AdjustSectorRow(digarr.GetID(),sector,row);
254       fTransformer->Sector2Slice(lslice,lrow,sector,row);
255 //      if(fSlice != lslice) continue;
256       if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
257       tempPt->fRow = lrow;
258       tempPt->fNDigit = ndigits[lrow];
259
260       Int_t localcount=0;
261       digarr.First();
262       do {
263         //dig=digarr.CurrentDigit();
264         //if (dig<=fParam->GetZeroSup()) continue;
265         time=digarr.CurrentRow();
266         pad=digarr.CurrentColumn();
267         dig = digarr.GetDigit(time,pad);
268         if (dig <= fParam->GetZeroSup()) continue;
269         if(time < fParam->GetMaxTBin()-1 && time > 0)
270           if(digarr.GetDigit(time-1,pad) <= fParam->GetZeroSup() &&
271              digarr.GetDigit(time+1,pad) <= fParam->GetZeroSup()) continue;
272
273         //Exclude data outside cone:
274         fTransformer->Raw2Local(xyz,sector,row,pad,time);
275         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
276           continue;
277
278         if(localcount >= ndigits[lrow])
279           LOG(AliL3Log::kFatal,"AliL3FileHandler::AliDigits2Binary","Memory")
280             <<AliL3Log::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
281             <<ndigits[lrow]<<ENDLOG;
282         
283         tempPt->fDigitData[localcount].fCharge=dig;
284         tempPt->fDigitData[localcount].fPad=pad;
285         tempPt->fDigitData[localcount].fTime=time;
286         localcount++;
287       } while (digarr.Next());
288
289       Byte_t *tmp = (Byte_t*)tempPt;
290       Int_t size = sizeof(AliL3DigitRowData)
291                                       + ndigits[lrow]*sizeof(AliL3DigitData);
292       tmp += size;
293       tempPt = (AliL3DigitRowData*)tmp;
294     }
295   t->Delete();
296   savedir->cd();
297     
298   return data;
299 }
300
301 void AliL3FileHandler::AliDigits2RootFile(AliL3DigitRowData *rowPt,Char_t *new_digitsfile)
302 {
303   //Write digits to a new alirootfile.
304
305   if(!fInAli)
306     {
307       printf("AliL3FileHandler::AliDigits2RootFile : No rootfile\n");
308       return;
309     }
310   if(!fParam)
311     {
312       printf("AliL3FileHandler::AliDigits2RootFile : No parameter object. Run on rootfile\n");
313       return;
314     }
315   if(!fTransformer)
316     {
317       printf("AliL3FileHandler::AliDigits2RootFile : No transform object\n");
318       return;
319     }
320   
321   //Get the original digitstree:
322   fInAli->cd();
323   AliTPCDigitsArray *old_array = new AliTPCDigitsArray();
324   old_array->Setup(fParam);
325   old_array->SetClass("AliSimDigits");
326   Bool_t ok = old_array->ConnectTree("TreeD_75x40_100x60_0");
327   if(!ok)
328     {
329       printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object\n");
330       return;
331     }
332   
333   Bool_t create=kFALSE;
334   TFile *digFile;
335   
336   digFile = TFile::Open(new_digitsfile,"NEW");
337   if(digFile->IsOpen())
338     {    
339       create = kTRUE;
340       fParam->Write(fParam->GetTitle());
341     }
342   else
343     {
344       LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
345         <<"Rootfile did already exist, so I will just open it for updates"<<ENDLOG;
346       digFile = TFile::Open(new_digitsfile,"UPDATE");
347       create=kFALSE;
348     }
349   if(!digFile->IsOpen())
350     {
351       LOG(AliL3Log::kError,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
352         <<"Error opening rootfile "<<new_digitsfile<<ENDLOG;
353       return;
354     }
355   
356   digFile->cd();
357     
358   //setup a new one, or connect it to the existing one:
359   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
360   arr->SetClass("AliSimDigits");
361   arr->Setup(fParam);
362   if(create)
363     arr->MakeTree();
364   else
365     {
366       Bool_t ok = arr->ConnectTree("TreeD_75x40_100x60_0");
367       if(!ok)
368         {
369           printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object in existing file\n");
370           return;
371         }
372     }
373   for(Int_t i=fRowMin; i<=fRowMax; i++)
374     {
375       
376       if(rowPt->fRow != i) printf("AliL3FileHandler::AliDigits2RootFile : Mismatching row numbering!!!\n");
377             
378       Int_t sector,row;
379       fTransformer->Slice2Sector(fSlice,i,sector,row);
380       AliDigits * dig = arr->CreateRow(sector,row);
381       AliDigits *old_dig = old_array->LoadRow(sector,row);
382       if(!old_dig)
383         printf("AliL3FileHandler::AliDigits2RootFile : No padrow %d %d\n",sector,row);
384       
385       AliL3DigitData *digPt = rowPt->fDigitData;
386       for(UInt_t j=0; j<rowPt->fNDigit; j++)
387         {
388           UShort_t charge = digPt[j].fCharge;
389           UChar_t pad = digPt[j].fPad;
390           UShort_t time = digPt[j].fTime;
391           
392           if(charge == 0) //Only write the digits that has not been removed
393             continue;
394           dig->SetDigitFast(old_dig->GetDigit(time,pad),time,pad);
395           ((AliSimDigits*)dig)->SetTrackIDFast(((AliSimDigits*)old_dig)->GetTrackID((Int_t)time,(Int_t)pad,0),time,pad,0);
396           ((AliSimDigits*)dig)->SetTrackIDFast(((AliSimDigits*)old_dig)->GetTrackID((Int_t)time,(Int_t)pad,1),time,pad,1);
397           ((AliSimDigits*)dig)->SetTrackIDFast(((AliSimDigits*)old_dig)->GetTrackID((Int_t)time,(Int_t)pad,2),time,pad,2);
398           
399         }
400       UpdateRowPointer(rowPt);
401       arr->StoreRow(sector,row);
402       arr->ClearRow(sector,row);  
403       old_array->ClearRow(sector,row);
404     }
405   digFile->cd();
406   char treeName[100];
407   sprintf(treeName,"TreeD_%s_0",fParam->GetTitle());
408   printf("Writing tree to file.....");
409   arr->GetTree()->Write(treeName,TObject::kOverwrite);
410   printf("done\n");
411   digFile->Close();
412   //arr->GetTree()->Delete();
413   //delete arr;
414 }
415
416 ///////////////////////////////////////// Point IO  
417 Bool_t AliL3FileHandler::AliPoints2Binary(){
418   Bool_t out = kTRUE;
419   UInt_t npoint;
420   AliL3SpacePointData *data = AliPoints2Memory(npoint);
421   out = Memory2Binary(npoint,data);
422   Free();
423   return out;
424 }
425
426 AliL3SpacePointData * AliL3FileHandler::AliPoints2Memory(UInt_t & npoint){
427   AliL3SpacePointData *data = 0;
428   npoint=0;
429   if(!fInAli){
430     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
431     <<"No Input avalible: no object TFile"<<ENDLOG;
432     return 0;
433   }
434   if(!fInAli->IsOpen()){
435     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
436     <<"No Input avalible: TFile not opend"<<ENDLOG;
437     return 0;
438   }
439   TDirectory *savedir = gDirectory;
440   fInAli->cd();
441   
442   Char_t cname[100];
443   Int_t eventn = 0;
444   sprintf(cname,"TreeC_TPC_%d",eventn);
445   AliTPCClustersArray carray;
446   carray.Setup(fParam);
447   carray.SetClusterType("AliTPCcluster");
448   Bool_t clusterok = carray.ConnectTree(cname);
449   if(!clusterok) return 0;
450
451   AliTPCClustersRow ** clusterrow = 
452                new AliTPCClustersRow*[ (int)carray.GetTree()->GetEntries()];
453   Int_t *rows = new int[ (int)carray.GetTree()->GetEntries()];
454   Int_t *sects = new int[  (int)carray.GetTree()->GetEntries()];
455   Int_t sum=0;
456
457   Int_t lslice,lrow;
458   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
459     AliSegmentID *s = carray.LoadEntry(i);
460     Int_t sector,row;
461     fParam->AdjustSectorRow(s->GetID(),sector,row);
462     rows[i] = row;
463     sects[i] = sector;
464     clusterrow[i] = 0;
465     fTransformer->Sector2Slice(lslice,lrow,sector,row);
466     if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
467     clusterrow[i] = carray.GetRow(sector,row);
468     if(clusterrow[i])
469       sum+=clusterrow[i]->GetArray()->GetEntriesFast();
470   }
471   UInt_t size = sum*sizeof(AliL3SpacePointData);
472
473   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliPoints2Memory","File")
474   <<AliL3Log::kDec<<"Found "<<sum<<" SpacePoints"<<ENDLOG;
475
476   data = (AliL3SpacePointData *) Allocate(size);
477   npoint = sum;
478   UInt_t n=0; 
479   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
480     if(!clusterrow[i]) continue;
481     Int_t row = rows[i];
482     Int_t sector = sects[i];
483     fTransformer->Sector2Slice(lslice,lrow,sector,row);
484     Int_t entries_in_row = clusterrow[i]->GetArray()->GetEntriesFast();
485     for(Int_t j = 0;j<entries_in_row;j++){
486       AliTPCcluster *c = (AliTPCcluster*)(*clusterrow[i])[j];
487       data[n].fZ = c->GetZ();
488       data[n].fY = c->GetY();
489       data[n].fX = fParam->GetPadRowRadii(sector,row);
490       data[n].fID = n+((fSlice&0x7f)<<25)+((fPatch&0x7)<<22);//uli
491       data[n].fPadRow = lrow;
492       data[n].fXYErr = c->GetSigmaY2();
493       data[n].fZErr = c->GetSigmaZ2();
494       if(fMC) fprintf(fMC,"%d %d\n",data[n].fID,c->GetLabel(0));
495       n++;
496     }
497   }
498   for(Int_t i=0;i<carray.GetTree()->GetEntries();i++){
499     Int_t row = rows[i];
500     Int_t sector = sects[i];
501     if(carray.GetRow(sector,row))
502       carray.ClearRow(sector,row);
503   }
504
505   delete [] clusterrow;
506   delete [] rows;
507   delete [] sects;
508   savedir->cd();   
509
510   return data;
511 }
512