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