#include <TBranch.h>
#include <TDirectory.h>
#include <TLinearFitter.h>
#include <TTree.h>
#include <TClonesArray.h>
#include <TTreeStream.h>
#include <TGeoMatrix.h>
#include <TGeoManager.h>
#include "AliLog.h"
#include "AliMathBase.h"
#include "AliESDEvent.h"
#include "AliGeomManager.h"
#include "AliRieman.h"
#include "AliTrackPointArray.h"
#include "AliTRDgeometry.h"
#include "AliTRDpadPlane.h"
#include "AliTRDcalibDB.h"
#include "AliTRDReconstructor.h"
#include "AliTRDCalibraFillHisto.h"
#include "AliTRDrecoParam.h"
#include "AliTRDcluster.h"
#include "AliTRDdigitsParam.h"
#include "AliTRDseedV1.h"
#include "AliTRDtrackV1.h"
#include "AliTRDtrackerV1.h"
#include "AliTRDtrackerDebug.h"
#include "AliTRDtrackingChamber.h"
#include "AliTRDchamberTimeBin.h"
ClassImp(AliTRDtrackerV1)
ClassImp(AliTRDtrackerV1::AliTRDLeastSquare)
ClassImp(AliTRDtrackerV1::AliTRDtrackFitterRieman)
AliTRDtrackerV1::ETRDtrackerV1BetheBloch AliTRDtrackerV1::fgBB = AliTRDtrackerV1::kGeant;
Double_t AliTRDtrackerV1::fgTopologicQA[kNConfigs] = {
0.5112, 0.5112, 0.5112, 0.0786, 0.0786,
0.0786, 0.0786, 0.0579, 0.0579, 0.0474,
0.0474, 0.0408, 0.0335, 0.0335, 0.0335
};
const Double_t AliTRDtrackerV1::fgkX0[kNPlanes] = {
300.2, 312.8, 325.4, 338.0, 350.6, 363.2};
Int_t AliTRDtrackerV1::fgNTimeBins = 0;
AliRieman* AliTRDtrackerV1::fgRieman = NULL;
TLinearFitter* AliTRDtrackerV1::fgTiltedRieman = NULL;
TLinearFitter* AliTRDtrackerV1::fgTiltedRiemanConstrained = NULL;
AliTRDtrackerV1::AliTRDtrackerV1(const AliTRDReconstructor *rec)
:AliTracker()
,fkReconstructor(NULL)
,fkRecoParam(NULL)
,fGeom(NULL)
,fClusters(NULL)
,fTracklets(NULL)
,fTracks(NULL)
,fTracksESD(NULL)
,fSieveSeeding(0)
,fEventInFile(-1)
{
SetReconstructor(rec);
if(!AliGeomManager::GetGeometry()){
AliFatal("Could not get geometry.");
}
fGeom = new AliTRDgeometry();
fGeom->CreateClusterMatrixArray();
TGeoHMatrix *matrix = NULL;
Double_t loc[] = {0., 0., 0.};
Double_t glb[] = {0., 0., 0.};
for(Int_t ily=kNPlanes; ily--;){
Int_t ism = 0;
while(!(matrix = fGeom->GetClusterMatrix(AliTRDgeometry::GetDetector(ily, 2, ism)))) ism++;
if(!matrix){
AliError(Form("Could not get transformation matrix for layer %d. Use default.", ily));
fR[ily] = fgkX0[ily];
continue;
}
matrix->LocalToMaster(loc, glb);
fR[ily] = glb[0]+ AliTRDgeometry::AnodePos()-.5*AliTRDgeometry::AmThick() - AliTRDgeometry::DrThick();
}
for (Int_t isector = 0; isector < AliTRDgeometry::kNsector; isector++) new(&fTrSec[isector]) AliTRDtrackingSector(fGeom, isector);
memset(fTrackQuality, 0, kMaxTracksStack*sizeof(Double_t));
memset(fSeedLayer, 0, kMaxTracksStack*sizeof(Int_t));
memset(fSeedTB, 0, kNSeedPlanes*sizeof(AliTRDchamberTimeBin*));
fTracksESD = new TClonesArray("AliESDtrack", 2*kMaxTracksStack);
fTracksESD->SetOwner();
}
AliTRDtrackerV1::~AliTRDtrackerV1()
{
if(fgRieman) delete fgRieman; fgRieman = NULL;
if(fgTiltedRieman) delete fgTiltedRieman; fgTiltedRieman = NULL;
if(fgTiltedRiemanConstrained) delete fgTiltedRiemanConstrained; fgTiltedRiemanConstrained = NULL;
for(Int_t isl =0; isl<kNSeedPlanes; isl++) if(fSeedTB[isl]) delete fSeedTB[isl];
if(fTracksESD){ fTracksESD->Delete(); delete fTracksESD; }
if(fTracks) {fTracks->Delete(); delete fTracks;}
if(fTracklets) {fTracklets->Delete(); delete fTracklets;}
if(IsClustersOwner() && fClusters) {
AliInfo(Form("tracker[%p] removing %d own clusters @ %p", (void*)this, fClusters->GetEntries(), (void*)fClusters));
fClusters->Delete(); delete fClusters;
}
if(fGeom) delete fGeom;
}
Int_t AliTRDtrackerV1::Clusters2Tracks(AliESDEvent *esd)
{
if(!fkRecoParam){
AliError("Reconstruction configuration not initialized. Call first AliTRDReconstructor::SetRecoParam().");
return 0;
}
Int_t ntracks = 0;
for(int ism=0; ism<AliTRDgeometry::kNsector; ism++){
ntracks += Clusters2TracksSM(ism, esd);
}
AliInfo(Form("Number of tracks: !TRDin[%d]", ntracks));
return ntracks;
}
Bool_t AliTRDtrackerV1::GetTrackPoint(Int_t index, AliTrackPoint &p) const
{
p.SetXYZ(0., 0., 0.);
AliTRDseedV1 *tracklet = GetTracklet(index);
if (!tracklet) return kFALSE;
Int_t det = tracklet->GetDetector();
Int_t sec = fGeom->GetSector(det);
Double_t alpha = (sec+.5)*AliTRDgeometry::GetAlpha(),
sinA = TMath::Sin(alpha),
cosA = TMath::Cos(alpha);
Double_t local[3];
local[0] = tracklet->GetX();
local[1] = tracklet->GetY();
local[2] = tracklet->GetZ();
Double_t global[3];
fGeom->RotateBack(det, local, global);
Double_t cov2D[3]; Float_t cov[6];
tracklet->GetCovAt(local[0], cov2D);
cov[0] = cov2D[0]*sinA*sinA;
cov[1] =-cov2D[0]*sinA*cosA;
cov[2] =-cov2D[1]*sinA;
cov[3] = cov2D[0]*cosA*cosA;
cov[4] = cov2D[1]*cosA;
cov[5] = cov2D[2];
p.SetXYZ(global[0],global[1],global[2], cov);
AliGeomManager::ELayerID iLayer = AliGeomManager::ELayerID(AliGeomManager::kTRD1+fGeom->GetLayer(det));
Int_t modId = fGeom->GetSector(det) * AliTRDgeometry::kNstack + fGeom->GetStack(det);
UShort_t volid = AliGeomManager::LayerToVolUID(iLayer, modId);
p.SetVolumeID(volid);
return kTRUE;
}
TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitter()
{
if(!fgTiltedRieman) fgTiltedRieman = new TLinearFitter(4, "hyp4");
return fgTiltedRieman;
}
TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitterConstraint()
{
if(!fgTiltedRiemanConstrained) fgTiltedRiemanConstrained = new TLinearFitter(2, "hyp2");
return fgTiltedRiemanConstrained;
}
AliRieman* AliTRDtrackerV1::GetRiemanFitter()
{
if(!fgRieman) fgRieman = new AliRieman(AliTRDseedV1::kNtb * AliTRDgeometry::kNlayer);
return fgRieman;
}
Int_t AliTRDtrackerV1::PropagateBack(AliESDEvent *event)
{
if(!fClusters || !fClusters->GetEntriesFast()){
AliInfo("No TRD clusters");
return 0;
}
AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
if (!calibra) AliInfo("Could not get Calibra instance");
if (!fgNTimeBins) fgNTimeBins = fkReconstructor->GetNTimeBins();
Int_t nFound = 0,
nBacked = 0,
nSeeds = 0,
nTRDseeds= 0,
nTPCseeds= 0;
Float_t foundMin = 20.0;
Float_t *quality = NULL;
Int_t *index = NULL;
fEventInFile = event->GetEventNumberInFile();
nSeeds = event->GetNumberOfTracks();
if(nSeeds){
quality = new Float_t[nSeeds];
index = new Int_t[4*nSeeds];
for (Int_t iSeed = nSeeds; iSeed--;) {
AliESDtrack *seed = event->GetTrack(iSeed);
Double_t covariance[15];
seed->GetExternalCovariance(covariance);
quality[iSeed] = covariance[0] + covariance[2];
}
TMath::Sort(nSeeds, quality, index,kFALSE);
}
Int_t expectedClr;
AliTRDtrackV1 track;
for (Int_t iSeed = 0; iSeed < nSeeds; iSeed++) {
AliESDtrack *seed = event->GetTrack(index[iSeed]);
Float_t p4 = seed->GetC(seed->GetBz());
ULong_t status = seed->GetStatus();
if ((status & AliESDtrack::kTRDout) != 0) continue;
if ((status & AliESDtrack::kTPCout)){
AliDebug(3, Form("Prolongate seed[%2d] which is TPC.", iSeed));
} else continue;
track.~AliTRDtrackV1();
new(&track) AliTRDtrackV1(*seed);
if(AliTRDgeometry::GetXtrdBeg() > (AliTRDReconstructor::GetMaxStep() + track.GetX()) && !PropagateToX(track, AliTRDgeometry::GetXtrdBeg(), AliTRDReconstructor::GetMaxStep())){
seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
continue;
}
if(!AdjustSector(&track)){
seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
continue;
}
if(TMath::Abs(track.GetSnp()) > AliTRDReconstructor::GetMaxSnp()) {
seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
continue;
}
nTPCseeds++;
AliDebug(2, Form("TRD propagate TPC seed[%d] = %d.", iSeed, index[iSeed]));
seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
track.SetReconstructor(fkReconstructor);
track.SetKink(Bool_t(seed->GetKinkIndex(0)));
track.SetPrimary(status & AliESDtrack::kTPCin);
expectedClr = FollowBackProlongation(track);
if(track.GetTrackIn()){
seed->UpdateTrackParams(&track, AliESDtrack::kTRDin);
nTRDseeds++;
}
if (expectedClr<0){
seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
continue;
} else {
nFound++;
track.CookPID();
track.CookLabel(1. - AliTRDReconstructor::GetLabelFraction());
if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(&track);
if (fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 0 || AliTRDReconstructor::GetStreamLevel()>0 ) {
AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(track);
calibTrack->SetOwner();
seed->AddCalibObject(calibTrack);
}
seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
track.UpdateESDtrack(seed);
}
if ((TMath::Abs(track.GetC(track.GetBz()) - p4) / TMath::Abs(p4) < 0.2) || (track.Pt() > 0.8)) {
Int_t foundClr = track.GetNumberOfClusters();
if (foundClr >= foundMin) {
if (track.GetChi2() / track.GetNumberOfClusters() < 4) {
}
Bool_t isGold = kFALSE;
if (track.GetChi2() / track.GetNumberOfClusters() < 5) {
if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
nBacked++;
isGold = kTRUE;
}
if ((!isGold) && (track.GetNCross() == 0) && (track.GetChi2() / track.GetNumberOfClusters() < 7)) {
if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
nBacked++;
isGold = kTRUE;
}
if ((!isGold) && (track.GetBackupTrack())) {
if ((track.GetBackupTrack()->GetNumberOfClusters() > foundMin) && ((track.GetBackupTrack()->GetChi2()/(track.GetBackupTrack()->GetNumberOfClusters()+1)) < 7)) {
seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
nBacked++;
isGold = kTRUE;
}
}
}
}
if(!(seed->GetStatus()&AliESDtrack::kTRDStop)) {
Int_t sm = track.GetSector();
Double_t xtof = 371.;
if(gGeoManager){
gGeoManager->cd(Form("/ALIC_1/B077_1/BSEGMO%d_1/BTOF%d_1", sm, sm));
TGeoHMatrix *m = NULL;
Double_t loc[]={0., 0., -.5*29.05}, glob[3];
if((m=gGeoManager->GetCurrentMatrix())){
m->LocalToMaster(loc, glob);
xtof = TMath::Sqrt(glob[0]*glob[0]+glob[1]*glob[1]);
}
}
if(xtof > (AliTRDReconstructor::GetMaxStep() + track.GetX()) && !PropagateToX(track, xtof, AliTRDReconstructor::GetMaxStep())){
seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
continue;
}
if(!AdjustSector(&track)){
seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
continue;
}
if(TMath::Abs(track.GetSnp()) > AliTRDReconstructor::GetMaxSnp()){
seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
continue;
}
seed->SetTRDQuality(track.StatusForTOF());
}
seed->SetTRDBudget(track.GetBudget(0));
}
if(index) delete [] index;
if(quality) delete [] quality;
AliInfo(Form("Number of seeds: TPCout[%d] TRDin[%d]", nTPCseeds, nTRDseeds));
AliInfo(Form("Number of tracks: TRDout[%d] TRDbackup[%d]", nFound, nBacked));
if (fkReconstructor->IsSeeding()) Clusters2Tracks(event);
return 0;
}
Int_t AliTRDtrackerV1::RefitInward(AliESDEvent *event)
{
Int_t nseed = 0;
Int_t found = 0;
if(!fClusters || !fClusters->GetEntriesFast()){
AliInfo("No TRD clusters");
return 0;
}
AliTRDtrackV1 track;
for (Int_t itrack = 0; itrack < event->GetNumberOfTracks(); itrack++) {
AliESDtrack *seed = event->GetTrack(itrack);
ULong_t status = seed->GetStatus();
new(&track) AliTRDtrackV1(*seed);
if (track.GetX() < 270.0) {
seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
continue;
}
if(!(status & AliESDtrack::kTRDout)) continue;
if(!(status & AliESDtrack::kTRDin)) continue;
nseed++;
track.ResetCovariance(50.0);
Bool_t kUPDATE = kFALSE;
Double_t xTPC = 250.0;
if(FollowProlongation(track)){
if (fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 0 || AliTRDReconstructor::GetStreamLevel()>0 ){
TObject *o = NULL; Int_t ic = 0;
AliTRDtrackV1 *calibTrack = NULL;
while((o = seed->GetCalibObject(ic++))){
if(!(calibTrack = dynamic_cast<AliTRDtrackV1*>(o))) continue;
calibTrack->SetTrackOut(&track);
}
}
if (PropagateToX(track, xTPC, AliTRDReconstructor::GetMaxStep())) {
seed->UpdateTrackParams(&track, AliESDtrack::kTRDrefit);
found++;
kUPDATE = kTRUE;
}
}
if(!kUPDATE) {
AliTRDtrackV1 tt(*seed);
if (PropagateToX(tt, xTPC, AliTRDReconstructor::GetMaxStep())) seed->UpdateTrackParams(&tt, AliESDtrack::kTRDbackup);
}
}
AliInfo(Form("Number of seeds: TRDout[%d]", nseed));
AliInfo(Form("Number of tracks: TRDrefit[%d]", found));
return 0;
}
Int_t AliTRDtrackerV1::FollowProlongation(AliTRDtrackV1 &t)
{
Int_t nClustersExpected = 0;
for (Int_t iplane = kNPlanes; iplane--;) {
Int_t index(-1);
AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
AliDebug(2, Form("Tracklet[%p] ly[%d] idx[%d]", (void*)tracklet, iplane, index));
if(!tracklet) continue;
if(!tracklet->IsOK()){
AliDebug(1, Form("Tracklet Det[%d] !OK", tracklet->GetDetector()));
continue;
}
Double_t x = tracklet->GetX();
if(x > t.GetX()+AliTRDReconstructor::GetMaxStep()) continue;
t.SetTracklet(tracklet, index);
if (x < (t.GetX()-AliTRDReconstructor::GetMaxStep()) && !PropagateToX(t, x+AliTRDReconstructor::GetMaxStep(), AliTRDReconstructor::GetMaxStep())) break;
if (!AdjustSector(&t)) break;
Double_t xyz0[3];
t.GetXYZ(xyz0);
Double_t alpha = t.GetAlpha(), y, z;
if (!t.GetProlongation(x,y,z)) break;
Double_t xyz1[3];
xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
xyz1[1] = x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
xyz1[2] = z;
Double_t length = TMath::Sqrt(
(xyz0[0]-xyz1[0])*(xyz0[0]-xyz1[0]) +
(xyz0[1]-xyz1[1])*(xyz0[1]-xyz1[1]) +
(xyz0[2]-xyz1[2])*(xyz0[2]-xyz1[2])
);
if(length>0.){
Double_t param[7];
if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param)<=0.) break;
Double_t xrho= param[0]*param[4];
Double_t xx0 = param[1];
t.PropagateTo(x, xx0, xrho);
if (!AdjustSector(&t)) break;
}
Double_t cov[3]; tracklet->GetCovAt(x, cov);
Double_t p[2] = { tracklet->GetY(), tracklet->GetZ()};
Double_t chi2 = ((AliExternalTrackParam)t).GetPredictedChi2(p, cov);
if(fkReconstructor->IsDebugStreaming()){
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
AliExternalTrackParam param0(t);
AliExternalTrackParam param1(t);
param1.Update(p, cov);
TVectorD vcov(3,cov);
TVectorD vpar(3,p);
cstreamer << "FollowProlongationInfo"
<< "EventNumber=" << eventNumber
<< "iplane="<<iplane
<< "vcov.="<<&vcov
<< "vpar.="<<&vpar
<< "tracklet.=" << tracklet
<< "chi2="<< chi2
<< "param0.=" << ¶m0
<< "param1.=" << ¶m1
<< "\n";
}
if (chi2 < fkRecoParam->GetChi2Cut() && ((AliExternalTrackParam&)t).Update(p, cov)){
t.SetNumberOfClusters();
t.UpdateChi2(chi2);
nClustersExpected += tracklet->GetN();
}
}
if(fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 1 || AliTRDReconstructor::GetStreamLevel()>1){
Int_t index;
for(int iplane=0; iplane<AliTRDgeometry::kNlayer; iplane++){
AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
if(!tracklet) continue;
t.SetTracklet(tracklet, index);
}
if(fkReconstructor->IsDebugStreaming()){
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
AliTRDtrackV1 track(t);
track.SetOwner();
cstreamer << "FollowProlongation"
<< "EventNumber=" << eventNumber
<< "ncl=" << nClustersExpected
<< "track.=" << &track
<< "\n";
}
}
return nClustersExpected;
}
Int_t AliTRDtrackerV1::FollowBackProlongation(AliTRDtrackV1 &t)
{
Int_t n = 0;
Double_t driftLength = .5*AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick();
AliTRDtrackingChamber *chamber = NULL;
Int_t debugLevel = fkReconstructor->IsDebugStreaming() ? fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) : 0;
if ( AliTRDReconstructor::GetStreamLevel()>0) debugLevel= AliTRDReconstructor::GetStreamLevel();
TTreeSRedirector *cstreamer = fkReconstructor->IsDebugStreaming() ? fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker) : 0x0;
Bool_t kStoreIn(kTRUE),
kStandAlone(kFALSE),
kUseTRD(fkRecoParam->IsOverPtThreshold(t.Pt()));
Int_t startLayer(0);
AliTRDseedV1 tracklet, *ptrTracklet = NULL;
AliTRDseedV1 *tracklets[kNPlanes];
memset(tracklets, 0, sizeof(AliTRDseedV1 *) * kNPlanes);
for(Int_t ip(kNPlanes); ip--;){
if(!(tracklets[ip] = t.GetTracklet(ip))) continue;
t.UnsetTracklet(ip);
if(tracklets[ip]->IsOK()) startLayer=ip;
kStandAlone = kTRUE;
kUseTRD = kTRUE;
}
AliDebug(4, Form("SA[%c] Start[%d]\n"
" [0]idx[%d] traklet[%p]\n"
" [1]idx[%d] traklet[%p]\n"
" [2]idx[%d] traklet[%p]\n"
" [3]idx[%d] traklet[%p]\n"
" [4]idx[%d] traklet[%p]\n"
" [5]idx[%d] traklet[%p]"
, kStandAlone?'y':'n', startLayer
, t.GetTrackletIndex(0), (void*)tracklets[0]
, t.GetTrackletIndex(1), (void*)tracklets[1]
, t.GetTrackletIndex(2), (void*)tracklets[2]
, t.GetTrackletIndex(3), (void*)tracklets[3]
, t.GetTrackletIndex(4), (void*)tracklets[4]
, t.GetTrackletIndex(5), (void*)tracklets[5]));
TGeoHMatrix *matrix = NULL;
Double_t x(0.), y(0.), z(0.);
for (Int_t ily=startLayer, sm=-1, stk=-1, det=-1; ily < AliTRDgeometry::kNlayer; ily++) {
AliDebug(2, Form("Propagate to x[%d] = %7.2f", ily, fR[ily]));
if (!t.GetProlongation(fR[ily], y, z)){
n=-1;
t.SetErrStat(AliTRDtrackV1::kProlongation);
AliDebug(4, Form("Failed Rough Prolongation to ly[%d] x[%7.2f] y[%7.2f] z[%7.2f]", ily, fR[ily], y, z));
break;
}
sm = t.GetSector();
stk = fGeom->GetStack(z, ily);
det = stk>=0 ? AliTRDgeometry::GetDetector(ily, stk, sm) : -1;
matrix = det>=0 ? fGeom->GetClusterMatrix(det) : NULL;
if( !fGeom->GetSMstatus(sm) ||
stk<0. ||
fGeom->IsHole(ily, stk, sm) ||
!matrix ){
AliDebug(4, Form("Missing Geometry ly[%d]. Guess radial position", ily));
if(fR[ily] > (AliTRDReconstructor::GetMaxStep() + t.GetX()) && !PropagateToX(t, fR[ily], AliTRDReconstructor::GetMaxStep())){
n=-1;
t.SetErrStat(AliTRDtrackV1::kPropagation);
AliDebug(4, "Failed Propagation [Missing Geometry]");
break;
}
if(!AdjustSector(&t)){
n=-1;
t.SetErrStat(AliTRDtrackV1::kAdjustSector);
AliDebug(4, "Failed Adjust Sector [Missing Geometry]");
break;
}
if(TMath::Abs(t.GetSnp()) > AliTRDReconstructor::GetMaxSnp()){
n=-1;
t.SetErrStat(AliTRDtrackV1::kSnp);
AliDebug(4, "Failed Max Snp [Missing Geometry]");
break;
}
t.SetErrStat(AliTRDtrackV1::kGeometry, ily);
continue;
}
Double_t loc[] = {AliTRDgeometry::AnodePos()- driftLength, 0., 0.};
Double_t glb[] = {0., 0., 0.};
matrix->LocalToMaster(loc, glb);
AliDebug(3, Form("Propagate to det[%3d] x_anode[%7.2f] (%f %f)", det, glb[0]+driftLength, glb[1], glb[2]));
x = glb[0] - AliTRDReconstructor::GetMaxStep();
if(x > (AliTRDReconstructor::GetMaxStep() + t.GetX()) && !PropagateToX(t, x, AliTRDReconstructor::GetMaxStep())){
n=-1;
t.SetErrStat(AliTRDtrackV1::kPropagation);
AliDebug(4, Form("Failed Initial Propagation to x[%7.2f]", x));
break;
}
if(!AdjustSector(&t)){
n=-1;
t.SetErrStat(AliTRDtrackV1::kAdjustSector);
AliDebug(4, "Failed Adjust Sector Start");
break;
}
if(TMath::Abs(t.GetSnp()) > AliTRDReconstructor::GetMaxSnp()) {
n=-1;
t.SetErrStat(AliTRDtrackV1::kSnp);
AliDebug(4, Form("Failed Max Snp[%f] MaxSnp[%f]", t.GetSnp(), AliTRDReconstructor::GetMaxSnp()));
break;
}
Bool_t doRecalculate = kFALSE;
if(sm != t.GetSector()){
sm = t.GetSector();
doRecalculate = kTRUE;
}
if(stk != fGeom->GetStack(z, ily)){
stk = fGeom->GetStack(z, ily);
doRecalculate = kTRUE;
}
if(doRecalculate){
det = AliTRDgeometry::GetDetector(ily, stk, sm);
if(!(matrix = fGeom->GetClusterMatrix(det))){
t.SetErrStat(AliTRDtrackV1::kGeometry, ily);
AliDebug(4, Form("Failed Geometry Matrix ly[%d]", ily));
continue;
}
matrix->LocalToMaster(loc, glb);
x = glb[0] - AliTRDReconstructor::GetMaxStep();
}
if (!t.GetProlongation(x+AliTRDReconstructor::GetMaxStep(), y, z)) {
n=-1;
t.SetErrStat(AliTRDtrackV1::kProlongation);
AliDebug(4, Form("Failed Prolongation to x[%7.2f] y[%7.2f] z[%7.2f]", x+AliTRDReconstructor::GetMaxStep(), y, z));
break;
}
if(fGeom->IsOnBoundary(det, y, z, .5)){
t.SetErrStat(AliTRDtrackV1::kBoundary, ily);
AliDebug(4, "Failed Track on Boundary");
continue;
}
Float_t prod(t.GetBz()*t.Charge());
ptrTracklet = tracklets[ily];
if(!ptrTracklet){
AliDebug(3, Form("Building tracklet det[%d]", det));
if(!fTrSec[sm].GetNChambers()){
t.SetErrStat(AliTRDtrackV1::kNoClusters, ily);
AliDebug(4, "Failed NoClusters");
continue;
}
if(fTrSec[sm].GetX(ily) < 1.){
t.SetErrStat(AliTRDtrackV1::kNoClusters, ily);
AliDebug(4, "Failed NoX");
continue;
}
if(!(chamber = fTrSec[sm].GetChamber(stk, ily))){
t.SetErrStat(AliTRDtrackV1::kNoClusters, ily);
AliDebug(4, "Failed No Detector");
continue;
}
if(chamber->GetNClusters() < fgNTimeBins*fkRecoParam ->GetFindableClusters()){
t.SetErrStat(AliTRDtrackV1::kNoClusters, ily);
AliDebug(4, "Failed Not Enough Clusters in Detector");
continue;
}
tracklet.~AliTRDseedV1();
ptrTracklet = new(&tracklet) AliTRDseedV1(det);
ptrTracklet->SetReconstructor(fkReconstructor);
ptrTracklet->SetKink(t.IsKink());
ptrTracklet->SetPrimary(t.IsPrimary());
ptrTracklet->SetPadPlane(fGeom->GetPadPlane(ily, stk));
ptrTracklet->SetX0(glb[0]+driftLength);
if(!ptrTracklet->Init(&t)){
n=-1;
t.SetErrStat(AliTRDtrackV1::kTrackletInit);
AliDebug(4, "Failed Tracklet Init");
break;
}
if(!ptrTracklet->AttachClusters(chamber, kTRUE, prod<0.?kTRUE:kFALSE, fEventInFile)){
t.SetErrStat(AliTRDtrackV1::kNoAttach, ily);
if(debugLevel>3){
AliTRDseedV1 trackletCp(*ptrTracklet);
UChar_t status(t.GetStatusTRD(ily));
(*cstreamer) << "FollowBackProlongation4"
<<"status=" << status
<<"tracklet.=" << &trackletCp
<< "\n";
}
AliDebug(4, "Failed Attach Clusters");
continue;
}
AliDebug(3, Form("Number of Clusters in Tracklet: %d", ptrTracklet->GetN()));
if(ptrTracklet->GetN() < fgNTimeBins*fkRecoParam->GetFindableClusters()){
t.SetErrStat(AliTRDtrackV1::kNoClustersTracklet, ily);
if(debugLevel>3){
AliTRDseedV1 trackletCp(*ptrTracklet);
UChar_t status(t.GetStatusTRD(ily));
(*cstreamer) << "FollowBackProlongation4"
<<"status=" << status
<<"tracklet.=" << &trackletCp
<< "\n";
}
AliDebug(4, "Failed N Clusters Attached");
continue;
}
ptrTracklet->UpdateUsed();
} else AliDebug(2, Form("Use external tracklet ly[%d]", ily));
if(!ptrTracklet->FitRobust(fGeom->GetPadPlane(ily, stk), matrix, t.GetBz(), t.Charge())){
t.SetErrStat(AliTRDtrackV1::kNoFit, ily);
AliDebug(4, "Failed Tracklet Fit");
continue;
}
ptrTracklet->SetXYZ(matrix);
x = ptrTracklet->GetX();
if(x > (AliTRDReconstructor::GetMaxStep() + t.GetX()) && !PropagateToX(t, x, AliTRDReconstructor::GetMaxStep())) {
n=-1;
t.SetErrStat(AliTRDtrackV1::kPropagation);
AliDebug(4, Form("Failed Propagation to Tracklet x[%7.2f]", x));
break;
}
if(!AdjustSector(&t)) {
n=-1;
t.SetErrStat(AliTRDtrackV1::kAdjustSector);
AliDebug(4, "Failed Adjust Sector");
break;
}
if(TMath::Abs(t.GetSnp()) > AliTRDReconstructor::GetMaxSnp()) {
n=-1;
t.SetErrStat(AliTRDtrackV1::kSnp);
AliDebug(4, Form("Failed Max Snp[%f] MaxSnp[%f]", t.GetSnp(), AliTRDReconstructor::GetMaxSnp()));
break;
}
Double_t cov[3]; ptrTracklet->GetCovAt(x, cov);
Double_t p[2] = { ptrTracklet->GetY(), ptrTracklet->GetZ()};
Double_t chi2 = ((AliExternalTrackParam)t).GetPredictedChi2(p, cov);
if(fkReconstructor->IsDebugStreaming()){
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
AliExternalTrackParam param0(t);
AliExternalTrackParam param1(t);
param1.Update(p, cov);
TVectorD vcov(3,cov);
TVectorD vpar(3,p);
(*cstreamer) << "FollowBackProlongationInfo"
<< "EventNumber=" << eventNumber
<< "chi2="<<chi2
<< "iplane="<<ily
<< "vcov.="<<&vcov
<< "vpar.="<<&vpar
<< "tracklet.=" << ptrTracklet
<< "param0.=" << ¶m0
<< "param1.=" << ¶m1
<< "\n";
}
if (chi2> fkRecoParam->GetChi2Cut()){
t.SetErrStat(AliTRDtrackV1::kChi2, ily);
if(debugLevel > 2){
UChar_t status(t.GetStatusTRD());
AliTRDseedV1 trackletCp(*ptrTracklet);
AliTRDtrackV1 trackCp(t);
trackCp.SetOwner();
(*cstreamer) << "FollowBackProlongation3"
<< "status=" << status
<< "tracklet.=" << &trackletCp
<< "track.=" << &trackCp
<< "\n";
}
AliDebug(4, Form("Failed Chi2[%f]", chi2));
continue;
}
if(kStoreIn){
t.SetTrackIn();
kStoreIn = kFALSE;
}
if(kUseTRD){
if(!((AliExternalTrackParam&)t).Update(p, cov)) {
n=-1;
t.SetErrStat(AliTRDtrackV1::kUpdate);
if(debugLevel > 2){
UChar_t status(t.GetStatusTRD());
AliTRDseedV1 trackletCp(*ptrTracklet);
AliTRDtrackV1 trackCp(t);
trackCp.SetOwner();
(*cstreamer) << "FollowBackProlongation3"
<< "status=" << status
<< "tracklet.=" << &trackletCp
<< "track.=" << &trackCp
<< "\n";
}
AliDebug(4, Form("Failed Track Update @ y[%7.2f] z[%7.2f] s2y[%f] s2z[%f] covyz[%f]", p[0], p[1], cov[0], cov[2], cov[1]));
break;
}
}
if(!kStandAlone) ptrTracklet->UseClusters();
AliTracker::FillResiduals(&t, p, cov, ptrTracklet->GetVolumeId());
ptrTracklet = SetTracklet(ptrTracklet);
Int_t index(fTracklets->GetEntriesFast()-1);
t.SetTracklet(ptrTracklet, index);
t.SetNumberOfClusters();
t.UpdateChi2(chi2);
n += ptrTracklet->GetN();
AliDebug(2, Form("Setting Tracklet[%d] @ Idx[%d]", ily, index));
Int_t failed(0);
if(!kStandAlone && (failed = t.MakeBackupTrack())) AliDebug(2, Form("Failed backup on cut[%d]", failed));
}
if(n && debugLevel > 1){
AliTRDtrackV1 track(t);
track.SetOwner();
(*cstreamer) << "FollowBackProlongation2"
<< "EventNumber=" << fEventInFile
<< "track.=" << &track
<< "\n";
}
return n;
}
Float_t AliTRDtrackerV1::FitRieman(AliTRDseedV1 *tracklets, Double_t *chi2, Int_t *const planes){
AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
fitter->Reset();
Int_t allplanes[] = {0, 1, 2, 3, 4, 5};
Int_t *ppl = &allplanes[0];
Int_t maxLayers = 6;
if(planes){
maxLayers = 4;
ppl = planes;
}
for(Int_t il = 0; il < maxLayers; il++){
if(!tracklets[ppl[il]].IsOK()) continue;
fitter->AddPoint(tracklets[ppl[il]].GetX0(), tracklets[ppl[il]].GetYfit(0), tracklets[ppl[il]].GetZfit(0),1,10);
}
fitter->Update();
memset(chi2, 0, sizeof(Double_t) * 2);
for(Int_t il = 0; il < maxLayers; il++){
tracklets[ppl[il]].Init(fitter);
if((!tracklets[ppl[il]].IsOK()) && (!planes)) continue;
chi2[0] += tracklets[ppl[il]].GetChi2Y();
chi2[1] += tracklets[ppl[il]].GetChi2Z();
}
return fitter->GetC();
}
void AliTRDtrackerV1::FitRieman(AliTRDcluster **seedcl, Double_t chi2[2])
{
AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
fitter->Reset();
for(Int_t i = 0; i < 4; i++){
fitter->AddPoint(seedcl[i]->GetX(), seedcl[i]->GetY(), seedcl[i]->GetZ(), 1., 10.);
}
fitter->Update();
chi2[0] = 0; chi2[1] = 0;
for(Int_t ipl = 0; ipl < kNSeedPlanes; ipl++){
chi2[0] += (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX()));
chi2[1] += (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX()));
}
}
Float_t AliTRDtrackerV1::FitTiltedRiemanConstraint(AliTRDseedV1 *tracklets, Double_t zVertex)
{
TLinearFitter *fitter = GetTiltedRiemanFitterConstraint();
fitter->StoreData(kTRUE);
fitter->ClearPoints();
AliTRDcluster *cl = NULL;
Float_t x, y, z, w, t, error, tilt;
Double_t uvt[2];
Int_t nPoints = 0;
for(Int_t ilr = 0; ilr < AliTRDgeometry::kNlayer; ilr++){
if(!tracklets[ilr].IsOK()) continue;
for(Int_t itb = 0; itb < AliTRDseedV1::kNclusters; itb++){
if(!tracklets[ilr].IsUsable(itb)) continue;
if(!(cl = tracklets[ilr].GetClusters(itb))) continue;
if(!cl->IsInChamber()) continue;
x = cl->GetX();
y = cl->GetY();
z = cl->GetZ();
tilt = tracklets[ilr].GetTilt();
t = 1./(x * x + y * y);
uvt[0] = 2. * x * t;
uvt[1] = 2. * x * t * tilt ;
w = 2. * (y + tilt * (z - zVertex)) * t;
error = 2. * TMath::Sqrt(cl->GetSigmaY2()+tilt*tilt*cl->GetSigmaZ2()) * t;
fitter->AddPoint(uvt, w, error);
nPoints++;
}
}
fitter->Eval();
Double_t a = fitter->GetParameter(0);
Double_t b = fitter->GetParameter(1);
Double_t curvature = a/TMath::Sqrt(b*b + 1);
Float_t chi2track = 0.0;
if (nPoints > 0) {
chi2track = fitter->GetChisquare()/Double_t(nPoints);
}
for(Int_t ip = 0; ip < AliTRDtrackerV1::kNPlanes; ip++)
tracklets[ip].SetC(curvature, 1);
if(AliLog::GetDebugLevel("TRD", "AliTRDtrackerV1")>3) printf("D-AliTRDtrackerV1::FitTiltedRiemanConstraint: Chi2[%f] C[%5.2e] pt[%8.3f]\n", chi2track, curvature, GetBz()*kB2C/curvature);
return chi2track;
}
Float_t AliTRDtrackerV1::FitTiltedRieman(AliTRDseedV1 *tracklets, Bool_t sigError)
{
TLinearFitter *fitter = GetTiltedRiemanFitter();
fitter->StoreData(kTRUE);
fitter->ClearPoints();
AliTRDLeastSquare zfitter;
AliTRDcluster *cl = NULL;
Double_t xref = CalculateReferenceX(tracklets);
Double_t x, y, z, t, tilt, dx, w, we, erry, errz;
Double_t uvt[4], sumPolY[5], sumPolZ[3];
memset(sumPolY, 0, sizeof(Double_t) * 5);
memset(sumPolZ, 0, sizeof(Double_t) * 3);
Int_t nPoints = 0;
for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
if(!tracklets[ipl].IsOK()) continue;
tilt = tracklets[ipl].GetTilt();
for(Int_t itb = 0; itb < AliTRDseedV1::kNclusters; itb++){
if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
if(!cl->IsInChamber()) continue;
if (!tracklets[ipl].IsUsable(itb)) continue;
x = cl->GetX();
y = cl->GetY();
z = cl->GetZ();
dx = x - xref;
t = 1./(x*x + y*y);
uvt[0] = 2. * x * t;
uvt[1] = t;
uvt[2] = 2. * tilt * t;
uvt[3] = 2. * tilt * dx * t;
w = 2. * (y + tilt*z) * t;
we = 2. * t;
we *= sigError ? TMath::Sqrt(cl->GetSigmaY2()+tilt*tilt*cl->GetSigmaZ2()) : 0.2;
fitter->AddPoint(uvt, w, we);
zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
erry = 1./(TMath::Sqrt(cl->GetSigmaY2()) + 0.1);
erry *= erry;
errz = 1./cl->GetSigmaZ2();
for(Int_t ipol = 0; ipol < 5; ipol++){
sumPolY[ipol] += erry;
erry *= x;
if(ipol < 3){
sumPolZ[ipol] += errz;
errz *= x;
}
}
nPoints++;
}
}
if (fitter->Eval()) return 1.e10;
zfitter.Eval();
Double_t offset = fitter->GetParameter(3);
Double_t slope = fitter->GetParameter(4);
Bool_t acceptablez = kTRUE;
Double_t zref = 0.0;
for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
if(!tracklets[iLayer].IsOK()) continue;
zref = offset + slope * (tracklets[iLayer].GetX0() - xref);
if (TMath::Abs(tracklets[iLayer].GetZfit(0) - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
acceptablez = kFALSE;
}
if (!acceptablez) {
Double_t dzmf = zfitter.GetFunctionParameter(1);
Double_t zmf = zfitter.GetFunctionValue(&xref);
fgTiltedRieman->FixParameter(3, zmf);
fgTiltedRieman->FixParameter(4, dzmf);
fitter->Eval();
fitter->ReleaseParameter(3);
fitter->ReleaseParameter(4);
offset = fitter->GetParameter(3);
slope = fitter->GetParameter(4);
}
Double_t a = fitter->GetParameter(0);
Double_t b = fitter->GetParameter(1);
Double_t c = fitter->GetParameter(2);
Double_t curvature = 1.0 + b*b - c*a;
if (curvature > 0.0) curvature = a / TMath::Sqrt(curvature);
Double_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
TMatrixD covarPolY(3,3);
covarPolY(0,0) = sumPolY[0]; covarPolY(1,1) = sumPolY[2]; covarPolY(2,2) = sumPolY[4];
covarPolY(0,1) = covarPolY(1,0) = sumPolY[1];
covarPolY(0,2) = covarPolY(2,0) = sumPolY[2];
covarPolY(2,1) = covarPolY(1,2) = sumPolY[3];
covarPolY.Invert();
TMatrixD covarPolZ(2,2);
covarPolZ(0,0) = sumPolZ[0]; covarPolZ(1,1) = sumPolZ[2];
covarPolZ(1,0) = covarPolZ(0,1) = sumPolZ[1];
covarPolZ.Invert();
Double_t dy, dz;
Double_t cov[15];
memset(cov, 0, sizeof(Double_t) * 15);
for(Int_t iLayer = 0; iLayer < AliTRDtrackerV1::kNPlanes; iLayer++) {
x = tracklets[iLayer].GetX0();
y = 0;
z = 0;
dy = 0;
dz = 0;
memset(cov, 0, sizeof(Double_t) * 3);
TMatrixD transform(3,3);
transform(0,0) = 1;
transform(0,1) = x;
transform(0,2) = x*x;
transform(1,1) = 1;
transform(1,2) = x;
transform(2,2) = 1;
TMatrixD covariance(transform, TMatrixD::kMult, covarPolY);
covariance *= transform.T();
TMatrixD transformZ(2,2);
transformZ(0,0) = transformZ(1,1) = 1;
transformZ(0,1) = x;
TMatrixD covarZ(transformZ, TMatrixD::kMult, covarPolZ);
covarZ *= transformZ.T();
Double_t res = (x * a + b);
res *= res;
res = 1.0 - c * a + b * b - res;
if (res >= 0) {
res = TMath::Sqrt(res);
y = (1.0 - res) / a;
}
cov[0] = covariance(0,0);
cov[2] = covarZ(0,0);
cov[1] = 0.;
Double_t x0 = -b / a;
if (-c * a + b * b + 1 > 0) {
if (1.0/(curvature * curvature) - (x - x0) * (x - x0) > 0.0) {
Double_t yderiv = (x - x0) / TMath::Sqrt(1.0/(curvature * curvature) - (x - x0) * (x - x0));
if (a < 0) yderiv *= -1.0;
dy = yderiv;
}
}
z = offset + slope * (x - xref);
dz = slope;
tracklets[iLayer].SetYref(0, y);
tracklets[iLayer].SetYref(1, dy);
tracklets[iLayer].SetZref(0, z);
tracklets[iLayer].SetZref(1, dz);
tracklets[iLayer].SetC(curvature);
tracklets[iLayer].SetCovRef(cov);
tracklets[iLayer].SetChi2(chi2track);
}
if(AliLog::GetDebugLevel("TRD", "AliTRDtrackerV1")>3) printf("D-AliTRDtrackerV1::FitTiltedRieman: Chi2[%f] C[%5.2e] pt[%8.3f]\n", chi2track, curvature, GetBz()*kB2C/curvature);
return chi2track;
}
Double_t AliTRDtrackerV1::FitLine(const AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t err, Int_t np, AliTrackPoint *points)
{
AliTRDLeastSquare yfitter, zfitter;
AliTRDcluster *cl = NULL;
AliTRDseedV1 work[kNPlanes], *tracklet = NULL;
if(!tracklets){
for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
if(!(tracklet = track->GetTracklet(ipl))) continue;
if(!tracklet->IsOK()) continue;
new(&work[ipl]) AliTRDseedV1(*tracklet);
}
tracklets = &work[0];
}
Double_t xref = CalculateReferenceX(tracklets);
Double_t x, y, z, dx, ye, yr, tilt;
for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
if(!tracklets[ipl].IsOK()) continue;
for(Int_t itb = 0; itb < fgNTimeBins; itb++){
if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
if (!tracklets[ipl].IsUsable(itb)) continue;
x = cl->GetX();
z = cl->GetZ();
dx = x - xref;
zfitter.AddPoint(&dx, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
}
}
zfitter.Eval();
Double_t z0 = zfitter.GetFunctionParameter(0);
Double_t dzdx = zfitter.GetFunctionParameter(1);
for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
if(!tracklets[ipl].IsOK()) continue;
for(Int_t itb = 0; itb < fgNTimeBins; itb++){
if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
if (!tracklets[ipl].IsUsable(itb)) continue;
x = cl->GetX();
y = cl->GetY();
z = cl->GetZ();
tilt = tracklets[ipl].GetTilt();
dx = x - xref;
yr = y + tilt*(z - z0 - dzdx*dx);
ye = tilt*TMath::Sqrt(cl->GetSigmaZ2());
ye += err ? tracklets[ipl].GetSigmaY() : 0.2;
yfitter.AddPoint(&dx, yr, ye);
}
}
yfitter.Eval();
Double_t y0 = yfitter.GetFunctionParameter(0);
Double_t dydx = yfitter.GetFunctionParameter(1);
Double_t chi2 = 0.;
if(np && points){
Float_t xyz[3];
for(int ip=0; ip<np; ip++){
points[ip].GetXYZ(xyz);
xyz[1] = y0 + dydx * (xyz[0] - xref);
xyz[2] = z0 + dzdx * (xyz[0] - xref);
points[ip].SetXYZ(xyz);
}
}
return chi2;
}
Double_t AliTRDtrackerV1::FitRiemanTilt(const AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t sigError, Int_t np, AliTrackPoint *points)
{
// R^{2} = (x-x_{0})^{2} + (y^{*}-y_{0})^{2}
// y^{*} = y - tg(h)(z - z_{t})
// z_{t} = z_{0}+dzdx*(x-x_{r})
// END_LATEX
// t = 1 / (x^{2} + y^{2})
// u = 2 * x * t
// v = 2 * tan(h) * t
// w = 2 * tan(h) * (x - x_{r}) * t
// END_LATEX
// a + b * u + c * t + d * v + e * w = 2 * (y + tg(h) * z) * t
// END_LATEX
// a = -1/y_{0}
// b = x_{0}/y_{0}
// c = (R^{2} -x_{0}^{2} - y_{0}^{2})/y_{0}
// d = z_{0}
// e = dz/dx
// END_LATEX
// #sigma = 2 * #sqrt{#sigma^{2}_{y} + (tilt corr ...) + tg^{2}(h) * #sigma^{2}_{z}} * t
// END_LATEX
// C = 1/R = a/(1 + b^{2} + c*a)
// END_LATEX
TLinearFitter *fitter = GetTiltedRiemanFitter();
fitter->StoreData(kTRUE);
fitter->ClearPoints();
AliTRDLeastSquare zfitter;
AliTRDcluster *cl = NULL;
AliTRDseedV1 work[kNPlanes], *tracklet = NULL;
if(!tracklets){
for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
if(!(tracklet = track->GetTracklet(ipl))) continue;
if(!tracklet->IsOK()) continue;
new(&work[ipl]) AliTRDseedV1(*tracklet);
}
tracklets = &work[0];
}
Double_t xref = CalculateReferenceX(tracklets);
if(AliLog::GetDebugLevel("TRD", "AliTRDtrackerV1")>3) printf("D-AliTRDtrackerV1::FitRiemanTilt:\nx0[(0)%6.2f (1)%6.2f (2)%6.2f (3)%6.2f (4)%6.2f (5)%6.2f] xref[%6.2f]", tracklets[0].GetX0(), tracklets[1].GetX0(), tracklets[2].GetX0(), tracklets[3].GetX0(), tracklets[4].GetX0(), tracklets[5].GetX0(), xref);
Double_t x, y, z, t, tilt, dx, w, we;
Double_t uvt[4];
Int_t nPoints = 0;
for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
if(!tracklets[ipl].IsOK()) continue;
for(Int_t itb = 0; itb < AliTRDseedV1::kNclusters; itb++){
if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
x = cl->GetX();
y = cl->GetY();
z = cl->GetZ();
tilt = tracklets[ipl].GetTilt();
dx = x - xref;
t = 1./(x*x + y*y);
uvt[0] = 2. * x * t;
uvt[1] = t;
uvt[2] = 2. * tilt * t;
uvt[3] = 2. * tilt * dx * t;
w = 2. * (y + tilt*z) * t;
we = 2. * t;
we *= sigError ? TMath::Sqrt(cl->GetSigmaY2()) : 0.2;
fitter->AddPoint(uvt, w, we);
zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
nPoints++;
}
}
if(fitter->Eval()) return 1.E10;
Double_t z0 = fitter->GetParameter(3);
Double_t dzdx = fitter->GetParameter(4);
Bool_t accept = kTRUE;
Double_t zref = 0.0;
for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
if(!tracklets[iLayer].IsOK()) continue;
zref = z0 + dzdx * (tracklets[iLayer].GetX0() - xref);
if (TMath::Abs(tracklets[iLayer].GetZfit(0) - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
accept = kFALSE;
}
if (!accept) {
zfitter.Eval();
Double_t dzmf = zfitter.GetFunctionParameter(1);
Double_t zmf = zfitter.GetFunctionValue(&xref);
fitter->FixParameter(3, zmf);
fitter->FixParameter(4, dzmf);
fitter->Eval();
fitter->ReleaseParameter(3);
fitter->ReleaseParameter(4);
z0 = fitter->GetParameter(3);
dzdx = fitter->GetParameter(4);
}
Double_t a = fitter->GetParameter(0);
Double_t b = fitter->GetParameter(1);
Double_t c = fitter->GetParameter(2);
Double_t y0 = 1. / a;
Double_t x0 = -b * y0;
Double_t tmp = y0*y0 + x0*x0 - c*y0;
if(tmp<=0.) return 1.E10;
Double_t radius = TMath::Sqrt(tmp);
Double_t curvature = 1.0 + b*b - c*a;
if (curvature > 0.0) curvature = a / TMath::Sqrt(curvature);
Double_t chi2 = fitter->GetChisquare()/Double_t(nPoints);
if(AliLog::GetDebugLevel("TRD", "AliTRDtrackerV1")>3) printf("D-AliTRDtrackerV1::FitRiemanTilt:x0[%6.2f] y0[%6.2f] R[%6.2f] chi2[%f]\n", x0, y0, radius, chi2);
if(!track){
for(Int_t ip = 0; ip < kNPlanes; ip++) {
x = tracklets[ip].GetX0();
tmp = radius*radius-(x-x0)*(x-x0);
if(tmp <= 0.) continue;
tmp = TMath::Sqrt(tmp);
tracklets[ip].SetYref(0, y0 - (y0>0.?1.:-1)*tmp);
tracklets[ip].SetYref(1, (x - x0) / tmp);
tracklets[ip].SetZref(0, z0 + dzdx * (x - xref));
tracklets[ip].SetZref(1, dzdx);
tracklets[ip].SetC(curvature);
tracklets[ip].SetChi2(chi2);
}
}
if(np && points){
Float_t xyz[3];
for(int ip=0; ip<np; ip++){
points[ip].GetXYZ(xyz);
xyz[1] = TMath::Abs(xyz[0] - x0) > radius ? 100. : y0 - (y0>0.?1.:-1.)*TMath::Sqrt((radius-(xyz[0]-x0))*(radius+(xyz[0]-x0)));
xyz[2] = z0 + dzdx * (xyz[0] - xref);
points[ip].SetXYZ(xyz);
}
}
return chi2;
}
Double_t AliTRDtrackerV1::FitKalman(AliTRDtrackV1 *track, AliTRDseedV1 * const tracklets, Bool_t up, Int_t np, AliTrackPoint *points)
{
Int_t ip = np ? 0 : 1;
while(ip<np){
if((up?-1:1) * (track->GetX() - points[ip].GetX()) > 0.) break;
ip++;
}
AliTRDseedV1 tracklet;
AliTRDseedV1 *ptrTracklet = NULL;
for (Int_t jplane = 0; jplane < kNPlanes; jplane++) {
Int_t iplane = up ? jplane : kNPlanes - 1 - jplane;
if(tracklets){
if(!(ptrTracklet = &tracklets[iplane])) continue;
}else{
if(!(ptrTracklet = track->GetTracklet(iplane))){
continue;
}
}
if(!ptrTracklet->IsOK()) continue;
Double_t x = ptrTracklet->GetX0();
while(ip < np){
if((up?-1:1) * (points[ip].GetX() - x) - AliTRDReconstructor::GetMaxStep() < 0) break;
if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), AliTRDReconstructor::GetMaxStep())) return -1.;
Double_t xyz[3];
track->GetXYZ(xyz);
track->Global2LocalPosition(xyz, track->GetAlpha());
points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
ip++;
}
if(((up?-1:1) * (x - track->GetX()) + AliTRDReconstructor::GetMaxStep() < 0) && !PropagateToX(*track, x + (up?-1:1)*AliTRDReconstructor::GetMaxStep(), AliTRDReconstructor::GetMaxStep())) return -1.;
if(!AdjustSector(track)) return -1;
if(TMath::Abs(track->GetSnp()) > AliTRDReconstructor::GetMaxSnp()) return -1;
if(!tracklets) track->SetTracklet(ptrTracklet, -1);
Double_t xyz0[3]; track->GetXYZ(xyz0);
Double_t alpha = track->GetAlpha();
Double_t xyz1[3], y, z;
if(!track->GetProlongation(x, y, z)) return -1;
xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
xyz1[2] = z;
if(TMath::Abs(xyz0[0] - xyz1[0]) < 1e-3 && TMath::Abs(xyz0[1] - xyz1[1]) < 1e-3) continue;
Double_t param[7];
if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param) <=0.) break;
Double_t xrho = param[0]*param[4];
Double_t xx0 = param[1];
track->PropagateTo(x, xx0, xrho);
if (!AdjustSector(track)) break;
Double_t cov[3]; ptrTracklet->GetCovAt(x, cov);
Double_t p[2] = { ptrTracklet->GetY(), ptrTracklet->GetZ()};
Double_t chi2 = ((AliExternalTrackParam*)track)->GetPredictedChi2(p, cov);
if(chi2<1e+10) ((AliExternalTrackParam*)track)->Update(p, cov);
if(!up) continue;
if(iplane>0 && track->GetTracklet(iplane-1) && ptrTracklet->GetN() + track->GetTracklet(iplane-1)->GetN() > 20) track->SetBudget(2, 0.);
}
while(ip < np){
if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), AliTRDReconstructor::GetMaxStep())) return -1.;
Double_t xyz[3];
track->GetXYZ(xyz);
track->Global2LocalPosition(xyz, track->GetAlpha());
points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
ip++;
}
return track->GetChi2();
}
Float_t AliTRDtrackerV1::CalculateChi2Z(const AliTRDseedV1 *tracklets, Double_t offset, Double_t slope, Double_t xref)
{
Float_t chi2Z = 0, nLayers = 0;
for (Int_t iLayer = 0; iLayer < AliTRDgeometry::kNlayer; iLayer++) {
if(!tracklets[iLayer].IsOK()) continue;
Double_t z = offset + slope * (tracklets[iLayer].GetX0() - xref);
chi2Z += TMath::Abs(tracklets[iLayer].GetZfit(0) - z);
nLayers++;
}
chi2Z /= TMath::Max((nLayers - 3.0),1.0);
return chi2Z;
}
Int_t AliTRDtrackerV1::PropagateToX(AliTRDtrackV1 &t, Double_t xToGo, Double_t maxStep)
{
Double_t xpos = t.GetX();
Double_t dir = (xpos < xToGo) ? 1.0 : -1.0;
while (((xToGo - xpos) * dir) > AliTRDReconstructor::GetEpsilon()) {
Double_t xyz0[3];
Double_t xyz1[3];
Double_t param[7];
Double_t x;
Double_t y;
Double_t z;
Double_t step = dir * TMath::Min(TMath::Abs(xToGo-xpos),maxStep);
t.GetXYZ(xyz0);
x = xpos + step;
if(t.GetProlongation(x,y,z)<0) return 0;
xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
xyz1[2] = z;
if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param)<=0.) return 0;
if (!t.PropagateTo(x, param[1], param[0]*param[4])) return 0;
if(!AdjustSector(&t)) return 0;
xpos = t.GetX();
}
return 1;
}
Bool_t AliTRDtrackerV1::ReadClusters(TTree *clusterTree)
{
Int_t nsize = Int_t(clusterTree->GetTotBytes() / (sizeof(AliTRDcluster)));
TObjArray *clusterArray = new TObjArray(nsize+1000);
TBranch *branch = clusterTree->GetBranch("TRDcluster");
if (!branch) {
AliError("Can't get the branch !");
return kFALSE;
}
branch->SetAddress(&clusterArray);
if(!fClusters){
Float_t nclusters = fkRecoParam->GetNClusters();
if(fkReconstructor->IsHLT()) nclusters /= AliTRDgeometry::kNsector;
fClusters = new TClonesArray("AliTRDcluster", Int_t(nclusters));
fClusters->SetOwner(kTRUE);
SetClustersOwner();
AliInfo(Form("Tracker owning clusters @ %p", (void*)fClusters));
}
Int_t nEntries = (Int_t) clusterTree->GetEntries();
Int_t nbytes = 0;
Int_t ncl = 0;
AliTRDcluster *c = NULL;
for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
nbytes += clusterTree->GetEvent(iEntry);
Int_t nCluster = clusterArray->GetEntriesFast();
for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
if(!(c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster))) continue;
new((*fClusters)[ncl++]) AliTRDcluster(*c);
delete (clusterArray->RemoveAt(iCluster));
}
}
delete clusterArray;
return kTRUE;
}
Int_t AliTRDtrackerV1::LoadClusters(TTree *cTree)
{
fkRecoParam = fkReconstructor->GetRecoParam();
if(!(fClusters = AliTRDReconstructor::GetClusters())){
AliWarning("Clusters unavailable from TRD reconstructor. Trying reading from tree ...");
} else {
if(!ReadClusters(cTree)) {
AliError("Reading clusters from tree failed.");
return 1;
}
}
if(!fClusters || !fClusters->GetEntriesFast()){
AliInfo("No TRD clusters");
return 1;
} else AliInfo(Form("Using :: clusters[%d] onl.tracklets[%d] onl.tracks[%d]",
fClusters?fClusters->GetEntriesFast():0,
AliTRDReconstructor::GetTracklets()?AliTRDReconstructor::GetTracklets()->GetEntriesFast():0,
AliTRDReconstructor::GetTracks()?AliTRDReconstructor::GetTracks()->GetEntriesFast():0));
BuildTrackingContainers();
return 0;
}
Int_t AliTRDtrackerV1::LoadClusters(TClonesArray * const clusters)
{
if(!clusters || !clusters->GetEntriesFast()){
AliInfo("No TRD clusters");
return 1;
} else AliInfo(Form("Using :: external.clusters[%d]", clusters->GetEntriesFast()));
fClusters = clusters;
fkRecoParam = fkReconstructor->GetRecoParam();
BuildTrackingContainers();
return 0;
}
Int_t AliTRDtrackerV1::BuildTrackingContainers()
{
Int_t nin(0), ncl(fClusters->GetEntriesFast());
while (ncl--) {
AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(ncl);
if(c->IsInChamber()) nin++;
if(fkReconstructor->IsHLT()) c->SetRPhiMethod(AliTRDcluster::kCOG);
Int_t detector = c->GetDetector();
Int_t sector = fGeom->GetSector(detector);
Int_t stack = fGeom->GetStack(detector);
Int_t layer = fGeom->GetLayer(detector);
fTrSec[sector].GetChamber(stack, layer, kTRUE)->InsertCluster(c, ncl);
}
for(int isector =0; isector<AliTRDgeometry::kNsector; isector++){
if(!fTrSec[isector].GetNChambers()) continue;
fTrSec[isector].Init(fkReconstructor);
}
return nin;
}
void AliTRDtrackerV1::UnloadClusters()
{
if(fTracks){
fTracks->Delete();
if(HasRemoveContainers()){delete fTracks; fTracks = NULL;}
}
if(fTracklets){
fTracklets->Delete();
if(HasRemoveContainers()){delete fTracklets; fTracklets = NULL;}
}
if(fClusters && IsClustersOwner()){
AliInfo(Form("tracker[%p] clearing %d own clusters @ %p", (void*)this, fClusters->GetEntries(), (void*)fClusters));
fClusters->Delete();
}
for (int i = 0; i < AliTRDgeometry::kNsector; i++) fTrSec[i].Clear();
AliTRDtrackerDebug::SetEventNumber(AliTRDtrackerDebug::GetEventNumber() + 1);
}
Bool_t AliTRDtrackerV1::AdjustSector(AliTRDtrackV1 *const track)
{
Double_t alpha = AliTRDgeometry::GetAlpha();
Double_t y = track->GetY();
Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
if (y > ymax) {
if (!track->Rotate( alpha)) {
return kFALSE;
}
}
else if (y < -ymax) {
if (!track->Rotate(-alpha)) {
return kFALSE;
}
}
return kTRUE;
}
AliTRDseedV1* AliTRDtrackerV1::GetTracklet(const AliTRDtrackV1 *const track, Int_t p, Int_t &idx)
{
idx = track->GetTrackletIndex(p);
AliTRDseedV1 *tracklet = (idx<0) ? NULL : (AliTRDseedV1*)fTracklets->UncheckedAt(idx);
return tracklet;
}
AliTRDseedV1* AliTRDtrackerV1::SetTracklet(const AliTRDseedV1 * const tracklet)
{
if(!fTracklets){
fTracklets = new TClonesArray("AliTRDseedV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
fTracklets->SetOwner(kTRUE);
}
Int_t nentries = fTracklets->GetEntriesFast();
return new ((*fTracklets)[nentries]) AliTRDseedV1(*tracklet);
}
AliTRDtrackV1* AliTRDtrackerV1::SetTrack(const AliTRDtrackV1 * const track)
{
if(!fTracks){
fTracks = new TClonesArray("AliTRDtrackV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
fTracks->SetOwner(kTRUE);
}
Int_t nentries = fTracks->GetEntriesFast();
return new ((*fTracks)[nentries]) AliTRDtrackV1(*track);
}
Int_t AliTRDtrackerV1::Clusters2TracksSM(Int_t sector, AliESDEvent *esd)
{
Int_t nTracks = 0;
Int_t nChambers = 0;
AliTRDtrackingChamber **stack = NULL, *chamber = NULL;
for(int istack = 0; istack<AliTRDgeometry::kNstack; istack++){
if(!(stack = fTrSec[sector].GetStack(istack))) continue;
nChambers = 0;
for(int ilayer=0; ilayer<AliTRDgeometry::kNlayer; ilayer++){
if(!(chamber = stack[ilayer])) continue;
if(chamber->GetNClusters() < fgNTimeBins * fkRecoParam->GetFindableClusters()) continue;
nChambers++;
}
if(nChambers < 4) continue;
nTracks += Clusters2TracksStack(stack, fTracksESD);
}
if(nTracks) AliDebug(2, Form("Number of tracks: SM_%02d[%d]", sector, nTracks));
for(int itrack=0; itrack<nTracks; itrack++){
AliESDtrack *esdTrack((AliESDtrack*)(fTracksESD->operator[](itrack)));
Int_t id = esd->AddTrack(esdTrack);
if (fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 0 || AliTRDReconstructor::GetStreamLevel()>0 ){
esdTrack=esd->GetTrack(id);
TObject *o(NULL); Int_t ic(0);
AliTRDtrackV1 *calibTrack(NULL);
while((o = esdTrack->GetCalibObject(ic++))){
if(!(calibTrack = dynamic_cast<AliTRDtrackV1*>(o))) continue;
calibTrack->SetESDid(esdTrack->GetID());
break;
}
}
}
AliTRDtrackerDebug::SetCandidateNumber(0);
AliTRDtrackerDebug::SetTrackNumber(0);
fTracksESD->Delete();
return nTracks;
}
Int_t AliTRDtrackerV1::Clusters2TracksStack(AliTRDtrackingChamber **stack, TClonesArray * const esdTrackList)
{
AliTRDtrackingChamber *chamber = NULL;
AliTRDtrackingChamber **ci = NULL;
AliTRDseedV1 sseed[kMaxTracksStack*6];
Int_t pars[4];
Int_t configs[kNConfigs];
ci = &stack[0];
for(Int_t ic = kNPlanes; ic--; ci++){
if(!(*ci)) continue;
(*ci)->Update();
}
Double_t quality = BuildSeedingConfigs(stack, configs);
if(fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 10 || AliTRDReconstructor::GetStreamLevel()>10){
AliInfo(Form("Plane config %d %d %d Quality %f"
, configs[0], configs[1], configs[2], quality));
}
Int_t ntracks,
ntracks1,
ntracks2 = 0;
fSieveSeeding = 0;
Int_t ic = 0; ci = &stack[0];
while(ic<kNPlanes && !(*ci)){ic++; ci++;}
if(!(*ci)) return ntracks2;
Int_t istack = fGeom->GetStack((*ci)->GetDetector());
do{
ntracks = 0; ntracks1 = 0;
for (Int_t iconf = 0; iconf<fkRecoParam->GetNumberOfSeedConfigs(); iconf++) {
pars[0] = configs[iconf];
pars[1] = ntracks;
pars[2] = istack;
ntracks = MakeSeeds(stack, &sseed[6*ntracks], pars);
if(ntracks == kMaxTracksStack) break;
}
AliDebug(2, Form("Candidate TRD tracks %d in iteration %d.", ntracks, fSieveSeeding));
if(!ntracks) break;
Int_t sort[kMaxTracksStack+1];
TMath::Sort(ntracks, fTrackQuality, sort, kTRUE);
if(AliLog::GetDebugLevel("TRD", "AliTRDtrackerV1") > 2){
AliDebug(3, "Track candidates classification:");
for (Int_t it(0); it < ntracks; it++) {
Int_t jt(sort[it]);
printf(" %2d idx[%d] Quality[%e]\n", it, jt, fTrackQuality[jt]);
}
}
Int_t ntracks0 = esdTrackList->GetEntriesFast();
Bool_t signedTrack[kMaxTracksStack];
Bool_t fakeTrack[kMaxTracksStack];
for (Int_t i=0; i<ntracks; i++){
signedTrack[i] = kFALSE;
fakeTrack[i] = kFALSE;
}
Int_t jSieve(0), rejectedCandidates(0);
do{
rejectedCandidates=0;
for (Int_t itrack = 0; itrack < ntracks; itrack++) {
Int_t trackIndex = sort[itrack];
if (signedTrack[trackIndex] || fakeTrack[trackIndex]) continue;
Int_t ncl = 0;
Int_t nused = 0;
Int_t nlayers = 0;
Int_t findable = 0;
for (Int_t jLayer = 0; jLayer < kNPlanes; jLayer++) {
Int_t jseed = kNPlanes*trackIndex+jLayer;
sseed[jseed].UpdateUsed();
if(!sseed[jseed].IsOK()) continue;
if (TMath::Abs(sseed[jseed].GetYref(0) / sseed[jseed].GetX0()) < 0.158) findable++;
ncl += sseed[jseed].GetN();
nused += sseed[jseed].GetNUsed();
nlayers++;
}
if (nused > 30){
AliDebug(4, Form("REJECTED : %d idx[%d] quality[%e] tracklets[%d] usedClusters[%d]", itrack, trackIndex, fTrackQuality[trackIndex], nlayers, nused));
fakeTrack[trackIndex] = kTRUE;
continue;
}
if (ncl>0 && Float_t(nused)/ncl >= .25){
AliDebug(4, Form("REJECTED : %d idx[%d] quality[%e] tracklets[%d] usedClusters[%d] used/ncl[%f]", itrack, trackIndex, fTrackQuality[trackIndex], nlayers, nused, Float_t(nused)/ncl));
fakeTrack[trackIndex] = kTRUE;
continue;
}
AliDebug(4, Form("Candidate[%d] Quality[%e] Tracklets[%d] Findable[%d] Ncl[%d] Nused[%d]", trackIndex, fTrackQuality[trackIndex], nlayers, findable, ncl, nused));
Bool_t skip = kFALSE;
switch(jSieve){
case 0:
if(nlayers > findable || nlayers < kNPlanes) {skip = kTRUE; break;}
if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
break;
case 1:
if(nlayers < findable){skip = kTRUE; break;}
if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -4.){skip = kTRUE; break;}
break;
case 2:
if(nlayers < kNPlanes) { skip = kTRUE; break;}
if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -6.0){skip = kTRUE; break;}
break;
case 3:
if (nlayers<4){skip = kTRUE; break;}
if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
break;
case 4:
if (nlayers<4){skip = kTRUE; break;}
break;
}
if(skip){
rejectedCandidates++;
AliDebug(4, Form("REJECTED : %d idx[%d] quality[%e] tracklets[%d] usedClusters[%d]", itrack, trackIndex, fTrackQuality[trackIndex], nlayers, nused));
continue;
} else AliDebug(4, Form("ACCEPTED : %d idx[%d] quality[%e] tracklets[%d] usedClusters[%d]", itrack, trackIndex, fTrackQuality[trackIndex], nlayers, nused));
signedTrack[trackIndex] = kTRUE;
AliTRDseedV1 *lseed =&sseed[trackIndex*kNPlanes];
AliTRDtrackV1 *track = MakeTrack(lseed);
if(!track){
AliDebug(1, "Track building failed.");
continue;
} else {
if(AliLog::GetDebugLevel("TRD", "AliTRDtrackerV1") > 1){
Int_t ich = 0; while(!(chamber = stack[ich])) ich++;
AliDebug(2, Form("Track pt=%7.2fGeV/c SM[%2d] Done.", track->Pt(), fGeom->GetSector(chamber->GetDetector())));
}
}
if(fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 1 && fkReconstructor->IsDebugStreaming()){
AliTRDseedV1 *dseed[6];
for(Int_t iseed = AliTRDgeometry::kNlayer; iseed--;) dseed[iseed] = new AliTRDseedV1(lseed[iseed]);
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t trackNumber = AliTRDtrackerDebug::GetTrackNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
cstreamer << "Clusters2TracksStack"
<< "EventNumber=" << eventNumber
<< "TrackNumber=" << trackNumber
<< "CandidateNumber=" << candidateNumber
<< "Iter=" << fSieveSeeding
<< "Like=" << fTrackQuality[trackIndex]
<< "S0.=" << dseed[0]
<< "S1.=" << dseed[1]
<< "S2.=" << dseed[2]
<< "S3.=" << dseed[3]
<< "S4.=" << dseed[4]
<< "S5.=" << dseed[5]
<< "Ncl=" << ncl
<< "NLayers=" << nlayers
<< "Findable=" << findable
<< "NUsed=" << nused
<< "\n";
}
AliESDtrack *esdTrack = new ((*esdTrackList)[ntracks0++]) AliESDtrack();
esdTrack->UpdateTrackParams(track, AliESDtrack::kTRDout);
esdTrack->SetLabel(track->GetLabel());
track->UpdateESDtrack(esdTrack);
if (fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 0 || AliTRDReconstructor::GetStreamLevel()>0 ){
AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(*track);
calibTrack->SetOwner();
esdTrack->AddCalibObject(calibTrack);
}
ntracks1++;
AliTRDtrackerDebug::SetTrackNumber(AliTRDtrackerDebug::GetTrackNumber() + 1);
}
jSieve++;
} while(jSieve<5 && rejectedCandidates);
if(!ntracks1) break;
ntracks2 += ntracks1;
if(fkReconstructor->IsHLT()) break;
fSieveSeeding++;
quality = BuildSeedingConfigs(stack, configs);
if(quality < 1.E-7) break;
for(Int_t ip = 0; ip < kNPlanes; ip++){
if(!(chamber = stack[ip])) continue;
chamber->Build(fGeom);
}
if(fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 10 || AliTRDReconstructor::GetStreamLevel()>10){
AliInfo(Form("Sieve level %d Plane config %d %d %d Quality %f", fSieveSeeding, configs[0], configs[1], configs[2], quality));
}
} while(fSieveSeeding<10);
return ntracks2;
}
Double_t AliTRDtrackerV1::BuildSeedingConfigs(AliTRDtrackingChamber **stack, Int_t *configs)
{
Double_t chamberQ[kNPlanes];memset(chamberQ, 0, kNPlanes*sizeof(Double_t));
AliTRDtrackingChamber *chamber = NULL;
for(int iplane=0; iplane<kNPlanes; iplane++){
if(!(chamber = stack[iplane])) continue;
chamberQ[iplane] = (chamber = stack[iplane]) ? chamber->GetQuality() : 0.;
}
Double_t tconfig[kNConfigs];memset(tconfig, 0, kNConfigs*sizeof(Double_t));
Int_t planes[] = {0, 0, 0, 0};
for(int iconf=0; iconf<kNConfigs; iconf++){
GetSeedingConfig(iconf, planes);
tconfig[iconf] = fgTopologicQA[iconf];
for(int iplane=0; iplane<4; iplane++) tconfig[iconf] *= chamberQ[planes[iplane]];
}
TMath::Sort((Int_t)kNConfigs, tconfig, configs, kTRUE);
return tconfig[configs[0]];
}
Int_t AliTRDtrackerV1::MakeSeeds(AliTRDtrackingChamber **stack, AliTRDseedV1 * const sseed, const Int_t * const ipar)
{
AliTRDtrackingChamber *chamber = NULL;
AliTRDcluster *c[kNSeedPlanes] = {NULL, NULL, NULL, NULL};
AliTRDseedV1 *cseed = &sseed[0];
Int_t ncl, mcl;
Int_t index[AliTRDchamberTimeBin::kMaxClustersLayer], jndex[AliTRDchamberTimeBin::kMaxClustersLayer];
Double_t chi2[4];
Int_t config = ipar[0];
Int_t ntracks = ipar[1];
Int_t istack = ipar[2];
Int_t planes[kNSeedPlanes]; GetSeedingConfig(config, planes);
Int_t planesExt[kNPlanes-kNSeedPlanes]; GetExtrapolationConfig(config, planesExt);
Double_t hL[kNPlanes];
Float_t padlength[kNPlanes];
Float_t padwidth[kNPlanes];
AliTRDpadPlane *pp = NULL;
for(int iplane=0; iplane<kNPlanes; iplane++){
pp = fGeom->GetPadPlane(iplane, istack);
hL[iplane] = TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle());
padlength[iplane] = pp->GetLengthIPad();
padwidth[iplane] = pp->GetWidthIPad();
}
Double_t x0[kNPlanes],
driftLength = .5*AliTRDgeometry::AmThick() - AliTRDgeometry::DrThick();
TGeoHMatrix *matrix = NULL;
Double_t loc[] = {AliTRDgeometry::AnodePos(), 0., 0.};
Double_t glb[] = {0., 0., 0.};
AliTRDtrackingChamber **cIter = &stack[0];
for(int iLayer=0; iLayer<kNPlanes; iLayer++,cIter++){
if(!(*cIter)) continue;
if(!(matrix = fGeom->GetClusterMatrix((*cIter)->GetDetector()))){
x0[iLayer] = fgkX0[iLayer];
continue;
}
matrix->LocalToMaster(loc, glb);
x0[iLayer] = glb[0];
}
AliDebug(2, Form("Making seeds Stack[%d] Config[%d] Tracks[%d]...", istack, config, ntracks));
ResetSeedTB();
Int_t nlayers = 0;
for(int isl=0; isl<kNSeedPlanes; isl++){
if(!(chamber = stack[planes[isl]])) continue;
if(!chamber->GetSeedingLayer(fSeedTB[isl], fGeom, fkReconstructor)) continue;
nlayers++;
}
if(nlayers < kNSeedPlanes) return ntracks;
Double_t cond0[4], cond1[4], cond2[4];
Int_t icl = 0;
while((c[3] = (*fSeedTB[3])[icl++])){
if(!c[3]) continue;
fSeedTB[0]->BuildCond(c[3], cond0, 0);
fSeedTB[0]->GetClusters(cond0, index, ncl);
Int_t jcl = 0;
while(jcl<ncl) {
c[0] = (*fSeedTB[0])[index[jcl++]];
if(!c[0]) continue;
Double_t dx = c[3]->GetX() - c[0]->GetX();
Double_t dzdx = (c[3]->GetZ() - c[0]->GetZ())/dx;
Double_t dydx = (c[3]->GetY() - c[0]->GetY())/dx;
fSeedTB[1]->BuildCond(c[0], cond1, 1, dzdx, dydx);
fSeedTB[1]->GetClusters(cond1, jndex, mcl);
Int_t kcl = 0;
while(kcl<mcl) {
c[1] = (*fSeedTB[1])[jndex[kcl++]];
if(!c[1]) continue;
fSeedTB[2]->BuildCond(c[1], cond2, 2, dzdx, dydx);
c[2] = fSeedTB[2]->GetNearestCluster(cond2);
if(!c[2]) continue;
AliDebug(3, Form("Seeding clusters\n 0[%6.3f %6.3f %6.3f]\n 1[%6.3f %6.3f %6.3f]\n 2[%6.3f %6.3f %6.3f]\n 3[%6.3f %6.3f %6.3f].",
c[0]->GetX(), c[0]->GetY(), c[0]->GetZ(),
c[1]->GetX(), c[1]->GetY(), c[1]->GetZ(),
c[2]->GetX(), c[2]->GetY(), c[2]->GetZ(),
c[3]->GetX(), c[3]->GetY(), c[3]->GetZ()));
for (Int_t il = 0; il < kNPlanes; il++) cseed[il].Reset();
FitRieman(c, chi2);
AliTRDseedV1 *tseed = &cseed[0];
cIter = &stack[0];
for(int iLayer=0; iLayer<kNPlanes; iLayer++, tseed++, cIter++){
Int_t det = (*cIter) ? (*cIter)->GetDetector() : -1;
tseed->SetDetector(det);
tseed->SetTilt(hL[iLayer]);
tseed->SetPadLength(padlength[iLayer]);
tseed->SetPadWidth(padwidth[iLayer]);
tseed->SetReconstructor(fkReconstructor);
tseed->SetX0(det<0 ? fR[iLayer]+driftLength : x0[iLayer]);
tseed->Init(GetRiemanFitter());
tseed->SetStandAlone(kTRUE);
}
Bool_t isFake = kFALSE;
if((fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) >= 2 && fkReconstructor->IsDebugStreaming())
||AliTRDReconstructor::GetStreamLevel()>=2 ){
if (c[0]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
if (c[1]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
if (c[2]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
Double_t xpos[4];
for(Int_t l = 0; l < kNSeedPlanes; l++) xpos[l] = fSeedTB[l]->GetX();
Float_t yref[4];
for(int il=0; il<4; il++) yref[il] = cseed[planes[il]].GetYref(0);
Int_t ll = c[3]->GetLabel(0);
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
AliRieman *rim = GetRiemanFitter();
TTreeSRedirector &cs0 = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
cs0 << "MakeSeeds0"
<<"EventNumber=" << eventNumber
<<"CandidateNumber=" << candidateNumber
<<"isFake=" << isFake
<<"config=" << config
<<"label=" << ll
<<"chi2z=" << chi2[0]
<<"chi2y=" << chi2[1]
<<"Y2exp=" << cond2[0]
<<"Z2exp=" << cond2[1]
<<"X0=" << xpos[0]
<<"X1=" << xpos[1]
<<"X2=" << xpos[2]
<<"X3=" << xpos[3]
<<"yref0=" << yref[0]
<<"yref1=" << yref[1]
<<"yref2=" << yref[2]
<<"yref3=" << yref[3]
<<"c0.=" << c[0]
<<"c1.=" << c[1]
<<"c2.=" << c[2]
<<"c3.=" << c[3]
<<"Seed0.=" << &cseed[planes[0]]
<<"Seed1.=" << &cseed[planes[1]]
<<"Seed2.=" << &cseed[planes[2]]
<<"Seed3.=" << &cseed[planes[3]]
<<"RiemanFitter.=" << rim
<<"\n";
}
if(chi2[0] > fkRecoParam->GetChi2Z()/*iter*/){
AliDebug(3, Form("Filter on chi2Z [%f].", chi2[0]));
AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
continue;
}
if(chi2[1] > fkRecoParam->GetChi2Y()/*iter*/){
AliDebug(3, Form("Filter on chi2Y [%f].", chi2[1]));
AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
continue;
}
Int_t mlayers = 0;
AliTRDcluster *cl = NULL;
for(int iLayer=0; iLayer<kNSeedPlanes; iLayer++){
Int_t jLayer = planes[iLayer];
Int_t nNotInChamber = 0;
if(!cseed[jLayer].AttachClusters(stack[jLayer], kTRUE)) continue;
if(kFALSE){
cseed[jLayer].UpdateUsed();
if(!cseed[jLayer].IsOK()) continue;
}else{
cseed[jLayer].Fit();
cseed[jLayer].UpdateUsed();
cseed[jLayer].ResetClusterIter();
while((cl = cseed[jLayer].NextCluster())){
if(!cl->IsInChamber()) nNotInChamber++;
}
if(cseed[jLayer].GetN() - (cseed[jLayer].GetNUsed() + nNotInChamber) < 5) continue;
}
mlayers++;
}
if(mlayers < kNSeedPlanes){
AliDebug(2, Form("Found only %d tracklets out of %d. Skip.", mlayers, kNSeedPlanes));
AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
continue;
}
if(fkReconstructor->IsHLT()){
for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
Int_t jLayer = planesExt[iLayer];
if(!(chamber = stack[jLayer])) continue;
if(!cseed[jLayer].AttachClusters(chamber, kTRUE)) continue;
cseed[jLayer].Fit();
}
fTrackQuality[ntracks] = 1.;
ntracks++;
if(ntracks == kMaxTracksStack) return ntracks;
cseed += 6;
continue;
}
Double_t chi2Vals[4];
chi2Vals[0] = FitTiltedRieman(&cseed[0], kTRUE);
for(int iLayer=0; iLayer<kNSeedPlanes; iLayer++){
Int_t jLayer = planes[iLayer];
cseed[jLayer].Fit(1);
}
Double_t like = CookLikelihood(&cseed[0], planes);
if (TMath::Log(1.E-9 + like) < fkRecoParam->GetTrackLikelihood()){
AliDebug(3, Form("Filter on likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
continue;
}
fSeedLayer[ntracks] = config;
Int_t elayers(0);
for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
Int_t jLayer = planesExt[iLayer];
if(!(chamber = stack[jLayer])) continue;
if ((jLayer == 0) && !(cseed[1].IsOK())) continue;
if ((jLayer == 5) && !(cseed[4].IsOK())) continue;
AliTRDseedV1 pseed = cseed[jLayer];
if(!pseed.AttachClusters(chamber, kTRUE)) continue;
pseed.Fit(1);
cseed[jLayer] = pseed;
chi2Vals[0] = FitTiltedRieman(cseed, kTRUE);
cseed[jLayer].Fit(1);
elayers++;
}
if((fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) >= 2 && fkReconstructor->IsDebugStreaming())
||AliTRDReconstructor::GetStreamLevel()>=2){
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
cstreamer << "MakeSeeds1"
<< "EventNumber=" << eventNumber
<< "CandidateNumber=" << candidateNumber
<< "S0.=" << &cseed[0]
<< "S1.=" << &cseed[1]
<< "S2.=" << &cseed[2]
<< "S3.=" << &cseed[3]
<< "S4.=" << &cseed[4]
<< "S5.=" << &cseed[5]
<< "FitterT.=" << tiltedRieman
<< "\n";
}
if(fkRecoParam->HasImproveTracklets()){
if(!ImproveSeedQuality(stack, cseed, chi2Vals[0])){
AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
AliDebug(3, "ImproveSeedQuality() failed.");
}
}
if(fkRecoParam->IsVertexConstrained()) chi2Vals[1] = FitTiltedRiemanConstraint(&cseed[0], GetZ());
else chi2Vals[1] = -1.;
chi2Vals[2] = GetChi2Z(&cseed[0]);
chi2Vals[3] = GetChi2Phi(&cseed[0]);
fTrackQuality[ntracks] = CalculateTrackLikelihood(&chi2Vals[0]);
if((fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) >= 2 && fkReconstructor->IsDebugStreaming())
||AliTRDReconstructor::GetStreamLevel()>=2){
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
TLinearFitter *fitterTC = GetTiltedRiemanFitterConstraint();
TLinearFitter *fitterT = GetTiltedRiemanFitter();
Int_t ncls = 0;
for(Int_t iseed = 0; iseed < kNPlanes; iseed++){
ncls += cseed[iseed].IsOK() ? cseed[iseed].GetN2() : 0;
}
cstreamer << "MakeSeeds2"
<< "EventNumber=" << eventNumber
<< "CandidateNumber=" << candidateNumber
<< "Chi2TR=" << chi2Vals[0]
<< "Chi2TC=" << chi2Vals[1]
<< "Nlayers=" << mlayers
<< "NClusters=" << ncls
<< "Like=" << like
<< "S0.=" << &cseed[0]
<< "S1.=" << &cseed[1]
<< "S2.=" << &cseed[2]
<< "S3.=" << &cseed[3]
<< "S4.=" << &cseed[4]
<< "S5.=" << &cseed[5]
<< "FitterT.=" << fitterT
<< "FitterTC.=" << fitterTC
<< "\n";
}
if(AliLog::GetDebugLevel("TRD", "AliTRDtrackerV1")){
Double_t pt[]={0., 0.};
for(Int_t il(0); il<kNPlanes; il++){
if(!cseed[il].IsOK()) continue;
pt[0] = GetBz()*kB2C/cseed[il].GetC();
pt[1] = GetBz()*kB2C/cseed[il].GetC(1);
break;
}
AliDebug(2, Form("Candidate[%2d] pt[%7.3f %7.3f] Q[%e]\n"
" [0] x[%6.2f] n[%2d] nu[%d] OK[%c]\n"
" [1] x[%6.2f] n[%2d] nu[%d] OK[%c]\n"
" [2] x[%6.2f] n[%2d] nu[%d] OK[%c]\n"
" [3] x[%6.2f] n[%2d] nu[%d] OK[%c]\n"
" [4] x[%6.2f] n[%2d] nu[%d] OK[%c]\n"
" [5] x[%6.2f] n[%2d] nu[%d] OK[%c]"
, ntracks, pt[0], pt[1], fTrackQuality[ntracks]
,cseed[0].GetX(), cseed[0].GetN(), cseed[0].GetNUsed(), cseed[0].IsOK()?'y':'n'
,cseed[1].GetX(), cseed[1].GetN(), cseed[1].GetNUsed(), cseed[1].IsOK()?'y':'n'
,cseed[2].GetX(), cseed[2].GetN(), cseed[2].GetNUsed(), cseed[2].IsOK()?'y':'n'
,cseed[3].GetX(), cseed[3].GetN(), cseed[3].GetNUsed(), cseed[3].IsOK()?'y':'n'
,cseed[4].GetX(), cseed[4].GetN(), cseed[4].GetNUsed(), cseed[4].IsOK()?'y':'n'
,cseed[5].GetX(), cseed[5].GetN(), cseed[5].GetNUsed(), cseed[5].IsOK()?'y':'n'));
}
ntracks++;
AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
if(ntracks == kMaxTracksStack){
AliWarning(Form("Number of seeds reached maximum allowed (%d) in stack.", kMaxTracksStack));
return ntracks;
}
cseed += 6;
}
}
}
return ntracks;
}
AliTRDtrackV1* AliTRDtrackerV1::MakeTrack(AliTRDseedV1 * const tracklet)
{
if(fkReconstructor->IsHLT()) FitTiltedRiemanConstraint(tracklet, 0);
Double_t alpha = AliTRDgeometry::GetAlpha();
Double_t shift = AliTRDgeometry::GetAlpha()/2.0;
Int_t idx(0); while(idx<kNPlanes && !tracklet[idx].IsOK()) idx++;
if(idx>2){ AliDebug(1, Form("Found suspect track start @ layer idx[%d]\n"
" %c[0] x0[%f] n[%d] nu[%d] OK[%c]\n"
" %c[1] x0[%f] n[%d] nu[%d] OK[%c]\n"
" %c[2] x0[%f] n[%d] nu[%d] OK[%c]\n"
" %c[3] x0[%f] n[%d] nu[%d] OK[%c]\n"
" %c[4] x0[%f] n[%d] nu[%d] OK[%c]\n"
" %c[5] x0[%f] n[%d] nu[%d] OK[%c]"
,idx
,idx==0?'*':' ', tracklet[0].GetX0(), tracklet[0].GetN(), tracklet[0].GetNUsed(), tracklet[0].IsOK()?'y':'n'
,idx==1?'*':' ', tracklet[1].GetX0(), tracklet[1].GetN(), tracklet[1].GetNUsed(), tracklet[1].IsOK()?'y':'n'
,idx==2?'*':' ', tracklet[2].GetX0(), tracklet[2].GetN(), tracklet[2].GetNUsed(), tracklet[2].IsOK()?'y':'n'
,idx==3?'*':' ', tracklet[3].GetX0(), tracklet[3].GetN(), tracklet[3].GetNUsed(), tracklet[3].IsOK()?'y':'n'
,idx==4?'*':' ', tracklet[4].GetX0(), tracklet[4].GetN(), tracklet[4].GetNUsed(), tracklet[4].IsOK()?'y':'n'
,idx==5?'*':' ', tracklet[5].GetX0(), tracklet[5].GetN(), tracklet[5].GetNUsed(), tracklet[5].IsOK()?'y':'n'));
return NULL;
}
Double_t dx(5.);
Double_t x(tracklet[idx].GetX0() - dx);
Double_t params[] = {
tracklet[idx].GetYref(0) - dx*tracklet[idx].GetYref(1)
,tracklet[idx].GetZref(0) - dx*tracklet[idx].GetZref(1)
,TMath::Sin(TMath::ATan(tracklet[idx].GetYref(1)))
,tracklet[idx].GetZref(1) / TMath::Sqrt(1. + tracklet[idx].GetYref(1) * tracklet[idx].GetYref(1))
,tracklet[idx].GetC(fkReconstructor->IsHLT()?1:0)
};
Int_t sector(fGeom->GetSector(tracklet[idx].GetDetector()));
Double_t c[15];
c[ 0] = 0.2;
c[ 1] = 0.0; c[ 2] = 2.0;
c[ 3] = 0.0; c[ 4] = 0.0; c[ 5] = 0.02;
c[ 6] = 0.0; c[ 7] = 0.0; c[ 8] = 0.0; c[ 9] = 0.1;
c[10] = 0.0; c[11] = 0.0; c[12] = 0.0; c[13] = 0.0; c[14] = params[4]*params[4]*0.01;
AliTRDtrackV1 track(tracklet, params, c, x, sector*alpha+shift);
AliTRDseedV1 *ptrTracklet = NULL;
if(kFALSE){
for (Int_t jLayer = 0; jLayer < AliTRDgeometry::kNlayer; jLayer++) {
track.UnsetTracklet(jLayer);
ptrTracklet = &tracklet[jLayer];
if(!ptrTracklet->IsOK()) continue;
if(TMath::Abs(ptrTracklet->GetYref(1) - ptrTracklet->GetYfit(1)) >= .2) continue;
ptrTracklet = SetTracklet(ptrTracklet);
ptrTracklet->UseClusters();
track.SetTracklet(ptrTracklet, fTracklets->GetEntriesFast()-1);
}
AliTRDtrackV1 *ptrTrack = SetTrack(&track);
ptrTrack->CookPID();
ptrTrack->CookLabel(.9);
ptrTrack->SetReconstructor(fkReconstructor);
return ptrTrack;
}
if(TMath::Abs(track.GetX()) + TMath::Abs(track.GetY()) + TMath::Abs(track.GetZ()) > 10000) return NULL;
track.ResetCovariance(1);
Int_t nc = TMath::Abs(FollowBackProlongation(track));
if((fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) > 5 && fkReconstructor->IsDebugStreaming())
||AliTRDReconstructor::GetStreamLevel()>5){
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
Double_t p[5];
track.GetExternalParameters(x, p);
TTreeSRedirector &cs = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
cs << "MakeTrack"
<< "EventNumber=" << eventNumber
<< "CandidateNumber=" << candidateNumber
<< "nc=" << nc
<< "X=" << x
<< "Y=" << p[0]
<< "Z=" << p[1]
<< "snp=" << p[2]
<< "tnd=" << p[3]
<< "crv=" << p[4]
<< "Yin=" << params[0]
<< "Zin=" << params[1]
<< "snpin=" << params[2]
<< "tndin=" << params[3]
<< "crvin=" << params[4]
<< "track.=" << &track
<< "\n";
}
if (nc < 30){
UnsetTrackletsTrack(&track);
return NULL;
}
AliTRDtrackV1 *ptrTrack = SetTrack(&track);
ptrTrack->SetReconstructor(fkReconstructor);
ptrTrack->CookLabel(.9);
for(Int_t il(kNPlanes); il--;){
if(!(ptrTracklet = ptrTrack->GetTracklet(il))) continue;
ptrTracklet->UseClusters();
}
ptrTrack->CookPID();
AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
if(!calibra){
AliInfo("Could not get Calibra instance.");
} else if(calibra->GetHisto2d()){
calibra->UpdateHistogramsV1(ptrTrack);
}
return ptrTrack;
}
Bool_t AliTRDtrackerV1::ImproveSeedQuality(AliTRDtrackingChamber **stack, AliTRDseedV1 *cseed, Double_t &chi2)
{
AliTRDtrackingChamber *chamber = NULL;
AliTRDseedV1 bseed[AliTRDgeometry::kNlayer];
Float_t quality(1.e3),
lQuality[AliTRDgeometry::kNlayer] = {1.e3, 1.e3, 1.e3, 1.e3, 1.e3, 1.e3};
Int_t rLayers(0);
for(Int_t jLayer=AliTRDgeometry::kNlayer; jLayer--;){
bseed[jLayer] = cseed[jLayer];
if(!bseed[jLayer].IsOK()) continue;
rLayers++;
lQuality[jLayer] = bseed[jLayer].GetQuality(kTRUE);
quality += lQuality[jLayer];
}
if (rLayers > 0) {
quality /= rLayers;
}
AliDebug(2, Form("Start N[%d] Q[%f] chi2[%f]", rLayers, quality, chi2));
for (Int_t iter = 0; iter < 4; iter++) {
Int_t nLayers(0); Float_t qualitynew(0.);
Int_t indexes[4*AliTRDgeometry::kNlayer];
TMath::Sort(Int_t(AliTRDgeometry::kNlayer), lQuality, indexes, kFALSE);
for(Int_t jLayer=AliTRDgeometry::kNlayer; jLayer--;) {
Int_t bLayer = indexes[jLayer];
bseed[bLayer].Reset("c");
if(!(chamber = stack[bLayer])) continue;
if(!bseed[bLayer].AttachClusters(chamber, kTRUE)) continue;
bseed[bLayer].Fit(1);
if(!bseed[bLayer].IsOK()) continue;
nLayers++;
lQuality[jLayer] = bseed[jLayer].GetQuality(kTRUE);
qualitynew += lQuality[jLayer];
}
if(rLayers > nLayers){
AliDebug(1, Form("Lost %d tracklets while improving.", rLayers-nLayers));
return iter>0?kTRUE:kFALSE;
} else rLayers=nLayers;
qualitynew /= rLayers;
if(qualitynew > quality){
AliDebug(4, Form("Quality[%f] worsen in iter[%d] to ref[%f].", qualitynew, iter, quality));
return iter>0?kTRUE:kFALSE;
} else quality = qualitynew;
Float_t chi2new = FitTiltedRieman(bseed, kTRUE);
if(chi2new > chi2){
AliDebug(4, Form("Chi2[%f] worsen in iter[%d] to ref[%f].", chi2new, iter, chi2));
return iter>0?kTRUE:kFALSE;
} else chi2 = chi2new;
for(Int_t jLayer=AliTRDgeometry::kNlayer; jLayer--;) cseed[jLayer]=bseed[jLayer];
AliDebug(2, Form("Iter[%d] Q[%f] chi2[%f]", iter, quality, chi2));
if((fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) >= 7 && fkReconstructor->IsDebugStreaming())
||AliTRDReconstructor::GetStreamLevel()>=7){
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
cstreamer << "ImproveSeedQuality"
<< "EventNumber=" << eventNumber
<< "CandidateNumber=" << candidateNumber
<< "Iteration=" << iter
<< "S0.=" << &cseed[0]
<< "S1.=" << &cseed[1]
<< "S2.=" << &cseed[2]
<< "S3.=" << &cseed[3]
<< "S4.=" << &cseed[4]
<< "S5.=" << &cseed[5]
<< "FitterT.=" << tiltedRieman
<< "\n";
}
}
return kTRUE;
}
Double_t AliTRDtrackerV1::CalculateTrackLikelihood(Double_t *chi2){
Double_t likeChi2TR = TMath::Exp(-chi2[0] * 0.0078);
Double_t likeChi2TC(1.);
if(chi2[1]>0.){
likeChi2TC = TMath::Exp(-chi2[1] * 0.677);
Double_t r = likeChi2TC/likeChi2TR;
if(r>1.e2){;}
else if(r<1.e2)
likeChi2TC =1.;
else{;}
}
Double_t likeChi2Z = TMath::Exp(-chi2[2] * 0.14);
Double_t likeChi2Phi= TMath::Exp(-chi2[3] * 3.23);
Double_t trackLikelihood = likeChi2Z * likeChi2TR * likeChi2TC * likeChi2Phi;
AliDebug(2, Form("Likelihood [%e]\n"
" Rieman : chi2[%f] likelihood[%6.2e]\n"
" Vertex : chi2[%f] likelihood[%6.2e]\n"
" Z : chi2[%f] likelihood[%6.2e]\n"
" Phi : chi2[%f] likelihood[%6.2e]"
, trackLikelihood
, chi2[0], likeChi2TR
, chi2[1], likeChi2TC
, chi2[2], likeChi2Z
, chi2[3], likeChi2Phi
));
if((fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) >= 2 && fkReconstructor->IsDebugStreaming())
||AliTRDReconstructor::GetStreamLevel()>=2){
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
cstreamer << "CalculateTrackLikelihood0"
<< "EventNumber=" << eventNumber
<< "CandidateNumber=" << candidateNumber
<< "LikeChi2Z=" << likeChi2Z
<< "LikeChi2TR=" << likeChi2TR
<< "LikeChi2TC=" << likeChi2TC
<< "LikeChi2Phi=" << likeChi2Phi
<< "TrackLikelihood=" << trackLikelihood
<< "\n";
}
return trackLikelihood;
}
Double_t AliTRDtrackerV1::CookLikelihood(AliTRDseedV1 *cseed, Int_t planes[4])
{
Double_t chi2y = GetChi2Y(&cseed[0]);
Double_t chi2z = GetChi2Z(&cseed[0]);
Float_t nclusters = 0.;
Double_t sumda = 0.;
for(UChar_t ilayer = 0; ilayer < 4; ilayer++){
Int_t jlayer = planes[ilayer];
nclusters += cseed[jlayer].GetN2();
sumda += TMath::Abs(cseed[jlayer].GetYfit(1) - cseed[jlayer].GetYref(1));
}
nclusters *= .25;
Double_t likea = TMath::Exp(-sumda * fkRecoParam->GetPhiSlope());
Double_t likechi2y = 0.0000000001;
if (fkReconstructor->IsCosmic() || chi2y < fkRecoParam->GetChi2YCut()) likechi2y += TMath::Exp(-TMath::Sqrt(chi2y) * fkRecoParam->GetChi2YSlope());
Double_t likechi2z = TMath::Exp(-chi2z * fkRecoParam->GetChi2ZSlope());
Double_t likeN = TMath::Exp(-(fkRecoParam->GetNMeanClusters() - nclusters) / fkRecoParam->GetNSigmaClusters());
Double_t like = likea * likechi2y * likechi2z * likeN;
if((fkRecoParam->GetStreamLevel(AliTRDrecoParam::kTracker) >= 2 && fkReconstructor->IsDebugStreaming())
||AliTRDReconstructor::GetStreamLevel()>=2){
Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
Int_t nTracklets = 0; Float_t meanNcls = 0;
for(Int_t iseed=0; iseed < kNPlanes; iseed++){
if(!cseed[iseed].IsOK()) continue;
nTracklets++;
meanNcls += cseed[iseed].GetN2();
}
if(nTracklets) meanNcls /= nTracklets;
TTreeSRedirector &cstreamer = *fkReconstructor->GetDebugStream(AliTRDrecoParam::kTracker);
cstreamer << "CookLikelihood"
<< "EventNumber=" << eventNumber
<< "CandidateNumber=" << candidateNumber
<< "tracklet0.=" << &cseed[0]
<< "tracklet1.=" << &cseed[1]
<< "tracklet2.=" << &cseed[2]
<< "tracklet3.=" << &cseed[3]
<< "tracklet4.=" << &cseed[4]
<< "tracklet5.=" << &cseed[5]
<< "sumda=" << sumda
<< "chi2y=" << chi2y
<< "chi2z=" << chi2z
<< "likea=" << likea
<< "likechi2y=" << likechi2y
<< "likechi2z=" << likechi2z
<< "nclusters=" << nclusters
<< "likeN=" << likeN
<< "like=" << like
<< "meanncls=" << meanNcls
<< "\n";
}
return like;
}
void AliTRDtrackerV1::GetSeedingConfig(Int_t iconfig, Int_t planes[4])
{
// <img src="gif/topologicQA.gif">
//End_Html
switch(iconfig){
case 0:
planes[0] = 2;
planes[1] = 3;
planes[2] = 4;
planes[3] = 5;
break;
case 1:
planes[0] = 1;
planes[1] = 2;
planes[2] = 3;
planes[3] = 4;
break;
case 2:
planes[0] = 0;
planes[1] = 1;
planes[2] = 2;
planes[3] = 3;
break;
case 3:
planes[0] = 1;
planes[1] = 2;
planes[2] = 3;
planes[3] = 5;
break;
case 4:
planes[0] = 0;
planes[1] = 1;
planes[2] = 2;
planes[3] = 4;
break;
case 5:
planes[0] = 1;
planes[1] = 3;
planes[2] = 4;
planes[3] = 5;
break;
case 6:
planes[0] = 0;
planes[1] = 2;
planes[2] = 3;
planes[3] = 4;
break;
case 7:
planes[0] = 0;
planes[1] = 3;
planes[2] = 4;
planes[3] = 5;
break;
case 8:
planes[0] = 0;
planes[1] = 1;
planes[2] = 2;
planes[3] = 5;
break;
case 9:
planes[0] = 1;
planes[1] = 2;
planes[2] = 4;
planes[3] = 5;
break;
case 10:
planes[0] = 0;
planes[1] = 1;
planes[2] = 3;
planes[3] = 4;
break;
case 11:
planes[0] = 0;
planes[1] = 1;
planes[2] = 4;
planes[3] = 5;
break;
case 12:
planes[0] = 0;
planes[1] = 2;
planes[2] = 4;
planes[3] = 5;
break;
case 13:
planes[0] = 0;
planes[1] = 2;
planes[2] = 3;
planes[3] = 5;
break;
case 14:
planes[0] = 0;
planes[1] = 1;
planes[2] = 3;
planes[3] = 5;
break;
}
}
void AliTRDtrackerV1::GetExtrapolationConfig(Int_t iconfig, Int_t planes[2])
{
switch(iconfig){
case 0:
planes[0] = 1;
planes[1] = 0;
break;
case 1:
planes[0] = 5;
planes[1] = 0;
break;
case 2:
planes[0] = 4;
planes[1] = 5;
break;
case 3:
planes[0] = 4;
planes[1] = 0;
break;
case 4:
planes[0] = 5;
planes[1] = 3;
break;
case 5:
planes[0] = 2;
planes[1] = 0;
break;
case 6:
planes[0] = 5;
planes[1] = 1;
break;
case 7:
planes[0] = 2;
planes[1] = 1;
break;
case 8:
planes[0] = 4;
planes[1] = 3;
break;
case 9:
planes[0] = 3;
planes[1] = 0;
break;
case 10:
planes[0] = 5;
planes[1] = 2;
break;
case 11:
planes[0] = 3;
planes[1] = 2;
break;
case 12:
planes[0] = 3;
planes[1] = 1;
break;
case 13:
planes[0] = 4;
planes[1] = 1;
break;
case 14:
planes[0] = 4;
planes[1] = 2;
break;
}
}
AliCluster* AliTRDtrackerV1::GetCluster(Int_t idx) const
{
if(!fClusters) return NULL;
Int_t ncls = fClusters->GetEntriesFast();
return idx >= 0 && idx < ncls ? (AliCluster*)fClusters->UncheckedAt(idx) : NULL;
}
AliTRDseedV1* AliTRDtrackerV1::GetTracklet(Int_t idx) const
{
if(!fTracklets) return NULL;
Int_t ntrklt = fTracklets->GetEntriesFast();
return idx >= 0 && idx < ntrklt ? (AliTRDseedV1*)fTracklets->UncheckedAt(idx) : NULL;
}
AliKalmanTrack* AliTRDtrackerV1::GetTrack(Int_t idx) const
{
if(!fTracks) return NULL;
Int_t ntrk = fTracks->GetEntriesFast();
return idx >= 0 && idx < ntrk ? (AliKalmanTrack*)fTracks->UncheckedAt(idx) : NULL;
}
void AliTRDtrackerV1::ResetSeedTB()
{
for(Int_t isl=0; isl<kNSeedPlanes; isl++){
if(!fSeedTB[isl]) fSeedTB[isl] = new AliTRDchamberTimeBin();
else fSeedTB[isl]->Clear();
}
}
Float_t AliTRDtrackerV1::GetChi2Y(const AliTRDseedV1 * const tracklets) const
{
Double_t chi2 = 0.; Int_t n = 0;
for(Int_t ipl = kNPlanes; ipl--;){
if(!tracklets[ipl].IsOK()) continue;
chi2 += tracklets[ipl].GetChi2Y();
n++;
}
return n ? chi2/n : 0.;
}
Float_t AliTRDtrackerV1::GetChi2Z(const AliTRDseedV1 *const tracklets) const
{
Double_t chi2 = 0; Int_t n = 0;
for(Int_t ipl = kNPlanes; ipl--;){
if(!tracklets[ipl].IsOK()) continue;
chi2 += tracklets[ipl].GetChi2Z();
n++;
}
return n ? chi2/n : 0.;
}
Float_t AliTRDtrackerV1::GetChi2Phi(const AliTRDseedV1 *const tracklets) const
{
Double_t chi2 = 0; Int_t n = 0;
for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
if(!tracklets[iLayer].IsOK()) continue;
chi2 += tracklets[iLayer].GetChi2Phi();
n++;
}
return n ? chi2/n: 0.;
}
Float_t AliTRDtrackerV1::CalculateReferenceX(const AliTRDseedV1 *const tracklets){
Int_t nDistances = 0;
Float_t meanDistance = 0.;
Int_t startIndex = 5;
for(Int_t il =5; il > 0; il--){
if(tracklets[il].IsOK() && tracklets[il -1].IsOK()){
Float_t xdiff = tracklets[il].GetX0() - tracklets[il -1].GetX0();
meanDistance += xdiff;
nDistances++;
}
if(tracklets[il].IsOK()) startIndex = il;
}
if(tracklets[0].IsOK()) startIndex = 0;
if(!nDistances){
Float_t xpos[2]; memset(xpos, 0, sizeof(Float_t) * 2);
Int_t iok = 0, idiff = 0;
for(Int_t il = 5; il >= 0; il--){
if(tracklets[il].IsOK()){
xpos[iok] = tracklets[il].GetX0();
iok++;
startIndex = il;
}
if(iok) idiff++;
if(iok > 1) break;
}
if(iok > 1){
meanDistance = (xpos[0] - xpos[1])/idiff;
}
else{
return 331.;
}
}
else{
meanDistance /= nDistances;
}
return tracklets[startIndex].GetX0() + (2.5 - startIndex) * meanDistance - 0.5 * (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
}
Double_t AliTRDtrackerV1::FitTiltedRiemanV1(AliTRDseedV1 *const tracklets){
AliTRDtrackFitterRieman fitter;
fitter.SetRiemanFitter(GetTiltedRiemanFitter());
fitter.Reset();
for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++) fitter.SetTracklet(il, &tracklets[il]);
Double_t chi2 = fitter.Eval();
Double_t cov[15]; Double_t x0;
memset(cov, 0, sizeof(Double_t) * 15);
for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
x0 = tracklets[il].GetX0();
tracklets[il].SetYref(0, fitter.GetYat(x0));
tracklets[il].SetZref(0, fitter.GetZat(x0));
tracklets[il].SetYref(1, fitter.GetDyDxAt(x0));
tracklets[il].SetZref(1, fitter.GetDzDx());
tracklets[il].SetC(fitter.GetCurvature());
fitter.GetCovAt(x0, cov);
tracklets[il].SetCovRef(cov);
tracklets[il].SetChi2(chi2);
}
return chi2;
}
void AliTRDtrackerV1::UnsetTrackletsTrack(const AliTRDtrackV1 * const track)
{
Int_t idx(-1);
for(Int_t il(0); il<kNPlanes; il++){
if((idx = track->GetTrackletIndex(il)) < 0) continue;
delete (fTracklets->RemoveAt(idx));
}
}
AliTRDtrackerV1::AliTRDLeastSquare::AliTRDLeastSquare(){
memset(fParams, 0, sizeof(Double_t) * 2);
memset(fSums, 0, sizeof(Double_t) * 6);
memset(fCovarianceMatrix, 0, sizeof(Double_t) * 3);
}
void AliTRDtrackerV1::AliTRDLeastSquare::AddPoint(const Double_t *const x, Double_t y, Double_t sigmaY){
Double_t weight = 1/(sigmaY > 1e-9 ? sigmaY : 1e-9);
weight *= weight;
const Double_t &xpt = *x;
fSums[0] += weight;
fSums[1] += weight * xpt;
fSums[2] += weight * y;
fSums[3] += weight * xpt * y;
fSums[4] += weight * xpt * xpt;
fSums[5] += weight * y * y;
}
void AliTRDtrackerV1::AliTRDLeastSquare::RemovePoint(const Double_t *const x, Double_t y, Double_t sigmaY){
Double_t weight = 1/(sigmaY > 1e-9 ? sigmaY : 1e-9);
weight *= weight;
const Double_t &xpt = *x;
fSums[0] -= weight;
fSums[1] -= weight * xpt;
fSums[2] -= weight * y;
fSums[3] -= weight * xpt * y;
fSums[4] -= weight * xpt * xpt;
fSums[5] -= weight * y * y;
}
Bool_t AliTRDtrackerV1::AliTRDLeastSquare::Eval(){
Double_t det = fSums[0] * fSums[4] - fSums[1] *fSums[1];
if(TMath::Abs(det)<1.e-30) return kFALSE;
fParams[0] = (fSums[2] * fSums[4] - fSums[1] * fSums[3])/det;
fParams[1] = (fSums[0] * fSums[3] - fSums[1] * fSums[2])/det;
Double_t den = fSums[0]*fSums[4] - fSums[1]*fSums[1];
fCovarianceMatrix[0] = fSums[4] / den;
fCovarianceMatrix[1] = fSums[0] / den;
fCovarianceMatrix[2] = -fSums[1] / den;
return kTRUE;
}
Double_t AliTRDtrackerV1::AliTRDLeastSquare::GetFunctionValue(const Double_t *const xpos) const {
return fParams[0] + fParams[1] * (*xpos);
}
void AliTRDtrackerV1::AliTRDLeastSquare::GetCovarianceMatrix(Double_t *storage) const {
memcpy(storage, fCovarianceMatrix, sizeof(Double_t) * 3);
}
void AliTRDtrackerV1::AliTRDLeastSquare::Reset(){
memset(fParams, 0, sizeof(Double_t) * 2);
memset(fCovarianceMatrix, 0, sizeof(Double_t) * 3);
memset(fSums, 0, sizeof(Double_t) * 6);
}
AliTRDtrackerV1::AliTRDtrackFitterRieman::AliTRDtrackFitterRieman():
fTrackFitter(NULL),
fZfitter(NULL),
fCovarPolY(NULL),
fCovarPolZ(NULL),
fXref(0.),
fSysClusterError(0.)
{
fZfitter = new AliTRDLeastSquare;
fCovarPolY = new TMatrixD(3,3);
fCovarPolZ = new TMatrixD(2,2);
memset(fTracklets, 0, sizeof(AliTRDseedV1 *) * 6);
memset(fParameters, 0, sizeof(Double_t) * 5);
memset(fSumPolY, 0, sizeof(Double_t) * 5);
memset(fSumPolZ, 0, sizeof(Double_t) * 2);
}
AliTRDtrackerV1::AliTRDtrackFitterRieman::~AliTRDtrackFitterRieman(){
if(fZfitter) delete fZfitter;
if(fCovarPolY) delete fCovarPolY;
if(fCovarPolZ) delete fCovarPolZ;
}
void AliTRDtrackerV1::AliTRDtrackFitterRieman::Reset(){
if(fTrackFitter){
fTrackFitter->StoreData(kTRUE);
fTrackFitter->ClearPoints();
}
if(fZfitter){
fZfitter->Reset();
}
fXref = 0.;
memset(fTracklets, 0, sizeof(AliTRDseedV1 *) * AliTRDgeometry::kNlayer);
memset(fParameters, 0, sizeof(Double_t) * 5);
memset(fSumPolY, 0, sizeof(Double_t) * 5);
memset(fSumPolZ, 0, sizeof(Double_t) * 2);
for(Int_t irow = 0; irow < fCovarPolY->GetNrows(); irow++)
for(Int_t icol = 0; icol < fCovarPolY->GetNcols(); icol++){
(*fCovarPolY)(irow, icol) = 0.;
if(irow < 2 && icol < 2)
(*fCovarPolZ)(irow, icol) = 0.;
}
}
void AliTRDtrackerV1::AliTRDtrackFitterRieman::SetTracklet(Int_t itr, AliTRDseedV1 *tracklet){
if(itr >= AliTRDgeometry::kNlayer) return;
fTracklets[itr] = tracklet;
}
Double_t AliTRDtrackerV1::AliTRDtrackFitterRieman::Eval(){
if(!fTrackFitter){
return 1e10;
}
fXref = CalculateReferenceX();
for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++) UpdateFitters(fTracklets[il]);
if(!fTrackFitter->GetNpoints()) return 1e10;
fTrackFitter->Eval();
fZfitter->Eval();
fParameters[3] = fTrackFitter->GetParameter(3);
fParameters[4] = fTrackFitter->GetParameter(4);
if(!CheckAcceptable(fParameters[3], fParameters[4])) {
fTrackFitter->FixParameter(3, fZfitter->GetFunctionValue(&fXref));
fTrackFitter->FixParameter(4, fZfitter->GetFunctionParameter(1));
fTrackFitter->Eval();
fTrackFitter->ReleaseParameter(3);
fTrackFitter->ReleaseParameter(4);
fParameters[3] = fTrackFitter->GetParameter(3);
fParameters[4] = fTrackFitter->GetParameter(4);
}
fParameters[0] = fTrackFitter->GetParameter(0);
fParameters[1] = fTrackFitter->GetParameter(1);
fParameters[2] = fTrackFitter->GetParameter(2);
(*fCovarPolY)(0,0) = fSumPolY[0]; (*fCovarPolY)(1,1) = fSumPolY[2]; (*fCovarPolY)(2,2) = fSumPolY[4];
(*fCovarPolY)(1,0) = (*fCovarPolY)(0,1) = fSumPolY[1];
(*fCovarPolY)(2,0) = (*fCovarPolY)(0,2) = fSumPolY[2];
(*fCovarPolY)(2,1) = (*fCovarPolY)(1,2) = fSumPolY[3];
fCovarPolY->Invert();
(*fCovarPolZ)(0,0) = fSumPolZ[0]; (*fCovarPolZ)(1,1) = fSumPolZ[2];
(*fCovarPolZ)(1,0) = (*fCovarPolZ)(0,1) = fSumPolZ[1];
fCovarPolZ->Invert();
return fTrackFitter->GetChisquare() / fTrackFitter->GetNpoints();
}
void AliTRDtrackerV1::AliTRDtrackFitterRieman::UpdateFitters(const AliTRDseedV1 * const tracklet){
AliTRDcluster *cl = NULL;
Double_t x, y, z, dx, t, w, we, yerr, zerr;
Double_t uvt[4];
if(!tracklet || !tracklet->IsOK()) return;
Double_t tilt = tracklet->GetTilt();
for(Int_t itb = 0; itb < AliTRDseedV1::kNclusters; itb++){
if(!(cl = tracklet->GetClusters(itb))) continue;
if(!cl->IsInChamber()) continue;
if (!tracklet->IsUsable(itb)) continue;
x = cl->GetX();
y = cl->GetY();
z = cl->GetZ();
dx = x - fXref;
t = 1./(x*x + y*y);
uvt[0] = 2. * x * t;
uvt[1] = t;
uvt[2] = 2. * tilt * t;
uvt[3] = 2. * tilt * dx * t;
w = 2. * (y + tilt*z) * t;
we = 2. * t;
we *= TMath::Sqrt(cl->GetSigmaY2()+tilt*tilt*cl->GetSigmaZ2());
yerr = 1./(TMath::Sqrt(cl->GetSigmaY2()) + fSysClusterError);
yerr *= yerr;
zerr = 1./cl->GetSigmaZ2();
for(Int_t ipol = 0; ipol < 5; ipol++){
fSumPolY[ipol] += yerr;
yerr *= x;
if(ipol < 3){
fSumPolZ[ipol] += zerr;
zerr *= x;
}
}
fTrackFitter->AddPoint(uvt, w, we);
fZfitter->AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
}
}
Bool_t AliTRDtrackerV1::AliTRDtrackFitterRieman::CheckAcceptable(Double_t offset, Double_t slope){
Bool_t acceptablez = kTRUE;
Double_t zref = 0.0;
for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
if(!fTracklets[iLayer]->IsOK()) continue;
zref = offset + slope * (fTracklets[iLayer]->GetX0() - fXref);
if (TMath::Abs(fTracklets[iLayer]->GetZfit(0) - zref) > fTracklets[iLayer]->GetPadLength() * 0.5 + 1.0)
acceptablez = kFALSE;
}
return acceptablez;
}
Double_t AliTRDtrackerV1::AliTRDtrackFitterRieman::GetYat(Double_t x) const {
Double_t y = 0;
Double_t disc = (x * fParameters[0] + fParameters[1]);
disc = 1 - fParameters[0]*fParameters[2] + fParameters[1]*fParameters[1] - disc*disc;
if (disc >= 0) {
disc = TMath::Sqrt(disc);
y = (1.0 - disc) / fParameters[0];
}
return y;
}
Double_t AliTRDtrackerV1::AliTRDtrackFitterRieman::GetZat(Double_t x) const {
return fParameters[3] + fParameters[4] * (x - fXref);
}
Double_t AliTRDtrackerV1::AliTRDtrackFitterRieman::GetDyDxAt(Double_t x) const {
Double_t x0 = -fParameters[1] / fParameters[0];
Double_t curvature = GetCurvature();
Double_t dy = 0;
if (-fParameters[2] * fParameters[0] + fParameters[1] * fParameters[1] + 1 > 0) {
if (1.0/(curvature * curvature) - (x - x0) * (x - x0) > 0.0) {
Double_t yderiv = (x - x0) / TMath::Sqrt(1.0/(curvature * curvature) - (x - x0) * (x - x0));
if (fParameters[0] < 0) yderiv *= -1.0;
dy = yderiv;
}
}
return dy;
}
Double_t AliTRDtrackerV1::AliTRDtrackFitterRieman::GetCurvature() const {
Double_t curvature = 1.0 + fParameters[1]*fParameters[1] - fParameters[2]*fParameters[0];
if (curvature > 0.0)
curvature = fParameters[0] / TMath::Sqrt(curvature);
return curvature;
}
void AliTRDtrackerV1::AliTRDtrackFitterRieman::GetCovAt(Double_t x, Double_t *cov) const {
TMatrixD transform(3,3);
transform(0,0) = transform(1,1) = transform(2,2) = 1;
transform(0,1) = transform(1,2) = x;
transform(0,2) = x*x;
TMatrixD covariance(transform, TMatrixD::kMult, *fCovarPolY);
covariance *= transform.T();
cov[0] = covariance(0,0);
TMatrixD transformZ(2,2);
transformZ(0,0) = transformZ(1,1) = 1;
transformZ(0,1) = x;
TMatrixD covarZ(transformZ, TMatrixD::kMult, *fCovarPolZ);
covarZ *= transformZ.T();
cov[1] = covarZ(0,0);
cov[2] = 0;
}
Double_t AliTRDtrackerV1::AliTRDtrackFitterRieman::CalculateReferenceX(){
Int_t nDistances = 0;
Float_t meanDistance = 0.;
Int_t startIndex = 5;
for(Int_t il =5; il > 0; il--){
if(fTracklets[il]->IsOK() && fTracklets[il -1]->IsOK()){
Float_t xdiff = fTracklets[il]->GetX0() - fTracklets[il -1]->GetX0();
meanDistance += xdiff;
nDistances++;
}
if(fTracklets[il]->IsOK()) startIndex = il;
}
if(fTracklets[0]->IsOK()) startIndex = 0;
if(!nDistances){
Float_t xpos[2]; memset(xpos, 0, sizeof(Float_t) * 2);
Int_t iok = 0, idiff = 0;
for(Int_t il = 5; il >= 0; il--){
if(fTracklets[il]->IsOK()){
xpos[iok] = fTracklets[il]->GetX0();
iok++;
startIndex = il;
}
if(iok) idiff++;
if(iok > 1) break;
}
if(iok > 1){
meanDistance = (xpos[0] - xpos[1])/idiff;
}
else{
return 331.;
}
}
else{
meanDistance /= nDistances;
}
return fTracklets[startIndex]->GetX0() + (2.5 - startIndex) * meanDistance - 0.5 * (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
}