Belle II Software development
Composite Class Reference

A class for composite particles, where the daughters must be ignored by the fitter. More...

#include <Composite.h>

Inheritance diagram for Composite:
ParticleBase RecoResonance

Public Types

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

Public Member Functions

 Composite (Belle2::Particle *bc, const ParticleBase *mother, const ConstraintConfiguration &config, bool massconstraint=false)
 constructor
 
virtual ~Composite ()
 destructor
 
virtual ErrCode initParticleWithMother (FitParams &fitparams) override
 init particle in case it has a mother
 
virtual ErrCode initMotherlessParticle (FitParams &fitparams) override
 init particle in case it has no mother
 
void updateParams ()
 update changed params
 
ErrCode projectComposite (const FitParams &fitparams, Projection &p) const
 project this particle constraint
 
virtual int dim () const override
 get dimension of constraint
 
int dimM () const
 get dimension of measurement
 
virtual ErrCode projectConstraint (Constraint::Type, const FitParams &, Projection &) const override
 project this constraint
 
virtual int type () const override
 get type
 
virtual int posIndex () const override
 get position index in statevectof x,y,z,tau,px,py,pz
 
virtual int tauIndex () const override
 get tau (lifetime) index in statevector
 
virtual int momIndex () const override
 get momentum index in statevector
 
virtual bool hasEnergy () const override
 return of this constraint/particle has an energy component
 
virtual bool hasPosition () const override
 return true FIXME
 
virtual void addToConstraintList (constraintlist &alist, int depth) const override
 add this to list
 
virtual ErrCode initCovariance (FitParams &) const
 init covariance matrix
 
virtual void updateIndex (int &offset)
 this sets the index for momentum, position, etc.
 
virtual std::string parname (int index) const
 get name of parameter i
 
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
 
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
 
virtual int nFinalChargedCandidates () const
 number of charged candidates
 

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

Eigen::Matrix< double, 7, 1 > m_params
 column vector to store the measurement
 
Eigen::Matrix< double, -1, -1, 0, 7, 7 > m_covariance
 only lower triangle filled!
 
bool m_hasEnergy
 flag
 
bool m_massconstraint
 flag
 
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

int m_index
 index
 
std::string m_name
 name
 

Detailed Description

A class for composite particles, where the daughters must be ignored by the fitter.

Currently used for bremsstrahlung correction and dummy invisibles

Definition at line 19 of file Composite.h.

Member Typedef Documentation

◆ constraintlist

typedef std::vector<Constraint> constraintlist
inherited

alias

Definition at line 48 of file ParticleBase.h.

◆ indexmap

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

alias

Definition at line 51 of file ParticleBase.h.

◆ ParticleContainer

typedef std::vector<ParticleBase*> ParticleContainer
protectedinherited

just an alias

Definition at line 174 of file ParticleBase.h.

Member Enumeration Documentation

◆ TFParticleType

enum TFParticleType
inherited

particle types

Definition at line 27 of file ParticleBase.h.

27 {kInteractionPoint,
28 kOrigin,
29 kComposite,
30 kRecoResonance,
31 kInternalParticle,
32 kRecoTrack,
33 kResonance,
34 kRecoNeutral,
35 kMissingParticle
36 };

Constructor & Destructor Documentation

◆ Composite()

Composite ( Belle2::Particle * bc,
const ParticleBase * mother,
const ConstraintConfiguration & config,
bool massconstraint = false )

constructor

Definition at line 20 of file Composite.cc.

22 : ParticleBase(particle, mother, &config), m_params(), m_hasEnergy(true), m_massconstraint(massConstraint)
23 {
24 updateParams();
25 }

◆ ~Composite()

virtual ~Composite ( )
inlinevirtual

destructor

Definition at line 27 of file Composite.h.

27{};

Member Function Documentation

◆ addDaughter()

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

add daughter

Definition at line 63 of file ParticleBase.cc.

64 {
65 auto newDaughter = ParticleBase::createParticle(cand, this, config, forceFitAll);
66 m_daughters.push_back(newDaughter);
67 return m_daughters.back();
68 }

◆ addToConstraintList()

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

add this to list

Implements ParticleBase.

Reimplemented in RecoResonance.

Definition at line 69 of file Composite.h.

70 {
71 alist.push_back(Constraint(this, Constraint::composite, depth, dimM())) ;
72 alist.push_back(Constraint(this, Constraint::geometric, depth, 3)) ;
73 if (m_massconstraint) {
74 alist.push_back(Constraint(this, Constraint::mass, depth, 1, 3));
75 }
76 }

◆ charge()

int charge ( ) const
inlineinherited

get charge

Definition at line 140 of file ParticleBase.h.

141 {
142 if (m_particle->getPDGCode()) {
143 double fltcharge = m_particle->getCharge();
144 return fltcharge < 0 ? int(fltcharge - 0.5) : int(fltcharge + 0.5);
145 } else {
146 return m_particle->getCharge() > 0 ? 1 : (m_particle->getCharge() < 0 ? -1 : 0);
147 }
148 }

◆ chiSquare()

double chiSquare ( const FitParams & fitparams) const
virtualinherited

get chi2

Definition at line 229 of file ParticleBase.cc.

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

◆ collectVertexDaughters()

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

get vertex daughters

Definition at line 154 of file ParticleBase.cc.

155 {
156 if (mother() && mother()->posIndex() == posindex) {
157 particles.push_back(this);
158 }
159
160 for (auto* daughter : m_daughters) {
161 daughter->collectVertexDaughters(particles, posindex);
162 }
163 }

◆ createOrigin()

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

create a custom origin particle or a beamspot

Definition at line 91 of file ParticleBase.cc.

96 {
97 return new Origin(daughter, config, forceFitAll);
98 }

◆ 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 100 of file ParticleBase.cc.

102 {
103 ParticleBase* rc = nullptr;
104
105 if (!mother) { // If there is no mother, this is the 'head of tree' particle (is never a resonance)
106 rc = new InternalParticle(particle, nullptr, config, forceFitAll);
107 } else if (particle->hasExtraInfo("bremsCorrected") // Has Bremsstrahlungs-recovery
108 && particle->getExtraInfo("bremsCorrected") != 0) { // and gammas are attached
109 rc = new Composite(particle, mother, config, true);
110 // if no gamma is attached, it is treated as a RecoTrack
111 } else if (particle->hasExtraInfo("treeFitterTreatMeAsInvisible")
112 && particle->getExtraInfo("treeFitterTreatMeAsInvisible") == 1) { // dummy particles with invisible flag
113 rc = new RecoResonance(particle, mother, config);
114 } else if (particle->getTrack()) { // external reconstructed track
115 rc = new RecoTrack(particle, mother);
116 } else if (particle->getParticleSource() == Belle2::Particle::EParticleSourceObject::c_ECLCluster
117 or particle->getParticleSource() == Belle2::Particle::EParticleSourceObject::c_KLMCluster) {
118 // external reconstructed neutral final-state particle
119 rc = new RecoNeutral(particle, mother);
120 } else if (particle->getMdstArrayIndex()) { // external composite e.g. V0
121 rc = new InternalParticle(particle, mother, config, forceFitAll);
122 } else { // 'internal' particles
123 if (isAResonance(particle)) {
124 rc = new Resonance(particle, mother, config, forceFitAll);
125 } else {
126 rc = new InternalParticle(particle, mother, config, forceFitAll);
127 }
128 }
129 return rc;
130 }

◆ dim()

virtual int dim ( ) const
inlineoverridevirtual

get dimension of constraint

Implements ParticleBase.

Reimplemented in RecoResonance.

Definition at line 42 of file Composite.h.

42{ return m_hasEnergy ? 8 : 7 ; }// (x,y,z,t,px,py,pz,(E))

◆ dimM()

int dimM ( ) const
inline

get dimension of measurement

Definition at line 45 of file Composite.h.

45{ return m_hasEnergy ? 7 : 6 ; }

◆ eneIndex()

int eneIndex ( ) const
inlineinherited

get energy index

Definition at line 134 of file ParticleBase.h.

134{ return hasEnergy() ? momIndex() + 3 : -1 ; }

◆ forceP4Sum()

virtual void forceP4Sum ( FitParams & ) const
inlinevirtualinherited

force p4 sum conservation all along the tree

Reimplemented in InternalParticle.

Definition at line 112 of file ParticleBase.h.

112{} ;

◆ hasEnergy()

virtual bool hasEnergy ( ) const
inlineoverridevirtual

return of this constraint/particle has an energy component

Reimplemented from ParticleBase.

Definition at line 63 of file Composite.h.

63{ return m_hasEnergy ; }

◆ hasPosition()

virtual bool hasPosition ( ) const
inlineoverridevirtual

return true FIXME

Reimplemented from ParticleBase.

Definition at line 66 of file Composite.h.

66{ return true ; }

◆ index()

int index ( ) const
inlineinherited

get index

Definition at line 91 of file ParticleBase.h.

91{ return m_index ; }

◆ initCovariance()

ErrCode initCovariance ( FitParams & fitparams) const
virtualinherited

init covariance matrix

Reimplemented in InternalParticle, Origin, RecoNeutral, and RecoTrack.

Definition at line 165 of file ParticleBase.cc.

166 {
167 // this is very sensitive and can heavily affect the efficiency of the fit
168 ErrCode status;
169
170 const int posindex = posIndex();
171 if (posindex >= 0) {
172 for (int i = 0; i < 3; ++i) {
173 fitparams.getCovariance()(posindex + i, posindex + i) = 1;
174 }
175 }
176
177 const int momindex = momIndex();
178 if (momindex >= 0) {
179 const int maxrow = hasEnergy() ? 4 : 3;
180 for (int i = 0; i < maxrow; ++i) {
181 fitparams.getCovariance()(momindex + i, momindex + i) = 10.;
182 }
183 }
184
185 const int tauindex = tauIndex();
186 if (tauindex >= 0) {
187 fitparams.getCovariance()(tauindex, tauindex) = 1.;
188 }
189 return status;
190 }

◆ initMotherlessParticle()

ErrCode initMotherlessParticle ( FitParams & fitparams)
overridevirtual

init particle in case it has no mother

Implements ParticleBase.

Reimplemented in RecoResonance.

Definition at line 32 of file Composite.cc.

33 {
34 const int posindex = posIndex();
35 const int momindex = momIndex();
36
37 fitparams.getStateVector().segment(posindex, 3) = m_params.segment(0, 3);
38
39 fitparams.getStateVector().segment(momindex, 4) = m_params.segment(3, 4);
40
41 return ErrCode(ErrCode::Status::success) ;
42 }

◆ initParticleWithMother()

ErrCode initParticleWithMother ( FitParams & fitparams)
overridevirtual

init particle in case it has a mother

Implements ParticleBase.

Reimplemented in RecoResonance.

Definition at line 27 of file Composite.cc.

28 {
29 return initTau(fitparams);
30 }

◆ initTau()

ErrCode initTau ( FitParams & par) const
protectedinherited

initialises tau as a length

Definition at line 428 of file ParticleBase.cc.

429 {
430 const int tauindex = tauIndex();
431 if (tauindex >= 0 && hasPosition()) {
432
433 const int posindex = posIndex();
434 const int mother_ps_index = mother()->posIndex();
435 const int dim = m_config->m_originDimension; // TODO can we configure this to be particle specific?
436
437 // tau has different meaning depending on the dimension of the constraint
438 // 2-> use x-y projection
439 const Eigen::Matrix < double, 1, -1, 1, 1, 3 > vertex_dist =
440 fitparams.getStateVector().segment(posindex, dim) - fitparams.getStateVector().segment(mother_ps_index, dim);
441 const Eigen::Matrix < double, 1, -1, 1, 1, 3 >
442 mom = fitparams.getStateVector().segment(posindex, dim);
443
444 // if an intermediate vertex is not well defined by a track or so it will be initialised with 0
445 // same for the momentum of for example B0, it might be initialised with 0
446 // in those cases use pdg value
447 const double mom_norm = mom.norm();
448 const double dot = std::abs(vertex_dist.dot(mom));
449 const double tau = dot / mom_norm;
450 if (0 == mom_norm || 0 == dot) {
451 const double mass = m_particle->getPDGMass();
452 if (mass > 0)
453 fitparams.getStateVector()(tauindex) = m_particle->getPDGLifetime() * 1e9 * Belle2::Const::speedOfLight / mass;
454 else
455 fitparams.getStateVector()(tauindex) = 0;
456 } else {
457 fitparams.getStateVector()(tauindex) = tau;
458 }
459 }
460
461 return ErrCode(ErrCode::Status::success);
462 }
static const double speedOfLight
[cm/ns]
Definition Const.h:695
T dot(GeneralVector< T > a, GeneralVector< T > b)
dot product of two general vectors

◆ 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 132 of file ParticleBase.cc.

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

◆ locate()

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

get particle base from basf2 particle

Definition at line 209 of file ParticleBase.cc.

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

◆ momIndex()

virtual int momIndex ( ) const
inlineoverridevirtual

get momentum index in statevector

Reimplemented from ParticleBase.

Reimplemented in RecoResonance.

Definition at line 60 of file Composite.h.

60{ return index() + 4 ; }

◆ mother()

const ParticleBase * mother ( ) const
inlineinherited

getMother() / hasMother()

Definition at line 94 of file ParticleBase.h.

94{ return m_mother; };

◆ nFinalChargedCandidates()

int nFinalChargedCandidates ( ) const
virtualinherited

number of charged candidates

Reimplemented in RecoTrack.

Definition at line 238 of file ParticleBase.cc.

239 {
240 int rc = 0;
241 for (auto* daughter : m_daughters) {
242 rc += daughter->nFinalChargedCandidates();
243 }
244 return rc;
245 }

◆ parname()

std::string parname ( int index) const
virtualinherited

get name of parameter i

Reimplemented in InternalParticle, RecoParticle, RecoResonance, and Resonance.

Definition at line 192 of file ParticleBase.cc.

193 {
194 std::string rc = m_name;
195 switch (thisindex) {
196 case 0: rc += "_x "; break;
197 case 1: rc += "_y "; break;
198 case 2: rc += "_z "; break;
199 case 3: rc += "_tau"; break;
200 case 4: rc += "_px "; break;
201 case 5: rc += "_py "; break;
202 case 6: rc += "_pz "; break;
203 case 7: rc += "_E "; break;
204 default:;
205 }
206 return rc;
207 }

◆ particle()

Belle2::Particle * particle ( ) const
inlineinherited

get basf2 particle

Definition at line 88 of file ParticleBase.h.

88{ return m_particle ; }

◆ posIndex()

virtual int posIndex ( ) const
inlineoverridevirtual

get position index in statevectof x,y,z,tau,px,py,pz

Reimplemented from ParticleBase.

Reimplemented in RecoResonance.

Definition at line 54 of file Composite.h.

54{ return index() ; }

◆ projectComposite()

ErrCode projectComposite ( const FitParams & fitparams,
Projection & p ) const

project this particle constraint

Definition at line 76 of file Composite.cc.

77 {
78 const int posindex = posIndex() ;
79 const int momindex = momIndex() ;
80
81 p.getResiduals().segment(0, 3) = m_params.segment(0, 3) - fitparams.getStateVector().segment(posindex, 3);
82
83 p.getResiduals().segment(3, 4) = m_params.segment(3, 4) - fitparams.getStateVector().segment(momindex, 4);
84
85 for (int row = 0; row < 3; ++row) {
86 p.getH()(row, posindex + row) = -1;
87 }
88
89 for (int row = 0; row < 4; ++row) {
90 p.getH()(3 + row, momindex + row) = -1;
91 }
92
93 p.getV().triangularView<Eigen::Lower>() = m_covariance.triangularView<Eigen::Lower>();
94
95 return ErrCode(ErrCode::Status::success) ;
96 }

◆ projectConstraint()

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

project this constraint

Reimplemented from ParticleBase.

Reimplemented in RecoResonance.

Definition at line 98 of file Composite.cc.

99 {
100 ErrCode status;
101 switch (type) {
102 case Constraint::mass:
103 status |= projectMassConstraint(fitparams, p);
104 break;
105 case Constraint::composite:
106 status |= projectComposite(fitparams, p);
107 break ;
108 case Constraint::geometric:
109 status |= projectGeoConstraint(fitparams, p);
110 break ;
111 default:
112 status |= ParticleBase::projectConstraint(type, fitparams, p);
113 }
114 return status ;
115 }

◆ 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 247 of file ParticleBase.cc.

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

◆ projectMassConstraint()

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

project mass constraint abstract

Definition at line 407 of file ParticleBase.cc.

409 {
410 assert(m_config);
411 if (m_config->m_massConstraintType == 0) {
412 return projectMassConstraintParticle(fitparams, p);
413 } else {
414 return projectMassConstraintDaughters(fitparams, p);
415 }
416 }

◆ 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 310 of file ParticleBase.cc.

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

◆ 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 double mass = 0;
379 if (particle()->hasExtraInfo("treeFitterMassConstraintValue")) mass = particle()->getExtraInfo("treeFitterMassConstraintValue");
380 else mass = particle()->getPDGMass();
381 const double mass2 = mass * mass;
382 const int momindex = momIndex();
383 Eigen::Vector4d state = fitparams.getStateVector().segment<4>(momindex);
384 Eigen::Vector3d p3 = state.head<3>();
385 const double px = state(0);
386 const double py = state(1);
387 const double pz = state(2);
388 const double E = state(3);
389
393 p.getResiduals()(0) = mass2 - E * E + px * px + py * py + pz * pz;
394
395 p.getH().block<1, 3>(0, momindex) = 2.0 * p3;
396 p.getH()(0, momindex + 3) = -2.0 * E;
397
398 // TODO 0 in most cases -> needs special treatment if width=0 to not crash chi2 calculation
399 // const double width = TDatabasePDG::Instance()->GetParticle(particle()->getPDGCode())->Width();
400 // transport measurement uncertainty into residual system
401 // f' = sigma_x^2 * (df/dx)^2
402 // p.getV()(0) = width * width * 4 * mass2;
403
404 return ErrCode(ErrCode::Status::success);
405 }

◆ removeDaughter()

void removeDaughter ( const ParticleBase * pb)
virtualinherited

remove daughter

Definition at line 71 of file ParticleBase.cc.

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

◆ retrieveIndexMap()

void retrieveIndexMap ( indexmap & anindexmap) const
virtualinherited

get index map

Definition at line 221 of file ParticleBase.cc.

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

◆ setIndex()

void setIndex ( int i)
inlineprotectedinherited

set Index (in statevector)

Definition at line 185 of file ParticleBase.h.

185{ m_index = i ; }

◆ setMother()

void setMother ( const ParticleBase * m)
inlineinherited

set mother

Definition at line 160 of file ParticleBase.h.

160{ m_mother = m ; }

◆ tauIndex()

virtual int tauIndex ( ) const
inlineoverridevirtual

get tau (lifetime) index in statevector

Reimplemented from ParticleBase.

Reimplemented in RecoResonance.

Definition at line 57 of file Composite.h.

57{ return index() + 3 ; }

◆ type()

virtual int type ( ) const
inlineoverridevirtual

get type

Implements ParticleBase.

Reimplemented in RecoResonance.

Definition at line 51 of file Composite.h.

51{ return kComposite ; }

◆ updateIndex()

void updateIndex ( int & offset)
virtualinherited

this sets the index for momentum, position, etc.

in the statevector

Definition at line 82 of file ParticleBase.cc.

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

◆ updateParams()

void updateParams ( )

update changed params

Definition at line 44 of file Composite.cc.

45 {
46 const ROOT::Math::XYZVector pos = particle()->getVertex();
47 const ROOT::Math::XYZVector mom = particle()->getMomentum();
48 const double energy = particle()->getEnergy();
49
50 m_params = Eigen::Matrix<double, 7, 1>::Zero(7, 1);
51 m_params(0) = pos.X();
52 m_params(1) = pos.Y();
53 m_params(2) = pos.Z();
54 m_params(3) = mom.X();
55 m_params(4) = mom.Y();
56 m_params(5) = mom.Z();
57
58 m_params(6) = energy;
59
60 m_covariance = Eigen::Matrix < double, -1, -1, 0, 7, 7 >::Zero(7, 7);
61 const TMatrixFSym cov7in = particle()->getMomentumVertexErrorMatrix(); //this is (p,E,x)
62
63 for (int row = 0; row < 4; ++row) { //first the p,E block
64 for (int col = 0; col <= row; ++col) {
65 m_covariance(3 + row, 3 + col) = cov7in[row][col];
66 }
67 }
68
69 for (int row = 0; row < 3; ++row) { //then the x block
70 for (int col = 0; col <= row; ++col) {
71 m_covariance(row, col) = cov7in[4 + row][4 + col];
72 }
73 }
74 }

Member Data Documentation

◆ m_config

const ConstraintConfiguration* m_config
protectedinherited

has all the constraint config

Definition at line 200 of file ParticleBase.h.

◆ m_covariance

Eigen::Matrix< double, -1, -1, 0, 7, 7 > m_covariance
protected

only lower triangle filled!

Definition at line 84 of file Composite.h.

◆ m_daughters

std::vector<ParticleBase*> m_daughters
protectedinherited

daughter container

Definition at line 194 of file ParticleBase.h.

◆ m_hasEnergy

bool m_hasEnergy
protected

flag

Definition at line 87 of file Composite.h.

◆ m_index

int m_index
privateinherited

index

Definition at line 204 of file ParticleBase.h.

◆ m_isStronglyDecayingResonance

bool m_isStronglyDecayingResonance
protectedinherited

decay length less than 1 micron

Definition at line 197 of file ParticleBase.h.

◆ m_massconstraint

bool m_massconstraint
protected

flag

Definition at line 90 of file Composite.h.

◆ m_mother

const ParticleBase* m_mother
protectedinherited

motherparticle

Definition at line 191 of file ParticleBase.h.

◆ m_name

std::string m_name
privateinherited

name

Definition at line 207 of file ParticleBase.h.

◆ m_params

Eigen::Matrix<double, 7, 1> m_params
protected

column vector to store the measurement

Definition at line 81 of file Composite.h.

◆ m_particle

Belle2::Particle* m_particle
protectedinherited

pointer to framework type

Definition at line 188 of file ParticleBase.h.


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