]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliESDpid.cxx
Changes for #82873: Module debugging broken (Christian)
[u/mrichter/AliRoot.git] / STEER / AliESDpid.cxx
CommitLineData
8c6a71ab 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
4f679a16 16/* $Id$ */
17
8c6a71ab 18//-----------------------------------------------------------------
19// Implementation of the combined PID class
4f679a16 20// For the Event Summary Data Class
21// produced by the reconstruction process
22// and containing information on the particle identification
8c6a71ab 23// Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
24//-----------------------------------------------------------------
25
f858b00e 26#include "TArrayI.h"
27#include "TArrayF.h"
28
10d100d4 29#include "AliLog.h"
30#include "AliPID.h"
f858b00e 31#include "AliTOFHeader.h"
8c6a71ab 32#include "AliESDpid.h"
af885e0f 33#include "AliESDEvent.h"
8c6a71ab 34#include "AliESDtrack.h"
35
36ClassImp(AliESDpid)
37
f858b00e 38Int_t AliESDpid::MakePID(AliESDEvent *event, Bool_t TPConly, Float_t timeZeroTOF) const {
10d100d4 39 //
40 // Calculate probabilities for all detectors, except if TPConly==kTRUE
41 // and combine PID
42 //
43 // Option TPConly==kTRUE is used during reconstruction,
44 // because ITS tracking uses TPC pid
45 // HMPID and TRD pid are done in detector reconstructors
46 //
47
48 /*
f858b00e 49 Float_t timeZeroTOF = 0;
10d100d4 50 if (subtractT0)
f858b00e 51 timeZeroTOF = event->GetT0();
10d100d4 52 */
53 Int_t nTrk=event->GetNumberOfTracks();
54 for (Int_t iTrk=0; iTrk<nTrk; iTrk++) {
55 AliESDtrack *track=event->GetTrack(iTrk);
56 MakeTPCPID(track);
57 if (!TPConly) {
58 MakeITSPID(track);
f858b00e 59 MakeTOFPID(track, timeZeroTOF);
10d100d4 60 //MakeHMPIDPID(track);
ce5d6908 61 //MakeTRDPID(track);
10d100d4 62 }
63 CombinePID(track);
64 }
65 return 0;
66}
8c6a71ab 67//_________________________________________________________________________
10d100d4 68void AliESDpid::MakeTPCPID(AliESDtrack *track) const
8c6a71ab 69{
4f679a16 70 //
10d100d4 71 // TPC pid using bethe-bloch and gaussian response
4f679a16 72 //
10d100d4 73 if ((track->GetStatus()&AliESDtrack::kTPCin )==0)
74 if ((track->GetStatus()&AliESDtrack::kTPCout)==0) return;
8c6a71ab 75
10d100d4 76 Double_t mom = track->GetP();
77 const AliExternalTrackParam *in=track->GetInnerParam();
78 if (in) mom = in->GetP();
8c6a71ab 79
10d100d4 80 Double_t p[AliPID::kSPECIES];
81 Double_t dedx=track->GetTPCsignal();
82 Bool_t mismatch=kTRUE, heavy=kTRUE;
8c6a71ab 83
10d100d4 84 for (Int_t j=0; j<AliPID::kSPECIES; j++) {
85 AliPID::EParticleType type=AliPID::EParticleType(j);
86 Double_t bethe=fTPCResponse.GetExpectedSignal(mom,type);
87 Double_t sigma=fTPCResponse.GetExpectedSigma(mom,track->GetTPCsignalN(),type);
88 if (TMath::Abs(dedx-bethe) > fRange*sigma) {
89 p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
90 } else {
91 p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
92 mismatch=kFALSE;
93 }
8c6a71ab 94
10d100d4 95 // Check for particles heavier than (AliPID::kSPECIES - 1)
96 if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
c630aafd 97
4a78b8c5 98 }
99
10d100d4 100 if (mismatch)
77638021 101 for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
10d100d4 102
103 track->SetTPCpid(p);
104
105 if (heavy) track->ResetStatus(AliESDtrack::kTPCpid);
106
107}
108//_________________________________________________________________________
109void AliESDpid::MakeITSPID(AliESDtrack *track) const
110{
111 //
112 // ITS PID
113 // Two options, depending on fITSPIDmethod:
114 // 1) Truncated mean method
115 // 2) Likelihood, using charges measured in all 4 layers and
116 // Landau+gaus response functions
117 //
118
119 if ((track->GetStatus()&AliESDtrack::kITSin)==0 &&
120 (track->GetStatus()&AliESDtrack::kITSout)==0) return;
121
122 Double_t mom=track->GetP();
123 if (fITSPIDmethod == kITSTruncMean) {
124 Double_t dedx=track->GetITSsignal();
15e979c9 125 Bool_t isSA=kTRUE;
126 Double_t momITS=mom;
127 ULong_t trStatus=track->GetStatus();
99daa709 128 if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
15e979c9 129 UChar_t clumap=track->GetITSClusterMap();
130 Int_t nPointsForPid=0;
131 for(Int_t i=2; i<6; i++){
132 if(clumap&(1<<i)) ++nPointsForPid;
133 }
134
135 if(nPointsForPid<3) { // track not to be used for combined PID purposes
136 track->ResetStatus(AliESDtrack::kITSpid);
137 return;
138 }
139
10d100d4 140 Double_t p[10];
15e979c9 141
10d100d4 142 Bool_t mismatch=kTRUE, heavy=kTRUE;
143 for (Int_t j=0; j<AliPID::kSPECIES; j++) {
144 Double_t mass=AliPID::ParticleMass(j);//GeV/c^2
15e979c9 145 Double_t bethe=fITSResponse.Bethe(momITS,mass);
146 Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA);
10d100d4 147 if (TMath::Abs(dedx-bethe) > fRange*sigma) {
148 p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
149 } else {
150 p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
151 mismatch=kFALSE;
152 }
153
154 // Check for particles heavier than (AliPID::kSPECIES - 1)
155 if (dedx < (bethe + fRange*sigma)) heavy=kFALSE;
156
4a78b8c5 157 }
c630aafd 158
10d100d4 159 if (mismatch)
160 for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
161
162 track->SetITSpid(p);
163
164 if (heavy) track->ResetStatus(AliESDtrack::kITSpid);
165 }
166 else { // Likelihood method
167 Double_t condprobfun[AliPID::kSPECIES];
168 Double_t qclu[4];
169 track->GetITSdEdxSamples(qclu);
170 fITSResponse.GetITSProbabilities(mom,qclu,condprobfun);
171 track->SetITSpid(condprobfun);
8c6a71ab 172 }
7a8614f3 173
10d100d4 174}
175//_________________________________________________________________________
f858b00e 176void AliESDpid::MakeTOFPID(AliESDtrack *track, Float_t /*timeZeroTOF*/) const
10d100d4 177{
178 //
179 // TOF PID using gaussian response
180 //
f858b00e 181
10d100d4 182 if ((track->GetStatus()&AliESDtrack::kTOFout)==0) return;
183 if ((track->GetStatus()&AliESDtrack::kTIME)==0) return;
a512bf97 184 if ((track->GetStatus()&AliESDtrack::kITSin)==0) return;
10d100d4 185
f858b00e 186 Int_t ibin = fTOFResponse.GetMomBin(track->GetP());
187 Float_t timezero = fTOFResponse.GetT0bin(ibin);
188
10d100d4 189 Double_t time[AliPID::kSPECIESN];
190 track->GetIntegratedTimes(time);
191
192 Double_t sigma[AliPID::kSPECIES];
193 for (Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++) {
194 sigma[iPart] = fTOFResponse.GetExpectedSigma(track->GetP(),time[iPart],AliPID::ParticleMass(iPart));
195 }
196
197 AliDebugGeneral("AliESDpid::MakeTOFPID",2,
198 Form("Expected TOF signals [ps]: %f %f %f %f %f",
199 time[AliPID::kElectron],
200 time[AliPID::kMuon],
201 time[AliPID::kPion],
202 time[AliPID::kKaon],
203 time[AliPID::kProton]));
204
205 AliDebugGeneral("AliESDpid::MakeTOFPID",2,
206 Form("Expected TOF std deviations [ps]: %f %f %f %f %f",
207 sigma[AliPID::kElectron],
208 sigma[AliPID::kMuon],
209 sigma[AliPID::kPion],
210 sigma[AliPID::kKaon],
211 sigma[AliPID::kProton]
212 ));
213
f858b00e 214 Double_t tof = track->GetTOFsignal() - timezero;
10d100d4 215
216 Double_t p[AliPID::kSPECIES];
217 Bool_t mismatch = kTRUE, heavy = kTRUE;
218 for (Int_t j=0; j<AliPID::kSPECIES; j++) {
219 Double_t sig = sigma[j];
f858b00e 220 if (TMath::Abs(tof-time[j]) > (fRange+2)*sig) {
221 p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
10d100d4 222 } else
223 p[j] = TMath::Exp(-0.5*(tof-time[j])*(tof-time[j])/(sig*sig))/sig;
224
225 // Check the mismatching
226 Double_t mass = AliPID::ParticleMass(j);
227 Double_t pm = fTOFResponse.GetMismatchProbability(track->GetP(),mass);
228 if (p[j]>pm) mismatch = kFALSE;
229
230 // Check for particles heavier than (AliPID::kSPECIES - 1)
231 if (tof < (time[j] + fRange*sig)) heavy=kFALSE;
232
233 }
234
235 if (mismatch)
77638021 236 for (Int_t j=0; j<AliPID::kSPECIES; j++) p[j]=1./AliPID::kSPECIES;
10d100d4 237
238 track->SetTOFpid(p);
239
240 if (heavy) track->ResetStatus(AliESDtrack::kTOFpid);
7170298c 241 if (!CheckTOFMatching(track)) track->SetStatus(AliESDtrack::kTOFmismatch);
f858b00e 242
10d100d4 243}
244//_________________________________________________________________________
b439f460 245void AliESDpid::MakeTRDPID(AliESDtrack *track) const
246{
247 //
248 // Method to recalculate the TRD PID probabilities
249 //
ea235c90 250 Double_t prob[AliPID::kSPECIES];
251 ComputeTRDProbability(track, AliPID::kSPECIES, prob);
b439f460 252 track->SetTRDpid(prob);
253}
254//_________________________________________________________________________
10d100d4 255void AliESDpid::CombinePID(AliESDtrack *track) const
256{
257 //
258 // Combine the information of various detectors
259 // to determine the Particle Identification
260 //
261 Int_t ns=AliPID::kSPECIES;
262 Double_t p[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
263
264 if (track->IsOn(AliESDtrack::kITSpid)) {
265 Double_t d[10];
266 track->GetITSpid(d);
267 for (Int_t j=0; j<ns; j++) p[j]*=d[j];
268 }
269
270 if (track->IsOn(AliESDtrack::kTPCpid)) {
271 Double_t d[10];
272 track->GetTPCpid(d);
273 for (Int_t j=0; j<ns; j++) p[j]*=d[j];
274 }
275
276 if (track->IsOn(AliESDtrack::kTRDpid)) {
277 Double_t d[10];
278 track->GetTRDpid(d);
279 for (Int_t j=0; j<ns; j++) p[j]*=d[j];
280 }
281
282 if (track->IsOn(AliESDtrack::kTOFpid)) {
283 Double_t d[10];
284 track->GetTOFpid(d);
285 for (Int_t j=0; j<ns; j++) p[j]*=d[j];
286 }
287
288 if (track->IsOn(AliESDtrack::kHMPIDpid)) {
289 Double_t d[10];
290 track->GetHMPIDpid(d);
291 for (Int_t j=0; j<ns; j++) p[j]*=d[j];
292 }
293
294 track->SetESDpid(p);
8c6a71ab 295}
f858b00e 296//_________________________________________________________________________
297Bool_t AliESDpid::CheckTOFMatching(AliESDtrack *track) const{
ea235c90 298 //
299 // Check pid matching of TOF with TPC as reference
300 //
f858b00e 301 Bool_t status = kFALSE;
302
303 Double_t exptimes[5];
304 track->GetIntegratedTimes(exptimes);
305
f858b00e 306 Float_t p = track->P();
307
7170298c 308 Float_t dedx = track->GetTPCsignal();
309 Float_t time = track->GetTOFsignal() - fTOFResponse.GetStartTime(p);
310
f858b00e 311 Double_t ptpc[3];
312 track->GetInnerPxPyPz(ptpc);
313 Float_t momtpc=TMath::Sqrt(ptpc[0]*ptpc[0] + ptpc[1]*ptpc[1] + ptpc[2]*ptpc[2]);
314
315 for(Int_t i=0;i < 5;i++){
316 AliPID::EParticleType type=AliPID::EParticleType(i);
317
318 Float_t resolutionTOF = fTOFResponse.GetExpectedSigma(p, exptimes[i], AliPID::ParticleMass(i));
319 if(TMath::Abs(exptimes[i] - time) < fRange * resolutionTOF){
320 Float_t dedxExp = fTPCResponse.GetExpectedSignal(momtpc,type);
321 Float_t resolutionTPC = fTPCResponse.GetExpectedSigma(momtpc,track->GetTPCsignalN(),type);
322
7170298c 323 if(TMath::Abs(dedx - dedxExp) < fRangeTOFMismatch * resolutionTPC){
f858b00e 324 status = kTRUE;
325 }
326 }
327 }
328
329 // for nuclei
c1252298 330 Float_t resolutionTOFpr = fTOFResponse.GetExpectedSigma(p, exptimes[4], AliPID::ParticleMass(4));
331 if(!status && (exptimes[4] + fRange*resolutionTOFpr < time)) status = kTRUE;
f858b00e 332
333
334 return status;
335}
336//_________________________________________________________________________
fca752cb 337void AliESDpid::SetTOFResponse(AliVEvent *vevent,EStartTimeType_t option){
f858b00e 338 //
339 // Set TOF response function
340 // Input option for event_time used
341 //
342
fca752cb 343 AliESDEvent *event=(AliESDEvent*)vevent;
344
f858b00e 345 Float_t t0spread = 0.; //event->GetEventTimeSpread();
346 if(t0spread < 10) t0spread = 80;
347
348 // T0 from TOF algorithm
f858b00e 349
350 Bool_t flagT0TOF=kFALSE;
351 Bool_t flagT0T0=kFALSE;
2f43a5ff 352 Float_t *startTime = new Float_t[fTOFResponse.GetNmomBins()];
353 Float_t *startTimeRes = new Float_t[fTOFResponse.GetNmomBins()];
d4d15b21 354 Int_t *startTimeMask = new Int_t[fTOFResponse.GetNmomBins()];
f858b00e 355
7eb68cb5 356 // T0-TOF arrays
2f43a5ff 357 Float_t *estimatedT0event = new Float_t[fTOFResponse.GetNmomBins()];
358 Float_t *estimatedT0resolution = new Float_t[fTOFResponse.GetNmomBins()];
f858b00e 359 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
360 estimatedT0event[i]=0.0;
361 estimatedT0resolution[i]=0.0;
d4d15b21 362 startTimeMask[i] = 0;
f858b00e 363 }
364
7eb68cb5 365 Float_t resT0A=75,resT0C=65,resT0AC=55;
366 if(event->GetT0TOF()){ // check if T0 detector information is available
f858b00e 367 flagT0T0=kTRUE;
368 }
369
370
371 AliTOFHeader *tofHeader =(AliTOFHeader*)event->GetTOFHeader();
372
7eb68cb5 373 if(tofHeader){ // read global info and T0-TOF info from ESD
4d197ca4 374 fTOFResponse.SetTimeResolution(tofHeader->GetTOFResolution());
7eb68cb5 375 t0spread = tofHeader->GetT0spread(); // read t0 sprad
4d197ca4 376 if(t0spread < 10) t0spread = 80;
f858b00e 377
378 flagT0TOF=kTRUE;
7eb68cb5 379 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){ // read T0-TOF default value
380 startTime[i]=tofHeader->GetDefaultEventTimeVal();
381 startTimeRes[i]=tofHeader->GetDefaultEventTimeRes();
edc91478 382 if(startTimeRes[i] < 1.e-5) startTimeRes[i] = t0spread;
f858b00e 383 }
384
385 TArrayI *ibin=tofHeader->GetNvalues();
386 TArrayF *t0Bin=tofHeader->GetEventTimeValues();
387 TArrayF *t0ResBin=tofHeader->GetEventTimeRes();
7eb68cb5 388 for(Int_t j=0;j < tofHeader->GetNbins();j++){ // fill T0-TOF in p-bins
f858b00e 389 Int_t icurrent = (Int_t)ibin->GetAt(j);
7eb68cb5 390 startTime[icurrent]=t0Bin->GetAt(j);
391 startTimeRes[icurrent]=t0ResBin->GetAt(j);
edc91478 392 if(startTimeRes[icurrent] < 1.e-5) startTimeRes[icurrent] = t0spread;
393 }
f858b00e 394 }
395
7eb68cb5 396 // for cut of 3 sigma on t0 spread
397 Float_t t0cut = 3 * t0spread;
398 if(t0cut < 500) t0cut = 500;
399
f858b00e 400 if(option == kFILL_T0){ // T0-FILL is used
401 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
402 estimatedT0event[i]=0.0;
403 estimatedT0resolution[i]=t0spread;
f858b00e 404 }
405 fTOFResponse.SetT0event(estimatedT0event);
406 fTOFResponse.SetT0resolution(estimatedT0resolution);
407 }
7eb68cb5 408
f858b00e 409 if(option == kTOF_T0){ // T0-TOF is used when available (T0-FILL otherwise) from ESD
410 if(flagT0TOF){
f858b00e 411 fTOFResponse.SetT0event(startTime);
412 fTOFResponse.SetT0resolution(startTimeRes);
d4d15b21 413 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
414 if(startTimeRes[i]<t0spread) startTimeMask[i]=1;
415 fTOFResponse.SetT0binMask(i,startTimeMask[i]);
416 }
f858b00e 417 }
418 else{
419 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
420 estimatedT0event[i]=0.0;
421 estimatedT0resolution[i]=t0spread;
d4d15b21 422 fTOFResponse.SetT0binMask(i,startTimeMask[i]);
f858b00e 423 }
424 fTOFResponse.SetT0event(estimatedT0event);
425 fTOFResponse.SetT0resolution(estimatedT0resolution);
426 }
427 }
d4d15b21 428 else if(option == kBest_T0){ // T0-T0 or T0-TOF are used when available (T0-FILL otherwise) from ESD
f858b00e 429 Float_t t0AC=-10000;
430 Float_t t0A=-10000;
431 Float_t t0C=-10000;
432 if(flagT0T0){
433 t0AC= event->GetT0TOF()[0];
434 t0A= event->GetT0TOF()[1];
435 t0C= event->GetT0TOF()[2];
436 }
437
438 Float_t t0t0Best = 0;
439 Float_t t0t0BestRes = 9999;
d4d15b21 440 Int_t t0used=0;
7eb68cb5 441 if(TMath::Abs(t0A) < t0cut && TMath::Abs(t0C) < t0cut && TMath::Abs(t0C-t0A) < 500){
f858b00e 442 t0t0Best = t0AC;
7eb68cb5 443 t0t0BestRes = resT0AC;
d4d15b21 444 t0used=6;
f858b00e 445 }
7eb68cb5 446 else if(TMath::Abs(t0C) < t0cut){
f858b00e 447 t0t0Best = t0C;
7eb68cb5 448 t0t0BestRes = resT0C;
d4d15b21 449 t0used=4;
7eb68cb5 450 }
451 else if(TMath::Abs(t0A) < t0cut){
452 t0t0Best = t0A;
453 t0t0BestRes = resT0A;
d4d15b21 454 t0used=2;
f858b00e 455 }
456
457 if(flagT0TOF){ // if T0-TOF info is available
458 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
7eb68cb5 459 if(t0t0BestRes < 999){
460 if(startTimeRes[i] < t0spread){
f858b00e 461 Double_t wtot = 1./startTimeRes[i]/startTimeRes[i] + 1./t0t0BestRes/t0t0BestRes;
462 Double_t t0best = startTime[i]/startTimeRes[i]/startTimeRes[i] + t0t0Best/t0t0BestRes/t0t0BestRes;
463 estimatedT0event[i]=t0best / wtot;
464 estimatedT0resolution[i]=1./TMath::Sqrt(wtot);
d4d15b21 465 startTimeMask[i] = t0used+1;
f858b00e 466 }
467 else {
468 estimatedT0event[i]=t0t0Best;
469 estimatedT0resolution[i]=t0t0BestRes;
d4d15b21 470 startTimeMask[i] = t0used;
f858b00e 471 }
472 }
473 else{
474 estimatedT0event[i]=startTime[i];
475 estimatedT0resolution[i]=startTimeRes[i];
d4d15b21 476 if(startTimeRes[i]<t0spread) startTimeMask[i]=1;
f858b00e 477 }
d4d15b21 478 fTOFResponse.SetT0binMask(i,startTimeMask[i]);
f858b00e 479 }
480 fTOFResponse.SetT0event(estimatedT0event);
481 fTOFResponse.SetT0resolution(estimatedT0resolution);
482 }
483 else{ // if no T0-TOF info is available
484 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
d4d15b21 485 fTOFResponse.SetT0binMask(i,t0used);
486 if(t0t0BestRes < 999){
487 estimatedT0event[i]=t0t0Best;
488 estimatedT0resolution[i]=t0t0BestRes;
489 }
490 else{
491 estimatedT0event[i]=0.0;
492 estimatedT0resolution[i]=t0spread;
493 }
f858b00e 494 }
495 fTOFResponse.SetT0event(estimatedT0event);
496 fTOFResponse.SetT0resolution(estimatedT0resolution);
497 }
498 }
7eb68cb5 499
d4d15b21 500 else if(option == kT0_T0){ // T0-T0 is used when available (T0-FILL otherwise) from ESD
f858b00e 501 Float_t t0AC=-10000;
502 Float_t t0A=-10000;
503 Float_t t0C=-10000;
504 if(flagT0T0){
505 t0AC= event->GetT0TOF()[0];
506 t0A= event->GetT0TOF()[1];
507 t0C= event->GetT0TOF()[2];
508 }
509
7eb68cb5 510 if(TMath::Abs(t0A) < t0cut && TMath::Abs(t0C) < t0cut && TMath::Abs(t0C-t0A) < 500){
f858b00e 511 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
512 estimatedT0event[i]=t0AC;
7eb68cb5 513 estimatedT0resolution[i]=resT0AC;
d4d15b21 514 fTOFResponse.SetT0binMask(i,6);
f858b00e 515 }
516 }
7eb68cb5 517 else if(TMath::Abs(t0C) < t0cut){
f858b00e 518 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
7eb68cb5 519 estimatedT0event[i]=t0C;
520 estimatedT0resolution[i]=resT0C;
d4d15b21 521 fTOFResponse.SetT0binMask(i,4);
f858b00e 522 }
523 }
7eb68cb5 524 else if(TMath::Abs(t0A) < t0cut){
f858b00e 525 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
7eb68cb5 526 estimatedT0event[i]=t0A;
527 estimatedT0resolution[i]=resT0A;
d4d15b21 528 fTOFResponse.SetT0binMask(i,2);
f858b00e 529 }
530 }
531 else{
532 for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
533 estimatedT0event[i]=0.0;
534 estimatedT0resolution[i]=t0spread;
d4d15b21 535 fTOFResponse.SetT0binMask(i,0);
f858b00e 536 }
537 }
538 fTOFResponse.SetT0event(estimatedT0event);
539 fTOFResponse.SetT0resolution(estimatedT0resolution);
540 }
2f43a5ff 541 delete [] startTime;
542 delete [] startTimeRes;
d4d15b21 543 delete [] startTimeMask;
2f43a5ff 544 delete [] estimatedT0event;
545 delete [] estimatedT0resolution;
f858b00e 546}