]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFT0.cxx
FindROOT + dictionary generating macro
[u/mrichter/AliRoot.git] / TOF / AliTOFT0.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //
20 // This is a TTask that made the calculation of the Time zero using TOF.
21 // Description: The algorithm used to calculate the time zero of
22 // interaction using TOF detector is the following.
23 // We select in the MonteCarlo some primary particles - or tracks in
24 // the following - that strike the TOF detector (the larger part are
25 // pions, kaons or protons).
26 // We choose a set of 10 selected tracks, for each track You have the
27 // length of the track when the TOF is reached (a standard TOF hit
28 // does not contain this additional information, this is the reason
29 // why we implemented a new time zero dedicated TOF hit class
30 // AliTOFhitT0; in order to store this type of hit You have to use the
31 // AliTOFv4T0 as TOF class in Your Config.C. In AliTOFv4T0 the
32 // StepManager was modified in order to fill the TOF hit branch with
33 // this type of hits; in fact the AliTOF::AddT0Hit is called rather
34 // that the usual AliTOF::AddHit), the momentum at generation (from
35 // TreeK) and the time of flight given by the TOF detector.
36 // (Observe that the ctor of the AliTOF class, when the AliTOFv4T0
37 // class is used, is called with the "tzero" option: it is in order
38 // create the fHits TClonesArray filled with AliTOFhitT0 objects,
39 // rather than with normal AliTOFhit)
40 // Then Momentum and time of flight for each track are smeared
41 // according to known experimental resolution (all sources of error
42 // have been token into account).
43 // Let consider now only one set of 10 tracks (the algorithm is the
44 // same for all sets).
45 // Assuming the (mass) hypothesis that each track can be AUT a pion,
46 // AUT a kaon, AUT a proton, we consider all the 3 at 10 possible
47 // cases.
48 // For each track in each (mass) configuration
49 // (a configuration can be
50 // e.g. pion/pion/kaon/proton/pion/proton/kaon/kaon/pion/pion)
51 // we calculate the time zero (we know in fact the velocity of the
52 // track after the assumption about its mass, the time of flight given
53 // by the TOF, and the corresponding path travelled till the TOF
54 // detector). Then for each mass configuration we have 10 time zero
55 // and we can calculate the ChiSquare for the current configuration
56 // using the weighted mean over all 10 time zero.
57 // We call the best assignment the mass configuration that gives the
58 // minimum value of the ChiSquare.
59 // We plot the weighted mean over all 10 time zero for the best
60 // assignment, the ChiSquare for the best assignment and the
61 // corresponding confidence level.
62 // The strong assumption is the MC selection of primary particles. It
63 // will be introduced in the future also some more realistic
64 // simulation about this point.
65 //
66 // Use case:
67 // root [0] AliTOFT0 * tzero = new AliTOFT0("galice.root")
68 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
69 // root [1] tzero->ExecuteTask()
70 // root [2] tzero->ExecuteTask("tim")
71 //             // available parameters:
72 //             tim - print benchmarking information
73 //             all - print usefull informations about the number of
74 //                   misidentified tracks and a comparison about the
75 //                   true configuration (known from MC) and the best
76 //                   assignment
77 //
78 //-- Author: F. Pierella
79 //
80 //_________________________________________________________________________
81
82 #include <TCanvas.h>
83 #include <TClonesArray.h>
84 #include <TFile.h>
85 //#include <TFolder.h>
86 #include <TFrame.h>
87 #include <TH1.h>
88 #include <TParticle.h>
89 #include <TBenchmark.h>
90 #include <TTask.h>
91 #include <TTree.h>
92 #include <TRandom.h>
93 #include <TROOT.h>
94
95 #include "AliLoader.h"
96 #include "AliMC.h"
97 #include "AliRun.h"
98
99 #include "AliTOFhitT0.h"
100 #include "AliTOFT0.h"
101 #include "AliTOF.h"
102
103 //extern TROOT *gROOT;
104 extern TRandom *gRandom;
105 extern TBenchmark *gBenchmark;
106
107 extern AliRun *gAlice;
108
109 ClassImp(AliTOFT0)
110
111 //____________________________________________________________________________ 
112 AliTOFT0::AliTOFT0():
113   TTask("AliTOFT0",""),
114   fNevents(0),
115   fTimeResolution(0),
116   fLowerMomBound(0),
117   fUpperMomBound(0),
118   fT0File(""),
119   fHeadersFile("")
120 {
121   // ctor
122 }
123            
124 //____________________________________________________________________________ 
125 AliTOFT0::AliTOFT0(char* headerFile, Int_t nEvents):
126   TTask("AliTOFT0",""), 
127   fNevents(nEvents),
128   fTimeResolution(1.2e-10),
129   fLowerMomBound(1.5),
130   fUpperMomBound(2.),
131   fT0File(""),
132   fHeadersFile(headerFile)
133 {
134   //
135   //
136   //
137
138   //fNevents=nEvents ; // Number of events for which calculate the T0, 
139                      // default 0: it means all evens in current file
140   //fLowerMomBound=1.5; // [GeV/c] default value
141   //fUpperMomBound=2. ; // [GeV/c] default value
142   //fTimeResolution   = 1.2e-10; // 120 ps by default   
143   //fHeadersFile = headerFile ;
144
145   TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
146
147   //File was not opened yet
148   if(file == 0){
149     if(fHeadersFile.Contains("rfio"))
150       file =    TFile::Open(fHeadersFile,"update") ;
151     else
152       file = new TFile(fHeadersFile.Data(),"update") ;
153     gAlice = (AliRun *) file->Get("gAlice") ;
154   }
155
156   // add Task to //root/Tasks folder
157   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
158   roottasks->Add(this) ; 
159 }
160
161 //____________________________________________________________________________ 
162 AliTOFT0::AliTOFT0(const AliTOFT0 & tzero):
163   TTask("AliTOFT0",""),
164   fNevents(0),
165   fTimeResolution(0),
166   fLowerMomBound(0),
167   fUpperMomBound(0),
168   fT0File(""),
169   fHeadersFile("")
170 {
171   // copy ctr
172
173 ( (AliTOFT0 &)tzero ).Copy(*this);
174 }
175
176 //____________________________________________________________________________ 
177 AliTOFT0& AliTOFT0::operator = (const AliTOFT0 & tzero)
178 {
179   if (this==&tzero) return *this;
180
181   fNevents=tzero.fNevents;
182   fTimeResolution=tzero.fTimeResolution;
183   fLowerMomBound=tzero.fLowerMomBound;
184   fUpperMomBound=tzero.fUpperMomBound;
185   fT0File=tzero.fT0File;
186   fHeadersFile=tzero.fHeadersFile;
187   return *this;
188
189 }
190
191 //____________________________________________________________________________ 
192   AliTOFT0::~AliTOFT0()
193 {
194   // dtor
195 }
196
197 //____________________________________________________________________________
198 void AliTOFT0::Exec(Option_t *option) 
199
200   //
201   // calculate T0 distribution for all events using chisquare 
202   //
203   Int_t ngood=0;
204   Int_t nmisidentified=0;
205   Int_t nmisidentified0=0;
206   Int_t nmisidentified1=0;
207   Int_t nmisidentified2=0;
208   Int_t nmisidentified3=0;
209   Int_t nmisidentified4=0;
210   Int_t nmisidentified5=0;
211   Int_t nmisidentified6=0;
212   Int_t nmisidentified7=0;
213   Int_t nmisidentified8=0;
214   Int_t nmisidentified9=0;
215   Int_t ipartold = -1;
216   Int_t ipart;
217   Int_t selected=0;
218   Int_t istop=0;
219   Float_t timeresolutioninns=fTimeResolution*(1.e+9); // convert in [ns]
220   const Int_t kUPDATE = 5; // for visual option
221   Int_t itimes=0;
222   TCanvas* c1=0;
223   TCanvas* c2=0;
224   TCanvas* c3=0;
225
226   if(strstr(option,"visual")){
227     // Create a new canvas.
228     //c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,500,500);
229     c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,370,370);
230     c1->SetFillColor(35);
231     c1->GetFrame()->SetFillColor(21);
232     c1->GetFrame()->SetBorderSize(6);
233     c1->GetFrame()->SetBorderMode(-1);
234
235     //c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",550,10,500,500);
236     c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",380,10,370,370);
237     c2->SetFillColor(35);
238     c2->GetFrame()->SetFillColor(21);
239     c2->GetFrame()->SetBorderSize(6);
240     c2->GetFrame()->SetBorderMode(-1);
241
242     //c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",280,550,500,500);
243     c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",760,10,370,370);
244     c3->SetFillColor(35);
245     c3->GetFrame()->SetFillColor(21);
246     c3->GetFrame()->SetBorderSize(6);
247     c3->GetFrame()->SetBorderMode(-1);
248   }
249
250   if(strstr(option,"tim") || strstr(option,"all"))
251     gBenchmark->Start("TOFT0");
252
253   TH1F *htzerobest= new TH1F("htzerobest","T0 for best assignment",200,-1.,1.);
254   TH1F* hchibest  = new TH1F("hchibest","ChiSquare Min Distribution",80,0.,40.);
255   TH1F* hchibestconflevel  = new TH1F("hchibestconflevel","ChiSquare Min Confidence Level",10,0.,1.);
256
257   // setting histo colors
258   if(strstr(option,"visual")){
259     htzerobest->SetFillColor(48);
260     hchibest->SetFillColor(50);
261     hchibestconflevel->SetFillColor(52);
262   }
263
264   Int_t   assparticle[10]={3,3,3,3,3,3,3,3,3,3};
265   Int_t   truparticle[10]={3,3,3,3,3,3,3,3,3,3};
266   Float_t t0best=999.;
267   Float_t timeofflight[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
268   Float_t momentum[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
269   Float_t timezero[10];
270   Float_t weightedtimezero[10];
271   Float_t beta[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
272   Float_t sqMomError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
273   Float_t sqTrackError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
274   Float_t massarray[3]={0.13957,0.493677,0.9382723};
275   Float_t dummychisquare=0.;
276   Float_t chisquare=999.;
277   Float_t tracktoflen[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
278
279   AliTOF *detTOF = (AliTOF *) gAlice->GetDetector ("TOF");
280
281   if (!detTOF) {
282     Error("AliTOFT0","TOF not found");
283     return;
284   }
285
286   if(strstr(option,"all")){
287     AliInfo(Form("Selecting primary tracks with momentum between %f GeV/c and %f GeV/c",  fLowerMomBound, fUpperMomBound));
288     AliInfo("Memorandum: 0 means PION | 1 means KAON | 2 means PROTON");
289   }
290
291   if (fNevents == 0) fNevents = (Int_t) AliRunLoader::Instance()->TreeE()->GetEntries();
292
293   for (Int_t ievent = 0; ievent < fNevents; ievent++) {
294     gAlice->GetEvent(ievent);
295     TTree *hitTree = detTOF->GetLoader()->TreeH ();
296     if (!hitTree)
297       return;
298     TParticle*    particle;
299     AliTOFhitT0*  tofHit;
300     TClonesArray* tofHits = detTOF->Hits();
301
302     Int_t lasttrack=-1;
303     Int_t nset=0;
304
305     hitTree->SetBranchStatus("*",0); // switch off all branches
306     hitTree->SetBranchStatus("TOF*",1); // switch on only TOF
307
308     // Start loop on primary tracks in the hits containers
309
310     Int_t ntracks = static_cast<Int_t>(hitTree->GetEntries());
311     for (Int_t track = 0; track < ntracks; track++)
312     {
313       if(nset>=5) break; // check on the number of set analyzed
314       
315       gAlice->GetMCApp()->ResetHits();
316       hitTree->GetEvent(track);
317
318       //AliMC *mcApplication = (AliMC*)gAlice->GetMCApp();
319
320       //particle = mcApplication->Particle(track);
321       Int_t nhits = tofHits->GetEntriesFast();
322
323       for (Int_t hit = 0; hit < nhits; hit++)
324       {
325         tofHit = (AliTOFhitT0 *) tofHits->UncheckedAt(hit);
326         ipart    = tofHit->GetTrack();
327         // check to discard the case when the same particle is selected more than one
328         // time 
329
330         if (ipart != ipartold){
331           
332           particle = (TParticle*)gAlice->GetMCApp()->Particle(ipart);
333           
334           Float_t idealtime=tofHit->GetTof();
335           //       Float_t time=idealtime;
336           Float_t time   = gRandom->Gaus(idealtime, fTimeResolution);
337           Float_t toflen=tofHit->GetLen();
338           toflen=toflen/100.; // toflen given in m
339           Int_t pdg   = particle->GetPdgCode();
340           Int_t abspdg   =TMath::Abs(pdg);
341           Float_t idealmom  = particle->P();
342           Float_t momres=idealmom*0.025; // 2.5% res token into account for all momenta
343           Float_t mom =gRandom->Gaus(idealmom,momres);
344
345           Bool_t isgoodpart=(abspdg==211 || abspdg==2212 || abspdg==321);
346
347           time*=1.E+9; // tof given in nanoseconds         
348           if (particle->GetFirstMother() < 0 && isgoodpart && mom<=fUpperMomBound && mom>=fLowerMomBound){
349             selected+=1;
350             istop=selected;
351             if(istop>10) break;
352             Int_t index=selected-1;
353             timeofflight[index]=time;
354             tracktoflen[index]=toflen;
355             momentum[index]=mom;
356             //AliInfo(Form(" %d  %d  %d ", timeofflight[index], tracktoflen[index], momentum[index]));
357             switch (abspdg) {
358             case 211:
359               truparticle[index]=0;
360               break ;
361             case 321:
362               truparticle[index]=1;
363               break ;
364             case 2212: 
365               truparticle[index]=2;
366               break ;
367             }
368             
369           }
370           ipartold = ipart;
371           
372           if(istop==10){ // start analysis on current set
373             nset+=1;
374             lasttrack=track;
375             istop=0;
376             selected=0;
377             //AliInfo("starting t0 calculation for current set");
378             for (Int_t i1=0; i1<3;i1++) {
379               beta[0]=momentum[0]/sqrt(massarray[i1]*massarray[i1]+momentum[0]*momentum[0]);
380               for (Int_t i2=0; i2<3;i2++) { 
381                 beta[1]=momentum[1]/sqrt(massarray[i2]*massarray[i2]+momentum[1]*momentum[1]);
382                 for (Int_t i3=0; i3<3;i3++) {
383                   beta[2]=momentum[2]/sqrt(massarray[i3]*massarray[i3]+momentum[2]*momentum[2]);
384                   for (Int_t i4=0; i4<3;i4++) {
385                     beta[3]=momentum[3]/sqrt(massarray[i4]*massarray[i4]+momentum[3]*momentum[3]);
386                     for (Int_t i5=0; i5<3;i5++) {
387                       beta[4]=momentum[4]/sqrt(massarray[i5]*massarray[i5]+momentum[4]*momentum[4]);
388                       for (Int_t i6=0; i6<3;i6++) {
389                         beta[5]=momentum[5]/sqrt(massarray[i6]*massarray[i6]+momentum[5]*momentum[5]);
390                         for (Int_t i7=0; i7<3;i7++) { 
391                           beta[6]=momentum[6]/sqrt(massarray[i7]*massarray[i7]+momentum[6]*momentum[6]);
392                           for (Int_t i8=0; i8<3;i8++) {
393                             beta[7]=momentum[7]/sqrt(massarray[i8]*massarray[i8]+momentum[7]*momentum[7]);
394                             for (Int_t i9=0; i9<3;i9++) {
395                               beta[8]=momentum[8]/sqrt(massarray[i9]*massarray[i9]+momentum[8]*momentum[8]);
396                               for (Int_t i10=0; i10<3;i10++) {  
397                                 beta[9]=momentum[9]/sqrt(massarray[i10]*massarray[i10]+momentum[9]*momentum[9]);
398                                 
399                                 Float_t meantzero=0.;
400                                 Float_t sumAllweights=0.;
401                                 for (Int_t itz=0; itz<10;itz++) {
402                                   sqMomError[itz]=((1.-beta[itz]*beta[itz])*0.025)*((1.-beta[itz]*beta[itz])*0.025)*(tracktoflen[itz]/(0.299792*beta[itz]))*(tracktoflen[itz]/(0.299792*beta[itz])); // this gives the square of the momentum error in nanoseconds
403                                   sqTrackError[itz]=(timeresolutioninns*timeresolutioninns+sqMomError[itz]); // total error for the current track
404                                   sumAllweights+=1./sqTrackError[itz];
405
406                                   timezero[itz]=(tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz];
407                                   weightedtimezero[itz]=((tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz])/sqTrackError[itz];// weighted time zero for current track
408                                   meantzero+=weightedtimezero[itz];
409                                 } // end loop for (Int_t itz=0; itz<10;itz++)
410                                 meantzero=meantzero/sumAllweights; // it is given in [ns]
411                                 
412                                 dummychisquare=0.;
413                                 // calculate the chisquare for the current assignment
414                                 for (Int_t icsq=0; icsq<10;icsq++) {
415                                   dummychisquare+=(timezero[icsq]-meantzero)*(timezero[icsq]-meantzero)/sqTrackError[icsq];
416                                 } // end loop for (Int_t icsq=0; icsq<10;icsq++) 
417
418                                 if(dummychisquare<=chisquare){
419                                   assparticle[0]=i1;
420                                   assparticle[1]=i2;
421                                   assparticle[2]=i3;
422                                   assparticle[3]=i4;
423                                   assparticle[4]=i5;
424                                   assparticle[5]=i6;
425                                   assparticle[6]=i7;
426                                   assparticle[7]=i8;
427                                   assparticle[8]=i9;
428                                   assparticle[9]=i10;
429                                   chisquare=dummychisquare;
430                                   t0best=meantzero;
431                                 } // close if(dummychisquare<=chisquare)
432                                 
433                               } // end loop on i10
434                             } // end loop on i9
435                           } // end loop on i8
436                         } // end loop on i7
437                       } // end loop on i6
438                     } // end loop on i5
439                   } // end loop on i4
440                 } // end loop on i3
441               } // end loop on i2
442             } // end loop on i1
443
444             if(truparticle[0]==assparticle[0] && truparticle[1]==assparticle[1] && truparticle[2]==assparticle[2]  && truparticle[3]==assparticle[3] && truparticle[4]==assparticle[4]&& truparticle[5]==assparticle[5] && truparticle[6]==assparticle[6] && truparticle[7]==assparticle[7]  && truparticle[8]==assparticle[8] && truparticle[9]==assparticle[9]) ngood+=1;
445             if(truparticle[0]!=assparticle[0]) nmisidentified0+=1;
446             if(truparticle[1]!=assparticle[1]) nmisidentified1+=1;
447             if(truparticle[2]!=assparticle[2]) nmisidentified2+=1;
448             if(truparticle[3]!=assparticle[3]) nmisidentified3+=1;
449             if(truparticle[4]!=assparticle[4]) nmisidentified4+=1;
450             if(truparticle[5]!=assparticle[5]) nmisidentified5+=1;
451             if(truparticle[6]!=assparticle[6]) nmisidentified6+=1;
452             if(truparticle[7]!=assparticle[7]) nmisidentified7+=1;
453             if(truparticle[8]!=assparticle[8]) nmisidentified8+=1;
454             if(truparticle[9]!=assparticle[9]) nmisidentified9+=1;
455             // filling histos
456             htzerobest->Fill(t0best);
457             hchibest->Fill(chisquare);
458             Double_t dblechisquare=(Double_t)chisquare;
459             Float_t confLevel=(Float_t)TMath::Prob(dblechisquare,9); // ndf 10-1=9
460             hchibestconflevel->Fill(confLevel);
461             itimes++;
462             if(strstr(option,"all")){
463               AliInfo(Form("True Assignment %d  %d  %d  %d  %d  %d  %d  %d  %d  %d", truparticle[0], truparticle[1], truparticle[2], truparticle[3], truparticle[4], truparticle[5], truparticle[6], truparticle[7], truparticle[8], truparticle[9]));
464               AliInfo(Form("Best Assignment %d  %d  %d  %d  %d  %d  %d  %d  %d  %d", assparticle[0], assparticle[1], assparticle[2], assparticle[3], assparticle[4], assparticle[5], assparticle[6], assparticle[7], assparticle[8], assparticle[9]));
465               AliInfo(Form("Minimum ChiSquare for current set   %f ", chisquare));
466               AliInfo(Form("Confidence Level (Minimum ChiSquare) %f", confLevel));
467             }
468             if (strstr(option,"visual") && itimes && (itimes%kUPDATE) == 0) {
469               if (itimes == kUPDATE){
470                 c1->cd();
471                 htzerobest->Draw();
472                 c2->cd();
473                 hchibest->Draw();
474                 c3->cd();
475                 hchibestconflevel->Draw();
476               }
477               c1->Modified();
478               c1->Update();
479               c2->Modified();
480               c2->Update();
481               c3->Modified();
482               c3->Update();
483               if (gSystem->ProcessEvents())
484                 break;
485             }
486             chisquare=999.;
487             t0best=999.;
488             
489           } // end for the current set. close if(istop==5)
490         } // end condition on ipartold
491       } // end loop on hits for the current track
492       if(istop>=10) break;
493     } // end loop on ntracks  
494   } //event loop
495   
496   if(strstr(option,"all")){
497     nmisidentified=(nmisidentified0+nmisidentified1+nmisidentified2+nmisidentified3+nmisidentified4+nmisidentified5+nmisidentified6+nmisidentified7+nmisidentified8+nmisidentified9);
498     AliInfo(Form("total number of tracks token into account  %i", 10*5*fNevents));
499     Float_t badPercentage=100.*(Float_t)nmisidentified/(10*5*fNevents);
500     AliInfo(Form("total misidentified                       %i (%f %%) ", nmisidentified, badPercentage));
501     AliInfo(Form("Total Number of set token into account     %i", 5*fNevents));
502     Float_t goodSetPercentage=100.*(Float_t)ngood/(5*fNevents);
503     AliInfo(Form("Number of set with no misidentified tracks %i (%f %%)", ngood, goodSetPercentage));
504   }
505
506   // free used memory for canvas
507   delete c1; c1=0;
508   delete c2; c2=0;
509   delete c3; c3=0;
510
511   // generating output filename only if not previously specified using SetTZeroFile
512
513   /*const*/ Int_t kSize = 70+fHeadersFile.Length()+1;
514   char outFileName[kSize];
515   strncpy(outFileName,"ht010tr120ps",kSize); // global time resolution has to be converted from Int_t to char
516                                       // in order to have in the output filename this parameter
517   strncat(outFileName,fHeadersFile,kSize);
518
519   if(fT0File.IsNull()) fT0File=outFileName;
520
521   TFile* houtfile = new TFile(fT0File,"recreate");
522   houtfile->cd();
523   htzerobest->Write(0,TObject::kOverwrite);
524   hchibest->Write(0,TObject::kOverwrite);
525   hchibestconflevel->Write(0,TObject::kOverwrite);
526   houtfile->Close();  
527   
528   
529   if(strstr(option,"tim") || strstr(option,"all")){
530     gBenchmark->Stop("TOFT0");
531     AliInfo("AliTOFT0:");
532     /*
533     cout << "   took " << gBenchmark->GetCpuTime("TOFT0") << " seconds in order to calculate T0 " 
534          << gBenchmark->GetCpuTime("TOFT0")/fNevents << " seconds per event " << endl ;
535     */
536     gBenchmark->Print("TOFT0");
537   }
538 }
539  
540 //__________________________________________________________________
541 void AliTOFT0::SetTZeroFile(char * file )
542 {
543   //
544   // Set T0 file name
545   //
546   printf("Destination file : %s \n", file) ;
547   fT0File=file;
548
549 }
550
551 //__________________________________________________________________
552 void AliTOFT0::Print(Option_t* /*option*/)const
553 {
554   //
555   // Print class content
556   //
557   printf("------------------- %s -------------\n", GetName()) ;
558   if(!fT0File.IsNull())
559     printf("  Writing T0 Distribution to file  %s \n",(char*) fT0File.Data());
560
561 }
562
563 //__________________________________________________________________
564 Bool_t AliTOFT0::operator==( AliTOFT0 const &tzero )const
565 {
566   //
567   // Equal operator
568   // 
569
570   if( (fTimeResolution==tzero.fTimeResolution)&&
571       (fLowerMomBound==tzero.fLowerMomBound)&&
572       (fUpperMomBound==tzero.fUpperMomBound))
573     return kTRUE ;
574   else
575     return kFALSE ;
576 }