10 #include <analysis/modules/ParticleVertexFitter/ParticleVertexFitterModule.h>
13 #include <framework/gearbox/Unit.h>
14 #include <framework/gearbox/Const.h>
15 #include <framework/logging/Logger.h>
16 #include <framework/particledb/EvtGenDatabasePDG.h>
19 #include <analysis/dataobjects/Particle.h>
20 #include <analysis/dataobjects/Btube.h>
22 #include <analysis/utility/CLHEPToROOT.h>
23 #include <analysis/utility/PCmsLabTransform.h>
24 #include <analysis/utility/ParticleCopy.h>
25 #include <analysis/utility/ROOTToCLHEP.h>
28 #include <framework/geometry/BFieldManager.h>
55 setDescription(
"Vertex fitter for modular analysis");
56 setPropertyFlags(c_ParallelProcessingCertified);
59 addParam(
"listName", m_listName,
"name of particle list",
string(
""));
60 addParam(
"confidenceLevel", m_confidenceLevel,
61 "Confidence level to accept the fit. Particle candidates with "
62 "p-value less than confidenceLevel are removed from the particle "
63 "list. If set to -1, all candidates are kept; if set to 0, "
64 "the candidates failing the fit are removed.",
66 addParam(
"vertexFitter", m_vertexFitter,
"KFit or Rave",
string(
"KFit"));
67 addParam(
"fitType", m_fitType,
"type of the kinematic fit (vertex, massvertex, mass)",
string(
"vertex"));
68 addParam(
"withConstraint", m_withConstraint,
69 "additional constraint on vertex: ipprofile, iptube, mother, iptubecut, pointing, btube",
71 addParam(
"decayString", m_decayString,
"specifies which daughter particles are included in the kinematic fit",
string(
""));
72 addParam(
"updateDaughters", m_updateDaughters,
"true: update the daughters after the vertex fit",
false);
73 addParam(
"smearing", m_smearing,
"smear IP tube width by given length", 0.002);
74 addParam(
"massConstraintList", m_massConstraintList,
75 "Type::[int]. List of daughter particles to mass constrain with int = pdg code. (only for MassFourCKFit)", {});
76 addParam(
"massConstraintListParticlename", m_massConstraintListParticlename,
77 "Type::[string]. List of daughter particles to mass constrain with string = particle name. (only for MassFourCKFit)", {});
80 void ParticleVertexFitterModule::initialize()
83 m_plist.isRequired(m_listName);
86 m_Bfield = BFieldManager::getField(TVector3(0, 0, 0)).Z() / Unit::T;
89 if (m_vertexFitter ==
"Rave")
90 analysis::RaveSetup::initialize(1, m_Bfield);
92 B2DEBUG(1,
"ParticleVertexFitterModule : magnetic field = " << m_Bfield);
95 if (m_decayString !=
"")
96 m_decaydescriptor.init(m_decayString);
98 if ((m_massConstraintList.size()) == 0 && (m_massConstraintListParticlename.size()) > 0) {
99 for (
auto& containedParticle : m_massConstraintListParticlename) {
100 TParticlePDG* particletemp = TDatabasePDG::Instance()->GetParticle((containedParticle).c_str());
101 m_massConstraintList.push_back(particletemp->PdgCode());
105 B2INFO(
"ParticleVertexFitter: Performing " << m_fitType <<
" fit on " << m_listName <<
" using " << m_vertexFitter);
106 if (m_decayString !=
"")
107 B2INFO(
"ParticleVertexFitter: Using specified decay string: " << m_decayString);
108 if (m_withConstraint !=
"")
109 B2INFO(
"ParticleVertexFitter: Additional " << m_withConstraint <<
" will be applied");
113 void ParticleVertexFitterModule::beginRun()
120 void ParticleVertexFitterModule::event()
122 if (m_vertexFitter ==
"Rave")
123 analysis::RaveSetup::initialize(1, m_Bfield);
125 m_BeamSpotCenter = m_beamSpotDB->getIPPosition();
126 m_beamSpotCov.ResizeTo(3, 3);
127 TMatrixDSym beamSpotCov(3);
128 if (m_withConstraint ==
"ipprofile") m_beamSpotCov = m_beamSpotDB->getCovVertex();
129 if (m_withConstraint ==
"iptube") {
130 if (m_smearing > 0 && m_vertexFitter ==
"KFit") {
131 ParticleVertexFitterModule::smearBeamSpot(m_smearing);
133 ParticleVertexFitterModule::findConstraintBoost(2.);
136 if (m_withConstraint ==
"iptubecut") {
137 m_BeamSpotCenter = TVector3(0.001, 0., .013);
138 findConstraintBoost(0.03);
140 if ((m_vertexFitter ==
"Rave") && (m_withConstraint ==
"ipprofile" || m_withConstraint ==
"iptube"
141 || m_withConstraint ==
"mother" || m_withConstraint ==
"iptubecut" || m_withConstraint ==
"btube"))
142 analysis::RaveSetup::getInstance()->setBeamSpot(m_BeamSpotCenter, m_beamSpotCov);
143 std::vector<unsigned int> toRemove;
144 unsigned int n = m_plist->getListSize();
145 for (
unsigned i = 0; i < n; i++) {
146 Particle* particle = m_plist->getParticle(i);
147 m_hasCovMatrix =
false;
148 if (m_updateDaughters ==
true) {
149 if (m_decayString.empty()) ParticleCopy::copyDaughters(particle);
150 else B2ERROR(
"Daughters update works only when all daughters are selected. Daughters will not be updated");
153 if (m_withConstraint ==
"mother") {
154 m_BeamSpotCenter = particle->getVertex();
155 m_beamSpotCov = particle->getVertexErrorMatrix();
158 TMatrixFSym mother_errMatrix(7);
159 mother_errMatrix = particle->getMomentumVertexErrorMatrix();
160 for (
int k = 0; k < 7; k++) {
161 for (
int j = 0; j < 7; j++) {
162 if (mother_errMatrix[k][j] > 0) {
163 m_hasCovMatrix =
true;
168 if (m_withConstraint ==
"btube") {
172 toRemove.push_back(particle->getArrayIndex());
180 ok = doVertexFit(particle);
183 particle->setPValue(-1);
184 if (particle->getPValue() < m_confidenceLevel)
185 toRemove.push_back(particle->getArrayIndex());
187 m_plist->removeParticles(toRemove);
190 if (m_vertexFitter ==
"Rave")
191 analysis::RaveSetup::getInstance()->reset();
194 bool ParticleVertexFitterModule::doVertexFit(
Particle* mother)
199 B2FATAL(
"ParticleVertexFitter: No magnetic field");
202 if (m_withConstraint !=
"ipprofile" &&
203 m_withConstraint !=
"iptube" &&
204 m_withConstraint !=
"mother" &&
205 m_withConstraint !=
"iptubecut" &&
206 m_withConstraint !=
"pointing" &&
207 m_withConstraint !=
"btube" &&
208 m_withConstraint !=
"")
209 B2FATAL(
"ParticleVertexFitter: " << m_withConstraint <<
" ***invalid Constraint ");
213 if (m_vertexFitter ==
"KFit") {
216 if (m_fitType ==
"vertex") {
217 if (m_withConstraint ==
"ipprofile") {
218 ok = doKVertexFit(mother,
true,
false);
219 }
else if (m_withConstraint ==
"iptube") {
220 ok = doKVertexFit(mother,
false,
true);
222 ok = doKVertexFit(mother,
false,
false);
227 if (m_fitType ==
"massvertex") {
228 if (m_withConstraint ==
"ipprofile" || m_withConstraint ==
"iptube" || m_withConstraint ==
"iptubecut") {
229 B2FATAL(
"ParticleVertexFitter: Invalid options - mass-constrained fit using KFit does not work with iptube or ipprofile constraint.");
230 }
else if (m_withConstraint ==
"pointing") {
231 ok = doKMassPointingVertexFit(mother);
233 ok = doKMassVertexFit(mother);
238 if (m_fitType ==
"mass") {
239 if (m_withConstraint ==
"ipprofile" || m_withConstraint ==
"iptube" || m_withConstraint ==
"iptubecut") {
240 B2FATAL(
"ParticleVertexFitter: Invalid options - mass fit using KFit does not work with iptube or ipprofile constraint.");
242 ok = doKMassFit(mother);
247 if (m_fitType ==
"fourC") {
248 if (m_withConstraint ==
"ipprofile" || m_withConstraint ==
"iptube" || m_withConstraint ==
"iptubecut") {
249 B2FATAL(
"ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
251 ok = doKFourCFit(mother);
256 if (m_fitType ==
"massfourC") {
257 if (m_withConstraint ==
"ipprofile" || m_withConstraint ==
"iptube" || m_withConstraint ==
"iptubecut") {
258 B2FATAL(
"ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
260 ok = doKMassFourCFit(mother);
265 if (m_fitType !=
"vertex"
266 && m_fitType !=
"massvertex"
267 && m_fitType !=
"mass"
268 && m_fitType !=
"fourC"
269 && m_fitType !=
"massfourC")
270 B2FATAL(
"ParticleVertexFitter: " << m_fitType <<
" ***invalid fit type for the vertex fitter ");
274 if (m_vertexFitter ==
"Rave") {
276 ok = doRaveFit(mother);
277 }
catch (
const rave::CheckedFloatException&) {
278 B2ERROR(
"Invalid inputs (nan/inf)?");
284 if (m_vertexFitter !=
"KFit" && m_vertexFitter !=
"Rave")
285 B2FATAL(
"ParticleVertexFitter: " << m_vertexFitter <<
" ***invalid vertex fitter ");
287 if (!ok)
return false;
296 bool ParticleVertexFitterModule::fillFitParticles(
const Particle* mother, std::vector<const Particle*>& fitChildren,
297 std::vector<const Particle*>& twoPhotonChildren)
299 if (m_decayString.empty()) {
301 for (
unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
302 const Particle* child = mother->getDaughter(ichild);
305 if (mother->getProperty() == Particle::PropertyFlags::c_IsUnspecified and child->getPValue() < 0) {
308 fitChildren.push_back(child);
311 fitChildren = m_decaydescriptor.getSelectionParticles(mother);
314 auto itr = fitChildren.begin();
315 while (itr != fitChildren.end()) {
318 if (child->getPValue() < 0) {
319 B2WARNING(
"Daughter with PDG code " << child->getPDGCode() <<
" does not have a valid error matrix.");
322 bool isTwoPhotonParticle =
false;
323 if (m_hasCovMatrix ==
false) {
324 if (child->getPDGCode() == Const::pi0.getPDGCode() or child->getPDGCode() == 221) {
325 if (child->getNDaughters() == 2) {
326 if (child->getDaughter(0)->getPDGCode() == Const::photon.getPDGCode()
327 && child->getDaughter(1)->getPDGCode() == Const::photon.getPDGCode()) {
328 isTwoPhotonParticle =
true;
333 if (isTwoPhotonParticle) {
335 twoPhotonChildren.push_back(child);
336 itr = fitChildren.erase(itr);
345 bool ParticleVertexFitterModule::redoTwoPhotonDaughterMassFit(
Particle* postFit,
const Particle* preFit,
362 TMatrixFSym errMatrix(3);
363 for (
int i = 0; i < 3; i++)
364 for (
int j = 0; j < 3; j++)
365 errMatrix(i, j) = posErrorMatrix[i][j];
367 g1ErrMatrix.SetSub(4, errMatrix);
368 g2ErrMatrix.SetSub(4, errMatrix);
384 int err = km.
doFit();
390 bool updateDaughters = m_updateDaughters;
391 m_updateDaughters =
false;
392 bool ok = makeKMassMother(km, postFit);
393 m_updateDaughters = updateDaughters;
398 bool ParticleVertexFitterModule::doKVertexFit(
Particle* mother,
bool ipProfileConstraint,
bool ipTubeConstraint)
400 if ((mother->getNDaughters() < 2 && !ipTubeConstraint) || mother->getNDaughters() < 1)
return false;
402 std::vector<const Particle*> fitChildren;
403 std::vector<const Particle*> twoPhotonChildren;
404 bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
409 if (twoPhotonChildren.size() > 1) {
410 B2FATAL(
"[ParticleVertexFitterModule::doKVertexFit] Vertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
413 if ((fitChildren.size() < 2 && !ipTubeConstraint) || fitChildren.size() < 1) {
414 B2WARNING(
"[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is too low.");
422 for (
auto& child : fitChildren)
425 if (ipProfileConstraint)
426 addIPProfileToKFit(kv);
428 if (ipTubeConstraint)
432 int err = kv.
doFit();
437 if (twoPhotonChildren.size() == 0)
439 ok = makeKVertexMother(kv, mother);
440 else if (twoPhotonChildren.size() == 1) {
446 const Particle* twoPhotonDaughter = twoPhotonChildren[0];
448 ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
456 for (
auto& child : fitChildren)
461 if (ipProfileConstraint)
462 addIPProfileToKFit(kv2);
469 ok = makeKVertexMother(kv2, mother);
475 bool ParticleVertexFitterModule::doKMassVertexFit(
Particle* mother)
477 if (mother->getNDaughters() < 2)
return false;
479 std::vector<const Particle*> fitChildren;
480 std::vector<const Particle*> twoPhotonChildren;
481 bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
486 if (twoPhotonChildren.size() > 1) {
487 B2FATAL(
"[ParticleVertexFitterModule::doKVertexFit] MassVertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
490 if (fitChildren.size() < 2) {
491 B2WARNING(
"[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
496 if (twoPhotonChildren.size() == 0) {
501 for (
auto child : fitChildren)
505 int err = kmv.
doFit();
510 ok = makeKMassVertexMother(kmv, mother);
511 }
else if (twoPhotonChildren.size() == 1) {
520 for (
auto child : fitChildren)
524 int err = kv.
doFit();
528 const Particle* twoPhotonDaughter = twoPhotonChildren[0];
530 ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
538 for (
auto child : fitChildren)
548 ok = makeKMassVertexMother(kmv2, mother);
555 bool ParticleVertexFitterModule::doKMassPointingVertexFit(
Particle* mother)
557 if (!(mother->hasExtraInfo(
"prodVertX") && mother->hasExtraInfo(
"prodVertY") && mother->hasExtraInfo(
"prodVertZ"))) {
561 if (mother->getNDaughters() < 2)
return false;
563 std::vector<const Particle*> fitChildren;
564 std::vector<const Particle*> twoPhotonChildren;
565 bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
570 if (twoPhotonChildren.size() > 0) {
571 B2FATAL(
"[ParticleVertexFitterModule::doKMassPointingVertexFit] MassPointingVertex fit using KFit does not support fit with two-photon daughters (yet).");
574 if (fitChildren.size() < 2) {
575 B2WARNING(
"[ParticleVertexFitterModule::doKMassPointingVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
584 for (
auto child : fitChildren)
588 HepPoint3D productionVertex(mother->getExtraInfo(
"prodVertX"),
589 mother->getExtraInfo(
"prodVertY"),
590 mother->getExtraInfo(
"prodVertZ"));
592 int err = kmpv.
doFit();
593 if (err != 0)
return false;
595 ok = makeKMassPointingVertexMother(kmpv, mother);
600 bool ParticleVertexFitterModule::doKMassFit(
Particle* mother)
602 if (mother->getNDaughters() < 2)
return false;
607 for (
unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
608 const Particle* child = mother->getDaughter(ichild);
610 if (child->getPValue() < 0)
return false;
618 int err = km.
doFit();
620 if (err != 0)
return false;
622 bool ok = makeKMassMother(km, mother);
627 bool ParticleVertexFitterModule::doKFourCFit(
Particle* mother)
629 if (mother->getNDaughters() < 2)
return false;
634 for (
unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
635 const Particle* child = mother->getDaughter(ichild);
637 if (child->getNDaughters() > 0) {
638 bool err = addChildofParticletoKFit(kf, child);
639 if (!err)
return false;
641 if (child->getPValue() < 0)
return false;
651 int err = kf.
doFit();
653 if (err != 0)
return false;
655 bool ok = makeKFourCMother(kf, mother);
660 bool ParticleVertexFitterModule::doKMassFourCFit(
Particle* mother)
662 if (mother->getNDaughters() < 2)
return false;
667 for (
unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
668 const Particle* child = mother->getDaughter(ichild);
670 if (child->getNDaughters() > 0) {
671 bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
672 std::abs(child->getPDGCode())) != m_massConstraintList.end();
673 std::vector<unsigned> childId;
674 bool err = addChildofParticletoMassKFit(kf, child, childId);
676 if (!err)
return false;
678 if (child->getPValue() < 0)
return false;
687 int err = kf.
doFit();
689 if (err != 0)
return false;
691 bool ok = makeMassKFourCMother(kf, mother);
701 if (fitError != analysis::KFitError::kNoError)
703 if (m_decayString.empty() && m_updateDaughters ==
true) {
707 std::vector<Particle*> daughters = mother->getDaughters();
710 if (daughters.size() != track_count)
713 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
714 daughters[iChild]->set4Vector(
716 daughters[iChild]->setVertex(
718 daughters[iChild]->setMomentumVertexErrorMatrix(
731 if (fitError != analysis::KFitError::kNoError)
733 if (m_decayString.empty() && m_updateDaughters ==
true) {
737 std::vector<Particle*> daughters = mother->getDaughters();
740 if (daughters.size() != track_count)
743 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
744 daughters[iChild]->set4Vector(
746 daughters[iChild]->setVertex(
748 daughters[iChild]->setMomentumVertexErrorMatrix(
761 if (fitError != analysis::KFitError::kNoError) {
765 if (m_decayString.empty() && m_updateDaughters ==
true) {
769 std::vector<Particle*> daughters = mother->getDaughters();
772 if (daughters.size() != track_count)
775 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
776 daughters[iChild]->set4Vector(
778 daughters[iChild]->setVertex(
780 daughters[iChild]->setMomentumVertexErrorMatrix(
794 if (fitError != analysis::KFitError::kNoError)
796 if (m_decayString.empty() && m_updateDaughters ==
true) {
800 std::vector<Particle*> daughters = mother->getDaughters();
803 if (daughters.size() != track_count)
806 for (
unsigned iChild = 0; iChild < track_count; iChild++) {
807 daughters[iChild]->set4Vector(
809 daughters[iChild]->setVertex(
811 daughters[iChild]->setMomentumVertexErrorMatrix(
825 if (fitError != analysis::KFitError::kNoError)
827 mother->addExtraInfo(
"FourCFitProb", kf.
getCHIsq());
828 mother->addExtraInfo(
"FourCFitChi2", kf.
getNDF());
829 if (m_decayString.empty() && m_updateDaughters ==
true) {
833 std::vector<Particle*> daughters = mother->getDaughters();
835 const unsigned nd = daughters.size();
837 std::vector<std::vector<unsigned>> pars;
838 std::vector<Particle*> allparticles;
839 for (
unsigned ichild = 0; ichild < nd; ichild++) {
840 const Particle* daughter = mother->getDaughter(ichild);
841 std::vector<unsigned> pard;
842 if (daughter->getNDaughters() > 0) {
843 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
844 pars.push_back(pard);
845 allparticles.push_back(daughters[ichild]);
848 pars.push_back(pard);
849 allparticles.push_back(daughters[ichild]);
855 if (l != track_count)
858 for (
unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
859 TLorentzVector childMoms;
861 TMatrixFSym childErrMatrixs(7);
862 for (
unsigned int iChild : pars[iDaug]) {
863 childMoms = childMoms +
864 CLHEPToROOT::getTLorentzVector(
866 childPoss = childPoss +
867 CLHEPToROOT::getTVector3(
869 TMatrixFSym childErrMatrix =
871 childErrMatrixs = childErrMatrixs + childErrMatrix;
873 allparticles[iDaug]->set4Vector(childMoms);
874 allparticles[iDaug]->setVertex(childPoss);
875 allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
886 if (fitError != analysis::KFitError::kNoError)
888 mother->addExtraInfo(
"MassFourCFitProb", TMath::Prob(kf.
getCHIsq(), kf.
getNDF()));
889 mother->addExtraInfo(
"MassFourCFitChi2", kf.
getCHIsq());
890 mother->addExtraInfo(
"MassFourCFitNDF", kf.
getNDF());
891 if (m_decayString.empty() && m_updateDaughters ==
true) {
895 std::vector<Particle*> daughters = mother->getDaughters();
897 const unsigned nd = daughters.size();
899 std::vector<std::vector<unsigned>> pars;
900 std::vector<Particle*> allparticles;
901 for (
unsigned ichild = 0; ichild < nd; ichild++) {
902 const Particle* daughter = mother->getDaughter(ichild);
903 std::vector<unsigned> pard;
904 if (daughter->getNDaughters() > 0) {
905 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
906 pars.push_back(pard);
907 allparticles.push_back(daughters[ichild]);
910 pars.push_back(pard);
911 allparticles.push_back(daughters[ichild]);
917 if (l != track_count)
920 for (
unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
921 TLorentzVector childMoms;
923 TMatrixFSym childErrMatrixs(7);
924 for (
unsigned int iChild : pars[iDaug]) {
925 childMoms = childMoms +
926 CLHEPToROOT::getTLorentzVector(
928 childPoss = childPoss +
929 CLHEPToROOT::getTVector3(
931 TMatrixFSym childErrMatrix =
933 childErrMatrixs = childErrMatrixs + childErrMatrix;
935 allparticles[iDaug]->set4Vector(childMoms);
936 allparticles[iDaug]->setVertex(childPoss);
937 allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
945 void ParticleVertexFitterModule::updateMapOfTrackAndDaughter(
unsigned& l, std::vector<std::vector<unsigned>>& pars,
946 std::vector<unsigned>& parm, std::vector<Particle*>& allparticles,
const Particle* daughter)
948 std::vector <Belle2::Particle*> childs = daughter->getDaughters();
949 for (
unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
950 const Particle* child = daughter->getDaughter(ichild);
951 std::vector<unsigned> pard;
952 if (child->getNDaughters() > 0) {
953 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
954 parm.insert(parm.end(), pard.begin(), pard.end());
955 pars.push_back(pard);
956 allparticles.push_back(childs[ichild]);
960 pars.push_back(pard);
961 allparticles.push_back(childs[ichild]);
968 bool ParticleVertexFitterModule::doRaveFit(
Particle* mother)
970 if ((m_decayString.empty() ||
971 (m_withConstraint ==
"" && m_fitType !=
"mass")) && mother->getNDaughters() < 2)
return false;
972 if (m_withConstraint ==
"") analysis::RaveSetup::getInstance()->unsetBeamSpot();
973 if (m_withConstraint ==
"ipprofile" || m_withConstraint ==
"iptube" || m_withConstraint ==
"mother"
974 || m_withConstraint ==
"iptubecut" || m_withConstraint ==
"btube")
975 analysis::RaveSetup::getInstance()->setBeamSpot(m_BeamSpotCenter, m_beamSpotCov);
978 if (m_fitType ==
"mass") rf.
setVertFit(
false);
980 if (m_decayString.empty()) {
983 std::vector<const Particle*> tracksVertex = m_decaydescriptor.getSelectionParticles(mother);
984 std::vector<std::string> tracksName = m_decaydescriptor.getSelectionNames();
986 if (allSelectedDaughters(mother, tracksVertex)) {
987 for (
auto& itrack : tracksVertex) {
988 if (itrack != mother) rf.
addTrack(itrack);
994 bool mothSel =
false;
996 for (
unsigned itrack = 0; itrack < tracksVertex.size(); itrack++) {
997 if (tracksVertex[itrack] != mother) {
999 B2DEBUG(1,
"ParticleVertexFitterModule: Adding particle " << tracksName[itrack] <<
" to vertex fit ");
1002 if (tracksVertex[itrack] == mother) mothSel =
true;
1007 bool mothIPfit =
false;
1008 if (tracksVertex.size() == 1 && mothSel ==
true && m_withConstraint !=
"" && nTrk == 0) {
1010 if (tracksVertex[0] != mother)
1011 B2FATAL(
"ParticleVertexFitterModule: FATAL Error in IP constrained mother fit");
1018 TMatrixDSym RerrMatrix(3);
1024 for (
auto& itrack : tracksVertex) {
1026 nvert = rsg.
fit(
"kalman");
1029 RerrMatrix = rsg.
getCov(0);
1031 TLorentzVector mom(mother->getMomentum(), mother->getEnergy());
1032 TMatrixDSym errMatrix(7);
1033 for (
int i = 0; i < 7; i++) {
1034 for (
int j = 0; j < 7; j++) {
1035 if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1036 else {errMatrix[i][j] = 0;}
1040 mother->writeExtraInfo(
"prodVertX", pos.X());
1041 mother->writeExtraInfo(
"prodVertY", pos.Y());
1042 mother->writeExtraInfo(
"prodVertZ", pos.Z());
1043 mother->writeExtraInfo(
"prodVertSxx", RerrMatrix[0][0]);
1044 mother->writeExtraInfo(
"prodVertSxy", RerrMatrix[0][1]);
1045 mother->writeExtraInfo(
"prodVertSxz", RerrMatrix[0][2]);
1046 mother->writeExtraInfo(
"prodVertSyx", RerrMatrix[1][0]);
1047 mother->writeExtraInfo(
"prodVertSyy", RerrMatrix[1][1]);
1048 mother->writeExtraInfo(
"prodVertSyz", RerrMatrix[1][2]);
1049 mother->writeExtraInfo(
"prodVertSzx", RerrMatrix[2][0]);
1050 mother->writeExtraInfo(
"prodVertSzy", RerrMatrix[2][1]);
1051 mother->writeExtraInfo(
"prodVertSzz", RerrMatrix[2][2]);
1053 mother->updateMomentum(mom, pos, errMatrix, prob);
1056 }
else {
return false;}
1066 TLorentzVector mom(mother->getMomentum(), mother->getEnergy());
1067 TMatrixDSym errMatrix(7);
1068 for (
int i = 0; i < 7; i++) {
1069 for (
int j = 0; j < 7; j++) {
1070 if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1071 else {errMatrix[i][j] = 0;}
1074 mother->updateMomentum(mom, pos, errMatrix, prob);
1075 }
else {
return false;}
1078 if (mothSel && nTrk > 1) {
1079 analysis::RaveSetup::getInstance()->setBeamSpot(pos, RerrMatrix);
1081 int nKfit = rf.
fit();
1083 analysis::RaveSetup::getInstance()->unsetBeamSpot();
1085 if (nKfit > 0) {
return true;}
1092 if (m_fitType ==
"vertex") {
1094 int nVert = rf.
fit();
1096 if (m_decayString.empty() && m_updateDaughters ==
true) rf.
updateDaughters();
1097 if (nVert != 1)
return false;
1099 if (m_fitType ==
"mass") {
1104 int nVert = rf.
fit();
1106 if (nVert != 1)
return false;
1108 if (m_fitType ==
"massvertex") {
1111 int nVert = rf.
fit();
1113 if (m_decayString.empty() && m_updateDaughters ==
true) rf.
updateDaughters();
1114 if (nVert != 1)
return false;
1117 B2FATAL(
"fitType : " << m_fitType <<
" ***invalid fit type ");
1123 bool ParticleVertexFitterModule::allSelectedDaughters(
const Particle* mother,
1124 const std::vector<const Particle*>& tracksVertex)
1127 if (mother->getNDaughters() == 0)
return false;
1129 int nNotIncluded = mother->getNDaughters();
1131 for (
unsigned i = 0; i < mother->getNDaughters(); i++) {
1133 for (
auto& vi : tracksVertex) {
1134 if (vi == mother->getDaughter(i)) {
1135 nNotIncluded = nNotIncluded - 1;
1140 if (allSelectedDaughters(mother->getDaughter(i), tracksVertex)) nNotIncluded--;
1143 if (nNotIncluded == 0) isAll =
true;
1149 for (
unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1150 const Particle* child = particle->getDaughter(ichild);
1151 if (child->getNDaughters() > 0) addChildofParticletoKFit(kf, child);
1153 if (child->getPValue() < 0)
return false;
1162 std::vector<unsigned>& particleId)
1164 for (
unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1165 const Particle* child = particle->getDaughter(ichild);
1166 if (child->getNDaughters() > 0) {
1167 bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
1168 std::abs(child->getPDGCode())) != m_massConstraintList.end();
1169 std::vector<unsigned> childId;
1170 addChildofParticletoMassKFit(kf, child, childId);
1172 particleId.insert(particleId.end(), childId.begin(), childId.end());
1174 if (child->getPValue() < 0)
return false;
1185 CLHEP::HepSymMatrix covMatrix(3, 0);
1187 for (
int i = 0; i < 3; i++) {
1188 pos[i] = m_BeamSpotCenter(i);
1189 for (
int j = 0; j < 3; j++) {
1190 covMatrix[i][j] = m_beamSpotCov(i, j);
1199 CLHEP::HepSymMatrix err(7, 0);
1201 for (
int i = 0; i < 3; i++) {
1202 for (
int j = 0; j < 3; j++) {
1203 err[i + 4][j + 4] = m_beamSpotCov(i, j);
1211 ROOTToCLHEP::getHepLorentzVector(iptube_mom),
1212 ROOTToCLHEP::getPoint3D(m_BeamSpotCenter),
1217 void ParticleVertexFitterModule::findConstraintBoost(
double cut)
1222 TVector3 boostDir = boost.Unit();
1224 TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1225 beamSpotCov(2, 2) = cut * cut;
1226 double thetab = boostDir.Theta();
1227 double phib = boostDir.Phi();
1229 double stb = TMath::Sin(thetab);
1230 double ctb = TMath::Cos(thetab);
1231 double spb = TMath::Sin(phib);
1232 double cpb = TMath::Cos(phib);
1235 TMatrix rz(3, 3); rz(2, 2) = 1;
1236 rz(0, 0) = cpb; rz(0, 1) = spb;
1237 rz(1, 0) = -1 * spb; rz(1, 1) = cpb;
1239 TMatrix ry(3, 3); ry(1, 1) = 1;
1240 ry(0, 0) = ctb; ry(0, 2) = -1 * stb;
1241 ry(2, 0) = stb; ry(2, 2) = ctb;
1243 TMatrix r(3, 3); r.Mult(rz, ry);
1244 TMatrix rt(3, 3); rt.Transpose(r);
1246 TMatrix TubePart(3, 3); TubePart.Mult(rt, beamSpotCov);
1247 TMatrix Tube(3, 3); Tube.Mult(TubePart, r);
1249 m_beamSpotCov(0, 0) = Tube(0, 0); m_beamSpotCov(0, 1) = Tube(0, 1); m_beamSpotCov(0, 2) = Tube(0, 2);
1250 m_beamSpotCov(1, 0) = Tube(1, 0); m_beamSpotCov(1, 1) = Tube(1, 1); m_beamSpotCov(1, 2) = Tube(1, 2);
1251 m_beamSpotCov(2, 0) = Tube(2, 0); m_beamSpotCov(2, 1) = Tube(2, 1); m_beamSpotCov(2, 2) = Tube(2, 2);
1254 void ParticleVertexFitterModule::smearBeamSpot(
double width)
1256 TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1257 for (
int i = 0; i < 3; i++)
1258 beamSpotCov(i, i) += width * width;
1260 m_beamSpotCov = beamSpotCov;
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.
Class to store reconstructed particles.
float getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
void updateMomentum(const TLorentzVector &p4, const TVector3 &vertex, const TMatrixFSym &errMatrix, float pValue)
Sets Lorentz vector, position, 7x7 error matrix and p-value.
int getPDGCode(void) const
Returns PDG code.
TLorentzVector get4Vector() const
Returns Lorentz vector.
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
FourCFitKFit is a derived class from KFitBase to perform 4 momentum-constraint kinematical fit.
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode setFourMomentum(const TLorentzVector &m)
Set an 4 Momentum for the FourC-constraint 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.
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.
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.
MassFourCFitKFit is a derived class from KFitBase to perform mass and 4 momentum-constraint kinematic...
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 setFourMomentum(const TLorentzVector &m)
Set an 4 Momentum for the mass-four-constraint 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.
MassVertexFitKFit is a derived class from KFitBase to perform mass-vertex-constraint kinematical 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.
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.
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
TVector3 getPos()
get the position of the fitted vertex.
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
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
TVector3 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.
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 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.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.