BmnRoot
Loading...
Searching...
No Matches
CbmStsKFTrackFitter.cxx
Go to the documentation of this file.
2
3#include "CbmKFMath.h"
4#include "CbmKFTrack.h"
5#include "CbmKFVertex.h"
6#include "CbmKFStsHit.h"
7
8#include "FairRootManager.h"
9#include "CbmStsTrack.h"
10#include "CbmVertex.h"
11#include "CbmMvdHit.h"
12
13#include "TClonesArray.h"
14#include "TMath.h"
15
16#include <iostream>
17#include "math.h"
18
19using std::cout;
20using std::endl;
21
23 fHits(),
24 fMvdHitsArray(0),
25 fStsHitsArray(0),
27 fCheckTrigSi(CheckTrigSi)
28{}
29
31 // Initialisation
32 FairRootManager *rootMgr = FairRootManager::Instance();
33 if(NULL == rootMgr) {
34 cout << "-E- CbmStsKFTrackFitter::Init(): "
35 << "ROOT manager is not instantiated!" << endl;
36 return;
37 }
38 if(!fCheckTrigSi){
39 fStsHitsArray = reinterpret_cast<TClonesArray*>(rootMgr->GetObject("StsHit"));
40 if ( !fStsHitsArray ) {
41 cout << "-W- CbmStsKFTrackFitter::Init: "
42 << "no STS hits array" << endl;
43 //return;
44 }
45 } else {
46
47 fStsHitsArray = reinterpret_cast<TClonesArray*>(rootMgr->GetObject("StsHitSi"));
48 if ( !fStsHitsArray ) {
49 cout << "-W- CbmStsKFTrackFitter::Init: "
50 << "no STS hits array" << endl;
51 //return;
52 }
53 }
54 fMvdHitsArray = reinterpret_cast<TClonesArray*>(rootMgr->GetObject("MvdHit"));
55 if( !fMvdHitsArray ) {
56 cout << "-W- CbmStsKFTrackFitter::Init: "
57 << "no MVD hits array" << endl;
58 //return;
59 }
60 fIsInitialised = 1;
61};
62
64
65 T.fHits.clear();
66
67 if( !fIsInitialised ) Init();
68
69 Int_t NStsHits = ( fStsHitsArray ) ?track->GetNStsHits() :0;
70 Int_t NMvdHits= ( fMvdHitsArray) ?track->GetNMvdHits() :0;
71//cout<<" hits arr size: "<<NStsHits<<endl;
72 fHits.resize( NMvdHits + NStsHits );
73 if ( NMvdHits > 0 ){
74 for (Int_t i=0; i<NMvdHits;i++){
75 Int_t j = track->GetMvdHitIndex(i);
76 fHits[i].Create( reinterpret_cast<CbmMvdHit*>(fMvdHitsArray->At(j)) );
77 T.fHits.push_back(&(fHits[i]));
78 }
79 }
80 if( NStsHits>0 && fStsHitsArray ){
81 for (Int_t i=0; i<NStsHits;i++){
82 Int_t j = track->GetStsHitIndex(i);
83 // cout<<" hits arr index: "<<j<<endl;
84 if(!bkg) fHits[NMvdHits+i].Create( reinterpret_cast<CbmStsHit*>(fStsHitsArray->At(j)) );
85 //cout<<" hit: "<<(CbmStsHit*)track->GetStsHitArr()->At(j)<<endl;
86 if(bkg) fHits[NMvdHits+i].Create( reinterpret_cast<CbmStsHit*>(track->GetStsHitArr()->At(i)) );
87 T.fHits.push_back(&(fHits[NMvdHits+i]));
88 }
89 }
90}
91
92Int_t CbmStsKFTrackFitter::DoFit( CbmStsTrack* track, Int_t pidHypo )
93{
94 track->SetPidHypo(pidHypo);
95
96 CbmKFTrack T;
97 T.SetPID(pidHypo);
98 SetKFHits( T, track);
99 for( Int_t i=0; i<6; i++ ) T.GetTrack()[i] = 0.; // no guess taken
100 T.Fit( 1 ); // fit downstream
101 CheckTrack( T );
102 T.Fit( 0 ); // fit upstream
103 CheckTrack( T );
104 Int_t err = T.Fit( 1 ); // fit downstream
105 Bool_t ok = (!err) && CheckTrack( T );
106 if( ok ){
107 T.GetTrackParam( *track->GetParamLast() ); // store fitted track & cov.matrix
108 err = T.Fit( 0 ); // fit upstream
109 ok = ok && (!err) && CheckTrack( T );
110 if( ok ) T.GetStsTrack( *track, 1 ); // store fitted track & cov.matrix & chi2 & NDF
111 }
112 if( !ok ){
113 Double_t *t = T.GetTrack();
114 Double_t *c = T.GetCovMatrix();
115 for( int i=0; i<6 ; i++) t[i] = 0;
116 for( int i=0; i<15; i++) c[i] = 0;
117 c[0] = c[2] = c[5] = c[9] = c[14] = 100.;
118 T.GetRefChi2() = 100.;
119 T.GetRefNDF() = 0;
120 T.GetStsTrack( *track, 0 );
121 T.GetStsTrack( *track, 1 );
122 track->SetFlag(1);
123 } else{
124 track->SetFlag(0);
125 }
126 return !ok;
127}
128
129
130void CbmStsKFTrackFitter::Extrapolate( FairTrackParam* track, Double_t z, FairTrackParam* e_track )
131{
132 if( !track ) return;
133 CbmKFTrack T;
134 T.SetTrackParam( *track );
135 T.Extrapolate( z );
136 if( e_track ) T.GetTrackParam( *e_track );
137}
138
139
140void CbmStsKFTrackFitter::Extrapolate( CbmStsTrack* track, Double_t z, FairTrackParam* e_track )
141{
142 if( !track ) return;
143 CbmKFTrack T;
144 T.SetPID( track->GetPidHypo() );
145 FairTrackParam *fpar = track->GetParamFirst(), *lpar = track->GetParamLast();
146 if( z<=fpar->GetZ() ){ // extrapolate first parameters
147 T.SetTrackParam( *fpar );
148 T.Extrapolate( z );
149 }else if( z<fpar->GetZ()+0.1 ){ // extrapolate first parameters
150 T.SetTrackParam( *fpar );
151 T.Propagate( z );
152 }else if( lpar->GetZ()<=z ){ // extrapolate last parameters
153 T.SetTrackParam( *lpar );
154 T.Extrapolate( z );
155 }else if( lpar->GetZ()-0.1<z ){ // extrapolate last parameters
156 T.SetTrackParam( *lpar );
157 T.Propagate( z );
158 }else { // refit with smoother
159 SetKFHits( T, track);
160 T.SetTrackParam( *fpar );
161 T.Smooth( z );
162 }
163 if( e_track ) T.GetTrackParam( *e_track );
164}
165
166
168{
169 if( !vtx ){
170 FairRootManager *fManger = FairRootManager::Instance();
171 vtx = reinterpret_cast<CbmVertex *>(fManger->GetObject("PrimaryVertex"));
172 if( !vtx ){
173 cout<< "-W- CbmStsKFTrackFitter::GetChiToVertex: No Primary Vertex found!"<<endl;
174 return 100.;
175 }
176 }
177 CbmKFTrack T;
178 T.SetStsTrack( *track, 1 );
179 T.Extrapolate( vtx->GetZ() );
180
181 TMatrixFSym tmp(3);
182 vtx->CovMatrix( tmp );
183 Double_t Cv[3] = { tmp(0,0), tmp(0,1), tmp(1,1) };
184
185 return CbmKFMath::getDeviation( T.GetTrack()[0], T.GetTrack()[1], T.GetCovMatrix(),
186 vtx->GetX(), vtx->GetY(), Cv );
187}
188
189
190Double_t CbmStsKFTrackFitter::FitToVertex( CbmStsTrack* track, CbmVertex *vtx, FairTrackParam *v_track )
191{
192 Double_t ret = 100.;
193 if( !track || !vtx || !v_track ) return ret;
194 CbmKFTrack T( *track );
195 CbmKFVertex V( *vtx );
196 T.Fit2Vertex( V );
197 T.GetTrackParam( *v_track );
198 if( T.GetRefNDF()>0 && T.GetRefChi2()>=0 ){
199 ret = T.GetRefChi2()/T.GetRefNDF();
200 if( finite(ret) ) ret = sqrt(ret);
201 }
202 return ret;
203}
204
205Bool_t CbmStsKFTrackFitter::CheckTrack( CbmKFTrack &T )
206{
207 Bool_t ok = 1;
208 Double_t *t = T.GetTrack();
209 Double_t *c = T.GetCovMatrix();
210 for( int i=0; i<6; i++) ok = ok && finite(t[i]) && TMath::Abs(t[i])<1.e5;
211 for( int i=0; i<15; i++) ok = ok && finite(c[i]);
212 ok = ok && finite( T.GetMass() ) && finite( T.GetRefChi2() );
213 if( ok ){
214 ok = ok && ( c[ 0] > 0 );
215 ok = ok && ( c[ 2] > 0 );
216 ok = ok && ( c[ 5] > 0 );
217 ok = ok && ( c[ 9] > 0 );
218 ok = ok && ( c[14] > 0 );
219 }
220 if( ok ){ // correct the cov matrix
221 Double_t c00 = TMath::Sqrt(c[ 0]);
222 Double_t c11 = TMath::Sqrt(c[ 2]);
223 Double_t c22 = TMath::Sqrt(c[ 5]);
224 Double_t c33 = TMath::Sqrt(c[ 9]);
225 Double_t c44 = TMath::Sqrt(c[14]);
226 Double_t a = c11*c00;
227 if( c[ 1] > a ) c[ 1] = a;
228 if( c[ 1] < -a ) c[ 1] = -a;
229 a = c22*c00;
230 if( c[ 3] > a ) c[ 3] = a;
231 if( c[ 3] < -a ) c[ 3] = -a;
232 a = c22*c11;
233 if( c[ 4] > a ) c[ 4] = a;
234 if( c[ 4] < -a ) c[ 4] = -a;
235 a = c33*c00;
236 if( c[ 6] > a ) c[ 6] = a;
237 if( c[ 6] < -a ) c[ 6] = -a;
238 a = c33*c11;
239 if( c[ 7] > a ) c[ 7] = a;
240 if( c[ 7] < -a ) c[ 7] = -a;
241 a = c33*c22;
242 if( c[ 8] > a ) c[ 8] = a;
243 if( c[ 8] < -a ) c[ 8] = -a;
244 a = c44*c00;
245 if( c[10] > a ) c[10] = a;
246 if( c[10] < -a ) c[10] = -a;
247 a = c44*c11;
248 if( c[11] > a ) c[11] = a;
249 if( c[11] < -a ) c[11] = -a;
250 a = c44*c22;
251 if( c[12] > a ) c[12] = a;
252 if( c[12] < -a ) c[12] = -a;
253 a = c44*c33;
254 if( c[13] > a ) c[13] = a;
255 if( c[13] < -a ) c[13] = -a;
256 }
257 if( !ok ){
258 for( int i=0; i<6 ; i++) t[i] = 0;
259 for( int i=0; i<15; i++) c[i] = 0;
260 c[0] = c[2] = c[5] = c[9] = c[14] = 100.;
261 T.GetRefChi2() = 100.;
262 T.GetRefNDF() = 0;
263 }
264 return ok;
265}
Bool_t fIsInitialised
kTRUE if Init() was called
friend F32vec4 sqrt(const F32vec4 &a)
Definition P4_F32vec4.h:34
int i
Definition P4_F32vec4.h:22
static Double_t getDeviation(Double_t x, Double_t y, Double_t C[], Double_t vx, Double_t vy, Double_t Cv[]=0)
Int_t Propagate(Double_t z_out, Double_t QP0, Bool_t line=false)
Int_t Extrapolate(Double_t z, Double_t *QP0=0, Bool_t line=false)
Access to i-th hit.
Int_t Fit(Bool_t downstream=1, Bool_t line=false)
void Fit2Vertex(CbmKFVertexInterface &vtx)
void GetStsTrack(CbmStsTrack &track, bool first=1)
std::vector< CbmKFHit * > fHits
Definition CbmKFTrack.h:31
Double_t GetMass()
Definition CbmKFTrack.h:55
void SetStsTrack(CbmStsTrack &track, bool first=1)
void SetTrackParam(FairTrackParam &track)
void GetTrackParam(FairTrackParam &track)
Double_t * GetTrack()
Is it electron.
Definition CbmKFTrack.h:51
Double_t & GetRefChi2()
array[15] of covariance matrix
Definition CbmKFTrack.h:53
Int_t & GetRefNDF()
Chi^2 after fit.
Definition CbmKFTrack.h:54
Double_t * GetCovMatrix()
array[6] of track parameters(x,y,tx,ty,qp,z)
Definition CbmKFTrack.h:52
void SetPID(Int_t pidHypo)
void SetKFHits(CbmKFTrack &T, CbmStsTrack *track)
CbmStsKFTrackFitter(Bool_t CheckTrigSi=false)
Double_t GetChiToVertex(CbmStsTrack *track, CbmVertex *vtx=0)
void Extrapolate(CbmStsTrack *track, Double_t z, FairTrackParam *e_track)
Double_t FitToVertex(CbmStsTrack *track, CbmVertex *vtx, FairTrackParam *v_track)
Int_t DoFit(CbmStsTrack *track, Int_t pidHypo=211)
FairTrackParam * GetParamLast()
Definition CbmStsTrack.h:70
Int_t GetNMvdHits() const
Definition CbmStsTrack.h:61
Int_t GetPidHypo() const
Definition CbmStsTrack.h:64
Int_t GetMvdHitIndex(Int_t iHit) const
Definition CbmStsTrack.h:63
TClonesArray * GetStsHitArr()
Definition CbmStsTrack.h:73
void SetFlag(Int_t flag)
Definition CbmStsTrack.h:85
Int_t GetNStsHits() const
Definition CbmStsTrack.h:60
Int_t GetStsHitIndex(Int_t iHit) const
Definition CbmStsTrack.h:62
FairTrackParam * GetParamFirst()
Definition CbmStsTrack.h:69
void SetPidHypo(Int_t pid)
Definition CbmStsTrack.h:82
Double_t GetZ() const
Definition CbmVertex.h:60
void CovMatrix(TMatrixFSym &covMat) const
Double_t GetX() const
Definition CbmVertex.h:58
Double_t GetY() const
Definition CbmVertex.h:59