Belle II Software development
RecoPhoton Class Reference

representation of the photon constraint More...

#include <RecoPhoton.h>

Inheritance diagram for RecoPhoton:
RecoParticle ParticleBase

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

 RecoPhoton (Belle2::Particle *bc, const ParticleBase *mother)
 constructor
 
virtual ~RecoPhoton ()
 destructor
 
virtual ErrCode initParticleWithMother (FitParams &fitparams) override
 init particle with mother
 
virtual ErrCode initMotherlessParticle (FitParams &fitparams) override
 init particle without mother
 
ErrCode initCovariance (FitParams &fitparams) const override
 init covariance
 
ErrCode initParams ()
 update or init params
 
ErrCode projectRecoConstraint (const FitParams &fitparams, Projection &p) const override
 project photon constraint
 
virtual int dimM () const override
 sets the size of the corresponding residual projection
 
virtual bool hasEnergy () const override
 how should the energy be calculated ?
 
virtual int dim () const override
 set the size of the particle in the statevector
 
virtual int type () const override
 type
 
virtual void addToConstraintList (constraintlist &alist, int depth) const override
 add to list
 
virtual std::string parname (int index) const override
 name
 
virtual int momIndex () const override
 get momentum index
 
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
 
virtual int nFinalChargedCandidates () const
 number of charged candidates
 

Static Public Member Functions

static bool useEnergy (const Belle2::Particle &cand)
 has energy in fit params?
 
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

const int m_dim
 dimension of residuals and 'width' of H
 
bool m_init
 was initialized*
 
bool m_useEnergy
 has energy ins statevector
 
Eigen::Matrix< double, 1, 4 > m_clusterPars
 constrains measured params (x_c, y_c, z_c, E_c)
 
Eigen::Matrix< double, 4, 4 > m_covariance
 covariance (x_c,y_c,z_c,E_c) of measured pars
 
int m_i1
 index with the highest momentum.
 
int m_i2
 random other index
 
int m_i3
 another random index
 
const float m_momentumScalingFactor
 scale the momentum / energy by this correction factor
 
int m_index
 index
 
std::string m_name
 name
 

Detailed Description

representation of the photon constraint

Definition at line 16 of file RecoPhoton.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

◆ RecoPhoton()

RecoPhoton ( Belle2::Particle * bc,
const ParticleBase * mother )

constructor

Definition at line 25 of file RecoPhoton.cc.

25 : RecoParticle(particle, mother),
26 m_dim(3),
27 m_init(false),
28 m_useEnergy(useEnergy(*particle)),
29 m_clusterPars(),
30 m_covariance(),
31 m_momentumScalingFactor(particle->getEffectiveMomentumScale())
32 {
33 initParams();
34 }

◆ ~RecoPhoton()

virtual ~RecoPhoton ( )
inlinevirtual

destructor

Definition at line 23 of file RecoPhoton.h.

23{};

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 to list

Implements ParticleBase.

Definition at line 53 of file RecoPhoton.h.

54 {
55 alist.push_back(Constraint(this, Constraint::photon, depth, dimM())) ;
56 }

◆ 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 }

◆ 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->getECLCluster()) { // external reconstructed photon
117 rc = new RecoPhoton(particle, mother);
118 } else if (particle->getKLMCluster()) { // external reconstructed klong
119 rc = new RecoKlong(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

set the size of the particle in the statevector

Reimplemented from RecoParticle.

Definition at line 47 of file RecoPhoton.h.

47{ return m_dim; }

◆ dimM()

virtual int dimM ( ) const
inlineoverridevirtual

sets the size of the corresponding residual projection

Implements RecoParticle.

Definition at line 41 of file RecoPhoton.h.

41{ return dim(); }

◆ eneIndex()

int eneIndex ( ) const
inlineinherited

get energy index

Definition at line 138 of file ParticleBase.h.

138{ 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 116 of file ParticleBase.h.

116{} ;

◆ hasEnergy()

virtual bool hasEnergy ( ) const
inlineoverridevirtual

how should the energy be calculated ?

from momentum or from E ?

Reimplemented from RecoParticle.

Definition at line 44 of file RecoPhoton.h.

44{ 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

Reimplemented from ParticleBase.

Definition at line 74 of file RecoPhoton.cc.

75 {
76 const int momindex = momIndex();
77 const int posindex = mother()->posIndex();
78
79 const double factorE = 1000 * m_covariance(3, 3);
80 const double factorX = 1000; // ~ 10cm error on initial vertex
81
82 fitparams.getCovariance().block<4, 4>(momindex, momindex) =
83 Eigen::Matrix<double, 4, 4>::Identity(4, 4) * factorE;
84
85 fitparams.getCovariance().block<3, 3>(posindex, posindex) =
86 Eigen::Matrix<double, 3, 3>::Identity(3, 3) * factorX;
87
88 return ErrCode(ErrCode::Status::success);
89 }

◆ initMotherlessParticle()

ErrCode initMotherlessParticle ( FitParams & fitparams)
overridevirtual

init particle without mother

Reimplemented from RecoParticle.

Definition at line 57 of file RecoPhoton.cc.

58 {
59 return ErrCode(ErrCode::Status::success);
60 }

◆ initParams()

ErrCode initParams ( )

update or init params

Definition at line 91 of file RecoPhoton.cc.

92 {
93 const Belle2::ECLCluster* cluster = particle()->getECLCluster();
94 const Belle2::ECLCluster::EHypothesisBit clusterhypo = particle()->getECLClusterEHypothesisBit();
95 const ROOT::Math::XYZVector centroid = cluster->getClusterPosition();
96 const double energy = cluster->getEnergy(clusterhypo);
97
98 m_init = true;
99 m_covariance = Eigen::Matrix<double, 4, 4>::Zero(4, 4);
100 Belle2::ClusterUtils C;
101
102 TMatrixDSym cov_EPhiTheta = C.GetCovarianceMatrix3x3FromCluster(cluster);
103
104 Eigen::Matrix<double, 2, 2> covPhiTheta = Eigen::Matrix<double, 2, 2>::Zero(2, 2);
105
106 for (int row = 0; row < 2; ++row) {
107 // we go through all elements here instead of selfadjoint view later
108 for (int col = 0; col < 2; ++col) {
109 covPhiTheta(row, col) = cov_EPhiTheta[row + 1][col + 1];
110 }
111 }
112
113 // the in going x-E correlations are 0 so we don't fill them
114 const double R = cluster->getR();
115 const double theta = cluster->getPhi();
116 const double phi = cluster->getTheta();
117
118 const double st = std::sin(theta);
119 const double ct = std::cos(theta);
120 const double sp = std::sin(phi);
121 const double cp = std::cos(phi);
122
123 Eigen::Matrix<double, 2, 3> polarToCartesian = Eigen::Matrix<double, 2, 3>::Zero(2, 3);
124
125 // polarToCartesian({phi,theta} -> {x,y,z} )
126 polarToCartesian(0, 0) = -1. * R * st * sp; // dx/dphi
127 polarToCartesian(0, 1) = R * st * cp; // dy/dphi
128 polarToCartesian(0, 2) = 0; // dz/dphi
129
130 polarToCartesian(1, 0) = R * ct * cp; // dx/dtheta
131 polarToCartesian(1, 1) = R * ct * sp; // dy/dtheta
132 polarToCartesian(1, 2) = -1. * R * st; // dz/dtheta
133
134 m_covariance.block<3, 3>(0, 0) = polarToCartesian.transpose() * covPhiTheta * polarToCartesian;
135
136 m_covariance(3, 3) = cov_EPhiTheta[0][0];
137 m_clusterPars(0) = centroid.X();
138 m_clusterPars(1) = centroid.Y();
139 m_clusterPars(2) = centroid.Z();
140 m_clusterPars(3) = energy;
141
142 auto p_vec = particle()->getMomentum();
143 // find highest momentum, eliminate dim with highest mom
144 if ((std::abs(p_vec.X()) >= std::abs(p_vec.Y())) && (std::abs(p_vec.X()) >= std::abs(p_vec.Z()))) {
145 m_i1 = 0;
146 m_i2 = 1;
147 m_i3 = 2;
148 } else if ((std::abs(p_vec.Y()) >= std::abs(p_vec.X())) && (std::abs(p_vec.Y()) >= std::abs(p_vec.Z()))) {
149 m_i1 = 1;
150 m_i2 = 0;
151 m_i3 = 2;
152 } else if ((std::abs(p_vec.Z()) >= std::abs(p_vec.Y())) && (std::abs(p_vec.Z()) >= std::abs(p_vec.X()))) {
153 m_i1 = 2;
154 m_i2 = 1;
155 m_i3 = 0;
156 } else {
157 B2ERROR("Could not estimate highest momentum for photon constraint. Aborting this fit.\n px: "
158 << p_vec.X() << " py: " << p_vec.Y() << " pz: " << p_vec.Z() << " calculated from Ec: " << m_clusterPars(3));
159 return ErrCode(ErrCode::Status::photondimerror);
160 }
161
162 return ErrCode(ErrCode::Status::success);
163 }
double R
typedef autogenerated by FFTW
const TMatrixDSym GetCovarianceMatrix3x3FromCluster(const ECLCluster *cluster, int particleHypo=Const::photon.getPDGCode())
Returns 3x3 covariance matrix (E, theta, phi)
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
Definition ECLCluster.h:31

◆ initParticleWithMother()

ErrCode initParticleWithMother ( FitParams & fitparams)
overridevirtual

init particle with mother

Implements ParticleBase.

Definition at line 36 of file RecoPhoton.cc.

37 {
38 const int posindexmother = mother()->posIndex();
39
40 Eigen::Matrix<double, 1, 3> vertexToCluster = Eigen::Matrix<double, 1, 3>::Zero(1, 3);
41 for (unsigned int i = 0; i < 3; ++i) {
42 vertexToCluster(i) = m_clusterPars(i) - fitparams.getStateVector()(posindexmother + i);
43 }
44
45 const double distanceToMother = vertexToCluster.norm();
46 const double energy = m_momentumScalingFactor * m_clusterPars(3); // apply scaling factor to correct energy bias
47 const int momindex = momIndex();
48
49 for (unsigned int i = 0; i < 3; ++i) {
50 // px = E dx/|dx|
51 fitparams.getStateVector()(momindex + i) = energy * vertexToCluster(i) / distanceToMother;
52 }
53
54 return ErrCode(ErrCode::Status::success);
55 }

◆ initTau()

ErrCode initTau ( FitParams & par) const
protectedinherited

initialises tau as a length

Definition at line 435 of file ParticleBase.cc.

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

get momentum index

Reimplemented from ParticleBase.

Definition at line 42 of file RecoParticle.h.

42{ return index(); }

◆ mother()

const ParticleBase * mother ( ) const
inlineinherited

getMother() / hasMother()

Definition at line 98 of file ParticleBase.h.

98{ 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
overridevirtualinherited

name

Reimplemented from ParticleBase.

Definition at line 23 of file RecoParticle.cc.

24 {
25 return ParticleBase::parname(index + 4);
26 }

◆ 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 }

◆ 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::Matrix < double, 1, -1, 1, 1, 3 > x_vec = fitparams.getStateVector().segment(posindex, dim);
259 Eigen::Matrix < double, 1, -1, 1, 1, 3 > x_m = fitparams.getStateVector().segment(posindexmother, dim);
260 Eigen::Matrix < double, 1, -1, 1, 1, 3 > 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 for (int row = 0; row < dim; ++row) {
295
296 double posxmother = x_m(row);
297 double posx = x_vec(row);
298 double momx = p_vec(row);
299
303 p.getResiduals()(row) = posxmother + tau * momx / mom - posx ;
304 p.getH()(row, posindexmother + row) = 1;
305 p.getH()(row, posindex + row) = -1;
306 p.getH()(row, tauindex) = momx / mom;
307 }
308
309 return ErrCode(ErrCode::Status::success);
310 }

◆ projectMassConstraint()

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

project mass constraint abstract

Definition at line 414 of file ParticleBase.cc.

416 {
417 assert(m_config);
418 if (m_config->m_massConstraintType == 0) {
419 return projectMassConstraintParticle(fitparams, p);
420 } else {
421 return projectMassConstraintDaughters(fitparams, p);
422 }
423 }

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

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

384 {
385 double mass = 0;
386 if (particle()->hasExtraInfo("treeFitterMassConstraintValue")) mass = particle()->getExtraInfo("treeFitterMassConstraintValue");
387 else mass = particle()->getPDGMass();
388 const double mass2 = mass * mass;
389 const int momindex = momIndex();
390 const double px = fitparams.getStateVector()(momindex);
391 const double py = fitparams.getStateVector()(momindex + 1);
392 const double pz = fitparams.getStateVector()(momindex + 2);
393 const double E = fitparams.getStateVector()(momindex + 3);
394
398 p.getResiduals()(0) = mass2 - E * E + px * px + py * py + pz * pz;
399
400 p.getH()(0, momindex) = 2.0 * px;
401 p.getH()(0, momindex + 1) = 2.0 * py;
402 p.getH()(0, momindex + 2) = 2.0 * pz;
403 p.getH()(0, momindex + 3) = -2.0 * E;
404
405 // TODO 0 in most cases -> needs special treatment if width=0 to not crash chi2 calculation
406 // const double width = TDatabasePDG::Instance()->GetParticle(particle()->getPDGCode())->Width();
407 // transport measurement uncertainty into residual system
408 // f' = sigma_x^2 * (df/dx)^2
409 // p.getV()(0) = width * width * 4 * mass2;
410
411 return ErrCode(ErrCode::Status::success);
412 }

◆ projectRecoConstraint()

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

project photon constraint

m : decay vertex mother p : momentum photon c : position cluster so: m + p = c thus (tau converts p to the correct units): 0 = c - m - tau * p we have 3 geometric equations and eliminate tau using the dimension with the highest momentum (because we have to divide by that momentum) only downside is we have to figure out which dimension this is the 4th equation is the energy which we keep as: 0 = E - |p| just to be sure, essentially this is always zero because p is build from E

Implements RecoParticle.

Definition at line 165 of file RecoPhoton.cc.

166 {
167 const int momindex = momIndex();
168 const int posindex = mother()->posIndex();
184
185 const Eigen::Matrix<double, 1, 3> x_vertex = fitparams.getStateVector().segment(posindex, 3);
186 const Eigen::Matrix<double, 1, 3> p_vec = fitparams.getStateVector().segment(momindex, 3);
187
188 if (0 == p_vec[m_i1]) {
189 return ErrCode(ErrCode::photondimerror);
190 }
191
192 // p_vec[m_i1] must not be 0
193 const double elim = (m_clusterPars[m_i1] - x_vertex[m_i1]) / p_vec[m_i1];
194 const double mom = p_vec.norm();
195
196 // r'
197 Eigen::Matrix<double, 3, 1> residual3 = Eigen::Matrix<double, 3, 1>::Zero(3, 1);
198 residual3(0) = m_clusterPars[m_i2] - x_vertex[m_i2] - p_vec[m_i2] * elim;
199 residual3(1) = m_clusterPars[m_i3] - x_vertex[m_i3] - p_vec[m_i3] * elim;
200 residual3(2) = m_momentumScalingFactor * m_clusterPars[3] - mom; // scale measured energy by scaling factor
201
202 // dr'/dm | m:={xc,yc,zc,Ec} the measured quantities
203 Eigen::Matrix<double, 3, 4> P = Eigen::Matrix<double, 3, 4>::Zero(3, 4);
204 // deriving by the cluster pars
205 P(0, m_i2) = 1;
206 P(0, m_i1) = -p_vec[m_i2] / p_vec[m_i1];
207
208 P(1, m_i3) = 1;
209 P(1, m_i1) = -p_vec[m_i3] / p_vec[m_i1];
210 P(2, 3) = 1; // dE/dEc
211
212 p.getResiduals().segment(0, 3) = residual3;
213
214 p.getV() = P * m_covariance.selfadjointView<Eigen::Lower>() * P.transpose();
215
216 // dr'/dm | m:={x,y,z,px,py,pz,E}
217 // x := x_vertex (decay vertex of mother)
218 p.getH()(0, posindex + m_i1) = p_vec[m_i2] / p_vec[m_i1];
219 p.getH()(0, posindex + m_i2) = -1.0;
220 p.getH()(0, posindex + m_i3) = 0;
221
222 p.getH()(1, posindex + m_i1) = p_vec[m_i3] / p_vec[m_i1];
223 p.getH()(1, posindex + m_i2) = 0;
224 p.getH()(1, posindex + m_i3) = -1.0;
225
226 // elim already divided by p_vec[m_i1]
227 p.getH()(0, momindex + m_i1) = p_vec[m_i2] * elim / p_vec[m_i1];
228 p.getH()(0, momindex + m_i2) = -1. * elim;
229 p.getH()(0, momindex + m_i3) = 0;
230
231 p.getH()(1, momindex + m_i1) = p_vec[m_i3] * elim / p_vec[m_i1];
232 p.getH()(1, momindex + m_i2) = 0;
233 p.getH()(1, momindex + m_i3) = -1. * elim;
234
235 p.getH()(2, momindex + m_i1) = -1. * p_vec[m_i1] / mom;
236 p.getH()(2, momindex + m_i2) = -1. * p_vec[m_i2] / mom;
237 p.getH()(2, momindex + m_i3) = -1. * p_vec[m_i3] / mom;
238 // the photon does not store an energy in the state vector
239 // so no p.getH()(2, momindex + 3) here
240
241 return ErrCode(ErrCode::Status::success);
242 }

◆ 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 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

Implements ParticleBase.

Definition at line 50 of file RecoPhoton.h.

50{ return kRecoPhoton ; }

◆ 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 }

◆ useEnergy()

bool useEnergy ( const Belle2::Particle & cand)
static

has energy in fit params?

Definition at line 62 of file RecoPhoton.cc.

63 {
64 bool rc = true;
65 const int pdg = particle.getPDGCode();
66 if (pdg &&
67 Belle2::Const::ParticleType(pdg) != Belle2::Const::photon &&
68 Belle2::Const::ParticleType(pdg) != Belle2::Const::pi0) {
69 rc = false;
70 }
71 return rc;
72 }
static const ParticleType pi0
neutral pion particle
Definition Const.h:674
static const ParticleType photon
photon particle
Definition Const.h:673

Member Data Documentation

◆ m_clusterPars

Eigen::Matrix<double, 1, 4> m_clusterPars
private

constrains measured params (x_c, y_c, z_c, E_c)

Definition at line 73 of file RecoPhoton.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, 4, 4> m_covariance
private

covariance (x_c,y_c,z_c,E_c) of measured pars

Definition at line 76 of file RecoPhoton.h.

◆ m_daughters

std::vector<ParticleBase*> m_daughters
protectedinherited

daughter container

Definition at line 198 of file ParticleBase.h.

◆ m_dim

const int m_dim
private

dimension of residuals and 'width' of H

Definition at line 64 of file RecoPhoton.h.

◆ m_i1

int m_i1
private

index with the highest momentum.

We have to make sure this does not change during the fit.

Definition at line 79 of file RecoPhoton.h.

◆ m_i2

int m_i2
private

random other index

Definition at line 81 of file RecoPhoton.h.

◆ m_i3

int m_i3
private

another random index

Definition at line 83 of file RecoPhoton.h.

◆ m_index

int m_index
privateinherited

index

Definition at line 208 of file ParticleBase.h.

◆ m_init

bool m_init
private

was initialized*

Definition at line 67 of file RecoPhoton.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 momentum / energy by this correction factor

Definition at line 86 of file RecoPhoton.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_particle

Belle2::Particle* m_particle
protectedinherited

pointer to framework type

Definition at line 192 of file ParticleBase.h.

◆ m_useEnergy

bool m_useEnergy
private

has energy ins statevector

Definition at line 70 of file RecoPhoton.h.


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