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