]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDseedV1.cxx
Complete implementation of pools, see #88914. From rev. 53550,53557,53568 (Ruben)
[u/mrichter/AliRoot.git] / TRD / AliTRDseedV1.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 /* $Id$ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 ////
20 //  The TRD offline tracklet
21 //
22 // The running horse of the TRD reconstruction. The following tasks are preformed:
23 //   1. Clusters attachment to tracks based on prior information stored at tracklet level (see AttachClusters)
24 //   2. Clusters position recalculation based on track information (see GetClusterXY and Fit)
25 //   3. Cluster error parametrization recalculation (see Fit)
26 //   4. Linear track approximation (Fit)
27 //   5. Optimal position (including z estimate for pad row cross tracklets) and covariance matrix of the track fit inside one TRD chamber (Fit)
28 //   6. Tilt pad correction and systematic effects (GetCovAt)
29 //   7. dEdx calculation (CookdEdx)
30 //   8. PID probabilities estimation (CookPID)
31 //
32 //  Authors:                                                              //
33 //    Alex Bercuci <A.Bercuci@gsi.de>                                     //
34 //    Markus Fasel <M.Fasel@gsi.de>                                       //
35 //                                                                        //
36 ////////////////////////////////////////////////////////////////////////////
37
38 #include "TMath.h"
39 #include "TTreeStream.h"
40 #include "TGraphErrors.h"
41
42 #include "AliLog.h"
43 #include "AliMathBase.h"
44 #include "AliRieman.h"
45 #include "AliCDBManager.h"
46
47 #include "AliTRDReconstructor.h"
48 #include "AliTRDpadPlane.h"
49 #include "AliTRDtransform.h"
50 #include "AliTRDcluster.h"
51 #include "AliTRDseedV1.h"
52 #include "AliTRDtrackV1.h"
53 #include "AliTRDcalibDB.h"
54 #include "AliTRDchamberTimeBin.h"
55 #include "AliTRDtrackingChamber.h"
56 #include "AliTRDtrackerV1.h"
57 #include "AliTRDrecoParam.h"
58 #include "AliTRDCommonParam.h"
59 #include "AliTRDtrackletOflHelper.h"
60
61 #include "Cal/AliTRDCalTrkAttach.h"
62 #include "Cal/AliTRDCalPID.h"
63 #include "Cal/AliTRDCalROC.h"
64 #include "Cal/AliTRDCalDet.h"
65
66 class AliTracker;
67
68 ClassImp(AliTRDseedV1)
69
70 //____________________________________________________________________
71 AliTRDseedV1::AliTRDseedV1(Int_t det) 
72   :AliTRDtrackletBase()
73   ,fkReconstructor(NULL)
74   ,fClusterIter(NULL)
75   ,fExB(0.)
76   ,fVD(0.)
77   ,fT0(0.)
78   ,fS2PRF(0.)
79   ,fDiffL(0.)
80   ,fDiffT(0.)
81   ,fClusterIdx(0)
82   ,fErrorMsg(0)
83   ,fN(0)
84   ,fDet(det)
85   ,fPt(0.)
86   ,fdX(0.)
87   ,fX0(0.)
88   ,fX(0.)
89   ,fY(0.)
90   ,fZ(0.)
91   ,fS2Y(0.)
92   ,fS2Z(0.)
93   ,fChi2(0.)
94 {
95   //
96   // Constructor
97   //
98   memset(fIndexes,0xFF,kNclusters*sizeof(fIndexes[0]));
99   memset(fClusters, 0, kNclusters*sizeof(AliTRDcluster*));
100   memset(fPad, 0, 4*sizeof(Float_t));
101   fYref[0] = 0.; fYref[1] = 0.; 
102   fZref[0] = 0.; fZref[1] = 0.; 
103   fYfit[0] = 0.; fYfit[1] = 0.; 
104   fZfit[0] = 0.; fZfit[1] = 0.; 
105   memset(fdEdx, 0, kNslices*sizeof(Float_t)); 
106   for(int ispec=0; ispec<AliPID::kSPECIES; ispec++) fProb[ispec]  = -1.;
107   fLabels[0]=-1; fLabels[1]=-1; // most freq MC labels
108   fLabels[2]=0;  // number of different labels for tracklet
109   memset(fRefCov, 0, 7*sizeof(Double_t));
110   // stand alone curvature
111   fC[0] = 0.; fC[1] = 0.; 
112   // covariance matrix [diagonal]
113   // default sy = 200um and sz = 2.3 cm 
114   fCov[0] = 4.e-4; fCov[1] = 0.; fCov[2] = 5.3; 
115   SetStandAlone(kFALSE);
116 }
117
118 //____________________________________________________________________
119 AliTRDseedV1::AliTRDseedV1(const AliTRDseedV1 &ref)
120   :AliTRDtrackletBase((AliTRDtrackletBase&)ref)
121   ,fkReconstructor(NULL)
122   ,fClusterIter(NULL)
123   ,fExB(0.)
124   ,fVD(0.)
125   ,fT0(0.)
126   ,fS2PRF(0.)
127   ,fDiffL(0.)
128   ,fDiffT(0.)
129   ,fClusterIdx(0)
130   ,fErrorMsg(0)
131   ,fN(0)
132   ,fDet(-1)
133   ,fPt(0.)
134   ,fdX(0.)
135   ,fX0(0.)
136   ,fX(0.)
137   ,fY(0.)
138   ,fZ(0.)
139   ,fS2Y(0.)
140   ,fS2Z(0.)
141   ,fChi2(0.)
142 {
143   //
144   // Copy Constructor performing a deep copy
145   //
146   if(this != &ref){
147     ref.Copy(*this);
148   }
149   SetBit(kOwner, kFALSE);
150   SetStandAlone(ref.IsStandAlone());
151 }
152
153
154 //____________________________________________________________________
155 AliTRDseedV1& AliTRDseedV1::operator=(const AliTRDseedV1 &ref)
156 {
157   //
158   // Assignment Operator using the copy function
159   //
160
161   if(this != &ref){
162     ref.Copy(*this);
163   }
164   SetBit(kOwner, kFALSE);
165
166   return *this;
167 }
168
169 //____________________________________________________________________
170 AliTRDseedV1::~AliTRDseedV1()
171 {
172   //
173   // Destructor. The RecoParam object belongs to the underlying tracker.
174   //
175
176   //printf("I-AliTRDseedV1::~AliTRDseedV1() : Owner[%s]\n", IsOwner()?"YES":"NO");
177
178   if(IsOwner()) {
179     for(int itb=0; itb<kNclusters; itb++){
180       if(!fClusters[itb]) continue; 
181       //AliInfo(Form("deleting c %p @ %d", fClusters[itb], itb));
182       delete fClusters[itb];
183       fClusters[itb] = NULL;
184     }
185   }
186 }
187
188 //____________________________________________________________________
189 void AliTRDseedV1::Clear(Option_t *)
190 {
191   // clean
192   if (IsOwner()) {
193     for(int itb=0; itb<kNclusters; itb++){
194       if(!fClusters[itb]) continue; 
195       //AliInfo(Form("deleting c %p @ %d", fClusters[itb], itb));
196       delete fClusters[itb];
197       fClusters[itb] = NULL;
198     }
199   }
200 }
201
202 //____________________________________________________________________
203 void AliTRDseedV1::Copy(TObject &ref) const
204 {
205   //
206   // Copy function
207   //
208
209   //AliInfo("");
210   AliTRDseedV1 &target = (AliTRDseedV1 &)ref; 
211
212   target.fkReconstructor = fkReconstructor;
213   target.fClusterIter   = NULL;
214   target.fExB           = fExB;
215   target.fVD            = fVD;
216   target.fT0            = fT0;
217   target.fS2PRF         = fS2PRF;
218   target.fDiffL         = fDiffL;
219   target.fDiffT         = fDiffT;
220   target.fClusterIdx    = 0;
221   target.fErrorMsg      = fErrorMsg;
222   target.fN             = fN;
223   target.fDet           = fDet;
224   target.fPt            = fPt;
225   target.fdX            = fdX;
226   target.fX0            = fX0;
227   target.fX             = fX;
228   target.fY             = fY;
229   target.fZ             = fZ;
230   target.fS2Y           = fS2Y;
231   target.fS2Z           = fS2Z;
232   target.fChi2          = fChi2;
233   
234   memcpy(target.fIndexes, fIndexes, kNclusters*sizeof(Int_t));
235   memcpy(target.fClusters, fClusters, kNclusters*sizeof(AliTRDcluster*));
236   memcpy(target.fPad, fPad, 4*sizeof(Float_t));
237   target.fYref[0] = fYref[0]; target.fYref[1] = fYref[1]; 
238   target.fZref[0] = fZref[0]; target.fZref[1] = fZref[1]; 
239   target.fYfit[0] = fYfit[0]; target.fYfit[1] = fYfit[1]; 
240   target.fZfit[0] = fZfit[0]; target.fZfit[1] = fZfit[1]; 
241   memcpy(target.fdEdx, fdEdx, kNslices*sizeof(Float_t)); 
242   memcpy(target.fProb, fProb, AliPID::kSPECIES*sizeof(Float_t)); 
243   memcpy(target.fLabels, fLabels, 3*sizeof(Int_t)); 
244   memcpy(target.fRefCov, fRefCov, 7*sizeof(Double_t)); 
245   target.fC[0] = fC[0]; target.fC[1] = fC[1];
246   memcpy(target.fCov, fCov, 3*sizeof(Double_t)); 
247   
248   TObject::Copy(ref);
249 }
250
251
252 //____________________________________________________________
253 void AliTRDseedV1::Init(const AliRieman *rieman)
254 {
255 // Initialize this tracklet using the riemann fit information
256
257
258   fZref[0] = rieman->GetZat(fX0);
259   fZref[1] = rieman->GetDZat(fX0);
260   fYref[0] = rieman->GetYat(fX0);
261   fYref[1] = rieman->GetDYat(fX0);
262   if(fkReconstructor && fkReconstructor->IsHLT()){
263     fRefCov[0] = 1;
264     fRefCov[2] = 10;
265   }else{
266     fRefCov[0] = rieman->GetErrY(fX0);
267     fRefCov[2] = rieman->GetErrZ(fX0);
268   }
269   fC[0]    = rieman->GetC(); 
270   fChi2    = rieman->GetChi2();
271 }
272
273
274 //____________________________________________________________
275 Bool_t AliTRDseedV1::Init(AliTRDtrackV1 *track)
276 {
277 // Initialize this tracklet using the track information
278 //
279 // Parameters:
280 //   track - the TRD track used to initialize the tracklet
281 // 
282 // Detailed description
283 // The function sets the starting point and direction of the
284 // tracklet according to the information from the TRD track.
285 // 
286 // Caution
287 // The TRD track has to be propagated to the beginning of the
288 // chamber where the tracklet will be constructed
289 //
290
291   Double_t y, z; 
292   if(!track->GetProlongation(fX0, y, z)) return kFALSE;
293   Update(track);
294   return kTRUE;
295 }
296
297
298 //_____________________________________________________________________________
299 void AliTRDseedV1::Reset(Option_t *opt)
300 {
301 //
302 // Reset seed. If option opt="c" is given only cluster arrays are cleared.
303 //
304   for(Int_t ic=kNclusters; ic--;) fIndexes[ic] = -1;
305   memset(fClusters, 0, kNclusters*sizeof(AliTRDcluster*));
306   fN=0; SetBit(kRowCross, kFALSE);
307   if(strcmp(opt, "c")==0) return;
308
309   fExB=0.;fVD=0.;fT0=0.;fS2PRF=0.;
310   fDiffL=0.;fDiffT=0.;
311   fClusterIdx=0;
312   fErrorMsg = 0;
313   fDet=-1;
314   fPt=0.;
315   fdX=0.;fX0=0.; fX=0.; fY=0.; fZ=0.;
316   fS2Y=0.; fS2Z=0.;
317   fC[0]=0.; fC[1]=0.; 
318   fChi2 = 0.;
319
320   memset(fPad, 0, 4*sizeof(Float_t));
321   fYref[0] = 0.; fYref[1] = 0.; 
322   fZref[0] = 0.; fZref[1] = 0.; 
323   fYfit[0] = 0.; fYfit[1] = 0.; 
324   fZfit[0] = 0.; fZfit[1] = 0.; 
325   memset(fdEdx, 0, kNslices*sizeof(Float_t)); 
326   for(int ispec=0; ispec<AliPID::kSPECIES; ispec++) fProb[ispec]  = -1.;
327   fLabels[0]=-1; fLabels[1]=-1; // most freq MC labels
328   fLabels[2]=0;  // number of different labels for tracklet
329   memset(fRefCov, 0, 7*sizeof(Double_t));
330   // covariance matrix [diagonal]
331   // default sy = 200um and sz = 2.3 cm 
332   fCov[0] = 4.e-4; fCov[1] = 0.; fCov[2] = 5.3; 
333 }
334
335 //____________________________________________________________________
336 void AliTRDseedV1::Update(const AliTRDtrackV1 *trk)
337
338   // update tracklet reference position from the TRD track
339
340   Double_t fSnp = trk->GetSnp();
341   Double_t fTgl = trk->GetTgl();
342   fPt = trk->Pt();
343   Double_t norm =1./TMath::Sqrt((1.-fSnp)*(1.+fSnp)); 
344   fYref[1] = fSnp*norm;
345   fZref[1] = fTgl*norm;
346   SetCovRef(trk->GetCovariance());
347
348   Double_t dx = trk->GetX() - fX0;
349   fYref[0] = trk->GetY() - dx*fYref[1];
350   fZref[0] = trk->GetZ() - dx*fZref[1];
351 }
352
353 //_____________________________________________________________________________
354 void AliTRDseedV1::UpdateUsed()
355 {
356   //
357   // Calculate number of used clusers in the tracklet
358   //
359
360   Int_t nused = 0, nshared = 0;
361   for (Int_t i = kNclusters; i--; ) {
362     if (!fClusters[i]) continue;
363     if(fClusters[i]->IsUsed()){ 
364       nused++;
365     } else if(fClusters[i]->IsShared()){
366       if(IsStandAlone()) nused++;
367       else nshared++;
368     }
369   }
370   SetNUsed(nused);
371   SetNShared(nshared);
372 }
373
374 //_____________________________________________________________________________
375 void AliTRDseedV1::UseClusters()
376 {
377   //
378   // Use clusters
379   //
380   // In stand alone mode:
381   // Clusters which are marked as used or shared from another track are
382   // removed from the tracklet
383   //
384   // In barrel mode:
385   // - Clusters which are used by another track become shared
386   // - Clusters which are attached to a kink track become shared
387   //
388   AliTRDcluster **c = &fClusters[0];
389   for (Int_t ic=kNclusters; ic--; c++) {
390     if(!(*c)) continue;
391     if(IsStandAlone()){
392       if((*c)->IsShared() || (*c)->IsUsed()){ 
393         if((*c)->IsShared()) SetNShared(GetNShared()-1);
394         else SetNUsed(GetNUsed()-1);
395         (*c) = NULL;
396         fIndexes[ic] = -1;
397         SetN(GetN()-1);
398         continue;
399       }
400     } else {
401       if((*c)->IsUsed() || IsKink()){
402         (*c)->SetShared();
403         continue;
404       }
405     }
406     (*c)->Use();
407   }
408 }
409
410
411
412 //____________________________________________________________________
413 void AliTRDseedV1::CookdEdx(Int_t nslices)
414 {
415 // Calculates average dE/dx for all slices and store them in the internal array fdEdx. 
416 //
417 // Parameters:
418 //  nslices : number of slices for which dE/dx should be calculated
419 // Output:
420 //  store results in the internal array fdEdx. This can be accessed with the method
421 //  AliTRDseedV1::GetdEdx()
422 //
423 // Detailed description
424 // Calculates average dE/dx for all slices. Depending on the PID methode 
425 // the number of slices can be 3 (LQ) or 8(NN). 
426 // The calculation of dQ/dl are done using the tracklet fit results (see AliTRDseedV1::GetdQdl(Int_t))
427 //
428 // The following effects are included in the calculation:
429 // 1. calibration values for t0 and vdrift (using x coordinate to calculate slice)
430 // 2. cluster sharing (optional see AliTRDrecoParam::SetClusterSharing())
431 // 3. cluster size
432 //
433
434   memset(fdEdx, 0, kNslices*sizeof(Float_t));
435   const Double_t kDriftLength = (.5 * AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
436
437   AliTRDcluster *c(NULL);
438   for(int ic=0; ic<AliTRDtrackerV1::GetNTimeBins(); ic++){
439     if(!(c = fClusters[ic]) && !(c = fClusters[ic+kNtb])) continue;
440     Float_t dx = TMath::Abs(fX0 - c->GetX());
441
442     // Filter clusters for dE/dx calculation
443
444     // 1.consider calibration effects for slice determination
445     Int_t slice;
446     if(dx<kDriftLength){ // TODO should be replaced by c->IsInChamber()
447       slice = Int_t(dx * nslices / kDriftLength);
448     } else slice = c->GetX() < fX0 ? nslices-1 : 0;
449
450
451     // 2. take sharing into account
452     Float_t w = /*c->IsShared() ? .5 :*/ 1.;
453
454     // 3. take into account large clusters TODO
455     //w *= c->GetNPads() > 3 ? .8 : 1.;
456
457     //CHECK !!!
458     fdEdx[slice]   += w * GetdQdl(ic); //fdQdl[ic];
459   } // End of loop over clusters
460 }
461
462 //_____________________________________________________________________________
463 void AliTRDseedV1::CookLabels()
464 {
465   //
466   // Cook 2 labels for seed
467   //
468
469   Int_t labels[200];
470   Int_t out[200];
471   Int_t nlab = 0;
472   for (Int_t i = 0; i < kNclusters; i++) {
473     if (!fClusters[i]) continue;
474     for (Int_t ilab = 0; ilab < 3; ilab++) {
475       if (fClusters[i]->GetLabel(ilab) >= 0) {
476         labels[nlab] = fClusters[i]->GetLabel(ilab);
477         nlab++;
478       }
479     }
480   }
481
482   fLabels[2] = AliMathBase::Freq(nlab,labels,out,kTRUE);
483   fLabels[0] = out[0];
484   if ((fLabels[2]  > 1) && (out[3] > 1)) fLabels[1] = out[2];
485 }
486
487 //____________________________________________________________
488 Float_t AliTRDseedV1::GetAnodeWireOffset(Float_t zt)
489 {
490 // Find position inside the amplification cell for reading drift velocity map
491
492   Float_t d = fPad[3] - zt;
493   if(d<0.){
494     AliError(Form("Fail AnodeWireOffset calculation z0[%+7.2f] zt[%+7.2f] d[%+7.2f].", fPad[3], zt, d));
495     return 0.125;
496   } 
497   d -= ((Int_t)(2 * d)) / 2.0;
498   if(d > 0.25) d = 0.5 - d;
499   return d;
500 }
501
502
503 //____________________________________________________________________
504 Float_t AliTRDseedV1::GetCharge(Bool_t useOutliers) const
505 {
506 // Computes total charge attached to tracklet. If "useOutliers" is set clusters 
507 // which are not in chamber are also used (default false)
508
509   AliTRDcluster *c(NULL); Float_t qt(0.);
510   for(int ic=0; ic<kNclusters; ic++){
511     if(!(c=fClusters[ic])) continue;
512     if(c->IsInChamber() && !useOutliers) continue;
513     qt += TMath::Abs(c->GetQ());
514   }
515   return qt;
516 }
517
518 //____________________________________________________________________
519 Int_t AliTRDseedV1::GetChargeGaps(Float_t sz[kNtb], Float_t pos[kNtb], Int_t isz[kNtb]) const
520 {
521 // Find number, size and position of charge gaps (consecutive missing time bins).
522 // Returns the number of gaps and fills their size in input array "sz" and position in array "pos"
523
524   Bool_t gap(kFALSE);
525   Int_t n(0);
526   Int_t ipos[kNtb]; memset(isz, 0, kNtb*sizeof(Int_t));memset(ipos, 0, kNtb*sizeof(Int_t));
527   for(int ic(0); ic<kNtb; ic++){
528     if(fClusters[ic] || fClusters[ic+kNtb]){
529       if(gap) n++;
530       continue;
531     }
532     gap = kTRUE;
533     isz[n]++;
534     ipos[n] = ic;
535   }
536   if(!n) return 0;
537
538   // write calibrated values
539   AliTRDcluster fake;
540   for(Int_t igap(0); igap<n; igap++){
541     sz[igap] = isz[igap]*fVD/AliTRDCommonParam::Instance()->GetSamplingFrequency();
542     fake.SetPadTime(ipos[igap]);
543     pos[igap] = fake.GetXloc(fT0, fVD);
544     if(isz[igap]>1){
545       fake.SetPadTime(ipos[igap]-isz[igap]+1);
546       pos[igap] += fake.GetXloc(fT0, fVD);
547       pos[igap] /= 2.;
548     }
549   }
550   return n;
551 }
552
553
554 //____________________________________________________________________
555 Bool_t AliTRDseedV1::GetEstimatedCrossPoint(Float_t &x, Float_t &z) const
556 {
557 // Algorithm to estimate cross point in the x-z plane for pad row cross tracklets.
558 // Returns true in case of success.
559   if(!IsRowCross()) return kFALSE;
560
561   x=0.; z=0.;
562   AliTRDcluster *c(NULL);
563   // Find radial range for first row
564   Float_t x1[] = {0., 1.e3};
565   for(int ic=0; ic<kNtb; ic++){
566     if(!(c=fClusters[ic])) continue;
567     if(!c->IsInChamber()) continue;
568     if(c->GetX() <= x1[1]) x1[1] = c->GetX();
569     if(c->GetX() >= x1[0]) x1[0] = c->GetX();
570     z=c->GetZ();
571   }
572   if((x1[0] - x1[1])<1.e-5) return kFALSE;
573
574   // Find radial range for second row
575   Bool_t kZ(kFALSE);
576   Float_t x2[] = {0., 1.e3};
577   for(int ic=kNtb; ic<kNclusters; ic++){
578     if(!(c=fClusters[ic])) continue;
579     if(!c->IsInChamber()) continue;
580     if(c->GetX() <= x2[1]) x2[1] = c->GetX();
581     if(c->GetX() >= x2[0]) x2[0] = c->GetX();
582     if(!kZ){
583       z+=c->GetZ();
584       z*=0.5;
585       kZ=kTRUE;
586     }
587   }
588   if((x2[0] - x2[1])<1.e-5) return kFALSE;
589
590   // Find intersection of the 2 radial regions
591   x = 0.5*((x1[0]+x1[1] > x2[0]+x2[1]) ? (x1[1]+x2[0]) : (x1[0]+x2[1]));
592   return kTRUE;
593 }
594
595 //____________________________________________________________________
596 Float_t AliTRDseedV1::GetQperTB(Int_t tb) const
597 {
598   //
599   // Charge of the clusters at timebin
600   //
601   Float_t Q = 0;
602   if(fClusters[tb] /*&& fClusters[tb]->IsInChamber()*/)
603     Q += TMath::Abs(fClusters[tb]->GetQ());
604   if(fClusters[tb+kNtb] /*&& fClusters[tb+kNtb]->IsInChamber()*/)
605     Q += TMath::Abs(fClusters[tb+kNtb]->GetQ());
606   return Q/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
607 }
608
609 //____________________________________________________________________
610 Float_t AliTRDseedV1::GetdQdl() const
611 {
612 // Calculate total charge / tracklet length for 1D PID
613 //
614   Float_t Q = GetCharge(kTRUE);
615   return Q/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
616 }
617
618 //____________________________________________________________________
619 Float_t AliTRDseedV1::GetdQdl(Int_t ic, Float_t *dl) const
620 {
621 // Using the linear approximation of the track inside one TRD chamber (TRD tracklet) 
622 // the charge per unit length can be written as:
623 // BEGIN_LATEX
624 // #frac{dq}{dl} = #frac{q_{c}}{dx * #sqrt{1 + #(){#frac{dy}{dx}}^{2}_{fit} + #(){#frac{dz}{dx}}^{2}_{ref}}}
625 // END_LATEX
626 // where qc is the total charge collected in the current time bin and dx is the length 
627 // of the time bin. 
628 // The following correction are applied :
629 //   - charge : pad row cross corrections
630 //              [diffusion and TRF assymetry] TODO
631 //   - dx     : anisochronity, track inclination - see Fit and AliTRDcluster::GetXloc() 
632 //              and AliTRDcluster::GetYloc() for the effects taken into account
633 // 
634 //Begin_Html
635 //<img src="TRD/trackletDQDT.gif">
636 //End_Html
637 // In the picture the energy loss measured on the tracklet as a function of drift time [left] and respectively 
638 // drift length [right] for different particle species is displayed.
639 // Author : Alex Bercuci <A.Bercuci@gsi.de>
640 //
641   Float_t dq = 0.;
642   // check whether both clusters are inside the chamber
643   Bool_t hasClusterInChamber = kFALSE;
644   if(fClusters[ic] && fClusters[ic]->IsInChamber()){
645     hasClusterInChamber = kTRUE;
646     dq += TMath::Abs(fClusters[ic]->GetQ());
647   }
648   if(fClusters[ic+kNtb] && fClusters[ic+kNtb]->IsInChamber()){
649     hasClusterInChamber = kTRUE;
650     dq += TMath::Abs(fClusters[ic+kNtb]->GetQ());
651   }
652   if(!hasClusterInChamber) return 0.;
653   if(dq<1.e-3) return 0.;
654
655   Double_t dx = fdX;
656   if(ic-1>=0 && ic+1<kNtb){
657     Float_t x2(0.), x1(0.);
658     // try to estimate upper radial position (find the cluster which is inside the chamber)
659     if(fClusters[ic-1] && fClusters[ic-1]->IsInChamber()) x2 = fClusters[ic-1]->GetX(); 
660     else if(fClusters[ic-1+kNtb] && fClusters[ic-1+kNtb]->IsInChamber()) x2 = fClusters[ic-1+kNtb]->GetX(); 
661     else if(fClusters[ic] && fClusters[ic]->IsInChamber()) x2 = fClusters[ic]->GetX()+fdX;
662     else x2 = fClusters[ic+kNtb]->GetX()+fdX;
663     // try to estimate lower radial position (find the cluster which is inside the chamber)
664     if(fClusters[ic+1] && fClusters[ic+1]->IsInChamber()) x1 = fClusters[ic+1]->GetX();
665     else if(fClusters[ic+1+kNtb] && fClusters[ic+1+kNtb]->IsInChamber()) x1 = fClusters[ic+1+kNtb]->GetX();
666     else if(fClusters[ic] && fClusters[ic]->IsInChamber()) x1 = fClusters[ic]->GetX()-fdX;
667     else x1 = fClusters[ic+kNtb]->GetX()-fdX;
668
669     dx = .5*(x2 - x1);
670   }
671   dx *= TMath::Sqrt(1. + fYfit[1]*fYfit[1] + fZref[1]*fZref[1]);
672   if(dl) (*dl) = dx;
673   if(dx>1.e-9) return dq/dx;
674   else return 0.;
675 }
676
677 //____________________________________________________________
678 Float_t AliTRDseedV1::GetMomentum(Float_t *err) const
679
680 // Returns momentum of the track after update with the current tracklet as:
681 // BEGIN_LATEX
682 // p=#frac{1}{1/p_{t}} #sqrt{1+tgl^{2}}
683 // END_LATEX
684 // and optionally the momentum error (if err is not null). 
685 // The estimated variance of the momentum is given by:
686 // BEGIN_LATEX
687 // #sigma_{p}^{2} = (#frac{dp}{dp_{t}})^{2} #sigma_{p_{t}}^{2}+(#frac{dp}{dtgl})^{2} #sigma_{tgl}^{2}+2#frac{dp}{dp_{t}}#frac{dp}{dtgl} cov(tgl,1/p_{t})
688 // END_LATEX
689 // which can be simplified to
690 // BEGIN_LATEX
691 // #sigma_{p}^{2} = p^{2}p_{t}^{4}tgl^{2}#sigma_{tgl}^{2}-2p^{2}p_{t}^{3}tgl cov(tgl,1/p_{t})+p^{2}p_{t}^{2}#sigma_{1/p_{t}}^{2}
692 // END_LATEX
693 //
694
695   Double_t p = fPt*TMath::Sqrt(1.+fZref[1]*fZref[1]);
696   Double_t p2 = p*p;
697   Double_t tgl2 = fZref[1]*fZref[1];
698   Double_t pt2 = fPt*fPt;
699   if(err){
700     Double_t s2 = 
701       p2*tgl2*pt2*pt2*fRefCov[4]
702      -2.*p2*fZref[1]*fPt*pt2*fRefCov[5]
703      +p2*pt2*fRefCov[6];
704     (*err) = TMath::Sqrt(s2);
705   }
706   return p;
707 }
708
709
710 //____________________________________________________________________
711 Int_t AliTRDseedV1::GetTBoccupancy() const
712 {
713 // Returns no. of TB occupied by clusters
714
715   Int_t n(0);
716   for(int ic(0); ic<kNtb; ic++){
717     if(!fClusters[ic] && !fClusters[ic+kNtb]) continue;
718     n++;
719   }
720   return n;
721 }
722
723 //____________________________________________________________________
724 Int_t AliTRDseedV1::GetTBcross() const
725 {
726 // Returns no. of TB occupied by 2 clusters for pad row cross tracklets
727
728   if(!IsRowCross()) return 0;
729   Int_t n(0);
730   for(int ic(0); ic<kNtb; ic++){
731     if(fClusters[ic] && fClusters[ic+kNtb]) n++;
732   }
733   return n;
734 }
735
736 //____________________________________________________________________
737 Float_t* AliTRDseedV1::GetProbability(Bool_t force)
738 {       
739   if(!force) return &fProb[0];
740   if(!CookPID()) return NULL;
741   return &fProb[0];
742 }
743
744 //____________________________________________________________
745 Bool_t AliTRDseedV1::CookPID()
746 {
747 // Fill probability array for tracklet from the DB.
748 //
749 // Parameters
750 //
751 // Output
752 //   returns pointer to the probability array and NULL if missing DB access 
753 //
754 // Retrieve PID probabilities for e+-, mu+-, K+-, pi+- and p+- from the DB according to tracklet information:
755 // - estimated momentum at tracklet reference point 
756 // - dE/dx measurements
757 // - tracklet length
758 // - TRD layer
759 // According to the steering settings specified in the reconstruction one of the following methods are used
760 // - Neural Network [default] - option "nn"  
761 // - 2D Likelihood - option "!nn"  
762
763   AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
764   if (!calibration) {
765     AliError("No access to calibration data");
766     return kFALSE;
767   }
768
769   if (!fkReconstructor) {
770     AliError("Reconstructor not set.");
771     return kFALSE;
772   }
773
774   // Retrieve the CDB container class with the parametric detector response
775   const AliTRDCalPID *pd = calibration->GetPIDObject(fkReconstructor->GetPIDMethod());
776   if (!pd) {
777     AliError("No access to AliTRDCalPID object");
778     return kFALSE;
779   }
780
781   // calculate tracklet length TO DO
782   Float_t length = (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick())/ TMath::Sqrt((1.0 - GetSnp()*GetSnp()) / (1.0 + GetTgl()*GetTgl()));
783   
784   //calculate dE/dx
785   CookdEdx(AliTRDCalPID::kNSlicesNN);
786   AliDebug(4, Form("p=%6.4f[GeV/c] dEdx{%7.2f %7.2f %7.2f %7.2f %7.2f %7.2f %7.2f %7.2f} l=%4.2f[cm]", GetMomentum(), fdEdx[0], fdEdx[1], fdEdx[2], fdEdx[3], fdEdx[4], fdEdx[5], fdEdx[6], fdEdx[7], length));
787
788   // Sets the a priori probabilities
789   Bool_t kPIDNN(fkReconstructor->GetPIDMethod()==AliTRDpidUtil::kNN);
790   for(int ispec=0; ispec<AliPID::kSPECIES; ispec++)
791     fProb[ispec] = pd->GetProbability(ispec, GetMomentum(), &fdEdx[0], length, kPIDNN?GetPlane():fkReconstructor->GetRecoParam()->GetPIDLQslices());
792   
793   return kTRUE;
794 }
795
796 //____________________________________________________________________
797 Float_t AliTRDseedV1::GetQuality(Bool_t kZcorr) const
798 {
799   //
800   // Returns a quality measurement of the current seed
801   //
802
803   Float_t zcorr = kZcorr ? GetTilt() * (fZfit[0] - fZref[0]) : 0.;
804   return 
805       .5 * TMath::Abs(18.0 - GetN())
806     + 10.* TMath::Abs(fYfit[1] - fYref[1])
807     + 5. * TMath::Abs(fYfit[0] - fYref[0] + zcorr)
808     + 2. * TMath::Abs(fZfit[0] - fZref[0]) / GetPadLength();
809 }
810
811 //____________________________________________________________________
812 void AliTRDseedV1::GetCovAt(Double_t x, Double_t *cov) const
813 {
814 // Computes covariance in the y-z plane at radial point x (in tracking coordinates) 
815 // and returns the results in the preallocated array cov[3] as :
816 //   cov[0] = Var(y)
817 //   cov[1] = Cov(yz)
818 //   cov[2] = Var(z)
819 //
820 // Details
821 //
822 // For the linear transformation
823 // BEGIN_LATEX
824 // Y = T_{x} X^{T}
825 // END_LATEX
826 //   The error propagation has the general form
827 // BEGIN_LATEX
828 // C_{Y} = T_{x} C_{X} T_{x}^{T} 
829 // END_LATEX
830 //  We apply this formula 2 times. First to calculate the covariance of the tracklet 
831 // at point x we consider: 
832 // BEGIN_LATEX
833 // T_{x} = (1 x); X=(y0 dy/dx); C_{X}=#(){#splitline{Var(y0) Cov(y0, dy/dx)}{Cov(y0, dy/dx) Var(dy/dx)}} 
834 // END_LATEX
835 // and secondly to take into account the tilt angle
836 // BEGIN_LATEX
837 // T_{#alpha} = #(){#splitline{cos(#alpha) __ sin(#alpha)}{-sin(#alpha) __ cos(#alpha)}}; X=(y z); C_{X}=#(){#splitline{Var(y)    0}{0   Var(z)}} 
838 // END_LATEX
839 //
840 // using simple trigonometrics one can write for this last case
841 // BEGIN_LATEX
842 // C_{Y}=#frac{1}{1+tg^{2}#alpha} #(){#splitline{(#sigma_{y}^{2}+tg^{2}#alpha#sigma_{z}^{2}) __ tg#alpha(#sigma_{z}^{2}-#sigma_{y}^{2})}{tg#alpha(#sigma_{z}^{2}-#sigma_{y}^{2}) __ (#sigma_{z}^{2}+tg^{2}#alpha#sigma_{y}^{2})}} 
843 // END_LATEX
844 // which can be aproximated for small alphas (2 deg) with
845 // BEGIN_LATEX
846 // C_{Y}=#(){#splitline{#sigma_{y}^{2} __ (#sigma_{z}^{2}-#sigma_{y}^{2})tg#alpha}{((#sigma_{z}^{2}-#sigma_{y}^{2})tg#alpha __ #sigma_{z}^{2}}} 
847 // END_LATEX
848 //
849 // before applying the tilt rotation we also apply systematic uncertainties to the tracklet 
850 // position which can be tunned from outside via the AliTRDrecoParam::SetSysCovMatrix(). They might 
851 // account for extra misalignment/miscalibration uncertainties. 
852 //
853 // Author :
854 // Alex Bercuci <A.Bercuci@gsi.de> 
855 // Date : Jan 8th 2009
856 //
857
858
859   Double_t xr     = fX0-x; 
860   Double_t sy2    = fCov[0] +2.*xr*fCov[1] + xr*xr*fCov[2];
861   Double_t sz2    = fS2Z;
862   //GetPadLength()*GetPadLength()/12.;
863
864   // insert systematic uncertainties
865   if(fkReconstructor){
866     Double_t sys[15]; memset(sys, 0, 15*sizeof(Double_t));
867     fkReconstructor->GetRecoParam()->GetSysCovMatrix(sys);
868     sy2 += sys[0];
869     sz2 += sys[1];
870   }
871
872   // rotate covariance matrix if no RC
873   if(!IsRowCross()){
874     Double_t t2 = GetTilt()*GetTilt();
875     Double_t correction = 1./(1. + t2);
876     cov[0] = (sy2+t2*sz2)*correction;
877     cov[1] = GetTilt()*(sz2 - sy2)*correction;
878     cov[2] = (t2*sy2+sz2)*correction;
879   } else {
880     cov[0] = sy2; cov[1] = 0.; cov[2] = sy2;
881   }
882
883   AliDebug(4, Form("C(%6.1f %+6.3f %6.1f)  RC[%c]", 1.e4*TMath::Sqrt(cov[0]), cov[1], 1.e4*TMath::Sqrt(cov[2]), IsRowCross()?'y':'n'));
884 }
885
886 //____________________________________________________________
887 Int_t AliTRDseedV1::GetCovSqrt(const Double_t * const c, Double_t *d)
888 {
889 // Helper function to calculate the square root of the covariance matrix. 
890 // The input matrix is stored in the vector c and the result in the vector d. 
891 // Both arrays have to be initialized by the user with at least 3 elements. Return negative in case of failure.
892 // 
893 // For calculating the square root of the symmetric matrix c
894 // the following relation is used:
895 // BEGIN_LATEX
896 // C^{1/2} = VD^{1/2}V^{-1}
897 // END_LATEX
898 // with V being the matrix with the n eigenvectors as columns. 
899 // In case C is symmetric the followings are true:
900 //   - matrix D is diagonal with the diagonal given by the eigenvalues of C
901 //   - V = V^{-1}
902 //
903 // Author A.Bercuci <A.Bercuci@gsi.de>
904 // Date   Mar 19 2009
905
906   const Double_t kZero(1.e-20);
907   Double_t l[2], // eigenvalues
908            v[3]; // eigenvectors
909   // the secular equation and its solution :
910   // (c[0]-L)(c[2]-L)-c[1]^2 = 0
911   // L^2 - L*Tr(c)+DET(c) = 0
912   // L12 = [Tr(c) +- sqrt(Tr(c)^2-4*DET(c))]/2
913   Double_t tr = c[0]+c[2],           // trace
914           det = c[0]*c[2]-c[1]*c[1]; // determinant
915   if(TMath::Abs(det)<kZero) return 1;
916   Double_t dd = TMath::Sqrt(tr*tr - 4*det);
917   l[0] = .5*(tr + dd*(c[0]>c[2]?-1.:1.));
918   l[1] = .5*(tr + dd*(c[0]>c[2]?1.:-1.));
919   if(l[0]<kZero || l[1]<kZero) return 2;
920   // the sym V matrix
921   // | v00   v10|
922   // | v10   v11|
923   Double_t den = (l[0]-c[0])*(l[0]-c[0])+c[1]*c[1];
924   if(den<kZero){ // almost diagonal
925     v[0] = TMath::Sign(0., c[1]);
926     v[1] = TMath::Sign(1., (l[0]-c[0]));
927     v[2] = TMath::Sign(0., c[1]*(l[0]-c[0])*(l[1]-c[2]));
928   } else {
929     Double_t tmp = 1./TMath::Sqrt(den);
930     v[0] = c[1]* tmp;
931     v[1] = (l[0]-c[0])*tmp;
932     if(TMath::Abs(l[1]-c[2])<kZero) v[2] = TMath::Sign(v[0]*(l[0]-c[0])/kZero, (l[1]-c[2]));
933     else v[2] = v[0]*(l[0]-c[0])/(l[1]-c[2]);
934   }
935   // the VD^{1/2}V is: 
936   l[0] = TMath::Sqrt(l[0]); l[1] = TMath::Sqrt(l[1]);
937   d[0] = v[0]*v[0]*l[0]+v[1]*v[1]*l[1];
938   d[1] = v[0]*v[1]*l[0]+v[1]*v[2]*l[1];
939   d[2] = v[1]*v[1]*l[0]+v[2]*v[2]*l[1];
940
941   return 0;
942 }
943
944 //____________________________________________________________
945 Double_t AliTRDseedV1::GetCovInv(const Double_t * const c, Double_t *d)
946 {
947 // Helper function to calculate the inverse of the covariance matrix.
948 // The input matrix is stored in the vector c and the result in the vector d. 
949 // Both arrays have to be initialized by the user with at least 3 elements
950 // The return value is the determinant or 0 in case of singularity.
951 //
952 // Author A.Bercuci <A.Bercuci@gsi.de>
953 // Date   Mar 19 2009
954
955   Double_t det = c[0]*c[2] - c[1]*c[1];
956   if(TMath::Abs(det)<1.e-20) return 0.;
957   Double_t invDet = 1./det;
958   d[0] = c[2]*invDet;
959   d[1] =-c[1]*invDet;
960   d[2] = c[0]*invDet;
961   return det;
962 }
963
964 //____________________________________________________________________
965 UShort_t AliTRDseedV1::GetVolumeId() const
966 {
967 // Returns geometry volume id by delegation 
968
969   for(Int_t ic(0);ic<kNclusters; ic++){
970     if(fClusters[ic]) return fClusters[ic]->GetVolumeId();
971   }
972   return 0;
973 }
974
975
976 //____________________________________________________________________
977 void AliTRDseedV1::Calibrate()
978 {
979 // Retrieve calibration and position parameters from OCDB. 
980 // The following information are used
981 //  - detector index
982 //  - column and row position of first attached cluster. If no clusters are attached 
983 // to the tracklet a random central chamber position (c=70, r=7) will be used.
984 //
985 // The following information is cached in the tracklet
986 //   t0 (trigger delay)
987 //   drift velocity
988 //   PRF width
989 //   omega*tau = tg(a_L)
990 //   diffusion coefficients (longitudinal and transversal)
991 //
992 // Author :
993 // Alex Bercuci <A.Bercuci@gsi.de> 
994 // Date : Jan 8th 2009
995 //
996
997   AliCDBManager *cdb = AliCDBManager::Instance();
998   if(cdb->GetRun() < 0){
999     AliError("OCDB manager not properly initialized");
1000     return;
1001   }
1002
1003   AliTRDcalibDB *calib = AliTRDcalibDB::Instance();
1004   AliTRDCalROC  *vdROC = calib->GetVdriftROC(fDet),
1005                 *t0ROC = calib->GetT0ROC(fDet);;
1006   const AliTRDCalDet *vdDet = calib->GetVdriftDet();
1007   const AliTRDCalDet *t0Det = calib->GetT0Det();
1008
1009   Int_t col = 70, row = 7;
1010   AliTRDcluster **c = &fClusters[0];
1011   if(GetN()){ 
1012     Int_t ic = 0;
1013     while (ic<kNclusters && !(*c)){ic++; c++;} 
1014     if(*c){
1015       col = (*c)->GetPadCol();
1016       row = (*c)->GetPadRow();
1017     }
1018   }
1019
1020   fT0    = (t0Det->GetValue(fDet) + t0ROC->GetValue(col,row)) / AliTRDCommonParam::Instance()->GetSamplingFrequency();
1021   fVD    = vdDet->GetValue(fDet) * vdROC->GetValue(col, row);
1022   fS2PRF = calib->GetPRFWidth(fDet, col, row); fS2PRF *= fS2PRF;
1023   fExB   = AliTRDCommonParam::Instance()->GetOmegaTau(fVD);
1024   AliTRDCommonParam::Instance()->GetDiffCoeff(fDiffL,
1025   fDiffT, fVD);
1026   AliDebug(4, Form("Calibration params for Det[%3d] Col[%3d] Row[%2d]\n  t0[%f]  vd[%f]  s2PRF[%f]  ExB[%f]  Dl[%f]  Dt[%f]", fDet, col, row, fT0, fVD, fS2PRF, fExB, fDiffL, fDiffT));
1027
1028
1029   SetBit(kCalib, kTRUE);
1030 }
1031
1032 //____________________________________________________________________
1033 void AliTRDseedV1::SetOwner()
1034 {
1035   //AliInfo(Form("own [%s] fOwner[%s]", own?"YES":"NO", fOwner?"YES":"NO"));
1036   
1037   if(TestBit(kOwner)) return;
1038   for(int ic=0; ic<kNclusters; ic++){
1039     if(!fClusters[ic]) continue;
1040     fClusters[ic] = new AliTRDcluster(*fClusters[ic]);
1041   }
1042   SetBit(kOwner);
1043 }
1044
1045 //____________________________________________________________
1046 void AliTRDseedV1::SetPadPlane(AliTRDpadPlane * const p)
1047 {
1048 // Shortcut method to initialize pad geometry.
1049   fPad[0] = p->GetLengthIPad();
1050   fPad[1] = p->GetWidthIPad();
1051   fPad[2] = TMath::Tan(TMath::DegToRad()*p->GetTiltingAngle());
1052   fPad[3] = p->GetRow0() + p->GetAnodeWireOffset();
1053 }
1054
1055
1056
1057 //____________________________________________________________________
1058 Bool_t  AliTRDseedV1::AttachClusters(AliTRDtrackingChamber *const chamber, Bool_t tilt, Bool_t chgPos, Int_t ev)
1059 {
1060 //
1061 // Projective algorithm to attach clusters to seeding tracklets. The following steps are performed :
1062 // 1. Collapse x coordinate for the full detector plane
1063 // 2. truncated mean on y (r-phi) direction
1064 // 3. purge clusters
1065 // 4. truncated mean on z direction
1066 // 5. purge clusters
1067 //
1068 // Parameters
1069 //  - chamber : pointer to tracking chamber container used to search the tracklet
1070 //  - tilt    : switch for tilt correction during road building [default true]
1071 //  - chgPos  : mark same[kFALSE] and opposite[kTRUE] sign tracks with respect to Bz field sign [default true]
1072 //  - ev      : event number for debug purposes [default = -1]
1073 // Output
1074 //  - true    : if tracklet found successfully. Failure can happend because of the following:
1075 //      -
1076 // Detailed description
1077 //  
1078 // We start up by defining the track direction in the xy plane and roads. The roads are calculated based
1079 // on tracking information (variance in the r-phi direction) and estimated variance of the standard 
1080 // clusters (see AliTRDcluster::SetSigmaY2()) corrected for tilt (see GetCovAt()). From this the road is
1081 // BEGIN_LATEX
1082 // r_{y} = 3*#sqrt{12*(#sigma^{2}_{Trk}(y) + #frac{#sigma^{2}_{cl}(y) + tg^{2}(#alpha_{L})#sigma^{2}_{cl}(z)}{1+tg^{2}(#alpha_{L})})}
1083 // r_{z} = 1.5*L_{pad}
1084 // END_LATEX
1085 // 
1086 // Author : Alexandru Bercuci <A.Bercuci@gsi.de>
1087 // Debug  : level = 2 for calibration
1088 //          level = 3 for visualization in the track SR
1089 //          level = 4 for full visualization including digit level
1090
1091   const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam(); //the dynamic cast in GetRecoParam is slow, so caching the pointer to it
1092
1093   if(!recoParam){
1094     AliError("Tracklets can not be used without a valid RecoParam.");
1095     return kFALSE;
1096   }
1097   AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1098   if (!calibration) {
1099     AliError("No access to calibration data");
1100     return kFALSE;
1101   }
1102   // Retrieve the CDB container class with the parametric likelihood
1103   const AliTRDCalTrkAttach *attach = calibration->GetAttachObject();
1104   if (!attach) {
1105     AliError("No usable AttachClusters calib object.");
1106     return kFALSE;
1107   }
1108
1109   // Initialize reco params for this tracklet
1110   // 1. first time bin in the drift region
1111   Int_t t0 = 14;
1112   Int_t kClmin = Int_t(recoParam->GetFindableClusters()*AliTRDtrackerV1::GetNTimeBins());
1113   Int_t kTBmin = 4;
1114
1115   Double_t sysCov[5]; recoParam->GetSysCovMatrix(sysCov); 
1116   Double_t s2yTrk= fRefCov[0], 
1117            s2yCl = 0., 
1118            s2zCl = GetPadLength()*GetPadLength()/12., 
1119            syRef = TMath::Sqrt(s2yTrk),
1120            t2    = GetTilt()*GetTilt();
1121   //define roads
1122   const Double_t kroady = 3.; //recoParam->GetRoad1y();
1123   const Double_t kroadz = GetPadLength() * recoParam->GetRoadzMultiplicator() + 1.;
1124   // define probing cluster (the perfect cluster) and default calibration
1125   Short_t sig[] = {0, 0, 10, 30, 10, 0,0};
1126   AliTRDcluster cp(fDet, 6, 75, 0, sig, 0);
1127   if(fkReconstructor->IsHLT()) cp.SetRPhiMethod(AliTRDcluster::kCOG);
1128   if(!IsCalibrated()) Calibrate();
1129
1130 /*  Int_t kroadyShift(0);
1131   Float_t bz(AliTrackerBase::GetBz());
1132   if(TMath::Abs(bz)>2.){
1133     if(bz<0.) kroadyShift = chgPos ? +1 : -1;
1134     else kroadyShift = chgPos ? -1 : +1;
1135   }*/
1136   AliDebug(4, Form("\n       syTrk[cm]=%4.2f dydxTrk[deg]=%+6.2f Chg[%c] rY[cm]=%4.2f rZ[cm]=%5.2f TC[%c]", syRef, TMath::ATan(fYref[1])*TMath::RadToDeg(), chgPos?'+':'-', kroady, kroadz, tilt?'y':'n'));
1137   Double_t phiTrk(TMath::ATan(fYref[1])),
1138            thtTrk(TMath::ATan(fZref[1]));
1139
1140   // working variables
1141   const Int_t kNrows = 16;
1142   const Int_t kNcls  = 3*kNclusters; // buffer size
1143   TObjArray clst[kNrows];
1144   Bool_t blst[kNrows][kNcls];
1145   Double_t cond[4],
1146            dx, dy, dz,
1147            yt, zt,
1148            zc[kNrows],
1149            xres[kNrows][kNcls], yres[kNrows][kNcls], zres[kNrows][kNcls], s2y[kNrows][kNcls];
1150   Int_t idxs[kNrows][kNcls], ncl[kNrows], ncls = 0;
1151   memset(ncl, 0, kNrows*sizeof(Int_t));
1152   memset(zc, 0, kNrows*sizeof(Double_t));
1153   memset(idxs, 0, kNrows*kNcls*sizeof(Int_t));
1154   memset(xres, 0, kNrows*kNcls*sizeof(Double_t));
1155   memset(yres, 0, kNrows*kNcls*sizeof(Double_t));
1156   memset(zres, 0, kNrows*kNcls*sizeof(Double_t));
1157   memset(s2y, 0, kNrows*kNcls*sizeof(Double_t));
1158   memset(blst, 0, kNrows*kNcls*sizeof(Bool_t));   //this is 8 times faster to memset than "memset(clst, 0, kNrows*kNcls*sizeof(AliTRDcluster*))"
1159
1160   Double_t roady(0.), s2Mean(0.), sMean(0.); Int_t ns2Mean(0);
1161
1162   // Do cluster projection and pick up cluster candidates
1163   AliTRDcluster *c(NULL);
1164   AliTRDchamberTimeBin *layer(NULL);
1165   Bool_t kBUFFER = kFALSE;
1166   for (Int_t it = 0; it < kNtb; it++) {
1167     if(!(layer = chamber->GetTB(it))) continue;
1168     if(!Int_t(*layer)) continue;
1169     // get track projection at layers position
1170     dx   = fX0 - layer->GetX();
1171     yt = fYref[0] - fYref[1] * dx;
1172     zt = fZref[0] - fZref[1] * dx;
1173     // get standard cluster error corrected for tilt if selected
1174     cp.SetLocalTimeBin(it);
1175     cp.SetSigmaY2(0.02, fDiffT, fExB, dx, -1./*zt*/, fYref[1]);
1176     s2yCl = cp.GetSigmaY2() + sysCov[0]; if(!tilt) s2yCl = (s2yCl + t2*s2zCl)/(1.+t2);
1177     if(TMath::Abs(it-12)<7){ s2Mean += cp.GetSigmaY2(); ns2Mean++;}
1178     // get estimated road in r-phi direction
1179     roady = TMath::Min(3.*TMath::Sqrt(12.*(s2yTrk + s2yCl)), kroady);
1180
1181     AliDebug(5, Form("\n"
1182       "  %2d xd[cm]=%6.3f yt[cm]=%7.2f zt[cm]=%8.2f\n"
1183       "      syTrk[um]=%6.2f syCl[um]=%6.2f syClTlt[um]=%6.2f\n"
1184       "      Ry[mm]=%f"
1185       , it, dx, yt, zt
1186       , 1.e4*TMath::Sqrt(s2yTrk), 1.e4*TMath::Sqrt(cp.GetSigmaY2()+sysCov[0]), 1.e4*TMath::Sqrt(s2yCl)
1187       , 1.e1*roady));
1188
1189     // get clusters from layer
1190     cond[0] = yt/*+0.5*kroadyShift*kroady*/; cond[2] = roady;
1191     cond[1] = zt; cond[3] = kroadz;
1192     Int_t n=0, idx[6]; layer->GetClusters(cond, idx, n, 6);
1193     for(Int_t ic = n; ic--;){
1194       c  = (*layer)[idx[ic]];
1195       dx = fX0 - c->GetX();
1196       yt = fYref[0] - fYref[1] * dx;
1197       zt = fZref[0] - fZref[1] * dx;
1198       dz = zt - c->GetZ();
1199       dy = yt - (c->GetY() + (tilt ? (GetTilt() * dz) : 0.));
1200       Int_t r = c->GetPadRow();
1201       clst[r].AddAtAndExpand(c, ncl[r]);
1202       blst[r][ncl[r]] = kTRUE;
1203       idxs[r][ncl[r]] = idx[ic];
1204       zres[r][ncl[r]] = dz/GetPadLength();
1205       yres[r][ncl[r]] = dy;
1206       xres[r][ncl[r]] = dx;
1207       zc[r]           = c->GetZ();
1208       // TODO temporary solution to avoid divercences in error parametrization
1209       s2y[r][ncl[r]]  = TMath::Min(c->GetSigmaY2()+sysCov[0], 0.025); 
1210       AliDebug(5, Form("   -> dy[cm]=%+7.4f yc[cm]=%7.2f row[%d] idx[%2d]", dy, c->GetY(), r, ncl[r]));
1211       ncl[r]++; ncls++;
1212
1213       if(ncl[r] >= kNcls) {
1214         AliWarning(Form("Cluster candidates row[%d] reached buffer limit[%d]. Some may be lost.", r, kNcls));
1215         kBUFFER = kTRUE;
1216         break;
1217       }
1218     }
1219     if(kBUFFER) break;
1220   }
1221   if(ncls<kClmin){ 
1222     AliDebug(1, Form("CLUSTERS FOUND %d LESS THAN THRESHOLD %d.", ncls, kClmin));
1223     SetErrorMsg(kAttachClFound);
1224     for(Int_t ir(kNrows);ir--;) clst[ir].Clear();
1225     return kFALSE;
1226   }
1227   if(ns2Mean<kTBmin){
1228     AliDebug(1, Form("CLUSTERS IN TimeBins %d LESS THAN THRESHOLD %d.", ns2Mean, kTBmin));
1229     SetErrorMsg(kAttachClFound);
1230     for(Int_t ir(kNrows);ir--;) clst[ir].Clear();
1231     return kFALSE;
1232   }
1233   s2Mean /= ns2Mean; sMean = TMath::Sqrt(s2Mean);
1234   //Double_t sRef(TMath::Sqrt(s2Mean+s2yTrk)); // reference error parameterization
1235
1236   // organize row candidates
1237   Int_t idxRow[kNrows], nrc(0); Double_t zresRow[kNrows];
1238   for(Int_t ir(0); ir<kNrows; ir++){
1239     idxRow[ir]=-1; zresRow[ir] = 999.;
1240     if(!ncl[ir]) continue;
1241     // get mean z resolution
1242     dz = 0.; for(Int_t ic = ncl[ir]; ic--;) dz += zres[ir][ic]; dz/=ncl[ir];
1243     // insert row
1244     idxRow[nrc] = ir; zresRow[nrc] = TMath::Abs(dz); nrc++;
1245   }
1246   AliDebug(4, Form("Found %d clusters in %d rows. Sorting ...", ncls, nrc));
1247
1248   // sort row candidates
1249   if(nrc>=2){
1250     if(nrc==2){
1251       if(zresRow[0]>zresRow[1]){ // swap
1252         Int_t itmp=idxRow[1]; idxRow[1] = idxRow[0]; idxRow[0] = itmp;
1253         Double_t dtmp=zresRow[1]; zresRow[1] = zresRow[0]; zresRow[0] = dtmp;
1254       }
1255       if(TMath::Abs(idxRow[1] - idxRow[0]) != 1){
1256         SetErrorMsg(kAttachRowGap);
1257         AliDebug(2, Form("Rows attached not continuous. Select first candidate.\n"
1258                     "       row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f",
1259                     idxRow[0], ncl[idxRow[0]], zresRow[0], idxRow[1], idxRow[1]<0?0:ncl[idxRow[1]], zresRow[1]));
1260         nrc=1; idxRow[1] = -1; zresRow[1] = 999.;
1261       }
1262     } else {
1263       Int_t idx0[kNrows];
1264       TMath::Sort(nrc, zresRow, idx0, kFALSE);
1265       nrc = 3; // select only maximum first 3 candidates
1266       Int_t iatmp[] = {-1, -1, -1}; Double_t datmp[] = {999., 999., 999.};
1267       for(Int_t irc(0); irc<nrc; irc++){
1268         iatmp[irc] = idxRow[idx0[irc]];
1269         datmp[irc] = zresRow[idx0[irc]];
1270       }
1271       idxRow[0] = iatmp[0]; zresRow[0] = datmp[0];
1272       idxRow[1] = iatmp[1]; zresRow[1] = datmp[1];
1273       idxRow[2] = iatmp[2]; zresRow[2] = datmp[2]; // temporary
1274       if(TMath::Abs(idxRow[1] - idxRow[0]) != 1){
1275         SetErrorMsg(kAttachRowGap);
1276         AliDebug(2, Form("Rows attached not continuous. Turn on selection.\n"
1277                     "row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f\n"
1278                     "row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f\n"
1279                     "row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f",
1280                     idxRow[0], ncl[idxRow[0]], zresRow[0],
1281                     idxRow[1], ncl[idxRow[1]], zresRow[1],
1282                     idxRow[2], ncl[idxRow[2]], zresRow[2]));
1283         if(TMath::Abs(idxRow[0] - idxRow[2]) == 1){ // select second candidate
1284           AliDebug(2, "Solved ! Remove second candidate.");
1285           nrc = 2;
1286           idxRow[1] = idxRow[2]; zresRow[1] = zresRow[2]; // swap
1287           idxRow[2] = -1; zresRow[2] = 999.;              // remove
1288         } else if(TMath::Abs(idxRow[1] - idxRow[2]) == 1){
1289           if(ncl[idxRow[1]]+ncl[idxRow[2]] > ncl[idxRow[0]]){
1290             AliDebug(2, "Solved ! Remove first candidate.");
1291             nrc = 2;
1292             idxRow[0] = idxRow[1]; zresRow[0] = zresRow[1]; // swap
1293             idxRow[1] = idxRow[2]; zresRow[1] = zresRow[2]; // swap
1294           } else {
1295             AliDebug(2, "Solved ! Remove second and third candidate.");
1296             nrc = 1;
1297             idxRow[1] = -1; zresRow[1] = 999.; // remove
1298             idxRow[2] = -1; zresRow[2] = 999.; // remove
1299           }
1300         } else {
1301           AliDebug(2, "Unsolved !!! Remove second and third candidate.");
1302           nrc = 1;
1303           idxRow[1] = -1; zresRow[1] = 999.; // remove
1304           idxRow[2] = -1; zresRow[2] = 999.; // remove
1305         }
1306       } else { // remove temporary candidate
1307         nrc = 2;
1308         idxRow[2] = -1; zresRow[2] = 999.;
1309       }
1310     }
1311   }
1312   AliDebug(4, Form("Sorted row candidates:\n"
1313       "  row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f"
1314       , idxRow[0], ncl[idxRow[0]], zresRow[0], idxRow[1], idxRow[1]<0?0:ncl[idxRow[1]], zresRow[1]));
1315
1316   // initialize debug streamer
1317   TTreeSRedirector *pstreamer(NULL);
1318   if(recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 3 && fkReconstructor->IsDebugStreaming()) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1319   if(pstreamer){
1320     // save config. for calibration
1321     TVectorD vdy[2], vdx[2], vs2[2];
1322     for(Int_t jr(0); jr<nrc; jr++){
1323       Int_t ir(idxRow[jr]);
1324       vdx[jr].ResizeTo(ncl[ir]); vdy[jr].ResizeTo(ncl[ir]); vs2[jr].ResizeTo(ncl[ir]);
1325       for(Int_t ic(ncl[ir]); ic--;){
1326         vdx[jr](ic) = xres[ir][ic];
1327         vdy[jr](ic) = yres[ir][ic];
1328         vs2[jr](ic) = s2y[ir][ic];
1329       }
1330     }
1331     (*pstreamer) << "AttachClusters4"
1332         << "r0="     << idxRow[0]
1333         << "dz0="    << zresRow[0]
1334         << "dx0="    << &vdx[0]
1335         << "dy0="    << &vdy[0]
1336         << "s20="    << &vs2[0]
1337         << "r1="     << idxRow[1]
1338         << "dz1="    << zresRow[1]
1339         << "dx1="    << &vdx[1]
1340         << "dy1="    << &vdy[1]
1341         << "s21="    << &vs2[1]
1342         << "\n";
1343     vdx[0].Clear(); vdy[0].Clear(); vs2[0].Clear();
1344     vdx[1].Clear(); vdy[1].Clear(); vs2[1].Clear();
1345     if(recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 4){    
1346       Int_t idx(idxRow[1]);
1347       if(idx<0){ 
1348         for(Int_t ir(0); ir<kNrows; ir++){ 
1349           if(clst[ir].GetEntries()>0) continue;
1350           idx = ir;
1351           break;
1352         }
1353       }
1354       (*pstreamer) << "AttachClusters5"
1355           << "c0.="    << &clst[idxRow[0]]
1356           << "c1.="    << &clst[idx]
1357           << "\n";
1358     }
1359   }
1360
1361 //=======================================================================================
1362   // Analyse cluster topology
1363   Double_t f[kNcls],     // likelihood factors for segments
1364            r[2][kNcls],  // d(dydx) of tracklet candidate with respect to track
1365            xm[2][kNcls], // mean <x>
1366            ym[2][kNcls], // mean <y>
1367            sm[2][kNcls], // mean <s_y>
1368            s[2][kNcls],  // sigma_y
1369            p[2][kNcls],  // prob of Gauss
1370            q[2][kNcls];  // charge/segment
1371   memset(f, 0, kNcls*sizeof(Double_t));
1372   Int_t index[2][kNcls], n[2][kNcls];
1373   memset(n, 0, 2*kNcls*sizeof(Int_t));
1374   Int_t mts(0), nts[2] = {0, 0};   // no of tracklet segments in row
1375   AliTRDpadPlane *pp(AliTRDtransform::Geometry().GetPadPlane(fDet));
1376   AliTRDtrackletOflHelper helper;
1377   Int_t lyDet(AliTRDgeometry::GetLayer(fDet));
1378   for(Int_t jr(0), n0(0); jr<nrc; jr++){
1379     Int_t ir(idxRow[jr]);
1380     // cluster segmentation
1381     Bool_t kInit(kFALSE);
1382     if(jr==0){ 
1383       n0 = helper.Init(pp, &clst[ir]); kInit = kTRUE;
1384       if(!n0 || (helper.ClassifyTopology() == AliTRDtrackletOflHelper::kNormal)){
1385         nts[jr] = 1; memset(index[jr], 0, ncl[ir]*sizeof(Int_t));
1386         n[jr][0] = ncl[ir];
1387       }
1388     }
1389     if(!n[jr][0]){
1390       nts[jr] = AliTRDtrackletOflHelper::Segmentation(ncl[ir], xres[ir], yres[ir], index[jr]);
1391       for(Int_t ic(ncl[ir]);ic--;) n[jr][index[jr][ic]]++;
1392     }
1393     mts += nts[jr];
1394     
1395     // tracklet segment processing
1396     for(Int_t its(0); its<nts[jr]; its++){
1397       if(n[jr][its]<=2) {   // don't touch small segments
1398         xm[jr][its] = 0.;ym[jr][its] = 0.;sm[jr][its] = 0.;
1399         for(Int_t ic(ncl[ir]); ic--;){
1400           if(its != index[jr][ic]) continue;
1401           ym[jr][its] += yres[ir][ic];
1402           xm[jr][its] += xres[ir][ic];
1403           sm[jr][its] += TMath::Sqrt(s2y[ir][ic]);
1404         }
1405         if(n[jr][its]==2){ xm[jr][its] *= 0.5; ym[jr][its] *= 0.5; sm[jr][its] *= 0.5;}
1406         xm[jr][its]= fX0 - xm[jr][its];
1407         r[jr][its] = 0.;
1408         s[jr][its] = 1.e-5;
1409         p[jr][its] = 1.;
1410         q[jr][its] = -1.;
1411         continue;
1412       }
1413       
1414       // for longer tracklet segments
1415       if(!kInit) n0 = helper.Init(pp, &clst[ir], index[jr], its);
1416       Int_t n1 = helper.GetRMS(r[jr][its], ym[jr][its], s[jr][its], fX0/*xm[jr][its]*/);
1417       p[jr][its]  = Double_t(n1)/n0;
1418       sm[jr][its] = helper.GetSyMean();
1419       q[jr][its]  = helper.GetQ()/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
1420       xm[jr][its] = fX0;
1421       Double_t dxm= fX0 - xm[jr][its];
1422       yt = fYref[0] - fYref[1]*dxm;
1423       zt = fZref[0] - fZref[1]*dxm;
1424       // correct tracklet fit for tilt
1425       ym[jr][its]+= GetTilt()*(zt - zc[ir]);
1426       r[jr][its] += GetTilt() * fZref[1];
1427       // correct tracklet fit for track position/inclination
1428       ym[jr][its] = yt - ym[jr][its];
1429       r[jr][its]  = (r[jr][its] - fYref[1])/(1+r[jr][its]*fYref[1]);
1430       // report inclination in radians
1431       r[jr][its] = TMath::ATan(r[jr][its]);
1432       if(jr) continue; // calculate only for first row likelihoods
1433         
1434       f[its] = attach->CookLikelihood(chgPos, lyDet, fPt, phiTrk, n[jr][its], ym[jr][its]/*sRef*/, r[jr][its]*TMath::RadToDeg(), s[jr][its]/sm[jr][its]);
1435     }
1436   }
1437   AliDebug(4, Form("   Tracklet candidates: row[%2d] = %2d row[%2d] = %2d:", idxRow[0], nts[0], idxRow[1], nts[1]));
1438   if(AliLog::GetDebugLevel("TRD", "AliTRDseedV1")>3){
1439     for(Int_t jr(0); jr<nrc; jr++){
1440       Int_t ir(idxRow[jr]);
1441       for(Int_t its(0); its<nts[jr]; its++){
1442         printf("  segId[%2d] row[%2d] Ncl[%2d] x[cm]=%7.2f dz[pu]=%4.2f dy[mm]=%+7.3f r[deg]=%+6.2f p[%%]=%6.2f s[um]=%7.2f\n",
1443             its, ir, n[jr][its], xm[jr][its], zresRow[jr], 1.e1*ym[jr][its], r[jr][its]*TMath::RadToDeg(), 100.*p[jr][its], 1.e4*s[jr][its]);
1444       }
1445     }
1446   }
1447   if(!pstreamer && recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 2 && fkReconstructor->IsDebugStreaming()) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1448   if(pstreamer){
1449     // save config. for calibration
1450     TVectorD vidx, vn, vx, vy, vr, vs, vsm, vp, vf;
1451     vidx.ResizeTo(ncl[idxRow[0]]+(idxRow[1]<0?0:ncl[idxRow[1]]));
1452     vn.ResizeTo(mts);
1453     vx.ResizeTo(mts);
1454     vy.ResizeTo(mts);
1455     vr.ResizeTo(mts);
1456     vs.ResizeTo(mts);
1457     vsm.ResizeTo(mts);
1458     vp.ResizeTo(mts);
1459     vf.ResizeTo(mts);
1460     for(Int_t jr(0), jts(0), jc(0); jr<nrc; jr++){
1461        Int_t ir(idxRow[jr]);
1462        for(Int_t its(0); its<nts[jr]; its++, jts++){
1463         vn[jts] = n[jr][its];
1464         vx[jts] = xm[jr][its];
1465         vy[jts] = ym[jr][its];
1466         vr[jts] = r[jr][its];
1467         vs[jts] = s[jr][its];
1468         vsm[jts]= sm[jr][its];
1469         vp[jts] = p[jr][its];
1470         vf[jts] = jr?-1.:f[its];
1471       }
1472       for(Int_t ic(0); ic<ncl[ir]; ic++, jc++) vidx[jc] = index[jr][ic];
1473     }
1474     (*pstreamer) << "AttachClusters3"
1475         << "idx="    << &vidx
1476         << "n="      << &vn
1477         << "x="      << &vx
1478         << "y="      << &vy
1479         << "r="      << &vr
1480         << "s="      << &vs
1481         << "sm="     << &vsm
1482         << "p="      << &vp
1483         << "f="      << &vf
1484         << "\n";
1485   }
1486
1487 //=========================================================
1488   // Get seed tracklet segment
1489   Int_t idx2[kNcls]; memset(idx2, 0, kNcls*sizeof(Int_t)); // seeding indexing
1490   if(nts[0]>1) TMath::Sort(nts[0], f, idx2);
1491   Int_t is(idx2[0]); // seed index
1492   Int_t     idxTrklt[kNcls],
1493             kts(0),
1494             nTrklt(n[0][is]);
1495   Double_t  fTrklt(f[is]),
1496             rTrklt(r[0][is]), 
1497             yTrklt(ym[0][is]), 
1498             sTrklt(s[0][is]), 
1499             smTrklt(sm[0][is]), 
1500             xTrklt(xm[0][is]), 
1501             pTrklt(p[0][is]),
1502             qTrklt(q[0][is]);
1503   memset(idxTrklt, 0, kNcls*sizeof(Int_t));
1504   // check seed idx2[0] exit if not found
1505   if(f[is]<1.e-2){
1506     AliDebug(1, Form("Seed   seg[%d] row[%2d] n[%2d] f[%f]<0.01.", is, idxRow[0], n[0][is], f[is]));
1507     SetErrorMsg(kAttachClAttach);
1508     if(!pstreamer && recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 1 && fkReconstructor->IsDebugStreaming()) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1509     if(pstreamer){
1510       UChar_t stat(0);
1511       if(IsKink()) SETBIT(stat, 1);
1512       if(IsStandAlone()) SETBIT(stat, 2);
1513       if(IsRowCross()) SETBIT(stat, 3);
1514       SETBIT(stat, 4); // set error bit
1515       TVectorD vidx; vidx.ResizeTo(1); vidx[0] = is;
1516       (*pstreamer) << "AttachClusters2"
1517           << "stat="   << stat
1518           << "ev="     << ev
1519           << "chg="    << chgPos
1520           << "det="    << fDet
1521           << "x0="     << fX0
1522           << "y0="     << fYref[0]
1523           << "z0="     << fZref[0]
1524           << "phi="    << phiTrk
1525           << "tht="    << thtTrk
1526           << "pt="     << fPt
1527           << "s2Trk="  << s2yTrk
1528           << "s2Cl="   << s2Mean
1529           << "idx="    << &vidx
1530           << "n="      << nTrklt
1531           << "f="      << fTrklt
1532           << "x="      << xTrklt
1533           << "y="      << yTrklt
1534           << "r="      << rTrklt
1535           << "s="      << sTrklt
1536           << "sm="     << smTrklt
1537           << "p="      << pTrklt
1538           << "q="      << qTrklt
1539           << "\n";
1540     }
1541     return kFALSE;
1542   }
1543   AliDebug(2, Form("Seed   seg[%d] row[%2d] n[%2d] dy[%f] r[%+5.2f] s[%+5.2f] f[%5.3f] q[%6.2f]", is, idxRow[0], n[0][is], ym[0][is], r[0][is]*TMath::RadToDeg(), s[0][is]/sm[0][is], f[is], q[0][is]));
1544
1545   // save seeding segment in the helper
1546   idxTrklt[kts++] = is;
1547   helper.Init(pp, &clst[idxRow[0]], index[0], is);
1548   AliTRDtrackletOflHelper test; // helper to test segment expantion
1549   Float_t rcLikelihood(0.); SetBit(kRowCross, kFALSE);
1550   Double_t dyRez[kNcls]; Int_t idx3[kNcls];
1551   
1552   //=========================================================
1553   // Define filter parameters from OCDB
1554   Int_t kNSgmDy[2]; attach->GetNsgmDy(kNSgmDy[0], kNSgmDy[1]);
1555   Float_t kLikeMinRelDecrease[2]; attach->GetLikeMinRelDecrease(kLikeMinRelDecrease[0], kLikeMinRelDecrease[1]);
1556   Float_t kRClikeLimit(attach->GetRClikeLimit());
1557
1558   //=========================================================
1559   // Try attaching next segments from first row (if any)
1560   if(nts[0]>1){
1561     Int_t jr(0), ir(idxRow[jr]);
1562     // organize  secondary sgms. in decreasing order of their distance from seed 
1563     memset(dyRez, 0, nts[jr]*sizeof(Double_t));
1564     for(Int_t jts(1); jts<nts[jr]; jts++) {
1565       Int_t its(idx2[jts]);
1566       Double_t rot(TMath::Tan(r[0][is]));
1567       dyRez[its] = TMath::Abs(ym[0][is] - ym[jr][its] + rot*(xm[0][is]-xm[jr][its]));
1568     }
1569     TMath::Sort(nts[jr], dyRez, idx3, kFALSE);
1570     for (Int_t jts(1); jts<nts[jr]; jts++) {
1571       Int_t its(idx3[jts]);
1572       if(dyRez[its] > kNSgmDy[jr]*smTrklt){
1573         AliDebug(2, Form("Reject seg[%d] row[%2d] n[%2d] dy[%f] > %d*s[%f].", its, idxRow[jr], n[jr][its], dyRez[its], kNSgmDy[jr], kNSgmDy[jr]*smTrklt));
1574         continue;
1575       }
1576       
1577       test = helper;
1578       Int_t n0 = test.Expand(&clst[ir], index[jr], its);
1579       Double_t rt, dyt, st, xt, smt, pt, qt, ft;
1580       Int_t n1 = test.GetRMS(rt, dyt, st, fX0/*xt*/);
1581       pt = Double_t(n1)/n0;
1582       smt = test.GetSyMean();
1583       qt  = test.GetQ()/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
1584       xt  = fX0;
1585       // correct position
1586       Double_t dxm= fX0 - xt;
1587       yt = fYref[0] - fYref[1]*dxm; 
1588       zt = fZref[0] - fZref[1]*dxm;
1589       // correct tracklet fit for tilt
1590       dyt+= GetTilt()*(zt - zc[idxRow[0]]);
1591       rt += GetTilt() * fZref[1];
1592       // correct tracklet fit for track position/inclination
1593       dyt  = yt - dyt;
1594       rt   = (rt - fYref[1])/(1+rt*fYref[1]);
1595       // report inclination in radians
1596       rt = TMath::ATan(rt);
1597         
1598       ft = (n0>=2) ? attach->CookLikelihood(chgPos, lyDet, fPt, phiTrk, n0,  dyt/*sRef*/, rt*TMath::RadToDeg(), st/smt) : 0.;
1599       Bool_t kAccept(ft>=fTrklt*(1.-kLikeMinRelDecrease[jr]));
1600       
1601       AliDebug(2, Form("%s seg[%d] row[%2d] n[%2d] dy[%f] r[%+5.2f] s[%+5.2f] f[%f] < %4.2f*F[%f].", 
1602         (kAccept?"Adding":"Reject"), its, idxRow[jr], n0, dyt, rt*TMath::RadToDeg(), st/smt, ft, 1.-kLikeMinRelDecrease[jr], fTrklt*(1.-kLikeMinRelDecrease[jr])));
1603       if(kAccept){
1604         idxTrklt[kts++] = its;
1605         nTrklt = n0;
1606         fTrklt = ft;
1607         rTrklt = rt;
1608         yTrklt = dyt;
1609         sTrklt = st;
1610         smTrklt= smt;
1611         xTrklt = xt;
1612         pTrklt = pt;
1613         qTrklt = qt;
1614         helper.Expand(&clst[ir], index[jr], its);
1615       }
1616     }
1617   }
1618   
1619   //=========================================================
1620   // Try attaching next segments from second row (if any)
1621   if(nts[1] && (rcLikelihood = zresRow[0]/zresRow[1]) > kRClikeLimit){
1622     // organize  secondaries in decreasing order of their distance from seed 
1623     Int_t jr(1), ir(idxRow[jr]);
1624     memset(dyRez, 0, nts[jr]*sizeof(Double_t));
1625     Double_t rot(TMath::Tan(r[0][is]));
1626     for(Int_t jts(0); jts<nts[jr]; jts++) {
1627       dyRez[jts] = TMath::Abs(ym[0][is] - ym[jr][jts] + rot*(xm[0][is]-xm[jr][jts]));
1628     }
1629     TMath::Sort(nts[jr], dyRez, idx3, kFALSE);
1630     for (Int_t jts(0); jts<nts[jr]; jts++) {
1631       Int_t its(idx3[jts]);
1632       if(dyRez[its] > kNSgmDy[jr]*smTrklt){
1633         AliDebug(2, Form("Reject seg[%d] row[%2d] n[%2d] dy[%f] > %d*s[%f].", its, idxRow[jr], n[jr][its], dyRez[its], kNSgmDy[jr], kNSgmDy[jr]*smTrklt));
1634         continue;
1635       }
1636       
1637       test = helper;
1638       Int_t n0 = test.Expand(&clst[ir], index[jr], its);
1639       Double_t rt, dyt, st, xt, smt, pt, qt, ft;
1640       Int_t n1 = test.GetRMS(rt, dyt, st, fX0/*xt*/);
1641       pt = Double_t(n1)/n0;
1642       smt = test.GetSyMean();
1643       qt  = test.GetQ()/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
1644       xt  = fX0;
1645       // correct position
1646       Double_t dxm= fX0 - xt;
1647       yt = fYref[0] - fYref[1]*dxm; 
1648       zt = fZref[0] - fZref[1]*dxm;
1649       // correct tracklet fit for tilt
1650       dyt+= GetTilt()*(zt - zc[idxRow[0]]);
1651       rt += GetTilt() * fZref[1];
1652       // correct tracklet fit for track position/inclination
1653       dyt  = yt - dyt;
1654       rt   = (rt - fYref[1])/(1+rt*fYref[1]);
1655       // report inclination in radians
1656       rt = TMath::ATan(rt);
1657         
1658       ft = (n0>=2) ? attach->CookLikelihood(chgPos, lyDet, fPt, phiTrk, n0,  dyt/*sRef*/, rt*TMath::RadToDeg(), st/smt) : 0.;
1659       Bool_t kAccept(ft>=fTrklt*(1.-kLikeMinRelDecrease[jr]));
1660       
1661       AliDebug(2, Form("%s seg[%d] row[%2d] n[%2d] dy[%f] r[%+5.2f] s[%+5.2f] f[%f] < %4.2f*F[%f].", 
1662         (kAccept?"Adding":"Reject"), its, idxRow[jr], n0, dyt, rt*TMath::RadToDeg(), st/smt, ft, 1.-kLikeMinRelDecrease[jr], fTrklt*(1.-kLikeMinRelDecrease[jr])));
1663       if(kAccept){
1664         idxTrklt[kts++] = its;
1665         nTrklt = n0;
1666         fTrklt = ft;
1667         rTrklt = rt;
1668         yTrklt = dyt;
1669         sTrklt = st;
1670         smTrklt= smt;
1671         xTrklt = xt;
1672         pTrklt = pt;
1673         qTrklt = qt;
1674         helper.Expand(&clst[ir], index[jr], its);
1675         SetBit(kRowCross, kTRUE); // mark pad row crossing
1676       }
1677     }
1678   }
1679   // clear local copy of clusters
1680   for(Int_t ir(0); ir<kNrows; ir++) clst[ir].Clear();
1681   
1682   if(!pstreamer && recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 1 && fkReconstructor->IsDebugStreaming()) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1683   if(pstreamer){
1684     UChar_t stat(0);
1685     if(IsKink()) SETBIT(stat, 1);
1686     if(IsStandAlone()) SETBIT(stat, 2);
1687     if(IsRowCross()) SETBIT(stat, 3);
1688     TVectorD vidx; vidx.ResizeTo(kts);
1689     for(Int_t its(0); its<kts; its++) vidx[its] = idxTrklt[its];
1690     (*pstreamer) << "AttachClusters2"
1691         << "stat="   << stat
1692         << "ev="     << ev
1693         << "chg="    << chgPos
1694         << "det="    << fDet
1695         << "x0="     << fX0
1696         << "y0="     << fYref[0]
1697         << "z0="     << fZref[0]
1698         << "phi="    << phiTrk
1699         << "tht="    << thtTrk
1700         << "pt="     << fPt
1701         << "s2Trk="  << s2yTrk
1702         << "s2Cl="   << s2Mean
1703         << "idx="    << &vidx
1704         << "n="      << nTrklt
1705         << "q="      << qTrklt
1706         << "f="      << fTrklt
1707         << "x="      << xTrklt
1708         << "y="      << yTrklt
1709         << "r="      << rTrklt
1710         << "s="      << sTrklt
1711         << "sm="     << smTrklt
1712         << "p="      << pTrklt
1713         << "\n";
1714   }
1715   
1716   
1717   //=========================================================
1718   // Store clusters
1719   Int_t nselected(0), nc(0);
1720   TObjArray *selected(helper.GetClusters());
1721   if(!selected || !(nselected = selected->GetEntriesFast())){
1722     AliError("Cluster candidates missing !!!");
1723     SetErrorMsg(kAttachClAttach);
1724     return kFALSE;
1725   }
1726   for(Int_t ic(0); ic<nselected; ic++){
1727     if(!(c = (AliTRDcluster*)selected->At(ic))) continue;
1728     Int_t it(c->GetPadTime()),
1729           jr(Int_t(helper.GetRow() != c->GetPadRow())),
1730           idx(it+kNtb*jr);
1731     if(fClusters[idx]){
1732       AliDebug(1, Form("Multiple clusters/tb for D[%03d] Tb[%02d] Row[%2d]", fDet, it, c->GetPadRow()));
1733       continue; // already booked
1734     }
1735     // TODO proper indexing of clusters !!
1736     fIndexes[idx]  = chamber->GetTB(it)->GetGlobalIndex(idxs[idxRow[jr]][ic]);
1737     fClusters[idx] = c;
1738     nc++;
1739   }
1740   AliDebug(2, Form("Clusters Found[%2d] Attached[%2d] RC[%c]", nselected, nc, IsRowCross()?'y':'n'));
1741
1742   // number of minimum numbers of clusters expected for the tracklet
1743   if (nc < kClmin){
1744     AliDebug(1, Form("NOT ENOUGH CLUSTERS %d ATTACHED TO THE TRACKLET [min %d] FROM FOUND %d.", nc, kClmin, ncls));
1745     SetErrorMsg(kAttachClAttach);
1746     return kFALSE;
1747   }
1748   SetN(nc);
1749
1750   // Load calibration parameters for this tracklet  
1751   //Calibrate();
1752
1753   // calculate dx for time bins in the drift region (calibration aware)
1754   Float_t x[2] = {0.,0.}; Int_t tb[2]={0,0};
1755   for (Int_t it = t0, irp=0; irp<2 && it < AliTRDtrackerV1::GetNTimeBins(); it++) {
1756     if(!fClusters[it]) continue;
1757     x[irp]  = fClusters[it]->GetX();
1758     tb[irp] = fClusters[it]->GetLocalTimeBin();
1759     irp++;
1760   }  
1761   Int_t dtb = tb[1] - tb[0];
1762   fdX = dtb ? (x[0] - x[1]) / dtb : 0.15;
1763   return kTRUE;
1764 }
1765
1766 //____________________________________________________________
1767 void AliTRDseedV1::Bootstrap(const AliTRDReconstructor *rec)
1768 {
1769 //   Fill in all derived information. It has to be called after recovery from file or HLT.
1770 //   The primitive data are
1771 //   - list of clusters
1772 //   - detector (as the detector will be removed from clusters)
1773 //   - position of anode wire (fX0) - temporary
1774 //   - track reference position and direction
1775 //   - momentum of the track
1776 //   - time bin length [cm]
1777 // 
1778 //   A.Bercuci <A.Bercuci@gsi.de> Oct 30th 2008
1779 //
1780   fkReconstructor = rec;
1781   AliTRDgeometry g;
1782   SetPadPlane(g.GetPadPlane(fDet));
1783
1784   //fSnp = fYref[1]/TMath::Sqrt(1+fYref[1]*fYref[1]);
1785   //fTgl = fZref[1];
1786   Int_t n = 0, nshare = 0, nused = 0;
1787   AliTRDcluster **cit = &fClusters[0];
1788   for(Int_t ic = kNclusters; ic--; cit++){
1789     if(!(*cit)) return;
1790     n++;
1791     if((*cit)->IsShared()) nshare++;
1792     if((*cit)->IsUsed()) nused++;
1793   }
1794   SetN(n); SetNUsed(nused); SetNShared(nshare);
1795   Fit();
1796   CookLabels();
1797   GetProbability();
1798 }
1799
1800
1801 //____________________________________________________________________
1802 Bool_t AliTRDseedV1::Fit(UChar_t opt)
1803 {
1804 //
1805 // Linear fit of the clusters attached to the tracklet
1806 //
1807 // Parameters :
1808 //   - opt : switch for tilt pad correction of cluster y position. Options are
1809 //           0 no correction [default]
1810 //           1 full tilt correction [dz/dx and z0]
1811 //           2 pseudo tilt correction [dz/dx from pad-chamber geometry]
1812 //
1813 // Output :
1814 //  True if successful
1815 //
1816 // Detailed description
1817 //
1818 //            Fit in the xy plane
1819 // 
1820 // The fit is performed to estimate the y position of the tracklet and the track 
1821 // angle in the bending plane. The clusters are represented in the chamber coordinate 
1822 // system (with respect to the anode wire - see AliTRDtrackerV1::FollowBackProlongation() 
1823 // on how this is set). The x and y position of the cluster and also their variances 
1824 // are known from clusterizer level (see AliTRDcluster::GetXloc(), AliTRDcluster::GetYloc(), 
1825 // AliTRDcluster::GetSX() and AliTRDcluster::GetSY()). 
1826 // If gaussian approximation is used to calculate y coordinate of the cluster the position 
1827 // is recalculated taking into account the track angle. The general formula to calculate the 
1828 // error of cluster position in the gaussian approximation taking into account diffusion and track
1829 // inclination is given for TRD by:
1830 // BEGIN_LATEX
1831 // #sigma^{2}_{y} = #sigma^{2}_{PRF} + #frac{x#delta_{t}^{2}}{(1+tg(#alpha_{L}))^{2}} + #frac{x^{2}tg^{2}(#phi-#alpha_{L})tg^{2}(#alpha_{L})}{12}
1832 // END_LATEX
1833 //
1834 // Since errors are calculated only in the y directions, radial errors (x direction) are mapped to y
1835 // by projection i.e.
1836 // BEGIN_LATEX
1837 // #sigma_{x|y} = tg(#phi) #sigma_{x}
1838 // END_LATEX
1839 // and also by the lorentz angle correction
1840 //
1841 //            Fit in the xz plane
1842 //
1843 // The "fit" is performed to estimate the radial position (x direction) where pad row cross happens. 
1844 // If no pad row crossing the z position is taken from geometry and radial position is taken from the xy 
1845 // fit (see below).
1846 // 
1847 // There are two methods to estimate the radial position of the pad row cross:
1848 //   1. leading cluster radial position : Here the lower part of the tracklet is considered and the last 
1849 // cluster registered (at radial x0) on this segment is chosen to mark the pad row crossing. The error 
1850 // of the z estimate is given by :
1851 // BEGIN_LATEX
1852 // #sigma_{z} = tg(#theta) #Delta x_{x_{0}}/12
1853 // END_LATEX
1854 // The systematic errors for this estimation are generated by the following sources:
1855 //   - no charge sharing between pad rows is considered (sharp cross)
1856 //   - missing cluster at row cross (noise peak-up, under-threshold signal etc.).
1857 // 
1858 //   2. charge fit over the crossing point : Here the full energy deposit along the tracklet is considered 
1859 // to estimate the position of the crossing by a fit in the qx plane. The errors in the q directions are 
1860 // parameterized as s_q = q^2. The systematic errors for this estimation are generated by the following sources:
1861 //   - no general model for the qx dependence
1862 //   - physical fluctuations of the charge deposit 
1863 //   - gain calibration dependence
1864 //
1865 //            Estimation of the radial position of the tracklet
1866 //
1867 // For pad row cross the radial position is taken from the xz fit (see above). Otherwise it is taken as the 
1868 // interpolation point of the tracklet i.e. the point where the error in y of the fit is minimum. The error
1869 // in the y direction of the tracklet is (see AliTRDseedV1::GetCovAt()):
1870 // BEGIN_LATEX
1871 // #sigma_{y} = #sigma^{2}_{y_{0}} + 2xcov(y_{0}, dy/dx) + #sigma^{2}_{dy/dx}
1872 // END_LATEX
1873 // and thus the radial position is:
1874 // BEGIN_LATEX
1875 // x = - cov(y_{0}, dy/dx)/#sigma^{2}_{dy/dx}
1876 // END_LATEX
1877 //
1878 //            Estimation of tracklet position error 
1879 //
1880 // The error in y direction is the error of the linear fit at the radial position of the tracklet while in the z 
1881 // direction is given by the cluster error or pad row cross error. In case of no pad row cross this is given by:
1882 // BEGIN_LATEX
1883 // #sigma_{y} = #sigma^{2}_{y_{0}} - 2cov^{2}(y_{0}, dy/dx)/#sigma^{2}_{dy/dx} + #sigma^{2}_{dy/dx}
1884 // #sigma_{z} = Pad_{length}/12
1885 // END_LATEX
1886 // For pad row cross the full error is calculated at the radial position of the crossing (see above) and the error 
1887 // in z by the width of the crossing region - being a matter of parameterization. 
1888 // BEGIN_LATEX
1889 // #sigma_{z} = tg(#theta) #Delta x_{x_{0}}/12
1890 // END_LATEX
1891 // In case of no tilt correction (default in the barrel tracking) the tilt is taken into account by the rotation of
1892 // the covariance matrix. See AliTRDseedV1::GetCovAt() for details.
1893 //
1894 // Author 
1895 // A.Bercuci <A.Bercuci@gsi.de>
1896
1897   if(!fkReconstructor){
1898     AliError("The tracklet needs the reconstruction setup. Please initialize by SetReconstructor().");
1899     return kFALSE;
1900   }
1901   if(!IsCalibrated()) Calibrate();
1902   if(opt>2){
1903     AliWarning(Form("Option [%d] outside range [0, 2]. Using default",opt));
1904     opt=0;
1905   }
1906
1907   const Int_t kClmin = 8;
1908   const Float_t kScalePulls = 10.; // factor to scale y pulls - NOT UNDERSTOOD
1909   // get track direction
1910   Double_t y0   = fYref[0];
1911   Double_t dydx = fYref[1]; 
1912   Double_t z0   = fZref[0];
1913   Double_t dzdx = fZref[1];
1914
1915   AliTRDtrackerV1::AliTRDLeastSquare fitterY;
1916   AliTRDtrackerV1::AliTRDLeastSquare fitterZ;
1917
1918   // book cluster information
1919   Double_t qc[kNclusters], xc[kNclusters], yc[kNclusters], zc[kNclusters], sy[kNclusters];
1920
1921   Bool_t tilt(opt==1)       // full tilt correction
1922         ,pseudo(opt==2)     // pseudo tilt correction
1923         ,rc(IsRowCross())   // row cross candidate 
1924         ,kDZDX(IsPrimary());// switch dzdx calculation for barrel primary tracks
1925   Int_t n(0);   // clusters used in fit 
1926   AliTRDcluster *c(NULL), *cc(NULL), **jc = &fClusters[0];
1927   const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam(); //the dynamic cast in GetRecoParam is slow, so caching the pointer to it
1928
1929   const Char_t *tcName[]={"NONE", "FULL", "HALF"};
1930   AliDebug(2, Form("Options : TC[%s] dzdx[%c]", tcName[opt], kDZDX?'Y':'N'));
1931
1932   
1933   for (Int_t ic=0; ic<kNclusters; ic++, ++jc) {
1934     xc[ic]  = -1.; yc[ic]  = 999.; zc[ic]  = 999.; sy[ic]  = 0.;
1935     if(!(c = (*jc))) continue;
1936     if(!c->IsInChamber()) continue;
1937     // compute pseudo tilt correction
1938     if(kDZDX){ 
1939       fZfit[0] = c->GetZ();
1940       if(rc){
1941         for(Int_t kc=AliTRDseedV1::kNtb; kc<AliTRDseedV1::kNclusters; kc++){
1942           if(!(cc=fClusters[kc])) continue;
1943           if(!cc->IsInChamber()) continue;
1944           fZfit[0] += cc->GetZ(); fZfit[0] *= 0.5;
1945           break;
1946         }
1947       }
1948       fZfit[1] = fZfit[0]/fX0;
1949       if(rc){
1950         fZfit[0] += fZfit[1]*0.5*AliTRDgeometry::CdrHght();
1951         fZfit[1] = fZfit[0]/fX0;
1952       }
1953       kDZDX=kFALSE;
1954     }
1955
1956     Float_t w = 1.;
1957     if(c->GetNPads()>4) w = .5;
1958     if(c->GetNPads()>5) w = .2;
1959
1960     // cluster charge
1961     qc[n]   = TMath::Abs(c->GetQ());
1962     // pad row of leading 
1963
1964     xc[n]   = fX0 - c->GetX();
1965
1966     // Recalculate cluster error based on tracking information
1967     c->SetSigmaY2(fS2PRF, fDiffT, fExB, xc[n], -1./*zcorr?zt:-1.*/, dydx);
1968     c->SetSigmaZ2(fPad[0]*fPad[0]/12.); // for HLT
1969     sy[n]  = TMath::Sqrt(c->GetSigmaY2());
1970
1971     yc[n]  = recoParam->UseGAUS() ? 
1972       c->GetYloc(y0, sy[n], GetPadWidth()): c->GetY();
1973     zc[n]   = c->GetZ();
1974
1975     //optional r-phi correction
1976     //printf("   n[%2d] yc[%7.5f] ", n, yc[n]);
1977     Float_t correction(0.);
1978     if(tilt) correction = fPad[2]*(xc[n]*dzdx + zc[n] - z0);
1979     else if(pseudo) correction = fPad[2]*(xc[n]*fZfit[1] + zc[n]-fZfit[0]);
1980     yc[n]-=correction;
1981     //printf("corr(%s%s)[%7.5f] yc1[%7.5f]\n", (tilt?"TC":""), (zcorr?"PC":""), correction, yc[n]);
1982
1983     AliDebug(5, Form("  tb[%2d] dx[%6.3f] y[%6.2f+-%6.3f]", c->GetLocalTimeBin(), xc[n], yc[n], sy[n]));
1984     fitterY.AddPoint(&xc[n], yc[n], sy[n]);
1985     if(rc) fitterZ.AddPoint(&xc[n], qc[n]*(ic<kNtb?1.:-1.), 1.);
1986     n++;
1987   }
1988
1989   // to few clusters
1990   if (n < kClmin){ 
1991     AliDebug(1, Form("Not enough clusters to fit. Clusters: Attached[%d] Fit[%d].", GetN(), n));
1992     SetErrorMsg(kFitCl);
1993     return kFALSE; 
1994   }
1995   // fit XY
1996   if(!fitterY.Eval()){
1997     AliDebug(1, "Fit Y failed.");
1998     SetErrorMsg(kFitFailedY);
1999     return kFALSE;
2000   }
2001   fYfit[0] = fitterY.GetFunctionParameter(0);
2002   fYfit[1] = -fitterY.GetFunctionParameter(1);
2003   // store covariance
2004   Double_t p[3];
2005   fitterY.GetCovarianceMatrix(p);
2006   fCov[0] = kScalePulls*p[1]; // variance of y0
2007   fCov[1] = kScalePulls*p[2]; // covariance of y0, dydx
2008   fCov[2] = kScalePulls*p[0]; // variance of dydx
2009   // the ref radial position is set at the minimum of 
2010   // the y variance of the tracklet
2011   fX   = -fCov[1]/fCov[2];
2012   fS2Y = fCov[0] +2.*fX*fCov[1] + fX*fX*fCov[2];
2013
2014   Float_t xs=fX+.5*AliTRDgeometry::CamHght();
2015   if(xs < 0. || xs > AliTRDgeometry::CamHght()+AliTRDgeometry::CdrHght()){
2016     AliDebug(1, Form("Ref radial position ouside chamber x[%5.2f].", fX));
2017     SetErrorMsg(kFitFailedY);
2018     return kFALSE;
2019   }
2020
2021 /*    // THE LEADING CLUSTER METHOD for z fit
2022     Float_t xMin = fX0;
2023     Int_t ic=n=kNclusters-1; jc = &fClusters[ic];
2024     AliTRDcluster *c0 =0x0, **kc = &fClusters[kNtb-1];
2025     for(; ic>kNtb; ic--, --jc, --kc){
2026       if((c0 = (*kc)) && c0->IsInChamber() && (xMin>c0->GetX())) xMin = c0->GetX();
2027       if(!(c = (*jc))) continue;
2028       if(!c->IsInChamber()) continue;
2029       zc[kNclusters-1] = c->GetZ(); 
2030       fX = fX0 - c->GetX();
2031     }
2032     fZfit[0] = .5*(zc[0]+zc[kNclusters-1]); fZfit[1] = 0.;
2033     // Error parameterization
2034     fS2Z     = fdX*fZref[1];
2035     fS2Z    *= fS2Z; fS2Z    *= 0.2887; //  1/sqrt(12)*/
2036
2037   // fit QZ
2038   if(opt!=1 && IsRowCross()){
2039     if(!fitterZ.Eval()) SetErrorMsg(kFitFailedZ);
2040     if(!HasError(kFitFailedZ) && TMath::Abs(fitterZ.GetFunctionParameter(1))>1.e-10){ 
2041       // TODO - one has to recalculate xy fit based on
2042       // better knowledge of z position
2043 //       Double_t x = -fitterZ.GetFunctionParameter(0)/fitterZ.GetFunctionParameter(1);
2044 //       Double_t z0 = .5*(zc[0]+zc[n-1]);
2045 //       fZfit[0] = z0 + fZfit[1]*x; 
2046 //       fZfit[1] = fZfit[0]/fX0; 
2047 //       redo fit on xy plane
2048     }
2049     // temporary external error parameterization
2050     fS2Z     = 0.05+0.4*TMath::Abs(fZref[1]); fS2Z *= fS2Z;
2051     // TODO correct formula
2052     //fS2Z     = sigma_x*TMath::Abs(fZref[1]);
2053   } else {
2054     //fZfit[0] = zc[0] + dzdx*0.5*AliTRDgeometry::CdrHght();
2055     fS2Z     = GetPadLength()*GetPadLength()/12.;
2056   }
2057   return kTRUE;
2058 }
2059
2060
2061 //____________________________________________________________________
2062 Bool_t AliTRDseedV1::FitRobust(Bool_t chg)
2063 {
2064 //
2065 // Linear fit of the clusters attached to the tracklet
2066 //
2067 // Author 
2068 // A.Bercuci <A.Bercuci@gsi.de>
2069
2070   TTreeSRedirector *pstreamer(NULL);
2071   const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam();   if(recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 3 && fkReconstructor->IsDebugStreaming()) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
2072
2073   // factor to scale y pulls.
2074   // ideally if error parametrization correct this is 1.
2075   //Float_t lyScaler = 1./(AliTRDgeometry::GetLayer(fDet)+1.);
2076   Float_t kScalePulls = 1.; 
2077   AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
2078   if(!calibration){ 
2079     AliWarning("No access to calibration data");
2080   } else {
2081     // Retrieve the CDB container class with the parametric likelihood
2082     const AliTRDCalTrkAttach *attach = calibration->GetAttachObject();
2083     if(!attach){ 
2084       AliWarning("No usable AttachClusters calib object.");
2085     } else { 
2086       kScalePulls = attach->GetScaleCov();//*lyScaler;
2087     }
2088     // Retrieve chamber status
2089     SetChmbGood(calibration->IsChamberGood(fDet));
2090     if(!IsChmbGood()) kScalePulls*=10.;
2091   }  
2092   Double_t xc[kNclusters], yc[kNclusters], sy[kNclusters];
2093   Int_t n(0),           // clusters used in fit 
2094         row[]={-1, 0}; // pad row spanned by the tracklet
2095   AliTRDcluster *c(NULL), **jc = &fClusters[0];
2096   for(Int_t ic=0; ic<kNtb; ic++, ++jc) {
2097     if(!(c = (*jc))) continue;
2098     if(!c->IsInChamber()) continue;
2099     if(row[0]<0){ 
2100       fZfit[0] = c->GetZ();
2101       fZfit[1] = 0.;
2102       row[0] = c->GetPadRow();
2103     }
2104     xc[n]  = c->GetX();
2105     yc[n]  = c->GetY();
2106     sy[n]  = c->GetSigmaY2()>0?(TMath::Min(TMath::Sqrt(c->GetSigmaY2()), 0.08)):0.08;
2107     n++;
2108   }
2109   Double_t corr = fPad[2]*fPad[0];
2110
2111   for(Int_t ic=kNtb; ic<kNclusters; ic++, ++jc) {
2112     if(!(c = (*jc))) continue;
2113     if(!c->IsInChamber()) continue;
2114     if(row[1]==0) row[1] = c->GetPadRow() - row[0];
2115     xc[n]  = c->GetX();
2116     yc[n]  = c->GetY() + corr*row[1];
2117     sy[n]  = c->GetSigmaY2()>0?(TMath::Min(TMath::Sqrt(c->GetSigmaY2()), 0.08)):0.08;
2118     n++;
2119   }
2120   UChar_t status(0);
2121   Double_t par[3] = {0.,0.,fX0}, cov[3];
2122   if(!AliTRDtrackletOflHelper::Fit(n, xc, yc, sy, par, 1.5, cov)){ 
2123     AliDebug(1, Form("Tracklet fit failed D[%03d].", fDet));
2124     SetErrorMsg(kFitCl);
2125     return kFALSE; 
2126   }
2127   fYfit[0] = par[0];
2128   fYfit[1] = par[1];
2129   // store covariance
2130   fCov[0] = kScalePulls*cov[0]; // variance of y0
2131   fCov[1] = kScalePulls*cov[2]; // covariance of y0, dydx
2132   fCov[2] = kScalePulls*cov[1]; // variance of dydx
2133   // the ref radial position is set at the minimum of 
2134   // the y variance of the tracklet
2135   fX   = 0.;//-fCov[1]/fCov[2];
2136   // check radial position
2137   Float_t xs=fX+.5*AliTRDgeometry::CamHght();
2138   if(xs < 0. || xs > AliTRDgeometry::CamHght()+AliTRDgeometry::CdrHght()){
2139     AliDebug(1, Form("Ref radial position x[%5.2f] ouside D[%3d].", fX, fDet));
2140     SetErrorMsg(kFitFailedY);
2141     return kFALSE;
2142   }
2143   fS2Y = fCov[0] + fX*fCov[1];
2144   fS2Z = fPad[0]*fPad[0]/12.;
2145   AliDebug(2, Form("[I]  x[cm]=%6.2f y[cm]=%+5.2f z[cm]=%+6.2f dydx[deg]=%+5.2f sy[um]=%6.2f sz[cm]=%6.2f", GetX(), GetY(), GetZ(), TMath::ATan(fYfit[1])*TMath::RadToDeg(), TMath::Sqrt(fS2Y)*1.e4, TMath::Sqrt(fS2Z)));
2146   if(IsRowCross()){
2147     Float_t x,z;
2148     if(!GetEstimatedCrossPoint(x,z)){
2149       AliDebug(2, Form("Failed(I) getting crossing point D[%03d].", fDet));
2150       SetErrorMsg(kFitFailedY);
2151       return kTRUE;
2152     }
2153     //if(IsPrimary()){ 
2154       fZfit[0] = fX0*z/x;
2155       fZfit[1] = z/x;
2156       fS2Z     = 0.05+0.4*TMath::Abs(fZfit[1]); fS2Z *= fS2Z;
2157     //}
2158     AliDebug(2, Form("s2y[%f] s2z[%f]", fS2Y, fS2Z));
2159     AliDebug(2, Form("[II] x[cm]=%6.2f y[cm]=%+5.2f z[cm]=%+6.2f dydx[deg]=%+5.2f sy[um]=%6.2f sz[um]=%6.2f dzdx[deg]=%+5.2f", GetX(), GetY(), GetZ(), TMath::ATan(fYfit[1])*TMath::RadToDeg(), TMath::Sqrt(fS2Y)*1.e4, TMath::Sqrt(fS2Z)*1.e4, TMath::ATan(fZfit[1])*TMath::RadToDeg()));
2160   }
2161   
2162   if(pstreamer){
2163     Float_t x= fX0 -fX,
2164             y = GetY(),
2165             yt = fYref[0]-fX*fYref[1];
2166     SETBIT(status, 2);
2167     TVectorD vcov(3); vcov[0]=cov[0];vcov[1]=cov[1];vcov[2]=cov[2];
2168     Double_t sm(0.), chi2(0.), tmp, dy[kNclusters];
2169     for(Int_t ic(0); ic<n; ic++){
2170       sm   += sy[ic];
2171       dy[ic] = yc[ic]-(fYfit[0]+(xc[ic]-fX0)*fYfit[1]); tmp = dy[ic]/sy[ic];
2172       chi2 += tmp*tmp;
2173     }
2174     sm /= n; chi2 = TMath::Sqrt(chi2);
2175     Double_t m(0.), s(0.);
2176     AliMathBase::EvaluateUni(n, dy, m, s, 0);
2177     (*pstreamer) << "FitRobust4"
2178       << "stat=" << status
2179       << "chg="  << chg
2180       << "ncl="  << n
2181       << "det="  << fDet
2182       << "x0="   << fX0
2183       << "y0="   << fYfit[0]
2184       << "x="    << x
2185       << "y="    << y
2186       << "dydx=" << fYfit[1]
2187       << "pt="   << fPt
2188       << "yt="   << yt
2189       << "dydxt="<< fYref[1]
2190       << "cov="  << &vcov
2191       << "chi2=" << chi2
2192       << "sm="   << sm
2193       << "ss="   << s
2194       << "\n";
2195   }
2196   return kTRUE;
2197 }
2198
2199 //___________________________________________________________________
2200 void AliTRDseedV1::Print(Option_t *o) const
2201 {
2202   //
2203   // Printing the seedstatus
2204   //
2205
2206   AliInfo(Form("Det[%3d] X0[%7.2f] Pad{L[%5.2f] W[%5.2f] Tilt[%+6.2f]}", fDet, fX0, GetPadLength(), GetPadWidth(), GetTilt()));
2207   AliInfo(Form("N[%2d] Nused[%2d] Nshared[%2d] [%d]", GetN(), GetNUsed(), GetNShared(), fN));
2208   AliInfo(Form("FLAGS : RC[%c] Kink[%c] SA[%c]", IsRowCross()?'y':'n', IsKink()?'y':'n', IsStandAlone()?'y':'n'));
2209   AliInfo(Form("CALIB PARAMS :  T0[%5.2f]  Vd[%5.2f]  s2PRF[%5.2f]  ExB[%5.2f]  Dl[%5.2f]  Dt[%5.2f]", fT0, fVD, fS2PRF, fExB, fDiffL, fDiffT));
2210
2211   Double_t cov[3], x=GetX();
2212   GetCovAt(x, cov);
2213   AliInfo("    |  x[cm]  |      y[cm]       |      z[cm]      |  dydx |  dzdx |");
2214   AliInfo(Form("Fit | %7.2f | %7.2f+-%7.2f | %7.2f+-%7.2f| %5.2f | ----- |", x, GetY(), TMath::Sqrt(cov[0]), GetZ(), TMath::Sqrt(cov[2]), fYfit[1]));
2215   AliInfo(Form("Ref | %7.2f | %7.2f+-%7.2f | %7.2f+-%7.2f| %5.2f | %5.2f |", x, fYref[0]-fX*fYref[1], TMath::Sqrt(fRefCov[0]), fZref[0]-fX*fYref[1], TMath::Sqrt(fRefCov[2]), fYref[1], fZref[1]));
2216   AliInfo(Form("P / Pt [GeV/c] = %f / %f", GetMomentum(), fPt));
2217   if(IsStandAlone()) AliInfo(Form("C Rieman / Vertex [1/cm] = %f / %f", fC[0], fC[1]));
2218   AliInfo(Form("dEdx [a.u.]    = %f / %f / %f / %f / %f/ %f / %f / %f", fdEdx[0], fdEdx[1], fdEdx[2], fdEdx[3], fdEdx[4], fdEdx[5], fdEdx[6], fdEdx[7]));
2219   AliInfo(Form("PID            = %5.3f / %5.3f / %5.3f / %5.3f / %5.3f", fProb[0], fProb[1], fProb[2], fProb[3], fProb[4]));
2220
2221   if(strcmp(o, "a")!=0) return;
2222
2223   AliTRDcluster* const* jc = &fClusters[0];
2224   for(int ic=0; ic<kNclusters; ic++, jc++) {
2225     if(!(*jc)) continue;
2226     (*jc)->Print(o);
2227   }
2228 }
2229
2230
2231 //___________________________________________________________________
2232 Bool_t AliTRDseedV1::IsEqual(const TObject *o) const
2233 {
2234   // Checks if current instance of the class has the same essential members
2235   // as the given one
2236
2237   if(!o) return kFALSE;
2238   const AliTRDseedV1 *inTracklet = dynamic_cast<const AliTRDseedV1*>(o);
2239   if(!inTracklet) return kFALSE;
2240
2241   for (Int_t i = 0; i < 2; i++){
2242     if ( fYref[i] != inTracklet->fYref[i] ) return kFALSE;
2243     if ( fZref[i] != inTracklet->fZref[i] ) return kFALSE;
2244   }
2245   
2246   if ( TMath::Abs(fS2Y - inTracklet->fS2Y)>1.e-10 ) return kFALSE;
2247   if ( TMath::Abs(GetTilt() - inTracklet->GetTilt())>1.e-10 ) return kFALSE;
2248   if ( TMath::Abs(GetPadLength() - inTracklet->GetPadLength())>1.e-10 ) return kFALSE;
2249   
2250   for (Int_t i = 0; i < kNclusters; i++){
2251 //     if ( fX[i] != inTracklet->GetX(i) ) return kFALSE;
2252 //     if ( fY[i] != inTracklet->GetY(i) ) return kFALSE;
2253 //     if ( fZ[i] != inTracklet->GetZ(i) ) return kFALSE;
2254     if ( fIndexes[i] != inTracklet->fIndexes[i] ) return kFALSE;
2255   }
2256 //   if ( fUsable != inTracklet->fUsable ) return kFALSE;
2257
2258   for (Int_t i=0; i < 2; i++){
2259     if ( fYfit[i] != inTracklet->fYfit[i] ) return kFALSE;
2260     if ( fZfit[i] != inTracklet->fZfit[i] ) return kFALSE;
2261     if ( fLabels[i] != inTracklet->fLabels[i] ) return kFALSE;
2262   }
2263   
2264 /*  if ( fMeanz != inTracklet->GetMeanz() ) return kFALSE;
2265   if ( fZProb != inTracklet->GetZProb() ) return kFALSE;*/
2266   if ( fN != inTracklet->fN ) return kFALSE;
2267   //if ( fNUsed != inTracklet->fNUsed ) return kFALSE;
2268   //if ( fFreq != inTracklet->GetFreq() ) return kFALSE;
2269   //if ( fNChange != inTracklet->GetNChange() ) return kFALSE;
2270    
2271   if ( TMath::Abs(fC[0] - inTracklet->fC[0])>1.e-10 ) return kFALSE;
2272   //if ( fCC != inTracklet->GetCC() ) return kFALSE;
2273   if ( TMath::Abs(fChi2 - inTracklet->fChi2)>1.e-10 ) return kFALSE;
2274   //  if ( fChi2Z != inTracklet->GetChi2Z() ) return kFALSE;
2275
2276   if ( fDet != inTracklet->fDet ) return kFALSE;
2277   if ( TMath::Abs(fPt - inTracklet->fPt)>1.e-10 ) return kFALSE;
2278   if ( TMath::Abs(fdX - inTracklet->fdX)>1.e-10 ) return kFALSE;
2279   
2280   for (Int_t iCluster = 0; iCluster < kNclusters; iCluster++){
2281     AliTRDcluster *curCluster = fClusters[iCluster];
2282     AliTRDcluster *inCluster = inTracklet->fClusters[iCluster];
2283     if (curCluster && inCluster){
2284       if (! curCluster->IsEqual(inCluster) ) {
2285         curCluster->Print();
2286         inCluster->Print();
2287         return kFALSE;
2288       }
2289     } else {
2290       // if one cluster exists, and corresponding 
2291       // in other tracklet doesn't - return kFALSE
2292       if(curCluster || inCluster) return kFALSE;
2293     }
2294   }
2295   return kTRUE;
2296 }
2297