10#include <analysis/modules/ParticleVertexFitter/ParticleVertexFitterModule.h>
13#include <framework/gearbox/Unit.h>
14#include <framework/gearbox/Const.h>
15#include <framework/geometry/B2Vector3.h>
16#include <framework/logging/Logger.h>
17#include <framework/particledb/EvtGenDatabasePDG.h>
20#include <analysis/dataobjects/Particle.h>
21#include <analysis/dataobjects/Btube.h>
22#include <mdst/dataobjects/V0.h>
25#include <analysis/utility/CLHEPToROOT.h>
26#include <analysis/utility/PCmsLabTransform.h>
27#include <analysis/utility/ParticleCopy.h>
28#include <analysis/utility/ROOTToCLHEP.h>
31#include <framework/geometry/BFieldManager.h>
34#include <analysis/VertexFitting/KFit/KFitConst.h>
64 "Confidence level to accept the fit. Particle candidates with "
65 "p-value less than confidenceLevel are removed from the particle "
66 "list. If set to -1, all candidates are kept; if set to 0, "
67 "the candidates failing the fit are removed.",
70 addParam(
"fitType",
m_fitType,
"type of the kinematic fit (vertex, massvertex, mass)",
string(
"vertex"));
72 "additional constraint on vertex: ipprofile, iptube, mother, iptubecut, pointing, btube",
74 addParam(
"decayString",
m_decayString,
"specifies which daughter particles are included in the kinematic fit",
string(
""));
79 "Type::[int]. List of daughter particles to mass constrain with int = pdg code. (only for MassFourCKFit)", {});
81 "Type::[string]. List of daughter particles to mass constrain with string = particle name. (only for MassFourCKFit)", {});
96 B2DEBUG(1,
"ParticleVertexFitterModule : magnetic field = " <<
m_Bfield);
104 TParticlePDG* particletemp = TDatabasePDG::Instance()->GetParticle((containedParticle).c_str());
111 B2INFO(
"ParticleVertexFitter: Using specified decay string: " <<
m_decayString);
113 B2INFO(
"ParticleVertexFitter: Additional " <<
m_withConstraint <<
" will be applied");
131 TMatrixDSym beamSpotCov(3);
148 std::vector<unsigned int> toRemove;
149 unsigned int nParticles =
m_plist->getListSize();
151 for (
unsigned iPart = 0; iPart < nParticles; iPart++) {
157 else B2ERROR(
"Daughters update works only when all daughters are selected. Daughters will not be updated");
165 TMatrixFSym mother_errMatrix(7);
166 mother_errMatrix = particle->getMomentumVertexErrorMatrix();
167 for (
int k = 0; k < 7; k++) {
168 for (
int j = 0; j < 7; j++) {
169 if (mother_errMatrix[k][j] > 0) {
180 toRemove.push_back(particle->getArrayIndex());
191 particle->setPValue(-1);
193 toRemove.push_back(particle->getArrayIndex());
196 m_plist->removeParticles(toRemove);
208 B2FATAL(
"ParticleVertexFitter: No magnetic field");
218 B2FATAL(
"ParticleVertexFitter: " <<
m_withConstraint <<
" ***invalid Constraint ");
225 B2FATAL(
"ParticleVertexFitter: KFit does not support yet selection of daughters via decay string except for vertex fit!");
241 B2FATAL(
"ParticleVertexFitter: Invalid options - mass-constrained fit using KFit does not work with iptube or ipprofile constraint.");
252 B2FATAL(
"ParticleVertexFitter: Invalid options - mass fit using KFit does not work with iptube or ipprofile constraint.");
261 B2FATAL(
"ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
270 B2FATAL(
"ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
279 B2FATAL(
"ParticleVertexFitter: Invalid options - recoil mass fit using KFit does not work with iptube or ipprofile constraint.");
292 B2FATAL(
"ParticleVertexFitter: " <<
m_fitType <<
" ***invalid fit type for the vertex fitter ");
299 }
catch (
const rave::CheckedFloatException&) {
300 B2ERROR(
"Invalid inputs (nan/inf)?");
307 B2FATAL(
"ParticleVertexFitter: " <<
m_vertexFitter <<
" ***invalid vertex fitter ");
309 if (!ok)
return false;
319 std::vector<const Particle*>& twoPhotonChildren)
323 for (
unsigned ichild = 0; ichild < mother->
getNDaughters(); ichild++) {
330 fitChildren.push_back(child);
336 auto itr = fitChildren.begin();
337 while (itr != fitChildren.end()) {
340 if (child->getPValue() < 0) {
341 B2WARNING(
"Daughter with PDG code " << child->getPDGCode() <<
" does not have a valid error matrix.");
344 bool isTwoPhotonParticle =
false;
347 if (child->getNDaughters() == 2) {
350 isTwoPhotonParticle =
true;
355 if (isTwoPhotonParticle) {
357 twoPhotonChildren.push_back(child);
358 itr = fitChildren.erase(itr);
368 const std::vector<const Particle*>& fitChildren)
370 if (fitChildren.empty())
371 B2WARNING(
"[ParticleVertexFitterModule::fillNotFitParticles] fitChildren is empty! Please call fillFitParticles firstly");
372 if (!notFitChildren.empty())
373 B2WARNING(
"[ParticleVertexFitterModule::fillNotFitParticles] notFitChildren is NOT empty!"
374 <<
" The function should be called only once");
380 std::function<bool(
const Particle*)> funcCheckInFit =
381 [&funcCheckInFit, ¬FitChildren, fitChildren](
const Particle * part) {
385 if (std::find(fitChildren.begin(), fitChildren.end(), part) != fitChildren.end())
389 if (part->getNDaughters() == 0)
394 bool isAnyChildrenInFit =
false;
395 vector<const Particle*> notFitChildren_tmp;
396 for (
unsigned ichild = 0; ichild < part->getNDaughters(); ichild++) {
398 const Particle* child = part->getDaughter(ichild);
399 bool isChildrenInFit = funcCheckInFit(child);
400 isAnyChildrenInFit = isChildrenInFit or isAnyChildrenInFit;
403 if (!isChildrenInFit)
404 notFitChildren_tmp.push_back(child);
408 if (isAnyChildrenInFit)
409 notFitChildren.insert(notFitChildren.end(), notFitChildren_tmp.begin(), notFitChildren_tmp.end());
413 return isAnyChildrenInFit;
418 for (
unsigned ichild = 0; ichild < mother->
getNDaughters(); ichild++) {
420 bool isGivenParticleOrAnyChildrenInFit = funcCheckInFit(child);
421 if (!isGivenParticleOrAnyChildrenInFit)
422 notFitChildren.push_back(child);
445 TMatrixFSym errMatrix(3);
446 for (
int i = 0; i < 3; i++)
447 for (
int j = 0; j < 3; j++)
448 errMatrix(i, j) = posErrorMatrix[i][j];
450 g1ErrMatrix.SetSub(4, errMatrix);
451 g2ErrMatrix.SetSub(4, errMatrix);
467 int err = km.
doFit();
485 std::vector<const Particle*> fitChildren;
486 std::vector<const Particle*> twoPhotonChildren;
487 bool validChildren =
fillFitParticles(mother, fitChildren, twoPhotonChildren);
492 std::vector<const Particle*> notFitChildren;
496 if (twoPhotonChildren.size() > 1) {
497 B2FATAL(
"[ParticleVertexFitterModule::doKVertexFit] Vertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
500 if ((fitChildren.size() < 2 && !ipTubeConstraint) || fitChildren.size() < 1) {
501 B2WARNING(
"[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is too low.");
509 if (mother->
getV0()) {
516 for (
auto& child : fitChildren)
519 if (ipProfileConstraint)
522 if (ipTubeConstraint)
526 int err = kv.
doFit();
536 if (twoPhotonChildren.size() == 0)
539 else if (twoPhotonChildren.size() == 1) {
545 const Particle* twoPhotonDaughter = twoPhotonChildren[0];
555 for (
auto& child : fitChildren)
560 if (ipProfileConstraint)
572 ROOT::Math::PxPyPzEVector total4Vector(mother->
get4Vector());
573 for (
auto& child : notFitChildren)
574 total4Vector += child->get4Vector();
584 std::vector<const Particle*> fitChildren;
585 std::vector<const Particle*> twoPhotonChildren;
586 bool validChildren =
fillFitParticles(mother, fitChildren, twoPhotonChildren);
591 if (twoPhotonChildren.size() > 1) {
592 B2FATAL(
"[ParticleVertexFitterModule::doKVertexFit] MassVertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
595 if (fitChildren.size() < 2) {
596 B2WARNING(
"[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
601 if (twoPhotonChildren.size() == 0) {
606 if (mother->
getV0()) {
613 for (
auto child : fitChildren)
617 int err = kmv.
doFit();
623 }
else if (twoPhotonChildren.size() == 1) {
632 for (
auto child : fitChildren)
636 int err = kv.
doFit();
640 const Particle* twoPhotonDaughter = twoPhotonChildren[0];
650 for (
auto child : fitChildren)
675 std::vector<const Particle*> fitChildren;
676 std::vector<const Particle*> twoPhotonChildren;
677 bool validChildren =
fillFitParticles(mother, fitChildren, twoPhotonChildren);
682 if (twoPhotonChildren.size() > 0) {
683 B2FATAL(
"[ParticleVertexFitterModule::doKMassPointingVertexFit] MassPointingVertex fit using KFit does not support fit with two-photon daughters (yet).");
686 if (fitChildren.size() < 2) {
687 B2WARNING(
"[ParticleVertexFitterModule::doKMassPointingVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
696 for (
auto child : fitChildren)
704 int err = kmpv.
doFit();
705 if (err != 0)
return false;
719 for (
unsigned ichild = 0; ichild < mother->
getNDaughters(); ichild++) {
722 if (child->getPValue() < 0)
return false;
730 int err = km.
doFit();
732 if (err != 0)
return false;
746 for (
unsigned ichild = 0; ichild < mother->
getNDaughters(); ichild++) {
749 if (child->getNDaughters() > 0) {
751 if (!err)
return false;
753 if (child->getPValue() < 0)
return false;
763 int err = kf.
doFit();
765 if (err != 0)
return false;
779 for (
unsigned ichild = 0; ichild < mother->
getNDaughters(); ichild++) {
782 if (child->getNDaughters() > 0) {
785 std::vector<unsigned> childId;
788 if (!err)
return false;
790 if (child->getPValue() < 0)
return false;
799 int err = kf.
doFit();
801 if (err != 0)
return false;
813 for (
unsigned ichild = 0; ichild < mother->
getNDaughters(); ichild++) {
816 if (child->getPValue() < 0)
return false;
828 int err = kf.
doFit();
830 if (err != 0)
return false;
848 std::vector<Particle*> daughters = mother->
getDaughters();
851 if (daughters.size() != track_count)
854 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
859 ROOT::Math::PxPyPzEVector i4Vector(kv.
getTrackMomentum(iChild).x() - a * dy,
863 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
865 daughters[iChild]->setVertex(
867 daughters[iChild]->setMomentumVertexErrorMatrix(
876 if (fitChildren.size() != track_count)
879 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
880 auto daughter =
const_cast<Particle*
>(fitChildren[iChild]);
886 ROOT::Math::PxPyPzEVector i4Vector(kv.
getTrackMomentum(iChild).x() - a * dy,
890 daughter->set4VectorDividingByMomentumScaling(i4Vector);
892 daughter->setVertex(CLHEPToROOT::getXYZVector(kv.
getTrackPosition(iChild)));
893 daughter->setMomentumVertexErrorMatrix(CLHEPToROOT::getTMatrixFSym(kv.
getTrackError(iChild)));
897 std::function<bool(
Particle*)> funcUpdateMomentum =
898 [&funcUpdateMomentum, fitChildren](
Particle * part) {
900 if (part->getNDaughters() == 0) {
902 if (std::find(fitChildren.begin(), fitChildren.end(), part) != fitChildren.end())
908 bool includeFitChildren =
false;
911 for (
auto daughter : part->getDaughters())
912 includeFitChildren = funcUpdateMomentum(daughter) || includeFitChildren;
914 if (includeFitChildren) {
916 ROOT::Math::PxPyPzEVector sum4Vector;
917 for (
auto daughter : part->getDaughters())
918 sum4Vector += daughter->get4Vector();
920 part->set4VectorDividingByMomentumScaling(sum4Vector);
923 return includeFitChildren;
928 funcUpdateMomentum(daughter);
947 std::vector<Particle*> daughters = mother->
getDaughters();
950 if (daughters.size() != track_count)
953 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
958 ROOT::Math::PxPyPzEVector i4Vector(kmv.
getTrackMomentum(iChild).x() - a * dy,
962 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
964 daughters[iChild]->setVertex(
966 daughters[iChild]->setMomentumVertexErrorMatrix(
987 std::vector<Particle*> daughters = mother->
getDaughters();
990 if (daughters.size() != track_count)
993 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
998 ROOT::Math::PxPyPzEVector i4Vector(kmpv.
getTrackMomentum(iChild).x() - a * dy,
1002 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1004 daughters[iChild]->setVertex(
1006 daughters[iChild]->setMomentumVertexErrorMatrix(
1026 std::vector<Particle*> daughters = mother->
getDaughters();
1029 if (daughters.size() != track_count)
1032 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
1037 ROOT::Math::PxPyPzEVector i4Vector(km.
getTrackMomentum(iChild).x() - a * dy,
1041 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1043 daughters[iChild]->setVertex(
1045 daughters[iChild]->setMomentumVertexErrorMatrix(
1067 std::vector<Particle*> daughters = mother->
getDaughters();
1069 const unsigned nd = daughters.size();
1071 std::vector<std::vector<unsigned>> pars;
1072 std::vector<Particle*> allparticles;
1073 for (
unsigned ichild = 0; ichild < nd; ichild++) {
1075 std::vector<unsigned> pard;
1076 if (daughter->getNDaughters() > 0) {
1078 pars.push_back(pard);
1079 allparticles.push_back(daughters[ichild]);
1082 pars.push_back(pard);
1083 allparticles.push_back(daughters[ichild]);
1089 if (l != track_count)
1092 for (
unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
1093 ROOT::Math::PxPyPzEVector childMoms;
1094 ROOT::Math::XYZVector childPoss;
1095 TMatrixFSym childErrMatrixs(7);
1096 for (
unsigned int iChild : pars[iDaug]) {
1097 childMoms = childMoms +
1098 CLHEPToROOT::getLorentzVector(
1100 childPoss = childPoss +
1101 CLHEPToROOT::getXYZVector(
1103 TMatrixFSym childErrMatrix =
1105 childErrMatrixs = childErrMatrixs + childErrMatrix;
1107 allparticles[iDaug]->set4Vector(childMoms);
1108 allparticles[iDaug]->setVertex(childPoss);
1109 allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
1129 std::vector<Particle*> daughters = mother->
getDaughters();
1131 const unsigned nd = daughters.size();
1133 std::vector<std::vector<unsigned>> pars;
1134 std::vector<Particle*> allparticles;
1135 for (
unsigned ichild = 0; ichild < nd; ichild++) {
1137 std::vector<unsigned> pard;
1138 if (daughter->getNDaughters() > 0) {
1140 pars.push_back(pard);
1141 allparticles.push_back(daughters[ichild]);
1144 pars.push_back(pard);
1145 allparticles.push_back(daughters[ichild]);
1151 if (l != track_count)
1154 for (
unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
1155 ROOT::Math::PxPyPzEVector childMoms;
1156 ROOT::Math::XYZVector childPoss;
1157 TMatrixFSym childErrMatrixs(7);
1158 for (
unsigned int iChild : pars[iDaug]) {
1159 childMoms = childMoms +
1160 CLHEPToROOT::getLorentzVector(
1162 childPoss = childPoss +
1163 CLHEPToROOT::getXYZVector(
1165 TMatrixFSym childErrMatrix =
1167 childErrMatrixs = childErrMatrixs + childErrMatrix;
1169 allparticles[iDaug]->set4Vector(childMoms);
1170 allparticles[iDaug]->setVertex(childPoss);
1171 allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
1191 std::vector<Particle*> daughters = mother->
getDaughters();
1194 if (daughters.size() != track_count)
1197 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
1202 ROOT::Math::PxPyPzEVector i4Vector(kf.
getTrackMomentum(iChild).x() - a * dy,
1206 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1208 daughters[iChild]->setVertex(
1210 daughters[iChild]->setMomentumVertexErrorMatrix(
1220 std::vector<unsigned>& parm, std::vector<Particle*>& allparticles,
const Particle* daughter)
1222 std::vector <Belle2::Particle*> daughters = daughter->getDaughters();
1223 for (
unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
1224 const Particle* child = daughter->getDaughter(ichild);
1225 std::vector<unsigned> pard;
1226 if (child->getNDaughters() > 0) {
1228 parm.insert(parm.end(), pard.begin(), pard.end());
1229 pars.push_back(pard);
1230 allparticles.push_back(daughters[ichild]);
1234 pars.push_back(pard);
1235 allparticles.push_back(daughters[ichild]);
1261 for (
auto& itrack : tracksVertex) {
1262 if (itrack != mother) rf.
addTrack(itrack);
1268 bool mothSel =
false;
1270 for (
unsigned itrack = 0; itrack < tracksVertex.size(); itrack++) {
1271 if (tracksVertex[itrack] != mother) {
1272 rsf.
addTrack(tracksVertex[itrack]);
1273 B2DEBUG(1,
"ParticleVertexFitterModule: Adding particle " << tracksName[itrack] <<
" to vertex fit ");
1276 if (tracksVertex[itrack] == mother) mothSel =
true;
1281 bool mothIPfit =
false;
1282 if (tracksVertex.size() == 1 && mothSel ==
true &&
m_withConstraint !=
"" && nTrk == 0) {
1284 if (tracksVertex[0] != mother)
1285 B2FATAL(
"ParticleVertexFitterModule: FATAL Error in IP constrained mother fit");
1291 ROOT::Math::XYZVector pos;
1292 TMatrixDSym RerrMatrix(3);
1298 for (
auto& itrack : tracksVertex) {
1300 nvert = rsg.
fit(
"kalman");
1303 RerrMatrix = rsg.
getCov(0);
1305 ROOT::Math::PxPyPzEVector mom(mother->
get4Vector());
1306 TMatrixDSym errMatrix(7);
1307 for (
int i = 0; i < 7; i++) {
1308 for (
int j = 0; j < 7; j++) {
1309 if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1310 else {errMatrix[i][j] = 0;}
1330 }
else {
return false;}
1340 ROOT::Math::PxPyPzEVector mom(mother->
get4Vector());
1341 TMatrixDSym errMatrix(7);
1342 for (
int i = 0; i < 7; i++) {
1343 for (
int j = 0; j < 7; j++) {
1344 if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1345 else {errMatrix[i][j] = 0;}
1349 }
else {
return false;}
1352 if (mothSel && nTrk > 1) {
1355 int nKfit = rf.
fit();
1359 if (nKfit > 0) {
return true;}
1368 int nVert = rf.
fit();
1371 if (nVert != 1)
return false;
1378 int nVert = rf.
fit();
1380 if (nVert != 1)
return false;
1385 int nVert = rf.
fit();
1388 if (nVert != 1)
return false;
1391 B2FATAL(
"fitType : " <<
m_fitType <<
" ***invalid fit type ");
1398 const std::vector<const Particle*>& tracksVertex)
1407 for (
auto& vi : tracksVertex) {
1409 nNotIncluded = nNotIncluded - 1;
1417 if (nNotIncluded == 0) isAll =
true;
1423 for (
unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1424 const Particle* child = particle->getDaughter(ichild);
1427 if (child->getPValue() < 0)
return false;
1436 std::vector<unsigned>& particleId)
1438 for (
unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1439 const Particle* child = particle->getDaughter(ichild);
1440 if (child->getNDaughters() > 0) {
1443 std::vector<unsigned> childId;
1446 particleId.insert(particleId.end(), childId.begin(), childId.end());
1448 if (child->getPValue() < 0)
return false;
1459 CLHEP::HepSymMatrix covMatrix(3, 0);
1461 for (
int i = 0; i < 3; i++) {
1463 for (
int j = 0; j < 3; j++) {
1473 CLHEP::HepSymMatrix err(7, 0);
1475 for (
int i = 0; i < 3; i++) {
1476 for (
int j = 0; j < 3; j++) {
1485 ROOTToCLHEP::getHepLorentzVector(iptube_mom),
1498 TMatrixDSym beamSpotCov =
m_beamSpotDB->getCovVertex();
1499 beamSpotCov(2, 2) = cut * cut;
1500 double thetab = boostDir.
Theta();
1501 double phib = boostDir.
Phi();
1503 double stb = TMath::Sin(thetab);
1504 double ctb = TMath::Cos(thetab);
1505 double spb = TMath::Sin(phib);
1506 double cpb = TMath::Cos(phib);
1509 TMatrix rz(3, 3); rz(2, 2) = 1;
1510 rz(0, 0) = cpb; rz(0, 1) = spb;
1511 rz(1, 0) = -1 * spb; rz(1, 1) = cpb;
1513 TMatrix ry(3, 3); ry(1, 1) = 1;
1514 ry(0, 0) = ctb; ry(0, 2) = -1 * stb;
1515 ry(2, 0) = stb; ry(2, 2) = ctb;
1517 TMatrix r(3, 3); r.Mult(rz, ry);
1518 TMatrix rt(3, 3); rt.Transpose(r);
1520 TMatrix TubePart(3, 3); TubePart.Mult(rt, beamSpotCov);
1521 TMatrix Tube(3, 3); Tube.Mult(TubePart, r);
1530 TMatrixDSym beamSpotCov =
m_beamSpotDB->getCovVertex();
1531 for (
int i = 0; i < 3; i++)
1532 beamSpotCov(i, i) += width * width;
1539 double chi2TrackL = 0;
1541 for (
int iTrack = 0; iTrack < kv.
getTrackCount(); iTrack++) {
1553 TVectorD boostD(0, 6, 0., 0., 0., 0., boost3.
X(), boost3.
Y(), boost3.
Z(),
"END");
1555 double dLBoost = dPos.
Dot(boost3);
1557 chi2TrackL += TMath::Power(dLBoost, 2) / err.Similarity(boostD);
DataType Phi() const
The azimuth angle.
DataType Z() const
access variable Z (= .at(2) without boundary check)
DataType Theta() const
The polar angle.
DataType X() const
access variable X (= .at(0) without boundary check)
DataType Y() const
access variable Y (= .at(1) without boundary check)
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
void SetXYZ(DataType x, DataType y, DataType z)
set all coordinates using data type
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
For each MCParticle with hits in the CDC, this class stores some summarising information on those hit...
TMatrixFSym getTubeMatrix() const
Returns Btube matrix.
Eigen::Matrix< double, 3, 1 > getTubeCenter() const
Returns Btube center.
int getPDGCode() const
PDG code.
static const ParticleType pi0
neutral pion particle
static const double speedOfLight
[cm/ns]
static const ParticleType photon
photon particle
bool init(const std::string &str)
Initialise the DecayDescriptor from given string.
std::vector< std::string > getSelectionNames()
Return list of human readable names of selected particles.
std::vector< const Particle * > getSelectionParticles(const Particle *particle)
Get a vector of pointers with selected daughters in the decay tree.
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
TMatrixDSym m_beamSpotCov
Beam spot covariance matrix.
bool makeKMassMother(analysis::MassFitKFit &kv, Particle *p)
Update mother particle after mass fit using KFit.
std::string m_withConstraint
additional constraint on vertex
bool doKMassPointingVertexFit(Particle *p)
Mass-constrained vertex fit with additional pointing constraint using KFit.
void smearBeamSpot(double width)
smear beam spot covariance
bool makeMassKFourCMother(analysis::MassFourCFitKFit &kv, Particle *p)
Update mother particle after MassFourC fit using KFit.
void updateMapOfTrackAndDaughter(unsigned &l, std::vector< std::vector< unsigned > > &pars, std::vector< unsigned > &pard, std::vector< Particle * > &allparticles, const Particle *daughter)
update the map of daughter and tracks, find out which tracks belong to each daughter.
bool doKVertexFit(Particle *p, bool ipProfileConstraint, bool ipTubeConstraint)
Unconstrained vertex fit using KFit.
bool doKMassVertexFit(Particle *p)
Mass-constrained vertex fit using KFit.
virtual void initialize() override
Initialize the Module.
void addIPTubeToKFit(analysis::VertexFitKFit &kv)
Adds IPTube constraint to the vertex fit using KFit.
bool addChildofParticletoMassKFit(analysis::MassFourCFitKFit &kf, const Particle *particle, std::vector< unsigned > &particleId)
Adds given particle's child to the MassFourCFitKFit.
virtual void event() override
Event processor.
bool makeKMassPointingVertexMother(analysis::MassPointingVertexFitKFit &kv, Particle *p)
Update mother particle after mass-constrained vertex fit with additional pointing constraint using KF...
bool m_updateDaughters
flag for daughters update
ParticleVertexFitterModule()
Constructor.
std::string m_decayString
daughter particles selection
std::vector< std::string > m_massConstraintListParticlename
Name of the particles to be mass constraint (massfourC)
std::string m_listName
particle list name
std::vector< int > m_massConstraintList
PDG codes of the particles to be mass constraint (massfourC)
bool m_hasCovMatrix
flag for mother covariance matrix (PseudoFitter)
bool makeKVertexMother(analysis::VertexFitKFit &kv, Particle *p)
Update mother particle after unconstrained vertex fit using KFit.
bool redoTwoPhotonDaughterMassFit(Particle *postFit, const Particle *preFit, const analysis::VertexFitKFit &kv)
Combines preFit particle and vertex information from vertex fit kv to create new postFit particle.
bool makeKMassVertexMother(analysis::MassVertexFitKFit &kv, Particle *p)
Update mother particle after mass-constrained vertex fit using KFit.
bool makeKFourCMother(analysis::FourCFitKFit &kv, Particle *p)
Update mother particle after FourC fit using KFit.
bool fillFitParticles(const Particle *mother, std::vector< const Particle * > &fitChildren, std::vector< const Particle * > &twoPhotonChildren)
Fills valid particle's children (with valid error matrix) in the vector of Particles that will enter ...
bool doKFourCFit(Particle *p)
FourC fit using KFit.
B2Vector3D m_BeamSpotCenter
Beam spot position.
virtual void beginRun() override
Called when entering a new run.
bool doKRecoilMassFit(Particle *p)
RecoilMass fit using KFit.
DBObjPtr< BeamSpot > m_beamSpotDB
Beam spot database object.
bool doVertexFit(Particle *p)
Main steering routine.
std::string m_vertexFitter
Vertex Fitter name.
bool doRaveFit(Particle *mother)
Fit using Rave.
bool makeKRecoilMassMother(analysis::RecoilMassKFit &kf, Particle *p)
Update mother particle after RecoilMass fit using KFit.
bool doKMassFit(Particle *p)
Mass fit using KFit.
double m_recoilMass
recoil mass for constraint
double m_confidenceLevel
required fit confidence level
bool doKMassFourCFit(Particle *p)
MassFourC fit using KFit.
DecayDescriptor m_decaydescriptor
Decay descriptor of decays to look for.
bool fillNotFitParticles(const Particle *mother, std::vector< const Particle * > ¬FitChildren, const std::vector< const Particle * > &fitChildren)
Fills valid particle's children (with valid error matrix) in the vector of Particles that will not en...
void findConstraintBoost(double cut)
calculate iptube constraint (quasi cylinder along boost direction) for RAVE fit
double m_Bfield
magnetic field from data base
bool addChildofParticletoKFit(analysis::FourCFitKFit &kv, const Particle *particle)
Adds given particle's child to the FourCFitKFit.
double getChi2TracksLBoost(const analysis::VertexFitKFit &kv)
calculate the chi2 using only lboost information of tracks
double m_smearing
smearing width applied to IP tube
std::string m_fitType
type of the kinematic fit
void addIPProfileToKFit(analysis::VertexFitKFit &kv)
Adds IPProfile constraint to the vertex fit using KFit.
StoreObjPtr< ParticleList > m_plist
particle list
bool allSelectedDaughters(const Particle *mother, const std::vector< const Particle * > &tracksVertex)
check if all the Daughters (o grand-daughters) are selected for the vertex fit
Class to store reconstructed particles.
void writeExtraInfo(const std::string &name, const double value)
Sets the user defined extraInfo.
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
const V0 * getV0() const
Returns the pointer to the V0 object that was used to create this Particle (if ParticleType == c_V0).
int getPDGCode(void) const
Returns PDG code.
int getProperty() const
Returns particle property as a bit pattern The values are defined in the PropertyFlags enum and descr...
unsigned getNDaughters(void) const
Returns number of daughter particles.
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
void addExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
void set4Vector(const ROOT::Math::PxPyPzEVector &p4)
Sets Lorentz vector.
void updateMomentum(const ROOT::Math::PxPyPzEVector &p4, const ROOT::Math::XYZVector &vertex, const TMatrixFSym &errMatrix, double pValue)
Sets Lorentz vector, position, 7x7 error matrix and p-value.
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
@ c_IsUnspecified
Ordinary particles.
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.
double getFittedVertexX() const
Get the X coordinate of the fitted vertex.
double getFittedVertexY() const
Get the Y coordinate of the fitted vertex.
double getFittedVertexZ() const
Get the Z coordinate of the fitted vertex.
FourCFitKFit is a derived class from KFitBase to perform 4 momentum-constraint kinematical fit.
enum KFitError::ECode setFourMomentum(const ROOT::Math::PxPyPzEVector &m)
Set an 4 Momentum for the FourC-constraint fit.
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a four momentum-constraint fit.
const CLHEP::HepSymMatrix getTrackError(const int id) const
Get an error matrix of the track.
const CLHEP::HepLorentzVector getTrackMomentum(const int id) const
Get a Lorentz vector of the track.
const HepPoint3D getTrackPosition(const int id) const
Get a position of the track.
virtual int getNDF(void) const
Get an NDF of the fit.
enum KFitError::ECode setMagneticField(const double mf)
Change a magnetic field from the default value KFitConst::kDefaultMagneticField.
const KFitTrack getTrack(const int id) const
Get a specified track object.
enum KFitError::ECode addParticle(const Particle *particle)
Add a particle to the fitter.
int getTrackCount(void) const
Get the number of added tracks.
ECode
ECode is a error code enumerate.
KFitTrack is a container of the track information (Lorentz vector, position, and error matrix),...
const CLHEP::HepSymMatrix getError(const int flag=KFitConst::kAfterFit) const
Get an error matrix of the track.
const HepPoint3D getPosition(const int flag=KFitConst::kAfterFit) const
Get a position of the track.
MassFitKFit is a derived class from KFitBase to perform mass-constraint kinematical fit.
enum KFitError::ECode setVertex(const HepPoint3D &v)
Set an initial vertex position for the mass-constraint fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a mass-constraint fit.
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-constraint fit.
enum KFitError::ECode setVertexError(const CLHEP::HepSymMatrix &e)
Set an initial vertex error matrix for the mass-constraint fit.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
MassFourCFitKFit is a derived class from KFitBase to perform mass and 4 momentum-constraint kinematic...
enum KFitError::ECode setFourMomentum(const ROOT::Math::PxPyPzEVector &m)
Set an 4 Momentum for the mass-four-constraint fit.
enum KFitError::ECode addMassConstraint(const double m, std::vector< unsigned > &childTrackId)
Set an invariant mass of daughter particle for the mass-four-momentum-constraint fit.
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a mass-four-momentum-constraint fit.
MassPointingVertexFitKFit is a derived class from KFitBase It performs a kinematical fit with three c...
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a mass-vertex-pointing constraint fit.
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-vertex-pointing constraint fit.
enum KFitError::ECode setProductionVertex(const HepPoint3D &v)
Set the production vertex of the particle.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
MassVertexFitKFit is a derived class from KFitBase to perform mass-vertex-constraint kinematical fit.
enum KFitError::ECode setInitialVertex(const HepPoint3D &v)
Set an initial vertex point for the mass-vertex constraint fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a mass-vertex-constraint fit.
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-vertex constraint fit.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
The RaveKinematicVertexFitter class is part of the RaveInterface together with RaveSetup.
void updateMother()
update the mother particle
void addMother(const Particle *aMotherParticlePtr)
All daughters of the argument of this function will be used as input for the vertex fit.
ROOT::Math::XYZVector getPos()
get the position of the fitted vertex.
void addTrack(const Particle *aParticlePtr)
add a track (in the format of a Belle2::Particle) to set of tracks that should be fitted to a vertex
void setMother(const Particle *aMotherParticlePtr)
Set Mother particle for Vertex/momentum update.
int fit()
do the kinematic vertex fit with all tracks previously added with the addTrack or addMother function.
void setVertFit(bool isVertFit=true)
Set vertex fit: set false in case of mass fit only.
void setMassConstFit(bool isConstFit=true)
Set mass constrained fit
double getPValue()
get the p value of the fitted vertex.
TMatrixDSym getVertexErrorMatrix()
get the covariance matrix (3x3) of the of the fitted vertex position.
void updateDaughters()
update the Daughters particles
void unsetBeamSpot()
unset beam spot constraint
static void initialize(int verbosity=1, double MagneticField=1.5)
Set everything up so everything needed for vertex fitting is there.
static RaveSetup * getInstance()
get the pointer to the instance to get/set any of options stored in RaveSetup
void setBeamSpot(const B2Vector3D &beamSpot, const TMatrixDSym &beamSpotCov)
The beam spot position and covariance is known you can set it here so that and a vertex in the beam s...
void reset()
frees memory allocated by initialize().
The RaveVertexFitter class is part of the RaveInterface together with RaveSetup.
TMatrixDSym getCov(VecSize vertexId=0) const
get the covariance matrix (3x3) of the of the fitted vertex position.
int fit(std::string options="default")
do the vertex fit with all tracks previously added with the addTrack or addMother function.
void addTrack(const Particle *const aParticlePtr)
add a track (in the format of a Belle2::Particle) to set of tracks that should be fitted to a vertex
B2Vector3D getPos(VecSize vertexId=0) const
get the position of the fitted vertex.
double getPValue(VecSize vertexId=0) const
get the p value of the fitted vertex.
RecoilMassKFit is a derived class from KFitBase to perform a kinematical fit with a recoil mass const...
enum KFitError::ECode setFourMomentum(const ROOT::Math::PxPyPzEVector &m)
Set a recoil mass .
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode setRecoilMass(const double m)
Set an invariant mass for the four momentum-constraint fit.
enum KFitError::ECode doFit(void)
Perform a recoil-mass constraint fit.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
VertexFitKFit is a derived class from KFitBase to perform vertex-constraint kinematical fit.
enum KFitError::ECode setIpTubeProfile(const CLHEP::HepLorentzVector &p, const HepPoint3D &x, const CLHEP::HepSymMatrix &e, const double q)
Set a virtual IP-tube track for the vertex constraint fit.
enum KFitError::ECode setInitialVertex(const HepPoint3D &v)
Set an initial vertex point for the vertex-vertex constraint fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
const CLHEP::HepSymMatrix getVertexError(void) const
Get a fitted vertex error matrix.
enum KFitError::ECode doFit(void)
Perform a vertex-constraint fit.
enum KFitError::ECode setIpProfile(const HepPoint3D &ip, const CLHEP::HepSymMatrix &ipe)
Set an IP-ellipsoid shape for the vertex constraint fit.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
B2Vector3< double > B2Vector3D
typedef for common usage with double
void copyDaughters(Particle *mother)
Function copies all (grand-)^n-daughter particles of the argument mother Particle.
Abstract base class for different kinds of events.
static const int kBeforeFit
Input parameter to specify before-fit when setting/getting a track attribute.