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