FEDRA emulsion software from the OPERA Collaboration
EdbTrackFitter Class Reference

#include <EdbTrackFitter.h>

Inheritance diagram for EdbTrackFitter:
Collaboration diagram for EdbTrackFitter:

Public Member Functions

float Chi2SegM (EdbSegP s1, EdbSegP s2, EdbSegP &s, EdbScanCond &cond1, EdbScanCond &cond2)
 
float Chi2SegMCS (const EdbSegP &s1, const EdbSegP &s2)
 
 EdbTrackFitter ()
 
int FitTrackLine (const EdbTrackP &tr, float &x, float &y, float &z, float &tx, float &ty, float &w)
 
int FitTrackLine (EdbTrackP &tr)
 
float PMS_KF (EdbTrackP &t, float p0=10., float probbest=0.5)
 
void Print ()
 
double ProbSegMCS (EdbSegP *s1, EdbSegP *s2)
 
void SetDefaultBrick ()
 
void SetNsegMax (int nseg)
 
bool SplitTrack (EdbTrackP &t, EdbTrackP &t1, int isplit)
 
int SplitTrackByKink (EdbTrackP *t, TObjArray &tracks, float maxkink)
 
virtual ~EdbTrackFitter ()
 

Static Public Member Functions

static float Chi2ACP (EdbSegP s1, EdbSegP s2, EdbScanCond &cond)
 
static float Chi2ASeg (EdbSegP &s1, EdbSegP &s2, EdbSegP &s, EdbScanCond &cond1, EdbScanCond &cond2)
 
static float Chi2ASegLL (EdbSegP &s1, EdbSegP &s2, EdbSegP &s, EdbScanCond &cond1, EdbScanCond &cond2)
 
static float Chi2PSeg (EdbSegP &s1, EdbSegP &s2, EdbSegP &seg, EdbScanCond &cond1, EdbScanCond &cond2)
 
static float Chi2Seg (EdbSegP *s1, EdbSegP *s2)
 
static float MaxChi2Seg (EdbTrackP &tr)
 
static float MaxKink (EdbTrackP &tr)
 
static float MeanChi2Seg (EdbTrackP &tr)
 
static float MeanKink (EdbTrackP &tr)
 
static int RMSprojXY (EdbTrackP &tr, float &ex, float &ey)
 
static float Theta (EdbSegP &s, EdbSegP &s1)
 

Public Attributes

bool eDE_correction
 take into account the energy loss or not More...
 
float eM
 mass of the particle (if negative - use the mass setted in the track) More...
 
float ePcut
 minimal momentum? More...
 
float ePdef
 default momentum More...
 
float eTPb
 ? More...
 
float eX0
 rad length of the media [microns] More...
 

Private Attributes

int eNsegMax
 max number of segments (for arrays allocation) More...
 

Constructor & Destructor Documentation

◆ EdbTrackFitter()

EdbTrackFitter::EdbTrackFitter ( )
38{
40}
void SetDefaultBrick()
Definition: EdbTrackFitter.cxx:43

◆ ~EdbTrackFitter()

virtual EdbTrackFitter::~EdbTrackFitter ( )
inlinevirtual
32{}

Member Function Documentation

◆ Chi2ACP()

float EdbTrackFitter::Chi2ACP ( EdbSegP  s1,
EdbSegP  s2,
EdbScanCond cond 
)
static

Exact copy of the function traditionally used for the linking of couples

fast estimation of chi2 in the special case when the position
errors of segments are negligible in respect to angular errors:
sigmaXY/dz << sigmaTXY
application: up/down linking, alignment (offset search)
All calculation are done in the track plane which remove the
dependency of the polar angle (phi)

126{
136
137 TVector3 v1,v2,v;
138 v1.SetXYZ( s1.TX(), s1.TY() , -1. );
139 v2.SetXYZ( s2.TX(), s2.TY() , -1. );
140 v.SetXYZ( -( s2.X() - s1.X() ),
141 -( s2.Y() - s1.Y() ),
142 -( s2.Z() - s1.Z() ) );
143
144 float phi = v.Phi();
145 v.RotateZ( -phi );
146 v1.RotateZ( -phi );
147 v2.RotateZ( -phi );
148
149 float dz = v.Z();
150 float tx = v.X()/dz;
151 float ty = v.Y()/dz;
152 float stx = cond.SigmaTX(tx);
153 float sty = cond.SigmaTY(ty);
154
155 float prob1=1., prob2=1.;
156 //if(iprob) {
157 prob1 = cond.ProbSeg(tx,ty,s1.W());
158 prob2 = cond.ProbSeg(tx,ty,s2.W());
159 // }
160 float dtx1 = (v1.X()-tx)*(v1.X()-tx)/stx/stx/prob1;
161 float dty1 = (v1.Y()-ty)*(v1.Y()-ty)/sty/sty/prob1;
162 float dtx2 = (v2.X()-tx)*(v2.X()-tx)/stx/stx/prob2;
163 float dty2 = (v2.Y()-ty)*(v2.Y()-ty)/sty/sty/prob2;
164
165 float chi2a=TMath::Sqrt(dtx1+dty1+dtx2+dty2)/2.;
166 return chi2a;
167}
brick dz
Definition: RecDispMC.C:107
float SigmaTX(float ax) const
Definition: EdbScanCond.h:106
float SigmaTY(float ay) const
Definition: EdbScanCond.h:107
float ProbSeg(float tx, float ty, float puls) const
Definition: EdbScanCond.cxx:119
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t X() const
Definition: EdbSegP.h:173
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174
Float_t W() const
Definition: EdbSegP.h:151
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ Chi2ASeg()

float EdbTrackFitter::Chi2ASeg ( EdbSegP s1,
EdbSegP s2,
EdbSegP seg,
EdbScanCond cond1,
EdbScanCond cond2 
)
static

fast estimation of chi2 in the special case when the position
errors of segments are negligible in respect to angular errors:
sigmaXY/dz << sigmaTXY
application: up/down linking,

All calculation are done in the track plane which remove the
dependency of the polar angle (phi)

218{
223 //
226
227 Float_t dz = s2.Z() - s1.Z();
228 seg.Set( 0,
229 0.5*(s1.X() + s2.X()),
230 0.5*(s1.Y() + s2.Y()),
231 (s2.X() - s1.X())/dz,
232 (s2.Y() - s1.Y())/dz,
233 s1.W() + s2.W(),
234 0);
235 seg.SetZ( 0.5*(s2.Z()+s1.Z()) );
236
237 Float_t phi = -seg.Phi();
238 Double_t s = TMath::Sin(phi);
239 Double_t c = TMath::Cos(phi);
240 Float_t tx = seg.TX()*c-seg.TY()*s;
241 Float_t ty = seg.TX()*s+seg.TY()*c;
242 Float_t tx1 = s1.TX()*c-s1.TY()*s;
243 Float_t ty1 = s1.TX()*s+s1.TY()*c;
244 Float_t tx2 = s2.TX()*c-s2.TY()*s;
245 Float_t ty2 = s2.TX()*s+s2.TY()*c;
246
247 Float_t stx = cond1.SigmaTX(tx);
248 Float_t sty = cond1.SigmaTY(ty);
249 Float_t prob1 = cond1.ProbSeg(tx,ty,s1.W());
250 Float_t prob2 = cond2.ProbSeg(tx,ty,s2.W());
251
252 Float_t dtx1 = (tx1-tx)*(tx1-tx)/stx/stx/prob1;
253 Float_t dty1 = (ty1-ty)*(ty1-ty)/sty/sty/prob1;
254 Float_t dtx2 = (tx2-tx)*(tx2-tx)/stx/stx/prob2;
255 Float_t dty2 = (ty2-ty)*(ty2-ty)/sty/sty/prob2;
256 Float_t chi2a = Sqrt(dtx1+dty1+dtx2+dty2)/2.;
257
258 seg.SetChi2(chi2a);
259 return chi2a;
260}
void SetZ(float z)
Definition: EdbSegP.h:125
void SetChi2(float chi2)
Definition: EdbSegP.h:135
Float_t Phi() const
Definition: EdbSegP.h:183
void Set(int id, float x, float y, float tx, float ty, float w, int flag)
Definition: EdbSegP.h:87
s
Definition: check_shower.C:55

◆ Chi2ASegLL()

float EdbTrackFitter::Chi2ASegLL ( EdbSegP s1,
EdbSegP s2,
EdbSegP s,
EdbScanCond cond1,
EdbScanCond cond2 
)
static

fast estimation of chi2 in the special case when the position
errors of segments are negligible in respect to angular errors:
sigmaXY/dz << sigmaTXY
application: up/down linking,

All calculation are done in the track plane which remove the
dependency of the polar angle (phi)
In this function the Likelihood estimated by Andrey and stored in s1.eChi2 and s2.eChi2 are used

171{
180
181 Float_t dz = s2.Z() - s1.Z();
182 seg.Set( 0,
183 0.5*(s1.X() + s2.X()),
184 0.5*(s1.Y() + s2.Y()),
185 (s2.X() - s1.X())/dz,
186 (s2.Y() - s1.Y())/dz,
187 s1.W() + s2.W(),
188 0);
189 seg.SetZ( 0.5*(s2.Z()+s1.Z()) );
190
191 Float_t phi = -seg.Phi();
192 Double_t s = TMath::Sin(phi);
193 Double_t c = TMath::Cos(phi);
194 Float_t tx = seg.TX()*c-seg.TY()*s;
195 Float_t ty = seg.TX()*s+seg.TY()*c;
196 Float_t tx1 = s1.TX()*c-s1.TY()*s;
197 Float_t ty1 = s1.TX()*s+s1.TY()*c;
198 Float_t tx2 = s2.TX()*c-s2.TY()*s;
199 Float_t ty2 = s2.TX()*s+s2.TY()*c;
200
201 Float_t stx = cond1.SigmaTX(tx);
202 Float_t sty = cond1.SigmaTY(ty);
203 Float_t prob1 = cond1.ProbLL(tx,ty,s1.Chi2());
204 Float_t prob2 = cond2.ProbLL(tx,ty,s2.Chi2());
205
206 Float_t dtx1 = (tx1-tx)*(tx1-tx)/stx/stx/prob1;
207 Float_t dty1 = (ty1-ty)*(ty1-ty)/sty/sty/prob1;
208 Float_t dtx2 = (tx2-tx)*(tx2-tx)/stx/stx/prob2;
209 Float_t dty2 = (ty2-ty)*(ty2-ty)/sty/sty/prob2;
210 Float_t chi2a = Sqrt(dtx1+dty1+dtx2+dty2)/2.;
211
212 seg.SetChi2(chi2a);
213 return chi2a;
214}
float ProbLL(float tx, float ty, float puls) const
Definition: EdbScanCond.cxx:145
Float_t Chi2() const
Definition: EdbSegP.h:157

◆ Chi2PSeg()

float EdbTrackFitter::Chi2PSeg ( EdbSegP s1,
EdbSegP s2,
EdbSegP seg,
EdbScanCond cond1,
EdbScanCond cond2 
)
static

fast estimation of chi2 in the special case when only the position errors used - no any angular informtion:
application: alignment

264{
268
269 seg.Set( 0,
270 0.5*(s1.X() + s2.X()),
271 0.5*(s1.Y() + s2.Y()),
272 0.5*(s1.TX() + s2.TX()),
273 0.5*(s1.TY() + s2.TY()),
274 s1.W() + s2.W(),
275 0);
276 seg.SetZ( 0.5*(s2.Z()+s1.Z()) );
277
278 Double_t dx = s2.X()-s1.X();
279 Double_t dy = s2.Y()-s1.Y();
280 Float_t sx = cond1.SigmaX(0);
281 Float_t sy = cond1.SigmaY(0);
282 Double_t chi2 = Sqrt( dx*dx/sx/sx + dy*dy/sy/sy );
283 seg.SetChi2(chi2);
284 return chi2;
285}
float SigmaX(float ax) const
Definition: EdbScanCond.h:102
float SigmaY(float ay) const
Definition: EdbScanCond.h:103
Float_t chi2
Definition: testBGReduction_By_ANN.C:14

◆ Chi2Seg()

float EdbTrackFitter::Chi2Seg ( EdbSegP s1,
EdbSegP s2 
)
static

Estimation of chi2 using the covariance matrix

Return value: Prob: is Chi2 probability (area of the tail of Chi2-distribution)
If we accept couples with Prob >= ProbMin then ProbMin is the
probability to reject the good couple

The mass and momentum of the tr are used for multiple scattering estimation

63{
71
72 double dz;
73 float prob;
74 VtVector par( (double)(tr->X()),
75 (double)(tr->Y()),
76 (double)(tr->TX()),
77 (double)(tr->TY()) );
78 VtSymMatrix cov(4); // covariance matrix for seg0 (measurements errors)
79 for(int k=0; k<4; k++)
80 for(int l=0; l<4; l++) cov(k,l) = (tr->COV())(k,l);
81
82 Double_t chi2=0.;
83 dz = s->Z()-tr->Z();
84 VtSqMatrix pred(4); //propagation matrix for track parameters (x,y,tx,ty)
85 pred.clear();
86 pred(0,0) = 1.;
87 pred(1,1) = 1.;
88 pred(2,2) = 1.;
89 pred(3,3) = 1.;
90 pred(0,2) = dz;
91 pred(1,3) = dz;
92 VtVector parpred(4); // prediction from seg0 to seg
93 parpred = pred*par;
94 VtSymMatrix covpred(4); // covariance matrix for prediction
95 covpred = pred*(cov*pred.T());
96
97 VtSymMatrix dmeas(4); // original covariance matrix for seg2
98 for(int k=0; k<4; k++)
99 for(int l=0; l<4; l++) dmeas(k,l) = (s->COV())(k,l);
100
101 covpred = covpred.dsinv();
102 dmeas = dmeas.dsinv();
103 cov = covpred + dmeas;
104 cov = cov.dsinv();
105
106 VtVector meas( (double)(s->X()),
107 (double)(s->Y()),
108 (double)(s->TX()),
109 (double)(s->TY()) );
110
111 par = cov*(covpred*parpred + dmeas*meas); // new parameters for seg
112 chi2 = (par-parpred)*(covpred*(par-parpred)) + (par-meas)*(dmeas*(par-meas));
113 prob = (float)TMath::Prob(chi2,4);
114
115 tr->Set(tr->ID(),(float)par(0),(float)par(1),(float)par(2),(float)par(3),tr->W(),tr->Flag());
116 tr->SetCOV( cov.array(), 4 );
117 tr->SetChi2((float)chi2);
118 tr->SetProb(prob);
119 tr->SetZ(s->Z());
120 tr->SetW(tr->W()+s->W());
121 return TMath::Sqrt(chi2/4.);
122}
T Prob(const T &rhs, int n)
Definition: Prob.hh:37
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
Definition: VtSqMatrix.hh:50
Definition: VtSymMatrix.hh:49
Definition: VtVector.hh:45

◆ Chi2SegM()

float EdbTrackFitter::Chi2SegM ( EdbSegP  s1,
EdbSegP  s2,
EdbSegP s,
EdbScanCond cond1,
EdbScanCond cond2 
)

full estimation of chi2 without covariance matrix - the result seems to be identical to Chi2Seg
VT: 19-Sep-2007
Input: 2 segments passed by value because them will be modified during calculations
Return value - the result of the fit - passed by value

1) calcualte the mean direction vector for the 2-seg group

316{
322
324
325 float dz = s2.Z() - s1.Z();
326 float tbx=0, tby=0, wbx=0, wby=0;
327 if(Abs(dz) > 0.1 ) {
328 tbx = (s2.X() - s1.X())/(s2.Z() - s1.Z());
329 tby = (s2.Y() - s1.Y())/(s2.Z() - s1.Z());
330 wbx = Sqrt( cond1.SigmaX(tbx)*cond1.SigmaX(tbx) + cond2.SigmaX(tbx)*cond2.SigmaX(tbx) ) / dz;
331 wby = Sqrt( cond1.SigmaY(tby)*cond1.SigmaY(tby) + cond2.SigmaY(tby)*cond2.SigmaY(tby) ) / dz;
332 }
333 float w1x = 1./cond1.SigmaTX(s1.TX()); w1x*=w1x;
334 float w1y = 1./cond1.SigmaTY(s1.TY()); w1y*=w1y;
335 float w2x = 1./cond2.SigmaTX(s2.TX()); w2x*=w2x;
336 float w2y = 1./cond2.SigmaTY(s2.TY()); w2y*=w2y;
337
338 float TX = (s1.TX()*w1x + s2.TX()*w2x + tbx*wbx)/(w1x+w2x+wbx);
339 float TY = (s1.TY()*w1y + s2.TY()*w2y + tby*wby)/(w1y+w2y+wby);
340
341 // 2) calcualte the COG of the 2-seg group
342
343 w1x = 1./cond1.SigmaX(TX); w1x *= w1x;
344 w1y = 1./cond1.SigmaY(TY); w1y *= w1y;
345 w2x = 1./cond2.SigmaX(TX); w2x *= w2x;
346 w2y = 1./cond2.SigmaY(TY); w2y *= w2y;
347 float Z = (s1.Z()*(w1x+w1y) + s2.Z()*(w2x+w2y))/(w1x+w1y+w2x+w2y);
348 float X = (s1.X()*(w1x+w1y) + s2.X()*(w2x+w2y))/(w1x+w1y+w2x+w2y);
349 float Y = (s1.Y()*(w1x+w1y) + s2.Y()*(w2x+w2y))/(w1x+w1y+w2x+w2y);
350 //printf("COG: x %f y %f z %f\n",X,Y,Z);
351
352 s.SetX(X);
353 s.SetY(Y);
354 s.SetTX(TX);
355 s.SetTY(TY);
356 s.SetW(s1.W()+s2.W());
357 s.SetZ(Z);
358
359 float PHI = ATan2(TY,TX); // angle of the 2-seg group plane
360 float T = Sqrt(TX*TX+TY*TY);
361
362 EdbAffine2D aff;
363 aff.Rotate(PHI);
364 //aff.ShiftX(X);
365 //aff.ShiftY(Y);
366 aff.Invert();
367 s1.Transform(&aff);
368 s2.Transform(&aff);
369 s.Transform(&aff);
370
371 float stx1 = cond1.SigmaTX(T), sty1 = cond1.SigmaTY(0);
372 float stx2 = cond2.SigmaTX(T), sty2 = cond2.SigmaTY(0);
373 w1x = 1./(stx1*stx1); w1y = 1./(sty1*sty1);
374 w2x = 1./(stx2*stx2); w2y = 1./(sty2*sty2);
375
376 //printf("w1x %f w2x %f w1y %f w2y %f\n",w1x,w2x,w1y,w2y);
377 float chi2t = ( (s1.TX()-T)*(s1.TX()-T)*w1x +
378 s1.TY()*s1.TY() *w1y +
379 (s2.TX()-T)*(s2.TX()-T)*w2x +
380 s2.TY()*s2.TY() *w2y )/4.;
381 //printf("angular component of chi2 = %f\n",chi2t);
382
383 float sx1 = cond1.SigmaX(T), sy1 = cond1.SigmaY(0);
384 float sx2 = cond2.SigmaX(T), sy2 = cond2.SigmaY(0);
385
386 float dx1 = s1.X()-(s.X()+(s1.Z()-s.Z())*s.TX());
387 float dy1 = s1.Y()-(s.Y()+(s1.Z()-s.Z())*s.TY());
388 float dx2 = s2.X()-(s.X()+(s2.Z()-s.Z())*s.TX());
389 float dy2 = s2.Y()-(s.Y()+(s2.Z()-s.Z())*s.TY());
390
391 float chi2pos = dx1*dx1/sx1/sx1+dy1*dy1/sy1/sy1+dx2*dx2/sx2/sx2+dy2*dy2/sy2/sy2/4.;
392 //printf("position component of chi2 = %f\n",chi2pos);
393 s.SetChi2( Sqrt(chi2t+chi2pos) );
394 aff.Invert();
395 s.Transform(&aff);
396 //s.Print();
397
398 return s.Chi2();
399}
Definition: EdbAffine.h:17
void Invert()
Definition: EdbAffine.cxx:103
void Rotate(float angle)
Definition: EdbAffine.cxx:354
virtual void Transform(const EdbAffine2D *a)
Double_t X
Definition: tlg2couples.C:76
Double_t Y
Definition: tlg2couples.C:76
Double_t TY
Definition: tlg2couples.C:78
Double_t TX
Definition: tlg2couples.C:78
Double_t Z
Definition: tlg2couples.C:104

◆ Chi2SegMCS()

float EdbTrackFitter::Chi2SegMCS ( const EdbSegP s1,
const EdbSegP s2 
)

double DZemul=s1.DZem()+s2.DZem();
Use hardcoded value because EdbSegP::DZem() contains sum of chi2 for microtracks.

293 {
294
297 double DZemul=88;
298
299 double dzem=0.5*(s1.DZ()+s2.DZ())+DZemul;
300 double dz=TMath::Abs(s1.Z()-s2.Z());
301 double dzpb=dz-dzem;
302 if(dzpb<=0)return 0;
303 double dist=EdbSegP::Distance(s1,s2);
304 dist*=dzpb/dz;
305 double mom=s1.P();
306 if(mom<=0)mom=ePdef;
307 double theta0=EdbPhysics::ThetaMCS(mom,eM,dist,eX0);
308
309 double theta=TMath::Sqrt2()*EdbSegP::Angle(s1,s2);
310 double chi=theta/theta0;
311 Log(5,"Chi2SegMCS",Form("dist=%6.4f (%6.4f), p=%6.2f, th=%4.3g(%4.3g) => chi2=%5.4g",dist,eX0,mom,theta,theta0,chi*chi));
312 return chi*chi;
313}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
static double ThetaMCS(float p, float mass, float dx, float x0)
Definition: EdbPhys.cxx:45
Float_t DZ() const
Definition: EdbSegP.h:154
Float_t P() const
Definition: EdbSegP.h:152
static Float_t Distance(const EdbSegP &s1, const EdbSegP &s2)
Definition: EdbSegP.cxx:437
static Float_t Angle(const EdbSegP &s1, const EdbSegP &s2)
Definition: EdbSegP.cxx:444
float eX0
rad length of the media [microns]
Definition: EdbTrackFitter.h:23
float ePdef
default momentum
Definition: EdbTrackFitter.h:25
float eM
mass of the particle (if negative - use the mass setted in the track)
Definition: EdbTrackFitter.h:24

◆ FitTrackLine() [1/2]

int EdbTrackFitter::FitTrackLine ( const EdbTrackP tr,
float &  x,
float &  y,
float &  z,
float &  tx,
float &  ty,
float &  w 
)

track fit by averaging of segments parameters and return the mean values

548{
550
551 int nseg=tr.N();
552 x=0; y=0; z=0; tx=0; ty=0; w=0;
553 EdbSegP *seg=0;
554 for(int i=0; i<nseg; i++) {
555 seg = tr.GetSegment(i);
556 x += seg->X();
557 y += seg->Y();
558 z += seg->Z();
559 tx += seg->TX();
560 ty += seg->TY();
561 w += seg->W();
562 }
563 x /= nseg;
564 y /= nseg;
565 z /= nseg;
566 tx /= nseg;
567 ty /= nseg;
568 return nseg;
569}
Definition: EdbSegP.h:21
void w(int rid=2, int nviews=2)
Definition: test.C:27

◆ FitTrackLine() [2/2]

int EdbTrackFitter::FitTrackLine ( EdbTrackP tr)

track fit by averaging of segments parameters and put them as the track parameters

521{
523 float x,y,z,tx,ty,w;
524 FitTrackLine(tr,x,y,z,tx,ty,w);
525 tr.Set(tr.ID(),x,y,tx,ty,w,tr.Flag());
526 tr.SetZ(z);
527
528 tr.ClearF();
529 int nseg=tr.N();
530 for(int i=0; i<nseg; i++) {
531 EdbSegP *s = tr.GetSegment(i);
532 float dz = s->Z()-z;
533 EdbSegP *segf = new EdbSegP(*s);
534
535 segf->Set( s->ID(), x + tx*dz, y+ty*dz, tx, ty, 1.,s->Flag());
536 segf->SetErrorP ( tr.SP() );
537 //segf.SetChi2(0.);
538 //segf.SetProb(1.);
539 segf->SetP( tr.P() );
540 tr.AddSegmentF(segf);
541 }
542
543 return tr.N();
544}
void SetP(float p)
Definition: EdbSegP.h:133
void SetErrorP(float sp2)
Definition: EdbSegP.h:94
int FitTrackLine(EdbTrackP &tr)
Definition: EdbTrackFitter.cxx:520

◆ MaxChi2Seg()

float EdbTrackFitter::MaxChi2Seg ( EdbTrackP tr)
static

return the maximal seg-to-seg chi2 along the track

403{
405 float chimax=0,chi=0;
406 if(t.N()<2) return chimax;
407 for(int i=0; i<t.N()-1; i++) {
408 EdbSegP s(*t.GetSegment(i));
409 chi = Chi2Seg(&s,t.GetSegment(i+1));
410 chimax = chi > chimax ? chi : chimax;
411 }
412 return chimax;
413}
static float Chi2Seg(EdbSegP *s1, EdbSegP *s2)
Definition: EdbTrackFitter.cxx:62
TTree * t
Definition: check_shower.C:4

◆ MaxKink()

float EdbTrackFitter::MaxKink ( EdbTrackP tr)
static

return the maximal seg-to-seg kink theta angle

430{
432 float kink=0,theta=0;
433 if(t.N()<2) return kink;
434 for(int i=0; i<t.N()-1; i++) {
435 theta = Theta(*t.GetSegment(i),*t.GetSegment(i+1));
436 kink = theta > kink ? theta : kink;
437 }
438 return kink;
439}
static float Theta(EdbSegP &s, EdbSegP &s1)
Definition: EdbTrackFitter.cxx:452

◆ MeanChi2Seg()

float EdbTrackFitter::MeanChi2Seg ( EdbTrackP tr)
static
417{
418 float meanchi=0;
419 if(t.N()<2) return meanchi;
420 for(int i=0; i<t.N()-1; i++) {
421 EdbSegP s(*t.GetSegment(i));
422 meanchi += Chi2Seg(&s,t.GetSegment(i+1));
423 }
424 return meanchi /= (t.N()-1);
425}

◆ MeanKink()

float EdbTrackFitter::MeanKink ( EdbTrackP tr)
static
443{
444 float meankink=0;
445 if(t.N()<2) return 0;
446 for(int i=0; i<t.N()-1; i++)
447 meankink += Theta(*t.GetSegment(i),*t.GetSegment(i+1));
448 return meankink /= (t.N()-1);
449}

◆ PMS_KF()

float EdbTrackFitter::PMS_KF ( EdbTrackP t,
float  p0 = 10.,
float  probbest = 0.5 
)

select track momentum in a way to have the given chi2-probablity calculated by KF

498{
500
501 if(t.N()<2) return 0;
502 float prob=0;
503 float pu=100., pl=0., p=p0;
504 int nstep = 0;
505 while( Abs(prob-probbest)>0.001 ) {
506 nstep++;
507 t.SetP(p);
508 t.FitTrackKFS(true);
509 prob = t.Prob();
510 if(prob<probbest) pu=p;
511 else pl=p;
512 if(nstep>30) break;
513 p = (pu+pl)/2.;
514 }
515 if(nstep>20) printf("Warning in EdbTrackFitter::PMS_KF: nstep=%d nseg=%d p=%f prob=%f\n",nstep,t.N(),p,prob);
516 return t.P();
517}
p
Definition: testBGReduction_AllMethods.C:8

◆ Print()

void EdbTrackFitter::Print ( )
54{
55 printf("EdbTrackFitter seetings:\n");
56 printf("eX0 = %f\n",eX0);
57 printf("eM = %f\n",eM);
58 printf("\n");
59}

◆ ProbSegMCS()

double EdbTrackFitter::ProbSegMCS ( EdbSegP s1,
EdbSegP s2 
)
288 {
289 if(s1->P()<=0)s1->SetP(ePdef);
290 return EdbPVRec::ProbeSeg(s1,s2,eX0,eM);
291}
static double ProbeSeg(const EdbTrackP *s1, EdbTrackP *s2, const float X0=5810.)
Definition: EdbPVRec.cxx:2658

◆ RMSprojXY()

int EdbTrackFitter::RMSprojXY ( EdbTrackP tr,
float &  ex,
float &  ey 
)
static

fit track by straight line using coordinates only and return RMS's

573{
575 int nseg=tr.N();
576 float x[100], y[100], z[100], w[100];
577 float x0,y0,z0,tx,ty;
578 for(int i=0; i<nseg; i++) {
579 EdbSegP *s = tr.GetSegment(i);
580 x[i]=s->X();
581 y[i]=s->Y();
582 z[i]=s->Z();
583 w[i]=s->W();
584 }
585 return EdbMath::LFIT3(x,y,z,w,nseg,x0,y0,z0,tx,ty,ex,ey);
586}
brick z0
Definition: RecDispMC.C:106
static int LFIT3(float *X, float *Y, float *Z, float *W, int L, float &X0, float &Y0, float &Z0, float &TX, float &TY, float &EX, float &EY)
Definition: EdbMath.cxx:313

◆ SetDefaultBrick()

void EdbTrackFitter::SetDefaultBrick ( )
44{
46 eTPb = 1000./1300.;
47 ePcut = 0.050;
48 eM = 0.13957;
49 eDE_correction = false;
50}
static float kX0_Pb()
Definition: EdbPhys.h:23
float eTPb
?
Definition: EdbTrackFitter.h:26
bool eDE_correction
take into account the energy loss or not
Definition: EdbTrackFitter.h:28
float ePcut
minimal momentum?
Definition: EdbTrackFitter.h:27

◆ SetNsegMax()

void EdbTrackFitter::SetNsegMax ( int  nseg)
inline
34{eNsegMax=nseg;}
int eNsegMax
max number of segments (for arrays allocation)
Definition: EdbTrackFitter.h:20

◆ SplitTrack()

bool EdbTrackFitter::SplitTrack ( EdbTrackP t,
EdbTrackP t1,
int  isplit 
)

split track t in 2 at the point isplit - will be the first segment of of t1

459{
461 if(t.N()<isplit) return false;
462 for(int i=t.N()-1; i>=isplit; i--) {
463 t1.AddSegment( t.GetSegment(i) );
464 t.RemoveSegment( t.GetSegment(i) );
465 }
466 t.SetCounters();
467 //t1.FitTrackKFS(true, X0, 0);
468 t1.SetCounters();
469 t1.SetM(t.M());
470 t1.SetP(t.P());
471 //t1.FitTrackKFS(true, X0, 0);
472
473 return true;
474}
void AddSegment(EdbSegP *s)
Definition: EdbPattern.h:214
void SetCounters()
Definition: EdbPattern.h:159
void SetM(float m)
Definition: EdbPattern.h:154

◆ SplitTrackByKink()

int EdbTrackFitter::SplitTrackByKink ( EdbTrackP t,
TObjArray &  tracks,
float  maxkink 
)

split track t in several tracks accourding to maxkink
return total number of tracks after splitting;
after splitting all new tracks added to the array "tracks"

478{
482
483 if(t->N()<1) return 0;
484 if(t->N()<2) return 1;
485 int nsplit=1;
486 for(int i=t->N()-1; i>0; i--)
487 if( Theta(*t->GetSegment(i),*t->GetSegment(i-1)) >= maxkink ) {
488 EdbTrackP *t1 = new EdbTrackP();
489 SplitTrack(*t,*t1, i);
490 tracks.Add(t1);
491 nsplit++;
492 }
493 return nsplit;;
494}
bool SplitTrack(EdbTrackP &t, EdbTrackP &t1, int isplit)
Definition: EdbTrackFitter.cxx:458
Definition: EdbPattern.h:113
TTree * tracks
Definition: check_tr.C:19

◆ Theta()

float EdbTrackFitter::Theta ( EdbSegP s,
EdbSegP s1 
)
static
453{
454 return Sqrt( (s.TX()-s1.TX())*(s.TX()-s1.TX()) + (s.TY()-s1.TY())*(s.TY()-s1.TY()) );
455}

Member Data Documentation

◆ eDE_correction

bool EdbTrackFitter::eDE_correction

take into account the energy loss or not

◆ eM

float EdbTrackFitter::eM

mass of the particle (if negative - use the mass setted in the track)

◆ eNsegMax

int EdbTrackFitter::eNsegMax
private

max number of segments (for arrays allocation)

◆ ePcut

float EdbTrackFitter::ePcut

minimal momentum?

◆ ePdef

float EdbTrackFitter::ePdef

default momentum

◆ eTPb

float EdbTrackFitter::eTPb

?

◆ eX0

float EdbTrackFitter::eX0

rad length of the media [microns]


The documentation for this class was generated from the following files: