1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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.
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
61 //-- Author: F. Pierella
62 //////////////////////////////////////////////////////////////////////////////
65 #include "AliTOFhitT0.h"
67 #include "AliTOFv4T0.h"
69 #include "AliDetector.h"
81 #include "TBenchmark.h"
82 #include "TParticle.h"
83 #include "TClonesArray.h"
85 #include <Riostream.h>
89 //____________________________________________________________________________
90 AliTOFT0::AliTOFT0():TTask("AliTOFT0","")
96 //____________________________________________________________________________
97 AliTOFT0::AliTOFT0(char* headerFile, Int_t nEvents):TTask("AliTOFT0","")
99 fNevents=nEvents ; // Number of events for which calculate the T0,
100 // default 0: it means all evens in current file
101 fLowerMomBound=1.5; // [GeV/c] default value
102 fUpperMomBound=2. ; // [GeV/c] default value
103 fTimeResolution = 1.2e-10; // 120 ps by default
104 fHeadersFile = headerFile ;
106 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
108 //File was not opened yet
110 if(fHeadersFile.Contains("rfio"))
111 file = TFile::Open(fHeadersFile,"update") ;
113 file = new TFile(fHeadersFile.Data(),"update") ;
114 gAlice = (AliRun *) file->Get("gAlice") ;
117 // add Task to //root/Tasks folder
118 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
119 roottasks->Add(this) ;
122 //____________________________________________________________________________
123 AliTOFT0::~AliTOFT0()
128 //____________________________________________________________________________
129 void AliTOFT0::Exec(Option_t *option)
132 // calculate T0 distribution for all events using chisquare
135 Int_t nmisidentified=0;
136 Int_t nmisidentified0=0;
137 Int_t nmisidentified1=0;
138 Int_t nmisidentified2=0;
139 Int_t nmisidentified3=0;
140 Int_t nmisidentified4=0;
141 Int_t nmisidentified5=0;
142 Int_t nmisidentified6=0;
143 Int_t nmisidentified7=0;
144 Int_t nmisidentified8=0;
145 Int_t nmisidentified9=0;
150 Float_t timeresolutioninns=fTimeResolution*(1.e+9); // convert in [ns]
151 const Int_t kUPDATE = 5; // for visual option
157 if(strstr(option,"visual")){
158 // Create a new canvas.
159 //c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,500,500);
160 c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,370,370);
161 c1->SetFillColor(35);
162 c1->GetFrame()->SetFillColor(21);
163 c1->GetFrame()->SetBorderSize(6);
164 c1->GetFrame()->SetBorderMode(-1);
166 //c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",550,10,500,500);
167 c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",380,10,370,370);
168 c2->SetFillColor(35);
169 c2->GetFrame()->SetFillColor(21);
170 c2->GetFrame()->SetBorderSize(6);
171 c2->GetFrame()->SetBorderMode(-1);
173 //c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",280,550,500,500);
174 c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",760,10,370,370);
175 c3->SetFillColor(35);
176 c3->GetFrame()->SetFillColor(21);
177 c3->GetFrame()->SetBorderSize(6);
178 c3->GetFrame()->SetBorderMode(-1);
181 if(strstr(option,"tim") || strstr(option,"all"))
182 gBenchmark->Start("TOFT0");
184 TH1F *htzerobest= new TH1F("htzerobest","T0 for best assignment",200,-1.,1.);
185 TH1F* hchibest = new TH1F("hchibest","ChiSquare Min Distribution",80,0.,40.);
186 TH1F* hchibestconflevel = new TH1F("hchibestconflevel","ChiSquare Min Confidence Level",10,0.,1.);
188 // setting histo colors
189 if(strstr(option,"visual")){
190 htzerobest->SetFillColor(48);
191 hchibest->SetFillColor(50);
192 hchibestconflevel->SetFillColor(52);
195 Int_t assparticle[10]={3,3,3,3,3,3,3,3,3,3};
196 Int_t truparticle[10]={3,3,3,3,3,3,3,3,3,3};
198 Float_t timeofflight[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
199 Float_t momentum[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
200 Float_t timezero[10];
201 Float_t weightedtimezero[10];
202 Float_t beta[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
203 Float_t sqMomError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
204 Float_t sqTrackError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
205 Float_t massarray[3]={0.13957,0.493677,0.9382723};
206 Float_t dummychisquare=0.;
207 Float_t chisquare=999.;
208 Float_t tracktoflen[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
210 AliTOF *TOF = (AliTOF *) gAlice->GetDetector ("TOF");
213 Error("AliTOFT0","TOF not found");
217 if(strstr(option,"all")){
218 cout << "Selecting primary tracks with momentum between " << fLowerMomBound << " GeV/c and " << fUpperMomBound << " GeV/c" << endl;
219 cout << "Memorandum: 0 means PION | 1 means KAON | 2 means PROTON" << endl;
222 if (fNevents == 0) fNevents = (Int_t) gAlice->TreeE()->GetEntries();
224 for (Int_t ievent = 0; ievent < fNevents; ievent++) {
225 gAlice->GetEvent(ievent);
226 TTree *TH = gAlice->TreeH ();
231 TClonesArray* TOFhits = TOF->Hits();
236 TH->SetBranchStatus("*",0); // switch off all branches
237 TH->SetBranchStatus("TOF*",1); // switch on only TOF
239 // Start loop on primary tracks in the hits containers
241 Int_t ntracks = static_cast<Int_t>(TH->GetEntries());
242 for (Int_t track = 0; track < ntracks; track++)
244 if(nset>=5) break; // check on the number of set analyzed
248 particle = gAlice->Particle(track);
249 Int_t nhits = TOFhits->GetEntriesFast();
251 for (Int_t hit = 0; hit < nhits; hit++)
253 tofHit = (AliTOFhitT0 *) TOFhits->UncheckedAt(hit);
254 ipart = tofHit->GetTrack();
255 // check to discard the case when the same particle is selected more than one
258 if (ipart != ipartold){
260 particle = (TParticle*)gAlice->Particle(ipart);
262 Float_t idealtime=tofHit->GetTof();
263 // Float_t time=idealtime;
264 Float_t time = gRandom->Gaus(idealtime, fTimeResolution);
265 Float_t toflen=tofHit->GetLen();
266 toflen=toflen/100.; // toflen given in m
267 Int_t pdg = particle->GetPdgCode();
268 Int_t abspdg =TMath::Abs(pdg);
269 Float_t idealmom = particle->P();
270 Float_t momres=idealmom*0.025; // 2.5% res token into account for all momenta
271 Float_t mom =gRandom->Gaus(idealmom,momres);
273 Bool_t isgoodpart=(abspdg==211 || abspdg==2212 || abspdg==321);
275 time*=1.E+9; // tof given in nanoseconds
276 if (particle->GetFirstMother() < 0 && isgoodpart && mom<=fUpperMomBound && mom>=fLowerMomBound){
280 Int_t index=selected-1;
281 timeofflight[index]=time;
282 tracktoflen[index]=toflen;
284 // cout << timeofflight[index] << " " << tracktoflen[index] << " " << momentum[index] << endl;
287 truparticle[index]=0;
290 truparticle[index]=1;
293 truparticle[index]=2;
300 if(istop==10){ // start analysis on current set
305 //cout << "starting t0 calculation for current set" << endl;
306 for (Int_t i1=0; i1<3;i1++) {
307 beta[0]=momentum[0]/sqrt(massarray[i1]*massarray[i1]+momentum[0]*momentum[0]);
308 for (Int_t i2=0; i2<3;i2++) {
309 beta[1]=momentum[1]/sqrt(massarray[i2]*massarray[i2]+momentum[1]*momentum[1]);
310 for (Int_t i3=0; i3<3;i3++) {
311 beta[2]=momentum[2]/sqrt(massarray[i3]*massarray[i3]+momentum[2]*momentum[2]);
312 for (Int_t i4=0; i4<3;i4++) {
313 beta[3]=momentum[3]/sqrt(massarray[i4]*massarray[i4]+momentum[3]*momentum[3]);
314 for (Int_t i5=0; i5<3;i5++) {
315 beta[4]=momentum[4]/sqrt(massarray[i5]*massarray[i5]+momentum[4]*momentum[4]);
316 for (Int_t i6=0; i6<3;i6++) {
317 beta[5]=momentum[5]/sqrt(massarray[i6]*massarray[i6]+momentum[5]*momentum[5]);
318 for (Int_t i7=0; i7<3;i7++) {
319 beta[6]=momentum[6]/sqrt(massarray[i7]*massarray[i7]+momentum[6]*momentum[6]);
320 for (Int_t i8=0; i8<3;i8++) {
321 beta[7]=momentum[7]/sqrt(massarray[i8]*massarray[i8]+momentum[7]*momentum[7]);
322 for (Int_t i9=0; i9<3;i9++) {
323 beta[8]=momentum[8]/sqrt(massarray[i9]*massarray[i9]+momentum[8]*momentum[8]);
324 for (Int_t i10=0; i10<3;i10++) {
325 beta[9]=momentum[9]/sqrt(massarray[i10]*massarray[i10]+momentum[9]*momentum[9]);
327 Float_t meantzero=0.;
328 Float_t sumAllweights=0.;
329 for (Int_t itz=0; itz<10;itz++) {
330 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
331 sqTrackError[itz]=(timeresolutioninns*timeresolutioninns+sqMomError[itz]); // total error for the current track
332 sumAllweights+=1./sqTrackError[itz];
334 timezero[itz]=(tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz];
335 weightedtimezero[itz]=((tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz])/sqTrackError[itz];// weighted time zero for current track
336 meantzero+=weightedtimezero[itz];
337 } // end loop for (Int_t itz=0; itz<10;itz++)
338 meantzero=meantzero/sumAllweights; // it is given in [ns]
341 // calculate the chisquare for the current assignment
342 for (Int_t icsq=0; icsq<10;icsq++) {
343 dummychisquare+=(timezero[icsq]-meantzero)*(timezero[icsq]-meantzero)/sqTrackError[icsq];
344 } // end loop for (Int_t icsq=0; icsq<10;icsq++)
346 if(dummychisquare<=chisquare){
357 chisquare=dummychisquare;
359 } // close if(dummychisquare<=chisquare)
372 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;
373 if(truparticle[0]!=assparticle[0]) nmisidentified0+=1;
374 if(truparticle[1]!=assparticle[1]) nmisidentified1+=1;
375 if(truparticle[2]!=assparticle[2]) nmisidentified2+=1;
376 if(truparticle[3]!=assparticle[3]) nmisidentified3+=1;
377 if(truparticle[4]!=assparticle[4]) nmisidentified4+=1;
378 if(truparticle[5]!=assparticle[5]) nmisidentified5+=1;
379 if(truparticle[6]!=assparticle[6]) nmisidentified6+=1;
380 if(truparticle[7]!=assparticle[7]) nmisidentified7+=1;
381 if(truparticle[8]!=assparticle[8]) nmisidentified8+=1;
382 if(truparticle[9]!=assparticle[9]) nmisidentified9+=1;
384 htzerobest->Fill(t0best);
385 hchibest->Fill(chisquare);
386 Double_t dblechisquare=(Double_t)chisquare;
387 Float_t confLevel=(Float_t)TMath::Prob(dblechisquare,9); // ndf 10-1=9
388 hchibestconflevel->Fill(confLevel);
390 if(strstr(option,"all")){
391 cout << "True Assignment " << truparticle[0] << truparticle[1] << truparticle[2] << truparticle[3] << truparticle[4] << truparticle[5] << truparticle[6] << truparticle[7] << truparticle[8] << truparticle[9] <<endl;
392 cout << "Best Assignment " << assparticle[0] << assparticle[1] << assparticle[2] << assparticle[3] << assparticle[4] << assparticle[5] << assparticle[6] << assparticle[7] << assparticle[8] << assparticle[9] << endl;
393 cout << "Minimum ChiSquare for current set " << chisquare << endl;
394 cout << "Confidence Level (Minimum ChiSquare) " << confLevel << endl;
396 if (strstr(option,"visual") && itimes && (itimes%kUPDATE) == 0) {
397 if (itimes == kUPDATE){
403 hchibestconflevel->Draw();
411 if (gSystem->ProcessEvents())
417 } // end for the current set. close if(istop==5)
418 } // end condition on ipartold
419 } // end loop on hits for the current track
421 } // end loop on ntracks
424 if(strstr(option,"all")){
425 nmisidentified=(nmisidentified0+nmisidentified1+nmisidentified2+nmisidentified3+nmisidentified4+nmisidentified5+nmisidentified6+nmisidentified7+nmisidentified8+nmisidentified9);
426 cout << "total number of tracks token into account " << 10*5*fNevents << endl;
427 Float_t badPercentage=100.*(Float_t)nmisidentified/(10*5*fNevents);
428 cout << "total misidentified " << nmisidentified << "("<< badPercentage << "%)" <<endl;
429 cout << "Total Number of set token into account " << 5*fNevents << endl;
430 Float_t goodSetPercentage=100.*(Float_t)ngood/(5*fNevents);
431 cout << "Number of set with no misidentified tracks " << ngood << "("<< goodSetPercentage << "%)" <<endl;
434 // free used memory for canvas
439 // generating output filename only if not previously specified using SetTZeroFile
440 char outFileName[70];
441 strcpy(outFileName,"ht010tr120ps"); // global time resolution has to be converted from Int_t to char
442 // in order to have in the output filename this parameter
443 strcat(outFileName,fHeadersFile);
445 if(fT0File.IsNull()) fT0File=outFileName;
447 TFile* houtfile = new TFile(fT0File,"recreate");
449 htzerobest->Write(0,TObject::kOverwrite);
450 hchibest->Write(0,TObject::kOverwrite);
451 hchibestconflevel->Write(0,TObject::kOverwrite);
455 if(strstr(option,"tim") || strstr(option,"all")){
456 gBenchmark->Stop("TOFT0");
457 cout << "AliTOFT0:" << endl ;
458 cout << " took " << gBenchmark->GetCpuTime("TOFT0") << " seconds in order to calculate T0 "
459 << gBenchmark->GetCpuTime("TOFT0")/fNevents << " seconds per event " << endl ;
464 //__________________________________________________________________
465 void AliTOFT0::SetTZeroFile(char * file ){
466 cout << "Destination file : " << file << endl ;
469 //__________________________________________________________________
470 void AliTOFT0::Print(Option_t* option)const
472 cout << "------------------- "<< GetName() << " -------------" << endl ;
473 if(!fT0File.IsNull())
474 cout << " Writing T0 Distribution to file " << (char*) fT0File.Data() << endl ;
477 //__________________________________________________________________
478 Bool_t AliTOFT0::operator==( AliTOFT0 const &tzero )const
483 if( (fTimeResolution==tzero.fTimeResolution)&&(fLowerMomBound==tzero.fLowerMomBound)&&(fUpperMomBound==tzero.fUpperMomBound))