Line data Source code
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 : // //
18 : // AliTOFtracker Class //
19 : // Task: Perform association of the ESD tracks to TOF Clusters //
20 : // and Update ESD track with associated TOF Cluster parameters //
21 : // //
22 : // -- Authors : S. Arcelli, C. Zampolli (Bologna University and INFN) //
23 : // -- Contacts: Annalisa.De.Caro@cern.ch //
24 : // -- : Chiara.Zampolli@bo.infn.it //
25 : // -- : Silvia.Arcelli@bo.infn.it //
26 : // //
27 : //--------------------------------------------------------------------//
28 :
29 : #include <Rtypes.h>
30 : #include <TROOT.h>
31 :
32 : #include <TSeqCollection.h>
33 : #include <TClonesArray.h>
34 : #include <TObjArray.h>
35 : #include <TGeoManager.h>
36 : #include <TTree.h>
37 : #include <TFile.h>
38 : #include <TH2F.h>
39 :
40 : #include "AliGeomManager.h"
41 : #include "AliESDtrack.h"
42 : #include "AliESDEvent.h"
43 : #include "AliESDpid.h"
44 : #include "AliLog.h"
45 : #include "AliTrackPointArray.h"
46 : #include "AliCDBManager.h"
47 :
48 : //#include "AliTOFpidESD.h"
49 : #include "AliTOFRecoParam.h"
50 : #include "AliTOFReconstructor.h"
51 : #include "AliTOFcluster.h"
52 : #include "AliTOFGeometry.h"
53 : #include "AliTOFtracker.h"
54 : #include "AliTOFtrack.h"
55 :
56 : extern TGeoManager *gGeoManager;
57 : //extern TROOT *gROOT;
58 :
59 :
60 26 : ClassImp(AliTOFtracker)
61 :
62 : //_____________________________________________________________________________
63 2 : AliTOFtracker::AliTOFtracker():
64 2 : fkRecoParam(0x0),
65 2 : fGeom(0x0),
66 2 : fN(0),
67 2 : fNseeds(0),
68 2 : fNseedsTOF(0),
69 2 : fngoodmatch(0),
70 2 : fnbadmatch(0),
71 2 : fnunmatch(0),
72 2 : fnmatch(0),
73 2 : fESDEv(0),
74 6 : fTracks(new TClonesArray("AliTOFtrack")),
75 6 : fSeeds(new TObjArray(100)),
76 6 : fTOFtrackPoints(new TObjArray(10)),
77 2 : fHDigClusMap(0x0),
78 2 : fHDigNClus(0x0),
79 2 : fHDigClusTime(0x0),
80 2 : fHDigClusToT(0x0),
81 2 : fHRecNClus(0x0),
82 2 : fHRecDist(0x0),
83 2 : fHRecSigYVsP(0x0),
84 2 : fHRecSigZVsP(0x0),
85 2 : fHRecSigYVsPWin(0x0),
86 2 : fHRecSigZVsPWin(0x0),
87 2 : fCalTree(0x0),
88 2 : fIch(-1),
89 2 : fToT(-1.),
90 2 : fTime(-1.),
91 2 : fExpTimePi(-1.),
92 2 : fExpTimeKa(-1.),
93 2 : fExpTimePr(-1.),
94 2 : fNTOFmatched(0)
95 10 : {
96 : //AliTOFtracker main Ctor
97 :
98 311112 : for (Int_t ii=0; ii<kMaxCluster; ii++) fClusters[ii]=0x0;
99 :
100 : // Getting the geometry
101 6 : fGeom = new AliTOFGeometry();
102 : /* RS?
103 : for(Int_t i=0; i< 20000;i++){
104 : fClusterESD[i] = NULL;
105 : fHit[i] = NULL;
106 : }
107 : */
108 2 : InitCheckHists();
109 :
110 4 : }
111 : //_____________________________________________________________________________
112 12 : AliTOFtracker::~AliTOFtracker() {
113 : //
114 : // Dtor
115 : //
116 :
117 2 : SaveCheckHists();
118 :
119 4 : if(!(AliCDBManager::Instance()->GetCacheFlag())){
120 0 : delete fkRecoParam;
121 : }
122 4 : delete fGeom;
123 4 : delete fHDigClusMap;
124 4 : delete fHDigNClus;
125 4 : delete fHDigClusTime;
126 4 : delete fHDigClusToT;
127 4 : delete fHRecNClus;
128 4 : delete fHRecDist;
129 4 : delete fHRecSigYVsP;
130 4 : delete fHRecSigZVsP;
131 4 : delete fHRecSigYVsPWin;
132 4 : delete fHRecSigZVsPWin;
133 4 : delete fCalTree;
134 2 : if (fTracks){
135 2 : fTracks->Delete();
136 4 : delete fTracks;
137 2 : fTracks=0x0;
138 2 : }
139 2 : if (fSeeds){
140 2 : fSeeds->Delete();
141 4 : delete fSeeds;
142 2 : fSeeds=0x0;
143 2 : }
144 2 : if (fTOFtrackPoints){
145 2 : fTOFtrackPoints->Delete();
146 4 : delete fTOFtrackPoints;
147 2 : fTOFtrackPoints=0x0;
148 2 : }
149 :
150 311112 : for (Int_t ii=0; ii<kMaxCluster; ii++)
151 155554 : if (fClusters[ii]) fClusters[ii]->Delete();
152 :
153 : /* RS?
154 : for(Int_t i=0; i< 20000;i++){
155 : if(fClusterESD[i]){
156 : delete fClusterESD[i];
157 : fClusterESD[i] = NULL;
158 : }
159 : if(fHit[i]){
160 : delete fHit[i];
161 : fHit[i] = NULL;
162 : }
163 : }
164 : */
165 :
166 6 : }
167 : //_____________________________________________________________________________
168 : void AliTOFtracker::GetPidSettings(AliESDpid *esdPID) {
169 : //
170 : // Sets TOF resolution from RecoParams
171 : //
172 0 : if (fkRecoParam)
173 0 : esdPID->GetTOFResponse().SetTimeResolution(fkRecoParam->GetTimeResolution());
174 : else
175 0 : AliWarning("fkRecoParam not yet set; cannot set PID settings");
176 0 : }
177 : //_____________________________________________________________________________
178 : Int_t AliTOFtracker::PropagateBack(AliESDEvent * const event) {
179 : //
180 : // Gets seeds from ESD event and Match with TOF Clusters
181 : //
182 16 : fESDEv = event;
183 : //
184 8 : if (fN==0) {
185 0 : AliInfo("No TOF recPoints to be matched with reconstructed tracks");
186 0 : return 0;
187 : }
188 :
189 : // initialize RecoParam for current event
190 24 : AliDebug(1,"Initializing params for TOF");
191 :
192 8 : fkRecoParam = AliTOFReconstructor::GetRecoParam(); // instantiate reco param from STEER...
193 :
194 8 : if (fkRecoParam == 0x0) {
195 0 : AliFatal("No Reco Param found for TOF!!!");
196 0 : }
197 : //fkRecoParam->Dump();
198 : //if(fkRecoParam->GetApplyPbPbCuts())fkRecoParam=fkRecoParam->GetPbPbparam();
199 : //fkRecoParam->PrintParameters();
200 :
201 : /* RS?
202 : // create clusters from hit
203 : for(Int_t i=0;i < fNTOFmatched;i++){
204 : fClusterESD[i] = new AliESDTOFCluster(fHit[i],event);
205 : fClusterESD[i]->SetStatus(fClusters[i]->GetStatus());
206 : }
207 : */
208 : //Initialise some counters
209 :
210 8 : fNseeds=0;
211 8 : fNseedsTOF=0;
212 8 : fngoodmatch=0;
213 8 : fnbadmatch=0;
214 8 : fnunmatch=0;
215 8 : fnmatch=0;
216 :
217 8 : Int_t ntrk=event->GetNumberOfTracks();
218 8 : fNseeds = ntrk;
219 :
220 :
221 : //Load ESD tracks into a local Array of ESD Seeds
222 320 : for (Int_t i=0; i<fNseeds; i++){
223 152 : fSeeds->AddLast(event->GetTrack(i));
224 : // event->GetTrack(i)->SetESDEvent(event); // RS: Why this is needed? The event is already set
225 : }
226 : //Prepare ESD tracks candidates for TOF Matching
227 8 : CollectESD();
228 :
229 16 : if (fNseeds==0 || fNseedsTOF==0) {
230 0 : AliInfo("No seeds to try TOF match");
231 0 : fSeeds->Clear();
232 0 : fTracks->Clear();
233 0 : return 0 ;
234 : }
235 :
236 : //First Step with Strict Matching Criterion
237 8 : MatchTracks(0);
238 :
239 : //Second Step with Looser Matching Criterion
240 8 : MatchTracks(1);
241 :
242 : //Third Step without kTOFout flag (just to update clusters)
243 8 : MatchTracks(2);
244 :
245 : //RS? event->SetTOFcluster(fNTOFmatched,fClusterESD);
246 :
247 8 : if (fN==0) {
248 0 : AliInfo("No TOF recPoints to be matched with reconstructed tracks");
249 0 : fSeeds->Clear();
250 0 : fTracks->Clear();
251 0 : return 0;
252 : }
253 :
254 8 : AliInfo(Form("Number of matched tracks = %d (good = %d, bad = %d)",fnmatch,fngoodmatch,fnbadmatch));
255 :
256 : //Update the matched ESD tracks
257 :
258 320 : for (Int_t i=0; i<ntrk; i++) {
259 : // RS: This is a bogus code since t and seed are the same object
260 : // AliESDtrack *t=event->GetTrack(i);
261 152 : AliESDtrack *seed =(AliESDtrack*)fSeeds->At(i);
262 152 : if ( (seed->GetStatus()&AliESDtrack::kTOFin)!=0 ) {
263 : //t->SetStatus(AliESDtrack::kTOFin);
264 : //if(seed->GetTOFsignal()>0){
265 103 : if ( (seed->GetStatus()&AliESDtrack::kTOFout)!=0 ) {
266 : //t->SetStatus(AliESDtrack::kTOFout);
267 : //t->SetTOFsignal(seed->GetTOFsignal());
268 : //t->SetTOFcluster(seed->GetTOFcluster());
269 : //t->SetTOFsignalToT(seed->GetTOFsignalToT());
270 : //t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
271 : //t->SetTOFsignalDz(seed->GetTOFsignalDz());
272 : //t->SetTOFsignalDx(seed->GetTOFsignalDx());
273 : //t->SetTOFDeltaBC(seed->GetTOFDeltaBC());
274 : //t->SetTOFL0L1(seed->GetTOFL0L1());
275 : //t->SetTOFCalChannel(seed->GetTOFCalChannel());
276 81 : Int_t tlab[3]; seed->GetTOFLabel(tlab);
277 : //t->SetTOFLabel(tlab);
278 :
279 81 : Double_t alphaA = (Double_t)seed->GetAlpha();
280 81 : Double_t xA = (Double_t)seed->GetX();
281 81 : Double_t yA = (Double_t)seed->GetY();
282 81 : Double_t zA = (Double_t)seed->GetZ();
283 81 : Double_t p1A = (Double_t)seed->GetSnp();
284 81 : Double_t p2A = (Double_t)seed->GetTgl();
285 81 : Double_t p3A = (Double_t)seed->GetSigned1Pt();
286 81 : const Double_t *covA = (Double_t*)seed->GetCovariance();
287 :
288 : // Make attention, please:
289 : // AliESDtrack::fTOFInfo array does not be stored in the AliESDs.root file
290 : // it is there only for a check during the reconstruction step.
291 81 : Float_t info[10]; seed->GetTOFInfo(info);
292 : //t->SetTOFInfo(info);
293 243 : AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
294 :
295 : // Check done:
296 : // by calling the AliESDtrack::UpdateTrackParams,
297 : // the current track parameters are changed
298 : // and it could cause refit problems.
299 : // We need to update only the following track parameters:
300 : // the track length and expected times.
301 : // Removed AliESDtrack::UpdateTrackParams call
302 : // Called AliESDtrack::SetIntegratedTimes(...) and
303 : // AliESDtrack::SetIntegratedLength() routines.
304 : /*
305 : AliTOFtrack *track = new AliTOFtrack(*seed);
306 : t->UpdateTrackParams(track,AliESDtrack::kTOFout); // to be checked - AdC
307 : delete track;
308 : Double_t time[AliPID::kSPECIESC]; t->GetIntegratedTimes(time);
309 : */
310 :
311 81 : Double_t time[AliPID::kSPECIESC]; seed->GetIntegratedTimes(time,AliPID::kSPECIESC);
312 : //t->SetIntegratedTimes(time);
313 :
314 : //Double_t length = seed->GetIntegratedLength();
315 : //t->SetIntegratedLength(length);
316 :
317 81 : Double_t alphaB = (Double_t)seed->GetAlpha();
318 81 : Double_t xB = (Double_t)seed->GetX();
319 81 : Double_t yB = (Double_t)seed->GetY();
320 81 : Double_t zB = (Double_t)seed->GetZ();
321 81 : Double_t p1B = (Double_t)seed->GetSnp();
322 81 : Double_t p2B = (Double_t)seed->GetTgl();
323 81 : Double_t p3B = (Double_t)seed->GetSigned1Pt();
324 81 : const Double_t *covB = (Double_t*)seed->GetCovariance();
325 243 : AliDebug(3,"Track params -now(before)-:");
326 243 : AliDebug(3,Form(" X: %f(%f), Y: %f(%f), Z: %f(%f) --- alpha: %f(%f)",
327 : xB,xA,
328 : yB,yA,
329 : zB,zA,
330 : alphaB,alphaA));
331 243 : AliDebug(3,Form(" p1: %f(%f), p2: %f(%f), p3: %f(%f)",
332 : p1B,p1A,
333 : p2B,p2A,
334 : p3B,p3A));
335 243 : AliDebug(3,Form(" cov1: %f(%f), cov2: %f(%f), cov3: %f(%f)"
336 : " cov4: %f(%f), cov5: %f(%f), cov6: %f(%f)"
337 : " cov7: %f(%f), cov8: %f(%f), cov9: %f(%f)"
338 : " cov10: %f(%f), cov11: %f(%f), cov12: %f(%f)"
339 : " cov13: %f(%f), cov14: %f(%f), cov15: %f(%f)",
340 : covB[0],covA[0],
341 : covB[1],covA[1],
342 : covB[2],covA[2],
343 : covB[3],covA[3],
344 : covB[4],covA[4],
345 : covB[5],covA[5],
346 : covB[6],covA[6],
347 : covB[7],covA[7],
348 : covB[8],covA[8],
349 : covB[9],covA[9],
350 : covB[10],covA[10],
351 : covB[11],covA[11],
352 : covB[12],covA[12],
353 : covB[13],covA[13],
354 : covB[14],covA[14]
355 : ));
356 243 : AliDebug(2,Form(" TOF params: %6d %f %f %f %f %f %6d %3d %f",
357 : i,
358 : seed->GetTOFsignalRaw(),
359 : seed->GetTOFsignal(),
360 : seed->GetTOFsignalToT(),
361 : seed->GetTOFsignalDz(),
362 : seed->GetTOFsignalDx(),
363 : seed->GetTOFCalChannel(),
364 : seed->GetTOFcluster(),
365 : seed->GetIntegratedLength()));
366 243 : AliDebug(2,Form(" %f %f %f %f %f %f %f %f %f",
367 : time[0], time[1], time[2], time[3], time[4], time[5], time[6], time[7], time[8]));
368 81 : }
369 : }
370 : }
371 : /* RS?
372 : if(fNTOFmatched){
373 : Int_t *matchmap = new Int_t[fNTOFmatched];
374 : event->SetTOFcluster(fNTOFmatched,fClusterESD,matchmap);
375 : for (Int_t i=0; i<ntrk; i++) { // remapping after TOF matching selection
376 : AliESDtrack *t=event->GetTrack(i);
377 : t->ReMapTOFcluster(fNTOFmatched,matchmap);
378 : }
379 :
380 : delete[] matchmap;
381 : }
382 : */
383 :
384 : //Make TOF PID
385 : // Now done in AliESDpid
386 : // fPid->MakePID(event,timeZero);
387 :
388 8 : fSeeds->Clear();
389 : //fTracks->Delete();
390 8 : fTracks->Clear();
391 8 : return 0;
392 :
393 8 : }
394 : //_________________________________________________________________________
395 : void AliTOFtracker::CollectESD() {
396 : //prepare the set of ESD tracks to be matched to clusters in TOF
397 :
398 : Int_t seedsTOF1=0;
399 : Int_t seedsTOF3=0;
400 : Int_t seedsTOF2=0;
401 :
402 16 : TClonesArray &aTOFTrack = *fTracks;
403 320 : for (Int_t i=0; i<fNseeds; i++) {
404 :
405 152 : AliESDtrack *t =(AliESDtrack*)fSeeds->At(i);
406 170 : if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
407 :
408 134 : AliTOFtrack *track = new AliTOFtrack(*t); // New
409 134 : Float_t x = (Float_t)track->GetX(); //New
410 :
411 : // TRD 'good' tracks
412 268 : if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) ) {
413 :
414 446 : AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
415 :
416 : // TRD 'good' tracks, already propagated at 371 cm
417 208 : if( x >= AliTOFGeometry::Rmin() ) {
418 :
419 106 : if ( track->PropagateToInnerTOF() ) {
420 :
421 3 : AliDebug(1,Form(" TRD propagated track till rho = %fcm."
422 : " And then the track has been propagated till rho = %fcm.",
423 : x, (Float_t)track->GetX()));
424 :
425 1 : track->SetSeedIndex(i);
426 1 : t->UpdateTrackParams(track,AliESDtrack::kTOFin);
427 1 : new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
428 1 : fNseedsTOF++;
429 1 : seedsTOF1++;
430 :
431 3 : AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
432 : }
433 4 : delete track;
434 :
435 : }
436 : else { // TRD 'good' tracks, propagated rho<371cm
437 :
438 102 : if ( track->PropagateToInnerTOF() ) {
439 :
440 306 : AliDebug(1,Form(" TRD propagated track till rho = %fcm."
441 : " And then the track has been propagated till rho = %fcm.",
442 : x, (Float_t)track->GetX()));
443 :
444 102 : track->SetSeedIndex(i);
445 102 : t->UpdateTrackParams(track,AliESDtrack::kTOFin);
446 102 : new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
447 102 : fNseedsTOF++;
448 102 : seedsTOF3++;
449 :
450 306 : AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
451 : }
452 204 : delete track;
453 :
454 : }
455 : //delete track;
456 : }
457 :
458 : else { // Propagate the rest of TPCbp
459 :
460 90 : AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
461 :
462 30 : if ( track->PropagateToInnerTOF() ) {
463 :
464 0 : AliDebug(1,Form(" TPC propagated track till rho = %fcm."
465 : " And then the track has been propagated till rho = %fcm.",
466 : x, (Float_t)track->GetX()));
467 :
468 0 : track->SetSeedIndex(i);
469 0 : t->UpdateTrackParams(track,AliESDtrack::kTOFin);
470 0 : new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
471 0 : fNseedsTOF++;
472 0 : seedsTOF2++;
473 :
474 0 : AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
475 : }
476 60 : delete track;
477 : }
478 134 : }
479 :
480 8 : AliInfo(Form("Number of TOF seeds = %d (kTRDout371 = %d, kTRDoutLess371 = %d, !kTRDout = %d)",fNseedsTOF,seedsTOF1,seedsTOF3,seedsTOF2));
481 :
482 : // Sort according uncertainties on track position
483 8 : fTracks->Sort();
484 :
485 8 : }
486 :
487 : //_________________________________________________________________________
488 : void AliTOFtracker::MatchTracks( Int_t mLastStep){
489 :
490 : // Parameters used/regulating the reconstruction
491 :
492 : //static Float_t corrLen=0.;//0.75;
493 : static Float_t detDepth=18.;
494 : static Float_t padDepth=0.5;
495 :
496 : const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
497 : const Float_t kTimeOffset = 0.; // time offset for tracking algorithm [ps]
498 :
499 48 : Float_t dY=AliTOFGeometry::XPad();
500 24 : Float_t dZ=AliTOFGeometry::ZPad();
501 :
502 24 : Float_t sensRadius = fkRecoParam->GetSensRadius();
503 24 : Float_t stepSize = fkRecoParam->GetStepSize();
504 24 : Float_t scaleFact = fkRecoParam->GetWindowScaleFact();
505 24 : Float_t dyMax=fkRecoParam->GetWindowSizeMaxY();
506 24 : Float_t dzMax=fkRecoParam->GetWindowSizeMaxZ();
507 24 : Float_t dCut=fkRecoParam->GetDistanceCut();
508 24 : if (dCut==3. && fNseedsTOF<=10) {
509 : dCut=10.;
510 0 : AliInfo(Form("Matching window=%f, since low multiplicity event (fNseedsTOF=%d)",
511 : dCut, fNseedsTOF));
512 0 : }
513 32 : if(mLastStep == 2) dCut=10.;
514 :
515 24 : if (AliTOFReconstructor::GetExtraTolerance()>0) {
516 0 : dCut += AliTOFReconstructor::GetExtraTolerance();
517 0 : AliInfoF("Extra %.2f tolerance on distance is added: dCut=%.2f",AliTOFReconstructor::GetExtraTolerance(),dCut);
518 0 : }
519 :
520 24 : Double_t maxChi2=fkRecoParam->GetMaxChi2TRD();
521 24 : Bool_t timeWalkCorr = fkRecoParam->GetTimeWalkCorr();
522 24 : if(!mLastStep){
523 24 : AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++");
524 24 : AliDebug(1,Form("TOF sens radius: %f",sensRadius));
525 24 : AliDebug(1,Form("TOF step size: %f",stepSize));
526 24 : AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
527 24 : AliDebug(1,Form("TOF Window max dy: %f",dyMax));
528 24 : AliDebug(1,Form("TOF Window max dz: %f",dzMax));
529 24 : AliDebug(1,Form("TOF distance Cut: %f",dCut));
530 24 : AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
531 24 : AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
532 : }
533 :
534 : //Match ESD tracks to clusters in TOF
535 :
536 : // Get the number of propagation steps
537 24 : Int_t nSteps=(Int_t)(detDepth/stepSize);
538 72 : AliDebug(1,Form(" Number of steps to be done %d",nSteps));
539 :
540 72 : AliDebug(1,"++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
541 :
542 : //PH Arrays (moved outside of the loop)
543 24 : Float_t * trackPos[4];
544 240 : for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
545 24 : Int_t * clind = new Int_t[fN];
546 :
547 : // Some init
548 : const Int_t kNclusterMax = 1000; // related to fN value
549 48024 : TGeoHMatrix global[kNclusterMax];
550 :
551 : //The matching loop
552 666 : for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
553 :
554 309 : fTOFtrackPoints->Delete();
555 :
556 618618 : for (Int_t ii=0; ii<kNclusterMax; ii++)
557 309000 : global[ii] = 0x0;
558 309 : AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
559 618 : AliESDtrack *t =(AliESDtrack*)fSeeds->At(track->GetSeedIndex());
560 : //if ( t->GetTOFsignal()>0. ) continue;
561 677 : if ( ((t->GetStatus()&AliESDtrack::kTOFout)!=0 ) && mLastStep < 2) continue;
562 500 : AliTOFtrack *trackTOFin = new AliTOFtrack(*track);
563 :
564 : // Determine a window around the track
565 250 : Double_t x,par[5];
566 250 : trackTOFin->GetExternalParameters(x,par);
567 250 : Double_t cov[15];
568 250 : trackTOFin->GetExternalCovariance(cov);
569 :
570 500 : if (cov[0]<0. || cov[2]<0.) {
571 0 : AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
572 : //delete trackTOFin;
573 : //continue;
574 : }
575 :
576 : Double_t dphi=
577 500 : scaleFact*
578 250 : ((5*TMath::Sqrt(TMath::Abs(cov[0])) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
579 : Double_t dz=
580 250 : scaleFact*
581 250 : (5*TMath::Sqrt(TMath::Abs(cov[2])) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
582 :
583 750 : Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
584 250 : if (phi<-TMath::Pi())phi+=2*TMath::Pi();
585 250 : if (phi>=TMath::Pi())phi-=2*TMath::Pi();
586 250 : Double_t z=par[1];
587 :
588 : //upper limit on window's size.
589 296 : if (dz> dzMax) dz=dzMax;
590 265 : if (dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
591 :
592 :
593 : // find the clusters in the window of the track
594 : Int_t nc=0;
595 3559 : for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
596 :
597 1650 : if (nc>=kNclusterMax) {
598 0 : AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
599 0 : break;
600 : }
601 :
602 1650 : AliTOFcluster *c=fClusters[k];
603 1897 : if (c->GetZ() > z+dz) break;
604 1467 : if (c->IsUsed() && mLastStep < 2) continue;
605 1339 : if (!c->GetStatus()) {
606 0 : AliDebug(1,"Cluster in channel declared bad!");
607 0 : continue; // skip bad channels as declared in OCDB
608 : }
609 :
610 1339 : Double_t dph=TMath::Abs(c->GetPhi()-phi);
611 1566 : if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
612 2244 : if (TMath::Abs(dph)>dphi) continue;
613 :
614 868 : Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
615 434 : Double_t p[2]={yc, c->GetZ()};
616 434 : Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
617 940 : if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
618 :
619 362 : clind[nc] = k;
620 362 : Char_t path[200];
621 362 : Int_t ind[5];
622 362 : ind[0]=c->GetDetInd(0);
623 362 : ind[1]=c->GetDetInd(1);
624 362 : ind[2]=c->GetDetInd(2);
625 362 : ind[3]=c->GetDetInd(3);
626 362 : ind[4]=c->GetDetInd(4);
627 362 : fGeom->GetVolumePath(ind,path);
628 362 : gGeoManager->cd(path);
629 724 : global[nc] = *gGeoManager->GetCurrentMatrix();
630 362 : nc++;
631 796 : }
632 :
633 500 : if (nc == 0 ) {
634 555 : AliDebug(1,Form("No available clusters for the track number %d",iseed));
635 61 : fnunmatch++;
636 122 : delete trackTOFin;
637 61 : continue;
638 : }
639 :
640 945 : AliDebug(1,Form(" Number of available TOF clusters for the track number %d: %d",iseed,nc));
641 :
642 : //start fine propagation
643 :
644 : Int_t nStepsDone = 0;
645 67926 : for( Int_t istep=0; istep<nSteps; istep++){
646 :
647 : // First of all, propagate the track...
648 33684 : Float_t xs = AliTOFGeometry::RinTOF()+istep*stepSize;
649 67371 : if (!(trackTOFin->PropagateTo(xs))) break;
650 :
651 : // ...and then, if necessary, rotate the track
652 33681 : Double_t ymax = xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
653 33681 : Double_t ysect = trackTOFin->GetY();
654 33681 : if (ysect > ymax) {
655 0 : if (!(trackTOFin->Rotate(AliTOFGeometry::GetAlpha()))) break;
656 33681 : } else if (ysect <-ymax) {
657 0 : if (!(trackTOFin->Rotate(-AliTOFGeometry::GetAlpha()))) break;
658 : }
659 :
660 33681 : nStepsDone++;
661 168405 : AliDebug(3,Form(" current step %d (%d) - nStepsDone=%d",istep,nSteps,nStepsDone));
662 :
663 : // store the running point (Globalrf) - fine propagation
664 :
665 33681 : Double_t r[3];
666 33681 : trackTOFin->GetXYZ(r);
667 33681 : trackPos[0][nStepsDone-1]= (Float_t) r[0];
668 33681 : trackPos[1][nStepsDone-1]= (Float_t) r[1];
669 33681 : trackPos[2][nStepsDone-1]= (Float_t) r[2];
670 67362 : trackPos[3][nStepsDone-1]= trackTOFin->GetIntegratedLength();
671 33681 : }
672 :
673 :
674 : #if 0
675 : /*****************/
676 : /**** OLD CODE ***/
677 : /*****************/
678 :
679 : Int_t nfound = 0;
680 : Bool_t accept = kFALSE;
681 : Bool_t isInside = kFALSE;
682 : for (Int_t istep=0; istep<nStepsDone; istep++) {
683 :
684 : Float_t ctrackPos[3];
685 : ctrackPos[0] = trackPos[0][istep];
686 : ctrackPos[1] = trackPos[1][istep];
687 : ctrackPos[2] = trackPos[2][istep];
688 :
689 : //now see whether the track matches any of the TOF clusters
690 :
691 : Float_t dist3d[3];
692 : accept = kFALSE;
693 : for (Int_t i=0; i<nc; i++) {
694 : isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
695 :
696 : if ( mLastStep ) {
697 : Float_t yLoc = dist3d[1];
698 : Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
699 : accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
700 : AliDebug(3," I am in the case mLastStep==kTRUE ");
701 : }
702 : else {
703 : accept = isInside;
704 : }
705 : if (accept) {
706 :
707 : fTOFtrackPoints->AddLast(new AliTOFtrackPoint(clind[i],
708 : TMath::Sqrt(dist3d[0]*dist3d[0] + dist3d[1]*dist3d[1] + dist3d[2]*dist3d[2]),
709 : dist3d[2],dist3d[0],dist3d[1],
710 : AliTOFGeometry::RinTOF()+istep*stepSize,trackPos[3][istep]));
711 :
712 : AliDebug(3,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
713 : nfound++;
714 : if(accept &&!mLastStep)break;
715 : }//end if accept
716 :
717 : } //end for on the clusters
718 : if(accept &&!mLastStep)break;
719 : } //end for on the steps
720 :
721 : /*****************/
722 : /**** OLD CODE ***/
723 : /*****************/
724 : #endif
725 :
726 378 : if ( nStepsDone == 0 ) {
727 189 : AliDebug(1,Form(" No track points for the track number %d",iseed));
728 0 : fnunmatch++;
729 0 : delete trackTOFin;
730 0 : continue;
731 : }
732 :
733 945 : AliDebug(2,Form(" Number of steps done for the track number %d: %d",iseed,nStepsDone));
734 :
735 : /*****************/
736 : /**** NEW CODE ***/
737 : /*****************/
738 :
739 : Int_t *isClusterMatchable = NULL;
740 189 : if(nc){
741 378 : isClusterMatchable = new Int_t[nc];
742 1102 : for (Int_t i=0; i<nc; i++) isClusterMatchable[i] = kFALSE;
743 189 : }
744 :
745 : Int_t nfound = 0;
746 : Bool_t accept = kFALSE;
747 : Bool_t isInside = kFALSE;
748 67740 : for (Int_t istep=0; istep<nStepsDone; istep++) {
749 :
750 : Bool_t gotInsideCluster = kFALSE;
751 : Int_t trackInsideCluster = -1;
752 :
753 33681 : Float_t ctrackPos[3];
754 33681 : ctrackPos[0] = trackPos[0][istep];
755 33681 : ctrackPos[1] = trackPos[1][istep];
756 33681 : ctrackPos[2] = trackPos[2][istep];
757 :
758 : //now see whether the track matches any of the TOF clusters
759 :
760 33681 : Float_t dist3d[3]={0.,0.,0.};
761 : accept = kFALSE;
762 197004 : for (Int_t i=0; i<nc; i++) {
763 :
764 : // ***** NEW *****
765 : /* check whether track was inside another cluster
766 : * and in case inhibit this cluster.
767 : * this will allow to only go on and add track points for
768 : * that cluster where the track got inside first */
769 65082 : if (gotInsideCluster && trackInsideCluster != i) {
770 1305 : AliDebug(3,Form(" A - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
771 : continue;
772 : }
773 322800 : AliDebug(3,Form(" B - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
774 :
775 : /* check whether track is inside this cluster */
776 516480 : for (Int_t hh=0; hh<3; hh++) dist3d[hh]=0.;
777 129120 : isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
778 :
779 : // ***** NEW *****
780 : /* if track is inside this cluster set flags which will then
781 : * inhibit to add track points for the other clusters */
782 64560 : if (isInside) {
783 : gotInsideCluster = kTRUE;
784 : trackInsideCluster = i;
785 587 : }
786 :
787 64560 : if ( mLastStep ) {
788 36541 : Float_t yLoc = dist3d[1];
789 36541 : Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
790 73966 : accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
791 182705 : AliDebug(3," I am in the case mLastStep==kTRUE ");
792 36541 : }
793 :
794 : //***** NEW *****
795 : /* add point everytime that:
796 : * - the track is inside the cluster
797 : * - the track got inside the cluster, even when it eventually exited the cluster
798 : * - the tracks is within dCut from the cluster
799 : */
800 191735 : if (accept || isInside || gotInsideCluster) {
801 :
802 4472 : fTOFtrackPoints->AddLast(new AliTOFtrackPoint(clind[i],
803 1118 : TMath::Sqrt(dist3d[0]*dist3d[0] + dist3d[1]*dist3d[1] + dist3d[2]*dist3d[2]),
804 1118 : dist3d[2],dist3d[0],dist3d[1],
805 1118 : AliTOFGeometry::RinTOF()+istep*stepSize,trackPos[3][istep]));
806 :
807 5590 : AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
808 1118 : nfound++;
809 :
810 5590 : AliDebug(3,Form(" C - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
811 :
812 : // store the match in the ESD
813 1756 : if (mLastStep==2 && !isClusterMatchable[i]) { // add TOF clusters to the track
814 : //
815 132 : isClusterMatchable[i] = kTRUE;
816 : //Tracking info
817 132 : Double_t mom=t->GetP();
818 660 : AliDebug(3,Form(" Momentum for track %d -> %f", iseed,mom));
819 132 : Double_t time[AliPID::kSPECIESC];
820 : // read from old structure (the one used by TPC in reco)
821 2640 : for(Int_t isp=0;isp<AliPID::kSPECIESC;isp++){
822 1188 : time[isp] = t->GetIntegratedTimesOld(isp); // in ps
823 2376 : Double_t mass=AliPID::ParticleMass(isp);
824 2376 : Double_t momz = mom*AliPID::ParticleCharge(isp);
825 1188 : time[isp]+=(trackPos[3][istep]-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(momz*momz+mass*mass)/momz;
826 : //time[isp]+=(trackPos[3][istep]-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
827 : }
828 : //
829 132 : AliESDTOFCluster* esdTOFCl = GetESDTOFCluster(clind[i]);
830 396 : if(!esdTOFCl->Update(t->GetID(),dist3d[0],dist3d[1],dist3d[2],trackPos[3][istep],time))//x,y,z -> tracking RF
831 102 : t->AddTOFcluster(esdTOFCl->GetESDID());
832 132 : }
833 :
834 : // ***** NEW *****
835 : /* do not break loop in any case
836 : * if the track got inside a cluster all other clusters
837 : * are inhibited */
838 : // if(accept &&!mLastStep)break;
839 :
840 : }//end if accept
841 :
842 : } //end for on the clusters
843 :
844 : // ***** NEW *****
845 : /* do not break loop in any case
846 : * if the track got inside a cluster all other clusters
847 : * are inhibited but we want to go on adding track points */
848 : // if(accept &&!mLastStep)break;
849 :
850 33681 : } //end for on the steps
851 567 : if(nc) delete[] isClusterMatchable;
852 :
853 378 : if (nfound == 0 ) {
854 319 : AliDebug(1,Form("No track points for the track number %d",iseed));
855 26 : fnunmatch++;
856 52 : delete trackTOFin;
857 26 : continue;
858 : }
859 :
860 815 : AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
861 :
862 : // now choose the cluster to be matched with the track.
863 :
864 : Int_t idclus=-1;
865 : Float_t recL = 0.;
866 : Float_t xpos=0.;
867 : Float_t mindist=1000.;
868 : Float_t mindistZ=0.;
869 : Float_t mindistY=0.;
870 : Float_t mindistX=stepSize;
871 2562 : for (Int_t iclus= 0; iclus<nfound;iclus++) {
872 2236 : AliTOFtrackPoint *matchableTOFcluster = (AliTOFtrackPoint*)fTOFtrackPoints->At(iclus);
873 : //if ( matchableTOFcluster->Distance()<mindist ) {
874 1362 : if ( TMath::Abs(matchableTOFcluster->DistanceX())<TMath::Abs(mindistX) &&
875 244 : TMath::Abs(matchableTOFcluster->DistanceX())<=stepSize ) {
876 244 : mindist = matchableTOFcluster->Distance();
877 244 : mindistZ = matchableTOFcluster->DistanceZ(); // Z distance in the
878 : // RF of the hit pad
879 : // closest to the
880 : // reconstructed
881 : // track
882 244 : mindistY = matchableTOFcluster->DistanceY(); // Y distance in the
883 : // RF of the hit pad
884 : // closest to the
885 : // reconstructed
886 : // track
887 244 : mindistX = matchableTOFcluster->DistanceX(); // X distance in the
888 : // RF of the hit pad
889 : // closest to the
890 : // reconstructed
891 : // track
892 244 : xpos = matchableTOFcluster->PropRadius();
893 244 : idclus = matchableTOFcluster->Index();
894 244 : recL = matchableTOFcluster->Length();// + corrLen*0.5;
895 :
896 1220 : AliDebug(2,Form(" %d(%d) --- %f (%f, %f, %f), step=%f -- idclus=%d --- seed=%d, trackId=%d, trackLab=%d", iclus,nfound,
897 : mindist,mindistX,mindistY,mindistZ,stepSize,idclus,iseed,track->GetSeedIndex(),track->GetLabel()));
898 :
899 : }
900 : } // loop on found TOF track points
901 :
902 163 : if (TMath::Abs(mindistX)>stepSize && idclus!=-1) {
903 0 : AliInfo(Form(" %d - not matched --- but idclus=%d, trackId=%d, trackLab=%d",iseed,
904 : idclus,track->GetSeedIndex(),track->GetLabel()));
905 : idclus=-1;
906 0 : }
907 :
908 326 : if (idclus==-1) {
909 163 : AliDebug(1,Form("Reconstructed track %d doesn't match any TOF cluster", iseed));
910 0 : fnunmatch++;
911 0 : delete trackTOFin;
912 0 : continue;
913 : }
914 :
915 815 : AliDebug(1,Form(" %d - matched",iseed));
916 :
917 163 : fnmatch++;
918 :
919 163 : AliTOFcluster *c=fClusters[idclus];
920 :
921 815 : AliDebug(3, Form("%7d %7d %10d %10d %10d %10d %7d",
922 : iseed,
923 : fnmatch-1,
924 : TMath::Abs(trackTOFin->GetLabel()),
925 : c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
926 : idclus)); // AdC
927 :
928 163 : c->Use();
929 :
930 : // Track length correction for matching Step 2
931 : /*
932 : if (mLastStep) {
933 : Float_t rc = TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
934 : Float_t rt = TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
935 : +trackPos[1][70]*trackPos[1][70]
936 : +trackPos[2][70]*trackPos[2][70]);
937 : Float_t dlt=rc-rt;
938 : recL=trackPos[3][70]+dlt;
939 : }
940 : */
941 : if (
942 420 : (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
943 163 : ||
944 188 : (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
945 94 : ||
946 188 : (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
947 : ) {
948 69 : fngoodmatch++;
949 :
950 345 : AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
951 :
952 : }
953 : else {
954 94 : fnbadmatch++;
955 :
956 470 : AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
957 :
958 : }
959 :
960 326 : delete trackTOFin;
961 :
962 : // Store quantities to be used in the TOF Calibration
963 163 : Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
964 : // t->SetTOFsignalToT(tToT);
965 163 : Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
966 : // t->SetTOFsignalRaw(rawTime);
967 : // t->SetTOFsignalDz(mindistZ);
968 : // t->SetTOFsignalDx(mindistY);
969 : // t->SetTOFDeltaBC(c->GetDeltaBC());
970 : // t->SetTOFL0L1(c->GetL0L1Latency());
971 :
972 163 : Float_t info[10] = {mindist,mindistY,mindistZ,
973 : 0.,0.,0.,0.,0.,0.,0.};
974 163 : t->SetTOFInfo(info);
975 815 : AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
976 :
977 :
978 163 : Int_t ind[5];
979 163 : ind[0]=c->GetDetInd(0);
980 163 : ind[1]=c->GetDetInd(1);
981 163 : ind[2]=c->GetDetInd(2);
982 163 : ind[3]=c->GetDetInd(3);
983 163 : ind[4]=c->GetDetInd(4);
984 163 : Int_t calindex = AliTOFGeometry::GetIndex(ind);
985 : // t->SetTOFCalChannel(calindex);
986 :
987 : // keep track of the track labels in the matched cluster
988 : Int_t tlab[3];
989 163 : tlab[0]=c->GetLabel(0);
990 163 : tlab[1]=c->GetLabel(1);
991 163 : tlab[2]=c->GetLabel(2);
992 815 : AliDebug(3,Form(" tdc time of the matched track %6d = ",c->GetTDC()));
993 163 : Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
994 815 : AliDebug(3,Form(" tof time of the matched track: %f = ",tof));
995 : Double_t tofcorr=tof;
996 163 : if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
997 815 : AliDebug(3,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
998 : //Set TOF time signal and pointer to the matched cluster
999 : // t->SetTOFsignal(tofcorr);
1000 163 : t->SetTOFcluster(idclus); // pointing to the recPoints tree
1001 :
1002 815 : AliDebug(3,Form(" Setting TOF raw time: %f, z distance: %f corrected time: %f ",rawTime,mindistZ,tofcorr));
1003 :
1004 : //Tracking info
1005 163 : Double_t time[AliPID::kSPECIESC];
1006 : // read from old structure (the one used by TPC in reco)
1007 3260 : for(Int_t isp=0;isp<AliPID::kSPECIESC;isp++){
1008 1467 : time[isp] = t->GetIntegratedTimesOld(isp); // in ps
1009 : }
1010 163 : Double_t mom=t->GetP();
1011 815 : AliDebug(3,Form(" Momentum for track %d -> %f", iseed,mom));
1012 3260 : for (Int_t j=0;j<AliPID::kSPECIESC;j++) {
1013 2934 : Double_t mass=AliPID::ParticleMass(j);
1014 2934 : Double_t momz = mom*AliPID::ParticleCharge(j);
1015 1467 : time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(momz*momz+mass*mass)/momz;
1016 : //time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
1017 : }
1018 :
1019 326 : AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
1020 326 : if (!(trackTOFout->PropagateTo(xpos))) {
1021 0 : delete trackTOFout;
1022 0 : continue;
1023 : }
1024 :
1025 : // If necessary, rotate the track
1026 163 : Double_t yATxposMax=xpos*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
1027 163 : Double_t yATxpos=trackTOFout->GetY();
1028 163 : if (yATxpos > yATxposMax) {
1029 0 : if (!(trackTOFout->Rotate(AliTOFGeometry::GetAlpha()))) {
1030 0 : delete trackTOFout;
1031 0 : continue;
1032 : }
1033 163 : } else if (yATxpos < -yATxposMax) {
1034 0 : if (!(trackTOFout->Rotate(-AliTOFGeometry::GetAlpha()))) {
1035 0 : delete trackTOFout;
1036 0 : continue;
1037 : }
1038 : }
1039 :
1040 : // Fill the track residual histograms and update track only if in the first two step (0 and 1)
1041 163 : if(mLastStep < 2){
1042 81 : FillResiduals(trackTOFout,c,kFALSE);
1043 :
1044 81 : t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
1045 :
1046 : // don't update old structure with TOF info
1047 : // t->SetIntegratedLength(recL);
1048 : // t->SetIntegratedTimes(time);
1049 : // t->SetTOFLabel(tlab);
1050 :
1051 : // add tof cluster to the track also for step 2
1052 81 : AliESDTOFCluster* esdTOFCl = GetESDTOFCluster(idclus);
1053 162 : esdTOFCl->Update(t->GetID(),mindistY,mindist,mindistZ,recL,time);
1054 162 : t->AddTOFcluster(esdTOFCl->GetESDID());
1055 : /* RS?
1056 : if(idclus < 20000){
1057 : fClusterESD[idclus]->Update(t->GetID(),mindistY,mindist,mindistZ,recL,time);//x,y,z -> tracking RF
1058 :
1059 : t->AddTOFcluster(idclus);
1060 : }
1061 : else{
1062 : AliInfo("Too many TOF clusters matched with tracks (> 20000)");
1063 : }
1064 : */
1065 81 : }
1066 : // Fill Reco-QA histos for Reconstruction
1067 163 : fHRecNClus->Fill(nc);
1068 163 : fHRecDist->Fill(mindist);
1069 326 : if (cov[0]>=0.)
1070 326 : fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
1071 : else
1072 0 : fHRecSigYVsP->Fill(mom,-TMath::Sqrt(-cov[0]));
1073 326 : if (cov[2]>=0.)
1074 326 : fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
1075 : else
1076 0 : fHRecSigZVsP->Fill(mom,-TMath::Sqrt(-cov[2]));
1077 163 : fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
1078 163 : fHRecSigZVsPWin->Fill(mom,dz);
1079 :
1080 : // Fill Tree for on-the-fly offline Calibration
1081 :
1082 326 : if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 ) ) {
1083 163 : fIch=calindex;
1084 163 : fToT=tToT;
1085 163 : fTime=rawTime;
1086 163 : fExpTimePi=time[2];
1087 163 : fExpTimeKa=time[3];
1088 163 : fExpTimePr=time[4];
1089 163 : fCalTree->Fill();
1090 : }
1091 326 : delete trackTOFout;
1092 576 : }
1093 :
1094 :
1095 432 : for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
1096 48 : delete [] clind;
1097 :
1098 24048 : }
1099 : //_________________________________________________________________________
1100 : Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
1101 : //--------------------------------------------------------------------
1102 : //This function loads the TOF clusters
1103 : //--------------------------------------------------------------------
1104 :
1105 16 : Int_t npadX = AliTOFGeometry::NpadX();
1106 8 : Int_t npadZ = AliTOFGeometry::NpadZ();
1107 8 : Int_t nStripA = AliTOFGeometry::NStripA();
1108 8 : Int_t nStripB = AliTOFGeometry::NStripB();
1109 8 : Int_t nStripC = AliTOFGeometry::NStripC();
1110 :
1111 8 : TBranch *branch=cTree->GetBranch("TOF");
1112 8 : if (!branch) {
1113 0 : AliError("can't get the branch with the TOF clusters !");
1114 0 : return 1;
1115 : }
1116 :
1117 14 : static TClonesArray dummy("AliTOFcluster",10000);
1118 8 : dummy.Clear();
1119 8 : TClonesArray *clusters=&dummy;
1120 8 : branch->SetAddress(&clusters);
1121 :
1122 8 : cTree->GetEvent(0);
1123 8 : Int_t nc=clusters->GetEntriesFast();
1124 8 : fHDigNClus->Fill(nc);
1125 :
1126 8 : AliInfo(Form("Number of clusters: %d",nc));
1127 :
1128 8 : fN = 0;
1129 8 : fNTOFmatched = 0;
1130 :
1131 416 : for (Int_t i=0; i<nc; i++) {
1132 200 : AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
1133 : //PH fClusters[i]=new AliTOFcluster(*c); fN++;
1134 :
1135 200 : if (!c->Misalign()) AliWarning("Can't misalign this cluster !"); // RS
1136 :
1137 200 : fClusters[i]=c; fN++;
1138 200 : c->SetESDID(-1);
1139 : // Fill Digits QA histos
1140 :
1141 200 : Int_t isector = c->GetDetInd(0);
1142 200 : Int_t iplate = c->GetDetInd(1);
1143 200 : Int_t istrip = c->GetDetInd(2);
1144 200 : Int_t ipadX = c->GetDetInd(4);
1145 200 : Int_t ipadZ = c->GetDetInd(3);
1146 :
1147 200 : Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
1148 200 : Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
1149 :
1150 : /* RS?
1151 : Int_t ind[5];
1152 : ind[0]=isector;
1153 : ind[1]=iplate;
1154 : ind[2]=istrip;
1155 : ind[3]=ipadZ;
1156 : ind[4]=ipadX;
1157 : Int_t calindex = AliTOFGeometry::GetIndex(ind);
1158 : Int_t tofLabels[3]={c->GetLabel(0),c->GetLabel(1),c->GetLabel(2)};
1159 : */
1160 : Int_t stripOffset = 0;
1161 200 : switch (iplate) {
1162 : case 0:
1163 : stripOffset = 0;
1164 18 : break;
1165 : case 1:
1166 : stripOffset = nStripC;
1167 26 : break;
1168 : case 2:
1169 138 : stripOffset = nStripC+nStripB;
1170 138 : break;
1171 : case 3:
1172 18 : stripOffset = nStripC+nStripB+nStripA;
1173 18 : break;
1174 : case 4:
1175 0 : stripOffset = nStripC+nStripB+nStripA+nStripB;
1176 0 : break;
1177 : default:
1178 0 : AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1179 0 : break;
1180 : };
1181 200 : Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
1182 200 : Int_t phiindex=npadX*isector+ipadX+1;
1183 200 : fHDigClusMap->Fill(zindex,phiindex);
1184 200 : fHDigClusTime->Fill(time);
1185 200 : fHDigClusToT->Fill(tot);
1186 :
1187 200 : fNTOFmatched++; // RS: Actually number of clusters
1188 : /* RS?
1189 : if(fNTOFmatched < 20000){
1190 : fHit[fNTOFmatched] = new AliESDTOFHit(AliTOFGeometry::TdcBinWidth()*c->GetTDC(),
1191 : AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW(),
1192 : AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3,
1193 : calindex,tofLabels,c->GetL0L1Latency(),
1194 : c->GetDeltaBC(),i,c->GetZ(),c->GetR(),c->GetPhi());
1195 : fNTOFmatched++;
1196 : }
1197 : */
1198 : }
1199 :
1200 : return 0;
1201 16 : }
1202 : //_________________________________________________________________________
1203 : void AliTOFtracker::UnloadClusters() {
1204 : //--------------------------------------------------------------------
1205 : //This function unloads TOF clusters
1206 : //--------------------------------------------------------------------
1207 424 : for (Int_t i=0; i<fN; i++) {
1208 : //PH delete fClusters[i];
1209 200 : fClusters[i] = 0x0;
1210 : }
1211 : /* RS
1212 : for(Int_t i=0; i< 20000;i++){
1213 : if(fClusterESD[i]){
1214 : delete fClusterESD[i];
1215 : fClusterESD[i] = NULL;
1216 : }
1217 : if(fHit[i]){
1218 : delete fHit[i];
1219 : fHit[i] = NULL;
1220 : }
1221 : }
1222 : */
1223 8 : fN=0;
1224 8 : fNTOFmatched = 0;
1225 8 : }
1226 :
1227 : //_________________________________________________________________________
1228 : Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
1229 : //--------------------------------------------------------------------
1230 : // This function returns the index of the nearest cluster
1231 : //--------------------------------------------------------------------
1232 500 : if (fN==0) return 0;
1233 292 : if (z <= fClusters[0]->GetZ()) return 0;
1234 208 : if (z > fClusters[fN-1]->GetZ()) return fN;
1235 208 : Int_t b=0, e=fN-1, m=(b+e)/2;
1236 2320 : for (; b<e; m=(b+e)/2) {
1237 1390 : if (z > fClusters[m]->GetZ()) b=m+1;
1238 : else e=m;
1239 : }
1240 : return m;
1241 250 : }
1242 :
1243 : //_________________________________________________________________________
1244 : Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
1245 : {
1246 : // Get track space point with index i
1247 : // Coordinates are in the global system
1248 82 : AliTOFcluster *cl = fClusters[index];
1249 : Float_t xyz[3];
1250 41 : xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
1251 41 : xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
1252 41 : xyz[2] = cl->GetZ();
1253 41 : Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
1254 41 : Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
1255 41 : Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
1256 41 : Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
1257 41 : Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
1258 41 : Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
1259 41 : Float_t cov[6];
1260 41 : cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
1261 41 : cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
1262 41 : cov[2] = -cosphi*sinth*costh*sigmaz2;
1263 41 : cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
1264 41 : cov[4] = -sinphi*sinth*costh*sigmaz2;
1265 41 : cov[5] = costh*costh*sigmaz2;
1266 41 : p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
1267 :
1268 : // Detector numbering scheme
1269 41 : Int_t nSector = AliTOFGeometry::NSectors();
1270 41 : Int_t nPlate = AliTOFGeometry::NPlates();
1271 41 : Int_t nStripA = AliTOFGeometry::NStripA();
1272 41 : Int_t nStripB = AliTOFGeometry::NStripB();
1273 41 : Int_t nStripC = AliTOFGeometry::NStripC();
1274 :
1275 41 : Int_t isector = cl->GetDetInd(0);
1276 41 : if (isector >= nSector)
1277 0 : AliError(Form("Wrong sector number in TOF (%d) !",isector));
1278 41 : Int_t iplate = cl->GetDetInd(1);
1279 41 : if (iplate >= nPlate)
1280 0 : AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1281 41 : Int_t istrip = cl->GetDetInd(2);
1282 :
1283 : Int_t stripOffset = 0;
1284 41 : switch (iplate) {
1285 : case 0:
1286 : stripOffset = 0;
1287 0 : break;
1288 : case 1:
1289 : stripOffset = nStripC;
1290 6 : break;
1291 : case 2:
1292 28 : stripOffset = nStripC+nStripB;
1293 28 : break;
1294 : case 3:
1295 7 : stripOffset = nStripC+nStripB+nStripA;
1296 7 : break;
1297 : case 4:
1298 0 : stripOffset = nStripC+nStripB+nStripA+nStripB;
1299 0 : break;
1300 : default:
1301 0 : AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1302 0 : break;
1303 : };
1304 :
1305 41 : Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
1306 41 : stripOffset +
1307 : istrip;
1308 41 : UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
1309 41 : p.SetVolumeID((UShort_t)volid);
1310 41 : return kTRUE;
1311 41 : }
1312 : //_________________________________________________________________________
1313 : void AliTOFtracker::InitCheckHists() {
1314 :
1315 : //Init histos for Digits/Reco QA and Calibration
1316 :
1317 :
1318 4 : TDirectory *dir = gDirectory;
1319 : TFile *logFileTOF = 0;
1320 :
1321 2 : TSeqCollection *list = gROOT->GetListOfFiles();
1322 2 : int n = list->GetEntries();
1323 : Bool_t isThere=kFALSE;
1324 16 : for(int i=0; i<n; i++) {
1325 5 : logFileTOF = (TFile*)list->At(i);
1326 5 : if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1327 : isThere=kTRUE;
1328 0 : break;
1329 : }
1330 : }
1331 :
1332 6 : if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
1333 2 : logFileTOF->cd();
1334 :
1335 4 : fCalTree = new TTree("CalTree", "Tree for TOF calibration");
1336 2 : fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
1337 2 : fCalTree->Branch("ToT",&fToT,"TOFToT/F");
1338 2 : fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
1339 2 : fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
1340 2 : fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
1341 2 : fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
1342 :
1343 : //Digits "QA"
1344 4 : fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
1345 4 : fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
1346 4 : fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
1347 4 : fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
1348 :
1349 : //Reco "QA"
1350 4 : fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
1351 4 : fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
1352 4 : fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
1353 4 : fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
1354 4 : fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
1355 4 : fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
1356 :
1357 2 : dir->cd();
1358 :
1359 2 : }
1360 :
1361 : //_________________________________________________________________________
1362 : void AliTOFtracker::SaveCheckHists() {
1363 :
1364 : //write histos for Digits/Reco QA and Calibration
1365 :
1366 4 : TDirectory *dir = gDirectory;
1367 : TFile *logFileTOF = 0;
1368 :
1369 2 : TSeqCollection *list = gROOT->GetListOfFiles();
1370 2 : int n = list->GetEntries();
1371 : Bool_t isThere=kFALSE;
1372 12 : for(int i=0; i<n; i++) {
1373 6 : logFileTOF = (TFile*)list->At(i);
1374 6 : if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1375 : isThere=kTRUE;
1376 2 : break;
1377 : }
1378 : }
1379 :
1380 2 : if(!isThere) {
1381 0 : AliError(Form("File TOFQA.root not found!! not wring histograms...."));
1382 0 : return;
1383 : }
1384 2 : logFileTOF->cd();
1385 2 : fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
1386 2 : fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
1387 2 : fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
1388 2 : fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
1389 2 : fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
1390 2 : fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
1391 2 : fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
1392 2 : fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
1393 2 : fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
1394 2 : fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
1395 2 : fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
1396 2 : logFileTOF->Flush();
1397 :
1398 2 : dir->cd();
1399 4 : }
1400 : //_________________________________________________________________________
1401 : Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) const {
1402 :
1403 : //dummy, for the moment
1404 : Float_t tofcorr=0.;
1405 0 : if(dist<AliTOFGeometry::ZPad()*0.5){
1406 : tofcorr=tof;
1407 : //place here the actual correction
1408 : }else{
1409 : tofcorr=tof;
1410 : }
1411 0 : return tofcorr;
1412 : }
1413 : //_________________________________________________________________________
1414 :
1415 : void AliTOFtracker::FillClusterArray(TObjArray* arr) const
1416 : {
1417 : //
1418 : // Returns the TOF cluster array
1419 : //
1420 :
1421 0 : if (fN==0)
1422 0 : arr = 0x0;
1423 : else
1424 0 : for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
1425 :
1426 0 : }
1427 :
1428 : //_________________________________________________________________________
1429 : AliESDTOFCluster* AliTOFtracker::GetESDTOFCluster(int clID)
1430 : {
1431 : // get ESDTOFcluster corresponding to fClusters[clID]. If the original cluster
1432 : // was not stored yet in the ESD, first do this
1433 426 : AliTOFcluster *c = fClusters[clID];
1434 : AliESDTOFCluster *clESD = 0;
1435 213 : int esdID = c->GetESDID(); // was this cluster already stored in the ESD clusters?
1436 213 : TClonesArray* esdTOFClArr = fESDEv->GetESDTOFClusters();
1437 426 : if (esdID<0) { // cluster was not stored yet, do this
1438 325 : esdID = esdTOFClArr->GetEntriesFast();
1439 112 : c->SetESDID(esdID);
1440 : // first store the hits of the cluster
1441 112 : TClonesArray* esdTOFHitArr = fESDEv->GetESDTOFHits();
1442 112 : int nh = esdTOFHitArr->GetEntriesFast();
1443 112 : Int_t tofLabels[3]={c->GetLabel(0),c->GetLabel(1),c->GetLabel(2)};
1444 112 : Int_t ind[5] = {c->GetDetInd(0), c->GetDetInd(1), c->GetDetInd(2), c->GetDetInd(3), c->GetDetInd(4) };
1445 112 : Int_t calindex = AliTOFGeometry::GetIndex(ind);
1446 : /*AliESDTOFHit* esdHit = */
1447 112 : new ( (*esdTOFHitArr)[nh] )
1448 224 : AliESDTOFHit( AliTOFGeometry::TdcBinWidth()*c->GetTDC(),
1449 112 : AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW(),
1450 112 : AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3,
1451 112 : calindex,tofLabels,c->GetL0L1Latency(),
1452 112 : c->GetDeltaBC(),esdID,c->GetZ(),c->GetR(),c->GetPhi());
1453 : //
1454 112 : clESD = new( (*esdTOFClArr)[esdID] ) AliESDTOFCluster( clID );
1455 112 : clESD->SetEvent(fESDEv);
1456 112 : clESD->SetStatus( c->GetStatus() );
1457 112 : clESD->SetESDID(esdID);
1458 : //
1459 : // register hits in the cluster
1460 112 : clESD->AddESDTOFHitIndex(nh);
1461 112 : }
1462 101 : else clESD = (AliESDTOFCluster*)esdTOFClArr->At(esdID); // cluster is aready stored in the ESD
1463 : //
1464 213 : return clESD;
1465 : //
1466 0 : }
|