1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliTRDseedV1.cxx 60233 2013-01-10 09:04:08Z abercuci $ */
18 ////////////////////////////////////////////////////////////////////////////
20 // The TRD offline tracklet
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)
33 // Alex Bercuci <A.Bercuci@gsi.de> //
34 // Markus Fasel <M.Fasel@gsi.de> //
36 ////////////////////////////////////////////////////////////////////////////
39 #include "TGeoManager.h"
40 #include "TTreeStream.h"
41 #include "TGraphErrors.h"
44 #include "AliMathBase.h"
45 #include "AliRieman.h"
46 #include "AliCDBManager.h"
48 #include "AliTRDReconstructor.h"
49 #include "AliTRDpadPlane.h"
50 #include "AliTRDtransform.h"
51 #include "AliTRDcluster.h"
52 #include "AliTRDseedV1.h"
53 #include "AliTRDtrackV1.h"
54 #include "AliTRDcalibDB.h"
55 #include "AliTRDchamberTimeBin.h"
56 #include "AliTRDtrackingChamber.h"
57 #include "AliTRDtrackerV1.h"
58 #include "AliTRDrecoParam.h"
59 #include "AliTRDCommonParam.h"
60 #include "AliTRDtrackletOflHelper.h"
62 #include "Cal/AliTRDCalTrkAttach.h"
63 #include "Cal/AliTRDCalPID.h"
64 #include "Cal/AliTRDCalROC.h"
65 #include "Cal/AliTRDCalDet.h"
69 ClassImp(AliTRDseedV1)
71 //____________________________________________________________________
72 AliTRDseedV1::AliTRDseedV1(Int_t det)
74 ,fkReconstructor(NULL)
99 memset(fIndexes,0xFF,kNclusters*sizeof(fIndexes[0]));
100 memset(fClusters, 0, kNclusters*sizeof(AliTRDcluster*));
101 memset(fPad, 0, 4*sizeof(Float_t));
102 fYref[0] = 0.; fYref[1] = 0.;
103 fZref[0] = 0.; fZref[1] = 0.;
104 fYfit[0] = 0.; fYfit[1] = 0.;
105 fZfit[0] = 0.; fZfit[1] = 0.;
106 memset(fdEdx, 0, kNdEdxSlices*sizeof(Float_t));
107 for(int ispec=0; ispec<AliPID::kSPECIES; ispec++) fProb[ispec] = -1.;
108 fLabels[0]=-1; fLabels[1]=-1; // most freq MC labels
109 fLabels[2]=0; // number of different labels for tracklet
110 memset(fRefCov, 0, 7*sizeof(Double_t));
111 // stand alone curvature
112 fC[0] = 0.; fC[1] = 0.;
113 // covariance matrix [diagonal]
114 // default sy = 200um and sz = 2.3 cm
115 fCov[0] = 4.e-4; fCov[1] = 0.; fCov[2] = 5.3;
116 SetStandAlone(kFALSE);
119 //____________________________________________________________________
120 AliTRDseedV1::AliTRDseedV1(const AliTRDseedV1 &ref)
121 :AliTRDtrackletBase((AliTRDtrackletBase&)ref)
122 ,fkReconstructor(NULL)
145 // Copy Constructor performing a deep copy
150 SetBit(kOwner, kFALSE);
151 SetStandAlone(ref.IsStandAlone());
155 //____________________________________________________________________
156 AliTRDseedV1& AliTRDseedV1::operator=(const AliTRDseedV1 &ref)
159 // Assignment Operator using the copy function
165 SetBit(kOwner, kFALSE);
170 //____________________________________________________________________
171 AliTRDseedV1::~AliTRDseedV1()
174 // Destructor. The RecoParam object belongs to the underlying tracker.
177 //printf("I-AliTRDseedV1::~AliTRDseedV1() : Owner[%s]\n", IsOwner()?"YES":"NO");
180 for(int itb=0; itb<kNclusters; itb++){
181 if(!fClusters[itb]) continue;
182 //AliInfo(Form("deleting c %p @ %d", fClusters[itb], itb));
183 delete fClusters[itb];
184 fClusters[itb] = NULL;
189 //____________________________________________________________________
190 void AliTRDseedV1::Copy(TObject &ref) const
197 AliTRDseedV1 &target = (AliTRDseedV1 &)ref;
199 target.fkReconstructor = fkReconstructor;
200 target.fClusterIter = NULL;
204 target.fS2PRF = fS2PRF;
205 target.fDiffL = fDiffL;
206 target.fDiffT = fDiffT;
207 target.fClusterIdx = 0;
208 target.fErrorMsg = fErrorMsg;
219 target.fChi2 = fChi2;
221 memcpy(target.fIndexes, fIndexes, kNclusters*sizeof(Int_t));
222 memcpy(target.fClusters, fClusters, kNclusters*sizeof(AliTRDcluster*));
223 memcpy(target.fPad, fPad, 4*sizeof(Float_t));
224 target.fYref[0] = fYref[0]; target.fYref[1] = fYref[1];
225 target.fZref[0] = fZref[0]; target.fZref[1] = fZref[1];
226 target.fYfit[0] = fYfit[0]; target.fYfit[1] = fYfit[1];
227 target.fZfit[0] = fZfit[0]; target.fZfit[1] = fZfit[1];
228 memcpy(target.fdEdx, fdEdx, kNdEdxSlices*sizeof(Float_t));
229 memcpy(target.fProb, fProb, AliPID::kSPECIES*sizeof(Float_t));
230 memcpy(target.fLabels, fLabels, 3*sizeof(Int_t));
231 memcpy(target.fRefCov, fRefCov, 7*sizeof(Double_t));
232 target.fC[0] = fC[0]; target.fC[1] = fC[1];
233 memcpy(target.fCov, fCov, 3*sizeof(Double_t));
239 //____________________________________________________________
240 void AliTRDseedV1::Init(const AliRieman *rieman)
242 // Initialize this tracklet using the riemann fit information
245 fZref[0] = rieman->GetZat(fX0);
246 fZref[1] = rieman->GetDZat(fX0);
247 fYref[0] = rieman->GetYat(fX0);
248 fYref[1] = rieman->GetDYat(fX0);
249 if(fkReconstructor && fkReconstructor->IsHLT()){
253 fRefCov[0] = rieman->GetErrY(fX0);
254 fRefCov[2] = rieman->GetErrZ(fX0);
256 fC[0] = rieman->GetC();
257 fChi2 = rieman->GetChi2();
261 //____________________________________________________________
262 Bool_t AliTRDseedV1::Init(const AliTRDtrackV1 *track)
264 // Initialize this tracklet using the track information
267 // track - the TRD track used to initialize the tracklet
269 // Detailed description
270 // The function sets the starting point and direction of the
271 // tracklet according to the information from the TRD track.
274 // The TRD track has to be propagated to the beginning of the
275 // chamber where the tracklet will be constructed
279 if(!track->GetProlongation(fX0, y, z)) return kFALSE;
285 //_____________________________________________________________________________
286 void AliTRDseedV1::Reset(Option_t *opt)
289 // Reset seed. If option opt="c" is given only cluster arrays are cleared.
291 for(Int_t ic=kNclusters; ic--;) fIndexes[ic] = -1;
292 memset(fClusters, 0, kNclusters*sizeof(AliTRDcluster*));
293 fN=0; SetBit(kRowCross, kFALSE);
294 if(strcmp(opt, "c")==0) return;
296 fExB=0.;fVD=0.;fT0=0.;fS2PRF=0.;
302 fdX=0.;fX0=0.; fX=0.; fY=0.; fZ=0.;
307 memset(fPad, 0, 4*sizeof(Float_t));
308 fYref[0] = 0.; fYref[1] = 0.;
309 fZref[0] = 0.; fZref[1] = 0.;
310 fYfit[0] = 0.; fYfit[1] = 0.;
311 fZfit[0] = 0.; fZfit[1] = 0.;
312 memset(fdEdx, 0, kNdEdxSlices*sizeof(Float_t));
313 for(int ispec=0; ispec<AliPID::kSPECIES; ispec++) fProb[ispec] = -1.;
314 fLabels[0]=-1; fLabels[1]=-1; // most freq MC labels
315 fLabels[2]=0; // number of different labels for tracklet
316 memset(fRefCov, 0, 7*sizeof(Double_t));
317 // covariance matrix [diagonal]
318 // default sy = 200um and sz = 2.3 cm
319 fCov[0] = 4.e-4; fCov[1] = 0.; fCov[2] = 5.3;
322 //____________________________________________________________________
323 void AliTRDseedV1::Update(const AliTRDtrackV1 *trk)
325 // update tracklet reference position from the TRD track
327 Double_t fSnp = trk->GetSnp();
328 Double_t fTgl = trk->GetTgl();
330 Double_t norm =1./TMath::Sqrt((1.-fSnp)*(1.+fSnp));
331 fYref[1] = fSnp*norm;
332 fZref[1] = fTgl*norm;
333 SetCovRef(trk->GetCovariance());
335 Double_t dx = trk->GetX() - fX0;
336 fYref[0] = trk->GetY() - dx*fYref[1];
337 fZref[0] = trk->GetZ() - dx*fZref[1];
340 //_____________________________________________________________________________
341 void AliTRDseedV1::UpdateUsed()
344 // Calculate number of used clusers in the tracklet
347 Int_t nused = 0, nshared = 0;
348 for (Int_t i = kNclusters; i--; ) {
349 if (!fClusters[i]) continue;
350 if(fClusters[i]->IsUsed()){
352 } else if(fClusters[i]->IsShared()){
353 if(IsStandAlone()) nused++;
361 //_____________________________________________________________________________
362 void AliTRDseedV1::UseClusters()
367 // In stand alone mode:
368 // Clusters which are marked as used or shared from another track are
369 // removed from the tracklet
372 // - Clusters which are used by another track become shared
373 // - Clusters which are attached to a kink track become shared
375 AliTRDcluster **c = &fClusters[0];
376 for (Int_t ic=kNclusters; ic--; c++) {
379 if((*c)->IsShared() || (*c)->IsUsed()){
380 if((*c)->IsShared()) SetNShared(GetNShared()-1);
381 else SetNUsed(GetNUsed()-1);
388 if((*c)->IsUsed() || IsKink()){
399 //____________________________________________________________________
400 void AliTRDseedV1::CookdEdx(Int_t nslices)
402 // Calculates average dE/dx for all slices and store them in the internal array fdEdx.
405 // nslices : number of slices for which dE/dx should be calculated
407 // store results in the internal array fdEdx. This can be accessed with the method
408 // AliTRDseedV1::GetdEdx()
410 // Detailed description
411 // Calculates average dE/dx for all slices. Depending on the PID methode
412 // the number of slices can be 3 (LQ) or 8(NN).
413 // The calculation of dQ/dl are done using the tracklet fit results (see AliTRDseedV1::GetdQdl(Int_t))
415 // The following effects are included in the calculation:
416 // 1. calibration values for t0 and vdrift (using x coordinate to calculate slice)
417 // 2. cluster sharing (optional see AliTRDrecoParam::SetClusterSharing())
421 memset(fdEdx, 0, kNdEdxSlices*sizeof(Float_t));
422 const Double_t kDriftLength = (.5 * AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
424 AliTRDcluster *c(NULL);
425 for(int ic=0; ic<AliTRDtrackerV1::GetNTimeBins(); ic++){
426 if(!(c = fClusters[ic]) && !(c = fClusters[ic+kNtb])) continue;
427 Float_t dx = TMath::Abs(fX0 - c->GetX());
429 // Filter clusters for dE/dx calculation
431 // 1.consider calibration effects for slice determination
433 if(dx<kDriftLength){ // TODO should be replaced by c->IsInChamber()
434 slice = Int_t(dx * nslices / kDriftLength);
435 } else slice = c->GetX() < fX0 ? nslices-1 : 0;
438 // 2. take sharing into account
439 Float_t w = /*c->IsShared() ? .5 :*/ 1.;
441 // 3. take into account large clusters TODO
442 //w *= c->GetNPads() > 3 ? .8 : 1.;
445 fdEdx[slice] += w * GetdQdl(ic); //fdQdl[ic];
446 } // End of loop over clusters
449 //_____________________________________________________________________________
450 void AliTRDseedV1::CookLabels()
453 // Cook 2 labels for seed
459 for (Int_t i = 0; i < kNclusters; i++) {
460 if (!fClusters[i]) continue;
461 for (Int_t ilab = 0; ilab < 3; ilab++) {
462 if (fClusters[i]->GetLabel(ilab) >= 0) {
463 labels[nlab] = fClusters[i]->GetLabel(ilab);
469 fLabels[2] = AliMathBase::Freq(nlab,labels,out,kTRUE);
471 if ((fLabels[2] > 1) && (out[3] > 1)) fLabels[1] = out[2];
474 //____________________________________________________________
475 Float_t AliTRDseedV1::GetAnodeWireOffset(Float_t zt)
477 // Find position inside the amplification cell for reading drift velocity map
479 Float_t d = fPad[3] - zt;
481 AliError(Form("Fail AnodeWireOffset calculation z0[%+7.2f] zt[%+7.2f] d[%+7.2f].", fPad[3], zt, d));
484 d -= ((Int_t)(2 * d)) / 2.0;
485 if(d > 0.25) d = 0.5 - d;
490 //____________________________________________________________________
491 Float_t AliTRDseedV1::GetCharge(Bool_t useOutliers) const
493 // Computes total charge attached to tracklet. If "useOutliers" is set clusters
494 // which are not in chamber are also used (default false)
496 AliTRDcluster *c(NULL); Float_t qt(0.);
497 for(int ic=0; ic<kNclusters; ic++){
498 if(!(c=fClusters[ic])) continue;
499 if(!c->IsInChamber() && !useOutliers) continue;
500 qt += TMath::Abs(c->GetQ());
505 //____________________________________________________________________
506 Int_t AliTRDseedV1::GetChargeGaps(Float_t sz[kNtb], Float_t pos[kNtb], Int_t isz[kNtb]) const
508 // Find number, size and position of charge gaps (consecutive missing time bins).
509 // Returns the number of gaps and fills their size in input array "sz" and position in array "pos"
513 Int_t ipos[kNtb]; memset(isz, 0, kNtb*sizeof(Int_t));memset(ipos, 0, kNtb*sizeof(Int_t));
514 for(int ic(0); ic<kNtb; ic++){
515 if(fClusters[ic] || fClusters[ic+kNtb]){
525 // write calibrated values
527 for(Int_t igap(0); igap<n; igap++){
528 sz[igap] = isz[igap]*fVD/AliTRDCommonParam::Instance()->GetSamplingFrequency();
529 fake.SetPadTime(ipos[igap]);
530 pos[igap] = fake.GetXloc(fT0, fVD);
532 fake.SetPadTime(ipos[igap]-isz[igap]+1);
533 pos[igap] += fake.GetXloc(fT0, fVD);
541 //____________________________________________________________________
542 Double_t AliTRDseedV1::EstimatedCrossPoint(AliTRDpadPlane *pp, Float_t bz)
544 // Algorithm to estimate cross point in the x-z plane for pad row cross tracklets or the z coordinate of pad row without pad row cross in the local chamber coordinates.
545 // Returns variance of the radial offset from anode wire in case of raw cross or 0 otherwise.
547 Int_t row[] = {-1, -1};
548 Double_t zoff(0.5 * (pp->GetRow0() + pp->GetRowEnd())), sx(0.), mean(0.5*pp->GetNrows()-0.5);
549 AliTRDcluster *c(NULL);
553 for(int ic=0; ic<kNtb; ic++){
554 if(!(c=fClusters[ic])) continue;
555 if(!c->IsInChamber()) continue;
556 row[0] = c->GetPadRow();
557 fZfit[0] = Int_t(mean-row[0])*pp->GetLengthIPad() +
558 0.5*(mean-row[0]>0.?1.:-1.)*(row[0]>0&&row[0]<pp->GetNrows()-1?pp->GetLengthIPad():pp->GetLengthOPad());
562 Float_t tbm[2] = {0.}; // mean value of time bin in rows
563 Int_t tb[kNtb]={0}, //array of time bins from first row
564 nc[2] = {0}, // no. of clusters in rows
565 mc(0); // no. of common clusters
566 Bool_t w[2] = {kFALSE, kFALSE}; // acceptance flag for rows
567 // Find radial range for first row
568 for(int ic(0); ic<kNtb; ic++){
570 if(!(c=fClusters[ic]) || !c->IsInChamber()) continue;
571 if(row[0]<0) row[0] = c->GetPadRow();
572 tb[nc[0]++] = ic; tbm[0] += ic;
578 // Find radial range for second row
579 for(int ic(kNtb), jc(0); ic<kNclusters; ic++, jc++){
580 if(!(c=fClusters[ic]) || !c->IsInChamber()) continue;
581 if(row[1]<0) row[1] = c->GetPadRow();
582 tbm[1] += jc; nc[1]++;
583 for(Int_t kc(0); kc<nc[0]; kc++)
585 tb[kc] += 100; // mark common cluster
594 //printf("0 : %f[%2d] 1 : %f[%2d] mc[%d]\n", tbm[0], nc[0], tbm[1], nc[1], mc);
596 AliError("Too few clusters to estimate tracklet.");
600 SetBit(kRowCross, kFALSE); // reset RC bit
601 if(w[1]) row[0] = row[1];
602 fZfit[0] = Int_t(mean-row[0])*pp->GetLengthIPad() +
603 0.5*(mean-row[0]>0.?1.:-1.)*(row[0]>0&&row[0]<pp->GetNrows()-1?pp->GetLengthIPad():pp->GetLengthOPad());
604 }else{ // find the best matching timebin
605 fZfit[0] = Int_t(mean-0.5*(row[0]+row[1]))*pp->GetLengthIPad();
606 Int_t itb(0), dtb(0);
607 if(!mc) { // no common range
608 itb = Int_t(0.5*(tbm[0] + tbm[1]));
609 dtb = Int_t(0.5*TMath::Abs(tbm[0] - tbm[1])); // simple parameterization of the cluster gap
611 Double_t rmax(100.); Int_t itbStart(-1), itbStop(0);
612 // compute distance from
613 for(Int_t jc(0); jc<nc[0]; jc++){
614 if(tb[jc] < 100) continue;
615 Int_t ltb(tb[jc]-100);
616 Double_t r = (1. - ltb/tbm[0])*(1. - ltb/tbm[1]);
617 //printf("tb[%2d] dr[%f %f %f] rmax[%f]\n", ltb, r, 1. - ltb/tbm[0], 1. - ltb/tbm[1], rmax);
618 if(TMath::Abs(r)<rmax){ rmax = TMath::Abs(r); itb = ltb; }
619 if(itbStart<0) itbStart = ltb;
622 dtb = itbStop-itbStart+1;
624 AliTRDCommonParam *cp = AliTRDCommonParam::Instance();
625 Double_t freq(cp?cp->GetSamplingFrequency():10.);
626 fS2Y = ((itb-0.5)/freq - fT0 - 0.189)*fVD; // xOff
627 sx = dtb*0.288675134594812921/freq; sx *= sx; sx += 1.56e-2; sx *= fVD*fVD;
632 Float_t dx(fX0-fS2Y);
633 fZfit[1] = (fZfit[0]+zoff)/dx;
635 // correct dzdx for the bias
636 UnbiasDZDX(IsRowCross(), bz);
638 // correct x_cross/sigma(x_cross) for the bias in dzdx
639 const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam();
641 fS2Y += recoParam->GetCorrDZDXxcross()*TMath::Abs(fZfit[1]);
642 sx += recoParam->GetCorrDZDXxcross()*recoParam->GetCorrDZDXxcross()*GetS2DZDX(fZfit[1]);
644 // correct sigma(x_cross) for the width of the crossing area
645 sx += GetS2XcrossDZDX(TMath::Abs(fZfit[1]));
647 // estimate z and error @ anode wire
648 fZfit[0] += fZfit[1]*fS2Y;
649 fS2Z = fZfit[1]*fZfit[1]*sx+fS2Y*fS2Y*GetS2DZDX(fZfit[1]);
654 //____________________________________________________________________
655 void AliTRDseedV1::UnbiasDZDX(Bool_t rc, Float_t bz)
657 // correct dzdx for the bias in z according to MC
658 const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam();
659 if(!recoParam) return;
660 fZfit[1] *= recoParam->GetCorrDZDX(rc)-(bz>0?0.01:0.);
661 if(rc) fZfit[1] += recoParam->GetCorrDZDXbiasRC(fZfit[1]<0);
664 //____________________________________________________________________
665 Double_t AliTRDseedV1::UnbiasY(Bool_t rc, Float_t bz)
667 // correct y coordinate for tail cancellation. This should be fixed by considering TC as a function of q/pt.
668 // rc : TRUE if tracklet crosses rows
669 // bz : magnetic field z component
671 const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam();
672 if(!recoParam) return 0.;
673 Double_t par[3]={0.};
674 Int_t idx(2*(rc?1:0)+Int_t(bz>0));
675 recoParam->GetYcorrTailCancel(idx, par);
676 return par[0]*TMath::Sin(par[1]*fYref[1])+par[2];
680 //____________________________________________________________________
681 Float_t AliTRDseedV1::GetQperTB(Int_t tb) const
684 // Charge of the clusters at timebin
687 if(fClusters[tb] /*&& fClusters[tb]->IsInChamber()*/)
688 q += TMath::Abs(fClusters[tb]->GetQ());
689 if(fClusters[tb+kNtb] /*&& fClusters[tb+kNtb]->IsInChamber()*/)
690 q += TMath::Abs(fClusters[tb+kNtb]->GetQ());
691 return q/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
694 //____________________________________________________________________
695 Float_t AliTRDseedV1::GetdQdl() const
697 // Calculate total charge / tracklet length for 1D PID
699 Float_t Q = GetCharge(kTRUE);
700 return Q/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
703 //____________________________________________________________________
704 Float_t AliTRDseedV1::GetdQdl(Int_t ic, Float_t *dl) const
706 // Using the linear approximation of the track inside one TRD chamber (TRD tracklet)
707 // the charge per unit length can be written as:
709 // #frac{dq}{dl} = #frac{q_{c}}{dx * #sqrt{1 + #(){#frac{dy}{dx}}^{2}_{fit} + #(){#frac{dz}{dx}}^{2}_{ref}}}
711 // where qc is the total charge collected in the current time bin and dx is the length
713 // The following correction are applied :
714 // - charge : pad row cross corrections
715 // [diffusion and TRF assymetry] TODO
716 // - dx : anisochronity, track inclination - see Fit and AliTRDcluster::GetXloc()
717 // and AliTRDcluster::GetYloc() for the effects taken into account
720 //<img src="TRD/trackletDQDT.gif">
722 // In the picture the energy loss measured on the tracklet as a function of drift time [left] and respectively
723 // drift length [right] for different particle species is displayed.
724 // Author : Alex Bercuci <A.Bercuci@gsi.de>
727 // check whether both clusters are inside the chamber
728 Bool_t hasClusterInChamber = kFALSE;
729 if(fClusters[ic] && fClusters[ic]->IsInChamber()){
730 hasClusterInChamber = kTRUE;
731 dq += TMath::Abs(fClusters[ic]->GetQ());
733 if(fClusters[ic+kNtb] && fClusters[ic+kNtb]->IsInChamber()){
734 hasClusterInChamber = kTRUE;
735 dq += TMath::Abs(fClusters[ic+kNtb]->GetQ());
737 if(!hasClusterInChamber) return 0.;
738 if(dq<1.e-3) return 0.;
741 if(ic-1>=0 && ic+1<kNtb){
742 Float_t x2(0.), x1(0.);
743 // try to estimate upper radial position (find the cluster which is inside the chamber)
744 if(fClusters[ic-1] && fClusters[ic-1]->IsInChamber()) x2 = fClusters[ic-1]->GetX();
745 else if(fClusters[ic-1+kNtb] && fClusters[ic-1+kNtb]->IsInChamber()) x2 = fClusters[ic-1+kNtb]->GetX();
746 else if(fClusters[ic] && fClusters[ic]->IsInChamber()) x2 = fClusters[ic]->GetX()+fdX;
747 else x2 = fClusters[ic+kNtb]->GetX()+fdX;
748 // try to estimate lower radial position (find the cluster which is inside the chamber)
749 if(fClusters[ic+1] && fClusters[ic+1]->IsInChamber()) x1 = fClusters[ic+1]->GetX();
750 else if(fClusters[ic+1+kNtb] && fClusters[ic+1+kNtb]->IsInChamber()) x1 = fClusters[ic+1+kNtb]->GetX();
751 else if(fClusters[ic] && fClusters[ic]->IsInChamber()) x1 = fClusters[ic]->GetX()-fdX;
752 else x1 = fClusters[ic+kNtb]->GetX()-fdX;
756 dx *= TMath::Sqrt(1. + fYfit[1]*fYfit[1] + fZref[1]*fZref[1]);
758 if(dx>1.e-9) return dq/dx;
762 //____________________________________________________________
763 Float_t AliTRDseedV1::GetMomentum(Float_t *err) const
765 // Returns momentum of the track after update with the current tracklet as:
767 // p=#frac{1}{1/p_{t}} #sqrt{1+tgl^{2}}
769 // and optionally the momentum error (if err is not null).
770 // The estimated variance of the momentum is given by:
772 // #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})
774 // which can be simplified to
776 // #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}
780 Double_t p = fPt*TMath::Sqrt(1.+fZref[1]*fZref[1]);
783 Double_t tgl2 = fZref[1]*fZref[1];
784 Double_t pt2 = fPt*fPt;
786 p2*tgl2*pt2*pt2*fRefCov[4]
787 -2.*p2*fZref[1]*fPt*pt2*fRefCov[5]
789 (*err) = TMath::Sqrt(s2);
795 //____________________________________________________________________
796 Int_t AliTRDseedV1::GetTBoccupancy() const
798 // Returns no. of TB occupied by clusters
801 for(int ic(0); ic<kNtb; ic++){
802 if(!fClusters[ic] && !fClusters[ic+kNtb]) continue;
808 //____________________________________________________________________
809 Int_t AliTRDseedV1::GetTBcross() const
811 // Returns no. of TB occupied by 2 clusters for pad row cross tracklets
813 if(!IsRowCross()) return 0;
815 for(int ic(0); ic<kNtb; ic++){
816 if(fClusters[ic] && fClusters[ic+kNtb]) n++;
821 //____________________________________________________________________
822 Float_t* AliTRDseedV1::GetProbability(Bool_t force)
824 if(!force) return &fProb[0];
825 if(!CookPID()) return NULL;
829 //____________________________________________________________
830 Bool_t AliTRDseedV1::CookPID()
832 // Fill probability array for tracklet from the DB.
837 // returns pointer to the probability array and NULL if missing DB access
839 // Retrieve PID probabilities for e+-, mu+-, K+-, pi+- and p+- from the DB according to tracklet information:
840 // - estimated momentum at tracklet reference point
841 // - dE/dx measurements
844 // According to the steering settings specified in the reconstruction one of the following methods are used
845 // - Neural Network [default] - option "nn"
846 // - 2D Likelihood - option "!nn"
848 AliWarning(Form("Obsolete function. Use AliTRDPIDResponse::GetResponse() instead."));
850 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
852 AliError("No access to calibration data");
856 if (!fkReconstructor) {
857 AliError("Reconstructor not set.");
861 // Retrieve the CDB container class with the parametric detector response
862 const AliTRDCalPID *pd = calibration->GetPIDObject(fkReconstructor->GetPIDMethod());
864 AliError("No access to AliTRDCalPID object");
868 // calculate tracklet length TO DO
869 Float_t length = (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick())/ TMath::Sqrt((1.0 - GetSnp()*GetSnp()) / (1.0 + GetTgl()*GetTgl()));
872 CookdEdx(AliTRDCalPID::kNSlicesNN);
873 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));
875 // Sets the a priori probabilities
876 Bool_t kPIDNN(fkReconstructor->GetPIDMethod()==AliTRDpidUtil::kNN);
877 for(int ispec=0; ispec<AliPID::kSPECIES; ispec++)
878 fProb[ispec] = pd->GetProbability(ispec, GetMomentum(), &fdEdx[0], length, kPIDNN?GetPlane():fkReconstructor->GetRecoParam()->GetPIDLQslices());
883 //____________________________________________________________________
884 Float_t AliTRDseedV1::GetQuality(Bool_t kZcorr) const
887 // Returns a quality measurement of the current seed
890 Float_t zcorr = kZcorr ? GetTilt() * (fZfit[0] - fZref[0]) : 0.;
892 .5 * TMath::Abs(18.0 - GetN())
893 + 10.* TMath::Abs(fYfit[1] - fYref[1])
894 + 5. * TMath::Abs(fYfit[0] - fYref[0] + zcorr)
895 + 2. * TMath::Abs(fZfit[0] - fZref[0]) / GetPadLength();
898 //____________________________________________________________________
899 void AliTRDseedV1::GetCovAt(Double_t /*x*/, Double_t *cov) const
901 // Computes covariance in the y-z plane at radial point x (in tracking coordinates)
902 // and returns the results in the preallocated array cov[3] as :
909 // For the linear transformation
913 // The error propagation has the general form
915 // C_{Y} = T_{x} C_{X} T_{x}^{T}
917 // We apply this formula 2 times. First to calculate the covariance of the tracklet
918 // at point x we consider:
920 // T_{x} = (1 x); X=(y0 dy/dx); C_{X}=#(){#splitline{Var(y0) Cov(y0, dy/dx)}{Cov(y0, dy/dx) Var(dy/dx)}}
922 // and secondly to take into account the tilt angle
924 // T_{#alpha} = #(){#splitline{cos(#alpha) __ sin(#alpha)}{-sin(#alpha) __ cos(#alpha)}}; X=(y z); C_{X}=#(){#splitline{Var(y) 0}{0 Var(z)}}
927 // using simple trigonometrics one can write for this last case
929 // 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})}}
931 // which can be aproximated for small alphas (2 deg) with
933 // 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}}}
936 // before applying the tilt rotation we also apply systematic uncertainties to the tracklet
937 // position which can be tunned from outside via the AliTRDrecoParam::SetSysCovMatrix(). They might
938 // account for extra misalignment/miscalibration uncertainties.
941 // Alex Bercuci <A.Bercuci@gsi.de>
942 // Date : Jan 8th 2009
946 //Double_t xr = fX0-x;
947 Double_t sy2 = fCov[0];// +2.*xr*fCov[1] + xr*xr*fCov[2];
949 //GetPadLength()*GetPadLength()/12.;
951 // insert systematic uncertainties
953 Double_t sys[15]; memset(sys, 0, 15*sizeof(Double_t));
954 fkReconstructor->GetRecoParam()->GetSysCovMatrix(sys);
959 // rotate covariance matrix if no RC
961 Double_t t2 = GetTilt()*GetTilt();
962 Double_t correction = 1./(1. + t2);
963 cov[0] = (sy2+t2*sz2)*correction;
964 cov[1] = GetTilt()*(sz2 - sy2)*correction;
965 cov[2] = (t2*sy2+sz2)*correction;
967 cov[0] = sy2; cov[1] = 0.; cov[2] = sz2;
970 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'));
973 //____________________________________________________________
974 Int_t AliTRDseedV1::GetCovSqrt(const Double_t * const c, Double_t *d)
976 // Helper function to calculate the square root of the covariance matrix.
977 // The input matrix is stored in the vector c and the result in the vector d.
978 // Both arrays have to be initialized by the user with at least 3 elements. Return negative in case of failure.
980 // For calculating the square root of the symmetric matrix c
981 // the following relation is used:
983 // C^{1/2} = VD^{1/2}V^{-1}
985 // with V being the matrix with the n eigenvectors as columns.
986 // In case C is symmetric the followings are true:
987 // - matrix D is diagonal with the diagonal given by the eigenvalues of C
990 // Author A.Bercuci <A.Bercuci@gsi.de>
993 const Double_t kZero(1.e-20);
994 Double_t l[2], // eigenvalues
995 v[3]; // eigenvectors
996 // the secular equation and its solution :
997 // (c[0]-L)(c[2]-L)-c[1]^2 = 0
998 // L^2 - L*Tr(c)+DET(c) = 0
999 // L12 = [Tr(c) +- sqrt(Tr(c)^2-4*DET(c))]/2
1000 Double_t tr = c[0]+c[2], // trace
1001 det = c[0]*c[2]-c[1]*c[1]; // determinant
1002 if(TMath::Abs(det)<kZero) return 1;
1003 Double_t dd = TMath::Sqrt(tr*tr - 4*det);
1004 l[0] = .5*(tr + dd*(c[0]>c[2]?-1.:1.));
1005 l[1] = .5*(tr + dd*(c[0]>c[2]?1.:-1.));
1006 if(l[0]<kZero || l[1]<kZero) return 2;
1010 Double_t den = (l[0]-c[0])*(l[0]-c[0])+c[1]*c[1];
1011 if(den<kZero){ // almost diagonal
1012 v[0] = TMath::Sign(0., c[1]);
1013 v[1] = TMath::Sign(1., (l[0]-c[0]));
1014 v[2] = TMath::Sign(0., c[1]*(l[0]-c[0])*(l[1]-c[2]));
1016 Double_t tmp = 1./TMath::Sqrt(den);
1018 v[1] = (l[0]-c[0])*tmp;
1019 if(TMath::Abs(l[1]-c[2])<kZero) v[2] = TMath::Sign(v[0]*(l[0]-c[0])/kZero, (l[1]-c[2]));
1020 else v[2] = v[0]*(l[0]-c[0])/(l[1]-c[2]);
1022 // the VD^{1/2}V is:
1023 l[0] = TMath::Sqrt(l[0]); l[1] = TMath::Sqrt(l[1]);
1024 d[0] = v[0]*v[0]*l[0]+v[1]*v[1]*l[1];
1025 d[1] = v[0]*v[1]*l[0]+v[1]*v[2]*l[1];
1026 d[2] = v[1]*v[1]*l[0]+v[2]*v[2]*l[1];
1031 //____________________________________________________________
1032 Double_t AliTRDseedV1::GetCovInv(const Double_t * const c, Double_t *d)
1034 // Helper function to calculate the inverse of the covariance matrix.
1035 // The input matrix is stored in the vector c and the result in the vector d.
1036 // Both arrays have to be initialized by the user with at least 3 elements
1037 // The return value is the determinant or 0 in case of singularity.
1039 // Author A.Bercuci <A.Bercuci@gsi.de>
1042 Double_t det = c[0]*c[2] - c[1]*c[1];
1043 if(TMath::Abs(det)<1.e-20) return 0.;
1044 Double_t invDet = 1./det;
1051 //____________________________________________________________________
1052 UShort_t AliTRDseedV1::GetVolumeId() const
1054 // Returns geometry volume id by delegation
1056 for(Int_t ic(0);ic<kNclusters; ic++){
1057 if(fClusters[ic]) return fClusters[ic]->GetVolumeId();
1063 //____________________________________________________________________
1064 void AliTRDseedV1::Calibrate()
1066 // Retrieve calibration and position parameters from OCDB.
1067 // The following information are used
1069 // - column and row position of first attached cluster. If no clusters are attached
1070 // to the tracklet a random central chamber position (c=70, r=7) will be used.
1072 // The following information is cached in the tracklet
1073 // t0 (trigger delay)
1076 // omega*tau = tg(a_L)
1077 // diffusion coefficients (longitudinal and transversal)
1080 // Alex Bercuci <A.Bercuci@gsi.de>
1081 // Date : Jan 8th 2009
1084 AliCDBManager *cdb = AliCDBManager::Instance();
1085 if(cdb->GetRun() < 0){
1086 AliError("OCDB manager not properly initialized");
1090 AliTRDcalibDB *calib = AliTRDcalibDB::Instance();
1091 AliTRDCalROC *vdROC = calib->GetVdriftROC(fDet),
1092 *t0ROC = calib->GetT0ROC(fDet);;
1093 const AliTRDCalDet *vdDet = calib->GetVdriftDet();
1094 const AliTRDCalDet *t0Det = calib->GetT0Det();
1096 Int_t col = 70, row = 7;
1097 AliTRDcluster **c = &fClusters[0];
1100 while (ic<kNclusters && !(*c)){ic++; c++;}
1102 col = (*c)->GetPadCol();
1103 row = (*c)->GetPadRow();
1107 fT0 = (t0Det->GetValue(fDet) + t0ROC->GetValue(col,row)) / AliTRDCommonParam::Instance()->GetSamplingFrequency();
1108 fVD = vdDet->GetValue(fDet) * vdROC->GetValue(col, row);
1109 fS2PRF = calib->GetPRFWidth(fDet, col, row); fS2PRF *= fS2PRF;
1110 fExB = AliTRDCommonParam::Instance()->GetOmegaTau(fVD);
1111 AliTRDCommonParam::Instance()->GetDiffCoeff(fDiffL,
1113 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));
1116 SetBit(kCalib, kTRUE);
1119 //____________________________________________________________________
1120 void AliTRDseedV1::SetOwner()
1122 //AliInfo(Form("own [%s] fOwner[%s]", own?"YES":"NO", fOwner?"YES":"NO"));
1124 if(TestBit(kOwner)) return;
1125 for(int ic=0; ic<kNclusters; ic++){
1126 if(!fClusters[ic]) continue;
1127 fClusters[ic] = new AliTRDcluster(*fClusters[ic]);
1132 //____________________________________________________________
1133 void AliTRDseedV1::SetPadPlane(AliTRDpadPlane * const p)
1135 // Shortcut method to initialize pad geometry.
1136 fPad[0] = p->GetLengthIPad();
1137 fPad[1] = p->GetWidthIPad();
1138 fPad[2] = TMath::Tan(TMath::DegToRad()*p->GetTiltingAngle());
1139 fPad[3] = p->GetRow0() + p->GetAnodeWireOffset();
1144 //____________________________________________________________________
1145 Bool_t AliTRDseedV1::AttachClusters(AliTRDtrackingChamber *const chamber, Bool_t tilt, Bool_t chgPos, Int_t ev)
1148 // Projective algorithm to attach clusters to seeding tracklets. The following steps are performed :
1149 // 1. Collapse x coordinate for the full detector plane
1150 // 2. truncated mean on y (r-phi) direction
1151 // 3. purge clusters
1152 // 4. truncated mean on z direction
1153 // 5. purge clusters
1156 // - chamber : pointer to tracking chamber container used to search the tracklet
1157 // - tilt : switch for tilt correction during road building [default true]
1158 // - chgPos : mark same[kFALSE] and opposite[kTRUE] sign tracks with respect to Bz field sign [default true]
1159 // - ev : event number for debug purposes [default = -1]
1161 // - true : if tracklet found successfully. Failure can happend because of the following:
1163 // Detailed description
1165 // We start up by defining the track direction in the xy plane and roads. The roads are calculated based
1166 // on tracking information (variance in the r-phi direction) and estimated variance of the standard
1167 // clusters (see AliTRDcluster::SetSigmaY2()) corrected for tilt (see GetCovAt()). From this the road is
1169 // 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})})}
1170 // r_{z} = 1.5*L_{pad}
1173 // Author : Alexandru Bercuci <A.Bercuci@gsi.de>
1174 // Debug : level = 2 for calibration
1175 // level = 3 for visualization in the track SR
1176 // level = 4 for full visualization including digit level
1178 const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam(); //the dynamic cast in GetRecoParam is slow, so caching the pointer to it
1181 AliError("Tracklets can not be used without a valid RecoParam.");
1184 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1186 AliError("No access to calibration data");
1189 // Retrieve the CDB container class with the parametric likelihood
1190 const AliTRDCalTrkAttach *attach = calibration->GetAttachObject();
1192 AliError("No usable AttachClusters calib object.");
1196 // Initialize reco params for this tracklet
1197 // 1. first time bin in the drift region
1199 Int_t kClmin = Int_t(recoParam->GetFindableClusters()*AliTRDtrackerV1::GetNTimeBins());
1202 Double_t sysCov[5]; recoParam->GetSysCovMatrix(sysCov);
1203 Double_t s2yTrk= fRefCov[0],
1205 s2zCl = GetPadLength()*GetPadLength()/12.,
1206 syRef = TMath::Sqrt(s2yTrk),
1207 t2 = GetTilt()*GetTilt();
1209 const Double_t kroady = 3.; //recoParam->GetRoad1y();
1210 const Double_t kroadz = GetPadLength() * recoParam->GetRoadzMultiplicator() + 1.;
1211 // define probing cluster (the perfect cluster) and default calibration
1212 Short_t sig[] = {0, 0, 10, 30, 10, 0,0};
1213 AliTRDcluster cp(fDet, 6, 75, 0, sig, 0);
1214 if(fkReconstructor->IsHLT()) cp.SetRPhiMethod(AliTRDcluster::kCOG);
1215 if(!IsCalibrated()) Calibrate();
1217 /* Int_t kroadyShift(0);
1218 Float_t bz(AliTrackerBase::GetBz());
1219 if(TMath::Abs(bz)>2.){
1220 if(bz<0.) kroadyShift = chgPos ? +1 : -1;
1221 else kroadyShift = chgPos ? -1 : +1;
1223 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'));
1224 Double_t phiTrk(TMath::ATan(fYref[1])),
1225 thtTrk(TMath::ATan(fZref[1]));
1227 // working variables
1228 const Int_t kNrows = 16;
1229 const Int_t kNcls = 3*kNclusters; // buffer size
1230 TObjArray clst[kNrows];
1231 Bool_t blst[kNrows][kNcls];
1236 xres[kNrows][kNcls], yres[kNrows][kNcls], zres[kNrows][kNcls], s2y[kNrows][kNcls];
1237 Int_t idxs[kNrows][kNcls], ncl[kNrows], ncls = 0;
1238 memset(ncl, 0, kNrows*sizeof(Int_t));
1239 memset(zc, 0, kNrows*sizeof(Double_t));
1240 memset(idxs, 0, kNrows*kNcls*sizeof(Int_t));
1241 memset(xres, 0, kNrows*kNcls*sizeof(Double_t));
1242 memset(yres, 0, kNrows*kNcls*sizeof(Double_t));
1243 memset(zres, 0, kNrows*kNcls*sizeof(Double_t));
1244 memset(s2y, 0, kNrows*kNcls*sizeof(Double_t));
1245 memset(blst, 0, kNrows*kNcls*sizeof(Bool_t)); //this is 8 times faster to memset than "memset(clst, 0, kNrows*kNcls*sizeof(AliTRDcluster*))"
1247 Double_t roady(0.), s2Mean(0.); Int_t ns2Mean(0);
1249 // Do cluster projection and pick up cluster candidates
1250 AliTRDcluster *c(NULL);
1251 AliTRDchamberTimeBin *layer(NULL);
1252 Bool_t kBUFFER = kFALSE;
1253 for (Int_t it = 0; it < kNtb; it++) {
1254 if(!(layer = chamber->GetTB(it))) continue;
1255 if(!Int_t(*layer)) continue;
1256 // get track projection at layers position
1257 dx = fX0 - layer->GetX();
1258 yt = fYref[0] - fYref[1] * dx;
1259 zt = fZref[0] - fZref[1] * dx;
1260 // get standard cluster error corrected for tilt if selected
1261 cp.SetLocalTimeBin(it);
1262 cp.SetSigmaY2(0.02, fDiffT, fExB, dx, -1./*zt*/, fYref[1]);
1263 s2yCl = cp.GetSigmaY2() + sysCov[0]; if(!tilt) s2yCl = (s2yCl + t2*s2zCl)/(1.+t2);
1264 if(TMath::Abs(it-12)<7){ s2Mean += cp.GetSigmaY2(); ns2Mean++;}
1265 // get estimated road in r-phi direction
1266 roady = TMath::Min(3.*TMath::Sqrt(12.*(s2yTrk + s2yCl)), kroady);
1268 AliDebug(5, Form("\n"
1269 " %2d xd[cm]=%6.3f yt[cm]=%7.2f zt[cm]=%8.2f\n"
1270 " syTrk[um]=%6.2f syCl[um]=%6.2f syClTlt[um]=%6.2f\n"
1273 , 1.e4*TMath::Sqrt(s2yTrk), 1.e4*TMath::Sqrt(cp.GetSigmaY2()+sysCov[0]), 1.e4*TMath::Sqrt(s2yCl)
1276 // get clusters from layer
1277 cond[0] = yt/*+0.5*kroadyShift*kroady*/; cond[2] = roady;
1278 cond[1] = zt; cond[3] = kroadz;
1279 Int_t n=0, idx[6]; layer->GetClusters(cond, idx, n, 6);
1280 for(Int_t ic = n; ic--;){
1281 c = (*layer)[idx[ic]];
1282 dx = fX0 - c->GetX();
1283 yt = fYref[0] - fYref[1] * dx;
1284 zt = fZref[0] - fZref[1] * dx;
1285 dz = zt - c->GetZ();
1286 dy = yt - (c->GetY() + (tilt ? (GetTilt() * dz) : 0.));
1287 Int_t r = c->GetPadRow();
1288 clst[r].AddAtAndExpand(c, ncl[r]);
1289 blst[r][ncl[r]] = kTRUE;
1290 idxs[r][ncl[r]] = idx[ic];
1291 zres[r][ncl[r]] = dz/GetPadLength();
1292 yres[r][ncl[r]] = dy;
1293 xres[r][ncl[r]] = dx;
1295 // TODO temporary solution to avoid divercences in error parametrization
1296 s2y[r][ncl[r]] = TMath::Min(c->GetSigmaY2()+sysCov[0], 0.025);
1297 AliDebug(5, Form(" -> dy[cm]=%+7.4f yc[cm]=%7.2f row[%d] idx[%2d]", dy, c->GetY(), r, ncl[r]));
1300 if(ncl[r] >= kNcls) {
1301 AliWarning(Form("Cluster candidates row[%d] reached buffer limit[%d]. Some may be lost.", r, kNcls));
1309 AliDebug(1, Form("CLUSTERS FOUND %d LESS THAN THRESHOLD %d.", ncls, kClmin));
1310 SetErrorMsg(kAttachClFound);
1311 for(Int_t ir(kNrows);ir--;) clst[ir].Clear();
1315 AliDebug(1, Form("CLUSTERS IN TimeBins %d LESS THAN THRESHOLD %d.", ns2Mean, kTBmin));
1316 SetErrorMsg(kAttachClFound);
1317 for(Int_t ir(kNrows);ir--;) clst[ir].Clear();
1320 s2Mean /= ns2Mean; //sMean = TMath::Sqrt(s2Mean);
1321 //Double_t sRef(TMath::Sqrt(s2Mean+s2yTrk)); // reference error parameterization
1323 // organize row candidates
1324 Int_t idxRow[kNrows], nrc(0); Double_t zresRow[kNrows];
1325 for(Int_t ir(0); ir<kNrows; ir++){
1326 idxRow[ir]=-1; zresRow[ir] = 999.;
1327 if(!ncl[ir]) continue;
1328 // get mean z resolution
1329 dz = 0.; for(Int_t ic = ncl[ir]; ic--;) dz += zres[ir][ic]; dz/=ncl[ir];
1331 idxRow[nrc] = ir; zresRow[nrc] = TMath::Abs(dz); nrc++;
1333 AliDebug(4, Form("Found %d clusters in %d rows. Sorting ...", ncls, nrc));
1335 // sort row candidates
1338 if(zresRow[0]>zresRow[1]){ // swap
1339 Int_t itmp=idxRow[1]; idxRow[1] = idxRow[0]; idxRow[0] = itmp;
1340 Double_t dtmp=zresRow[1]; zresRow[1] = zresRow[0]; zresRow[0] = dtmp;
1342 if(TMath::Abs(idxRow[1] - idxRow[0]) != 1){
1343 SetErrorMsg(kAttachRowGap);
1344 AliDebug(2, Form("Rows attached not continuous. Select first candidate.\n"
1345 " row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f",
1346 idxRow[0], ncl[idxRow[0]], zresRow[0], idxRow[1], idxRow[1]<0?0:ncl[idxRow[1]], zresRow[1]));
1347 nrc=1; idxRow[1] = -1; zresRow[1] = 999.;
1351 TMath::Sort(nrc, zresRow, idx0, kFALSE);
1352 nrc = 3; // select only maximum first 3 candidates
1353 Int_t iatmp[] = {-1, -1, -1}; Double_t datmp[] = {999., 999., 999.};
1354 for(Int_t irc(0); irc<nrc; irc++){
1355 iatmp[irc] = idxRow[idx0[irc]];
1356 datmp[irc] = zresRow[idx0[irc]];
1358 idxRow[0] = iatmp[0]; zresRow[0] = datmp[0];
1359 idxRow[1] = iatmp[1]; zresRow[1] = datmp[1];
1360 idxRow[2] = iatmp[2]; zresRow[2] = datmp[2]; // temporary
1361 if(TMath::Abs(idxRow[1] - idxRow[0]) != 1){
1362 SetErrorMsg(kAttachRowGap);
1363 AliDebug(2, Form("Rows attached not continuous. Turn on selection.\n"
1364 "row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f\n"
1365 "row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f\n"
1366 "row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f",
1367 idxRow[0], ncl[idxRow[0]], zresRow[0],
1368 idxRow[1], ncl[idxRow[1]], zresRow[1],
1369 idxRow[2], ncl[idxRow[2]], zresRow[2]));
1370 if(TMath::Abs(idxRow[0] - idxRow[2]) == 1){ // select second candidate
1371 AliDebug(2, "Solved ! Remove second candidate.");
1373 idxRow[1] = idxRow[2]; zresRow[1] = zresRow[2]; // swap
1374 idxRow[2] = -1; zresRow[2] = 999.; // remove
1375 } else if(TMath::Abs(idxRow[1] - idxRow[2]) == 1){
1376 if(ncl[idxRow[1]]+ncl[idxRow[2]] > ncl[idxRow[0]]){
1377 AliDebug(2, "Solved ! Remove first candidate.");
1379 idxRow[0] = idxRow[1]; zresRow[0] = zresRow[1]; // swap
1380 idxRow[1] = idxRow[2]; zresRow[1] = zresRow[2]; // swap
1382 AliDebug(2, "Solved ! Remove second and third candidate.");
1384 idxRow[1] = -1; zresRow[1] = 999.; // remove
1385 idxRow[2] = -1; zresRow[2] = 999.; // remove
1388 AliDebug(2, "Unsolved !!! Remove second and third candidate.");
1390 idxRow[1] = -1; zresRow[1] = 999.; // remove
1391 idxRow[2] = -1; zresRow[2] = 999.; // remove
1393 } else { // remove temporary candidate
1395 idxRow[2] = -1; zresRow[2] = 999.;
1399 AliDebug(4, Form("Sorted row candidates:\n"
1400 " row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f row[%2d] Ncl[%2d] <dz>[cm]=%+8.2f"
1401 , idxRow[0], ncl[idxRow[0]], zresRow[0], idxRow[1], idxRow[1]<0?0:ncl[idxRow[1]], zresRow[1]));
1403 // initialize debug streamer
1404 TTreeSRedirector *pstreamer(NULL);
1405 if((recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 3 && fkReconstructor->IsDebugStreaming())||
1406 AliTRDReconstructor::GetStreamLevel()>30) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1408 // save config. for calibration
1409 TVectorD vdy[2], vdx[2], vs2[2];
1410 for(Int_t jr(0); jr<nrc; jr++){
1411 Int_t ir(idxRow[jr]);
1412 vdx[jr].ResizeTo(ncl[ir]); vdy[jr].ResizeTo(ncl[ir]); vs2[jr].ResizeTo(ncl[ir]);
1413 for(Int_t ic(ncl[ir]); ic--;){
1414 vdx[jr](ic) = xres[ir][ic];
1415 vdy[jr](ic) = yres[ir][ic];
1416 vs2[jr](ic) = s2y[ir][ic];
1419 (*pstreamer) << "AttachClusters4"
1420 << "r0=" << idxRow[0]
1421 << "dz0=" << zresRow[0]
1422 << "dx0=" << &vdx[0]
1423 << "dy0=" << &vdy[0]
1424 << "s20=" << &vs2[0]
1425 << "r1=" << idxRow[1]
1426 << "dz1=" << zresRow[1]
1427 << "dx1=" << &vdx[1]
1428 << "dy1=" << &vdy[1]
1429 << "s21=" << &vs2[1]
1431 vdx[0].Clear(); vdy[0].Clear(); vs2[0].Clear();
1432 vdx[1].Clear(); vdy[1].Clear(); vs2[1].Clear();
1433 if(recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 4 ||AliTRDReconstructor::GetStreamLevel()>4){
1434 Int_t idx(idxRow[1]);
1436 for(Int_t ir(0); ir<kNrows; ir++){
1437 if(clst[ir].GetEntries()>0) continue;
1442 (*pstreamer) << "AttachClusters5"
1443 << "c0.=" << &clst[idxRow[0]]
1444 << "c1.=" << &clst[idx]
1449 //=======================================================================================
1450 // Analyse cluster topology
1451 Double_t f[kNcls], // likelihood factors for segments
1452 r[2][kNcls], // d(dydx) of tracklet candidate with respect to track
1453 xm[2][kNcls], // mean <x>
1454 ym[2][kNcls], // mean <y>
1455 sm[2][kNcls], // mean <s_y>
1456 s[2][kNcls], // sigma_y
1457 p[2][kNcls], // prob of Gauss
1458 q[2][kNcls]; // charge/segment
1459 memset(f, 0, kNcls*sizeof(Double_t));
1460 Int_t index[2][kNcls], n[2][kNcls];
1461 memset(n, 0, 2*kNcls*sizeof(Int_t));
1462 Int_t mts(0), nts[2] = {0, 0}; // no of tracklet segments in row
1463 AliTRDpadPlane *pp(AliTRDtransform::Geometry().GetPadPlane(fDet));
1464 AliTRDtrackletOflHelper helper;
1465 Int_t lyDet(AliTRDgeometry::GetLayer(fDet));
1466 for(Int_t jr(0), n0(0); jr<nrc; jr++){
1467 Int_t ir(idxRow[jr]);
1468 // cluster segmentation
1469 Bool_t kInit(kFALSE);
1471 n0 = helper.Init(pp, &clst[ir]); kInit = kTRUE;
1472 if(!n0 || (helper.ClassifyTopology() == AliTRDtrackletOflHelper::kNormal)){
1473 nts[jr] = 1; memset(index[jr], 0, ncl[ir]*sizeof(Int_t));
1478 nts[jr] = AliTRDtrackletOflHelper::Segmentation(ncl[ir], xres[ir], yres[ir], index[jr]);
1479 for(Int_t ic(ncl[ir]);ic--;) n[jr][index[jr][ic]]++;
1483 // tracklet segment processing
1484 for(Int_t its(0); its<nts[jr]; its++){
1485 if(n[jr][its]<=2) { // don't touch small segments
1486 xm[jr][its] = 0.;ym[jr][its] = 0.;sm[jr][its] = 0.;
1487 for(Int_t ic(ncl[ir]); ic--;){
1488 if(its != index[jr][ic]) continue;
1489 ym[jr][its] += yres[ir][ic];
1490 xm[jr][its] += xres[ir][ic];
1491 sm[jr][its] += TMath::Sqrt(s2y[ir][ic]);
1493 if(n[jr][its]==2){ xm[jr][its] *= 0.5; ym[jr][its] *= 0.5; sm[jr][its] *= 0.5;}
1494 xm[jr][its]= fX0 - xm[jr][its];
1502 // for longer tracklet segments
1503 if(!kInit) n0 = helper.Init(pp, &clst[ir], index[jr], its);
1504 Int_t n1 = helper.GetRMS(r[jr][its], ym[jr][its], s[jr][its], fX0/*xm[jr][its]*/);
1505 p[jr][its] = Double_t(n1)/n0;
1506 sm[jr][its] = helper.GetSyMean();
1507 q[jr][its] = helper.GetQ()/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
1509 Double_t dxm= fX0 - xm[jr][its];
1510 yt = fYref[0] - fYref[1]*dxm;
1511 zt = fZref[0] - fZref[1]*dxm;
1512 // correct tracklet fit for tilt
1513 ym[jr][its]+= GetTilt()*(zt - zc[ir]);
1514 r[jr][its] += GetTilt() * fZref[1];
1515 // correct tracklet fit for track position/inclination
1516 ym[jr][its] = yt - ym[jr][its];
1517 r[jr][its] = (r[jr][its] - fYref[1])/(1+r[jr][its]*fYref[1]);
1518 // report inclination in radians
1519 r[jr][its] = TMath::ATan(r[jr][its]);
1520 if(jr) continue; // calculate only for first row likelihoods
1522 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]);
1525 AliDebug(4, Form(" Tracklet candidates: row[%2d] = %2d row[%2d] = %2d:", idxRow[0], nts[0], idxRow[1], nts[1]));
1526 if(AliLog::GetDebugLevel("TRD", "AliTRDseedV1")>3){
1527 for(Int_t jr(0); jr<nrc; jr++){
1528 Int_t ir(idxRow[jr]);
1529 for(Int_t its(0); its<nts[jr]; its++){
1530 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",
1531 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]);
1536 ( (recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 2 && fkReconstructor->IsDebugStreaming()) ||
1537 AliTRDReconstructor::GetStreamLevel()>2 )
1538 ) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1540 // save config. for calibration
1541 TVectorD vidx, vn, vx, vy, vr, vs, vsm, vp, vf;
1542 vidx.ResizeTo(ncl[idxRow[0]]+(idxRow[1]<0?0:ncl[idxRow[1]]));
1551 for(Int_t jr(0), jts(0), jc(0); jr<nrc; jr++){
1552 Int_t ir(idxRow[jr]);
1553 for(Int_t its(0); its<nts[jr]; its++, jts++){
1554 vn[jts] = n[jr][its];
1555 vx[jts] = xm[jr][its];
1556 vy[jts] = ym[jr][its];
1557 vr[jts] = r[jr][its];
1558 vs[jts] = s[jr][its];
1559 vsm[jts]= sm[jr][its];
1560 vp[jts] = p[jr][its];
1561 vf[jts] = jr?-1.:f[its];
1563 for(Int_t ic(0); ic<ncl[ir]; ic++, jc++) vidx[jc] = index[jr][ic];
1565 (*pstreamer) << "AttachClusters3"
1578 //=========================================================
1579 // Get seed tracklet segment
1580 Int_t idx2[kNcls]; memset(idx2, 0, kNcls*sizeof(Int_t)); // seeding indexing
1581 if(nts[0]>1) TMath::Sort(nts[0], f, idx2);
1582 Int_t is(idx2[0]); // seed index
1583 Int_t idxTrklt[kNcls],
1586 Double_t fTrklt(f[is]),
1594 memset(idxTrklt, 0, kNcls*sizeof(Int_t));
1595 // check seed idx2[0] exit if not found
1597 AliDebug(1, Form("Seed seg[%d] row[%2d] n[%2d] f[%f]<0.01.", is, idxRow[0], n[0][is], f[is]));
1598 SetErrorMsg(kAttachClAttach);
1600 ( (recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 1 && fkReconstructor->IsDebugStreaming()) ||
1601 AliTRDReconstructor::GetStreamLevel()>1 )
1602 ) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1605 if(IsKink()) SETBIT(stat, 1);
1606 if(IsStandAlone()) SETBIT(stat, 2);
1607 if(IsRowCross()) SETBIT(stat, 3);
1608 SETBIT(stat, 4); // set error bit
1609 TVectorD vidx; vidx.ResizeTo(1); vidx[0] = is;
1610 (*pstreamer) << "AttachClusters2"
1616 << "y0=" << fYref[0]
1617 << "z0=" << fZref[0]
1621 << "s2Trk=" << s2yTrk
1622 << "s2Cl=" << s2Mean
1637 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]));
1639 // save seeding segment in the helper
1640 idxTrklt[kts++] = is;
1641 helper.Init(pp, &clst[idxRow[0]], index[0], is);
1642 AliTRDtrackletOflHelper test; // helper to test segment expantion
1643 Float_t rcLikelihood(0.); SetBit(kRowCross, kFALSE);
1644 Double_t dyRez[kNcls]; Int_t idx3[kNcls];
1646 //=========================================================
1647 // Define filter parameters from OCDB
1648 Int_t kNSgmDy[2]; attach->GetNsgmDy(kNSgmDy[0], kNSgmDy[1]);
1649 Float_t kLikeMinRelDecrease[2]; attach->GetLikeMinRelDecrease(kLikeMinRelDecrease[0], kLikeMinRelDecrease[1]);
1650 Float_t kRClikeLimit(attach->GetRClikeLimit());
1652 //=========================================================
1653 // Try attaching next segments from first row (if any)
1655 Int_t jr(0), ir(idxRow[jr]);
1656 // organize secondary sgms. in decreasing order of their distance from seed
1657 memset(dyRez, 0, nts[jr]*sizeof(Double_t));
1658 for(Int_t jts(1); jts<nts[jr]; jts++) {
1659 Int_t its(idx2[jts]);
1660 Double_t rot(TMath::Tan(r[0][is]));
1661 dyRez[its] = TMath::Abs(ym[0][is] - ym[jr][its] + rot*(xm[0][is]-xm[jr][its]));
1663 TMath::Sort(nts[jr], dyRez, idx3, kFALSE);
1664 for (Int_t jts(1); jts<nts[jr]; jts++) {
1665 Int_t its(idx3[jts]);
1666 if(dyRez[its] > kNSgmDy[jr]*smTrklt){
1667 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));
1672 Int_t n0 = test.Expand(&clst[ir], index[jr], its);
1673 Double_t rt, dyt, st, xt, smt, pt, qt, ft;
1674 Int_t n1 = test.GetRMS(rt, dyt, st, fX0/*xt*/);
1675 pt = Double_t(n1)/n0;
1676 smt = test.GetSyMean();
1677 qt = test.GetQ()/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
1680 Double_t dxm= fX0 - xt;
1681 yt = fYref[0] - fYref[1]*dxm;
1682 zt = fZref[0] - fZref[1]*dxm;
1683 // correct tracklet fit for tilt
1684 dyt+= GetTilt()*(zt - zc[idxRow[0]]);
1685 rt += GetTilt() * fZref[1];
1686 // correct tracklet fit for track position/inclination
1688 rt = (rt - fYref[1])/(1+rt*fYref[1]);
1689 // report inclination in radians
1690 rt = TMath::ATan(rt);
1692 ft = (n0>=2) ? attach->CookLikelihood(chgPos, lyDet, fPt, phiTrk, n0, dyt/*sRef*/, rt*TMath::RadToDeg(), st/smt) : 0.;
1693 Bool_t kAccept(ft>=fTrklt*(1.-kLikeMinRelDecrease[jr]));
1695 AliDebug(2, Form("%s seg[%d] row[%2d] n[%2d] dy[%f] r[%+5.2f] s[%+5.2f] f[%f] < %4.2f*F[%f].",
1696 (kAccept?"Adding":"Reject"), its, idxRow[jr], n0, dyt, rt*TMath::RadToDeg(), st/smt, ft, 1.-kLikeMinRelDecrease[jr], fTrklt*(1.-kLikeMinRelDecrease[jr])));
1698 idxTrklt[kts++] = its;
1708 helper.Expand(&clst[ir], index[jr], its);
1713 //=========================================================
1714 // Try attaching next segments from second row (if any)
1715 if(nts[1] && (rcLikelihood = zresRow[0]/zresRow[1]) > kRClikeLimit){
1716 // organize secondaries in decreasing order of their distance from seed
1717 Int_t jr(1), ir(idxRow[jr]);
1718 memset(dyRez, 0, nts[jr]*sizeof(Double_t));
1719 Double_t rot(TMath::Tan(r[0][is]));
1720 for(Int_t jts(0); jts<nts[jr]; jts++) {
1721 dyRez[jts] = TMath::Abs(ym[0][is] - ym[jr][jts] + rot*(xm[0][is]-xm[jr][jts]));
1723 TMath::Sort(nts[jr], dyRez, idx3, kFALSE);
1724 for (Int_t jts(0); jts<nts[jr]; jts++) {
1725 Int_t its(idx3[jts]);
1726 if(dyRez[its] > kNSgmDy[jr]*smTrklt){
1727 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));
1732 Int_t n0 = test.Expand(&clst[ir], index[jr], its);
1733 Double_t rt, dyt, st, xt, smt, pt, qt, ft;
1734 Int_t n1 = test.GetRMS(rt, dyt, st, fX0/*xt*/);
1735 pt = Double_t(n1)/n0;
1736 smt = test.GetSyMean();
1737 qt = test.GetQ()/TMath::Sqrt(1. + fYref[1]*fYref[1] + fZref[1]*fZref[1]);
1740 Double_t dxm= fX0 - xt;
1741 yt = fYref[0] - fYref[1]*dxm;
1742 zt = fZref[0] - fZref[1]*dxm;
1743 // correct tracklet fit for tilt
1744 dyt+= GetTilt()*(zt - zc[idxRow[0]]);
1745 rt += GetTilt() * fZref[1];
1746 // correct tracklet fit for track position/inclination
1748 rt = (rt - fYref[1])/(1+rt*fYref[1]);
1749 // report inclination in radians
1750 rt = TMath::ATan(rt);
1752 ft = (n0>=2) ? attach->CookLikelihood(chgPos, lyDet, fPt, phiTrk, n0, dyt/*sRef*/, rt*TMath::RadToDeg(), st/smt) : 0.;
1753 Bool_t kAccept(ft>=fTrklt*(1.-kLikeMinRelDecrease[jr]));
1755 AliDebug(2, Form("%s seg[%d] row[%2d] n[%2d] dy[%f] r[%+5.2f] s[%+5.2f] f[%f] < %4.2f*F[%f].",
1756 (kAccept?"Adding":"Reject"), its, idxRow[jr], n0, dyt, rt*TMath::RadToDeg(), st/smt, ft, 1.-kLikeMinRelDecrease[jr], fTrklt*(1.-kLikeMinRelDecrease[jr])));
1758 idxTrklt[kts++] = its;
1768 helper.Expand(&clst[ir], index[jr], its);
1769 SetBit(kRowCross, kTRUE); // mark pad row crossing
1773 // clear local copy of clusters
1774 for(Int_t ir(0); ir<kNrows; ir++) clst[ir].Clear();
1777 ((recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 1 && fkReconstructor->IsDebugStreaming()) ||
1778 AliTRDReconstructor::GetStreamLevel()>1 )
1779 ) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
1782 if(IsKink()) SETBIT(stat, 1);
1783 if(IsStandAlone()) SETBIT(stat, 2);
1784 if(IsRowCross()) SETBIT(stat, 3);
1785 TVectorD vidx; vidx.ResizeTo(kts);
1786 for(Int_t its(0); its<kts; its++) vidx[its] = idxTrklt[its];
1787 (*pstreamer) << "AttachClusters2"
1793 << "y0=" << fYref[0]
1794 << "z0=" << fZref[0]
1798 << "s2Trk=" << s2yTrk
1799 << "s2Cl=" << s2Mean
1814 //=========================================================
1816 Int_t nselected(0), nc(0);
1817 TObjArray *selected(helper.GetClusters());
1818 if(!selected || !(nselected = selected->GetEntriesFast())){
1819 AliError("Cluster candidates missing !!!");
1820 SetErrorMsg(kAttachClAttach);
1823 for(Int_t ic(0); ic<nselected; ic++){
1824 if(!(c = (AliTRDcluster*)selected->At(ic))) continue;
1825 Int_t it(c->GetPadTime()),
1826 jr(Int_t(helper.GetRow() != c->GetPadRow())),
1829 AliDebug(1, Form("Multiple clusters/tb for D[%03d] Tb[%02d] Row[%2d]", fDet, it, c->GetPadRow()));
1830 continue; // already booked
1832 // TODO proper indexing of clusters !!
1833 fIndexes[idx] = chamber->GetTB(it)->GetGlobalIndex(idxs[idxRow[jr]][ic]);
1837 AliDebug(2, Form("Clusters Found[%2d] Attached[%2d] RC[%c]", nselected, nc, IsRowCross()?'y':'n'));
1839 // number of minimum numbers of clusters expected for the tracklet
1841 AliDebug(1, Form("NOT ENOUGH CLUSTERS %d ATTACHED TO THE TRACKLET [min %d] FROM FOUND %d.", nc, kClmin, ncls));
1842 SetErrorMsg(kAttachClAttach);
1847 // Load calibration parameters for this tracklet
1850 // calculate dx for time bins in the drift region (calibration aware)
1851 Float_t x[2] = {0.,0.}; Int_t tb[2]={0,0};
1852 for (Int_t it = t0, irp=0; irp<2 && it < AliTRDtrackerV1::GetNTimeBins(); it++) {
1853 if(!fClusters[it]) continue;
1854 x[irp] = fClusters[it]->GetX();
1855 tb[irp] = fClusters[it]->GetLocalTimeBin();
1858 Int_t dtb = tb[1] - tb[0];
1859 fdX = dtb ? (x[0] - x[1]) / dtb : 0.15;
1863 //____________________________________________________________
1864 void AliTRDseedV1::Bootstrap(const AliTRDReconstructor *rec)
1866 // Fill in all derived information. It has to be called after recovery from file or HLT.
1867 // The primitive data are
1868 // - list of clusters
1869 // - detector (as the detector will be removed from clusters)
1870 // - position of anode wire (fX0) - temporary
1871 // - track reference position and direction
1872 // - momentum of the track
1873 // - time bin length [cm]
1875 // A.Bercuci <A.Bercuci@gsi.de> Oct 30th 2008
1877 fkReconstructor = rec;
1879 SetPadPlane(g.GetPadPlane(fDet));
1881 //fSnp = fYref[1]/TMath::Sqrt(1+fYref[1]*fYref[1]);
1883 Int_t n = 0, nshare = 0, nused = 0;
1884 AliTRDcluster **cit = &fClusters[0];
1885 for(Int_t ic = kNclusters; ic--; cit++){
1888 if((*cit)->IsShared()) nshare++;
1889 if((*cit)->IsUsed()) nused++;
1891 SetN(n); SetNUsed(nused); SetNShared(nshare);
1898 //____________________________________________________________________
1899 Bool_t AliTRDseedV1::Fit(UChar_t opt)
1902 // Linear fit of the clusters attached to the tracklet
1905 // - opt : switch for tilt pad correction of cluster y position. Options are
1906 // 0 no correction [default]
1907 // 1 full tilt correction [dz/dx and z0]
1908 // 2 pseudo tilt correction [dz/dx from pad-chamber geometry]
1911 // True if successful
1913 // Detailed description
1915 // Fit in the xy plane
1917 // The fit is performed to estimate the y position of the tracklet and the track
1918 // angle in the bending plane. The clusters are represented in the chamber coordinate
1919 // system (with respect to the anode wire - see AliTRDtrackerV1::FollowBackProlongation()
1920 // on how this is set). The x and y position of the cluster and also their variances
1921 // are known from clusterizer level (see AliTRDcluster::GetXloc(), AliTRDcluster::GetYloc(),
1922 // AliTRDcluster::GetSX() and AliTRDcluster::GetSY()).
1923 // If gaussian approximation is used to calculate y coordinate of the cluster the position
1924 // is recalculated taking into account the track angle. The general formula to calculate the
1925 // error of cluster position in the gaussian approximation taking into account diffusion and track
1926 // inclination is given for TRD by:
1928 // #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}
1931 // Since errors are calculated only in the y directions, radial errors (x direction) are mapped to y
1932 // by projection i.e.
1934 // #sigma_{x|y} = tg(#phi) #sigma_{x}
1936 // and also by the lorentz angle correction
1938 // Fit in the xz plane
1940 // The "fit" is performed to estimate the radial position (x direction) where pad row cross happens.
1941 // If no pad row crossing the z position is taken from geometry and radial position is taken from the xy
1944 // There are two methods to estimate the radial position of the pad row cross:
1945 // 1. leading cluster radial position : Here the lower part of the tracklet is considered and the last
1946 // cluster registered (at radial x0) on this segment is chosen to mark the pad row crossing. The error
1947 // of the z estimate is given by :
1949 // #sigma_{z} = tg(#theta) #Delta x_{x_{0}}/12
1951 // The systematic errors for this estimation are generated by the following sources:
1952 // - no charge sharing between pad rows is considered (sharp cross)
1953 // - missing cluster at row cross (noise peak-up, under-threshold signal etc.).
1955 // 2. charge fit over the crossing point : Here the full energy deposit along the tracklet is considered
1956 // to estimate the position of the crossing by a fit in the qx plane. The errors in the q directions are
1957 // parameterized as s_q = q^2. The systematic errors for this estimation are generated by the following sources:
1958 // - no general model for the qx dependence
1959 // - physical fluctuations of the charge deposit
1960 // - gain calibration dependence
1962 // Estimation of the radial position of the tracklet
1964 // For pad row cross the radial position is taken from the xz fit (see above). Otherwise it is taken as the
1965 // interpolation point of the tracklet i.e. the point where the error in y of the fit is minimum. The error
1966 // in the y direction of the tracklet is (see AliTRDseedV1::GetCovAt()):
1968 // #sigma_{y} = #sigma^{2}_{y_{0}} + 2xcov(y_{0}, dy/dx) + #sigma^{2}_{dy/dx}
1970 // and thus the radial position is:
1972 // x = - cov(y_{0}, dy/dx)/#sigma^{2}_{dy/dx}
1975 // Estimation of tracklet position error
1977 // The error in y direction is the error of the linear fit at the radial position of the tracklet while in the z
1978 // direction is given by the cluster error or pad row cross error. In case of no pad row cross this is given by:
1980 // #sigma_{y} = #sigma^{2}_{y_{0}} - 2cov^{2}(y_{0}, dy/dx)/#sigma^{2}_{dy/dx} + #sigma^{2}_{dy/dx}
1981 // #sigma_{z} = Pad_{length}/12
1983 // For pad row cross the full error is calculated at the radial position of the crossing (see above) and the error
1984 // in z by the width of the crossing region - being a matter of parameterization.
1986 // #sigma_{z} = tg(#theta) #Delta x_{x_{0}}/12
1988 // In case of no tilt correction (default in the barrel tracking) the tilt is taken into account by the rotation of
1989 // the covariance matrix. See AliTRDseedV1::GetCovAt() for details.
1992 // A.Bercuci <A.Bercuci@gsi.de>
1994 if(!fkReconstructor){
1995 AliError("The tracklet needs the reconstruction setup. Please initialize by SetReconstructor().");
1998 if(!IsCalibrated()) Calibrate();
2000 AliWarning(Form("Option [%d] outside range [0, 2]. Using default",opt));
2004 const Int_t kClmin = 8;
2005 const Float_t kScalePulls = 10.; // factor to scale y pulls - NOT UNDERSTOOD
2006 // get track direction
2007 Double_t y0 = fYref[0];
2008 Double_t dydx = fYref[1];
2009 Double_t z0 = fZref[0];
2010 Double_t dzdx = fZref[1];
2012 AliTRDtrackerV1::AliTRDLeastSquare fitterY;
2013 AliTRDtrackerV1::AliTRDLeastSquare fitterZ;
2015 // book cluster information
2016 Double_t qc[kNclusters], xc[kNclusters], yc[kNclusters], zc[kNclusters], sy[kNclusters];
2018 Bool_t tilt(opt==1) // full tilt correction
2019 ,pseudo(opt==2) // pseudo tilt correction
2020 ,rc(IsRowCross()) // row cross candidate
2021 ,kDZDX(IsPrimary());// switch dzdx calculation for barrel primary tracks
2022 Int_t n(0); // clusters used in fit
2023 AliTRDcluster *c(NULL), *cc(NULL), **jc = &fClusters[0];
2024 const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam(); //the dynamic cast in GetRecoParam is slow, so caching the pointer to it
2026 const Char_t *tcName[]={"NONE", "FULL", "HALF"};
2027 AliDebug(2, Form("Options : TC[%s] dzdx[%c]", tcName[opt], kDZDX?'Y':'N'));
2030 for (Int_t ic=0; ic<kNclusters; ic++, ++jc) {
2031 xc[ic] = -1.; yc[ic] = 999.; zc[ic] = 999.; sy[ic] = 0.;
2032 if(!(c = (*jc))) continue;
2033 if(!c->IsInChamber()) continue;
2034 // compute pseudo tilt correction
2036 fZfit[0] = c->GetZ();
2038 for(Int_t kc=AliTRDseedV1::kNtb; kc<AliTRDseedV1::kNclusters; kc++){
2039 if(!(cc=fClusters[kc])) continue;
2040 if(!cc->IsInChamber()) continue;
2041 fZfit[0] += cc->GetZ(); fZfit[0] *= 0.5;
2045 fZfit[1] = fZfit[0]/fX0;
2047 fZfit[0] += fZfit[1]*0.5*AliTRDgeometry::CdrHght();
2048 fZfit[1] = fZfit[0]/fX0;
2053 // TODO use this information to adjust cluster error parameterization
2055 // if(c->GetNPads()>4) w = .5;
2056 // if(c->GetNPads()>5) w = .2;
2059 qc[n] = TMath::Abs(c->GetQ());
2060 // pad row of leading
2062 xc[n] = fX0 - c->GetX();
2064 // Recalculate cluster error based on tracking information
2065 c->SetSigmaY2(fS2PRF, fDiffT, fExB, xc[n], -1./*zcorr?zt:-1.*/, dydx);
2066 c->SetSigmaZ2(fPad[0]*fPad[0]/12.); // for HLT
2067 sy[n] = TMath::Sqrt(c->GetSigmaY2());
2069 yc[n] = recoParam->UseGAUS() ?
2070 c->GetYloc(y0, sy[n], GetPadWidth()): c->GetY();
2073 //optional r-phi correction
2074 //printf(" n[%2d] yc[%7.5f] ", n, yc[n]);
2075 Float_t correction(0.);
2076 if(tilt) correction = fPad[2]*(xc[n]*dzdx + zc[n] - z0);
2077 else if(pseudo) correction = fPad[2]*(xc[n]*fZfit[1] + zc[n]-fZfit[0]);
2079 //printf("corr(%s%s)[%7.5f] yc1[%7.5f]\n", (tilt?"TC":""), (zcorr?"PC":""), correction, yc[n]);
2081 AliDebug(5, Form(" tb[%2d] dx[%6.3f] y[%6.2f+-%6.3f]", c->GetLocalTimeBin(), xc[n], yc[n], sy[n]));
2082 fitterY.AddPoint(&xc[n], yc[n], sy[n]);
2083 if(rc) fitterZ.AddPoint(&xc[n], qc[n]*(ic<kNtb?1.:-1.), 1.);
2089 AliDebug(1, Form("Not enough clusters to fit. Clusters: Attached[%d] Fit[%d].", GetN(), n));
2090 SetErrorMsg(kFitCl);
2094 if(!fitterY.Eval()){
2095 AliDebug(1, "Fit Y failed.");
2096 SetErrorMsg(kFitFailedY);
2099 fYfit[0] = fitterY.GetFunctionParameter(0);
2100 fYfit[1] = -fitterY.GetFunctionParameter(1);
2103 fitterY.GetCovarianceMatrix(p);
2104 fCov[0] = kScalePulls*p[1]; // variance of y0
2105 fCov[1] = kScalePulls*p[2]; // covariance of y0, dydx
2106 fCov[2] = kScalePulls*p[0]; // variance of dydx
2107 // the ref radial position is set at the minimum of
2108 // the y variance of the tracklet
2109 fX = -fCov[1]/fCov[2];
2110 fS2Y = fCov[0] +2.*fX*fCov[1] + fX*fX*fCov[2];
2112 Float_t xs=fX+.5*AliTRDgeometry::CamHght();
2113 if(xs < 0. || xs > AliTRDgeometry::CamHght()+AliTRDgeometry::CdrHght()){
2114 AliDebug(1, Form("Ref radial position ouside chamber x[%5.2f].", fX));
2115 SetErrorMsg(kFitFailedY);
2119 /* // THE LEADING CLUSTER METHOD for z fit
2121 Int_t ic=n=kNclusters-1; jc = &fClusters[ic];
2122 AliTRDcluster *c0 =0x0, **kc = &fClusters[kNtb-1];
2123 for(; ic>kNtb; ic--, --jc, --kc){
2124 if((c0 = (*kc)) && c0->IsInChamber() && (xMin>c0->GetX())) xMin = c0->GetX();
2125 if(!(c = (*jc))) continue;
2126 if(!c->IsInChamber()) continue;
2127 zc[kNclusters-1] = c->GetZ();
2128 fX = fX0 - c->GetX();
2130 fZfit[0] = .5*(zc[0]+zc[kNclusters-1]); fZfit[1] = 0.;
2131 // Error parameterization
2132 fS2Z = fdX*fZref[1];
2133 fS2Z *= fS2Z; fS2Z *= 0.2887; // 1/sqrt(12)*/
2136 if(opt!=1 && IsRowCross()){
2137 if(!fitterZ.Eval()) SetErrorMsg(kFitFailedZ);
2138 if(!HasError(kFitFailedZ) && TMath::Abs(fitterZ.GetFunctionParameter(1))>1.e-10){
2139 // TODO - one has to recalculate xy fit based on
2140 // better knowledge of z position
2141 // Double_t x = -fitterZ.GetFunctionParameter(0)/fitterZ.GetFunctionParameter(1);
2142 // Double_t z0 = .5*(zc[0]+zc[n-1]);
2143 // fZfit[0] = z0 + fZfit[1]*x;
2144 // fZfit[1] = fZfit[0]/fX0;
2145 // redo fit on xy plane
2147 // temporary external error parameterization
2148 fS2Z = 0.05+0.4*TMath::Abs(fZref[1]); fS2Z *= fS2Z;
2149 // TODO correct formula
2150 //fS2Z = sigma_x*TMath::Abs(fZref[1]);
2152 //fZfit[0] = zc[0] + dzdx*0.5*AliTRDgeometry::CdrHght();
2153 fS2Z = GetPadLength()*GetPadLength()/12.;
2159 //____________________________________________________________________
2160 Bool_t AliTRDseedV1::FitRobust(AliTRDpadPlane *pp, TGeoHMatrix *mDet, Float_t bz, Int_t chg, Int_t opt)
2163 // Linear fit of the clusters attached to the tracklet
2164 // The fit is performed in local chamber coordinates (27.11.2013) to take into account correctly the misalignment
2165 // Also the pad row cross is checked here and some background is removed
2168 // A.Bercuci <A.Bercuci@gsi.de>
2170 TTreeSRedirector *pstreamer(NULL);
2171 const AliTRDrecoParam* const recoParam = fkReconstructor->GetRecoParam();
2172 if( (recoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 3 && fkReconstructor->IsDebugStreaming()) ||
2173 AliTRDReconstructor::GetStreamLevel()>3 ) pstreamer = fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
2175 // factor to scale y pulls.
2176 // ideally if error parametrization correct this is 1.
2177 //Float_t lyScaler = 1./(AliTRDgeometry::GetLayer(fDet)+1.);
2178 Float_t kScalePulls = 1.;
2179 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
2181 AliWarning("No access to calibration data");
2183 // Retrieve the CDB container class with the parametric likelihood
2184 const AliTRDCalTrkAttach *attach = calibration->GetAttachObject();
2186 AliWarning("No usable AttachClusters calib object.");
2188 //kScalePulls = attach->GetScaleCov();//*lyScaler;
2190 // Retrieve chamber status
2191 SetChmbGood(calibration->IsChamberGood(fDet));
2192 if(!IsChmbGood()) kScalePulls*=10.;
2194 AliTRDCommonParam *cp = AliTRDCommonParam::Instance();
2195 Double_t freq(cp?cp->GetSamplingFrequency():10.);
2197 // evaluate locally z and dzdx from TRD only information
2198 if(EstimatedCrossPoint(pp, bz)<0.) return kFALSE;
2200 //printf("D%03d RC[%c] dzdx[%f %f] opt[%d]\n", fDet, IsRowCross()?'y':'n', fZref[1], fZfit[1], opt);
2201 Double_t //xchmb = 0.5 * AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick(),
2202 //zchmb = 0.5 * (pp->GetRow0() + pp->GetRowEnd()),
2203 z0(0.5 * (pp->GetRow0() + pp->GetRowEnd()) + fZfit[0]),
2204 DZ(pp->GetRow0() - pp->GetRowEnd() - pp->GetAnodeWireOffset() + fZfit[0]),
2206 Double_t xc[kNclusters], yc[kNclusters], dz(0.), dzdx(0.),
2207 s2dz(0.), s2dzdx(0.), sy[kNclusters],
2208 s2x((8.33e-2/freq/freq+1.56e-2)*fVD*fVD), // error of 1tb + error of mean time (TRF)
2209 t2(fPad[2]*fPad[2]), loc[3]={0.};
2210 Int_t n(0), // clusters used in fit
2211 row[]={-1, -1};// pad row spanned by the tracklet
2212 Double_t ycorr(UnbiasY(IsRowCross(), bz)),
2213 kS2Ycorr(recoParam->GetS2Ycorr(IsRowCross(), chg>0));
2215 AliTRDcluster *c(NULL), **jc = &fClusters[0];
2216 for(Int_t ic=0; ic<kNtb; ic++, ++jc) {
2217 if(!(c = (*jc))) continue;
2218 if(!c->IsInChamber()) continue;
2220 row[0] = c->GetPadRow();
2221 z = pp->GetRowPos(row[0]) - 0.5*pp->GetRowSize(row[0]);
2223 case 0: // no dz correction (only for RC tracklet) and dzdx from chamber position assuming primary
2224 dzdx = IsRowCross()?fZfit[1]:0.;
2225 s2dzdx= IsRowCross()?GetS2DZDX(dzdx):0.;
2226 dz = IsRowCross()?(z - z0):0.;
2227 s2dz = IsRowCross()?fS2Z:0.;
2229 case 1: // dz correction only for RC tracklet and dzdx from reference
2231 dz = IsRowCross()?(z - z0):0.;
2233 case 2: // full z correction (z0 & dzdx from reference)
2235 dz = c->GetZ()-fZref[0];
2238 AliError(Form("Wrong option fit %d !", opt));
2242 //Use local cluster coordinates
2243 //A.Bercuci 27.11.13/30.06.14
2244 Double_t trk[] = {c->GetX(), c->GetY(), c->GetZ()};
2245 mDet->MasterToLocal(trk, loc);
2246 xc[n] = AliTRDgeometry::AnodePos()-loc[0]; //c->GetXloc(fT0, fVD); // c->GetX();
2247 yc[n] = loc[1]; //c->GetYloc(pp->GetColPos(col) + .5*cs, fS2PRF, cs) - xc[n]*fExB; //c->GetY();
2248 yc[n]-= fPad[2]*(dz+xc[n]*dzdx);
2250 if(IsRowCross()){ // estimate closest distance to anode wire
2252 d -= ((Int_t)(2 * d)) / 2.0;
2253 if (d > 0.25) d = 0.5 - d;
2255 // recalculate cluster error from knowledge of the track inclination in the bending plane
2256 // and eventually distance to anode wire
2257 c->SetSigmaY2(fS2PRF, fDiffT, fExB, xc[n], d, fYref[1]);
2258 s2x = c->GetSX(c->GetLocalTimeBin(), d); s2x*=s2x;
2259 sy[n] = c->GetSigmaY2()>0?(TMath::Min(Double_t(c->GetSigmaY2()), 6.4e-3)):6.4e-3;
2260 sy[n]+= t2*(s2dz+xc[n]*xc[n]*s2dzdx+dzdx*dzdx*s2x);
2261 sy[n] = TMath::Sqrt(sy[n]);
2264 for(Int_t ic=kNtb; ic<kNclusters; ic++, ++jc) {
2265 if(!(c = (*jc))) continue;
2266 if(!c->IsInChamber()) continue;
2268 row[1] = c->GetPadRow();
2269 z = pp->GetRowPos(row[1]) - 0.5*pp->GetRowSize(row[1]);
2271 case 0: // no dz correction (only for RC tracklet) and dzdx from chamber position assuming primary
2275 case 1: // dz correction only for RC tracklet and dzdx from reference
2279 case 2: // full z correction (z0 & dzdx from reference)
2281 dz = c->GetZ()-fZref[0];
2284 AliError(Form("Wrong option fit %d !", opt));
2289 //Use local cluster coordinates - the code should be identical with AliTRDtransform::Transform() !!!
2290 //A.Bercuci 27.11.13
2291 Double_t trk[] = {c->GetX(), c->GetY(), c->GetZ()};
2292 mDet->MasterToLocal(trk, loc);
2293 xc[n] = AliTRDgeometry::AnodePos()-loc[0]; //c->GetXloc(fT0, fVD); // c->GetX();
2294 yc[n] = loc[1]; //c->GetYloc(pp->GetColPos(col) + .5*cs, fS2PRF, cs) - xc[n]*fExB; //c->GetY();
2295 yc[n]-= fPad[2]*(dz+xc[n]*dzdx);
2299 d -= ((Int_t)(2 * d)) / 2.0;
2300 if (d > 0.25) d = 0.5 - d;
2301 c->SetSigmaY2(fS2PRF, fDiffT, fExB, xc[n], d, fYref[1]);
2302 s2x = c->GetSX(c->GetLocalTimeBin(), d); s2x*=s2x;
2303 sy[n] = c->GetSigmaY2()>0?(TMath::Min(Double_t(c->GetSigmaY2()), 6.4e-3)):6.4e-3;
2304 sy[n]+= t2*(s2dz+xc[n]*xc[n]*s2dzdx+dzdx*dzdx*s2x);
2305 sy[n] = TMath::Sqrt(sy[n]);
2310 // the ref radial position is set close to the minimum of
2311 // the y variance of the tracklet
2312 fX = 0.;//set reference to anode wire
2313 Double_t par[3] = {0.,0.,fX}, cov[3];
2314 if(!AliTRDtrackletOflHelper::Fit(n, xc, yc, sy, par, 1.5, cov)){
2315 AliDebug(1, Form("Tracklet fit failed D[%03d].", fDet));
2316 SetErrorMsg(kFitCl);
2319 fYfit[0] = par[0] - fX * par[1];
2321 //printf(" yfit: %f [%f] x[%e] dydx[%f]\n", fYfit[0], par[0], fX, par[1]);
2323 fCov[0] = kS2Ycorr*cov[0]; // variance of y0
2324 fCov[1] = kScalePulls*cov[2]; // covariance of y0, dydx
2325 fCov[2] = kScalePulls*cov[1]; // variance of dydx
2326 // check radial position
2327 Float_t xs=fX+.5*AliTRDgeometry::CamHght();
2328 if(xs < 0. || xs > AliTRDgeometry::CamHght()+AliTRDgeometry::CdrHght()){
2329 AliDebug(1, Form("Ref radial position x[%5.2f] ouside D[%3d].", fX, fDet));
2330 SetErrorMsg(kFitFailedY);
2334 Double_t padEffLength(fPad[0] - TMath::Abs(dzdx));
2335 fS2Z = padEffLength*padEffLength/12.;
2337 AliDebug(2, Form("[I] x[cm]=%6.2f y[cm]=%+5.2f z[cm]=%+6.2f dydx[deg]=%+5.2f", GetX(), GetY(), GetZ(), TMath::ATan(fYfit[1])*TMath::RadToDeg()));
2342 yt = fYref[0]-fX*fYref[1];
2344 TVectorD vcov(3); vcov[0]=cov[0];vcov[1]=cov[1];vcov[2]=cov[2];
2345 Double_t sm(0.), chi2(0.), tmp, dy[kNclusters];
2346 for(Int_t ic(0); ic<n; ic++){
2348 dy[ic] = yc[ic]-(fYfit[0]+(xc[ic]-fX0)*fYfit[1]); tmp = dy[ic]/sy[ic];
2351 sm /= n; chi2 = TMath::Sqrt(chi2);
2352 Double_t m(0.), s(0.);
2353 AliMathBase::EvaluateUni(n, dy, m, s, 0);
2354 (*pstreamer) << "FitRobust4"
2355 << "stat=" << status
2360 << "y0=" << fYfit[0]
2363 << "dydx=" << fYfit[1]
2366 << "dydxt="<< fYref[1]
2376 //___________________________________________________________________
2377 void AliTRDseedV1::SetXYZ(TGeoHMatrix *mDet)
2379 // Apply alignment to the local position of tracklet
2380 // A.Bercuci @ 27.11.2013
2382 Double_t loc[] = {AliTRDgeometry::AnodePos(), GetLocalY(), fZfit[0]}, trk[3]={0.};
2383 mDet->LocalToMaster(loc, trk);
2388 // if(!IsRowCross()){/*fZfit[1] *= 1.09;*/ return;}
2389 // // recalculate local z coordinate assuming primary track for row cross tracklets
2390 // Double_t zoff(fZ-fZfit[0]); // no alignment aware !
2391 // //printf("SetXYZ : zoff[%f] zpp[%f]\n", zoff, zpp);
2392 // fZfit[0] = fX0*fZfit[1] - zoff;
2393 // // recalculate tracking coordinates based on the new z coordinate
2394 // loc[2] = fZfit[0];
2395 // mDet->LocalToMaster(loc, trk);
2398 // fZ = trk[2];//-zcorr[stk];
2399 //fZfit[1] = /*(IsRowCross()?1.05:1.09)**/fZ/(fX0-fS2Y);
2403 //___________________________________________________________________
2404 void AliTRDseedV1::Print(Option_t *o) const
2407 // Printing the seedstatus
2410 AliInfo(Form("Det[%3d] X0[%7.2f] Pad{L[%5.2f] W[%5.2f] Tilt[%+6.2f]}", fDet, fX0, GetPadLength(), GetPadWidth(), GetTilt()));
2411 AliInfo(Form("N[%2d] Nused[%2d] Nshared[%2d] [%d]", GetN(), GetNUsed(), GetNShared(), fN));
2412 AliInfo(Form("FLAGS : RC[%c] Kink[%c] SA[%c]", IsRowCross()?'y':'n', IsKink()?'y':'n', IsStandAlone()?'y':'n'));
2413 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));
2415 Double_t cov[3], x=GetX();
2417 AliInfo(" | x[cm] | y[cm] | z[cm] | dydx | dzdx |");
2418 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]));
2419 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]));
2420 AliInfo(Form("P / Pt [GeV/c] = %f / %f", GetMomentum(), fPt));
2421 if(IsStandAlone()) AliInfo(Form("C Rieman / Vertex [1/cm] = %f / %f", fC[0], fC[1]));
2422 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]));
2423 AliInfo(Form("PID = %5.3f / %5.3f / %5.3f / %5.3f / %5.3f", fProb[0], fProb[1], fProb[2], fProb[3], fProb[4]));
2425 if(strcmp(o, "a")!=0) return;
2427 AliTRDcluster* const* jc = &fClusters[0];
2428 for(int ic=0; ic<kNclusters; ic++, jc++) {
2429 if(!(*jc)) continue;
2435 //___________________________________________________________________
2436 Bool_t AliTRDseedV1::IsEqual(const TObject *o) const
2438 // Checks if current instance of the class has the same essential members
2441 if(!o) return kFALSE;
2442 const AliTRDseedV1 *inTracklet = dynamic_cast<const AliTRDseedV1*>(o);
2443 if(!inTracklet) return kFALSE;
2445 for (Int_t i = 0; i < 2; i++){
2446 if ( fYref[i] != inTracklet->fYref[i] ) return kFALSE;
2447 if ( fZref[i] != inTracklet->fZref[i] ) return kFALSE;
2450 if ( TMath::Abs(fS2Y - inTracklet->fS2Y)>1.e-10 ) return kFALSE;
2451 if ( TMath::Abs(GetTilt() - inTracklet->GetTilt())>1.e-10 ) return kFALSE;
2452 if ( TMath::Abs(GetPadLength() - inTracklet->GetPadLength())>1.e-10 ) return kFALSE;
2454 for (Int_t i = 0; i < kNclusters; i++){
2455 // if ( fX[i] != inTracklet->GetX(i) ) return kFALSE;
2456 // if ( fY[i] != inTracklet->GetY(i) ) return kFALSE;
2457 // if ( fZ[i] != inTracklet->GetZ(i) ) return kFALSE;
2458 if ( fIndexes[i] != inTracklet->fIndexes[i] ) return kFALSE;
2460 // if ( fUsable != inTracklet->fUsable ) return kFALSE;
2462 for (Int_t i=0; i < 2; i++){
2463 if ( fYfit[i] != inTracklet->fYfit[i] ) return kFALSE;
2464 if ( fZfit[i] != inTracklet->fZfit[i] ) return kFALSE;
2465 if ( fLabels[i] != inTracklet->fLabels[i] ) return kFALSE;
2468 /* if ( fMeanz != inTracklet->GetMeanz() ) return kFALSE;
2469 if ( fZProb != inTracklet->GetZProb() ) return kFALSE;*/
2470 if ( fN != inTracklet->fN ) return kFALSE;
2471 //if ( fNUsed != inTracklet->fNUsed ) return kFALSE;
2472 //if ( fFreq != inTracklet->GetFreq() ) return kFALSE;
2473 //if ( fNChange != inTracklet->GetNChange() ) return kFALSE;
2475 if ( TMath::Abs(fC[0] - inTracklet->fC[0])>1.e-10 ) return kFALSE;
2476 //if ( fCC != inTracklet->GetCC() ) return kFALSE;
2477 if ( TMath::Abs(fChi2 - inTracklet->fChi2)>1.e-10 ) return kFALSE;
2478 // if ( fChi2Z != inTracklet->GetChi2Z() ) return kFALSE;
2480 if ( fDet != inTracklet->fDet ) return kFALSE;
2481 if ( TMath::Abs(fPt - inTracklet->fPt)>1.e-10 ) return kFALSE;
2482 if ( TMath::Abs(fdX - inTracklet->fdX)>1.e-10 ) return kFALSE;
2484 for (Int_t iCluster = 0; iCluster < kNclusters; iCluster++){
2485 AliTRDcluster *curCluster = fClusters[iCluster];
2486 AliTRDcluster *inCluster = inTracklet->fClusters[iCluster];
2487 if (curCluster && inCluster){
2488 if (! curCluster->IsEqual(inCluster) ) {
2489 curCluster->Print();
2494 // if one cluster exists, and corresponding
2495 // in other tracklet doesn't - return kFALSE
2496 if(curCluster || inCluster) return kFALSE;