10 #include <Math/Vector4D.h> 
   12 #include <analysis/dataobjects/Particle.h> 
   14 #include <framework/logging/Logger.h> 
   15 #include <framework/gearbox/Const.h> 
   17 #include <analysis/VertexFitting/TreeFitter/FitManager.h> 
   18 #include <analysis/VertexFitting/TreeFitter/FitParams.h> 
   19 #include <analysis/VertexFitting/TreeFitter/DecayChain.h> 
   20 #include <analysis/VertexFitting/TreeFitter/ParticleBase.h> 
   22 namespace TreeFitter {
 
   28                          const bool useReferencing
 
   31     m_decaychain(nullptr),
 
   35     m_updateDaugthers(updateDaughters),
 
   38     m_useReferencing(useReferencing),
 
   53     const int nitermax = 100;
 
   54     const int maxndiverging = 3;
 
   55     const double dChisqConv = 
m_prec;
 
   59     if (
m_status == VertexStatus::UnFitted) {
 
   68       bool finished = 
false;
 
   70       for (niter = 0; niter < nitermax && !finished; ++niter) {
 
   87             if ((std::abs(deltachisq) / 
m_chiSquare < dChisqConv)) {
 
   92             } 
else if (deltachisq > 0 && ++ndiverging >= maxndiverging) {
 
   94               m_status = VertexStatus::NonConverged;
 
  105       if (niter == nitermax && 
m_status != VertexStatus::Success) {
 
  107         m_status = VertexStatus::NonConverged;
 
  115     if (
m_status == VertexStatus::Success) {
 
  124     return (
m_status == VertexStatus::Success);
 
  130     Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> cov = 
m_fitparams->
getCovariance().selfadjointView<Eigen::Lower>();
 
  133     if (posindex < 0 && pb->mother()) {
 
  140       for (
int row = 0; row < 4; ++row) {
 
  141         for (
int col = 0; col < 4; ++col) {
 
  142           returncov(row, col) = cov(momindex + row, momindex + col);
 
  146       for (
int row = 0; row < 3; ++row) {
 
  147         for (
int col = 0; col < 3; ++col) {
 
  148           returncov(row + 4, col + 4) = cov(posindex + row, posindex + col);
 
  153       Eigen::Matrix<double, 6, 6> cov6 =
 
  154         Eigen::Matrix<double, 6, 6>::Zero(6, 6);
 
  156       for (
int row = 0; row < 3; ++row) {
 
  157         for (
int col = 0; col < 3; ++col) {
 
  158           cov6(row, col) = cov(momindex + row, momindex + col);
 
  159           cov6(row + 3, col + 3) = cov(posindex + row, posindex + col);
 
  164       Eigen::Matrix<double, 3, 1> momVec =
 
  167       double energy2 = momVec.transpose() * momVec;
 
  168       energy2 += mass * mass;
 
  169       double energy = 
sqrt(energy2);
 
  171       Eigen::Matrix<double, 7, 6> jacobian =
 
  172         Eigen::Matrix<double, 7, 6>::Zero(7, 6);
 
  174       for (
int col = 0; col < 3; ++col) {
 
  175         jacobian(col, col) = 1;
 
  177         jacobian(col + 4, col + 3) = 1;
 
  180       Eigen::Matrix<double, 7, 7> cov7
 
  181         = jacobian * cov6.selfadjointView<Eigen::Lower>() * jacobian.transpose();
 
  183       for (
int row = 0; row < 7; ++row) {
 
  184         for (
int col = 0; col < 7; ++col) {
 
  185           returncov(row, col) = cov7(row, col);
 
  199     return pb != 
nullptr;
 
  206     if (posindex < 0 && pb.
mother()) {
 
  220           cand.
writeExtraInfo(
"modifiedPValue", TMath::Prob(fitparchi2, 3));
 
  225           if (motherPosIndex >= 0) {
 
  235       ROOT::Math::PxPyPzEVector p;
 
  247       TMatrixFSym cov7b2(7);
 
  254       const std::tuple<double, double>life = 
getLifeTime(cand);
 
  264     const bool updateableMother = 
updateCand(cand, isTreeHead);
 
  266     if (updateableMother and not cand.
hasExtraInfo(
"bremsCorrected") and
 
  267         not(cand.
hasExtraInfo(
"treeFitterTreatMeAsInvisible") and cand.
getExtraInfo(
"treeFitterTreatMeAsInvisible") == 1)) {
 
  269       for (
int i = 0; i < ndaughters; i++) {
 
  281       const int momindex = pb->
momIndex();
 
  282       const int tauIndex = pb->
tauIndex();
 
  286       Eigen::Matrix<double, 4, 4> comb_cov =  Eigen::Matrix<double, 4, 4>::Zero(4, 4);
 
  290       const double lenErr = std::get<1>(lenTuple);
 
  291       comb_cov(0, 0) = lenErr * lenErr;
 
  296       comb_cov.block<3, 3>(1, 1) = mom_cov;
 
  300       const double mom = mom_vec.norm();
 
  301       const double mom3 = mom * mom * mom;
 
  303       const double len = std::get<0>(lenTuple);
 
  304       const double t = len / mom * mBYc;
 
  306       Eigen::Matrix<double, 1, 4> jac = Eigen::Matrix<double, 1, 4>::Zero();
 
  307       jac(0) =  1. / mom * mBYc;
 
  308       jac(1) = -1. * len * mom_vec(0) / mom3 * mBYc;
 
  309       jac(2) = -1. * len * mom_vec(1) / mom3 * mBYc;
 
  310       jac(3) = -1. * len * mom_vec(2) / mom3 * mBYc;
 
  312       const double tErr2 = jac * comb_cov.selfadjointView<Eigen::Lower>() * jac.transpose();
 
  314       return std::make_tuple(t, 
std::sqrt(tErr2));
 
  316     return std::make_tuple(-999, -999);
 
  328       const int tauindex = pb->
tauIndex();
 
  330       const double lenErr2 = fitparams.
getCovariance()(tauindex, tauindex);
 
  331       return std::make_tuple(len, 
std::sqrt(lenErr2));
 
  333     return std::make_tuple(-999, -999);
 
  338     std::tuple<double, double> rc = std::make_tuple(-999, -999);
 
static const double speedOfLight
[cm/ns]
Class to store reconstructed particles.
std::string getName() const override
Return name of this particle.
void writeExtraInfo(const std::string &name, const double value)
Sets the user defined extraInfo.
void setVertex(const ROOT::Math::XYZVector &vertex)
Sets position (decay vertex)
void set4VectorDividingByMomentumScaling(const ROOT::Math::PxPyPzEVector &p4)
Sets Lorentz vector dividing by the momentum scaling factor.
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
unsigned getNDaughters(void) const
Returns number of daughter particles.
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
void setMomentumVertexErrorMatrix(const TMatrixFSym &errMatrix)
Sets 7x7 error matrix.
void setPValue(double pValue)
Sets chi^2 probability of fit.
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
double getExtraInfo(const std::string &name) const
Return given value if set.
const std::vector< int > m_massConstraintListPDG
list of pdg codes to mass constrain
this class does a lot of stuff: Build decaytree structure allowing to index particles and handle the ...
ErrCode initialize(FitParams &par)
initialize the chain
const ParticleBase * cand()
get candidate
const ParticleBase * locate(Belle2::Particle *bc) const
convert Belle2::particle into particle base(fitter base particle)
int dim() const
get dimension
ErrCode filterWithReference(FitParams &par, const FitParams &ref)
filter with respect to a previous iteration for better stability
ErrCode filter(FitParams &par)
filter down the chain
abstract errorocode be aware that the default is success
bool failure() const
returns true if errorcode is error
void reset()
reset the errorcode to default (success!)
const ConstraintConfiguration m_config
config container
DecayChain * m_decaychain
the decay tree
ErrCode m_errCode
errorcode
void getCovFromPB(const ParticleBase *pb, TMatrixFSym &returncov) const
extract cov from particle base
void updateTree(Belle2::Particle &particle, const bool isTreeHead) const
update the Belle2::Particles with the fit results
FitParams * m_fitparams
parameters to be fitted
const bool m_updateDaugthers
if this is set all daughters will be updated otherwise only the head of the tree
~FitManager()
destructor does stuff
int m_ndf
number of degrees of freedom for this topology
double m_prec
precision that is needed for status:converged (delta chi2)
bool fit()
main fit function that uses the kalman filter
int m_status
status of the current iteration
Belle2::Particle * particle()
getter for the head of the tree
Belle2::Particle * m_particle
head of the tree
double m_chiSquare
chi2 of the current iteration
bool m_useReferencing
use referencing
VertexStatus
status flag of the fit-itereation (the step in the newton method)
std::tuple< double, double > getLifeTime(Belle2::Particle &cand) const
get lifetime
std::tuple< double, double > getDecayLength(const ParticleBase *pb) const
get decay length
bool updateCand(Belle2::Particle &particle, const bool isTreeHead) const
update particles parameters with the fit results
Class to store and manage fitparams (statevector)
Eigen::Matrix< double, -1, -1, 0, MAX_MATRIX_SIZE, MAX_MATRIX_SIZE > & getCovariance()
getter for the states covariance
bool testCovariance() const
test if the covariance makes sense
Eigen::Matrix< double, -1, 1, 0, MAX_MATRIX_SIZE, 1 > & getStateVector()
getter for the fit parameters/statevector
double chiSquare() const
get chi2 of statevector
int nDof() const
get numer of degrees of freedom
base class for all particles
Belle2::Particle * particle() const
get basf2 particle
virtual int dim() const =0
get dimension of constraint
virtual void forceP4Sum(FitParams &) const
force p4 sum conservation all along the tree
virtual int posIndex() const
get vertex index (in statevector!)
virtual int momIndex() const
get momentum index
virtual bool hasEnergy() const
get momentum dimension
virtual int tauIndex() const
get tau index
const ParticleBase * mother() const
getMother() / hasMother()
double sqrt(double a)
sqrt for double