ATLAS Offline Software
Public Member Functions | Protected Attributes | Private Attributes | List of all members
TrigPrimaryVertexTrack Class Reference

#include <TrigPrimaryVertexTrack.h>

Inheritance diagram for TrigPrimaryVertexTrack:
Collaboration diagram for TrigPrimaryVertexTrack:

Public Member Functions

 TrigPrimaryVertexTrack (const TrigInDetTrack *)
 constructor for L2 tracks More...
 
 TrigPrimaryVertexTrack (const Trk::Track *)
 constructor for EF (offline) tracks More...
 
 ~TrigPrimaryVertexTrack ()
 
const TrigInDetTrackgetTrigTrack ()
 getter for L2 tracks More...
 
const Trk::TrackgetTrkTrack ()
 getter for EF (offline) tracks More...
 
virtual double getChi2Distance (TrigL2Vertex *)
 implementation of abstract method from the base class More...
 
virtual void updateVertex (TrigL2Vertex *)
 implementation of abstract method from the base class More...
 
virtual MsgStream & report (MsgStream &) const
 
void setIndex (int)
 to be used by TrigVertexingTool More...
 
int getIndex () const
 to be used by TrigVertexingTool More...
 
int getTrackType ()
 0: L2 track, 1: EF(offline) track More...
 
bool isActive ()
 if true this track will be used in the vertex fit otherwise it will be masked More...
 
void activate ()
 sets m_isActive to true More...
 
void mask ()
 sets m_isActive to false More...
 
const double * Perigee () const
 track parameters at the perigee More...
 
double PerigeeCovariance (int, int) const
 covariance of track parameters at the perigee More...
 
double getChi2Contribution ()
 chi2-contribution to the vertex fit More...
 

Protected Attributes

double m_resid [2]
 
double m_V [2][2]
 
double m_D [2][MAX_SIZE_VERT_COVM]
 

Private Attributes

const TrigInDetTrackm_pTrigTrack
 
const Trk::Trackm_pTrkTrack
 
int m_nTrackType
 
int m_index = 0
 
double m_Vqq [3][3] {}
 
double m_Vuq [2][3] {}
 
double m_Vuu [2][2] {}
 
double m_u [2] {}
 
double m_q [3] {}
 
double m_Perigee [5] {}
 
double m_PerigeeCovariance [5][5] {}
 
double m_A [2][3] {}
 
double m_B [2][3] {}
 
double m_dChi2
 
bool m_active
 

Detailed Description

Definition at line 14 of file TrigPrimaryVertexTrack.h.

Constructor & Destructor Documentation

◆ TrigPrimaryVertexTrack() [1/2]

TrigPrimaryVertexTrack::TrigPrimaryVertexTrack ( const TrigInDetTrack pT)

constructor for L2 tracks

Definition at line 11 of file TrigPrimaryVertexTrack.cxx.

12 {
13  m_nTrackType=1;m_active=true;
14  m_pTrigTrack=pT;m_pTrkTrack=NULL;m_dChi2=-100.0;
15  double Ck[5][5];
16  int i,j;
17 
18  const TrigInDetTrackFitPar* p=pT->param();
19 
20  double a=-2.0*exp(-p->eta())/(1.0+exp(-2.0*p->eta()));
21  m_u[0]=p->a0();
22  m_u[1]=p->z0();
23  m_q[0]=p->phi0();
24  m_q[1]=2.0*atan(exp(-p->eta()));
25  m_q[2]=p->pT()*1e-3;
26 
27 
28  Ck[0][0]=(*p->cov())[0];Ck[0][1]=Ck[1][0]=(*p->cov())[2];
29  Ck[0][2]=Ck[2][0]=(*p->cov())[1];Ck[0][3]=Ck[3][0]=(*p->cov())[3];
30  Ck[0][4]=Ck[4][0]=(*p->cov())[4];Ck[1][1]=(*p->cov())[9];
31  Ck[1][2]=Ck[2][1]=(*p->cov())[6];Ck[1][3]=Ck[3][1]=(*p->cov())[10];
32  Ck[1][4]=Ck[4][1]=(*p->cov())[11];Ck[2][2]=(*p->cov())[5];
33  Ck[2][3]=Ck[3][2]=(*p->cov())[7];Ck[2][4]=Ck[4][2]=(*p->cov())[8];
34  Ck[3][3]=(*p->cov())[12];Ck[3][4]=Ck[4][3]=(*p->cov())[13];
35  Ck[4][4]=(*p->cov())[14];
36  for(i=0;i<5;i++)
37  {
38  Ck[3][i]=a*Ck[3][i];Ck[i][3]=Ck[3][i];
39  }
40  Ck[3][3]*=a;
41  for(i=0;i<5;i++)
42  {
43  Ck[4][i]=Ck[4][i]*1e-3; Ck[i][4]=Ck[4][i];
44  }
45  Ck[4][4]*=1e-3;
46 
47 
48  for(i=0;i<2;i++) for(j=0;j<2;j++) m_Vuu[i][j]=Ck[i][j];
49  for(i=0;i<2;i++) for(j=0;j<3;j++) m_Vuq[i][j]=Ck[i][j+2];
50  for(i=0;i<3;i++) for(j=0;j<3;j++) m_Vqq[i][j]=Ck[i+2][j+2];
51 
52  m_Perigee[0]=m_u[0];m_Perigee[1]=m_u[1];m_Perigee[2]=m_q[0];m_Perigee[3]=m_q[1];m_Perigee[4]=m_q[2];
53 }

◆ TrigPrimaryVertexTrack() [2/2]

TrigPrimaryVertexTrack::TrigPrimaryVertexTrack ( const Trk::Track pT)

constructor for EF (offline) tracks

Definition at line 55 of file TrigPrimaryVertexTrack.cxx.

56 {
57  m_nTrackType=2;m_active=true;
58  m_pTrkTrack=pT;m_pTrigTrack=NULL;m_dChi2=-100.0;
59 
60  const Trk::Perigee* pP= pT->perigeeParameters();
61 
62  if(pP!=NULL) {
63  m_u[0]=pP->parameters()[Trk::d0];
64  m_u[1]=pP->parameters()[Trk::z0];
65  m_q[0]=pP->parameters()[Trk::phi0];
66  m_q[1]=pP->parameters()[Trk::theta];
67 
68  const double ptC=1000.0*pP->parameters()[Trk::qOverP];
69  const double inv_ptC = 1. / ptC;
70 
71  m_q[2]=sin(pP->parameters()[Trk::theta])*inv_ptC;
72 
73  const AmgSymMatrix(5)& TC= (*pP->covariance());
74 
75  const double a=cos(pP->parameters()[Trk::theta])*inv_ptC;
76  const double b=-sin(pP->parameters()[Trk::theta])/(pP->parameters()[Trk::qOverP]*ptC);
77 
78 
79  const double Ck[5][5] = {{TC(0,0), TC(1,0), TC(2,0), TC(3,0), a*TC(3,0) + b*TC(4,0)},
80  {TC(1,0), TC(1,1), TC(2,1), TC(3,1), a*TC(3,1) + b*TC(4,1)},
81  {TC(2,0), TC(2,1), TC(2,2), TC(3,2), a*TC(3,2) + b*TC(4,2)},
82  {TC(3,0), TC(3,1), TC(3,2), TC(3,3), a*TC(3,3) + b*TC(4,3)},
83  {a*TC(3,0) + b*TC(4,0), a*TC(3,1) + b*TC(4,1), a*TC(3,2) + b*TC(4,2), a*TC(3,3) + b*TC(4,3), a*a*TC(3,3) + 2*a*b*TC(4,3) + b*b*TC(4,4)}};
84 
85  for(int i=0;i<2;i++) for(int j=0;j<2;j++) m_Vuu[i][j]=Ck[i][j];
86  for(int i=0;i<3;i++) for(int j=0;j<3;j++) m_Vqq[i][j]=Ck[i+2][j+2];
87  for(int i=0;i<2;i++) for(int j=0;j<3;j++) m_Vuq[i][j]=Ck[i][j+2];
88 
89  m_Perigee[0]=m_u[0];m_Perigee[1]=m_u[1];m_Perigee[2]=m_q[0];m_Perigee[3]=m_q[1];m_Perigee[4]=m_q[2];
90  }
91 }

◆ ~TrigPrimaryVertexTrack()

TrigPrimaryVertexTrack::~TrigPrimaryVertexTrack ( )

Definition at line 93 of file TrigPrimaryVertexTrack.cxx.

94 {
95 
96 }

Member Function Documentation

◆ activate()

void TrigPrimaryVertexTrack::activate ( )

sets m_isActive to true

Definition at line 103 of file TrigPrimaryVertexTrack.cxx.

104 {
105  m_active=true;
106 }

◆ getChi2Contribution()

double TrigPrimaryVertexTrack::getChi2Contribution ( )

chi2-contribution to the vertex fit

Definition at line 243 of file TrigPrimaryVertexTrack.cxx.

244 {
245  return m_dChi2;
246 }

◆ getChi2Distance()

double TrigPrimaryVertexTrack::getChi2Distance ( TrigL2Vertex pV)
virtual

implementation of abstract method from the base class

Implements TrigVertexFittingNode.

Definition at line 149 of file TrigPrimaryVertexTrack.cxx.

150 {
151  const double C=0.02997;
152  const double B=20.84;
153 
154  double Sk[2][2],detr,chi2;
155  double AC[2][3],BV[2][3],h[2];
156  int i,j,k;
157  double psi,sinPsi,xv,yv,zv,
158  cosPhi0,sinPhi0,ctt,sint,phi0,theta0,P0;
159  const double alpha = C*B*1e-3;
160  const double inv_alpha = 1. / alpha;
161 
162  xv=pV->getParametersVector()[0];
163  yv=pV->getParametersVector()[1];
164  zv=pV->getParametersVector()[2];
165 
166  phi0=m_q[0];
167  theta0=m_q[1];
168  P0=m_q[2];
169 
170  cosPhi0=cos(phi0);sinPhi0=sin(phi0);
171  sinPsi=-alpha*(xv*cosPhi0+yv*sinPhi0)/P0;
172  if(fabs(sinPsi)>1.0) return -999.9;
173  const double cosPsi=sqrt(1.0-sinPsi*sinPsi);
174  const double inv_cosPsi = 1. / cosPsi;
175  psi=asin(sinPsi);
176  sint=sin(theta0);
177  ctt=cos(theta0)/sint;
178 
179  m_A[0][0]=-sin(phi0+psi)*inv_cosPsi;
180  m_A[0][1]= cos(phi0+psi)*inv_cosPsi;
181  m_A[0][2]=0.0;
182 
183  m_A[1][0]=-ctt*cosPhi0*inv_cosPsi;
184  m_A[1][1]=-ctt*sinPhi0*inv_cosPsi;
185  m_A[1][2]=1.0;
186 
187  m_B[0][0]=-xv*m_A[0][1]+yv*m_A[0][0];
188  m_B[0][1]=0.0;
189  m_B[0][2]=(1.0-inv_cosPsi)*inv_alpha;
190 
191  m_B[1][0]=-xv*m_A[1][1]+yv*m_A[1][0];
192  m_B[1][1]=-P0*psi/(alpha*sint*sint);
193  m_B[1][2]=ctt*(psi-sinPsi*inv_cosPsi)*inv_alpha;
194 
195  h[0]=yv*cosPhi0-xv*sinPhi0+P0*(1-cosPsi)*inv_alpha;
196  h[1]=zv+P0*ctt*psi*inv_alpha;
197 
198  m_resid[0]=m_u[0]-h[0];
199  m_resid[1]=m_u[1]-h[1];
200 
201  for(i=0;i<2;i++) for(j=0;j<2;j++) Sk[i][j]=m_Vuu[i][j];
202  for(i=0;i<2;i++) for(j=0;j<3;j++)
203  {
204  AC[i][j]=0.0;
205  for(k=0;k<3;k++) AC[i][j]+=m_A[i][k]*pV->m_Gk[k][j];
206  }
207  for(i=0;i<2;i++) for(j=0;j<2;j++)
208  {
209  for(k=0;k<3;k++) Sk[i][j]+=AC[i][k]*m_A[j][k];
210  }
211  for(i=0;i<2;i++)
212  for(j=0;j<3;j++)
213  {
214  BV[i][j]=0.0;
215  for(k=0;k<3;k++) BV[i][j]+=m_B[i][k]*m_Vqq[k][j];
216  }
217  for(i=0;i<2;i++)
218  for(j=0;j<2;j++)
219  {
220  for(k=0;k<3;k++) Sk[i][j]+=BV[i][k]*m_B[j][k];
221  }
222  Sk[0][0]-=2.0*(m_Vuq[0][0]*m_B[0][0]+m_Vuq[0][1]*m_B[0][1]+m_Vuq[0][2]*m_B[0][2]);
223  Sk[1][1]-=2.0*(m_Vuq[1][0]*m_B[1][0]+m_Vuq[1][1]*m_B[1][1]+m_Vuq[1][2]*m_B[1][2]);
224  Sk[0][1]-=m_Vuq[1][0]*m_B[0][0]+m_Vuq[1][1]*m_B[0][1]+m_Vuq[1][2]*m_B[0][2]+
225  m_Vuq[0][0]*m_B[1][0]+m_Vuq[0][1]*m_B[1][1]+m_Vuq[0][2]*m_B[1][2];
226  Sk[1][0]=Sk[0][1];
227 
228  detr=1.0/(Sk[0][0]*Sk[1][1]-Sk[0][1]*Sk[1][0]);
229  m_V[0][0]=Sk[1][1]*detr;
230  m_V[1][1]=Sk[0][0]*detr;
231  m_V[0][1]=m_V[1][0]=-Sk[0][1]*detr;
232 
233  chi2=m_V[0][0]*m_resid[0]*m_resid[0]+m_V[1][1]*m_resid[1]*m_resid[1]+
234  2.0*m_V[0][1]*m_resid[1]*m_resid[0];
235 
236  for(i=0;i<2;i++)
237  for(j=0;j<3;j++) m_D[i][j]=AC[i][j];
238 
239  m_dChi2=chi2;
240  return chi2;
241 }

◆ getIndex()

int TrigPrimaryVertexTrack::getIndex ( ) const

to be used by TrigVertexingTool

Definition at line 134 of file TrigPrimaryVertexTrack.cxx.

135 {
136  return m_index;
137 }

◆ getTrackType()

int TrigPrimaryVertexTrack::getTrackType ( )

0: L2 track, 1: EF(offline) track

Definition at line 113 of file TrigPrimaryVertexTrack.cxx.

114 {
115  return m_nTrackType;
116 }

◆ getTrigTrack()

const TrigInDetTrack * TrigPrimaryVertexTrack::getTrigTrack ( )

getter for L2 tracks

Definition at line 119 of file TrigPrimaryVertexTrack.cxx.

120 {
121  return m_pTrigTrack;
122 }

◆ getTrkTrack()

const Trk::Track * TrigPrimaryVertexTrack::getTrkTrack ( )

getter for EF (offline) tracks

Definition at line 124 of file TrigPrimaryVertexTrack.cxx.

125 {
126  return m_pTrkTrack;
127 }

◆ isActive()

bool TrigPrimaryVertexTrack::isActive ( )

if true this track will be used in the vertex fit otherwise it will be masked

Definition at line 108 of file TrigPrimaryVertexTrack.cxx.

109 {
110  return m_active;
111 }

◆ mask()

void TrigPrimaryVertexTrack::mask ( )

sets m_isActive to false

Definition at line 98 of file TrigPrimaryVertexTrack.cxx.

99 {
100  m_active=false;
101 }

◆ Perigee()

const double * TrigPrimaryVertexTrack::Perigee ( ) const

track parameters at the perigee

Definition at line 139 of file TrigPrimaryVertexTrack.cxx.

140 {
141  return &m_Perigee[0];
142 }

◆ PerigeeCovariance()

double TrigPrimaryVertexTrack::PerigeeCovariance ( int  i,
int  j 
) const

covariance of track parameters at the perigee

Definition at line 144 of file TrigPrimaryVertexTrack.cxx.

145 {
146  return m_PerigeeCovariance[i][j];
147 }

◆ report()

MsgStream & TrigPrimaryVertexTrack::report ( MsgStream &  out) const
virtual

Implements TrigVertexFittingNode.

Definition at line 273 of file TrigPrimaryVertexTrack.cxx.

274 {
275  int i;
276 
277  out<<"Primary track "<<m_index<<endmsg;
278  for(i=0;i<2;i++)
279  {
280  out<<" u"<<i<<" = "<<m_u[i]<<" "<<m_Vuu[i][0]<<" "<<m_Vuu[i][1]<<endmsg;
281  }
282  for(i=0;i<3;i++)
283  {
284  out<<" q"<<i<<" = "<<m_q[i]<<" "<<m_Vqq[i][0]<<" "<<m_Vqq[i][1]<<" "<<m_Vqq[i][2]<<endmsg;
285  }
286  return out;
287 }

◆ setIndex()

void TrigPrimaryVertexTrack::setIndex ( int  i)

to be used by TrigVertexingTool

Definition at line 129 of file TrigPrimaryVertexTrack.cxx.

130 {
131  m_index=i;
132 }

◆ updateVertex()

void TrigPrimaryVertexTrack::updateVertex ( TrigL2Vertex pV)
virtual

implementation of abstract method from the base class

Implements TrigVertexFittingNode.

Definition at line 248 of file TrigPrimaryVertexTrack.cxx.

249 {
250  int i,j,k;
251 
252  double K[2][3];
253 
254  for(i=0;i<2;i++)
255  {
256  for(j=0;j<3;j++)
257  {
258  K[i][j]=0.0;
259  for(k=0;k<2;k++) K[i][j]+=m_D[k][j]*m_V[k][i];
260  }
261  }
262  for(i=0;i<3;i++)
263  {
264  pV->getParametersVector()[i]+=K[0][i]*m_resid[0]+K[1][i]*m_resid[1];
265  for(j=i;j<3;j++)
266  {
267  pV->m_Gk[i][j]-=K[0][i]*m_D[0][j]+K[1][i]*m_D[1][j];
268  pV->m_Gk[j][i]=pV->m_Gk[i][j];
269  }
270  }
271 }

Member Data Documentation

◆ m_A

double TrigPrimaryVertexTrack::m_A[2][3] {}
private

Definition at line 47 of file TrigPrimaryVertexTrack.h.

◆ m_active

bool TrigPrimaryVertexTrack::m_active
private

Definition at line 50 of file TrigPrimaryVertexTrack.h.

◆ m_B

double TrigPrimaryVertexTrack::m_B[2][3] {}
private

Definition at line 48 of file TrigPrimaryVertexTrack.h.

◆ m_D

double TrigVertexFittingNode::m_D[2][MAX_SIZE_VERT_COVM]
protectedinherited

Definition at line 48 of file TrigL2Vertex.h.

◆ m_dChi2

double TrigPrimaryVertexTrack::m_dChi2
private

Definition at line 49 of file TrigPrimaryVertexTrack.h.

◆ m_index

int TrigPrimaryVertexTrack::m_index = 0
private

Definition at line 39 of file TrigPrimaryVertexTrack.h.

◆ m_nTrackType

int TrigPrimaryVertexTrack::m_nTrackType
private

Definition at line 38 of file TrigPrimaryVertexTrack.h.

◆ m_Perigee

double TrigPrimaryVertexTrack::m_Perigee[5] {}
private

Definition at line 45 of file TrigPrimaryVertexTrack.h.

◆ m_PerigeeCovariance

double TrigPrimaryVertexTrack::m_PerigeeCovariance[5][5] {}
private

Definition at line 46 of file TrigPrimaryVertexTrack.h.

◆ m_pTrigTrack

const TrigInDetTrack* TrigPrimaryVertexTrack::m_pTrigTrack
private

Definition at line 36 of file TrigPrimaryVertexTrack.h.

◆ m_pTrkTrack

const Trk::Track* TrigPrimaryVertexTrack::m_pTrkTrack
private

Definition at line 37 of file TrigPrimaryVertexTrack.h.

◆ m_q

double TrigPrimaryVertexTrack::m_q[3] {}
private

Definition at line 44 of file TrigPrimaryVertexTrack.h.

◆ m_resid

double TrigVertexFittingNode::m_resid[2]
protectedinherited

Definition at line 46 of file TrigL2Vertex.h.

◆ m_u

double TrigPrimaryVertexTrack::m_u[2] {}
private

Definition at line 43 of file TrigPrimaryVertexTrack.h.

◆ m_V

double TrigVertexFittingNode::m_V[2][2]
protectedinherited

Definition at line 47 of file TrigL2Vertex.h.

◆ m_Vqq

double TrigPrimaryVertexTrack::m_Vqq[3][3] {}
private

Definition at line 40 of file TrigPrimaryVertexTrack.h.

◆ m_Vuq

double TrigPrimaryVertexTrack::m_Vuq[2][3] {}
private

Definition at line 41 of file TrigPrimaryVertexTrack.h.

◆ m_Vuu

double TrigPrimaryVertexTrack::m_Vuu[2][2] {}
private

Definition at line 42 of file TrigPrimaryVertexTrack.h.


The documentation for this class was generated from the following files:
CalculateHighPtTerm.pT
pT
Definition: ICHEP2016/CalculateHighPtTerm.py:57
TrigPrimaryVertexTrack::m_A
double m_A[2][3]
Definition: TrigPrimaryVertexTrack.h:47
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
TrigPrimaryVertexTrack::m_pTrkTrack
const Trk::Track * m_pTrkTrack
Definition: TrigPrimaryVertexTrack.h:37
TrigPrimaryVertexTrack::m_Vqq
double m_Vqq[3][3]
Definition: TrigPrimaryVertexTrack.h:40
TrigPrimaryVertexTrack::m_Perigee
double m_Perigee[5]
Definition: TrigPrimaryVertexTrack.h:45
TrigVertexFittingNode::m_D
double m_D[2][MAX_SIZE_VERT_COVM]
Definition: TrigL2Vertex.h:48
TrigPrimaryVertexTrack::m_nTrackType
int m_nTrackType
Definition: TrigPrimaryVertexTrack.h:38
Trk::ParametersT
Dummy class used to allow special convertors to be called for surfaces owned by a detector element.
Definition: EMErrorDetail.h:25
TrigInDetTrackFitPar
Definition: TrigInDetTrackFitPar.h:67
InDetAccessor::phi0
@ phi0
Definition: InDetAccessor.h:33
DMTest::C
C_v1 C
Definition: C.h:26
python.AthDsoLogger.out
out
Definition: AthDsoLogger.py:71
TRT_PAI_gasdata::AC
const float AC
Definition: TRT_PAI_gasdata.h:27
Trk::z0
@ z0
Definition: ParamDefs.h:70
TrigVertexFittingNode::m_V
double m_V[2][2]
Definition: TrigL2Vertex.h:47
drawFromPickle.cos
cos
Definition: drawFromPickle.py:36
drawFromPickle.exp
exp
Definition: drawFromPickle.py:36
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:52
drawFromPickle.atan
atan
Definition: drawFromPickle.py:36
TrigPrimaryVertexTrack::m_q
double m_q[3]
Definition: TrigPrimaryVertexTrack.h:44
TrigVertexFittingNode::m_resid
double m_resid[2]
Definition: TrigL2Vertex.h:46
TrigPrimaryVertexTrack::m_index
int m_index
Definition: TrigPrimaryVertexTrack.h:39
lumiFormat.i
int i
Definition: lumiFormat.py:92
Trk::theta
@ theta
Definition: ParamDefs.h:72
endmsg
#define endmsg
Definition: AnalysisConfig_Ntuple.cxx:63
chi2
double chi2(TH1 *h0, TH1 *h1)
Definition: comparitor.cxx:522
TrigPrimaryVertexTrack::m_PerigeeCovariance
double m_PerigeeCovariance[5][5]
Definition: TrigPrimaryVertexTrack.h:46
Trk::d0
@ d0
Definition: ParamDefs.h:69
plotBeamSpotMon.b
b
Definition: plotBeamSpotMon.py:77
TrigL2Vertex::m_Gk
double m_Gk[MAX_SIZE_VERT_COVM][MAX_SIZE_VERT_COVM]
Definition: TrigL2Vertex.h:182
dqt_zlumi_alleff_HIST.B
B
Definition: dqt_zlumi_alleff_HIST.py:110
TrigPrimaryVertexTrack::m_pTrigTrack
const TrigInDetTrack * m_pTrigTrack
Definition: TrigPrimaryVertexTrack.h:36
TrigPrimaryVertexTrack::m_Vuu
double m_Vuu[2][2]
Definition: TrigPrimaryVertexTrack.h:42
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
a
TList * a
Definition: liststreamerinfos.cxx:10
h
TrigPrimaryVertexTrack::m_Vuq
double m_Vuq[2][3]
Definition: TrigPrimaryVertexTrack.h:41
TrigL2Vertex::getParametersVector
double * getParametersVector()
returns vector of vertex fit parameters: vertex position + refitted track momenta at-perigee (sic !...
Definition: TrigL2Vertex.cxx:650
Trk::qOverP
@ qOverP
perigee
Definition: ParamDefs.h:73
TrigPrimaryVertexTrack::m_B
double m_B[2][3]
Definition: TrigPrimaryVertexTrack.h:48
drawFromPickle.sin
sin
Definition: drawFromPickle.py:36
TrigPrimaryVertexTrack::m_u
double m_u[2]
Definition: TrigPrimaryVertexTrack.h:43
TrigPrimaryVertexTrack::m_active
bool m_active
Definition: TrigPrimaryVertexTrack.h:50
TrigPrimaryVertexTrack::m_dChi2
double m_dChi2
Definition: TrigPrimaryVertexTrack.h:49
Trk::phi0
@ phi0
Definition: ParamDefs.h:71
fitman.k
k
Definition: fitman.py:528