Belle II Software  release-08-01-10
RecoPhoton Class Reference

representation of the photon constraint More...

#include <RecoPhoton.h>

Inheritance diagram for RecoPhoton:
Collaboration diagram for RecoPhoton:

Public Types

enum  TFParticleType {
  kInteractionPoint ,
  kOrigin ,
  kComposite ,
  kRecoResonance ,
  kInternalParticle ,
  kRecoTrack ,
  kResonance ,
  kRecoPhoton ,
  kRecoKlong ,
  kMissingParticle
}
 particle types

 
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 More...
 
virtual int dimM () const override
 sets the size of the corresponding residual projection
 
virtual bool hasEnergy () const override
 how should the energy be calculated ? from momentum or from E ?

 
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. More...
 
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 More...
 
virtual ErrCode projectMassConstraintParticle (const FitParams &, Projection &) const
 project mass constraint using the particles parameters More...
 
virtual ErrCode projectMassConstraintDaughters (const FitParams &, Projection &) const
 project mass constraint using the parameters of the daughters More...
 
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. More...
 

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. More...
 
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 Function Documentation

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

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

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

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

◆ 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 162 of file RecoPhoton.cc.

163  {
164  const int momindex = momIndex() ;
165  const int posindex = mother()->posIndex();
182  const Eigen::Matrix<double, 1, 3> x_vertex = fitparams.getStateVector().segment(posindex, 3);
183  const Eigen::Matrix<double, 1, 3> p_vec = fitparams.getStateVector().segment(momindex, 3);
184 
185  if (0 == p_vec[m_i1]) { return ErrCode(ErrCode::photondimerror); }
186 
187  // p_vec[m_i1] must not be 0
188  const double elim = (m_clusterPars[m_i1] - x_vertex[m_i1]) / p_vec[m_i1];
189  const double mom = p_vec.norm();
190 
191  // r'
192  Eigen::Matrix<double, 3, 1> residual3 = Eigen::Matrix<double, 3, 1>::Zero(3, 1);
193  residual3(0) = m_clusterPars[m_i2] - x_vertex[m_i2] - p_vec[m_i2] * elim;
194  residual3(1) = m_clusterPars[m_i3] - x_vertex[m_i3] - p_vec[m_i3] * elim;
195  residual3(2) = m_momentumScalingFactor * m_clusterPars[3] - mom; // scale measured energy by scaling factor
196 
197  // dr'/dm | m:={xc,yc,zc,Ec} the measured quantities
198  Eigen::Matrix<double, 3, 4> P = Eigen::Matrix<double, 3, 4>::Zero(3, 4);
199  // deriving by the cluster pars
200  P(0, m_i2) = 1;
201  P(0, m_i1) = - p_vec[m_i2] / p_vec[m_i1];
202 
203  P(1, m_i3) = 1;
204  P(1, m_i1) = - p_vec[m_i3] / p_vec[m_i1];
205  P(2, 3) = 1; // dE/dEc
206 
207  p.getResiduals().segment(0, 3) = residual3;
208 
209  p.getV() = P * m_covariance.selfadjointView<Eigen::Lower>() * P.transpose();
210 
211  // dr'/dm | m:={x,y,z,px,py,pz,E}
212  // x := x_vertex (decay vertex of mother)
213  p.getH()(0, posindex + m_i1) = p_vec[m_i2] / p_vec[m_i1];
214  p.getH()(0, posindex + m_i2) = -1.0;
215  p.getH()(0, posindex + m_i3) = 0;
216 
217  p.getH()(1, posindex + m_i1) = p_vec[m_i3] / p_vec[m_i1];
218  p.getH()(1, posindex + m_i2) = 0;
219  p.getH()(1, posindex + m_i3) = -1.0;
220 
221  // elim already divided by p_vec[m_i1]
222  p.getH()(0, momindex + m_i1) = p_vec[m_i2] * elim / p_vec[m_i1];
223  p.getH()(0, momindex + m_i2) = -1. * elim;
224  p.getH()(0, momindex + m_i3) = 0;
225 
226  p.getH()(1, momindex + m_i1) = p_vec[m_i3] * elim / p_vec[m_i1];
227  p.getH()(1, momindex + m_i2) = 0;
228  p.getH()(1, momindex + m_i3) = -1. * elim;
229 
230  p.getH()(2, momindex + m_i1) = -1. * p_vec[m_i1] / mom;
231  p.getH()(2, momindex + m_i2) = -1. * p_vec[m_i2] / mom;
232  p.getH()(2, momindex + m_i3) = -1. * p_vec[m_i3] / mom;
233  // the photon does not store an energy in the state vector
234  // so no p.getH()(2, momindex + 3) here
235 
236  return ErrCode(ErrCode::Status::success);
237  }
virtual int posIndex() const
get vertex index (in statevector!)
Definition: ParticleBase.h:122
const ParticleBase * mother() const
getMother() / hasMother()
Definition: ParticleBase.h:98
virtual int momIndex() const override
get momentum index
Definition: RecoParticle.h:42
Eigen::Matrix< double, 4, 4 > m_covariance
covariance (x_c,y_c,z_c,E_c) of measured pars
Definition: RecoPhoton.h:76
const float m_momentumScalingFactor
scale the momentum / energy by this correction factor
Definition: RecoPhoton.h:86
int m_i3
another random index
Definition: RecoPhoton.h:83
int m_i1
index with the highest momentum.
Definition: RecoPhoton.h:79
int m_i2
random other index
Definition: RecoPhoton.h:81
Eigen::Matrix< double, 1, 4 > m_clusterPars
constrains measured params (x_c, y_c, z_c, E_c)
Definition: RecoPhoton.h:73

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

Member Data Documentation

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


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