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 **************************************************************************/
18 //_________________________________________________________________________
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
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.
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
78 //-- Author: F. Pierella
80 //_________________________________________________________________________
83 #include "TClonesArray.h"
88 #include "TParticle.h"
89 #include "TBenchmark.h"
96 #include "AliTOFhitT0.h"
100 extern TBenchmark *gBenchmark;
101 extern TSystem *gSystem;
102 extern TRandom *gRandom;
105 extern AliRun *gAlice;
110 //____________________________________________________________________________
111 AliTOFT0::AliTOFT0():TTask("AliTOFT0","")
117 //____________________________________________________________________________
118 AliTOFT0::AliTOFT0(char* headerFile, Int_t nEvents):TTask("AliTOFT0","")
124 fNevents=nEvents ; // Number of events for which calculate the T0,
125 // default 0: it means all evens in current file
126 fLowerMomBound=1.5; // [GeV/c] default value
127 fUpperMomBound=2. ; // [GeV/c] default value
128 fTimeResolution = 1.2e-10; // 120 ps by default
129 fHeadersFile = headerFile ;
131 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
133 //File was not opened yet
135 if(fHeadersFile.Contains("rfio"))
136 file = TFile::Open(fHeadersFile,"update") ;
138 file = new TFile(fHeadersFile.Data(),"update") ;
139 gAlice = (AliRun *) file->Get("gAlice") ;
142 // add Task to //root/Tasks folder
143 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
144 roottasks->Add(this) ;
147 //____________________________________________________________________________
148 AliTOFT0::AliTOFT0(const AliTOFT0 & tzero):TTask("AliTOFT0","")
152 ( (AliTOFT0 &)tzero ).Copy(*this);
155 //____________________________________________________________________________
156 AliTOFT0::~AliTOFT0()
161 //____________________________________________________________________________
162 void AliTOFT0::Exec(Option_t *option)
165 // calculate T0 distribution for all events using chisquare
168 Int_t nmisidentified=0;
169 Int_t nmisidentified0=0;
170 Int_t nmisidentified1=0;
171 Int_t nmisidentified2=0;
172 Int_t nmisidentified3=0;
173 Int_t nmisidentified4=0;
174 Int_t nmisidentified5=0;
175 Int_t nmisidentified6=0;
176 Int_t nmisidentified7=0;
177 Int_t nmisidentified8=0;
178 Int_t nmisidentified9=0;
183 Float_t timeresolutioninns=fTimeResolution*(1.e+9); // convert in [ns]
184 const Int_t kUPDATE = 5; // for visual option
190 if(strstr(option,"visual")){
191 // Create a new canvas.
192 //c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,500,500);
193 c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,370,370);
194 c1->SetFillColor(35);
195 c1->GetFrame()->SetFillColor(21);
196 c1->GetFrame()->SetBorderSize(6);
197 c1->GetFrame()->SetBorderMode(-1);
199 //c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",550,10,500,500);
200 c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",380,10,370,370);
201 c2->SetFillColor(35);
202 c2->GetFrame()->SetFillColor(21);
203 c2->GetFrame()->SetBorderSize(6);
204 c2->GetFrame()->SetBorderMode(-1);
206 //c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",280,550,500,500);
207 c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",760,10,370,370);
208 c3->SetFillColor(35);
209 c3->GetFrame()->SetFillColor(21);
210 c3->GetFrame()->SetBorderSize(6);
211 c3->GetFrame()->SetBorderMode(-1);
214 if(strstr(option,"tim") || strstr(option,"all"))
215 gBenchmark->Start("TOFT0");
217 TH1F *htzerobest= new TH1F("htzerobest","T0 for best assignment",200,-1.,1.);
218 TH1F* hchibest = new TH1F("hchibest","ChiSquare Min Distribution",80,0.,40.);
219 TH1F* hchibestconflevel = new TH1F("hchibestconflevel","ChiSquare Min Confidence Level",10,0.,1.);
221 // setting histo colors
222 if(strstr(option,"visual")){
223 htzerobest->SetFillColor(48);
224 hchibest->SetFillColor(50);
225 hchibestconflevel->SetFillColor(52);
228 Int_t assparticle[10]={3,3,3,3,3,3,3,3,3,3};
229 Int_t truparticle[10]={3,3,3,3,3,3,3,3,3,3};
231 Float_t timeofflight[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
232 Float_t momentum[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
233 Float_t timezero[10];
234 Float_t weightedtimezero[10];
235 Float_t beta[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
236 Float_t sqMomError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
237 Float_t sqTrackError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
238 Float_t massarray[3]={0.13957,0.493677,0.9382723};
239 Float_t dummychisquare=0.;
240 Float_t chisquare=999.;
241 Float_t tracktoflen[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
243 AliTOF *detTOF = (AliTOF *) gAlice->GetDetector ("TOF");
246 Error("AliTOFT0","TOF not found");
250 if(strstr(option,"all")){
251 AliInfo(Form("Selecting primary tracks with momentum between %d GeV/c and %d GeV/c", fLowerMomBound, fUpperMomBound));
252 AliInfo("Memorandum: 0 means PION | 1 means KAON | 2 means PROTON")
255 if (fNevents == 0) fNevents = (Int_t) gAlice->TreeE()->GetEntries();
257 for (Int_t ievent = 0; ievent < fNevents; ievent++) {
258 gAlice->GetEvent(ievent);
259 TTree *hitTree = detTOF->TreeH ();
264 TClonesArray* tofHits = detTOF->Hits();
269 hitTree->SetBranchStatus("*",0); // switch off all branches
270 hitTree->SetBranchStatus("TOF*",1); // switch on only TOF
272 // Start loop on primary tracks in the hits containers
274 Int_t ntracks = static_cast<Int_t>(hitTree->GetEntries());
275 for (Int_t track = 0; track < ntracks; track++)
277 if(nset>=5) break; // check on the number of set analyzed
280 hitTree->GetEvent(track);
282 AliMC *mcApplication = (AliMC*)gAlice->GetMCApp();
284 particle = mcApplication->Particle(track);
285 Int_t nhits = tofHits->GetEntriesFast();
287 for (Int_t hit = 0; hit < nhits; hit++)
289 tofHit = (AliTOFhitT0 *) tofHits->UncheckedAt(hit);
290 ipart = tofHit->GetTrack();
291 // check to discard the case when the same particle is selected more than one
294 if (ipart != ipartold){
296 particle = (TParticle*)gAlice->GetMCApp()->Particle(ipart);
298 Float_t idealtime=tofHit->GetTof();
299 // Float_t time=idealtime;
300 Float_t time = gRandom->Gaus(idealtime, fTimeResolution);
301 Float_t toflen=tofHit->GetLen();
302 toflen=toflen/100.; // toflen given in m
303 Int_t pdg = particle->GetPdgCode();
304 Int_t abspdg =TMath::Abs(pdg);
305 Float_t idealmom = particle->P();
306 Float_t momres=idealmom*0.025; // 2.5% res token into account for all momenta
307 Float_t mom =gRandom->Gaus(idealmom,momres);
309 Bool_t isgoodpart=(abspdg==211 || abspdg==2212 || abspdg==321);
311 time*=1.E+9; // tof given in nanoseconds
312 if (particle->GetFirstMother() < 0 && isgoodpart && mom<=fUpperMomBound && mom>=fLowerMomBound){
316 Int_t index=selected-1;
317 timeofflight[index]=time;
318 tracktoflen[index]=toflen;
320 //AliInfo(Form(" %d %d %d ", timeofflight[index], tracktoflen[index], momentum[index]));
323 truparticle[index]=0;
326 truparticle[index]=1;
329 truparticle[index]=2;
336 if(istop==10){ // start analysis on current set
341 //AliInfo("starting t0 calculation for current set");
342 for (Int_t i1=0; i1<3;i1++) {
343 beta[0]=momentum[0]/sqrt(massarray[i1]*massarray[i1]+momentum[0]*momentum[0]);
344 for (Int_t i2=0; i2<3;i2++) {
345 beta[1]=momentum[1]/sqrt(massarray[i2]*massarray[i2]+momentum[1]*momentum[1]);
346 for (Int_t i3=0; i3<3;i3++) {
347 beta[2]=momentum[2]/sqrt(massarray[i3]*massarray[i3]+momentum[2]*momentum[2]);
348 for (Int_t i4=0; i4<3;i4++) {
349 beta[3]=momentum[3]/sqrt(massarray[i4]*massarray[i4]+momentum[3]*momentum[3]);
350 for (Int_t i5=0; i5<3;i5++) {
351 beta[4]=momentum[4]/sqrt(massarray[i5]*massarray[i5]+momentum[4]*momentum[4]);
352 for (Int_t i6=0; i6<3;i6++) {
353 beta[5]=momentum[5]/sqrt(massarray[i6]*massarray[i6]+momentum[5]*momentum[5]);
354 for (Int_t i7=0; i7<3;i7++) {
355 beta[6]=momentum[6]/sqrt(massarray[i7]*massarray[i7]+momentum[6]*momentum[6]);
356 for (Int_t i8=0; i8<3;i8++) {
357 beta[7]=momentum[7]/sqrt(massarray[i8]*massarray[i8]+momentum[7]*momentum[7]);
358 for (Int_t i9=0; i9<3;i9++) {
359 beta[8]=momentum[8]/sqrt(massarray[i9]*massarray[i9]+momentum[8]*momentum[8]);
360 for (Int_t i10=0; i10<3;i10++) {
361 beta[9]=momentum[9]/sqrt(massarray[i10]*massarray[i10]+momentum[9]*momentum[9]);
363 Float_t meantzero=0.;
364 Float_t sumAllweights=0.;
365 for (Int_t itz=0; itz<10;itz++) {
366 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
367 sqTrackError[itz]=(timeresolutioninns*timeresolutioninns+sqMomError[itz]); // total error for the current track
368 sumAllweights+=1./sqTrackError[itz];
370 timezero[itz]=(tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz];
371 weightedtimezero[itz]=((tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz])/sqTrackError[itz];// weighted time zero for current track
372 meantzero+=weightedtimezero[itz];
373 } // end loop for (Int_t itz=0; itz<10;itz++)
374 meantzero=meantzero/sumAllweights; // it is given in [ns]
377 // calculate the chisquare for the current assignment
378 for (Int_t icsq=0; icsq<10;icsq++) {
379 dummychisquare+=(timezero[icsq]-meantzero)*(timezero[icsq]-meantzero)/sqTrackError[icsq];
380 } // end loop for (Int_t icsq=0; icsq<10;icsq++)
382 if(dummychisquare<=chisquare){
393 chisquare=dummychisquare;
395 } // close if(dummychisquare<=chisquare)
408 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;
409 if(truparticle[0]!=assparticle[0]) nmisidentified0+=1;
410 if(truparticle[1]!=assparticle[1]) nmisidentified1+=1;
411 if(truparticle[2]!=assparticle[2]) nmisidentified2+=1;
412 if(truparticle[3]!=assparticle[3]) nmisidentified3+=1;
413 if(truparticle[4]!=assparticle[4]) nmisidentified4+=1;
414 if(truparticle[5]!=assparticle[5]) nmisidentified5+=1;
415 if(truparticle[6]!=assparticle[6]) nmisidentified6+=1;
416 if(truparticle[7]!=assparticle[7]) nmisidentified7+=1;
417 if(truparticle[8]!=assparticle[8]) nmisidentified8+=1;
418 if(truparticle[9]!=assparticle[9]) nmisidentified9+=1;
420 htzerobest->Fill(t0best);
421 hchibest->Fill(chisquare);
422 Double_t dblechisquare=(Double_t)chisquare;
423 Float_t confLevel=(Float_t)TMath::Prob(dblechisquare,9); // ndf 10-1=9
424 hchibestconflevel->Fill(confLevel);
426 if(strstr(option,"all")){
427 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]));
428 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]));
429 AliInfo(Form("Minimum ChiSquare for current set %d ", chisquare));
430 AliInfo(Form("Confidence Level (Minimum ChiSquare) %d", confLevel));
432 if (strstr(option,"visual") && itimes && (itimes%kUPDATE) == 0) {
433 if (itimes == kUPDATE){
439 hchibestconflevel->Draw();
447 if (gSystem->ProcessEvents())
453 } // end for the current set. close if(istop==5)
454 } // end condition on ipartold
455 } // end loop on hits for the current track
457 } // end loop on ntracks
460 if(strstr(option,"all")){
461 nmisidentified=(nmisidentified0+nmisidentified1+nmisidentified2+nmisidentified3+nmisidentified4+nmisidentified5+nmisidentified6+nmisidentified7+nmisidentified8+nmisidentified9);
462 AliInfo(Form("total number of tracks token into account %i", 10*5*fNevents));
463 Float_t badPercentage=100.*(Float_t)nmisidentified/(10*5*fNevents);
464 AliInfo(Form("total misidentified %i (%d %) ", nmisidentified, badPercentage));
465 AliInfo(Form("Total Number of set token into account %i", 5*fNevents));
466 Float_t goodSetPercentage=100.*(Float_t)ngood/(5*fNevents);
467 AliInfo(Form("Number of set with no misidentified tracks %i (%d %)", ngood, goodSetPercentage));
470 // free used memory for canvas
475 // generating output filename only if not previously specified using SetTZeroFile
476 char outFileName[70];
477 strcpy(outFileName,"ht010tr120ps"); // global time resolution has to be converted from Int_t to char
478 // in order to have in the output filename this parameter
479 strcat(outFileName,fHeadersFile);
481 if(fT0File.IsNull()) fT0File=outFileName;
483 TFile* houtfile = new TFile(fT0File,"recreate");
485 htzerobest->Write(0,TObject::kOverwrite);
486 hchibest->Write(0,TObject::kOverwrite);
487 hchibestconflevel->Write(0,TObject::kOverwrite);
491 if(strstr(option,"tim") || strstr(option,"all")){
492 gBenchmark->Stop("TOFT0");
493 AliInfo("AliTOFT0:");
495 cout << " took " << gBenchmark->GetCpuTime("TOFT0") << " seconds in order to calculate T0 "
496 << gBenchmark->GetCpuTime("TOFT0")/fNevents << " seconds per event " << endl ;
498 gBenchmark->Print("TOFT0");
502 //__________________________________________________________________
503 void AliTOFT0::SetTZeroFile(char * file )
508 printf("Destination file : %s \n", file) ;
513 //__________________________________________________________________
514 void AliTOFT0::Print(Option_t* /*option*/)const
519 printf("------------------- %s -------------\n", GetName()) ;
520 if(!fT0File.IsNull())
521 printf(" Writing T0 Distribution to file %s \n",(char*) fT0File.Data());
525 //__________________________________________________________________
526 Bool_t AliTOFT0::operator==( AliTOFT0 const &tzero )const
532 if( (fTimeResolution==tzero.fTimeResolution)&&(fLowerMomBound==tzero.fLowerMomBound)&&(fUpperMomBound==tzero.fUpperMomBound))