Belle II Software light-2406-ragdoll
RecoTrack Class Reference

representation of all charged final states FIXME rename since this name is taken in tracking
More...

#include <RecoTrack.h>

Inheritance diagram for RecoTrack:
Collaboration diagram for RecoTrack:

Public Types

enum  TFParticleType {
  kInteractionPoint ,
  kOrigin ,
  kComposite ,
  kRecoResonance ,
  kInternalParticle ,
  kRecoTrack ,
  kResonance ,
  kRecoPhoton ,
  kRecoKlong ,
  kMissingParticle
}
 particle types
More...
 
typedef std::vector< Constraintconstraintlist
 alias
 
typedef std::vector< std::pair< const ParticleBase *, int > > indexmap
 alias
 

Public Member Functions

 RecoTrack (Belle2::Particle *bc, const ParticleBase *mother)
 constructor
 
virtual ~RecoTrack ()
 destructor
 
virtual ErrCode initParticleWithMother (FitParams &fitparams) override
 init with mother particle (replacing initPar2)

 
virtual ErrCode initMotherlessParticle (FitParams &fitparams) override
 init without mother particle

 
ErrCode initCovariance (FitParams &fitparams) const override
 init covariance matrix of this particle constraint
 
ErrCode updFltToMotherCopy (const FitParams *fitparams)
 update m_flt
 
virtual ErrCode projectRecoConstraint (const FitParams &, Projection &) const override
 project the constraint (calculate residuals)

 
ErrCode updateParams (double flt)
 updated the cached parameters
 
virtual int dimM () const override
 dimension (5)
 
virtual int type () const override
 type of the constraint

 
virtual int nFinalChargedCandidates () const override
 number of final charged candidates
 
virtual void addToConstraintList (constraintlist &alist, int depth) const override
 add to the list of constraints
 
ErrCode updFltToMother (const FitParams &fitparams)
 update flight length to mother
 
void setFlightLength (double flt)
 setter for the flight length
 
virtual std::string parname (int index) const override
 name
 
virtual int dim () const override
 this here sets the size in the state vector we are only interested in the momenta of photons and tracks as the position were the track ends or the cluster is, is not relevant for physics

 
virtual int momIndex () const override
 get momentum index
 
virtual bool hasEnergy () const override
 has an energy in the statevector?
 
virtual ErrCode projectConstraint (Constraint::Type, const FitParams &, Projection &) const override
 abstract abstract projection
 
virtual void updateIndex (int &offset)
 this sets the index for momentum, position, etc.
 
const ParticleBaselocate (Belle2::Particle *particle) const
 get particle base from basf2 particle
 
Belle2::Particleparticle () const
 get basf2 particle

 
int index () const
 get index

 
const ParticleBasemother () const
 getMother() / hasMother()
 
virtual ErrCode projectGeoConstraint (const FitParams &, Projection &) const
 project geometrical constraint
 
virtual ErrCode projectMassConstraintParticle (const FitParams &, Projection &) const
 project mass constraint using the particles parameters
 
virtual ErrCode projectMassConstraintDaughters (const FitParams &, Projection &) const
 project mass constraint using the parameters of the daughters
 
virtual ErrCode projectMassConstraint (const FitParams &, Projection &) const
 project mass constraint abstract
 
virtual void forceP4Sum (FitParams &) const
 force p4 sum conservation all along the tree
 
virtual int posIndex () const
 get vertex index (in statevector!)
 
virtual int tauIndex () const
 get tau index
 
virtual bool hasPosition () const
 get false

 
int eneIndex () const
 get energy index

 
virtual double chiSquare (const FitParams &) const
 get chi2
 
int charge () const
 get charge
 
virtual ParticleBaseaddDaughter (Belle2::Particle *, const ConstraintConfiguration &config, bool forceFitAll=false)
 add daughter

 
virtual void removeDaughter (const ParticleBase *pb)
 remove daughter
 
virtual void retrieveIndexMap (indexmap &anindexmap) const
 get index map

 
void setMother (const ParticleBase *m)
 set mother

 
void collectVertexDaughters (std::vector< ParticleBase * > &particles, int posindex)
 get vertex daughters
 

Static Public Member Functions

static ParticleBasecreateParticle (Belle2::Particle *particle, const ParticleBase *mother, const ConstraintConfiguration &config, bool forceFitAll=false)
 create the according treeFitter particle obj for a basf2 particle type

 
static ParticleBasecreateOrigin (Belle2::Particle *daughter, const ConstraintConfiguration &config, bool forceFitAll)
 create a custom origin particle or a beamspot
 

Protected Types

typedef std::vector< ParticleBase * > ParticleContainer
 just an alias
 

Protected Member Functions

ErrCode initTau (FitParams &par) const
 initialises tau as a length

 
void setIndex (int i)
 set Index (in statevector)
 

Static Protected Member Functions

static bool isAResonance (Belle2::Particle *particle)
 controls if a particle is treated as a resonance(lifetime=0) or a particle that has a finite lifetime.
 

Protected Attributes

Belle2::Particlem_particle
 pointer to framework type

 
const ParticleBasem_mother
 motherparticle
 
std::vector< ParticleBase * > m_daughters
 daughter container

 
bool m_isStronglyDecayingResonance
 decay length less than 1 micron

 
const ConstraintConfigurationm_config
 has all the constraint config
 

Private Attributes

double m_bfield
 B field along z

 
const Belle2::TrackFitResultm_trackfit
 trackfit result from reconstruction

 
bool m_cached
 flag to mark the particle as initialised

 
double m_flt
 helix arc length at vertex
 
Eigen::Matrix< double, 1, 5 > m_params
 column vector to store the measurement
 
Eigen::Matrix< double, 5, 5 > m_covariance
 only lower triangle filled!
 
const float m_momentumScalingFactor
 scale the momenta by this correction factor
 
int m_index
 index
 
std::string m_name
 name

 

Detailed Description

representation of all charged final states FIXME rename since this name is taken in tracking

Definition at line 18 of file RecoTrack.h.

Member Typedef Documentation

◆ constraintlist

typedef std::vector<Constraint> constraintlist
inherited

alias

Definition at line 52 of file ParticleBase.h.

◆ indexmap

typedef std::vector< std::pair<const ParticleBase*, int> > indexmap
inherited

alias

Definition at line 55 of file ParticleBase.h.

◆ ParticleContainer

typedef std::vector<ParticleBase*> ParticleContainer
protectedinherited

just an alias

Definition at line 178 of file ParticleBase.h.

Member Enumeration Documentation

◆ TFParticleType

enum TFParticleType
inherited

particle types

Definition at line 30 of file ParticleBase.h.

30 {kInteractionPoint,
31 kOrigin,
32 kComposite,
33 kRecoResonance,
34 kInternalParticle,
35 kRecoTrack,
36 kResonance,
37 kRecoPhoton,
38 kRecoKlong,
39 kMissingParticle
40 };

Constructor & Destructor Documentation

◆ RecoTrack()

RecoTrack ( Belle2::Particle bc,
const ParticleBase mother 
)

constructor

Definition at line 23 of file RecoTrack.cc.

23 :
25 m_bfield(0),
27 m_cached(false),
28 m_flt(0),
29 m_params(5),
30 m_covariance(5, 5),
32 {
33 m_bfield = Belle2::BFieldManager::getFieldInTesla(ROOT::Math::XYZVector(0, 0, 0)).Z(); //Bz in Tesla
34 m_covariance = Eigen::Matrix<double, 5, 5>::Zero(5, 5);
35 }
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:61
double getEffectiveMomentumScale() const
Returns effective momentum scale which is the product of the momentum scaling and smearing factors.
Definition: Particle.h:614
const TrackFitResult * getTrackFitResult() const
Returns the pointer to the TrackFitResult that was used to create this Particle (ParticleType == c_Tr...
Definition: Particle.cc:854
Belle2::Particle * particle() const
get basf2 particle
Definition: ParticleBase.h:92
const ParticleBase * mother() const
getMother() / hasMother()
Definition: ParticleBase.h:98
RecoParticle(Belle2::Particle *bc, const ParticleBase *mother)
constructor
Definition: RecoParticle.cc:15
const Belle2::TrackFitResult * m_trackfit
trackfit result from reconstruction
Definition: RecoTrack.h:72
Eigen::Matrix< double, 1, 5 > m_params
column vector to store the measurement
Definition: RecoTrack.h:81
const float m_momentumScalingFactor
scale the momenta by this correction factor
Definition: RecoTrack.h:87
double m_bfield
B field along z
Definition: RecoTrack.h:69
double m_flt
helix arc length at vertex
Definition: RecoTrack.h:78
Eigen::Matrix< double, 5, 5 > m_covariance
only lower triangle filled!
Definition: RecoTrack.h:84
bool m_cached
flag to mark the particle as initialised
Definition: RecoTrack.h:75

◆ ~RecoTrack()

virtual ~RecoTrack ( )
inlinevirtual

destructor

Definition at line 25 of file RecoTrack.h.

25{};

Member Function Documentation

◆ addDaughter()

ParticleBase * addDaughter ( Belle2::Particle cand,
const ConstraintConfiguration config,
bool  forceFitAll = false 
)
virtualinherited

add daughter

Definition at line 65 of file ParticleBase.cc.

66 {
67 auto newDaughter = ParticleBase::createParticle(cand, this, config, forceFitAll);
68 m_daughters.push_back(newDaughter);
69 return m_daughters.back();
70 }
static ParticleBase * createParticle(Belle2::Particle *particle, const ParticleBase *mother, const ConstraintConfiguration &config, bool forceFitAll=false)
create the according treeFitter particle obj for a basf2 particle type
std::vector< ParticleBase * > m_daughters
daughter container
Definition: ParticleBase.h:198

◆ addToConstraintList()

virtual void addToConstraintList ( constraintlist alist,
int  depth 
) const
inlineoverridevirtual

add to the list of constraints

Implements ParticleBase.

Definition at line 55 of file RecoTrack.h.

56 {
57 alist.push_back(Constraint(this, Constraint::track, depth, dimM(), 1)) ;
58 }
virtual int dimM() const override
dimension (5)
Definition: RecoTrack.h:46

◆ charge()

int charge ( ) const
inlineinherited

get charge

Definition at line 144 of file ParticleBase.h.

145 {
146 if (m_particle->getPDGCode()) {
147 double fltcharge = m_particle->getCharge();
148 return fltcharge < 0 ? int(fltcharge - 0.5) : int(fltcharge + 0.5);
149 } else {
150 return m_particle->getCharge() > 0 ? 1 : (m_particle->getCharge() < 0 ? -1 : 0);
151 }
152 }
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:454
double getCharge(void) const
Returns particle charge.
Definition: Particle.cc:622
Belle2::Particle * m_particle
pointer to framework type
Definition: ParticleBase.h:192

◆ chiSquare()

double chiSquare ( const FitParams fitparams) const
virtualinherited

get chi2

Definition at line 231 of file ParticleBase.cc.

232 {
233 double rc = 0;
234 for (auto* daughter : m_daughters) {
235 rc += daughter->chiSquare(fitparams);
236 }
237 return rc;
238 }

◆ collectVertexDaughters()

void collectVertexDaughters ( std::vector< ParticleBase * > &  particles,
int  posindex 
)
inherited

get vertex daughters

Definition at line 156 of file ParticleBase.cc.

157 {
158 if (mother() && mother()->posIndex() == posindex) {
159 particles.push_back(this);
160 }
161
162 for (auto* daughter : m_daughters) {
163 daughter->collectVertexDaughters(particles, posindex);
164 }
165 }
virtual int posIndex() const
get vertex index (in statevector!)
Definition: ParticleBase.h:122

◆ createOrigin()

ParticleBase * createOrigin ( Belle2::Particle daughter,
const ConstraintConfiguration config,
bool  forceFitAll 
)
staticinherited

create a custom origin particle or a beamspot

Definition at line 93 of file ParticleBase.cc.

98 {
99 return new Origin(daughter, config, forceFitAll);
100 }

◆ createParticle()

ParticleBase * createParticle ( Belle2::Particle particle,
const ParticleBase mother,
const ConstraintConfiguration config,
bool  forceFitAll = false 
)
staticinherited

create the according treeFitter particle obj for a basf2 particle type

Definition at line 102 of file ParticleBase.cc.

104 {
105 ParticleBase* rc = nullptr;
106
107 if (!mother) { // If there is no mother, this is the 'head of tree' particle (is never a resonance)
108 rc = new InternalParticle(particle, nullptr, config, forceFitAll);
109 } else if (particle->hasExtraInfo("bremsCorrected") // Has Bremsstrahlungs-recovery
110 && particle->getExtraInfo("bremsCorrected") != 0) { // and gammas are attached
111 rc = new Composite(particle, mother, config, true);
112 // if no gamma is attached, it is treated as a RecoTrack
113 } else if (particle->hasExtraInfo("treeFitterTreatMeAsInvisible")
114 && particle->getExtraInfo("treeFitterTreatMeAsInvisible") == 1) { // dummy particles with invisible flag
115 rc = new RecoResonance(particle, mother, config);
116 } else if (particle->getTrack()) { // external reconstructed track
117 rc = new RecoTrack(particle, mother);
118 } else if (particle->getECLCluster()) { // external reconstructed photon
119 rc = new RecoPhoton(particle, mother);
120 } else if (particle->getKLMCluster()) { // external reconstructed klong
121 rc = new RecoKlong(particle, mother);
122 } else if (particle->getMdstArrayIndex()) { // external composite e.g. V0
123 rc = new InternalParticle(particle, mother, config, forceFitAll);
124 } else { // 'internal' particles
125 if (isAResonance(particle)) {
126 rc = new Resonance(particle, mother, config, forceFitAll);
127 } else {
128 rc = new InternalParticle(particle, mother, config, forceFitAll);
129 }
130 }
131 return rc;
132 }
const KLMCluster * getKLMCluster() const
Returns the pointer to the KLMCluster object that was used to create this Particle (ParticleType == c...
Definition: Particle.cc:926
const Track * getTrack() const
Returns the pointer to the Track object that was used to create this Particle (ParticleType == c_Trac...
Definition: Particle.cc:845
const ECLCluster * getECLCluster() const
Returns the pointer to the ECLCluster object that was used to create this Particle (if ParticleType =...
Definition: Particle.cc:891
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
Definition: Particle.cc:1266
unsigned getMdstArrayIndex(void) const
Returns 0-based index of MDST store array object (0 for composite particles)
Definition: Particle.h:487
double getExtraInfo(const std::string &name) const
Return given value if set.
Definition: Particle.cc:1289
static bool isAResonance(Belle2::Particle *particle)
controls if a particle is treated as a resonance(lifetime=0) or a particle that has a finite lifetime...
ParticleBase(Belle2::Particle *particle, const ParticleBase *mother, const ConstraintConfiguration *config=nullptr)
default constructor
Definition: ParticleBase.cc:30

◆ dim()

virtual int dim ( ) const
inlineoverridevirtualinherited

this here sets the size in the state vector we are only interested in the momenta of photons and tracks as the position were the track ends or the cluster is, is not relevant for physics

Implements ParticleBase.

Reimplemented in RecoKlong, and RecoPhoton.

Definition at line 39 of file RecoParticle.h.

39{ return 3; }

◆ dimM()

virtual int dimM ( ) const
inlineoverridevirtual

dimension (5)

Implements RecoParticle.

Definition at line 46 of file RecoTrack.h.

46{ return 5 ; }

◆ eneIndex()

int eneIndex ( ) const
inlineinherited

get energy index

Definition at line 138 of file ParticleBase.h.

138{ return hasEnergy() ? momIndex() + 3 : -1 ; }
virtual int momIndex() const
get momentum index
Definition: ParticleBase.h:128
virtual bool hasEnergy() const
get momentum dimension
Definition: ParticleBase.h:132

◆ forceP4Sum()

virtual void forceP4Sum ( FitParams ) const
inlinevirtualinherited

force p4 sum conservation all along the tree

Reimplemented in InternalParticle.

Definition at line 116 of file ParticleBase.h.

116{} ;

◆ hasEnergy()

virtual bool hasEnergy ( ) const
inlineoverridevirtualinherited

has an energy in the statevector?

Reimplemented from ParticleBase.

Reimplemented in RecoKlong, and RecoPhoton.

Definition at line 45 of file RecoParticle.h.

45{ return false; }

◆ hasPosition()

virtual bool hasPosition ( ) const
inlinevirtualinherited

get false

Reimplemented in Composite, InternalParticle, and Resonance.

Definition at line 135 of file ParticleBase.h.

135{ return false ; }

◆ index()

int index ( ) const
inlineinherited

get index

Definition at line 95 of file ParticleBase.h.

95{ return m_index ; }

◆ initCovariance()

ErrCode initCovariance ( FitParams fitparams) const
overridevirtual

init covariance matrix of this particle constraint

Reimplemented from ParticleBase.

Definition at line 56 of file RecoTrack.cc.

57 {
58 // we only need a rough estimate of the covariance
59 TMatrixFSym p4Err = particle()->getMomentumErrorMatrix();
60 const int momindex = momIndex();
61
62 for (int row = 0; row < 3; ++row) {
63 fitparams.getCovariance()(momindex + row, momindex + row) = 1000 * p4Err[row][row];
64 }
65
66 return ErrCode(ErrCode::Status::success);
67 }
TMatrixFSym getMomentumErrorMatrix() const
Returns the 4x4 momentum error matrix.
Definition: Particle.cc:435
virtual int momIndex() const override
get momentum index
Definition: RecoParticle.h:42

◆ initMotherlessParticle()

ErrCode initMotherlessParticle ( FitParams fitparams)
overridevirtual

init without mother particle

Reimplemented from RecoParticle.

Definition at line 51 of file RecoTrack.cc.

52 {
53 return ErrCode(ErrCode::Status::success);
54 }

◆ initParticleWithMother()

ErrCode initParticleWithMother ( FitParams fitparams)
overridevirtual

init with mother particle (replacing initPar2)

Implements ParticleBase.

Definition at line 37 of file RecoTrack.cc.

38 {
39 //initPar2
40 if (m_flt == 0) {
41 const_cast<RecoTrack*>(this)->updFltToMother(fitparams);
42 }
44 const int momindex = momIndex();
45 fitparams.getStateVector()(momindex) = recoP.X();
46 fitparams.getStateVector()(momindex + 1) = recoP.Y();
47 fitparams.getStateVector()(momindex + 2) = recoP.Z();
48 return ErrCode(ErrCode::Status::success);
49 }
ROOT::Math::XYZVector getMomentumAtArcLength2D(const double &arcLength2D, const double &bz) const
Calculates the momentum vector at the given two dimensional arc length.
Definition: Helix.cc:270
Helix getHelix() const
Conversion to framework Helix (without covariance).
ErrCode updFltToMother(const FitParams &fitparams)
update flight length to mother
Definition: RecoTrack.cc:70
RecoTrack(Belle2::Particle *bc, const ParticleBase *mother)
constructor
Definition: RecoTrack.cc:23

◆ initTau()

ErrCode initTau ( FitParams par) const
protectedinherited

initialises tau as a length

Definition at line 426 of file ParticleBase.cc.

427 {
428 const int tauindex = tauIndex();
429 if (tauindex >= 0 && hasPosition()) {
430
431 const int posindex = posIndex();
432 const int mother_ps_index = mother()->posIndex();
433 const int dim = m_config->m_originDimension; // TODO can we configure this to be particle specific?
434
435 // tau has different meaning depending on the dimension of the constraint
436 // 2-> use x-y projection
437 const Eigen::Matrix < double, 1, -1, 1, 1, 3 > vertex_dist =
438 fitparams.getStateVector().segment(posindex, dim) - fitparams.getStateVector().segment(mother_ps_index, dim);
439 const Eigen::Matrix < double, 1, -1, 1, 1, 3 >
440 mom = fitparams.getStateVector().segment(posindex, dim) - fitparams.getStateVector().segment(mother_ps_index, dim);
441
442 // if an intermediate vertex is not well defined by a track or so it will be initialised with 0
443 // same for the momentum of for example B0, it might be initialised with 0
444 // in those cases use pdg value
445 const double mom_norm = mom.norm();
446 const double dot = std::abs(vertex_dist.dot(mom));
447 const double tau = dot / mom_norm;
448 if (0 == mom_norm || 0 == dot) {
449 const double mass = m_particle->getPDGMass();
450 if (mass > 0)
451 fitparams.getStateVector()(tauindex) = m_particle->getPDGLifetime() * 1e9 * Belle2::Const::speedOfLight / mass;
452 else
453 fitparams.getStateVector()(tauindex) = 0;
454 } else {
455 fitparams.getStateVector()(tauindex) = tau;
456 }
457 }
458
459 return ErrCode(ErrCode::Status::success);
460 }
static const double speedOfLight
[cm/ns]
Definition: Const.h:695
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
Definition: Particle.cc:604
double getPDGLifetime() const
Returns particle nominal lifetime.
Definition: Particle.cc:613
const int m_originDimension
dimension of the origin constraint and ALL geometric gcosntraints
virtual int dim() const =0
get dimension of constraint
const ConstraintConfiguration * m_config
has all the constraint config
Definition: ParticleBase.h:204
virtual int tauIndex() const
get tau index
Definition: ParticleBase.h:125
virtual bool hasPosition() const
get false
Definition: ParticleBase.h:135

◆ isAResonance()

bool isAResonance ( Belle2::Particle particle)
staticprotectedinherited

controls if a particle is treated as a resonance(lifetime=0) or a particle that has a finite lifetime.

A finite life time means it will register a geo constraint for this particle

Definition at line 134 of file ParticleBase.cc.

135 {
136 bool rc = false ;
137 const int pdgcode = std::abs(particle->getPDGCode());
138
139 if (pdgcode && !(particle->getMdstArrayIndex())) {
140 switch (pdgcode) {
141 case 22: //photon conversion
142 rc = false;
143 break ;
144
145 case -11: //bremsstrahlung
146 case 11:
147 rc = true ;
148 break ;
149 default: //everything with boosted flight length less than 1 micrometer
150 rc = (pdgcode && particle->getPDGLifetime() < 1e-14);
151 }
152 }
153 return rc ;
154 }

◆ locate()

const ParticleBase * locate ( Belle2::Particle particle) const
inherited

get particle base from basf2 particle

Definition at line 211 of file ParticleBase.cc.

212 {
213 const ParticleBase* rc = (m_particle == particle) ? this : nullptr;
214 if (!rc) {
215 for (auto* daughter : m_daughters) {
216 rc = daughter->locate(particle);
217 if (rc) {break;}
218 }
219 }
220 return rc;
221 }

◆ momIndex()

virtual int momIndex ( ) const
inlineoverridevirtualinherited

get momentum index

Reimplemented from ParticleBase.

Definition at line 42 of file RecoParticle.h.

42{ return index(); }
int index() const
get index
Definition: ParticleBase.h:95

◆ mother()

const ParticleBase * mother ( ) const
inlineinherited

getMother() / hasMother()

Definition at line 98 of file ParticleBase.h.

98{ return m_mother; };
const ParticleBase * m_mother
motherparticle
Definition: ParticleBase.h:195

◆ nFinalChargedCandidates()

virtual int nFinalChargedCandidates ( ) const
inlineoverridevirtual

number of final charged candidates

Reimplemented from ParticleBase.

Definition at line 52 of file RecoTrack.h.

52{ return 1 ; }

◆ parname()

std::string parname ( int  index) const
overridevirtualinherited

name

Reimplemented from ParticleBase.

Definition at line 23 of file RecoParticle.cc.

24 {
25 return ParticleBase::parname(index + 4);
26 }
virtual std::string parname(int index) const
get name of parameter i

◆ particle()

Belle2::Particle * particle ( ) const
inlineinherited

get basf2 particle

Definition at line 92 of file ParticleBase.h.

92{ return m_particle ; }

◆ posIndex()

virtual int posIndex ( ) const
inlinevirtualinherited

get vertex index (in statevector!)

Reimplemented in Composite, InternalParticle, Origin, RecoResonance, and Resonance.

Definition at line 122 of file ParticleBase.h.

122{ return -1 ; }

◆ projectConstraint()

ErrCode projectConstraint ( Constraint::Type  type,
const FitParams fitparams,
Projection p 
) const
overridevirtualinherited

abstract abstract projection

Reimplemented from ParticleBase.

Definition at line 28 of file RecoParticle.cc.

29 {
30 ErrCode status ;
31 switch (type) {
32 case Constraint::track :
33 case Constraint::photon :
34 case Constraint::klong :
35 status |= projectRecoConstraint(fitparams, p);
36 break ;
37 default:
38 status |= ParticleBase::projectConstraint(type, fitparams, p);
39 }
40 return status;
41 }
virtual ErrCode projectConstraint(Constraint::Type, const FitParams &, Projection &) const
project constraint.
virtual int type() const =0
get particle type
virtual ErrCode projectRecoConstraint(const FitParams &fitparams, Projection &p) const =0
abstract projection

◆ projectGeoConstraint()

ErrCode projectGeoConstraint ( const FitParams fitparams,
Projection p 
) const
virtualinherited

project geometrical constraint

the direction of the momentum is very well known from the kinematic constraints that is why we do not extract the distance as a vector here

Definition at line 249 of file ParticleBase.cc.

250 {
251 assert(m_config);
252 // only allow 2d for head of tree particles that are beam constrained
253 const int dim = m_config->m_originDimension == 2 && std::abs(m_particle->getPDGCode()) == m_config->m_headOfTreePDG ? 2 : 3;
254 const int posindexmother = mother()->posIndex();
255 const int posindex = posIndex();
256 const int tauindex = tauIndex();
257 const int momindex = momIndex();
258
259 const double tau = fitparams.getStateVector()(tauindex);
260 Eigen::Matrix < double, 1, -1, 1, 1, 3 > x_vec = fitparams.getStateVector().segment(posindex, dim);
261 Eigen::Matrix < double, 1, -1, 1, 1, 3 > x_m = fitparams.getStateVector().segment(posindexmother, dim);
262 Eigen::Matrix < double, 1, -1, 1, 1, 3 > p_vec = fitparams.getStateVector().segment(momindex, dim);
263 const double mom = p_vec.norm();
264 const double mom3 = mom * mom * mom;
265
266 if (3 == dim) {
267 // we can already set these
268 //diagonal momentum
269 p.getH()(0, momindex) = tau * (p_vec(1) * p_vec(1) + p_vec(2) * p_vec(2)) / mom3 ;
270 p.getH()(1, momindex + 1) = tau * (p_vec(0) * p_vec(0) + p_vec(2) * p_vec(2)) / mom3 ;
271 p.getH()(2, momindex + 2) = tau * (p_vec(0) * p_vec(0) + p_vec(1) * p_vec(1)) / mom3 ;
272
273 //offdiagonal momentum
274 p.getH()(0, momindex + 1) = - tau * p_vec(0) * p_vec(1) / mom3 ;
275 p.getH()(0, momindex + 2) = - tau * p_vec(0) * p_vec(2) / mom3 ;
276
277 p.getH()(1, momindex + 0) = - tau * p_vec(1) * p_vec(0) / mom3 ;
278 p.getH()(1, momindex + 2) = - tau * p_vec(1) * p_vec(2) / mom3 ;
279
280 p.getH()(2, momindex + 0) = - tau * p_vec(2) * p_vec(0) / mom3 ;
281 p.getH()(2, momindex + 1) = - tau * p_vec(2) * p_vec(1) / mom3 ;
282
283 } else if (2 == dim) {
284
285 // NOTE THAT THESE ARE DIFFERENT IN 2d
286 p.getH()(0, momindex) = tau * (p_vec(1) * p_vec(1)) / mom3 ;
287 p.getH()(1, momindex + 1) = tau * (p_vec(0) * p_vec(0)) / mom3 ;
288
289 //offdiagonal momentum
290 p.getH()(0, momindex + 1) = - tau * p_vec(0) * p_vec(1) / mom3 ;
291 p.getH()(1, momindex + 0) = - tau * p_vec(1) * p_vec(0) / mom3 ;
292 } else {
293 B2FATAL("Dimension of Geometric constraint is not 2 or 3. This will crash many things. You should feel bad.");
294 }
295
296 for (int row = 0; row < dim; ++row) {
297
298 double posxmother = x_m(row);
299 double posx = x_vec(row);
300 double momx = p_vec(row);
301
305 p.getResiduals()(row) = posxmother + tau * momx / mom - posx ;
306 p.getH()(row, posindexmother + row) = 1;
307 p.getH()(row, posindex + row) = -1;
308 p.getH()(row, tauindex) = momx / mom;
309 }
310
311 return ErrCode(ErrCode::Status::success);
312 }
int m_headOfTreePDG
PDG code of the head particle.

◆ projectMassConstraint()

ErrCode projectMassConstraint ( const FitParams fitparams,
Projection p 
) const
virtualinherited

project mass constraint abstract

Definition at line 405 of file ParticleBase.cc.

407 {
408 assert(m_config);
409 if (m_config->m_massConstraintType == 0) {
410 return projectMassConstraintParticle(fitparams, p);
411 } else {
412 return projectMassConstraintDaughters(fitparams, p);
413 }
414 }
const bool m_massConstraintType
const flag for the type of the mass constraint
virtual ErrCode projectMassConstraintParticle(const FitParams &, Projection &) const
project mass constraint using the particles parameters
virtual ErrCode projectMassConstraintDaughters(const FitParams &, Projection &) const
project mass constraint using the parameters of the daughters

◆ projectMassConstraintDaughters()

ErrCode projectMassConstraintDaughters ( const FitParams fitparams,
Projection p 
) const
virtualinherited

project mass constraint using the parameters of the daughters

be aware that the signs here are important E-|p|-m extracts a negative mass and messes with the momentum !

Definition at line 314 of file ParticleBase.cc.

316 {
317 const double mass = particle()->getPDGMass();
318 const double mass2 = mass * mass;
319 double px = 0;
320 double py = 0;
321 double pz = 0;
322 double E = 0;
323
324 // the parameters of the daughters must be used otherwise the mass constraint does not have an effect on the extracted daughter momenta
325 for (const auto* daughter : m_daughters) {
326 const int momindex = daughter->momIndex();
327 // in most cases the daughters will be final states so we cache the value to use it in the energy column
328 const double px_daughter = fitparams.getStateVector()(momindex);
329 const double py_daughter = fitparams.getStateVector()(momindex + 1);
330 const double pz_daughter = fitparams.getStateVector()(momindex + 2);
331
332 px += px_daughter;
333 py += py_daughter;
334 pz += pz_daughter;
335 if (daughter->hasEnergy()) {
336 E += fitparams.getStateVector()(momindex + 3);
337 } else {
338 // final states dont have an energy index
339 const double m = daughter->particle()->getPDGMass();
340 E += std::sqrt(m * m + px_daughter * px_daughter + py_daughter * py_daughter + pz_daughter * pz_daughter);
341 }
342 }
343
347 p.getResiduals()(0) = mass2 - E * E + px * px + py * py + pz * pz;
348
349 for (const auto* daughter : m_daughters) {
350 //dr/dx = d/dx m2-{E1+E2+...}^2+{p1+p2+...}^2 = 2*x (x= E or p)
351 const int momindex = daughter->momIndex();
352 p.getH()(0, momindex) = 2.0 * px;
353 p.getH()(0, momindex + 1) = 2.0 * py;
354 p.getH()(0, momindex + 2) = 2.0 * pz;
355
356 if (daughter->hasEnergy()) {
357 p.getH()(0, momindex + 3) = -2.0 * E;
358 } else {
359 const double px_daughter = fitparams.getStateVector()(momindex);
360 const double py_daughter = fitparams.getStateVector()(momindex + 1);
361 const double pz_daughter = fitparams.getStateVector()(momindex + 2);
362 const double m = daughter->particle()->getPDGMass();
363
364 const double E_daughter = std::sqrt(m * m + px_daughter * px_daughter + py_daughter * py_daughter + pz_daughter * pz_daughter);
365 const double E_by_E_daughter = E / E_daughter;
366 p.getH()(0, momindex) -= 2.0 * E_by_E_daughter * px_daughter;
367 p.getH()(0, momindex + 1) -= 2.0 * E_by_E_daughter * py_daughter;
368 p.getH()(0, momindex + 2) -= 2.0 * E_by_E_daughter * pz_daughter;
369 }
370
371 }
372 return ErrCode(ErrCode::Status::success);
373 }

◆ projectMassConstraintParticle()

ErrCode projectMassConstraintParticle ( const FitParams fitparams,
Projection p 
) const
virtualinherited

project mass constraint using the particles parameters

be aware that the signs here are important E-|p|-m extracts a negative mass and messes with the momentum !

Definition at line 375 of file ParticleBase.cc.

377 {
378 const double mass = particle()->getPDGMass();
379 const double mass2 = mass * mass;
380 const int momindex = momIndex();
381 const double px = fitparams.getStateVector()(momindex);
382 const double py = fitparams.getStateVector()(momindex + 1);
383 const double pz = fitparams.getStateVector()(momindex + 2);
384 const double E = fitparams.getStateVector()(momindex + 3);
385
389 p.getResiduals()(0) = mass2 - E * E + px * px + py * py + pz * pz;
390
391 p.getH()(0, momindex) = 2.0 * px;
392 p.getH()(0, momindex + 1) = 2.0 * py;
393 p.getH()(0, momindex + 2) = 2.0 * pz;
394 p.getH()(0, momindex + 3) = -2.0 * E;
395
396 // TODO 0 in most cases -> needs special treatment if width=0 to not crash chi2 calculation
397 // const double width = TDatabasePDG::Instance()->GetParticle(particle()->getPDGCode())->Width();
398 // transport measurement uncertainty into residual system
399 // f' = sigma_x^2 * (df/dx)^2
400 // p.getV()(0) = width * width * 4 * mass2;
401
402 return ErrCode(ErrCode::Status::success);
403 }

◆ projectRecoConstraint()

ErrCode projectRecoConstraint ( const FitParams fitparams,
Projection p 
) const
overridevirtual

project the constraint (calculate residuals)

Implements RecoParticle.

Definition at line 98 of file RecoTrack.cc.

99 {
100 ErrCode status;
101 const int posindexmother = mother()->posIndex();
102 const int momindex = momIndex();
103
104 Eigen::Matrix<double, 1, 6> positionAndMom = Eigen::Matrix<double, 1, 6>::Zero(1, 6);
105 positionAndMom.segment(0, 3) = fitparams.getStateVector().segment(posindexmother, 3);
106 positionAndMom.segment(3, 3) = fitparams.getStateVector().segment(momindex, 3);
107 Eigen::Matrix<double, 5, 6> jacobian = Eigen::Matrix<double, 5, 6>::Zero(5, 6);
108
109 const Belle2::Helix helix = Belle2::Helix(ROOT::Math::XYZVector(positionAndMom(0), positionAndMom(1), positionAndMom(2)),
110 ROOT::Math::XYZVector(positionAndMom(3), positionAndMom(4), positionAndMom(5)),
111 charge(), m_bfield);
112
114 positionAndMom(0),
115 positionAndMom(1),
116 positionAndMom(2),
117 positionAndMom(3),
118 positionAndMom(4),
119 positionAndMom(5),
120 m_bfield,
121 charge()
122 );
123 if (!m_cached) {
124 auto* nonconst = const_cast<RecoTrack*>(this);
125 if (m_flt == 0) { nonconst->updFltToMother(fitparams); }
126 nonconst->updateParams(m_flt);
127 }
128
129 Eigen::Matrix<double, 1, 5> helixpars(5);
130 helixpars(0) = helix.getD0();
131 helixpars(1) = helix.getPhi0();
132 helixpars(2) = helix.getOmega();
133 helixpars(3) = helix.getZ0();
134 helixpars(4) = helix.getTanLambda();
135
136 p.getResiduals().segment(0, 5) = m_params - helixpars;
137
138 //account for periodic boundary in phi residual
139 double phiResidual = p.getResiduals().segment(0, 5)(1);
140 phiResidual = std::fmod(phiResidual + pi, twoPi);
141 if (phiResidual < 0) phiResidual += twoPi;
142 phiResidual -= pi;
143 p.getResiduals().segment(0, 5)(1) = phiResidual;
144
145 p.getV().triangularView<Eigen::Lower>() = m_covariance.triangularView<Eigen::Lower>();
146
147 p.getH().block<5, 3>(0, posindexmother) = -1.0 * jacobian.block<5, 3>(0, 0);
148 p.getH().block<5, 3>(0, momindex) = -1.0 * jacobian.block<5, 3>(0, 3);
149
150 return status;
151 }
This class represents an ideal helix in perigee parameterization.
Definition: Helix.h:59
double getOmega() const
Getter for omega, which is a signed curvature measure of the track.
Definition: Helix.h:387
double getD0() const
Getter for d0, which is the signed distance to the perigee in the r-phi plane.
Definition: Helix.h:372
double getTanLambda() const
Getter for tan lambda, which is the z over two dimensional arc length slope of the track.
Definition: Helix.h:393
double getZ0() const
Getter for z0, which is the z coordinate of the perigee.
Definition: Helix.h:390
double getPhi0() const
Getter for phi0, which is the azimuth angle of the transverse momentum at the perigee.
Definition: Helix.h:378
static void getJacobianToCartesianFrameworkHelix(Eigen::Matrix< double, 5, 6 > &jacobian, const double x, const double y, const double z, const double px, const double py, const double pz, const double bfield, const double charge)
get the jacobian dh={helix pars}/dx={x,y,z,px,py,pz} for the implementation of the framework helix.
Definition: HelixUtils.cc:407
int charge() const
get charge
Definition: ParticleBase.h:144

◆ removeDaughter()

void removeDaughter ( const ParticleBase pb)
virtualinherited

remove daughter

Definition at line 73 of file ParticleBase.cc.

74 {
75 auto iter = std::find(m_daughters.begin(), m_daughters.end(), pb);
76 if (iter != m_daughters.end()) {
77 delete *iter;
78 m_daughters.erase(iter);
79 } else {
80 B2ERROR("Cannot remove particle, because not found ...");
81 }
82 }

◆ retrieveIndexMap()

void retrieveIndexMap ( indexmap anindexmap) const
virtualinherited

get index map

Definition at line 223 of file ParticleBase.cc.

224 {
225 map.push_back(std::pair<const ParticleBase*, int>(this, index()));
226 for (auto* daughter : m_daughters) {
227 daughter->retrieveIndexMap(map);
228 }
229 }

◆ setFlightLength()

void setFlightLength ( double  flt)
inline

setter for the flight length

Definition at line 64 of file RecoTrack.h.

64{ m_flt = flt ; }

◆ setIndex()

void setIndex ( int  i)
inlineprotectedinherited

set Index (in statevector)

Definition at line 189 of file ParticleBase.h.

189{ m_index = i ; }

◆ setMother()

void setMother ( const ParticleBase m)
inlineinherited

set mother

Definition at line 164 of file ParticleBase.h.

164{ m_mother = m ; }

◆ tauIndex()

virtual int tauIndex ( ) const
inlinevirtualinherited

get tau index

Reimplemented in Composite, InternalParticle, Origin, RecoResonance, and Resonance.

Definition at line 125 of file ParticleBase.h.

125{ return -1 ; }

◆ type()

virtual int type ( ) const
inlineoverridevirtual

type of the constraint

Implements ParticleBase.

Definition at line 49 of file RecoTrack.h.

49{ return kRecoTrack ; }

◆ updateIndex()

void updateIndex ( int &  offset)
virtualinherited

this sets the index for momentum, position, etc.

in the statevector

Definition at line 84 of file ParticleBase.cc.

85 {
86 for (auto* daughter : m_daughters) {
87 daughter->updateIndex(offset);
88 }
89 m_index = offset;
90 offset += dim();
91 }

◆ updateParams()

ErrCode updateParams ( double  flt)

updated the cached parameters

Definition at line 79 of file RecoTrack.cc.

80 {
81 m_flt = flt;
82 std::vector<float> trackParameter = m_trackfit->getTau();
83 for (unsigned int i = 0; i < trackParameter.size(); ++i) {
84 m_params(i) = trackParameter[i];
85 if (i == 2) m_params(2) /= m_momentumScalingFactor; // index 2 is the curvature of the track
86 }
87 TMatrixDSym cov = m_trackfit->getCovariance5();
88 for (int row = 0; row < 5; ++row) {
89 for (int col = 0; col <= row; ++col) {
90 m_covariance(row, col) = cov[row][col];
91 }
92 }
93
94 m_cached = true;
95 return ErrCode(ErrCode::Status::success);
96 }
TMatrixDSym getCovariance5() const
Getter for covariance matrix of perigee parameters in matrix form.
std::vector< float > getTau() const
Getter for all perigee parameters.

◆ updFltToMother()

ErrCode updFltToMother ( const FitParams fitparams)

update flight length to mother

Definition at line 70 of file RecoTrack.cc.

71 {
72 const int posindexmother = mother()->posIndex();
74 fitparams.getStateVector()(posindexmother),
75 fitparams.getStateVector()(posindexmother + 1));
76 return ErrCode(ErrCode::Status::success);
77 } ;
double getArcLength2DAtXY(const double &x, const double &y) const
Calculates the two dimensional arc length at which the circle in the xy projection is closest to the ...
Definition: Helix.cc:141

Member Data Documentation

◆ m_bfield

double m_bfield
private

B field along z

Definition at line 69 of file RecoTrack.h.

◆ m_cached

bool m_cached
private

flag to mark the particle as initialised

Definition at line 75 of file RecoTrack.h.

◆ m_config

const ConstraintConfiguration* m_config
protectedinherited

has all the constraint config

Definition at line 204 of file ParticleBase.h.

◆ m_covariance

Eigen::Matrix<double, 5, 5> m_covariance
private

only lower triangle filled!

Definition at line 84 of file RecoTrack.h.

◆ m_daughters

std::vector<ParticleBase*> m_daughters
protectedinherited

daughter container

Definition at line 198 of file ParticleBase.h.

◆ m_flt

double m_flt
private

helix arc length at vertex

Definition at line 78 of file RecoTrack.h.

◆ m_index

int m_index
privateinherited

index

Definition at line 208 of file ParticleBase.h.

◆ m_isStronglyDecayingResonance

bool m_isStronglyDecayingResonance
protectedinherited

decay length less than 1 micron

Definition at line 201 of file ParticleBase.h.

◆ m_momentumScalingFactor

const float m_momentumScalingFactor
private

scale the momenta by this correction factor

Definition at line 87 of file RecoTrack.h.

◆ m_mother

const ParticleBase* m_mother
protectedinherited

motherparticle

Definition at line 195 of file ParticleBase.h.

◆ m_name

std::string m_name
privateinherited

name

Definition at line 211 of file ParticleBase.h.

◆ m_params

Eigen::Matrix<double, 1, 5> m_params
private

column vector to store the measurement

Definition at line 81 of file RecoTrack.h.

◆ m_particle

Belle2::Particle* m_particle
protectedinherited

pointer to framework type

Definition at line 192 of file ParticleBase.h.

◆ m_trackfit

const Belle2::TrackFitResult* m_trackfit
private

trackfit result from reconstruction

Definition at line 72 of file RecoTrack.h.


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