]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUTrackerGlo.cxx
Updated hydro macro
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUTrackerGlo.cxx
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   TO CHECK
18   GetBz usage
19
20  */
21
22 //-------------------------------------------------------------------------
23 //  Implementation of the ITS Upgrade tracker mother class.              
24 //-------------------------------------------------------------------------
25 #include <TTree.h>
26 #include <Riostream.h> 
27 #include <TMath.h>
28 #include <TFile.h>
29
30 #include "AliITSUTrackerGlo.h"
31 #include "AliESDEvent.h"
32 #include "AliESDtrack.h"
33 #include "AliITSURecoDet.h"
34 #include "AliITSURecoSens.h"
35 #include "AliITSUReconstructor.h"
36 #include "AliITSReconstructor.h"
37 #include "AliITSUSeed.h"
38 #include "AliITSUAux.h"
39 #include "AliITSUClusterPix.h"
40 using namespace AliITSUAux;
41 using namespace TMath;
42
43 //----------------- tmp stuff -----------------
44
45 ClassImp(AliITSUTrackerGlo)
46
47 const Double_t AliITSUTrackerGlo::fgkToler =  1e-6;// tolerance for layer on-surface check
48
49
50 //_________________________________________________________________________
51 AliITSUTrackerGlo::AliITSUTrackerGlo(AliITSUReconstructor* rec)
52 :  fReconstructor(rec)
53   ,fITS(0)
54   ,fCurrESDtrack(0)
55   ,fCurrESDtrMClb(kDummyLabel)
56   ,fCurrMass(kPionMass)
57   ,fCountProlongationTrials(0)
58   ,fCountITSin(0)
59   ,fCountITSout(0)
60   ,fCountITSrefit(0)
61   ,fHypStore(100)
62   ,fNBranchesAdded(0)
63   ,fNCandidatesAdded(0)
64   ,fCurrHyp(0)
65   ,fSeedsPool("AliITSUSeed",0)
66   ,fFreeSeedsID(0)
67   ,fNFreeSeeds(0)
68   ,fLastSeedID(0)
69   ,fDefTrackConds(0)
70   ,fCurrTrackCond(0)
71   ,fTrackPhase(-1)
72   ,fClInfo(0)
73 #ifdef  _FILL_CONTROL_HISTOS_
74   ,fCHistoArr(0)
75 #endif
76 {
77   // Default constructor
78   if (rec) Init(rec);
79 }
80
81 //_________________________________________________________________________
82 AliITSUTrackerGlo::~AliITSUTrackerGlo()
83 {
84  // Default destructor
85  //  
86   delete[] fClInfo;
87   //
88 #ifdef  _FILL_CONTROL_HISTOS_
89   if (fCHistoArr) {
90     TFile* ctrOut = TFile::Open("itsuControlHistos.root","recreate");
91     ctrOut->cd();
92     AliInfo("Storing control histos");
93     fCHistoArr->Print();
94     //    ctrOut->WriteObject(fCHistoArr,"controlH","kSingleKey");
95     fCHistoArr->Write();
96     ctrOut->Close();
97     delete ctrOut;
98     fCHistoArr = 0;
99   }
100 #endif 
101   //
102 }
103
104 //_________________________________________________________________________
105 void AliITSUTrackerGlo::Init(AliITSUReconstructor* rec)
106 {
107   // init with external reconstructor
108   //
109   fITS = rec->GetITSInterface();
110   int nLr = fITS->GetNLayersActive();
111   fClInfo = new Int_t[nLr<<1];
112   //
113   fSeedsPool.ExpandCreateFast(1000); // RS TOCHECK
114   fFreeSeedsID.Set(1000);
115   //
116 }
117
118 //_________________________________________________________________________
119 void AliITSUTrackerGlo::CreateDefaultTrackCond()
120 {
121   // creates default tracking conditions to be used when recoparam does not provide them
122   int nLr = fITS->GetNLayersActive();
123   fClInfo = new Int_t[nLr<<1];
124   //
125   AliITSUTrackCond* cond = new AliITSUTrackCond();
126   //
127   cond->SetNLayers(fITS->GetNLayersActive());
128   cond->AddNewCondition(nLr);
129   cond->AddGroupPattern( 0xffff ); // require all layers hit
130   //
131   fDefTrackConds.AddLast(cond);
132   //
133   AliInfo("Created conditions: ");
134   for (int i=0;i<fDefTrackConds.GetEntriesFast();i++) fDefTrackConds[i]->Print();
135   //
136 }
137
138
139 //_________________________________________________________________________
140 Int_t AliITSUTrackerGlo::Clusters2Tracks(AliESDEvent *esdEv)
141 {
142   //
143   SetTrackingPhase(kClus2Tracks);
144   //
145 #ifdef  _FILL_CONTROL_HISTOS_
146   if (!fCHistoArr) BookControlHistos();
147 #endif
148
149   TObjArray *trackConds = 0;
150   //
151   fCountProlongationTrials = 0;
152   fCountITSin = 0;
153   fCountITSout = 0;
154   fCountITSrefit = 0;
155   //
156   ResetSeedsPool();
157   int nTrESD = esdEv->GetNumberOfTracks();
158   AliInfo(Form("Will try to find prolongations for %d tracks",nTrESD));
159   int nTrackCond = AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();
160   if (nTrackCond<1) {
161     if (!fDefTrackConds.GetEntriesFast()) {
162       AliInfo("No tracking conditions found in recoparams, creating default one requesting all layers hit");
163       CreateDefaultTrackCond();
164     }
165     trackConds = &fDefTrackConds;
166     nTrackCond = trackConds->GetEntriesFast();
167   } 
168   else trackConds = AliITSUReconstructor::GetRecoParam()->GetTrackingConditions();
169   //
170   static Bool_t first = kTRUE;
171   if (first) {
172     AliITSUReconstructor::GetRecoParam()->Print();
173     first = kFALSE;
174   }
175   fHypStore.Delete();
176   if (fHypStore.GetSize()<nTrESD) fHypStore.Expand(nTrESD+100);
177   //
178   fITS->ProcessClusters();
179   //
180   FlagSplitClusters(); // tmp RS
181   //
182   for (int icnd=0;icnd<nTrackCond;icnd++) {
183     fCurrTrackCond = (AliITSUTrackCond*)trackConds->UncheckedAt(icnd);
184     // select ESD tracks to propagate
185     for (int itr=0;itr<nTrESD;itr++) {
186       fCurrESDtrack = esdEv->GetTrack(itr);
187       fCurrESDtrMClb = fCurrESDtrack->GetLabel();
188       //
189       if (!NeedToProlong(fCurrESDtrack)) continue;  // are we interested in this track?
190       AliDebug(+1,Form("Processing track %d | M=%.3f Pt=%.3f | MCLabel: %d",itr,fCurrESDtrack->GetMass(kTRUE),fCurrESDtrack->Pt(),fCurrESDtrMClb));//RS
191       FindTrack(fCurrESDtrack, itr);
192     }   
193     //
194     if (AliDebugLevelClass()>+2) {
195       AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
196       fHypStore.Print();
197     }
198     FinalizeHypotheses();
199   }
200   //
201   AliInfo(Form("%d ITSin for %d tried TPC seeds out of %d ESD tracks\n",fCountITSin,fCountProlongationTrials,nTrESD));
202   return 0;
203 }
204
205 //_________________________________________________________________________
206 Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
207 {
208   //
209   // Do outward fits in ITS
210   //
211   SetTrackingPhase(kPropBack);
212   int nTrESD = esdEv->GetNumberOfTracks();
213   AliDebug(1,Form("Will propagate back %d tracks",nTrESD));
214   //
215   double bz0 = GetBz();
216   Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
217   AliITSUTrackHyp dummyTr,*currTr=0;
218   const double kWatchStep=10.; // for larger steps watch arc vs segment difference
219   Double_t times[AliPID::kSPECIES];
220   //
221   //
222   for (int itr=0;itr<nTrESD;itr++) {
223     fCurrESDtrack = esdEv->GetTrack(itr);
224     fCurrESDtrMClb = fCurrESDtrack->GetLabel();
225     // Start time integral and add distance from current position to vertex 
226     if (fCurrESDtrack->IsOn(AliESDtrack::kITSout)) continue; // already done
227     //
228     fCurrESDtrack->GetXYZ(xyzTrk); 
229     Double_t dst = 0.;     // set initial track lenght, tof
230     {
231       double dxs = xyzTrk[0] - xyzVtx[0];
232       double dys = xyzTrk[1] - xyzVtx[1];
233       double dzs = xyzTrk[2] - xyzVtx[2];
234       // RS: for large segment steps d use approximation of cicrular arc s by
235       // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
236       // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
237       // Hence s^2/d^2 = (1+1/6 p^2)^2
238       dst = dxs*dxs + dys*dys;
239       if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
240         double crv = Abs(fCurrESDtrack->GetC(bz0));
241         double fcarc = 1.+crv*crv*dst/6.;
242         dst *= fcarc*fcarc;
243       }
244       dst += dzs*dzs;
245       dst = Sqrt(dst); 
246     }
247     //    
248     fCurrESDtrack->SetStatus(AliESDtrack::kTIME);
249     //
250     if (!fCurrESDtrack->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
251       dummyTr.AliExternalTrackParam::operator=(*fCurrESDtrack);
252       dummyTr.StartTimeIntegral();
253       dummyTr.AddTimeStep(dst);
254       dummyTr.GetIntegratedTimes(times); 
255       fCurrESDtrack->SetIntegratedTimes(times);
256       fCurrESDtrack->SetIntegratedLength(dummyTr.GetIntegratedLength());
257       continue;
258     }
259     //
260     currTr = GetTrackHyp(itr);
261     currTr->StartTimeIntegral();
262     currTr->AddTimeStep(dst);
263     //    printf("Before resetCov: "); currTr->AliExternalTrackParam::Print();
264     currTr->ResetCovariance(10000);
265     if (RefitTrack(currTr,fITS->GetRMax())) { // propagate to exit from the ITS/TPC screen
266       UpdateESDTrack(currTr,AliESDtrack::kITSout);
267       fCountITSout++;
268     }
269     else {
270       AliDebug(2,Form("Refit Failed for track %d | ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
271       //currTr->AliExternalTrackParam::Print();
272       //currTr->GetWinner()->Print();
273     }
274     //
275   }
276   //
277   AliInfo(Form("%d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
278                fCountITSout,fCountITSin,fCountProlongationTrials,nTrESD));
279   //
280   return 0;
281 }
282
283 //_________________________________________________________________________
284 Int_t AliITSUTrackerGlo::RefitInward(AliESDEvent *esdEv)
285 {
286   //
287   // refit the tracks inward, using current cov.matrix
288   //
289   SetTrackingPhase(kRefitInw);
290   Int_t nTrESD = esdEv->GetNumberOfTracks();
291   //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
292
293   AliDebug(1,Form("Will refit inward %d tracks",nTrESD));
294   AliITSUTrackHyp *currTr=0;
295   //
296   for (int itr=0;itr<nTrESD;itr++) {
297     fCurrESDtrack = esdEv->GetTrack(itr);
298     fCurrESDtrMClb = fCurrESDtrack->GetLabel();
299     // Start time integral and add distance from current position to vertex 
300     UInt_t trStat = fCurrESDtrack->GetStatus();
301     if ( !(trStat & AliESDtrack::kITSout) ) continue;
302     if (   trStat & AliESDtrack::kITSrefit ) continue; // already done
303     if (  (trStat & AliESDtrack::kTPCout) && !(trStat & AliESDtrack::kTPCrefit) ) continue;
304     //
305     currTr = GetTrackHyp(itr);
306     currTr->AliExternalTrackParam::operator=(*fCurrESDtrack);  // fetch current ESDtrack kinematics
307     if (RefitTrack(currTr,fITS->GetRMin())) { // propagate up to inside radius of the beam pipe
308       UpdateESDTrack(currTr,AliESDtrack::kITSrefit);
309       fCountITSrefit++;
310     }
311     else {
312       AliDebug(2,Form("Refit Failed for track %d |ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
313     }
314   }    
315   //
316   AliInfo(Form("%d ITSrefit in %d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
317                fCountITSrefit,fCountITSout,fCountITSin,fCountProlongationTrials,nTrESD));
318   //
319   //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
320   return 0;
321 }
322
323 //_________________________________________________________________________
324 Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
325 {
326   // read from tree (if pointer provided) or directly from the ITS reco interface
327   //
328   return fReconstructor->LoadClusters(treeRP);
329
330
331 //_________________________________________________________________________
332 void AliITSUTrackerGlo::UnloadClusters()
333 {
334   //
335   // To be implemented 
336   //
337   
338   Info("UnloadClusters","To be implemented");
339
340 //_________________________________________________________________________
341 AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
342 {
343   //
344   // To be implemented 
345   //
346   
347   Info("GetCluster","To be implemented");
348   return 0x0;
349
350
351 //_________________________________________________________________________
352 Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr)
353 {
354   // do we need to match this track to ITS?
355   //
356   static double bz = GetBz();
357   if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
358       esdTr->IsOn(AliESDtrack::kTPCout) ||
359       esdTr->IsOn(AliESDtrack::kITSin)  ||
360       esdTr->GetKinkIndex(0)>0) return kFALSE;
361   //
362   if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
363   //
364   float dtz[2];
365   esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz); 
366   // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
367   if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation()) 
368        && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
369            Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
370   //
371   return kTRUE;
372 }
373
374 //_________________________________________________________________________
375 void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
376 {
377   // find prolongaion candidates finding for single seed
378   //
379   AliITSUSeed seedUC;  // copy of the seed from the upper layer
380   AliITSUSeed seedT;   // transient seed between the seedUC and new prolongation hypothesis
381   //
382   if (!InitHypothesis(esdTr,esdID)) return;  // initialize prolongations hypotheses tree
383   //
384   AliITSURecoSens *hitSens[AliITSURecoSens::kNNeighbors+1];
385   //
386   TObjArray clArr; // container for transfer of clusters matching to seed
387   //
388   int nLrActive = fITS->GetNLayersActive();
389   for (int ila=nLrActive;ila--;) {
390     int ilaUp = ila+1;                         // prolong seeds from layer above
391     //
392     // for the outermost layer the seed is created from the ESD track
393     int nSeedsUp = (ilaUp==nLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
394     int maxNBranches   = fCurrTrackCond->GetMaxBranches(ila);
395     int maxNCandidates = fCurrTrackCond->GetMaxCandidates(ila);
396     //
397     for (int isd=0;isd<nSeedsUp;isd++) {
398       AliITSUSeed* seedU;
399       if (ilaUp==nLrActive) {
400         seedU = 0;
401         seedUC.InitFromESDTrack(esdTr);
402       }
403       else {
404         seedU = fCurrHyp->GetSeed(ilaUp,isd);  // seed on prev.active layer to prolong  
405         seedUC = *seedU;                       // its copy will be prolonged
406         seedUC.SetParent(seedU);        
407       }
408       seedUC.ResetFMatrix();                    // reset the matrix for propagation to next layer
409       // go till next active layer
410       AliDebug(2,Form("working on Lr:%d Seed:%d of %d for esdID=%d (MClb:%d) | pT=%.3f",ila,isd,nSeedsUp,esdID,fCurrESDtrMClb,seedUC.Pt()));
411       if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) {
412         //
413         AliDebug(2,Form("Transport failed | esdID=%d (MClb:%d)",esdID,fCurrESDtrMClb));
414         // Check if the seed satisfies to track definition
415         if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
416         continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
417       }
418       AliITSURecoLayer* lrA = fITS->GetLayerActive(ila);
419       if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
420         if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
421         continue;
422       }
423       int nsens = lrA->FindSensors(&fTrImpData[kTrPhi0], hitSens);  // find detectors which may be hit by the track
424       AliDebug(2,Form("Will check %d sensors on lr:%d | esdID=%d (MClb:%d)",nsens,ila,esdID,fCurrESDtrMClb));
425       //
426       seedUC.SetLr(ila);
427       //
428       double bz = GetBz();
429       for (int isn=nsens;isn--;) {
430         AliITSURecoSens* sens = hitSens[isn];
431         int ncl = sens->GetNClusters();
432         if (!ncl) continue;
433         seedT = seedUC;
434         //
435         // We need to propagate the seed to sensor on lrA staying in the frame of the sensor from prev.layer,
436         // since the transport matrix should be defined in this frame.
437         double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
438         if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) {
439           if (AliDebugLevelClass()>2) {
440             printf("Failed on GetTrackingXAtXAlpha: X=%.4f alp=%.3f\n",sens->GetXTF(),sens->GetPhiTF());
441             seedT.Print("etp");
442           }
443           continue;
444         }
445         if (xs<seedT.GetX()) {
446           if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
447         }
448         else { // some low precision tracks may hit the sensor plane outside of the layer radius
449           if (AliDebugLevelClass()>2) {
450             if (!seedT.ContainsFake()) {
451               printf("WRONG PROP on L%d, sens%d of %d: %.4f > %.4f\n",ila,isn,nsens,xs,seedT.GetX());
452               sens->Print();
453               seedT.Print("etp");           
454             }
455           }
456           if (!seedT.PropagateParamOnlyTo(xs,bz)) continue;
457         }
458         //      if (!seedT.PropagateToX(xs,bz)) continue;
459         //      if (!seedT.Rotate(sens->GetPhiTF())) continue;
460         if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
461         //
462         int clID0 = sens->GetFirstClusterId();
463         for (int icl=ncl;icl--;) {
464           //      AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
465           int res = CheckCluster(&seedT,ila,clID0+icl);
466           //      AliLog::SetClassDebugLevel("AliITSUTrackerGlo", 0);
467           //
468           if (res==kStopSearchOnSensor) break;     // stop looking on this sensor
469           if (res==kClusterNotMatching) continue;  // cluster does not match
470           // cluster is matching and it was added to the hypotheses tree
471         }
472       }
473       // cluster search is done. Do we need to have a version of this seed skipping current layer
474       if (!NeedToKill(&seedUC,kMissingCluster)) {
475         AliITSUSeed* seedSkp = NewSeedFromPool(&seedUC);
476         double penalty = -AliITSUReconstructor::GetRecoParam()->GetMissPenalty(ila);
477         // to do: make penalty to account for probability to miss the cluster for good reason
478         seedSkp->SetChi2Cl(penalty);
479         AddProlongationHypothesis(seedSkp,ila);      
480       }
481       // transfer the new branches of the seed to the hypothesis container
482       if (fNBranchesAdded) ValidateAllowedBranches(maxNBranches);
483       //
484     }
485     if (fNCandidatesAdded) ValidateAllowedCandidates(ila,maxNCandidates);
486     //    ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
487     if (AliDebugLevelClass()>2) { //RS
488       printf(">>> All hypotheses on lr %d: \n",ila);
489       for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();}
490     }
491     //
492     /*
493     for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
494       if (ila!=0) continue;
495       double vecL[5] = {0};
496       double matL[15] = {0};
497       AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
498       while(sp->GetParent()) {
499         sp->Smooth(vecL,matL);
500         if (sp->GetLayerID()>=fITS->GetNLayersActive()-1) break;
501         sp = (AliITSUSeed*)sp->GetParent();
502       }
503       */
504   }
505   //
506   SaveCurrentTrackHypotheses();
507   //
508 }
509
510 //_________________________________________________________________________
511 Bool_t AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
512 {
513   // init prolongaion candidates finding for single seed
514   fCurrHyp = GetTrackHyp(esdID);
515   if (fCurrHyp) return kTRUE;
516   //
517   fCountProlongationTrials++;
518   //
519   fCurrMass = esdTr->GetMass();
520   if (fCurrMass<kPionMass*0.9) fCurrMass = kPionMass; // don't trust to mu, e identification from TPCin
521   //
522   fCurrHyp = new AliITSUTrackHyp(fITS->GetNLayersActive());
523   fCurrHyp->SetESDTrack(esdTr);
524   fCurrHyp->SetUniqueID(esdID);
525   fCurrHyp->SetMass(fCurrMass);
526   SetTrackHyp(fCurrHyp,esdID);
527   //
528   return kTRUE;
529   // TO DO
530 }
531
532 //_________________________________________________________________________
533 Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo)
534 {
535   // transport seed from layerFrom to the entrance (along the direction of the track) of layerTo
536   //  
537   //
538   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
539   //
540   int dir = lTo > lFrom ? 1:-1;
541   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
542   Bool_t checkFirst = kTRUE;
543   while(lFrom!=lTo) {
544     if (lrFr) {
545       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) {
546         //printf("FailHere0 Dir=%d\n",dir);
547         //seed->Print("etp");
548         return kFALSE; // go till the end of current layer
549       }
550       checkFirst = kFALSE;
551     }
552     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
553     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
554     //
555     // go the entrance of the layer, assuming no materials in between
556     double xToGo = lrTo->GetR(-dir);
557     //    double xts = xToGo;
558     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
559       //printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
560       //seed->Print("etp");
561       //
562       
563       return kFALSE;
564     }
565     AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
566     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
567       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
568       //seed->Print("etp");
569       return kFALSE;
570     }
571     lrFr = lrTo;
572   }
573   return kTRUE;
574   //
575 }
576
577 //_________________________________________________________________________
578 Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo)
579 {
580   // transport track from layerFrom to the entrance of layerTo
581   //  
582   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
583   //
584   int dir = lTo > lFrom ? 1:-1;
585   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
586   Bool_t checkFirst = kTRUE;
587   while(lFrom!=lTo) {
588     if (lrFr) {
589       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
590       checkFirst = kFALSE;
591     }
592     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
593     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
594     //
595     // go the entrance of the layer, assuming no materials in between
596     double xToGo = lrTo->GetR(-dir);
597     //    double xts = xToGo;
598     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
599       //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
600       //      seed->Print("etp");
601       return kFALSE;
602     }
603     AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
604     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
605       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
606       //seed->Print("etp");
607       return kFALSE;
608     }
609     lrFr = lrTo;
610   }
611   return kTRUE;
612   //
613 }
614
615 //_________________________________________________________________________
616 Bool_t AliITSUTrackerGlo::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
617 {
618   // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
619   //  
620   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
621   //
622   int dir = lTo > lFrom ? 1:-1;
623   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
624   Bool_t checkFirst = kTRUE;
625   while(lFrom!=lTo) {
626     if (lrFr) {
627       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
628       checkFirst = kFALSE;
629     }
630     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
631     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
632     //
633     // go the entrance of the layer, assuming no materials in between
634     double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
635     //
636     //    double xts = xToGo;
637     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
638       //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
639       //      seed->Print("etp");
640       return kFALSE;
641     }
642     if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
643     //
644     AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
645     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
646       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
647       //seed->Print("etp");
648       return kFALSE;
649     }
650     lrFr = lrTo;
651   }
652   return kTRUE;
653   //
654 }
655
656 //_________________________________________________________________________
657 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
658 {
659   // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
660   // If check is requested, do this only provided the track has not exited the layer already
661   double xToGo = lr->GetR(dir);
662   if (check) { // do we need to track till the surface of the current layer ?
663     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
664     AliDebug(3,Form(" dir:%d Cur: %e Tgt: %e",dir,Sqrt(curR2),xToGo));
665     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
666     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
667   }
668   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
669   // go via layer to its boundary, applying material correction.
670   if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
671   //
672   return kTRUE;
673   //
674 }
675
676 //_________________________________________________________________________
677 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
678 {
679   // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
680   // If check is requested, do this only provided the track has not exited the layer already
681   double xToGo = lr->GetR(dir);
682   if (check) { // do we need to track till the surface of the current layer ?
683     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
684     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
685     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
686   }
687   AliDebug(2,Form(" dir=%d : from R=%.4f -> R=%.4f",dir,Sqrt(seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY()), xToGo));
688   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
689   // go via layer to its boundary, applying material correction.
690   if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
691   //
692   return kTRUE;
693   //
694 }
695
696
697 //_________________________________________________________________________
698 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
699 {
700   // go to the entrance of lr in direction dir, w/o applying material corrections.
701   // If check is requested, do this only provided the track did not reach the layer already
702   double xToGo = lr->GetR(-dir);
703   if (check) { // do we need to track till the surface of the current layer ?
704     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
705     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
706     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
707   }
708   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
709   // go via layer to its boundary, applying material correction.
710   if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
711   return kTRUE;
712   //
713 }
714
715 //_________________________________________________________________________
716 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
717 {
718   // go to the entrance of lr in direction dir, w/o applying material corrections.
719   // If check is requested, do this only provided the track did not reach the layer already
720   double xToGo = lr->GetR(-dir);
721   if (check) { // do we need to track till the surface of the current layer ?
722     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
723     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
724     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
725   }
726   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
727   // go via layer to its boundary, applying material correction.
728   if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
729   return kTRUE;
730   //
731 }
732
733 //_________________________________________________________________________
734 Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
735 {
736   // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
737   // as well as some aux info
738   double bz = GetBz();
739   AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
740   seed->GetXYZ(&fTrImpData[kTrXIn]);    // lab position at the entrance from above
741   static AliExternalTrackParam sc;   // seed copy for manipulations
742   sc = *seed;
743   //
744   fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
745   if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
746   double dr  = lrA->GetDR();                              // approximate X dist at the inner radius
747   if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
748     // special case: track does not reach inner radius, might be tangential
749     double r = sc.GetD(0,0,bz);
750     double x;
751     if (!sc.GetXatLabR(r,x,bz,-1)) {
752       sc.Print();
753       AliFatal(Form("This should not happen: r=%f",r));
754     }
755     dr = Abs(sc.GetX() - x);
756     if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) {
757       sc.Print();
758       AliFatal(Form("This should not happen: x=%f",x));
759     }
760   }
761   //
762   fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
763   double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
764   double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
765   sgy = Sqrt(sgy)*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadY();
766   sgz = Sqrt(sgz)*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadZ();
767   fTrImpData[kTrPhi0] = 0.5*(fTrImpData[kTrPhiOut]+fTrImpData[kTrPhiIn]);
768   fTrImpData[kTrZ0]   = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
769   fTrImpData[kTrDPhi] = 0.5*Abs(fTrImpData[kTrPhiOut]-fTrImpData[kTrPhiIn]) + sgy/lrA->GetR();
770   fTrImpData[kTrDZ]   = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn])   + sgz;
771   //  
772   return kTRUE;
773 }
774
775 //________________________________________
776 void AliITSUTrackerGlo::ResetSeedsPool()
777 {
778   // mark all seeds in the pool as unused
779   AliDebug(-1,Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
780   fNFreeSeeds = 0;
781   fSeedsPool.Clear(); // seeds don't allocate memory
782 }
783
784
785 //________________________________________
786 void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd) 
787 {
788   // account that this seed is "deleted" 
789   int id = sd->GetPoolID();
790   if (id<0) {
791     AliError(Form("Freeing of seed %p NOT from the pool is requested",sd)); 
792     return;
793   }
794   //  AliInfo(Form("%d %p",id, seed));
795   fSeedsPool.RemoveAt(id);
796   if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
797   fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
798 }
799
800 //_________________________________________________________________________
801 Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID) 
802 {
803   // Check if the cluster (in tracking frame!) is matching to track. 
804   // The track must be already propagated to sensor tracking frame.
805   // Returns:  kStopSearchOnSensor if the search on given sensor should be stopped, 
806   //           kClusterMatching    if the cluster is matching
807   //           kClusterMatching    otherwise
808   //
809   const double kTolerX = 5e-4;
810   AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
811   //
812   Bool_t goodCl = kFALSE;
813   int currLabel = Abs(fCurrESDtrMClb);
814   //
815   if (cl->GetLabel(0)>=0) {
816     for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
817   }
818   Bool_t goodSeed = !track->ContainsFake();
819   //
820   if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
821     if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
822       if (goodCl&&goodSeed && AliDebugLevelClass()>2 ) {
823         AliDebug(2,Form("Lost good cl on L:%d failed propagation. |ESDtrack#%d (MClb:%d)",lr,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
824         track->Print("etp");
825         cl->Print();
826       }
827       return kStopSearchOnSensor; // propagation failed, seedT is intact
828     }
829   }
830   double dy = cl->GetY()-track->GetY();
831   double dz = cl->GetZ()-track->GetZ();
832   //
833 #ifdef  _FILL_CONTROL_HISTOS_
834   int hcOffs = fTrackPhase*kHistosPhase + lr;
835   double htrPt=-1;
836   if (goodCl&&goodSeed && fCHistoArr /*&& track->GetChi2Penalty()<1e-5*/) {
837     htrPt = track->Pt();
838     ((TH2*)fCHistoArr->At(kHResY+hcOffs))->Fill(htrPt,dy);
839     ((TH2*)fCHistoArr->At(kHResZ+hcOffs))->Fill(htrPt,dz);
840     double errY = track->GetSigmaY2();
841     double errZ = track->GetSigmaZ2();
842     if (errY>0) ((TH2*)fCHistoArr->At(kHResYP+hcOffs))->Fill(htrPt,dy/Sqrt(errY));
843     if (errZ>0) ((TH2*)fCHistoArr->At(kHResZP+hcOffs))->Fill(htrPt,dz/Sqrt(errZ));
844   }
845 #endif  
846   //
847   double dy2 = dy*dy;
848   double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
849     AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadY()*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadY(); // RS TOOPTIMIZE
850   if (dy2>tol2) {                          // the clusters are sorted in Z(col) then in Y(row). 
851     if (goodCl&&goodSeed && AliDebugLevelClass()>2) {    
852       AliDebug(2,Form("Lost good cl: dy2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",dy2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
853       track->Print("etp");
854       cl->Print();
855       PrintSeedClusters(track);
856     }
857     if (dy>0) return kStopSearchOnSensor;  // No chance that other cluster of this sensor will match (all Y's will be even larger)
858     else      return kClusterNotMatching;   // Other clusters may match
859   }
860   double dz2 = dz*dz;
861   tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
862     AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadZ()*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadZ(); // RS TOOPTIMIZE
863   if (dz2>tol2) {
864     if (goodCl&&goodSeed && AliDebugLevelClass()>2) {
865       AliDebug(2,Form("Lost good cl on L:%d : dz2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",lr,dz2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
866       track->Print("etp");
867       cl->Print();
868       PrintSeedClusters(track);
869     }
870     return kClusterNotMatching; // Other clusters may match
871   }
872   //
873   // check chi2
874   Double_t p[2]={cl->GetY(), cl->GetZ()};
875   Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
876   double chi2 = track->GetPredictedChi2(p,cov);
877   //
878 #ifdef  _FILL_CONTROL_HISTOS_
879   if (htrPt>0) {
880     ((TH2*)fCHistoArr->At(kHChi2Cl+hcOffs))->Fill(htrPt,chi2);
881   }
882 #endif
883   //
884   if (chi2>AliITSUReconstructor::GetRecoParam()->GetMaxTr2ClChi2(lr)) {
885     if (goodCl&&goodSeed && AliDebugLevelClass()>2) {
886       AliDebug(2,Form("Lost good cl on L:%d : Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
887                       lr,chi2,AliITSUReconstructor::GetRecoParam()->GetMaxTr2ClChi2(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
888       track->Print("etp");
889       cl->Print("");
890       PrintSeedClusters(track);
891     }
892     return kClusterNotMatching;
893   }
894   //
895   track = NewSeedFromPool(track);  // input track will be reused, use its clone for updates
896   if (!track->Update()) {
897     if (goodCl&&goodSeed && AliDebugLevelClass()>2) {
898       AliDebug(2,Form("Lost good cluster on L:%d : failed to update",lr));
899       track->Print("etp");
900       cl->Print("");
901       PrintSeedClusters(track);
902     }
903     MarkSeedFree(track);
904     return kClusterNotMatching;
905   }
906   track->SetChi2Cl(chi2);
907   track->SetLrClusterID(lr,clID);
908   cl->IncreaseClusterUsage();
909   //
910   track->SetFake(!goodCl);
911   //
912   AliDebug(2,Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)| ",
913                goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
914   //
915   AddSeedBranch(track);
916   //
917   return kClusterMatching;
918 }
919
920 //_________________________________________________________________________
921 Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
922 {
923   // check if the seed should not be discarded
924   const UShort_t kMask = 0xffff;
925   if (flag==kMissingCluster) {
926     int lastChecked = seed->GetLayerID();
927     UShort_t patt = seed->GetHitsPattern();
928     if (lastChecked) patt |= ~(kMask<<lastChecked); // not all layers were checked, complete unchecked once by potential hits
929     Bool_t seedOK = fCurrTrackCond->CheckPattern(patt);
930     return !seedOK;
931   }
932   return kTRUE;
933 }
934
935 //______________________________________________________________________________
936 Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr) 
937 {
938   // propagate seed to given x applying material correction if requested
939   const Double_t kEpsilon = 1e-5;
940   Double_t xpos     = seed->GetX();
941   Int_t dir         = (xpos<xToGo) ? 1:-1;
942   Double_t xyz0[3],xyz1[3],param[7];
943   //
944   Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
945   if (matCorr || updTime) seed->GetXYZ(xyz1);   //starting global position
946   while ( (xToGo-xpos)*dir > kEpsilon){
947     Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
948     Double_t x    = xpos+step;
949     Double_t bz=GetBz();   // getting the local Bz
950     if (!seed->PropagateToX(x,bz)) return kFALSE;
951     double ds = 0;
952     if (matCorr || updTime) {
953       xyz0[0]=xyz1[0]; // global pos at the beginning of step
954       xyz0[1]=xyz1[1];
955       xyz0[2]=xyz1[2];
956       seed->GetXYZ(xyz1);    //  // global pos at the end of step
957       if (matCorr) {
958         MeanMaterialBudget(xyz0,xyz1,param);    
959         Double_t xrho=param[0]*param[4], xx0=param[1];
960         if (dir>0) xrho = -xrho; // outward should be negative
961         if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) {AliInfo("Fail2"); seed->Print("etp"); return kFALSE;}
962         ds = param[4];
963       }
964        else { // matCorr is not requested but time integral is
965         double d0 = xyz1[0]-xyz0[0];
966         double d1 = xyz1[1]-xyz0[1];
967         double d2 = xyz1[2]-xyz0[2];    
968         ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
969       }     
970     }
971     if (updTime) seed->AddTimeStep(ds);
972     xpos = seed->GetX();
973   }
974   return kTRUE;
975 }
976
977 //______________________________________________________________________________
978 Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr) 
979 {
980   // propagate seed to given x applying material correction if requested
981   const Double_t kEpsilon = 1e-5;
982   Double_t xpos     = seed->GetX();
983   Int_t dir         = (xpos<xToGo) ? 1:-1;
984   Double_t xyz0[3],xyz1[3],param[7];
985   //
986   if (AliDebugLevelClass()>1) {
987     AliDebug(2,Form("Before Propagate to X=%f with M=%.3f MaxStep=%.4f MatCorr=%d",xToGo,mass,maxStep,matCorr));
988     seed->AliExternalTrackParam::Print();
989   }
990   Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
991   if (matCorr || updTime) seed->GetXYZ(xyz1);   //starting global position
992   while ( (xToGo-xpos)*dir > kEpsilon){
993     Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
994     Double_t x    = xpos+step;
995     Double_t bz=GetBz();   // getting the local Bz
996     if (!seed->PropagateTo(x,bz))  return kFALSE;
997     double ds = 0;
998     if (matCorr || updTime) {
999       xyz0[0]=xyz1[0]; // global pos at the beginning of step
1000       xyz0[1]=xyz1[1];
1001       xyz0[2]=xyz1[2];
1002       seed->GetXYZ(xyz1);    //  // global pos at the end of step
1003       //
1004       if (matCorr) {
1005         MeanMaterialBudget(xyz0,xyz1,param);    
1006         Double_t xrho=param[0]*param[4], xx0=param[1];
1007         if (dir>0) xrho = -xrho; // outward should be negative
1008         if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
1009         ds = param[4];
1010       }
1011       else { // matCorr is not requested but time integral is
1012         double d0 = xyz1[0]-xyz0[0];
1013         double d1 = xyz1[1]-xyz0[1];
1014         double d2 = xyz1[2]-xyz0[2];    
1015         ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1016       }
1017     }
1018     if (updTime) seed->AddTimeStep(ds);
1019     //
1020     xpos = seed->GetX();
1021   }
1022   if (AliDebugLevelClass()>1) {
1023     AliDebug(2,Form("After Propagate to X=%f",xToGo));
1024     seed->AliExternalTrackParam::Print();
1025   }
1026   return kTRUE;
1027 }
1028
1029 //______________________________________________________________________________
1030 void AliITSUTrackerGlo::SaveCurrentTrackHypotheses()
1031 {
1032   // RS: shall we clean up killed seeds?
1033   ((TObjArray*)fCurrHyp->GetLayerSeeds(0))->Sort();
1034   fCurrHyp = 0;
1035   // TODO
1036   
1037 }
1038
1039 //______________________________________________________________________________
1040 void AliITSUTrackerGlo::FinalizeHypotheses()
1041 {
1042   // select winner for each hypothesis, remove cl. sharing conflicts
1043   AliDebug(2,"TODO");
1044   //
1045   int nh = fHypStore.GetEntriesFast();
1046   for (int ih=0;ih<nh;ih++) {
1047     AliITSUTrackHyp* hyp = (AliITSUTrackHyp*) fHypStore.UncheckedAt(ih); 
1048     if (!hyp || !hyp->DefineWinner()) continue; // TODO
1049     if (hyp->GetESDTrack()->IsOn(AliESDtrack::kITSin)) continue;
1050     fCountITSin++;
1051     CookMCLabel(hyp);
1052     UpdateESDTrack(hyp,AliESDtrack::kITSin);
1053   }
1054
1055 }
1056
1057 //______________________________________________________________________________
1058 void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp,Int_t flag)
1059 {
1060   // update ESD track with current best hypothesis
1061   AliESDtrack* esdTr = hyp->GetESDTrack();
1062   if (!esdTr) return;
1063   AliITSUSeed* win = hyp->GetWinner();
1064   if (!win) return;
1065   switch (flag) {
1066   case AliESDtrack::kITSin: 
1067     esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1068     // TODO: set cluster info
1069     break;
1070     //
1071   case AliESDtrack::kITSout: 
1072     esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1073     // TODO: avoid setting friend
1074     break;
1075     //
1076   case AliESDtrack::kITSrefit: 
1077     esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1078     // TODO: avoid setting cluster info
1079     break;
1080   default:
1081     AliFatal(Form("Unknown flag %d",flag));
1082   }
1083   //
1084   // transfer module indices
1085   // TODO
1086 }
1087
1088 //______________________________________________________________________________
1089 Bool_t AliITSUTrackerGlo::RefitTrack(AliITSUTrackHyp* trc, Double_t rDest, Bool_t stopAtLastCl)
1090 {
1091   // refit track till radius rDest
1092   AliITSUTrackHyp tmpTr;
1093   //
1094   double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
1095   int dir,lrStart,lrStop;
1096   //
1097   dir = rCurr<rDest ? 1 : -1;
1098   lrStart = fITS->FindFirstLayerID(rCurr,dir);
1099   lrStop  = fITS->FindLastLayerID(rDest,dir);
1100   if (AliDebugLevelClass()>2) {
1101     printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
1102     printf("Before refit: "); trc->AliExternalTrackParam::Print();
1103   }
1104   if (lrStop<0 || lrStart<0) AliFatal(Form("Failed to find start(%d) or last(%d) layers. Track from %.3f to %.3f",lrStart,lrStop,rCurr,rDest));
1105   //
1106   int ncl = trc->FetchClusterInfo(fClInfo);
1107   fCurrMass = trc->GetMass();
1108   tmpTr.AliKalmanTrack::operator=(*trc);
1109   tmpTr.SetChi2(0);
1110   int iclLr[2],nclLr,clCount=0;
1111   //
1112   for (int ilr=lrStart;ilr!=lrStop;ilr+=dir) {
1113     AliITSURecoLayer* lr = fITS->GetLayer(ilr);
1114     if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
1115     int ilrA2,ilrA = lr->GetActiveID();
1116     // passive layer or active w/o hits will be traversed on the way to next cluster
1117     if (!lr->IsActive() || fClInfo[ilrA2=(ilrA<<1)]<0) continue; 
1118     //
1119     // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
1120     nclLr=0;
1121     if (dir>0) { // clusters are stored in increasing radius order
1122       iclLr[nclLr++]=fClInfo[ilrA2++];
1123       if (fClInfo[ilrA2]>=0) iclLr[nclLr++]=fClInfo[ilrA2];
1124     }
1125     else {
1126       if ( fClInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=fClInfo[ilrA2+1];
1127       iclLr[nclLr++]=fClInfo[ilrA2];
1128     }
1129     //
1130     Bool_t transportedToLayer = kFALSE;
1131     for (int icl=0;icl<nclLr;icl++) {
1132       AliITSUClusterPix* clus =  (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
1133       AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
1134       if (!tmpTr.Rotate(sens->GetPhiTF())) {
1135         AliDebug(2,Form("Failed on rotate to %f | ESDtrack#%d (MClb:%d)",sens->GetPhiTF(),fCurrESDtrack->GetID(),fCurrESDtrMClb));
1136         return kFALSE;
1137       }
1138       //
1139       double xClus = sens->GetXTF()+clus->GetX();
1140       if (!transportedToLayer) {
1141         if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
1142           AliDebug(2,Form("Failed to transport %d -> %d | ESDtrack#%d (MClb:%d)\n",lrStart,ilr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1143           //tmpTr.AliExternalTrackParam::Print(); trc->GetWinner()->Print("etp");
1144           return kFALSE; // go to the entrance to the layer
1145         }
1146         lrStart = ilr;
1147         transportedToLayer = kTRUE;
1148       }
1149       //
1150       if (AliDebugLevelClass()>1) {
1151         AliDebug(2,Form("Propagate to cl:%d on lr %d Need to go %.4f -> %.4f",icl,ilrA,tmpTr.GetX(),xClus));
1152         //      printf("Before: "); tmpTr.AliExternalTrackParam::Print();
1153       }
1154       //
1155       if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) {
1156         AliDebug(2,Form("Failed on propagate to %f (dir=%d) | ESDtrack#%d (MClb:%d)",xClus,dir,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1157         //tmpTr.AliExternalTrackParam::Print(""); trc->GetWinner()->Print("etp");
1158         return kFALSE;
1159       }
1160       //
1161 #ifdef  _FILL_CONTROL_HISTOS_
1162       int hcOffs = fTrackPhase*kHistosPhase + ilrA;
1163       double htrPt=-1;
1164       if (fCHistoArr && trc->GetLabel()>=0/* && trc->Charge()>0*/) {
1165         htrPt = tmpTr.Pt();
1166         double dy = clus->GetY()-tmpTr.GetY();
1167         double dz = clus->GetZ()-tmpTr.GetZ();
1168         ((TH2*)fCHistoArr->At(kHResY+hcOffs))->Fill(htrPt,dy);
1169         ((TH2*)fCHistoArr->At(kHResZ+hcOffs))->Fill(htrPt,dz);
1170         double errY = tmpTr.GetSigmaY2();
1171         double errZ = tmpTr.GetSigmaZ2();
1172         if (errY>0) ((TH2*)fCHistoArr->At(kHResYP+hcOffs))->Fill(htrPt,dy/Sqrt(errY));
1173         if (errZ>0) ((TH2*)fCHistoArr->At(kHResZP+hcOffs))->Fill(htrPt,dz/Sqrt(errZ));
1174       }
1175 #endif  
1176       //
1177       double chi2;
1178       if ( (chi2=tmpTr.Update(clus))<0 ) {
1179         AliDebug(2,Form("Failed on Update | ESDtrack#%d (MClb:%d)",fCurrESDtrack->GetID(),fCurrESDtrMClb));             
1180         return kFALSE;
1181       }
1182 #ifdef  _FILL_CONTROL_HISTOS_
1183       if (htrPt>0) {
1184         ((TH2*)fCHistoArr->At(kHChi2Cl+hcOffs))->Fill(htrPt,chi2);
1185       }
1186 #endif 
1187       if (AliDebugLevelClass()>1) {
1188         printf("AfterRefit: "); tmpTr.AliExternalTrackParam::Print();
1189       }
1190       if (stopAtLastCl && ++clCount==ncl) return kTRUE; // it was requested to not propagate after last update
1191     }
1192     //
1193   }
1194   // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
1195   // Still, try to go as close as possible to rDest.
1196   //
1197   if (lrStart!=lrStop) {
1198     //printf("Going to last layer %d -> %d\n",lrStart,lrStop);
1199     if (!TransportToLayer(&tmpTr,lrStart,lrStop)) {
1200       AliDebug(1,Form("Failed on TransportToLayer %d->%d | ESDtrack#%d (MClb:%d)",lrStart,lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1201       //tmpTr.AliExternalTrackParam::Print();
1202       //trc->GetWinner()->Print("etp");
1203       return kTRUE;
1204     }    
1205     if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) {
1206       AliDebug(3,Form("Failed on GoToExitFromLayer %d | ESDtrack#%d (MClb:%d)",lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1207       return kTRUE; // go till the exit from layer
1208     }
1209     //
1210     //printf("On exit from last layer\n");
1211     //tmpTr.AliExternalTrackParam::Print();
1212     // go to the destination radius
1213     if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),dir)) return kTRUE;
1214     if (!PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) return kTRUE;
1215   }
1216   trc->AliKalmanTrack::operator=(tmpTr);
1217   if (AliDebugLevelClass()>2) {
1218     printf("After refit (now at lr %d): ",lrStop); trc->AliExternalTrackParam::Print();
1219   }
1220   return kTRUE;
1221 }
1222
1223 //__________________________________________________________________
1224 void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp) 
1225 {
1226   // build MC label
1227   //
1228   const int kMaxLbPerCl = 3;
1229   int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
1230   Int_t lr,nLab=0,nCl=0;
1231   AliITSUSeed *seed = hyp->GetWinner();
1232   while(seed) {
1233     int clID = seed->GetLrCluster(lr);
1234     if (clID>=0) {
1235       AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1236       nCl++;
1237       for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
1238         int trLb = cl->GetLabel(imc);
1239         if (trLb<0) break;
1240         // search this mc track in already accounted ones
1241         int iLab;
1242         for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
1243         if (iLab<nLab) lbStat[iLab]++;
1244         else {
1245           lbID[nLab] = trLb;
1246           lbStat[nLab++] = 1;
1247         }
1248       } // loop over given cluster's labels
1249     }
1250     seed = (AliITSUSeed*)seed->GetParent();
1251   } // loop over clusters
1252   // 
1253   if (nCl && nLab) {
1254     int maxLab=0,nTPCok=0;
1255     AliESDtrack* esdTr = hyp->GetESDTrack();
1256     int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
1257     for (int ilb=nLab;ilb--;) {
1258       int st = lbStat[ilb];
1259       if (lbStat[maxLab]<st) maxLab=ilb;
1260       if (tpcLab==lbID[ilb]) nTPCok += st;
1261     }
1262     hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
1263     hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
1264     hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbStat[maxLab] : -lbStat[maxLab]); // winner label
1265     return;
1266   }
1267   //
1268   hyp->SetFakeRatio(-1.);
1269   hyp->SetLabel( kDummyLabel );
1270   hyp->SetITSLabel( kDummyLabel );
1271 }
1272
1273 //__________________________________________________________________
1274 Bool_t AliITSUTrackerGlo::AddSeedBranch(AliITSUSeed* seed)
1275 {
1276   // add new prolongation branch for the seed to temporary array. The seeds are sorted according to Compare f-n
1277   if (fNCandidatesAdded+fNBranchesAdded>=AliITSUTrackCond::kMaxCandidates ) {
1278     AliError(Form("Number of candidates for current seed reached limit of AliITSUTrackCond::kMaxCandidates=%d, increase it",AliITSUTrackCond::kMaxCandidates));
1279     return kFALSE;
1280   }
1281   AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded]; // note: fNCandidatesAdded is incremented after adding all branches of current seed
1282   int slot=fNBranchesAdded++;
1283   for (int slotF=slot;slotF--;) { // slotF is always slot-1
1284     AliITSUSeed* si = branches[slotF];
1285     if (si->Compare(seed)<0) break; // found the last seed with better quality
1286     // otherwise, shift the worse seed to the next slot
1287     branches[slot] = si;
1288     slot = slotF; // slot should be slotF+1
1289   }
1290   // if needed, move worse seeds
1291   branches[slot] = seed;
1292   return kTRUE;
1293   //
1294 }
1295
1296 //__________________________________________________________________
1297 void AliITSUTrackerGlo::ValidateAllowedBranches(Int_t acceptMax)
1298 {
1299   // keep allowed number of branches for current seed and disable extras
1300   int nb = Min(fNBranchesAdded,acceptMax);
1301   //  if (nb<fNBranchesAdded) printf("ValidateAllowedBranches: %d of %d (%d) on lr %d\n",nb,fNBranchesAdded,fNCandidatesAdded,ilr);
1302   // disable unused branches
1303   AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded];
1304   for (int ib=nb;ib<fNBranchesAdded;ib++) {
1305     /*
1306     if (!branches[ib]->ContainsFake()) {
1307       printf("Suppress good branch as %d of %d |",ib,fNBranchesAdded); branches[ib]->Print();
1308       printf("Survivors : \n");
1309       for (int j=0;j<nb;j++) branches[j]->Print();
1310     }
1311     */
1312     MarkSeedFree(branches[ib]);
1313   }
1314   fNCandidatesAdded += nb; // update total candidates counter
1315   fNBranchesAdded = 0; // reset branches counter
1316   //
1317 }
1318
1319 //__________________________________________________________________
1320 void AliITSUTrackerGlo::ValidateAllowedCandidates(Int_t ilr, Int_t acceptMax)
1321 {
1322   // transfer allowed number of branches to hypothesis container
1323   //
1324   // sort candidates in increasing order of chi2
1325   Int_t index[AliITSUTrackCond::kMaxCandidates];
1326   Float_t chi2[AliITSUTrackCond::kMaxCandidates];
1327   for (int i=fNCandidatesAdded;i--;) chi2[i] = fLayerCandidates[i]->GetChi2GloNrm();
1328   Sort(fNCandidatesAdded,chi2,index,kFALSE);
1329   //
1330   int nb = Min(fNCandidatesAdded,acceptMax);
1331   //  if (nb<fNCandidatesAdded) printf("ValidateAllowedCandidates: %d of %d on lr %d\n",nb,fNCandidatesAdded,ilr);
1332   //
1333   for (int ib=0;ib<nb;ib++) AddProlongationHypothesis(fLayerCandidates[index[ib]],ilr);
1334   // disable unused candidates
1335   for (int ib=nb;ib<fNCandidatesAdded;ib++) {
1336     /*
1337     if (!fLayerCandidates[index[ib]]->ContainsFake()) {
1338       printf("Suppress good candidate as %d of %d |",index[ib],fNCandidatesAdded); fLayerCandidates[index[ib]]->Print();
1339     }
1340     */
1341     MarkSeedFree(fLayerCandidates[index[ib]]);    
1342   }
1343   fNCandidatesAdded = 0; // reset candidates counter
1344   //
1345 }
1346
1347 //__________________________________________________________________
1348 void AliITSUTrackerGlo::FlagSplitClusters()
1349 {
1350   // set special bit on split clusters using MC info
1351   for (int ilr=fITS->GetNLayersActive();ilr--;) {
1352     int nsplit=0;
1353     AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1354     for (int isn=lr->GetNSensors();isn--;) {
1355       AliITSURecoSens* sens = lr->GetSensor(isn);
1356       int nCl = sens->GetNClusters();
1357       if (!nCl) continue;
1358       int cl0 = sens->GetFirstClusterId();
1359       for (int icl=nCl;icl--;) {
1360         AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(cl0+icl);
1361         for (int icl1=icl;icl1--;) {
1362           AliITSUClusterPix *cl1 = (AliITSUClusterPix*)lr->GetCluster(cl0+icl1);
1363           if (cl->HasCommonTrack(cl1)) {
1364             if (!cl->IsSplit())  nsplit++;
1365             if (!cl1->IsSplit()) nsplit++;
1366             cl->SetSplit();
1367             cl1->SetSplit();
1368           }
1369         }
1370       }
1371     }
1372     AliInfo(Form("%4d out of %4d clusters are split",nsplit,lr->GetNClusters()));
1373   }
1374   //
1375 }
1376
1377 //__________________________________________________________________
1378 Bool_t AliITSUTrackerGlo::ContainsSplitCluster(const AliITSUSeed* seed, Int_t maxSize)
1379 {
1380   // check if the seed contains split cluster with size < maxSize
1381   int lrID,clID;
1382   if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1383     AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1384     if (cl->IsSplit() && cl->GetNPix()<maxSize ) return kTRUE;
1385   }
1386   return seed->GetParent() ? ContainsSplitCluster((AliITSUSeed*)seed->GetParent(),maxSize) : kFALSE;
1387   //
1388 }
1389
1390 //__________________________________________________________________
1391 void AliITSUTrackerGlo::PrintSeedClusters(const AliITSUSeed* seed, Option_t* option)
1392 {
1393   // print seeds clusters
1394   int lrID,clID;
1395   if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1396     AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1397     cl->Print(option);
1398   }
1399   if (seed->GetParent()) PrintSeedClusters((AliITSUSeed*)seed->GetParent(), option);
1400   //
1401 }
1402
1403 #ifdef  _FILL_CONTROL_HISTOS_
1404 //__________________________________________________________________
1405 void AliITSUTrackerGlo::BookControlHistos()
1406 {
1407   // book special control histos
1408   if (!fCHistoArr) { // create control histos
1409     const int kNResDef=7;
1410     const double kResDef[kNResDef]={0.05,0.05,0.3, 0.05,1,0.5,1.5};
1411     fCHistoArr = new TObjArray();
1412     fCHistoArr->SetOwner(kTRUE);
1413     const double ptMax=10;
1414     const double plMax=10;
1415     const double chiMax=100;
1416     const int nptbins=50;
1417     const int nresbins=400;
1418     const int nplbins=50;
1419     const int nchbins=200;
1420     int nblr = fITS->GetNLayersActive();
1421     TString ttl;
1422     for (int stp=0;stp<kNTrackingPhases;stp++) {
1423       for (int ilr=0;ilr<nblr;ilr++) {
1424         int hoffs = stp*kHistosPhase + ilr;
1425         double mxdf = ilr>=kNResDef ? kResDef[kNResDef-1] : kResDef[ilr];
1426         ttl = Form("S%d_residY%d",stp,ilr);
1427         TH2F* hdy = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1428         fCHistoArr->AddAtAndExpand(hdy,hoffs + kHResY);
1429         hdy->SetDirectory(0);
1430         //
1431         ttl = Form("S%d_residYPull%d",stp,ilr); 
1432         TH2F* hdyp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1433         fCHistoArr->AddAtAndExpand(hdyp,hoffs + kHResYP);
1434         hdyp->SetDirectory(0);
1435         //
1436         ttl = Form("S%d_residZ%d",stp,ilr);     
1437         TH2F* hdz = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1438         fCHistoArr->AddAtAndExpand(hdz,hoffs + kHResZ);
1439         hdz->SetDirectory(0);
1440         //
1441         ttl = Form("S%d_residZPull%d",stp,ilr);         
1442         TH2F* hdzp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1443         hdzp->SetDirectory(0);
1444         fCHistoArr->AddAtAndExpand(hdzp,hoffs + kHResZP);
1445         //
1446         ttl = Form("S%d_chi2Cl%d",stp,ilr);             
1447         TH2F* hchi = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1448         hchi->SetDirectory(0);
1449         fCHistoArr->AddAtAndExpand(hchi,hoffs + kHChi2Cl);
1450       }
1451     }
1452   }
1453   //  
1454 }
1455 #endif