ATLAS Offline Software
Loading...
Searching...
No Matches
TrigVertexFitInputTrack Class Reference

#include <TrigL2Vertex.h>

Inheritance diagram for TrigVertexFitInputTrack:
Collaboration diagram for TrigVertexFitInputTrack:

Public Member Functions

 TrigVertexFitInputTrack (const TrigInDetTrack *, double)
 constructor for L2 tracks
 TrigVertexFitInputTrack (const Trk::Track *, double)
 constructor for EF (offline) tracks
 ~TrigVertexFitInputTrack ()
const TrigInDetTrackgetTrigTrack ()
 getter for L2 tracks
const Trk::TrackgetTrkTrack ()
 getter for EF (offline) tracks
void initializeVertex (class TrigL2Vertex *)
 resets and fills its part of the covariance and parameter vector
bool linearize (class TrigL2Vertex *)
 re-calculates linearization of the measurement model in the vicinity of parameters provided by input TrigL2Vertex
virtual double getChi2Distance (class TrigL2Vertex *)
 implementation of abstract method from the base class
virtual void updateVertex (class TrigL2Vertex *)
 implementation of abstract method from the base class
virtual MsgStream & report (MsgStream &) const
void setIndex (int)
 to be used by TrigVertexingTool
int getIndex () const
 to be used by TrigVertexingTool
int getTrackType ()
 0: L2 track, 1: EF(offline) track
bool isActive ()
 if true this track will be used in the vertex fit otherwise it will be masked
void activate ()
 sets isActive to true
void mask ()
 sets isActive to false
void setMass (double)
 sets a mass of this particle
double getMass () const
 gets a mass of this particle
const double * Perigee () const
 track parameters at the perigee
double PerigeeCovariance (int, int) const
 covariance of track parameters at the perigee

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 {}
double m_mass {}
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_h [2] {}
bool m_active {true}

Detailed Description

Definition at line 61 of file TrigL2Vertex.h.

Constructor & Destructor Documentation

◆ TrigVertexFitInputTrack() [1/2]

TrigVertexFitInputTrack::TrigVertexFitInputTrack ( const TrigInDetTrack * pT,
double mass = 0.0 )

constructor for L2 tracks

Definition at line 9 of file TrigL2Vertex.cxx.

9 :
10 m_pTrigTrack(pT),
11 m_nTrackType(1),
12 m_mass(mass)
13{
14 double Ck[5][5];
15 int i{0},j{0};
16
17 const TrigInDetTrackFitPar* p=pT->param();
18
19 double a=-2.0*exp(-p->eta())/(1.0+exp(-2.0*p->eta()));
20 m_u[0]=p->a0();
21 m_u[1]=p->z0();
22 m_q[0]=p->phi0();
23 m_q[1]=2.0*atan(exp(-p->eta()));
24 m_q[2]=p->pT()/1000.0;
25
26
27 Ck[0][0]=(*p->cov())[0];Ck[0][1]=Ck[1][0]=(*p->cov())[2];
28 Ck[0][2]=Ck[2][0]=(*p->cov())[1];Ck[0][3]=Ck[3][0]=(*p->cov())[3];
29 Ck[0][4]=Ck[4][0]=(*p->cov())[4];Ck[1][1]=(*p->cov())[9];
30 Ck[1][2]=Ck[2][1]=(*p->cov())[6];Ck[1][3]=Ck[3][1]=(*p->cov())[10];
31 Ck[1][4]=Ck[4][1]=(*p->cov())[11];Ck[2][2]=(*p->cov())[5];
32 Ck[2][3]=Ck[3][2]=(*p->cov())[7];Ck[2][4]=Ck[4][2]=(*p->cov())[8];
33 Ck[3][3]=(*p->cov())[12];Ck[3][4]=Ck[4][3]=(*p->cov())[13];
34 Ck[4][4]=(*p->cov())[14];
35 for(i=0;i<5;i++)
36 {
37 Ck[3][i]=a*Ck[3][i];Ck[i][3]=Ck[3][i];
38 }
39 Ck[3][3]*=a;
40 for(i=0;i<5;i++)
41 {
42 Ck[4][i]=Ck[4][i]/1000.0;Ck[i][4]=Ck[4][i];
43 }
44 Ck[4][4]/=1000.0;
45
46
47 for(i=0;i<2;i++) for(j=0;j<2;j++) m_Vuu[i][j]=Ck[i][j];
48 for(i=0;i<2;i++) for(j=0;j<3;j++) m_Vuq[i][j]=Ck[i][j+2];
49 for(i=0;i<3;i++) for(j=0;j<3;j++) m_Vqq[i][j]=Ck[i+2][j+2];
50
51 /*
52 m_Vuu[0][0]=(*p->cov())[0];
53 m_Vuu[0][1]=m_Vuu[1][0]=(*p->cov())[2];
54 m_Vuu[1][1]=(*p->cov())[9];
55
56 m_Vuq[0][0]=(*p->cov())[1];m_Vuq[0][1]=(*p->cov())[3]*a;
57 m_Vuq[0][2]=(*p->cov())[4]*1e-3;m_Vuq[1][0]=(*p->cov())[6];m_Vuq[1][1]=(*p->cov())[10]*a;
58 m_Vuq[1][2]=(*p->cov())[11]*1e-3;
59
60 m_Vqq[0][0]=(*p->cov())[5];m_Vqq[0][1]=m_Vqq[1][0]=(*p->cov())[7]*a;
61 m_Vqq[0][2]=m_Vqq[2][0]=(*p->cov())[8]*1e-3;
62 m_Vqq[1][1]=(*p->cov())[12]*a*a;m_Vqq[1][2]=m_Vqq[2][1]=(*p->cov())[13]*a*1e-3;
63 m_Vqq[2][2]=(*p->cov())[14]*1e-6;
64 */
65 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];
66}
static Double_t a
const TrigInDetTrack * m_pTrigTrack

◆ TrigVertexFitInputTrack() [2/2]

TrigVertexFitInputTrack::TrigVertexFitInputTrack ( const Trk::Track * pT,
double mass = 0.0 )

constructor for EF (offline) tracks

Definition at line 68 of file TrigL2Vertex.cxx.

68 :
69 m_pTrkTrack(pT),
70 m_nTrackType(2),
71 m_mass(mass)
72{
73 double Ck[5][5];
74
75 const Trk::Perigee* pP= pT->perigeeParameters();
76
77 m_u[0]=pP->parameters()[Trk::d0];
78 m_u[1]=pP->parameters()[Trk::z0];
79 m_q[0]=pP->parameters()[Trk::phi0];
80 m_q[1]=pP->parameters()[Trk::theta];
81
82 double ptC=1000.0*pP->parameters()[Trk::qOverP];
83
84 m_q[2]=sin(pP->parameters()[Trk::theta])/ptC;
85
86 if( pP->covariance() ){
87 // has covariance matrix
88 const AmgSymMatrix(5)* TC = pP->covariance();
89
90 Ck[0][0]=(*TC)(Trk::d0,Trk::d0);
91 Ck[0][1]=Ck[1][0]=(*TC)(Trk::d0,Trk::z0);
92 Ck[0][2]=Ck[2][0]=(*TC)(Trk::d0,Trk::phi0);
93 Ck[0][3]=Ck[3][0]=(*TC)(Trk::d0,Trk::theta);
94 Ck[0][4]=Ck[4][0]=(*TC)(Trk::d0,Trk::qOverP);
95 Ck[1][1]=(*TC)(Trk::z0,Trk::z0);
96 Ck[1][2]=Ck[2][1]=(*TC)(Trk::z0,Trk::phi0);
97 Ck[1][3]=Ck[3][1]=(*TC)(Trk::z0,Trk::theta);
98 Ck[1][4]=Ck[4][1]=(*TC)(Trk::z0,Trk::qOverP);
99 Ck[2][2]=(*TC)(Trk::phi0,Trk::phi0);
100 Ck[2][3]=Ck[3][2]=(*TC)(Trk::phi0,Trk::theta);
101 Ck[2][4]=Ck[4][2]=(*TC)(Trk::phi0,Trk::qOverP);
102 Ck[3][3]=(*TC)(Trk::theta,Trk::theta);
103 Ck[3][4]=Ck[4][3]=(*TC)(Trk::theta,Trk::qOverP);
104 Ck[4][4]=(*TC)(Trk::qOverP,Trk::qOverP);
105
106 const double a = cos(pP->parameters()[Trk::theta])/ptC;
107 const double b = -sin(pP->parameters()[Trk::theta])/(pP->parameters()[Trk::qOverP]*ptC);
108
109 Ck[3][3]=Ck[3][3]+2.0*a*Ck[3][4]+a*a*Ck[4][4];
110 Ck[3][4]=Ck[4][3]=b*Ck[3][4]+a*b*Ck[4][4];
111 Ck[4][4]=b*b*Ck[4][4];
112 Ck[0][3]=Ck[3][0]=Ck[0][3]+a*Ck[0][4];Ck[0][4]*=b;Ck[4][0]=Ck[0][4];
113 Ck[1][3]=Ck[3][1]=Ck[1][3]+a*Ck[1][4];Ck[1][4]*=b;Ck[4][1]=Ck[1][4];
114 Ck[2][3]=Ck[3][2]=Ck[2][3]+a*Ck[2][4];Ck[2][4]*=b;Ck[4][2]=Ck[2][4];
115 for(int i=0;i<2;i++) for(int j=0;j<2;j++) m_Vuu[i][j]=Ck[i][j];
116 for(int i=0;i<3;i++) for(int j=0;j<3;j++) m_Vqq[i][j]=Ck[i+2][j+2];
117 for(int i=0;i<2;i++) for(int j=0;j<3;j++) m_Vuq[i][j]=Ck[i][j+2];
118
119 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];
120 }
121 else{
122 m_active=false;
123 // no covariance
124 }
125
126}
#define AmgSymMatrix(dim)
const Trk::Track * m_pTrkTrack
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ phi0
Definition ParamDefs.h:65
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64

◆ ~TrigVertexFitInputTrack()

TrigVertexFitInputTrack::~TrigVertexFitInputTrack ( )

Definition at line 128 of file TrigL2Vertex.cxx.

129{
130
131}

Member Function Documentation

◆ activate()

void TrigVertexFitInputTrack::activate ( )

sets isActive to true

Definition at line 143 of file TrigL2Vertex.cxx.

144{
145 m_active=true;
146}

◆ getChi2Distance()

double TrigVertexFitInputTrack::getChi2Distance ( class TrigL2Vertex * pV)
virtual

implementation of abstract method from the base class

Implements TrigVertexFittingNode.

Definition at line 196 of file TrigL2Vertex.cxx.

197{
198 const double C=0.02997;
199 const double B=20.84;
200 const double alpha=C*B/1000.0;
201
202 double Sk[2][2],detr,chi2;
203 double AC[2][3],BV[2][3];
204 int i{0},j{0},k{0};
205
206 const int shift=m_index*3;
207
208 k=3+shift;
209 for(i=0;i<3;i++)
210 for(j=0;j<3;j++)
211 pV->m_Gk[i+k][j+k]=m_Vqq[i][j];
212
214 const double xv=pV->getParametersVector()[0];
215 const double yv=pV->getParametersVector()[1];
216 const double zv=pV->getParametersVector()[2];
217 const double phi0=pV->getParametersVector()[3+shift];
218 const double theta0=pV->getParametersVector()[4+shift];
219 const double P0=pV->getParametersVector()[5+shift];
220
221 const double cosPhi0=cos(phi0);
222 const double sinPhi0=sin(phi0);
223 const double sinPsi=-alpha*(xv*cosPhi0+yv*sinPhi0)/P0;
224 if(fabs(sinPsi)>1.0) return -999.9;
225 const double cosPsi=sqrt(1.0-sinPsi*sinPsi);
226 const double psi=asin(sinPsi);
227 const double sint=sin(theta0);
228 const double ctt=cos(theta0)/sint;
229
230 m_A[0][0]=-sin(phi0+psi)/cosPsi;
231 m_A[0][1]= cos(phi0+psi)/cosPsi;
232 m_A[0][2]=0.0;
233
234 m_A[1][0]=-ctt*cosPhi0/cosPsi;
235 m_A[1][1]=-ctt*sinPhi0/cosPsi;
236 m_A[1][2]=1.0;
237
238 m_B[0][0]=-xv*m_A[0][1]+yv*m_A[0][0];
239 m_B[0][1]=0.0;
240 m_B[0][2]=(1.0-1.0/cosPsi)/alpha;
241
242 m_B[1][0]=-xv*m_A[1][1]+yv*m_A[1][0];
243 m_B[1][1]=-P0*psi/(alpha*sint*sint);
244 m_B[1][2]=ctt*(psi-sinPsi/cosPsi)/alpha;
245
246 m_h[0]=yv*cosPhi0-xv*sinPhi0+P0*(1-cosPsi)/alpha;
247 m_h[1]=zv+P0*ctt*psi/alpha;
248
249 m_resid[0]=m_u[0]-m_h[0];
250 m_resid[1]=m_u[1]-m_h[1];
251
252 for(i=0;i<2;i++) for(j=0;j<2;j++) Sk[i][j]=m_Vuu[i][j];
253 for(i=0;i<2;i++) for(j=0;j<3;j++)
254 {
255 AC[i][j]=0.0;
256 for(k=0;k<3;k++) AC[i][j]+=m_A[i][k]*pV->m_Gk[k][j];
257 }
258 for(i=0;i<2;i++) for(j=0;j<2;j++)
259 {
260 for(k=0;k<3;k++) Sk[i][j]+=AC[i][k]*m_A[j][k];
261 }
262 for(i=0;i<2;i++)
263 for(j=0;j<3;j++)
264 {
265 BV[i][j]=0.0;
266 for(k=0;k<3;k++) BV[i][j]+=m_B[i][k]*m_Vqq[k][j];
267 }
268 for(i=0;i<2;i++)
269 for(j=0;j<2;j++)
270 {
271 for(k=0;k<3;k++) Sk[i][j]+=BV[i][k]*m_B[j][k];
272 }
273 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]);
274 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]);
275 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]+
276 m_Vuq[0][0]*m_B[1][0]+m_Vuq[0][1]*m_B[1][1]+m_Vuq[0][2]*m_B[1][2];
277 Sk[1][0]=Sk[0][1];
278
279 detr=1.0/(Sk[0][0]*Sk[1][1]-Sk[0][1]*Sk[1][0]);
280 m_V[0][0]=Sk[1][1]*detr;
281 m_V[1][1]=Sk[0][0]*detr;
282 m_V[0][1]=m_V[1][0]=-Sk[0][1]*detr;
283
284 chi2=m_V[0][0]*m_resid[0]*m_resid[0]+m_V[1][1]*m_resid[1]*m_resid[1]+
285 2.0*m_V[0][1]*m_resid[1]*m_resid[0];
286
287 for(i=0;i<2;i++)
288 {
289 for(j=0;j<3;j++) m_D[i][j]=AC[i][j];
290 for(j=3;j<3+shift;j++)
291 {
292 m_D[i][j]=0.0;
293 for(k=0;k<3;k++)
294 m_D[i][j]+=m_A[i][k]*pV->m_Gk[k][j];
295 }
296 for(j=0;j<3;j++) m_D[i][j+3+shift]=BV[i][j]-m_Vuq[i][j];
297 }
298
299 return chi2;
300}
void initializeVertex(class TrigL2Vertex *)
resets and fills its part of the covariance and parameter vector
double m_D[2][MAX_SIZE_VERT_COVM]
double chi2(TH1 *h0, TH1 *h1)
struct color C

◆ getIndex()

int TrigVertexFitInputTrack::getIndex ( ) const

to be used by TrigVertexingTool

Definition at line 174 of file TrigL2Vertex.cxx.

175{
176 return m_index;
177}

◆ getMass()

double TrigVertexFitInputTrack::getMass ( ) const
inline

gets a mass of this particle

Definition at line 81 of file TrigL2Vertex.h.

◆ getTrackType()

int TrigVertexFitInputTrack::getTrackType ( )

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

Definition at line 153 of file TrigL2Vertex.cxx.

154{
155 return m_nTrackType;
156}

◆ getTrigTrack()

const TrigInDetTrack * TrigVertexFitInputTrack::getTrigTrack ( )

getter for L2 tracks

Definition at line 159 of file TrigL2Vertex.cxx.

160{
161 return m_pTrigTrack;
162}

◆ getTrkTrack()

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

getter for EF (offline) tracks

Definition at line 164 of file TrigL2Vertex.cxx.

165{
166 return m_pTrkTrack;
167}

◆ initializeVertex()

void TrigVertexFitInputTrack::initializeVertex ( class TrigL2Vertex * pV)

resets and fills its part of the covariance and parameter vector

Definition at line 189 of file TrigL2Vertex.cxx.

190{
191 const int i = 3+m_index*3;
192 for(int j=0;j<3;j++)
193 pV->getParametersVector()[i+j]=m_q[j];
194}

◆ isActive()

bool TrigVertexFitInputTrack::isActive ( )

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

Definition at line 148 of file TrigL2Vertex.cxx.

149{
150 return m_active;
151}

◆ linearize()

bool TrigVertexFitInputTrack::linearize ( class TrigL2Vertex * )

re-calculates linearization of the measurement model in the vicinity of parameters provided by input TrigL2Vertex

◆ mask()

void TrigVertexFitInputTrack::mask ( )

sets isActive to false

Definition at line 138 of file TrigL2Vertex.cxx.

139{
140 m_active=false;
141}

◆ Perigee()

const double * TrigVertexFitInputTrack::Perigee ( ) const

track parameters at the perigee

Definition at line 179 of file TrigL2Vertex.cxx.

180{
181 return &m_Perigee[0];
182}

◆ PerigeeCovariance()

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

covariance of track parameters at the perigee

Definition at line 184 of file TrigL2Vertex.cxx.

185{
186 return m_PerigeeCovariance[i][j];
187}
double m_PerigeeCovariance[5][5]

◆ report()

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

Implements TrigVertexFittingNode.

Definition at line 329 of file TrigL2Vertex.cxx.

330{
331 out<<"Track "<<m_index<<" : mass="<<m_mass<<endmsg;
332 for(int i=0;i<2;i++)
333 {
334 out<<" u"<<i<<" = "<<m_u[i]<<" "<<m_Vuu[i][0]<<" "<<m_Vuu[i][1]<<endmsg;
335 }
336 for(int i=0;i<3;i++)
337 {
338 out<<" q"<<i<<" = "<<m_q[i]<<" "<<m_Vqq[i][0]<<" "<<m_Vqq[i][1]<<" "<<m_Vqq[i][2]<<endmsg;
339 }
340 return out;
341}
#define endmsg

◆ setIndex()

void TrigVertexFitInputTrack::setIndex ( int i)

to be used by TrigVertexingTool

Definition at line 169 of file TrigL2Vertex.cxx.

170{
171 m_index=i;
172}

◆ setMass()

void TrigVertexFitInputTrack::setMass ( double m)

sets a mass of this particle

Definition at line 133 of file TrigL2Vertex.cxx.

134{
135 m_mass=m;
136}

◆ updateVertex()

void TrigVertexFitInputTrack::updateVertex ( class TrigL2Vertex * pV)
virtual

implementation of abstract method from the base class

Implements TrigVertexFittingNode.

Definition at line 303 of file TrigL2Vertex.cxx.

304{
305 int i{0},j{0},k{0};
306 const int nSize=6+3*m_index;
307
308 double K[2][MAX_SIZE_VERT_COVM];
309
310 for(i=0;i<2;i++)
311 {
312 for(j=0;j<nSize;j++)
313 {
314 K[i][j]=0.0;
315 for(k=0;k<2;k++) K[i][j]+=m_D[k][j]*m_V[k][i];
316 }
317 }
318 for(i=0;i<nSize;i++)
319 {
320 pV->getParametersVector()[i]+=K[0][i]*m_resid[0]+K[1][i]*m_resid[1];
321 for(j=i;j<nSize;j++)
322 {
323 pV->m_Gk[i][j]-=K[0][i]*m_D[0][j]+K[1][i]*m_D[1][j];
324 pV->m_Gk[j][i]=pV->m_Gk[i][j];
325 }
326 }
327}
#define MAX_SIZE_VERT_COVM

Member Data Documentation

◆ m_A

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

Definition at line 97 of file TrigL2Vertex.h.

97{};

◆ m_active

bool TrigVertexFitInputTrack::m_active {true}
private

Definition at line 101 of file TrigL2Vertex.h.

101{true};

◆ m_B

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

Definition at line 98 of file TrigL2Vertex.h.

98{};

◆ m_D

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

Definition at line 48 of file TrigL2Vertex.h.

48{};

◆ m_h

double TrigVertexFitInputTrack::m_h[2] {}
private

Definition at line 99 of file TrigL2Vertex.h.

99{};

◆ m_index

int TrigVertexFitInputTrack::m_index {}
private

Definition at line 88 of file TrigL2Vertex.h.

88{};

◆ m_mass

double TrigVertexFitInputTrack::m_mass {}
private

Definition at line 89 of file TrigL2Vertex.h.

89{};

◆ m_nTrackType

int TrigVertexFitInputTrack::m_nTrackType {}
private

Definition at line 87 of file TrigL2Vertex.h.

87{};

◆ m_Perigee

double TrigVertexFitInputTrack::m_Perigee[5] {}
private

Definition at line 95 of file TrigL2Vertex.h.

95{};

◆ m_PerigeeCovariance

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

Definition at line 96 of file TrigL2Vertex.h.

96{};

◆ m_pTrigTrack

const TrigInDetTrack* TrigVertexFitInputTrack::m_pTrigTrack {}
private

Definition at line 85 of file TrigL2Vertex.h.

85{};

◆ m_pTrkTrack

const Trk::Track* TrigVertexFitInputTrack::m_pTrkTrack {}
private

Definition at line 86 of file TrigL2Vertex.h.

86{};

◆ m_q

double TrigVertexFitInputTrack::m_q[3] {}
private

Definition at line 94 of file TrigL2Vertex.h.

94{};

◆ m_resid

double TrigVertexFittingNode::m_resid[2] {}
protectedinherited

Definition at line 46 of file TrigL2Vertex.h.

46{};

◆ m_u

double TrigVertexFitInputTrack::m_u[2] {}
private

Definition at line 93 of file TrigL2Vertex.h.

93{};

◆ m_V

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

Definition at line 47 of file TrigL2Vertex.h.

47{};

◆ m_Vqq

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

Definition at line 90 of file TrigL2Vertex.h.

90{};

◆ m_Vuq

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

Definition at line 91 of file TrigL2Vertex.h.

91{};

◆ m_Vuu

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

Definition at line 92 of file TrigL2Vertex.h.

92{};

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