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

#include <TrigL2Vertex.h>

Inheritance diagram for TrigVertexFitConstraint:
Collaboration diagram for TrigVertexFitConstraint:

Public Member Functions

 TrigVertexFitConstraint (double, const TrigVertexFitInputTrack *, const TrigVertexFitInputTrack *)
 two-track mass constraint
 TrigVertexFitConstraint (double, const TrigVertexFitInputTrack *, const TrigVertexFitInputTrack *, const TrigVertexFitInputTrack *)
 three-track mass constraint
 ~TrigVertexFitConstraint ()
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
double getValue ()
 returns a mass of the constraint

Protected Attributes

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

Private Member Functions

double calculateInvariantMass (TrigL2Vertex *pV)

Private Attributes

std::list< const TrigVertexFitInputTrack * > m_trackList
double m_value

Detailed Description

Definition at line 104 of file TrigL2Vertex.h.

Constructor & Destructor Documentation

◆ TrigVertexFitConstraint() [1/2]

TrigVertexFitConstraint::TrigVertexFitConstraint ( double m,
const TrigVertexFitInputTrack * pT1,
const TrigVertexFitInputTrack * pT2 )

two-track mass constraint

Definition at line 343 of file TrigL2Vertex.cxx.

345 : m_value(m)
346{
347 m_trackList.clear();
348 m_trackList.push_back(pT1);
349 m_trackList.push_back(pT2);
350}
std::list< const TrigVertexFitInputTrack * > m_trackList

◆ TrigVertexFitConstraint() [2/2]

TrigVertexFitConstraint::TrigVertexFitConstraint ( double m,
const TrigVertexFitInputTrack * pT1,
const TrigVertexFitInputTrack * pT2,
const TrigVertexFitInputTrack * pT3 )

three-track mass constraint

Definition at line 352 of file TrigL2Vertex.cxx.

355 : m_value(m)
356{
357 m_trackList.clear();
358 m_trackList.push_back(pT1);m_trackList.push_back(pT2);m_trackList.push_back(pT3);
359}

◆ ~TrigVertexFitConstraint()

TrigVertexFitConstraint::~TrigVertexFitConstraint ( )

Definition at line 362 of file TrigL2Vertex.cxx.

363{
364 m_trackList.clear();
365}

Member Function Documentation

◆ calculateInvariantMass()

double TrigVertexFitConstraint::calculateInvariantMass ( TrigL2Vertex * pV)
private

Definition at line 367 of file TrigL2Vertex.cxx.

368{
369 const double C=0.029997;
370 const double B=20.84;
371
372 double invMass=0.0,alpha=C*B/1000.0;
373 double P[3];double E=0.0;
374
375 double* Rk = pV->getParametersVector();
376
377 int offset=0;
378 for(int i=0;i<3;i++) P[i]=0.0;
379
380 for(std::list<const TrigVertexFitInputTrack*>::iterator it=m_trackList.begin();it!=m_trackList.end();++it)
381 {
382 offset=3+3*(*it)->getIndex();
383 double mass=(*it)->getMass()/1000.0;
384 double pT=fabs(Rk[offset+2]);
385 double p=pT/sin(Rk[offset+1]);
386
387 double psi=-asin(alpha*(Rk[0]*cos(Rk[offset])+Rk[1]*sin(Rk[offset]))/Rk[offset+2]);
388 double phiV=Rk[offset]+psi;
389 P[0]+=pT*cos(phiV);
390 P[1]+=pT*sin(phiV);
391 P[2]+=pT*cos(Rk[offset+1])/sin(Rk[offset+1]);
392 E+=sqrt(mass*mass+p*p);
393 }
394 invMass=sqrt(E*E-P[0]*P[0]-P[1]*P[1]-P[2]*P[2]);
395
396 return invMass;
397}
static Double_t P(Double_t *tt, Double_t *par)
double * getParametersVector()
returns vector of vertex fit parameters: vertex position + refitted track momenta at-perigee (sic !...
struct color C
double invMass(const I4Momentum &pA, const I4Momentum &pB)
invariant mass from two I4momentum references
Definition P4Helpers.h:252

◆ getChi2Distance()

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

implementation of abstract method from the base class

Implements TrigVertexFittingNode.

Definition at line 400 of file TrigL2Vertex.cxx.

401{
402 const double C=0.029997;
403 const double B=20.84;
404
405 double invMass=0.0,alpha=C*B/1000.0;
406 double P[3];
407 double E=0.0;
408 int i{0},j{0};
409 bool linFailed=false;
410 double* Rk = pV->getParametersVector();
411
412 int offset=0;
413 for(i=0;i<3;i++) P[i]=0.0;
414
415 int nSize=3+3*pV->getTracks()->size();
416
417 for(std::list<const TrigVertexFitInputTrack*>::iterator it=m_trackList.begin();it!=m_trackList.end();++it)
418 {
419 offset=3+3*(*it)->getIndex();
420 double mass=(*it)->getMass()/1000.0;
421 double pT=fabs(Rk[offset+2]);
422 double p=pT/sin(Rk[offset+1]);
423
424 double psi=-asin(alpha*(Rk[0]*cos(Rk[offset])+Rk[1]*sin(Rk[offset]))/Rk[offset+2]);
425 double phiV=Rk[offset]+psi;
426 P[0]+=pT*cos(phiV);
427 P[1]+=pT*sin(phiV);
428 P[2]+=pT*cos(Rk[offset+1])/sin(Rk[offset+1]);
429 E+=sqrt(mass*mass+p*p);
430 }
431 invMass=sqrt(E*E-P[0]*P[0]-P[1]*P[1]-P[2]*P[2]);
432
433 m_resid[0]=m_value/1000.0-invMass;
434 m_D[0][0]=0.0;m_D[0][1]=0.0;m_D[0][2]=0.0;
435
436 for(std::list<const TrigVertexFitInputTrack*>::iterator it=m_trackList.begin();it!=m_trackList.end();++it)
437 {
438 offset=3+3*(*it)->getIndex();
439
440 double mass=(*it)->getMass()/1000.0;
441 double Ck=(Rk[offset+2]<0.0)?-1.0:1.0;
442 double sinT=sin(Rk[offset+1]);
443 double cosT=cos(Rk[offset+1]);
444 double pT=fabs(Rk[offset+2]);
445 double p=pT/sinT;
446 double e=sqrt(p*p+mass*mass);
447 double sinF=sin(Rk[offset]);
448 double cosF=cos(Rk[offset]);
449
450 double sinPsi=-alpha*(Rk[0]*cosF+Rk[1]*sinF)/Rk[offset+2];
451 if(fabs(sinPsi)>1.0)
452 {
453 linFailed=true;
454 break;
455 }
456 double psi=asin(sinPsi);
457 double cosPsi=sqrt(1.0-sinPsi*sinPsi);
458 double phiV=Rk[offset]+psi;
459
460 double aCos=alpha*Ck/cosPsi;
461 double dP=P[1]*cos(phiV)-P[0]*sin(phiV);
462 double eE=E*Rk[offset+2]/(e*sinT);
463
464 m_D[0][0]+=dP*cosF*aCos;
465 m_D[0][1]+=dP*sinF*aCos;
466 m_D[0][offset]=-dP*Ck*(Rk[offset+2]-Ck*aCos*(Rk[1]*cosF-Rk[0]*sinF));
467 m_D[0][offset+1]=(Rk[offset+2]/(sinT*sinT))*(P[2]*Ck-eE*cosT);
468 m_D[0][offset+2]=eE/sinT-Ck*(P[0]*cos(phiV)+P[1]*sin(phiV)+P[2]*cosT/sinT)+dP*Ck*sinPsi/cosPsi;
469 }
470 for(i=0;i<nSize;i++) m_D[0][i]/=invMass;
471 if(linFailed) return -999.9;
472 /*
473 for(i=0;i<nSize;i++) printf("%E ",m_D[0][i]);
474 printf("\n");
475 printf("----- numerical Jacobian -----\n");
476
477 double oldPar, newPar, Delta, newMass, der;
478 for(i=0;i<nSize;i++)
479 {
480 oldPar=Rk[i];
481 Delta=0.00001*oldPar;
482 newPar=oldPar+Delta;
483 Rk[i]=newPar;
484 newMass=calculateInvariantMass(pV);
485 der=(newMass-invMass)/Delta;
486 m_D[1][i]=der;
487 Rk[i]=oldPar;
488 }
489 for(i=0;i<nSize;i++) printf("%E ",m_D[1][i]);
490 printf("\n");
491 printf("------------------------------\n");
492 */
493
494 double covM=1e-12;
495
496 for(i=0;i<nSize;i++)
497 {
498 for(j=i;j<nSize;j++)
499 {
500 double dCov=pV->m_Gk[i][j]*m_D[0][i]*m_D[0][j];
501 if(i!=j) dCov*=2.0;
502 covM+=dCov;
503 }
504 }
505 m_V[0][0]=1.0/covM;
506 return ((m_resid[0]*m_resid[0])*m_V[0][0]);
507}
double m_D[2][MAX_SIZE_VERT_COVM]

◆ getValue()

double TrigVertexFitConstraint::getValue ( )

returns a mass of the constraint

Definition at line 540 of file TrigL2Vertex.cxx.

541{
542 return m_value;
543}

◆ report()

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

Implements TrigVertexFittingNode.

Definition at line 534 of file TrigL2Vertex.cxx.

535{
536 out<<"Mass constraint with "<<m_trackList.size()<<" : mass="<<m_value<<endmsg;
537 return out;
538}
#define endmsg

◆ updateVertex()

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

implementation of abstract method from the base class

Implements TrigVertexFittingNode.

Definition at line 509 of file TrigL2Vertex.cxx.

510{
511 double K[MAX_SIZE_VERT_COVM];
512 int i{0},j{0};
513 const int nSize=3+3*pV->getTracks()->size();
514 double gain{0};
515
516 for(i=0;i<nSize;i++)
517 {
518 gain=0.0;
519 for(j=0;j<nSize;j++) gain+=pV->m_Gk[i][j]*m_D[0][j];
520 m_D[1][i]=gain;
521 K[i]=gain*m_V[0][0];
522 }
523 for(i=0;i<nSize;i++)
524 {
525 pV->getParametersVector()[i]+=K[i]*m_resid[0];
526 for(j=i;j<nSize;j++)
527 {
528 pV->m_Gk[i][j]-=K[i]*m_D[1][j];
529 pV->m_Gk[j][i]=pV->m_Gk[i][j];
530 }
531 }
532}
#define MAX_SIZE_VERT_COVM

Member Data Documentation

◆ m_D

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

Definition at line 48 of file TrigL2Vertex.h.

48{};

◆ m_resid

double TrigVertexFittingNode::m_resid[2] {}
protectedinherited

Definition at line 46 of file TrigL2Vertex.h.

46{};

◆ m_trackList

std::list<const TrigVertexFitInputTrack*> TrigVertexFitConstraint::m_trackList
private

Definition at line 117 of file TrigL2Vertex.h.

◆ m_V

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

Definition at line 47 of file TrigL2Vertex.h.

47{};

◆ m_value

double TrigVertexFitConstraint::m_value
private

Definition at line 118 of file TrigL2Vertex.h.


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