3 // Author: Uli Frankenfeld <mailto:franken@fi.uib.no>, Anders Vestbo <mailto:vestbo$fi.uib.no>
4 //*-- Copyright © Uli
9 #include "AliL3Transform.h"
10 #include "AliL3Logging.h"
11 #include "AliL3MemHandler.h"
12 #include "AliL3FileHandler.h"
14 #include "AliTPCDigitsArray.h"
15 #include "AliTPCClustersArray.h"
16 #include "AliTPCcluster.h"
17 #include "AliTPCClustersRow.h"
19 #include "AliL3DigitData.h"
20 #include "AliL3TrackSegmentData.h"
21 #include "AliL3SpacePointData.h"
22 #include "AliL3TrackArray.h"
24 //_____________________________________________________________
27 // The HLT ROOT <-> binary files handling class
29 // This class provides the interface between AliROOT files,
30 // and HLT binary files. It should be used for converting
31 // TPC data stored in AliROOT format (outputfile from a simulation),
32 // into the data format currently used by in the HLT framework.
33 // This enables the possibility to always use the same data format,
34 // whether you are using a binary file as an input, or a AliROOT file.
36 // For example on how to create binary files from a AliROOT simulation,
37 // see example macro exa/Binary.C.
39 // For reading a AliROOT file into HLT format in memory, do the following:
41 // AliL3FileHandler file;
42 // file.SetAliInput("galice.root");
43 // AliL3DigitRowData *dataPt = (AliL3DigitRowData*)file.AliDigits2Memory(nrows,eventnr);
45 // All the data are then stored in memory and accessible via the pointer dataPt.
46 // Accesing the data is then identical to the example 1) showed in AliL3MemHandler class.
48 // For converting the data back, and writing it to a new AliROOT file do:
50 // AliL3FileHandler file;
51 // file.SetAliInput("galice.root");
52 // file.Init(slice,patch,NumberOfRowsInPatch);
53 // file.AliDigits2RootFile(dataPt,"new_galice.root");
54 // file.CloseAliInput();
56 ClassImp(AliL3FileHandler)
58 AliL3FileHandler::AliL3FileHandler()
69 AliL3FileHandler::~AliL3FileHandler()
72 if(fMC) CloseMCOutput();
73 if(fDigitsTree) delete fDigitsTree;
74 if(fInAli) CloseAliInput();
78 Bool_t AliL3FileHandler::SetMCOutput(char *name)
80 fMC = fopen(name,"w");
82 LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
83 <<"Pointer to File = 0x0 "<<ENDLOG;
89 Bool_t AliL3FileHandler::SetMCOutput(FILE *file)
93 LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
94 <<"Pointer to File = 0x0 "<<ENDLOG;
100 void AliL3FileHandler::CloseMCOutput()
103 LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseMCOutPut","File Close")
104 <<"Nothing to Close"<<ENDLOG;
111 Bool_t AliL3FileHandler::SetAliInput()
113 if(!fInAli->IsOpen()){
114 LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
115 <<"Ali File "<<fInAli->GetName()<<" does not exist"<<ENDLOG;
118 fParam = (AliTPCParam*)fInAli->Get("75x40_100x60");
120 LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
121 <<"No AliTPCParam 75x40_100x60 in File "<<fInAli->GetName()<<ENDLOG;
127 Bool_t AliL3FileHandler::SetAliInput(char *name)
129 //Open the AliROOT file with name.
131 fInAli= new TFile(name,"READ");
133 LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
134 <<"Pointer to TFile = 0x0 "<<ENDLOG;
137 return SetAliInput();
140 Bool_t AliL3FileHandler::SetAliInput(TFile *file)
142 //Specify already opened AliROOT file to use as an input.
146 LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
147 <<"Pointer to TFile = 0x0 "<<ENDLOG;
150 return SetAliInput();
153 void AliL3FileHandler::CloseAliInput()
156 LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseAliInput","File Close")
157 <<"Nothing to Close"<<ENDLOG;
160 if(fInAli->IsOpen()) fInAli->Close();
166 Bool_t AliL3FileHandler::IsDigit()
168 //Check if there is a TPC digit tree in the current file.
169 //Return kTRUE if tree was found, and kFALSE if not found.
172 LOG(AliL3Log::kWarning,"AliL3FileHandler::IsDigit","File")
173 <<"Pointer to TFile = 0x0 "<<ENDLOG;
174 return kTRUE; //may you are use binary input which is Digits!!
176 TTree *t=(TTree*)fInAli->Get("TreeD_75x40_100x60_0");
178 LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
179 <<"Found Digit Tree -> Use Fast Cluster Finder"<<ENDLOG;
183 LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
184 <<"No Digit Tree -> Use Cluster Tree"<<ENDLOG;
189 ///////////////////////////////////////// Digit IO
190 Bool_t AliL3FileHandler::AliDigits2Binary(Int_t event)
195 AliL3DigitRowData* data = AliDigits2Memory(nrow,event);
196 out = Memory2Binary(nrow,data);
202 Bool_t AliL3FileHandler::AliDigits2CompBinary(Int_t event)
204 //Convert AliROOT TPC data, into HLT data format.
205 //event specifies the event you want in the aliroot file.
209 AliL3DigitRowData* digits = AliDigits2Memory(ndigits,event);
210 out = Memory2CompBinary(ndigits,digits);
216 AliL3DigitRowData * AliL3FileHandler::AliDigits2Memory(UInt_t & nrow,Int_t event)
218 //Read data from AliROOT file into memory, and store it in the HLT data format.
219 //Returns a pointer to the data.
221 AliL3DigitRowData *data = 0;
224 LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
225 <<"No Input avalible: no object TFile"<<ENDLOG;
228 if(!fInAli->IsOpen()){
229 LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
230 <<"No Input avalible: TFile not opend"<<ENDLOG;
235 LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","Transformer")
236 <<"No transformer object"<<ENDLOG;
241 GetDigitsTree(event);
244 Int_t time,pad,sector,row;
247 Int_t entries = (Int_t)fDigitsTree->GetEntries();
248 Int_t ndigits[entries];
251 for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
253 fDigitsTree->GetEvent(n);
254 fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
255 fTransformer->Sector2Slice(lslice,lrow,sector,row);
256 //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
257 if(lslice < fSlice) continue;
258 if(lslice != fSlice) break;
259 if(lrow < fRowMin) continue;
260 if(lrow > fRowMax) break;
266 time=fDigits->CurrentRow();
267 pad=fDigits->CurrentColumn();
268 dig = fDigits->GetDigit(time,pad);
269 if(dig<=fParam->GetZeroSup()) continue;
270 if(time < fParam->GetMaxTBin()-1 && time > 0)
271 if(fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()
272 && fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup())
275 fTransformer->Raw2Local(xyz,sector,row,pad,time);
276 if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
279 ndigits[lrow]++; //for this row only
280 ndigitcount++; //total number of digits to be published
282 } while (fDigits->Next());
286 Int_t size = sizeof(AliL3DigitData)*ndigitcount
287 + nrows*sizeof(AliL3DigitRowData);
289 LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2Memory","Digits")
290 <<AliL3Log::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
292 data=(AliL3DigitRowData*) Allocate(size);
293 nrow = (UInt_t)nrows;
294 AliL3DigitRowData *tempPt = data;
295 for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
297 fDigitsTree->GetEvent(n);
299 fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
300 fTransformer->Sector2Slice(lslice,lrow,sector,row);
301 //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
302 if(lslice < fSlice) continue;
303 if(lslice != fSlice) break;
304 if(lrow < fRowMin) continue;
305 if(lrow > fRowMax) break;
308 tempPt->fNDigit = ndigits[lrow];
313 //dig=fDigits->CurrentDigit();
314 //if (dig<=fParam->GetZeroSup()) continue;
315 time=fDigits->CurrentRow();
316 pad=fDigits->CurrentColumn();
317 dig = fDigits->GetDigit(time,pad);
318 if (dig <= fParam->GetZeroSup()) continue;
319 if(time < fParam->GetMaxTBin()-1 && time > 0)
320 if(fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup() &&
321 fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()) continue;
323 //Exclude data outside cone:
324 fTransformer->Raw2Local(xyz,sector,row,pad,time);
325 if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
328 if(localcount >= ndigits[lrow])
329 LOG(AliL3Log::kFatal,"AliL3FileHandler::AliDigits2Binary","Memory")
330 <<AliL3Log::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
331 <<ndigits[lrow]<<ENDLOG;
333 tempPt->fDigitData[localcount].fCharge=dig;
334 tempPt->fDigitData[localcount].fPad=pad;
335 tempPt->fDigitData[localcount].fTime=time;
337 tempPt->fDigitData[localcount].fTrackID[0] = fDigits->GetTrackID(time,pad,0);
338 tempPt->fDigitData[localcount].fTrackID[1] = fDigits->GetTrackID(time,pad,1);
339 tempPt->fDigitData[localcount].fTrackID[2] = fDigits->GetTrackID(time,pad,2);
342 } while (fDigits->Next());
344 Byte_t *tmp = (Byte_t*)tempPt;
345 Int_t size = sizeof(AliL3DigitRowData)
346 + ndigits[lrow]*sizeof(AliL3DigitData);
348 tempPt = (AliL3DigitRowData*)tmp;
355 Bool_t AliL3FileHandler::GetDigitsTree(Int_t event)
357 //Connects to the TPC digit tree in the AliROOT file.
361 sprintf(dname,"TreeD_75x40_100x60_%d",event);
362 fDigitsTree = (TTree*)fInAli->Get("TreeD_75x40_100x60_0");
365 LOG(AliL3Log::kError,"AliL3FileHandler::GetDigitsTree","Digits Tree")
366 <<AliL3Log::kHex<<"Error getting digitstree "<<(Int_t)fDigitsTree<<ENDLOG;
369 fDigitsTree->GetBranch("Segment")->SetAddress(&fDigits);
373 void AliL3FileHandler::AliDigits2RootFile(AliL3DigitRowData *rowPt,Char_t *new_digitsfile)
375 //Write the data stored in rowPt, into a new AliROOT file.
376 //The data is stored in the AliROOT format
377 //This is specially a nice thing if you have modified data, and wants to run it
378 //through the offline reconstruction chain.
379 //The arguments is a pointer to the data, and the name of the new AliROOT file.
380 //Remember to pass the original AliROOT file (the one that contains the original
381 //simulated data) to this object, in order to retrieve the MC id's of the digits.
385 printf("AliL3FileHandler::AliDigits2RootFile : No rootfile\n");
390 printf("AliL3FileHandler::AliDigits2RootFile : No parameter object. Run on rootfile\n");
395 printf("AliL3FileHandler::AliDigits2RootFile : No transform object\n");
399 //Get the original digitstree:
401 AliTPCDigitsArray *old_array = new AliTPCDigitsArray();
402 old_array->Setup(fParam);
403 old_array->SetClass("AliSimDigits");
404 Bool_t ok = old_array->ConnectTree("TreeD_75x40_100x60_0");
407 printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object\n");
411 Bool_t create=kFALSE;
414 digFile = TFile::Open(new_digitsfile,"NEW");
415 if(digFile->IsOpen())
418 fParam->Write(fParam->GetTitle());
422 LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
423 <<"Rootfile did already exist, so I will just open it for updates"<<ENDLOG;
424 digFile = TFile::Open(new_digitsfile,"UPDATE");
427 if(!digFile->IsOpen())
429 LOG(AliL3Log::kError,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
430 <<"Error opening rootfile "<<new_digitsfile<<ENDLOG;
436 //setup a new one, or connect it to the existing one:
437 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
438 arr->SetClass("AliSimDigits");
444 Bool_t ok = arr->ConnectTree("TreeD_75x40_100x60_0");
447 printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object in existing file\n");
453 for(Int_t i=fRowMin; i<=fRowMax; i++)
456 if((Int_t)rowPt->fRow != i) printf("AliL3FileHandler::AliDigits2RootFile : Mismatching row numbering!!!\n");
459 fTransformer->Slice2Sector(fSlice,i,sector,row);
460 AliSimDigits * dig = (AliSimDigits*)arr->CreateRow(sector,row);
461 AliSimDigits *old_dig = (AliSimDigits*)old_array->LoadRow(sector,row);
463 printf("AliL3FileHandler::AliDigits2RootFile : No padrow %d %d\n",sector,row);
465 AliL3DigitData *digPt = rowPt->fDigitData;
467 for(UInt_t j=0; j<rowPt->fNDigit; j++)
469 Int_t charge = (Int_t)digPt[j].fCharge;
470 Int_t pad = (Int_t)digPt[j].fPad;
471 Int_t time = (Int_t)digPt[j].fTime;
473 if(charge == 0) //Only write the digits that has not been removed
476 dig->SetDigitFast(charge,time,pad);
478 Int_t trackID[3] = {old_dig->GetTrackID(time,pad,0),old_dig->GetTrackID(time,pad,1),old_dig->GetTrackID(time,pad,2)};
480 Int_t s_time = time - 1;
481 while(trackID[0] < 0)
483 if(s_time >= 0 && s_time < fTransformer->GetNTimeBins() && s_pad >= 0 && s_pad < fTransformer->GetNPads(i))
485 if(old_dig->GetTrackID(s_time,s_pad,0) > 0)
487 trackID[0]=old_dig->GetTrackID(s_time,s_pad,0);
488 trackID[1]=old_dig->GetTrackID(s_time,s_pad,1);
489 trackID[2]=old_dig->GetTrackID(s_time,s_pad,2);
492 if(s_pad == pad && s_time == time - 1)
494 else if(s_pad == pad && s_time == time + 1)
495 {s_pad = pad - 1; s_time = time;}
496 else if(s_pad == pad - 1 && s_time == time)
498 else if(s_pad == pad - 1 && s_time == time - 1)
500 else if(s_pad == pad - 1 && s_time == time + 1)
501 {s_pad = pad + 1; s_time = time;}
502 else if(s_pad == pad + 1 && s_time == time)
504 else if(s_pad == pad + 1 && s_time == time - 1)
510 dig->SetTrackIDFast(trackID[0],time,pad,0);
511 dig->SetTrackIDFast(trackID[1],time,pad,1);
512 dig->SetTrackIDFast(trackID[2],time,pad,2);
515 //cout<<"Wrote "<<digcounter<<" on row "<<i<<endl;
516 UpdateRowPointer(rowPt);
517 arr->StoreRow(sector,row);
518 arr->ClearRow(sector,row);
519 old_array->ClearRow(sector,row);
523 sprintf(treeName,"TreeD_%s_0",fParam->GetTitle());
524 printf("Writing tree to file.....");
525 arr->GetTree()->Write(treeName,TObject::kOverwrite);
528 //arr->GetTree()->Delete();
532 ///////////////////////////////////////// Point IO
533 Bool_t AliL3FileHandler::AliPoints2Binary(){
536 AliL3SpacePointData *data = AliPoints2Memory(npoint);
537 out = Memory2Binary(npoint,data);
542 AliL3SpacePointData * AliL3FileHandler::AliPoints2Memory(UInt_t & npoint){
543 AliL3SpacePointData *data = 0;
546 LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
547 <<"No Input avalible: no object TFile"<<ENDLOG;
550 if(!fInAli->IsOpen()){
551 LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
552 <<"No Input avalible: TFile not opend"<<ENDLOG;
557 LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","Transformer")
558 <<"No transformer object"<<ENDLOG;
561 TDirectory *savedir = gDirectory;
566 sprintf(cname,"TreeC_TPC_%d",eventn);
567 AliTPCClustersArray carray;
568 carray.Setup(fParam);
569 carray.SetClusterType("AliTPCcluster");
570 Bool_t clusterok = carray.ConnectTree(cname);
571 if(!clusterok) return 0;
573 AliTPCClustersRow ** clusterrow =
574 new AliTPCClustersRow*[ (int)carray.GetTree()->GetEntries()];
575 Int_t *rows = new int[ (int)carray.GetTree()->GetEntries()];
576 Int_t *sects = new int[ (int)carray.GetTree()->GetEntries()];
580 for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
581 AliSegmentID *s = carray.LoadEntry(i);
583 fParam->AdjustSectorRow(s->GetID(),sector,row);
587 fTransformer->Sector2Slice(lslice,lrow,sector,row);
588 if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
589 clusterrow[i] = carray.GetRow(sector,row);
591 sum+=clusterrow[i]->GetArray()->GetEntriesFast();
593 UInt_t size = sum*sizeof(AliL3SpacePointData);
595 LOG(AliL3Log::kDebug,"AliL3FileHandler::AliPoints2Memory","File")
596 <<AliL3Log::kDec<<"Found "<<sum<<" SpacePoints"<<ENDLOG;
598 data = (AliL3SpacePointData *) Allocate(size);
601 for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
602 if(!clusterrow[i]) continue;
604 Int_t sector = sects[i];
605 fTransformer->Sector2Slice(lslice,lrow,sector,row);
606 Int_t entries_in_row = clusterrow[i]->GetArray()->GetEntriesFast();
607 for(Int_t j = 0;j<entries_in_row;j++){
608 AliTPCcluster *c = (AliTPCcluster*)(*clusterrow[i])[j];
609 data[n].fZ = c->GetZ();
610 data[n].fY = c->GetY();
611 data[n].fX = fParam->GetPadRowRadii(sector,row);
612 data[n].fID = n+((fSlice&0x7f)<<25)+((fPatch&0x7)<<22);//uli
613 data[n].fPadRow = lrow;
614 data[n].fXYErr = c->GetSigmaY2();
615 data[n].fZErr = c->GetSigmaZ2();
616 if(fMC) fprintf(fMC,"%d %d\n",data[n].fID,c->GetLabel(0));
620 for(Int_t i=0;i<carray.GetTree()->GetEntries();i++){
622 Int_t sector = sects[i];
623 if(carray.GetRow(sector,row))
624 carray.ClearRow(sector,row);
627 delete [] clusterrow;