]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFT0.cxx
Suppression of some warnings
[u/mrichter/AliRoot.git] / TOF / AliTOFT0.cxx
CommitLineData
d599d913 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
88cb7938 16/* $Id$ */
17
d599d913 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
88cb7938 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 "AliRun.h"
85#include "AliTOF.h"
d599d913 86#include "AliTOFT0.h"
87#include "AliTOFhitT0.h"
d599d913 88#include "AliTOFv4T0.h"
d599d913 89
90ClassImp(AliTOFT0)
91
92//____________________________________________________________________________
93 AliTOFT0::AliTOFT0():TTask("AliTOFT0","")
94{
95 // ctor
ea7a588a 96 fNevents = 0 ;
d599d913 97}
98
99//____________________________________________________________________________
100 AliTOFT0::AliTOFT0(char* headerFile, Int_t nEvents):TTask("AliTOFT0","")
101{
102 fNevents=nEvents ; // Number of events for which calculate the T0,
103 // default 0: it means all evens in current file
104 fLowerMomBound=1.5; // [GeV/c] default value
105 fUpperMomBound=2. ; // [GeV/c] default value
106 fTimeResolution = 1.2e-10; // 120 ps by default
107 fHeadersFile = headerFile ;
d599d913 108
109 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
110
111 //File was not opened yet
112 if(file == 0){
113 if(fHeadersFile.Contains("rfio"))
114 file = TFile::Open(fHeadersFile,"update") ;
115 else
116 file = new TFile(fHeadersFile.Data(),"update") ;
117 gAlice = (AliRun *) file->Get("gAlice") ;
118 }
119
120 // add Task to //root/Tasks folder
121 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
122 roottasks->Add(this) ;
123}
124
5c016a7b 125//____________________________________________________________________________
126 AliTOFT0::AliTOFT0(const AliTOFT0 & tzero):TTask("AliTOFT0","")
127{
128( (AliTOFT0 &)tzero ).Copy(*this);
129}
130
d599d913 131//____________________________________________________________________________
132 AliTOFT0::~AliTOFT0()
133{
134 // dtor
135}
136
137//____________________________________________________________________________
138void AliTOFT0::Exec(Option_t *option)
139{
140 //
141 // calculate T0 distribution for all events using chisquare
142 //
143 Int_t ngood=0;
144 Int_t nmisidentified=0;
145 Int_t nmisidentified0=0;
146 Int_t nmisidentified1=0;
147 Int_t nmisidentified2=0;
148 Int_t nmisidentified3=0;
149 Int_t nmisidentified4=0;
150 Int_t nmisidentified5=0;
151 Int_t nmisidentified6=0;
152 Int_t nmisidentified7=0;
153 Int_t nmisidentified8=0;
154 Int_t nmisidentified9=0;
155 Int_t ipartold = -1;
156 Int_t ipart;
157 Int_t selected=0;
158 Int_t istop=0;
159 Float_t timeresolutioninns=fTimeResolution*(1.e+9); // convert in [ns]
160 const Int_t kUPDATE = 5; // for visual option
161 Int_t itimes=0;
162 TCanvas* c1=0;
163 TCanvas* c2=0;
164 TCanvas* c3=0;
165
166 if(strstr(option,"visual")){
167 // Create a new canvas.
168 //c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,500,500);
169 c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,370,370);
170 c1->SetFillColor(35);
171 c1->GetFrame()->SetFillColor(21);
172 c1->GetFrame()->SetBorderSize(6);
173 c1->GetFrame()->SetBorderMode(-1);
174
175 //c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",550,10,500,500);
176 c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",380,10,370,370);
177 c2->SetFillColor(35);
178 c2->GetFrame()->SetFillColor(21);
179 c2->GetFrame()->SetBorderSize(6);
180 c2->GetFrame()->SetBorderMode(-1);
181
182 //c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",280,550,500,500);
183 c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",760,10,370,370);
184 c3->SetFillColor(35);
185 c3->GetFrame()->SetFillColor(21);
186 c3->GetFrame()->SetBorderSize(6);
187 c3->GetFrame()->SetBorderMode(-1);
188 }
189
190 if(strstr(option,"tim") || strstr(option,"all"))
191 gBenchmark->Start("TOFT0");
192
193 TH1F *htzerobest= new TH1F("htzerobest","T0 for best assignment",200,-1.,1.);
194 TH1F* hchibest = new TH1F("hchibest","ChiSquare Min Distribution",80,0.,40.);
195 TH1F* hchibestconflevel = new TH1F("hchibestconflevel","ChiSquare Min Confidence Level",10,0.,1.);
196
197 // setting histo colors
198 if(strstr(option,"visual")){
199 htzerobest->SetFillColor(48);
200 hchibest->SetFillColor(50);
201 hchibestconflevel->SetFillColor(52);
202 }
203
204 Int_t assparticle[10]={3,3,3,3,3,3,3,3,3,3};
205 Int_t truparticle[10]={3,3,3,3,3,3,3,3,3,3};
206 Float_t t0best=999.;
207 Float_t timeofflight[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
208 Float_t momentum[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
209 Float_t timezero[10];
210 Float_t weightedtimezero[10];
211 Float_t beta[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
212 Float_t sqMomError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
213 Float_t sqTrackError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
214 Float_t massarray[3]={0.13957,0.493677,0.9382723};
215 Float_t dummychisquare=0.;
216 Float_t chisquare=999.;
217 Float_t tracktoflen[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
218
219 AliTOF *TOF = (AliTOF *) gAlice->GetDetector ("TOF");
220
221 if (!TOF) {
222 Error("AliTOFT0","TOF not found");
223 return;
224 }
225
226 if(strstr(option,"all")){
227 cout << "Selecting primary tracks with momentum between " << fLowerMomBound << " GeV/c and " << fUpperMomBound << " GeV/c" << endl;
228 cout << "Memorandum: 0 means PION | 1 means KAON | 2 means PROTON" << endl;
229 }
230
231 if (fNevents == 0) fNevents = (Int_t) gAlice->TreeE()->GetEntries();
232
233 for (Int_t ievent = 0; ievent < fNevents; ievent++) {
234 gAlice->GetEvent(ievent);
88cb7938 235 TTree *TH = TOF->TreeH ();
d599d913 236 if (!TH)
237 return;
238 TParticle* particle;
239 AliTOFhitT0* tofHit;
240 TClonesArray* TOFhits = TOF->Hits();
241
242 Int_t lasttrack=-1;
243 Int_t nset=0;
5fff655e 244
245 TH->SetBranchStatus("*",0); // switch off all branches
246 TH->SetBranchStatus("TOF*",1); // switch on only TOF
247
d599d913 248 // Start loop on primary tracks in the hits containers
249
250 Int_t ntracks = static_cast<Int_t>(TH->GetEntries());
251 for (Int_t track = 0; track < ntracks; track++)
252 {
253 if(nset>=5) break; // check on the number of set analyzed
254
255 gAlice->ResetHits();
256 TH->GetEvent(track);
257 particle = gAlice->Particle(track);
258 Int_t nhits = TOFhits->GetEntriesFast();
259
260 for (Int_t hit = 0; hit < nhits; hit++)
261 {
262 tofHit = (AliTOFhitT0 *) TOFhits->UncheckedAt(hit);
263 ipart = tofHit->GetTrack();
264 // check to discard the case when the same particle is selected more than one
265 // time
266
267 if (ipart != ipartold){
268
269 particle = (TParticle*)gAlice->Particle(ipart);
270
271 Float_t idealtime=tofHit->GetTof();
272 // Float_t time=idealtime;
273 Float_t time = gRandom->Gaus(idealtime, fTimeResolution);
274 Float_t toflen=tofHit->GetLen();
275 toflen=toflen/100.; // toflen given in m
276 Int_t pdg = particle->GetPdgCode();
277 Int_t abspdg =TMath::Abs(pdg);
278 Float_t idealmom = particle->P();
279 Float_t momres=idealmom*0.025; // 2.5% res token into account for all momenta
280 Float_t mom =gRandom->Gaus(idealmom,momres);
281
282 Bool_t isgoodpart=(abspdg==211 || abspdg==2212 || abspdg==321);
283
284 time*=1.E+9; // tof given in nanoseconds
285 if (particle->GetFirstMother() < 0 && isgoodpart && mom<=fUpperMomBound && mom>=fLowerMomBound){
286 selected+=1;
287 istop=selected;
288 if(istop>10) break;
289 Int_t index=selected-1;
290 timeofflight[index]=time;
291 tracktoflen[index]=toflen;
292 momentum[index]=mom;
293 // cout << timeofflight[index] << " " << tracktoflen[index] << " " << momentum[index] << endl;
294 switch (abspdg) {
295 case 211:
296 truparticle[index]=0;
297 break ;
298 case 321:
299 truparticle[index]=1;
300 break ;
301 case 2212:
302 truparticle[index]=2;
303 break ;
304 }
305
306 }
307 ipartold = ipart;
308
309 if(istop==10){ // start analysis on current set
310 nset+=1;
311 lasttrack=track;
312 istop=0;
313 selected=0;
314 //cout << "starting t0 calculation for current set" << endl;
315 for (Int_t i1=0; i1<3;i1++) {
ea7a588a 316 beta[0]=momentum[0]/sqrt(massarray[i1]*massarray[i1]+momentum[0]*momentum[0]);
d599d913 317 for (Int_t i2=0; i2<3;i2++) {
ea7a588a 318 beta[1]=momentum[1]/sqrt(massarray[i2]*massarray[i2]+momentum[1]*momentum[1]);
d599d913 319 for (Int_t i3=0; i3<3;i3++) {
ea7a588a 320 beta[2]=momentum[2]/sqrt(massarray[i3]*massarray[i3]+momentum[2]*momentum[2]);
d599d913 321 for (Int_t i4=0; i4<3;i4++) {
ea7a588a 322 beta[3]=momentum[3]/sqrt(massarray[i4]*massarray[i4]+momentum[3]*momentum[3]);
d599d913 323 for (Int_t i5=0; i5<3;i5++) {
ea7a588a 324 beta[4]=momentum[4]/sqrt(massarray[i5]*massarray[i5]+momentum[4]*momentum[4]);
d599d913 325 for (Int_t i6=0; i6<3;i6++) {
ea7a588a 326 beta[5]=momentum[5]/sqrt(massarray[i6]*massarray[i6]+momentum[5]*momentum[5]);
d599d913 327 for (Int_t i7=0; i7<3;i7++) {
ea7a588a 328 beta[6]=momentum[6]/sqrt(massarray[i7]*massarray[i7]+momentum[6]*momentum[6]);
d599d913 329 for (Int_t i8=0; i8<3;i8++) {
ea7a588a 330 beta[7]=momentum[7]/sqrt(massarray[i8]*massarray[i8]+momentum[7]*momentum[7]);
d599d913 331 for (Int_t i9=0; i9<3;i9++) {
ea7a588a 332 beta[8]=momentum[8]/sqrt(massarray[i9]*massarray[i9]+momentum[8]*momentum[8]);
d599d913 333 for (Int_t i10=0; i10<3;i10++) {
d599d913 334 beta[9]=momentum[9]/sqrt(massarray[i10]*massarray[i10]+momentum[9]*momentum[9]);
335
d599d913 336 Float_t meantzero=0.;
337 Float_t sumAllweights=0.;
338 for (Int_t itz=0; itz<10;itz++) {
339 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
340 sqTrackError[itz]=(timeresolutioninns*timeresolutioninns+sqMomError[itz]); // total error for the current track
341 sumAllweights+=1./sqTrackError[itz];
ea7a588a 342
343 timezero[itz]=(tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz];
344 weightedtimezero[itz]=((tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz])/sqTrackError[itz];// weighted time zero for current track
d599d913 345 meantzero+=weightedtimezero[itz];
346 } // end loop for (Int_t itz=0; itz<10;itz++)
347 meantzero=meantzero/sumAllweights; // it is given in [ns]
348
349 dummychisquare=0.;
350 // calculate the chisquare for the current assignment
351 for (Int_t icsq=0; icsq<10;icsq++) {
352 dummychisquare+=(timezero[icsq]-meantzero)*(timezero[icsq]-meantzero)/sqTrackError[icsq];
353 } // end loop for (Int_t icsq=0; icsq<10;icsq++)
354
355 if(dummychisquare<=chisquare){
356 assparticle[0]=i1;
357 assparticle[1]=i2;
358 assparticle[2]=i3;
359 assparticle[3]=i4;
360 assparticle[4]=i5;
361 assparticle[5]=i6;
362 assparticle[6]=i7;
363 assparticle[7]=i8;
364 assparticle[8]=i9;
365 assparticle[9]=i10;
366 chisquare=dummychisquare;
367 t0best=meantzero;
368 } // close if(dummychisquare<=chisquare)
369
370 } // end loop on i10
371 } // end loop on i9
372 } // end loop on i8
373 } // end loop on i7
374 } // end loop on i6
375 } // end loop on i5
376 } // end loop on i4
377 } // end loop on i3
378 } // end loop on i2
379 } // end loop on i1
380
381 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;
382 if(truparticle[0]!=assparticle[0]) nmisidentified0+=1;
383 if(truparticle[1]!=assparticle[1]) nmisidentified1+=1;
384 if(truparticle[2]!=assparticle[2]) nmisidentified2+=1;
385 if(truparticle[3]!=assparticle[3]) nmisidentified3+=1;
386 if(truparticle[4]!=assparticle[4]) nmisidentified4+=1;
387 if(truparticle[5]!=assparticle[5]) nmisidentified5+=1;
388 if(truparticle[6]!=assparticle[6]) nmisidentified6+=1;
389 if(truparticle[7]!=assparticle[7]) nmisidentified7+=1;
390 if(truparticle[8]!=assparticle[8]) nmisidentified8+=1;
391 if(truparticle[9]!=assparticle[9]) nmisidentified9+=1;
392 // filling histos
393 htzerobest->Fill(t0best);
394 hchibest->Fill(chisquare);
395 Double_t dblechisquare=(Double_t)chisquare;
396 Float_t confLevel=(Float_t)TMath::Prob(dblechisquare,9); // ndf 10-1=9
397 hchibestconflevel->Fill(confLevel);
398 itimes++;
399 if(strstr(option,"all")){
400 cout << "True Assignment " << truparticle[0] << truparticle[1] << truparticle[2] << truparticle[3] << truparticle[4] << truparticle[5] << truparticle[6] << truparticle[7] << truparticle[8] << truparticle[9] <<endl;
401 cout << "Best Assignment " << assparticle[0] << assparticle[1] << assparticle[2] << assparticle[3] << assparticle[4] << assparticle[5] << assparticle[6] << assparticle[7] << assparticle[8] << assparticle[9] << endl;
402 cout << "Minimum ChiSquare for current set " << chisquare << endl;
403 cout << "Confidence Level (Minimum ChiSquare) " << confLevel << endl;
404 }
405 if (strstr(option,"visual") && itimes && (itimes%kUPDATE) == 0) {
406 if (itimes == kUPDATE){
407 c1->cd();
408 htzerobest->Draw();
409 c2->cd();
410 hchibest->Draw();
411 c3->cd();
412 hchibestconflevel->Draw();
413 }
414 c1->Modified();
415 c1->Update();
416 c2->Modified();
417 c2->Update();
418 c3->Modified();
419 c3->Update();
420 if (gSystem->ProcessEvents())
421 break;
422 }
423 chisquare=999.;
424 t0best=999.;
425
426 } // end for the current set. close if(istop==5)
427 } // end condition on ipartold
428 } // end loop on hits for the current track
429 if(istop>=10) break;
430 } // end loop on ntracks
431 } //event loop
432
433 if(strstr(option,"all")){
434 nmisidentified=(nmisidentified0+nmisidentified1+nmisidentified2+nmisidentified3+nmisidentified4+nmisidentified5+nmisidentified6+nmisidentified7+nmisidentified8+nmisidentified9);
435 cout << "total number of tracks token into account " << 10*5*fNevents << endl;
436 Float_t badPercentage=100.*(Float_t)nmisidentified/(10*5*fNevents);
437 cout << "total misidentified " << nmisidentified << "("<< badPercentage << "%)" <<endl;
438 cout << "Total Number of set token into account " << 5*fNevents << endl;
439 Float_t goodSetPercentage=100.*(Float_t)ngood/(5*fNevents);
440 cout << "Number of set with no misidentified tracks " << ngood << "("<< goodSetPercentage << "%)" <<endl;
441 }
442
443 // free used memory for canvas
444 delete c1; c1=0;
445 delete c2; c2=0;
446 delete c3; c3=0;
447
448 // generating output filename only if not previously specified using SetTZeroFile
449 char outFileName[70];
450 strcpy(outFileName,"ht010tr120ps"); // global time resolution has to be converted from Int_t to char
451 // in order to have in the output filename this parameter
452 strcat(outFileName,fHeadersFile);
453
454 if(fT0File.IsNull()) fT0File=outFileName;
455
456 TFile* houtfile = new TFile(fT0File,"recreate");
457 houtfile->cd();
458 htzerobest->Write(0,TObject::kOverwrite);
459 hchibest->Write(0,TObject::kOverwrite);
460 hchibestconflevel->Write(0,TObject::kOverwrite);
461 houtfile->Close();
462
463
464 if(strstr(option,"tim") || strstr(option,"all")){
465 gBenchmark->Stop("TOFT0");
466 cout << "AliTOFT0:" << endl ;
467 cout << " took " << gBenchmark->GetCpuTime("TOFT0") << " seconds in order to calculate T0 "
468 << gBenchmark->GetCpuTime("TOFT0")/fNevents << " seconds per event " << endl ;
469 cout << endl ;
470 }
471}
472
473//__________________________________________________________________
474void AliTOFT0::SetTZeroFile(char * file ){
475 cout << "Destination file : " << file << endl ;
476 fT0File=file;
477}
478//__________________________________________________________________
5c016a7b 479void AliTOFT0::Print(Option_t* /*option*/)const
d599d913 480{
481 cout << "------------------- "<< GetName() << " -------------" << endl ;
482 if(!fT0File.IsNull())
483 cout << " Writing T0 Distribution to file " << (char*) fT0File.Data() << endl ;
484}
485
486//__________________________________________________________________
487Bool_t AliTOFT0::operator==( AliTOFT0 const &tzero )const
488{
489 // Equal operator.
490 //
491
492 if( (fTimeResolution==tzero.fTimeResolution)&&(fLowerMomBound==tzero.fLowerMomBound)&&(fUpperMomBound==tzero.fUpperMomBound))
493 return kTRUE ;
494 else
495 return kFALSE ;
496}
497