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