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