FEDRA emulsion software from the OPERA Collaboration
EdbVertex Class Reference

#include <EdbVertex.h>

Inheritance diagram for EdbVertex:
Collaboration diagram for EdbVertex:

Public Member Functions

void AddVTA (EdbVTA *vta)
 
Int_t CheckDiscardedTracks ()
 
float CheckImp (const EdbTrackP *tr)
 
EdbVTACheckImp (const EdbTrackP *tr, float ImpMax, int zpos, float dist)
 
float CheckImpGeom (const EdbTrackP *tr)
 
Float_t Chi2Track (EdbTrackP *tr, int zpos, float X0=0.)
 
void Clear ()
 
void ClearNeighborhood ()
 
void ClearV ()
 
Int_t Compare (const TObject *o) const
 
Float_t DistSeg (EdbSegP *seg, float X0=0.)
 
Float_t DistTrack (EdbTrackP *tr, int zpos, float X0=0.)
 
void Edb2Vt (const EdbSegP &s, VERTEX::Track &t, float X0=0., float m=0.139)
 
void Edb2Vt (const EdbTrackP &tr, VERTEX::Track &t, float X0=0., float m=0.139)
 
 EdbVertex ()
 
Int_t EstimateVertexFlag ()
 
Bool_t EstimateVertexMath (float &xv, float &yv, float &zv, float &d)
 
Int_t Flag () const
 
EdbVertexGetConnectedVertex (int nv)
 
EdbVertexGetConnectedVertexForTrack (int it)
 
EdbVTAGetMaxImpVTA ()
 
EdbTrackPGetTrack (int i)
 
EdbTrackPGetTrackN (int i)
 
EdbSegPGetTrackV (int i, bool usesegpar=false)
 
EdbVTAGetVTa (int i)
 
EdbVTAGetVTn (int i)
 
ULong_t Hash () const
 
Int_t ID () const
 
Float_t Impact (int i)
 
Float_t ImpTrack (int i)
 
Bool_t IsEqual (const TObject *o) const
 
Bool_t IsSortable () const
 
Float_t MaxAperture ()
 
Float_t MaxImpact ()
 
Int_t MCEvt () const
 
EdbTrackPMeanTrack ()
 
Float_t MinDist ()
 
Int_t N () const
 
Int_t Nn () const
 
Int_t Nv ()
 
void Print ()
 
Float_t Quality ()
 
void RemoveVTA (EdbVTA *vta)
 
void ResetTracks ()
 
void SetFlag (int flag=0)
 
void SetID (int ID=0)
 
void SetMC (int mEvt=0)
 
void SetQuality (float q=0)
 
void SetV (VERTEX::Vertex *v)
 
void SetXYZ (float x, float y, float z)
 
Bool_t TrackInVertex (EdbTrackP *t)
 
VERTEX::VertexV () const
 
Float_t Volume ()
 
TList * VTa ()
 
TList * VTn ()
 
Float_t VX () const
 
Float_t VY () const
 
Float_t VZ () const
 
Float_t X () const
 
Float_t Y () const
 
Float_t Z () const
 
Int_t Zpos (int i)
 
virtual ~EdbVertex ()
 

Private Attributes

Int_t eFlag
 
Int_t eID
 
Int_t eMCEvt
 
Float_t eQuality
 Probability/(vsigmax**2+vsigmay**2) More...
 
VERTEX::VertexeV
 pointer to VtVertex object More...
 
TList eVTa
 attached tracks More...
 
TList eVTn
 vertex neighborhood tracks and segments More...
 
Float_t eX
 for generated vertexes - the real vertex position More...
 
Float_t eY
 for reconstructed ones - average of track connection
More...
 
Float_t eZ
 to avoid accuracy problem More...
 

Constructor & Destructor Documentation

◆ EdbVertex()

EdbVertex::EdbVertex ( )
97{
98 eID= 0;
99 eV = 0;
100 eX = 0.;
101 eY = 0.;
102 eZ = 0.;
103 eFlag = 0;
104 eQuality=0.;
105 eMCEvt=-999;
106}
Float_t eY
for reconstructed ones - average of track connection
Definition: EdbVertex.h:77
Int_t eID
Definition: EdbVertex.h:88
Float_t eX
for generated vertexes - the real vertex position
Definition: EdbVertex.h:76
VERTEX::Vertex * eV
pointer to VtVertex object
Definition: EdbVertex.h:91
Float_t eQuality
Probability/(vsigmax**2+vsigmay**2)
Definition: EdbVertex.h:89
Float_t eZ
to avoid accuracy problem
Definition: EdbVertex.h:78
Int_t eFlag
Definition: EdbVertex.h:81
Int_t eMCEvt
Definition: EdbVertex.h:87

◆ ~EdbVertex()

EdbVertex::~EdbVertex ( )
virtual
110{
111 Clear();
112}
void Clear()
Definition: EdbVertex.cxx:125

Member Function Documentation

◆ AddVTA()

void EdbVertex::AddVTA ( EdbVTA vta)
356{
357 if(vta->Flag()!=2) eVTn.Add(vta);
358 else eVTa.Add(vta);
359}
Int_t Flag() const
Definition: EdbVertex.h:48
TList eVTa
attached tracks
Definition: EdbVertex.h:74
TList eVTn
vertex neighborhood tracks and segments
Definition: EdbVertex.h:73

◆ CheckDiscardedTracks()

Int_t EdbVertex::CheckDiscardedTracks ( )
169{
170 int ndsc = 0;
171 EdbVTA *vta = 0;
172 for (int i=0; i<N(); i++)
173 {
174 vta = GetVTa(i);
175 if (vta) {
176 if( GetTrack(i)->Vertex(vta->Zpos()) != this ) ndsc++;
177 }
178 else ndsc++;
179 }
180 return ndsc;
181}
Definition: EdbVertex.h:26
Int_t Zpos() const
Definition: EdbVertex.h:47
EdbVTA * GetVTa(int i)
Definition: EdbVertex.h:139
Int_t N() const
Definition: EdbVertex.h:121
EdbTrackP * GetTrack(int i)
Definition: EdbVertex.h:141
Definition: VtVertex.hh:88

◆ CheckImp() [1/2]

float EdbVertex::CheckImp ( const EdbTrackP tr)
442{
443 Track *t = new Track();
444 Vertex *v = this->V();
445 Edb2Vt(*tr, *t);
446 return distance(*t,*v);
447}
VERTEX::Vertex * V() const
Definition: EdbVertex.h:154
void Edb2Vt(const EdbTrackP &tr, VERTEX::Track &t, float X0=0., float m=0.139)
Definition: EdbVertex.cxx:568
Definition: VtTrack.hh:64
TTree * t
Definition: check_shower.C:4
@ Track
Definition: tlg2couples.C:52
double distance(const Track &t, const Vertex &v)
spatial distance track - vertex
Definition: VtDistance.C:49

◆ CheckImp() [2/2]

EdbVTA * EdbVertex::CheckImp ( const EdbTrackP tr,
float  ImpMax,
int  zpos,
float  dist 
)
451{
452 EdbVTA *vta = 0;
453 EdbTrackP *tr1 = (EdbTrackP *)tr;
454 if (!tr) return vta;
455 Track *t = new Track();
456 Vertex *v = this->V();
457 Edb2Vt(*tr, *t);
458 float imp = distance(*t,*v);
459 if (imp > ImpMax) { delete t; return vta;}
460 vta = new EdbVTA(tr1, this);
461 vta->SetZpos(zpos);
462 vta->SetFlag(0);
463 vta->SetImp(imp);
464 vta->SetDist(dist);
465 AddVTA(vta);
466 delete t;
467 return vta;
468}
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
Definition: EdbPattern.h:113
void SetImp(float imp)
Definition: EdbVertex.h:57
void SetDist(float dist)
Definition: EdbVertex.h:58
void SetZpos(int zpos)
Definition: EdbVertex.h:55
void SetFlag(int flag)
Definition: EdbVertex.h:56
void AddVTA(EdbVTA *vta)
Definition: EdbVertex.cxx:355
float ImpMax
Definition: check_vertex.C:30

◆ CheckImpGeom()

float EdbVertex::CheckImpGeom ( const EdbTrackP tr)
431{
432 float pv[3] = {VX(), VY(), VZ()};
433 float p1[3] = { tr->X(), tr->Y(), tr->Z() };
434 float p2[3] = { tr->X() + tr->TX()*1000., tr->Y() + tr->TY()*1000., tr->Z()+1000. };
435 bool inside=0;
436 float imp = EdbMath::DistancePointLine3(pv, p1,p2, &inside);
437 return imp;
438}
static double DistancePointLine3(float Point[3], float LineStart[3], float LineEnd[3], bool inside)
Definition: EdbMath.cxx:61
Float_t VX() const
Definition: EdbVertex.h:133
Float_t VY() const
Definition: EdbVertex.h:134
Float_t VZ() const
Definition: EdbVertex.h:135

◆ Chi2Track()

float EdbVertex::Chi2Track ( EdbTrackP tr,
int  zpos,
float  X0 = 0. 
)

Chi2-distance from track to already existing vertex

636{
638
639 if (!track) return 0.;
640 double distchi2 = -2.;
641 EdbSegP *seg = track->TrackExtremity( zpos); // usesegpar??
642 if (eV)
643 {
644 if (track->NF() <= 0) return -1.;
645 if (track->P() <= 0.) track->SetP(1.);
646 if (track->M() <= 0.) track->SetM(.1395);
647 if (track->SP() <= 0.) track->SetErrorP(1.);
648 Track *t=new Track();
649 Edb2Vt( *seg, *t, X0, track->M() );
650 distchi2 = -3.;
651 if (eV->valid())
652 {
653 t->rm(track->M());
654 distchi2 = eV->distance(*t);
655 }
656 SafeDelete(t);
657 }
658 return (float)distchi2;
659}
brick X0
Definition: RecDispMC.C:112
Definition: EdbSegP.h:21
double distance(double x, double y, double z) const
$\chi^2$ distance to space point, $ndf = 3$
Definition: VtVertex.C:803
bool valid() const
is vertex valid?
Definition: bitview.h:14

◆ Clear()

void EdbVertex::Clear ( )
126{
127 eVTa.Delete("slow");
128 eVTn.Delete("slow");
129 ClearV();
130 eX = 0.;
131 eY = 0.;
132 eZ = 0.;
133 eID = 0;
134 eFlag = 0;
135 eQuality = 0.;
136 eMCEvt=-999;
137}
void ClearV()
Definition: EdbVertex.cxx:115

◆ ClearNeighborhood()

void EdbVertex::ClearNeighborhood ( )
276{
277 int nn = Nn();
278 if (nn>0) eVTn.Clear("nodelete");
279 for(int i=0; i<nn; i++) delete GetVTn(i);
280}
EdbVTA * GetVTn(int i)
Definition: EdbVertex.h:140
Int_t Nn() const
Definition: EdbVertex.h:122

◆ ClearV()

void EdbVertex::ClearV ( )

clear VtVertex object

116{
118 if (eV) {
119 eV->clear(); // should delete also tracks ownered by vertex
120 SafeDelete(eV);
121 }
122}
void clear()
clear all track-vertex relations, makes vertex invalid
Definition: VtVertex.C:196

◆ Compare()

int EdbVertex::Compare ( const TObject *  o) const
283{
284 /*printf("Inside compare\n");*/
285 if ( eQuality > ((EdbVertex *)o)->eQuality ) return -1;
286 else if ( eQuality == ((EdbVertex *)o)->eQuality ) return 0;
287 else return 1;
288}
Definition: EdbVertex.h:69
AcqOdyssey * o
Definition: hwinit.C:2

◆ DistSeg()

float EdbVertex::DistSeg ( EdbSegP seg,
float  X0 = 0. 
)

distance from segment to already fitted vertex

673{
675 if(seg)
676 if(eV)
677 if(eV->valid()) {
678 Track t;
679 Edb2Vt( *seg, t, X0, 0. );
680 return (float)distance(t, *eV);
681 }
682 return 100000.;
683}

◆ DistTrack()

float EdbVertex::DistTrack ( EdbTrackP tr,
int  zpos,
float  X0 = 0. 
)

distance from track to already fitted vertex

663{
665 if (!track) return 0.;
666 EdbSegP *seg = track->TrackExtremity( zpos); // usesegpar??
667 if (!seg) return 0.;
668 return DistSeg(seg,X0);
669}
Float_t DistSeg(EdbSegP *seg, float X0=0.)
Definition: EdbVertex.cxx:672

◆ Edb2Vt() [1/2]

void EdbVertex::Edb2Vt ( const EdbSegP s,
VERTEX::Track t,
float  X0 = 0.,
float  m = 0.139 
)

Input: EdbSegP tr - track parameters near vertex
X0 - rad length for ms estimation
m - mass of the particle
if X0 or m are negative - ignore multiple scattering
Output: VERTEX:Track t - object propagated to the vertex position with the estimated errors matrix
Used: eX,eY,eZ - the reference point of this vertex

575{
582
583 double dz = eZ - tr.Z();
584 double tx = (double)tr.TX();
585 double ty = (double)tr.TY();
586 double x = (double)tr.X() + tx*dz - eX;
587 double y = (double)tr.Y() + ty*dz - eY;
588 double z = 0.;
589 float p = tr.P();
590
591 VtSymMatrix dms(4); // multiple scattering matrix
592 dms.clear();
593 if ( X0 > 0. && m > 0.)
594 {
595 double dPb = dz*TMath::Sqrt(1.+tx*tx+ty*ty); // thickness of the Pb+emulsion cell in microns
596 double theta0sq = EdbPhysics::ThetaMS2( p, m, dPb, X0 );
597 //printf("( p, m, dPb, X0 dz) = %f %f %f %f %f theta0 = %g\n", p, m, dPb, X0, dz, TMath::Sqrt(theta0sq));
598 dms(0,0) = theta0sq*dz*dz/3.;
599 dms(1,1) = dms(0,0);
600 dms(2,2) = theta0sq;
601 dms(3,3) = dms(2,2);
602 dms(2,0) = theta0sq*dz/2.;
603 dms(3,1) = dms(2,0);
604 dms(0,2) = dms(2,0);
605 dms(1,3) = dms(2,0);
606 }
607
608 VtSqMatrix pred(4); //propagation matrix for track parameters (x,y,tx,ty)
609 pred.clear();
610 pred(0,0) = 1.;
611 pred(1,1) = 1.;
612 pred(2,2) = 1.;
613 pred(3,3) = 1.;
614 pred(0,2) = dz;
615 pred(1,3) = dz;
616
617 VtSymMatrix cov(4); // covariance matrix for seg0
618 for(int k=0; k<4; k++)
619 for(int l=0; l<4; l++) cov(k,l) = (tr.COV())(k,l);
620
621 VtSymMatrix covpred(4); // covariance matrix for prediction
622 covpred = pred*(cov*(pred.T()))+dms;
623
624 CMatrix covp; // covariance matrix for the track
625 covp.clear();
626 for(int k=0; k<4; k++)
627 for(int l=0; l<4; l++) covp(k,l) = covpred(k,l);
628 covp(4,4) = (tr.COV())(4,4);
629
630 t.set(x, y, z, tx, ty, (double)p, covp);
631 t.rm((double)m);
632}
brick dz
Definition: RecDispMC.C:107
static double ThetaMS2(float p, float mass, float dx, float X0)
Definition: EdbPhys.cxx:51
Definition: CMatrix.hh:63
void clear(void)
set matrix elements to 0
Definition: VtMatrix.C:394
Definition: VtSqMatrix.hh:50
Definition: VtSymMatrix.hh:49
p
Definition: testBGReduction_AllMethods.C:8

◆ Edb2Vt() [2/2]

void EdbVertex::Edb2Vt ( const EdbTrackP tr,
VERTEX::Track t,
float  X0 = 0.,
float  m = 0.139 
)
569{
570 Edb2Vt( *((EdbSegP*)&tr), t, X0, m);
571}

◆ EstimateVertexFlag()

Int_t EdbVertex::EstimateVertexFlag ( )
185{
186 int flag=-1;
187 int n0 = 0, n1=0, nn=0;
188 EdbVTA *vta = 0;
189 for (int i=0; i<N(); i++)
190 {
191 vta = GetVTa(i);
192 if ( vta->Zpos()==0 ) n0++;
193 else if( vta->Zpos()==1 ) n1++;
194 else nn++;
195 }
196
197 if ( n0>0 && n1 >0) flag=1; // end & start
198 else if ( n0==0 && n1>0) flag=0; // start & start
199 else if ( n0>0 && n1==0) flag=2; // end & end
200 else flag= -1;
201 if(nn) flag= -1; // flag -1 should newer happened: if so - should be debugged
202
203 SetFlag(flag);
204 Log(3,"EdbVertex::EstimateVertexFlag","%d with n0=%d n1=%d, nn=%d", flag, n0,n1,nn );
205 return flag;
206}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
void SetFlag(int flag=0)
Definition: EdbVertex.h:158

◆ EstimateVertexMath()

bool EdbVertex::EstimateVertexMath ( float &  xv,
float &  yv,
float &  zv,
float &  d 
)
472{
473 int nt = N();
474 if(!nt) return false;
475
476 double tx_sum = 0.;
477 double x_sum = 0.;
478 double xw_sum = 0.;
479 double xtx_sum = 0.;
480 double tx2_sum = 0.;
481 double ty_sum = 0.;
482 double y_sum = 0.;
483 double yw_sum = 0.;
484 double yty_sum = 0.;
485 double ty2_sum = 0.;
486 const EdbTrackP *tr = 0;
487 const EdbSegP *seg = 0;
488
489 double x,y,tx,ty,xweight,yweight,xweight2,yweight2;
490
491 // fill cumulants
492 for( int i = 0; i < nt; i++ ) {
493
494 tr = GetTrack(i);
495
496 seg = tr->TrackExtremity( Zpos(i)); // usesegpar??
497
498 x = seg->X();
499 y = seg->Y();
500 tx = seg->TX();
501 ty = seg->TY();
502 xweight = 1./(seg->COV())(0,0);
503 yweight = 1./(seg->COV())(1,1);
504 xweight2 = xweight*xweight;
505 yweight2 = yweight*yweight;
506
507 tx_sum += tx * xweight;
508 x_sum += x * xweight;
509 xw_sum += xweight;
510 xtx_sum += x * tx * xweight2;
511 tx2_sum += tx * tx * xweight2;
512
513 ty_sum += ty * yweight;
514 y_sum += y * yweight;
515 yw_sum += yweight;
516 yty_sum += y * ty * yweight2;
517 ty2_sum += ty * ty * yweight2;
518
519 } // for track
520
521 double det = -tx2_sum - ty2_sum + tx_sum*tx_sum/xw_sum + ty_sum*ty_sum/yw_sum;
522
523 if(det == 0.) {
524 return false;
525 }
526
527 zv = ( xtx_sum + yty_sum - tx_sum*x_sum/xw_sum - ty_sum*y_sum/yw_sum ) / det;
528 xv = ( x_sum + tx_sum * zv ) / xw_sum;
529 yv = ( y_sum + ty_sum * zv ) / yw_sum;
530
531
532// float zTolerance=300.;
533
534// for( int i = 0; i < nt; i++ ) {
535// tr = GetTrack(i);
536// if (Zpos(i)) seg = tr->TrackZmin();
537// else seg = tr->TrackZmax();
538// if( zv > (seg->Z() + zTolerance) ) return false;
539// if( zv < (seg->Z() - zTolerance) ) return false;
540// }
541
542 double drx;
543 double dry;
544 double drz;
545 double drt;
546 double drms = 0.;
547
548 for( int i = 0; i < nt; i++ ) {
549
550 tr = GetTrack(i);
551
552 seg = tr->TrackExtremity( Zpos(i)); // usesegpar??
553
554 drx = seg->X() - xv;
555 dry = seg->Y() - yv;
556 drz = seg->Z() - zv;
557 drt = (drx*seg->TX() + dry*seg->TY() + drz);
558 drms += drx*drx + dry*dry + drz*drz -
559 (drt*drt)/(1.+seg->TX()*seg->TX()+seg->TY()*seg->TY());
560 }
561
562 d = TMath::Sqrt(drms/nt);
563
564 return true;
565}
void d()
Definition: RecDispEX.C:381
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t X() const
Definition: EdbSegP.h:173
TMatrixD & COV() const
Definition: EdbSegP.h:123
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
Int_t Zpos(int i)
Definition: EdbVertex.h:127

◆ Flag()

Int_t EdbVertex::Flag ( ) const
inline
124{return eFlag;}

◆ GetConnectedVertex()

EdbVertex * EdbVertex::GetConnectedVertex ( int  nv)
406{
407 EdbTrackP *tr = 0;
408 EdbVertex *vc = 0;
409 int n = 0;
410 for (int i=0; i<N(); i++)
411 {
412 if ((tr = GetTrack(i)))
413 {
414 if ((Zpos(i) == 1) && (vc = tr->VertexE()))
415 {
416 if (n == nv) return vc;
417 n++;
418 }
419 else if ((Zpos(i) == 0) && (vc = tr->VertexS()))
420 {
421 if (n == nv) return vc;
422 n++;
423 }
424 }
425 }
426 return 0;
427}

◆ GetConnectedVertexForTrack()

EdbVertex * EdbVertex::GetConnectedVertexForTrack ( int  it)
385{
386 EdbTrackP *tr = 0;
387 EdbVertex *vc = 0;
388 if (it<N())
389 {
390 if ((tr = GetTrack(it)))
391 {
392 if ((Zpos(it) == 1) && (vc = tr->VertexE()))
393 {
394 return vc;
395 }
396 else if ((Zpos(it) == 0) && (vc = tr->VertexS()))
397 {
398 return vc;
399 }
400 }
401 }
402 return 0;
403}

◆ GetMaxImpVTA()

EdbVTA * EdbVertex::GetMaxImpVTA ( )

return vta with a biggest impact par

141{
143 int ntr = N();
144 EdbVTA *maxvta=0;
145 if(ntr<1) maxvta=0;
146 else if(ntr==1) maxvta=GetVTa(0);
147 else {
148 float maximp=0;
149 for(int i=0; i<ntr; i++) {
150 EdbVTA *vta = GetVTa(i);
151 float imp = vta->Imp();
152 if( imp > maximp ) { maximp=imp; maxvta=vta; }
153 }
154 }
155 return maxvta;
156}
Float_t Imp() const
Definition: EdbVertex.h:49

◆ GetTrack()

EdbTrackP * EdbVertex::GetTrack ( int  i)
inline
141{return GetVTa(i)->GetTrack();}
EdbTrackP * GetTrack() const
Definition: EdbVertex.h:51

◆ GetTrackN()

EdbTrackP * EdbVertex::GetTrackN ( int  i)
inline
142{return GetVTn(i)->GetTrack();}

◆ GetTrackV()

EdbSegP * EdbVertex::GetTrackV ( int  i,
bool  usesegpar = false 
)
160{
161 EdbTrackP *t = GetTrack(i); if(!t) return 0;
162 EdbSegP *s = t->TrackExtremity(Zpos(i), usesegpar);
163 if( s->P()>=0 && (s->P() != t->P()) ) Log(1,"GetTrackV","Warning! segment momentum=%f is not equal to the track momentum=%f",s->P(), t->P() );
164 return s;
165}
s
Definition: check_shower.C:55

◆ GetVTa()

EdbVTA * EdbVertex::GetVTa ( int  i)
inline
139{return (EdbVTA*)(eVTa.At(i));}

◆ GetVTn()

EdbVTA * EdbVertex::GetVTn ( int  i)
inline
140{return (EdbVTA*)(eVTn.At(i));}

◆ Hash()

ULong_t EdbVertex::Hash ( ) const
inline
128{return eID;}

◆ ID()

Int_t EdbVertex::ID ( ) const
inline
126{return eID;}

◆ Impact()

Float_t EdbVertex::Impact ( int  i)
inline
149{ return GetVTa(i)? GetVTa(i)->Imp(): 1000000.; }

◆ ImpTrack()

Float_t EdbVertex::ImpTrack ( int  i)
inline
152{ return DistSeg( GetTrackV(i) ); }
EdbSegP * GetTrackV(int i, bool usesegpar=false)
Definition: EdbVertex.cxx:159

◆ IsEqual()

bool EdbVertex::IsEqual ( const TObject *  o) const
291{
292 /*printf("Inside isequal\n");*/
293 if ( eID != ((EdbVertex *)o)->eID ) return false;
294 if ( eQuality != ((EdbVertex *)o)->eQuality ) return false;
295 else return true;
296}

◆ IsSortable()

Bool_t EdbVertex::IsSortable ( ) const
inline
129{return kTRUE;}

◆ MaxAperture()

float EdbVertex::MaxAperture ( )
253{
254 float aper=0.;
255 int ntr = N();
256 if(ntr<2) return aper;
257 EdbTrackP *t1=0;
258 EdbTrackP *t2=0;
259
260 float tx=0,ty=0,a=0;
261 for (int i=0; i<ntr-1; i++) {
262 t1 = GetTrack(i);
263 for (int j=i+1; j<ntr; j++) {
264 t2 = GetTrack(j);
265
266 tx= t1->TX() - t2->TX();
267 ty= t1->TY() - t2->TY();
268 a = TMath::Sqrt( tx*tx+ty*ty );
269 if( a>aper) aper=a;
270 }
271 }
272 return aper;
273}
void a()
Definition: check_aligned.C:59

◆ MaxImpact()

Float_t EdbVertex::MaxImpact ( )
inline
116{ EdbVTA *vta=GetMaxImpVTA(); return vta? vta->Imp(): 0; }
EdbVTA * GetMaxImpVTA()
Definition: EdbVertex.cxx:140

◆ MCEvt()

Int_t EdbVertex::MCEvt ( ) const
inline
125{return eMCEvt;}

◆ MeanTrack()

EdbTrackP * EdbVertex::MeanTrack ( )

calculate mean track trajectory waited with momentum

300{
302 int ntr = N();
303 EdbTrackP *mean = new EdbTrackP();
304 float psum=0;
305 for(int i=0; i<ntr; i++) {
306 EdbTrackP *t = GetTrack(i);
307 mean->Set( 0 , mean->X() + t->X()*t->P()
308 , mean->Y() + t->Y()*t->P()
309 , mean->TX() + t->TX()*t->P()
310 , mean->TY() + t->TY()*t->P()
311 , mean->W() + t->W()*t->P()
312 , 0 );
313 mean->SetZ( mean->Z() + t->Z()*t->P() );
314 psum += t->P();
315 }
316 mean->Set( 0 , mean->X() / psum
317 , mean->Y() / psum
318 , mean->TX() / psum
319 , mean->TY() / psum
320 , mean->W() / psum
321 , 0 );
322 mean->SetZ( mean->Z() / psum );
323 mean->SetP(psum);
324 return mean;
325}
void SetZ(float z)
Definition: EdbSegP.h:125
void SetP(float p)
Definition: EdbSegP.h:133
Float_t W() const
Definition: EdbSegP.h:151
void Set(int id, float x, float y, float tx, float ty, float w, int flag)
Definition: EdbSegP.h:87

◆ MinDist()

Float_t EdbVertex::MinDist ( )
240{
241 float mind=999999999999.;
242 int ntr = N();
243 if(ntr<2) return mind;
244 for (int i=0; i<ntr; i++) {
245 float d = GetVTa(i)->Dist();
246 if(Abs(d)<mind) mind = Abs(d);
247 }
248 return mind;
249}
Float_t Dist() const
Definition: EdbVertex.h:50

◆ N()

Int_t EdbVertex::N ( ) const
inline
121{return eVTa.GetSize();}

◆ Nn()

Int_t EdbVertex::Nn ( ) const
inline
122{return eVTn.GetSize();}

◆ Nv()

int EdbVertex::Nv ( )

return the number of linked vertex: tracks attached by the other edge to the different vertex

368{
370 EdbTrackP *tr = 0;
371 EdbVertex *vc = 0;
372 int nv = 0;
373 for (int i=0; i<N(); i++)
374 {
375 if ((tr = GetTrack(i))) {
376 if ((Zpos(i) == 1) && (vc = tr->VertexE())) nv++;
377 else if ((Zpos(i) == 0) && (vc = tr->VertexS())) nv++;
378 }
379 }
380 return nv;
381}

◆ Print()

void EdbVertex::Print ( )
329{
330 int ntr = N();
331 printf( "\n********************** Vertex %d with flag %d and %d tracks **************************\n", ID(), Flag(), ntr );
332 printf( "Fit quality : probability = %f Chi2 = %f\n", V()->prob(), V()->chi2() );
333 printf( "Vertex Position : %12.3f %12.3f %12.3f \n", VX(), VY(), VZ() );
334 printf( "Position Errors : %12.3f %12.3f %12.3f\n", V()->vxerr(), V()->vyerr(), V()->vzerr() );
335 //printf("---------------------------------------------------------\n");
336 printf( "Track ID Nseg Mass P Chi2/ndf Prob Chi2Contrib Impact Ztr\n");
337 for(int i=0; i<ntr; i++) {
338 EdbTrackP *tr = GetTrack(i);
339 float ztr = tr->TrackExtremity( Zpos(i) )->Z();
340 printf("%4d %4d %4d %7.4f %7.2f %5.2f %7.4f %7.3f %8.2f %10.2f\n",
341 i, tr->ID(), tr->N(), tr->M(), tr->P(),
342 tr->Chi2()/tr->N(), tr->Prob(), V()->track_chi2(i), Impact(i), ztr );
343 }
344 /*
345 printf("------- mean \"jet\" direction: ---------\n");
346 EdbTrackP *mean = MeanTrack();
347 printf(" X Y Z TX TY P\n");
348 printf("%12.1f %12.1f %12.1f %8.4f %8.4f %10.2f\n",
349 mean->X(), mean->Y(), mean->Z(), mean->TX(), mean->TY(), mean->P() );
350 SafeDelete(mean);
351 */
352 printf("****************************************************************************************\n");
353}
Int_t ID() const
Definition: EdbVertex.h:126
Float_t Impact(int i)
Definition: EdbVertex.h:149
Int_t Flag() const
Definition: EdbVertex.h:124
Float_t chi2
Definition: testBGReduction_By_ANN.C:14

◆ Quality()

Float_t EdbVertex::Quality ( )
inline
136{return eQuality;}

◆ RemoveVTA()

void EdbVertex::RemoveVTA ( EdbVTA vta)
362{
363 if(vta->Flag()!=2) eVTn.Remove(vta);
364 else eVTa.Remove(vta);
365}

◆ ResetTracks()

void EdbVertex::ResetTracks ( )

Assign the eVTAS or eVTAE of the tracks to the current vertex

210{
212 EdbVTA *vta = 0;
213 for (int i=0; i<N(); i++)
214 {
215 vta = GetVTa(i);
216 if (vta) GetTrack(i)->AddVTA(vta);
217 }
218}
void AddVTA(EdbVTA *vta)
Definition: EdbPattern.cxx:450

◆ SetFlag()

void EdbVertex::SetFlag ( int  flag = 0)
inline
158{eFlag = flag;}

◆ SetID()

void EdbVertex::SetID ( int  ID = 0)
inline
156{eID = ID;}

◆ SetMC()

void EdbVertex::SetMC ( int  mEvt = 0)
inline
159{eMCEvt=mEvt;}

◆ SetQuality()

void EdbVertex::SetQuality ( float  q = 0)
inline
161{eQuality = q;}
q
Definition: testBGReduction_AllMethods.C:55

◆ SetV()

void EdbVertex::SetV ( VERTEX::Vertex v)
inline
160{eV=v;}

◆ SetXYZ()

void EdbVertex::SetXYZ ( float  x,
float  y,
float  z 
)
inline
157{eX=x; eY=y; eZ=z;}

◆ TrackInVertex()

Bool_t EdbVertex::TrackInVertex ( EdbTrackP t)
222{
223 int ntr = N();
224 if(!ntr) return 0;
225 for (int i=0; i<ntr; i++) if(t == GetTrack(i)) return 1;
226 return 0;
227}

◆ V()

VERTEX::Vertex * EdbVertex::V ( ) const
inline
154{return eV;}

◆ Volume()

Float_t EdbVertex::Volume ( )
inline
114{ if(V()) return V()->vxerr()*V()->vyerr()*V()->vzerr(); else return 0; }
double vyerr() const
$\sqrt{\sigma_{vy}^2}$ vertex $y$-error
double vzerr() const
$\sqrt{\sigma_{vz}^2}$ vertex $z$-error
double vxerr() const
$\sqrt{\sigma_{vx}^2}$ vertex $x$-error

◆ VTa()

TList * EdbVertex::VTa ( )
inline
137{return &eVTa;}

◆ VTn()

TList * EdbVertex::VTn ( )
inline
138{return &eVTn;}

◆ VX()

Float_t EdbVertex::VX ( ) const
inline
133{return eV ? (eV->vx() + eX) : 1000000.;}
float vx() const
$x$ of vertex
Definition: VtVertex.C:233

◆ VY()

Float_t EdbVertex::VY ( ) const
inline
134{return eV ? (eV->vy() + eY) : 1000000.;}
float vy() const
$y$ of vertex
Definition: VtVertex.C:234

◆ VZ()

Float_t EdbVertex::VZ ( ) const
inline
135{return eV ? (eV->vz() + eZ) : 1000000.;}
float vz() const
$z$ of vertex
Definition: VtVertex.C:235

◆ X()

Float_t EdbVertex::X ( ) const
inline
130{return eX;}

◆ Y()

Float_t EdbVertex::Y ( ) const
inline
131{return eY;}

◆ Z()

Float_t EdbVertex::Z ( ) const
inline
132{return eZ;}

◆ Zpos()

Int_t EdbVertex::Zpos ( int  i)
inline
127{return GetVTa(i)->Zpos();}

Member Data Documentation

◆ eFlag

Int_t EdbVertex::eFlag
private

0 - neutral (tracks starts attached only) 1 - charge (tracks ends&starts attached) 2 - back neutral (tracks ends attached only) 3 - neutral, linked (has common track with other vertex) 4 - charge, linked 5 - back neutral, linked

◆ eID

Int_t EdbVertex::eID
private

◆ eMCEvt

Int_t EdbVertex::eMCEvt
private

◆ eQuality

Float_t EdbVertex::eQuality
private

Probability/(vsigmax**2+vsigmay**2)

◆ eV

VERTEX::Vertex* EdbVertex::eV
private

pointer to VtVertex object

◆ eVTa

TList EdbVertex::eVTa
private

attached tracks

◆ eVTn

TList EdbVertex::eVTn
private

vertex neighborhood tracks and segments

◆ eX

Float_t EdbVertex::eX
private

for generated vertexes - the real vertex position

◆ eY

Float_t EdbVertex::eY
private

for reconstructed ones - average of track connection

◆ eZ

Float_t EdbVertex::eZ
private

to avoid accuracy problem

points, used as local coordiantes origin (0,0,0)


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