9 #define SSTR( x ) dynamic_cast< std::ostringstream & >(( std::ostringstream() << std::dec << x ) ).str()
12 #include <analysis/modules/ParticleKinematicFitter/ParticleKinematicFitterModule.h>
13 #include <analysis/OrcaKinFit/BaseFitObject.h>
14 #include <analysis/OrcaKinFit/OPALFitterGSL.h>
15 #include <analysis/OrcaKinFit/JetFitObject.h>
16 #include <analysis/OrcaKinFit/NewtonFitterGSL.h>
17 #include <analysis/OrcaKinFit/NewFitterGSL.h>
18 #include <analysis/OrcaKinFit/PxPyPzMFitObject.h>
20 #include <mdst/dataobjects/ECLCluster.h>
23 #include <framework/gearbox/Const.h>
24 #include <framework/logging/Logger.h>
27 #include <analysis/dataobjects/Particle.h>
30 #include <analysis/utility/PCmsLabTransform.h>
31 #include <analysis/utility/ParticleCopy.h>
33 using namespace CLHEP;
41 namespace OrcaKinFit {
58 setDescription(
"Kinematic fitter for modular analysis");
59 setPropertyFlags(c_ParallelProcessingCertified);
62 addParam(
"listName", m_listName,
"Name of particle list.",
string(
""));
63 addParam(
"kinematicFitter", m_kinematicFitter,
"Available: OrcaKinFit.",
string(
"OrcaKinFit"));
64 addParam(
"orcaFitterEngine", m_orcaFitterEngine,
"OrcaKinFit engine: NewFitterGSL, NewtonFitterGSL, OPALFitterGSL.",
65 string(
"OPALFitterGSL"));
66 addParam(
"orcaTracer", m_orcaTracer,
"OrcaKinFit tracer: None, Text.",
string(
"None"));
67 addParam(
"orcaConstraint", m_orcaConstraint,
"OrcaKinFit constraint: HardBeam, RecoilMass.",
string(
"HardBeam"));
68 addParam(
"debugFitter", m_debugFitter,
"Switch on/off internal debugging output if available.",
false);
69 addParam(
"debugFitterLevel", m_debugFitterLevel,
"Internal debugging output level if available.", 10);
70 addParam(
"addUnmeasuredPhoton", m_addUnmeasuredPhoton,
"Add one unmeasured photon (-3C).",
false);
71 addParam(
"add3CPhoton", m_add3CPhoton,
"Add one photon with unmeasured energy (-1C).",
false);
72 addParam(
"decayString", m_decayString,
"Specifies which daughter particles are included in the kinematic fit.",
string(
""));
73 addParam(
"updateMother", m_updateMother,
"Update the mother kinematics.",
true);
74 addParam(
"updateDaughters", m_updateDaughters,
"Update the daughter kinematics.",
false);
75 addParam(
"recoilMass", m_recoilMass,
"Recoil mass in GeV. RecoilMass constraint only.", 0.0);
76 addParam(
"invMass", m_invMass,
"Invariant mass in GeV. Mass constraint only.", 0.0);
80 void ParticleKinematicFitterModule::initialize()
82 m_eventextrainfo.registerInDataStore();
84 if (m_decayString !=
"") {
85 m_decaydescriptor.init(m_decayString);
86 B2INFO(
"ParticleKinematicFitter: Using specified decay string: " << m_decayString);
89 m_plist.isRequired(m_listName);
92 void ParticleKinematicFitterModule::terminate()
94 B2INFO(
"ParticleKinematicFitterModule::terminate");
97 void ParticleKinematicFitterModule::event()
99 B2DEBUG(17,
"ParticleKinematicFitterModule::event");
101 unsigned int n = m_plist->getListSize();
103 for (
unsigned i = 0; i < n; i++) {
104 Particle* particle = m_plist->getParticle(i);
106 if (m_updateDaughters ==
true) {
107 if (m_decayString.empty()) ParticleCopy::copyDaughters(particle);
108 else B2ERROR(
"Daughters update works only when all daughters are selected. Daughters will not be updated");
111 bool ok = doKinematicFit(particle);
113 if (!ok) particle->setPValue(-1.);
118 bool ParticleKinematicFitterModule::doKinematicFit(
Particle* mother)
121 B2DEBUG(17,
"ParticleKinematicFitterModule::doKinematicFit");
126 if (m_kinematicFitter ==
"OrcaKinFit") {
129 if (m_decayString !=
"") {
130 B2FATAL(
"ParticleKinematicFitterModule: OrcaKinFit does not support yet selection of daughters via decay string!") ;
134 if (!(m_orcaFitterEngine ==
"OPALFitterGSL" or
135 m_orcaFitterEngine ==
"NewtonFitterGSL" or
136 m_orcaFitterEngine ==
"NewFitterGSL")) {
137 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaFitterEngine <<
" is an invalid OrcaKinFit fitter engine!");
141 if (!(m_orcaConstraint ==
"HardBeam" or
142 m_orcaConstraint ==
"RecoilMass" or
143 m_orcaConstraint ==
"Mass")) {
144 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaConstraint <<
" is an invalid OrcaKinFit constraint!");
148 ok = doOrcaKinFitFit(mother);
153 B2FATAL(
"ParticleKinematicFitter: " << m_kinematicFitter <<
" is an invalid kinematic fitter!");
156 if (!ok)
return false;
163 bool ParticleKinematicFitterModule::doOrcaKinFitFit(
Particle* mother)
165 if (mother->getNDaughters() <= 1) {
166 B2WARNING(
"ParticleKinematicFitterModule: Cannot fit with " << mother->getNDaughters() <<
" daughters.");
171 std::vector<Particle*> particleChildren;
172 bool validChildren = fillFitParticles(mother, particleChildren);
174 if (!validChildren) {
175 B2WARNING(
"ParticleKinematicFitterModule: Cannot find valid children for the fit.");
184 if (m_debugFitter) debugfitter = m_debugFitterLevel;
187 if (m_orcaFitterEngine ==
"OPALFitterGSL") {
189 }
else if (m_orcaFitterEngine ==
"NewtonFitterGSL") {
192 }
else if (m_orcaFitterEngine ==
"NewFitterGSL") {
194 (
dynamic_cast<NewFitterGSL*
>(pfitter))->setDebug(debugfitter);
196 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaFitterEngine <<
" is an invalid OrcaKinFit fitter engine!");
200 if (!pfitter)
return false;
210 for (
unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
211 addParticleToOrcaKinFit(fitter, particleChildren[iChild], iChild);
215 if (m_addUnmeasuredPhoton) addUnmeasuredGammaToOrcaKinFit(fitter);
218 addConstraintsToFitter(fitter);
221 addTracerToFitter(fitter);
224 storeOrcaKinFitParticles(
"Measured", fitter, particleChildren, mother);
226 double prob = fitter.fit();
227 double chi2 = fitter.getChi2();
228 int niter = fitter.getIterations();
229 int ndof = fitter.getDoF();
230 int errorcode = fitter.getError();
233 B2DEBUG(17,
"ParticleKinematicFitterModule: -------------------------------------------");
234 B2DEBUG(17,
"ParticleKinematicFitterModule: Fit result of OrcaKinFit using " << m_orcaFitterEngine);
235 B2DEBUG(17,
"ParticleKinematicFitterModule: prob " << prob);
236 B2DEBUG(17,
"ParticleKinematicFitterModule: chi2 " << chi2);
237 B2DEBUG(17,
"ParticleKinematicFitterModule: iterations " << niter);
238 B2DEBUG(17,
"ParticleKinematicFitterModule: ndf " << ndof);
239 B2DEBUG(17,
"ParticleKinematicFitterModule: errorcode " << errorcode);
240 B2DEBUG(17,
"ParticleKinematicFitterModule: -------------------------------------------");
243 if (m_updateMother) updateOrcaKinFitMother(fitter, particleChildren, mother);
246 if (m_updateDaughters) updateOrcaKinFitDaughters(fitter, mother);
249 storeOrcaKinFitParticles(
"Fitted", fitter, particleChildren, mother);
252 mother->addExtraInfo(
"OrcaKinFitProb", prob);
253 mother->setPValue(prob);
254 mother->addExtraInfo(
"OrcaKinFitChi2", chi2);
255 mother->addExtraInfo(
"OrcaKinFitErrorCode", errorcode);
258 if (m_addUnmeasuredPhoton) {
259 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
260 for (
auto fo : *fitObjectContainer) {
261 const std::string name = fo->getName();
262 if (name ==
"unmeasured") {
264 TLorentzVector tlv = getTLorentzVector(fitobject);
265 mother->addExtraInfo(
"OrcaKinFitUnmeasuredTheta", tlv.Theta());
266 mother->addExtraInfo(
"OrcaKinFitUnmeasuredPhi", tlv.Phi());
267 mother->addExtraInfo(
"OrcaKinFitUnmeasuredE", tlv.E());
271 mother->addExtraInfo(
"OrcaKinFitUnmeasuredErrorTheta", getFitObjectError(fitobject, 1));
272 mother->addExtraInfo(
"OrcaKinFitUnmeasuredErrorPhi", getFitObjectError(fitobject, 2));
273 mother->addExtraInfo(
"OrcaKinFitUnmeasuredErrorE", getFitObjectError(fitobject, 0));
284 bool ParticleKinematicFitterModule::fillFitParticles(
Particle* mother, std::vector<Particle*>& particleChildren)
286 for (
unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
287 auto* child =
const_cast<Particle*
>(mother->getDaughter(ichild));
289 if (child->getNDaughters() > 0) {
290 bool err = fillFitParticles(child, particleChildren);
292 B2WARNING(
"ParticleKinematicFitterModule: Cannot find valid children for the fit.");
295 }
else if (child->getPValue() > 0) {
296 particleChildren.push_back(child);
298 B2ERROR(
"Daughter with PDG code " << child->getPDGCode() <<
" does not have a valid p-value: p=" << child->getPValue() <<
", E=" <<
299 child->getEnergy() <<
" GeV");
306 bool ParticleKinematicFitterModule::AddFour(
Particle* mother)
308 TMatrixFSym MomentumVertexErrorMatrix(7);
309 for (
unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
310 auto* child =
const_cast<Particle*
>(mother->getDaughter(ichild));
312 if (child->getPValue() > 0) {
313 MomentumVertexErrorMatrix += child->getMomentumVertexErrorMatrix();
314 }
else if (child->getNDaughters() > 0) {
316 MomentumVertexErrorMatrix += child->getMomentumVertexErrorMatrix();
318 B2ERROR(
"Daughter with PDG code " << child->getPDGCode() <<
" does not have a valid p-value: p=" << child->getPValue() <<
", E=" <<
319 child->getEnergy() <<
" GeV");
323 mother->setMomentumVertexErrorMatrix(MomentumVertexErrorMatrix);
328 void ParticleKinematicFitterModule::addParticleToOrcaKinFit(
BaseFitter& fitter,
Particle* particle,
const int index)
330 B2DEBUG(17,
"ParticleKinematicFitterModule: adding a particle to the fitter!");
332 if (m_add3CPhoton && index == 0) {
333 if (particle -> getPDGCode() != Const::photon.getPDGCode()) {
334 B2ERROR(
"In 3C Kinematic fit, the first daughter should be the Unmeasured Photon!");
337 double startingE = particle -> getECLCluster() -> getEnergy(particle -> getECLClusterEHypothesisBit());
338 double startingPhi = particle -> getECLCluster() -> getPhi();
339 double startingTheta = particle -> getECLCluster() -> getTheta();
341 double startingeE = particle->getECLCluster() -> getUncertaintyEnergy();
342 double startingePhi = particle->getECLCluster() -> getUncertaintyPhi();
343 double startingeTheta = particle->getECLCluster() -> getUncertaintyTheta();
345 B2DEBUG(17, startingPhi <<
" " << startingTheta <<
" " << startingePhi <<
" " << startingeTheta);
348 pfitobject =
new JetFitObject(startingE, startingTheta, startingPhi, startingeE, startingeTheta, startingePhi, 0.);
349 pfitobject->
setParam(0, startingE,
false,
false);
350 pfitobject->
setParam(1, startingTheta,
true,
false);
351 pfitobject->
setParam(2, startingPhi,
true,
false);
353 std::string fitObjectName =
"unmeasured";
354 pfitobject->
setName(fitObjectName.c_str());
358 addFitObjectToConstraints(fitobject);
361 fitter.addFitObject(fitobject);
365 CLHEP::HepLorentzVector clheplorentzvector = getCLHEPLorentzVector(particle);
368 CLHEP::HepSymMatrix clhepmomentumerrormatrix = getCLHEPMomentumErrorMatrix(particle);
372 pfitobject =
new PxPyPzMFitObject(clheplorentzvector, clhepmomentumerrormatrix);
373 std::string fitObjectName =
"particle_" + SSTR(index);
374 pfitobject->
setName(fitObjectName.c_str());
378 addFitObjectToConstraints(fitobject);
381 fitter.addFitObject(fitobject);
389 CLHEP::HepSymMatrix ParticleKinematicFitterModule::getCLHEPMomentumErrorMatrix(
Particle* particle)
391 CLHEP::HepSymMatrix covMatrix(4);
392 TMatrixFSym errMatrix = particle->getMomentumErrorMatrix();
394 for (
int i = 0; i < 4; i++) {
395 for (
int j = 0; j < 4; j++) {
396 covMatrix[i][j] = errMatrix[i][j];
404 CLHEP::HepSymMatrix ParticleKinematicFitterModule::getCLHEPMomentumVertexErrorMatrix(
Particle* particle)
406 CLHEP::HepSymMatrix covMatrix(7);
407 TMatrixFSym errMatrix = particle->getMomentumVertexErrorMatrix();
409 for (
int i = 0; i < 7; i++) {
410 for (
int j = 0; j < 7; j++) {
411 covMatrix[i][j] = errMatrix[i][j];
419 CLHEP::HepLorentzVector ParticleKinematicFitterModule::getCLHEPLorentzVector(
Particle* particle)
421 CLHEP::HepLorentzVector mom(particle->getPx(), particle->getPy(), particle->getPz(), particle->getEnergy());
428 TLorentzVector mom(fitobject->
getPx(), fitobject->
getPy(), fitobject->
getPz(), fitobject->
getE());
434 TMatrixFSym ParticleKinematicFitterModule::getTMatrixFSymMomentumErrorMatrix()
436 TMatrixFSym errMatrix(4);
438 for (
int i = 0; i < 4; i++) {
439 for (
int j = i; j < 4; j++) {
440 errMatrix[i][j] = 0.0;
448 TMatrixFSym ParticleKinematicFitterModule::getTMatrixFSymMomentumVertexErrorMatrix()
450 TMatrixFSym errMatrix(7);
452 for (
int i = 0; i < 7; i++) {
453 for (
int j = i; j < 7; j++) {
454 errMatrix[i][j] = 0.0;
461 TLorentzVector ParticleKinematicFitterModule::getTLorentzVectorConstraints()
464 if (m_orcaConstraint ==
"HardBeam") {
465 TLorentzVector constraints4vector(m_hardConstraintPx.getValue(),
466 m_hardConstraintPy.getValue(),
467 m_hardConstraintPz.getValue(),
468 m_hardConstraintE.getValue());
469 return constraints4vector;
471 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaConstraint <<
" is an invalid OrcaKinFit constraint!");
475 return TLorentzVector(0., 0., 0., 0.);
479 void ParticleKinematicFitterModule::setConstraints()
482 if (m_orcaConstraint ==
"HardBeam") {
491 m_hardConstraintPx.resetFOList();
492 m_hardConstraintPy.resetFOList();
493 m_hardConstraintPz.resetFOList();
494 m_hardConstraintE.resetFOList();
496 m_hardConstraintPx.setName(
"Sum(p_x) [hard]");
497 m_hardConstraintPy.setName(
"Sum(p_y) [hard]");
498 m_hardConstraintPz.setName(
"Sum(p_z) [hard]");
499 m_hardConstraintE.setName(
"Sum(E) [hard]");
501 }
else if (m_orcaConstraint ==
"RecoilMass") {
505 m_hardConstraintRecoilMass =
RecoilMassConstraint(m_recoilMass, boost.Px(), boost.Py(), boost.Pz(), boost.E());
507 m_hardConstraintRecoilMass.resetFOList();
508 m_hardConstraintRecoilMass.setName(
"Recoil Mass [hard]");
510 }
else if (m_orcaConstraint ==
"Mass") {
513 m_hardConstraintMass.resetFOList();
514 m_hardConstraintMass.setName(
"Mass [hard]");
516 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaConstraint <<
" is an invalid OrcaKinFit constraint!");
521 void ParticleKinematicFitterModule::resetFitter(
BaseFitter& fitter)
523 B2DEBUG(17,
"ParticleKinematicFitterModule: Resetting the fitter");
531 if (m_orcaConstraint ==
"HardBeam") {
532 m_hardConstraintPx.addToFOList(fitobject);
533 m_hardConstraintPy.addToFOList(fitobject);
534 m_hardConstraintPz.addToFOList(fitobject);
535 m_hardConstraintE.addToFOList(fitobject);
536 }
else if (m_orcaConstraint ==
"RecoilMass") {
537 m_hardConstraintRecoilMass.addToFOList(fitobject);
538 }
else if (m_orcaConstraint ==
"Mass") {
539 m_hardConstraintMass.addToFOList(fitobject);
541 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaConstraint <<
" is an invalid OrcaKinFit constraint!");
547 void ParticleKinematicFitterModule::addConstraintsToFitter(
BaseFitter& fitter)
549 if (m_orcaConstraint ==
"HardBeam") {
550 fitter.addConstraint(m_hardConstraintPx);
551 fitter.addConstraint(m_hardConstraintPy);
552 fitter.addConstraint(m_hardConstraintPz);
553 fitter.addConstraint(m_hardConstraintE);
554 }
else if (m_orcaConstraint ==
"RecoilMass") {
555 fitter.addConstraint(m_hardConstraintRecoilMass);
556 }
else if (m_orcaConstraint ==
"Mass") {
557 fitter.addConstraint(m_hardConstraintMass);
561 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaConstraint <<
" is an invalid OrcaKinFit constraint!");
566 void ParticleKinematicFitterModule::addTracerToFitter(
BaseFitter& fitter)
568 if (m_orcaTracer ==
"Text") {
570 fitter.setTracer(m_textTracer);
571 }
else if (m_orcaTracer !=
"None") {
572 B2FATAL(
"ParticleKinematicFitterModule: " << m_orcaTracer <<
" is an invalid OrcaKinFit tracer!");
576 void ParticleKinematicFitterModule::addUnmeasuredGammaToOrcaKinFit(
BaseFitter& fitter)
579 TLorentzVector tlv = getTLorentzVectorConstraints();
580 double startingE = tlv.E();
581 double startingPhi = tlv.Phi();
582 double startingTheta = tlv.Theta();
586 pfitobject =
new JetFitObject(startingE, startingTheta, startingPhi, 0.0, 0.0, 0.0, 0.);
587 pfitobject->
setParam(0, startingE,
false,
false);
588 pfitobject->
setParam(1, startingTheta,
false,
false);
589 pfitobject->
setParam(2, startingPhi,
false,
false);
591 std::string fitObjectName =
"unmeasured";
592 pfitobject->
setName(fitObjectName.c_str());
596 addFitObjectToConstraints(fitobject);
599 fitter.addFitObject(fitobject);
606 std::vector <Belle2::Particle*> bDau = mother->getDaughters();
607 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
609 const unsigned nd = bDau.size();
611 std::vector<std::vector<unsigned>> pars;
612 std::vector<Particle*> allparticles;
613 for (
unsigned ichild = 0; ichild < nd; ichild++) {
614 const Particle* daughter = mother->getDaughter(ichild);
615 std::vector<unsigned> pard;
616 if (daughter->getNDaughters() > 0) {
617 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
618 pars.push_back(pard);
619 allparticles.push_back(bDau[ichild]);
622 pars.push_back(pard);
623 allparticles.push_back(bDau[ichild]);
628 if (l == fitObjectContainer->size() - m_addUnmeasuredPhoton) {
630 if (fitter.getError() == 0) {
631 for (
unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
633 TMatrixFSym errMatrixU(7);
634 if (pars[iDaug].size() > 0) {
635 for (
unsigned int iChild : pars[iDaug]) {
638 TLorentzVector tlv_sub = getTLorentzVector(fitobject);
639 TMatrixFSym errMatrixU_sub = getCovMat7(fitobject);
641 errMatrixU = errMatrixU + errMatrixU_sub;
644 B2FATAL(
"ParticleKinematicFitterModule: no fitObject could be used to update the daughter!");
646 TVector3 pos = allparticles[iDaug]->getVertex();
647 TMatrixFSym errMatrix = allparticles[iDaug]->getMomentumVertexErrorMatrix();
648 TMatrixFSym errMatrixMom = allparticles[iDaug]->getMomentumErrorMatrix();
649 TMatrixFSym errMatrixVer = allparticles[iDaug]->getVertexErrorMatrix();
651 for (
int i = 0; i < 3; i++) {
652 for (
int j = i; j < 3; j++) {
653 errMatrixU[i + 4][j + 4] = errMatrixVer[i][j];
657 allparticles[iDaug]->set4Vector(tlv);
658 allparticles[iDaug]->setVertex(pos);
659 allparticles[iDaug]->setMomentumVertexErrorMatrix(errMatrixU);
665 B2ERROR(
"updateOrcaKinFitDaughters: Cannot update daughters, mismatch between number of daughters and number of fitobjects!");
671 void ParticleKinematicFitterModule::updateMapOfTrackAndDaughter(
unsigned& l, std::vector<std::vector<unsigned>>& pars,
672 std::vector<unsigned>& parm, std::vector<Particle*>& allparticles,
const Particle* daughter)
674 std::vector <Belle2::Particle*> dDau = daughter->getDaughters();
675 for (
unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
676 const Particle* child = daughter->getDaughter(ichild);
677 std::vector<unsigned> pard;
678 if (child->getNDaughters() > 0) {
679 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
680 parm.insert(parm.end(), pard.begin(), pard.end());
681 pars.push_back(pard);
682 allparticles.push_back(dDau[ichild]);
686 pars.push_back(pard);
687 allparticles.push_back(dDau[ichild]);
694 void ParticleKinematicFitterModule::updateOrcaKinFitMother(
BaseFitter& fitter, std::vector<Particle*>& particleChildren,
698 TVector3 pos = mother->getVertex();
699 TMatrixFSym errMatrix = mother->getMomentumVertexErrorMatrix();
700 float pvalue = mother->getPValue();
703 TLorentzVector momnew(0., 0., 0., 0.);
705 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
706 for (
unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
709 TLorentzVector tlv = getTLorentzVector(fitobject);
716 mother->updateMomentum(momnew, pos, errMatrix, pvalue);
720 bool ParticleKinematicFitterModule::storeOrcaKinFitParticles(
const std::string& prefix,
BaseFitter& fitter,
721 std::vector<Particle*>& particleChildren,
Particle* mother)
723 bool updated =
false;
724 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
726 for (
unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
729 TLorentzVector tlv = getTLorentzVector(fitobject);
732 std::string extraVariableParticlePx =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_Px";
733 std::string extraVariableParticlePy =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_Py";
734 std::string extraVariableParticlePz =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_Pz";
735 std::string extraVariableParticleE =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_E";
736 std::string extraVariableParticlePxErr =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_PxErr";
737 std::string extraVariableParticlePyErr =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_PyErr";
738 std::string extraVariableParticlePzErr =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_PzErr";
739 std::string extraVariableParticleEErr =
"OrcaKinFit" + prefix +
"_" + SSTR(iChild) +
"_EErr";
741 mother->addExtraInfo(extraVariableParticlePx, tlv.Px());
742 mother->addExtraInfo(extraVariableParticlePy, tlv.Py());
743 mother->addExtraInfo(extraVariableParticlePz, tlv.Pz());
744 mother->addExtraInfo(extraVariableParticleE, tlv.E());
745 mother->addExtraInfo(extraVariableParticlePxErr, getFitObjectError(fitobject, 0));
746 mother->addExtraInfo(extraVariableParticlePyErr, getFitObjectError(fitobject, 1));
747 mother->addExtraInfo(extraVariableParticlePzErr, getFitObjectError(fitobject, 2));
748 mother->addExtraInfo(extraVariableParticleEErr, -1.0);
755 float ParticleKinematicFitterModule::getFitObjectError(
ParticleFitObject* fitobject,
int ilocal)
759 if (pxpypzmfitobject) {
762 B2FATAL(
"ParticleKinematicFitterModule: not implemented yet");
772 if (pxpypzmfitobject) {
774 TMatrixFSym errMatrix(3);
777 for (
int i = 0; i < 3; i++) {
778 for (
int j = 0; j < 3; j++) {
779 errMatrix[i][j] = pxpypzmfitobject->getCov(i, j);
785 B2FATAL(
"ParticleKinematicFitterModule: not implemented yet");
793 TMatrixFSym fitCovMatrix(3);
797 auto* jetfitObject =
static_cast<JetFitObject*
>(fitobject);
800 fitCovMatrix = getFitObjectCovMat(fitobject);
801 TLorentzVector tlv = getTLorentzVector(fitobject);
803 const double energy = tlv.E();
804 const double theta = tlv.Theta();
805 const double phi = tlv.Phi();
807 const double st = sin(theta);
808 const double ct = cos(theta);
809 const double sp = sin(phi);
810 const double cp = cos(phi);
815 A(0, 1) = energy * cp * ct ;
816 A(0, 2) = -energy * sp * st ;
818 A(1, 1) = energy * sp * ct ;
819 A(1, 2) = energy * cp * st ;
821 A(2, 1) = -energy * st ;
827 TMatrixFSym D = fitCovMatrix.Similarity(A);
831 B2FATAL(
"ParticleKinematicFitterModule: not implemented yet");
836 if (pxpypzmfitobject) {
838 fitCovMatrix = getFitObjectCovMat(fitobject);
841 TLorentzVector tlv = getTLorentzVector(fitobject);
853 A[3][0] = tlv.Px() / tlv.E();
854 A[3][1] = tlv.Py() / tlv.E();
855 A[3][2] = tlv.Pz() / tlv.E();
858 TMatrixFSym D = fitCovMatrix.Similarity(A);
862 B2FATAL(
"ParticleKinematicFitterModule: not implemented yet");
In the store you can park objects that have to be accessed by various modules.
virtual const char * getParamName(int ilocal) const =0
Get name of parameter ilocal.
virtual double getError(int ilocal) const
Get error of parameter ilocal.
virtual bool setParam(int ilocal, double par_, bool measured_, bool fixed_=false)
Set value and measured flag of parameter i; return: significant change.
void setName(const char *name_)
Set object's name.
Abstract base class for fitting engines of kinematic fits.
Class for jets with (E, eta, phi) in kinematic fits.
Implements constraint 0 = mass1 - mass2 - m.
Implements a constraint of the form efact*sum(E)+pxfact*sum(px)+pyfact*sum(py)+pzfact*sum(pz)=value.
A kinematic fitter using the Newton-Raphson method to solve the equations.
Description of the fit algorithm and interface:
virtual double getPx() const
Return px.
virtual double getPz() const
Return pz.
virtual double getPy() const
Return py.
virtual double getE() const
Return E.
Class to store reconstructed particles.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.