Belle II Software light-2411-aldebaran
ParticleKinematicFitterModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8#include <sstream>
9#define SSTR( x ) dynamic_cast< std::ostringstream && >(( std::ostringstream() << std::dec << x ) ).str()
10
11// kinfitter
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>
19
20#include <mdst/dataobjects/ECLCluster.h>
21
22// framework utilities
23#include <framework/gearbox/Const.h>
24#include <framework/logging/Logger.h>
25
26// analysis dataobjects
27#include <analysis/dataobjects/Particle.h>
28
29// analysis utilities
30#include <analysis/utility/PCmsLabTransform.h>
31#include <analysis/utility/ParticleCopy.h>
32#include <analysis/ClusterUtility/ClusterUtils.h>
33
34#include <cmath>
35
36using namespace CLHEP;
37using namespace std;
38using namespace Belle2;
39using namespace Belle2::OrcaKinFit;
40
41//-----------------------------------------------------------------
42// Register module
43//-----------------------------------------------------------------
44
45REG_MODULE(ParticleKinematicFitter);
46
47//-----------------------------------------------------------------
48// Implementation
49//-----------------------------------------------------------------
50
51ParticleKinematicFitterModule::ParticleKinematicFitterModule() : Module(), m_textTracer(nullptr), m_eventextrainfo("",
52 DataStore::c_Event)
53{
54 setDescription("Kinematic fitter for modular analysis");
56
57 // Add parameters
58 addParam("listName", m_listName, "Name of particle list.", string(""));
59 addParam("kinematicFitter", m_kinematicFitter, "Available: OrcaKinFit.", string("OrcaKinFit"));
60 addParam("orcaFitterEngine", m_orcaFitterEngine, "OrcaKinFit engine: NewFitterGSL, NewtonFitterGSL, OPALFitterGSL.",
61 string("OPALFitterGSL"));
62 addParam("orcaTracer", m_orcaTracer, "OrcaKinFit tracer: None, Text.", string("None"));
63 addParam("orcaConstraint", m_orcaConstraint, "OrcaKinFit constraint: HardBeam, RecoilMass.", string("HardBeam"));
64 addParam("debugFitter", m_debugFitter, "Switch on/off internal debugging output if available.", false);
65 addParam("debugFitterLevel", m_debugFitterLevel, "Internal debugging output level if available.", 10);
66 addParam("addUnmeasuredPhoton", m_addUnmeasuredPhoton, "Add one unmeasured photon (-3C).", false);
67 addParam("fixUnmeasuredToHER", m_fixUnmeasuredPhotonToHER, "fix the momentum of the unmeasured photon to HER (+2C).", false);
68 addParam("fixUnmeasuredToLER", m_fixUnmeasuredPhotonToLER, "fix the momentum of the unmeasured photon to LER (+2C).", false);
69 addParam("add3CPhoton", m_add3CPhoton, "Add one photon with unmeasured energy (-1C).", false);
70 addParam("liftPhotonTheta", m_liftPhotonTheta, "Lift theta constraint of 3CPhoton. Valid when add3CPhoton is true.", false);
71 addParam("decayString", m_decayString, "Specifies which daughter particles are included in the kinematic fit.", string(""));
72 addParam("updateMother", m_updateMother, "Update the mother kinematics.", true);
73 addParam("updateDaughters", m_updateDaughters, "Update the daughter kinematics.", false);
74 addParam("recoilMass", m_recoilMass, "Recoil mass in GeV. RecoilMass constraint only.", 0.0);
75 addParam("invMass", m_invMass, "Invariant mass in GeV. Mass constraint only.", 0.0);
76 addParam("variablePrefix", m_prefix, "Prefix attached to extra info variables.", string(""));
77
78}
79
81{
82 m_eventextrainfo.registerInDataStore();
83
84 if (m_decayString != "") {
86 B2INFO("ParticleKinematicFitter: Using specified decay string: " << m_decayString);
87 }
88
89 m_plist.isRequired(m_listName);
90}
91
93{
94 B2INFO("ParticleKinematicFitterModule::terminate");
95}
96
98{
99 B2DEBUG(17, "ParticleKinematicFitterModule::event");
100
101 unsigned int n = m_plist->getListSize();
102
103 for (unsigned i = 0; i < n; i++) {
104 Particle* particle = m_plist->getParticle(i);
105
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");
109 }
110
111 bool ok = doKinematicFit(particle);
112
113 if (!ok) particle->setPValue(-1.);
114 }
115}
116
118{
119 B2DEBUG(17, "ParticleKinematicFitterModule::doKinematicFit");
120
121 bool ok = false;
122
123 // fitting with OrcaKinFit
124 if (m_kinematicFitter == "OrcaKinFit") {
125
126 // select subset of particles for the fit
127 if (m_decayString != "") {
128 B2FATAL("ParticleKinematicFitterModule: OrcaKinFit does not support yet selection of daughters via decay string!") ;
129 }
130
131 // check requested fit engine
132 if (!(m_orcaFitterEngine == "OPALFitterGSL" or
133 m_orcaFitterEngine == "NewtonFitterGSL" or
134 m_orcaFitterEngine == "NewFitterGSL")) {
135 B2FATAL("ParticleKinematicFitterModule: " << m_orcaFitterEngine << " is an invalid OrcaKinFit fitter engine!");
136 }
137
138 // check requested constraint
139 if (!(m_orcaConstraint == "HardBeam" or
140 m_orcaConstraint == "RecoilMass" or
141 m_orcaConstraint == "Mass")) {
142 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
143 }
144
145 // basic check is good, go to fitting routine
146 ok = doOrcaKinFitFit(mother);
147 } else { // invalid fitter
148 B2FATAL("ParticleKinematicFitter: " << m_kinematicFitter << " is an invalid kinematic fitter!");
149 }
150
151 if (!ok) return false;
152
153 return true;
154
155}
156
158{
159 if (mother->getNDaughters() <= 1) {
160 B2WARNING("ParticleKinematicFitterModule: Cannot fit with " << mother->getNDaughters() << " daughters.");
161 return false;
162 }
163
164 // fill particles
165 std::vector<Particle*> particleChildren;
166 bool validChildren = fillFitParticles(mother, particleChildren);
167
168 if (!validChildren) {
169 B2WARNING("ParticleKinematicFitterModule: Cannot find valid children for the fit.");
170 return false;
171 }
172
173 // set fit engine
174 BaseFitter* pfitter;
175
176 // internal debugger
177 int debugfitter = 0;
178 if (m_debugFitter) debugfitter = m_debugFitterLevel;
179
180 // choose minimization
181 if (m_orcaFitterEngine == "OPALFitterGSL") {
182 pfitter = new OPALFitterGSL(); // OPAL fitter has no debugger
183 } else if (m_orcaFitterEngine == "NewtonFitterGSL") {
184 pfitter = new NewtonFitterGSL();
185 (dynamic_cast<NewtonFitterGSL*>(pfitter))->setDebug(debugfitter);
186 } else if (m_orcaFitterEngine == "NewFitterGSL") {
187 pfitter = new NewFitterGSL();
188 (dynamic_cast<NewFitterGSL*>(pfitter))->setDebug(debugfitter);
189 } else {
190 B2FATAL("ParticleKinematicFitterModule: " << m_orcaFitterEngine << " is an invalid OrcaKinFit fitter engine!");
191 return false;
192 }
193
194 if (!pfitter) return false;
195 BaseFitter& fitter(*pfitter);
196
197 // reset fitter
198 resetFitter(fitter);
199
200 // set constraints (not connected to a fitter or particles at this point!)
202
203 // add fit particles from particle list to the fitter and to all constraints
204 for (unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
205 addParticleToOrcaKinFit(fitter, particleChildren[iChild], iChild);
206 }
207
208 // add unmeasured photon to the fitter and to all constraints
210
211 // add constraints to the fitter
213
214 // add tracers to the fitter
215 addTracerToFitter(fitter);
216
217 //store information before the fit
218 storeOrcaKinFitParticles("Measured", fitter, particleChildren, mother);
219
220 double prob = fitter.fit();
221 double chi2 = fitter.getChi2();
222 int niter = fitter.getIterations();
223 int ndof = fitter.getDoF();
224 int errorcode = fitter.getError();
225
226 B2DEBUG(17, "ParticleKinematicFitterModule: -------------------------------------------");
227 B2DEBUG(17, "ParticleKinematicFitterModule: Fit result of OrcaKinFit using " << m_orcaFitterEngine);
228 B2DEBUG(17, "ParticleKinematicFitterModule: prob " << prob);
229 B2DEBUG(17, "ParticleKinematicFitterModule: chi2 " << chi2);
230 B2DEBUG(17, "ParticleKinematicFitterModule: iterations " << niter);
231 B2DEBUG(17, "ParticleKinematicFitterModule: ndf " << ndof);
232 B2DEBUG(17, "ParticleKinematicFitterModule: errorcode " << errorcode);
233 B2DEBUG(17, "ParticleKinematicFitterModule: -------------------------------------------");
234
235 // default update mother information
236 if (m_updateMother) updateOrcaKinFitMother(fitter, particleChildren, mother);
237
238 // update daughter information if that is requested
239 if (m_updateDaughters) updateOrcaKinFitDaughters(fitter, mother);
240
241 // store information after the fit
242 storeOrcaKinFitParticles("Fitted", fitter, particleChildren, mother);
243
244 //store general fit results
245 mother->addExtraInfo(m_prefix + "OrcaKinFitProb", prob);
246 mother->setPValue(prob);
247 mother->addExtraInfo(m_prefix + "OrcaKinFitChi2", chi2);
248 mother->addExtraInfo(m_prefix + "OrcaKinFitErrorCode", errorcode);
249
250 // if we added an unmeasured photon, add the kinematics to the mother - at some point we may want to create a particle list from this?
251 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
252 for (auto fo : *fitObjectContainer) {
254 const std::string name = fo->getName();
255 if (name.find("Unmeasured") != std::string::npos) {
256 auto* fitobject = static_cast<ParticleFitObject*>(fo);
257 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
258 mother->addExtraInfo(m_prefix + "OrcaKinFit" + name + "Theta", tlv.Theta());
259 mother->addExtraInfo(m_prefix + "OrcaKinFit" + name + "Phi", tlv.Phi());
260 mother->addExtraInfo(m_prefix + "OrcaKinFit" + name + "E", tlv.E());
261
262 // Uncertainty
263 // const double err0 = getFitObjectError(fitobject, 0);
264 mother->addExtraInfo(m_prefix + "OrcaKinFit" + name + "ErrorTheta", getFitObjectError(fitobject, 1));
265 mother->addExtraInfo(m_prefix + "OrcaKinFit" + name + "ErrorPhi", getFitObjectError(fitobject, 2));
266 mother->addExtraInfo(m_prefix + "OrcaKinFit" + name + "ErrorE", getFitObjectError(fitobject, 0));
267 }
268 }
269 delete fo;
270 }
271
272 delete pfitter;
273 delete m_textTracer;
274 return true;
275}
276
277bool ParticleKinematicFitterModule::fillFitParticles(Particle* mother, std::vector<Particle*>& particleChildren)
278{
279 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
280 auto* child = const_cast<Particle*>(mother->getDaughter(ichild));
281
282 if (child->getNDaughters() > 0) {
283 bool err = fillFitParticles(child, particleChildren);
284 if (!err) {
285 B2WARNING("ParticleKinematicFitterModule: Cannot find valid children for the fit.");
286 return false;
287 }
288 } else if (child->getPValue() > 0) {
289 particleChildren.push_back(child);
290 } else {
291 B2ERROR("Daughter with PDG code " << child->getPDGCode() << " does not have a valid p-value: p=" << child->getPValue() << ", E=" <<
292 child->getEnergy() << " GeV");
293 return false; // error matrix not valid
294 }
295 }
296 return true;
297}
298
300{
301 TMatrixFSym MomentumVertexErrorMatrix(7);
302 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
303 auto* child = const_cast<Particle*>(mother->getDaughter(ichild));
304
305 if (child->getPValue() > 0) {
306 MomentumVertexErrorMatrix += child->getMomentumVertexErrorMatrix();
307 } else if (child->getNDaughters() > 0) {
308 AddFour(child);
309 MomentumVertexErrorMatrix += child->getMomentumVertexErrorMatrix();
310 } else {
311 B2ERROR("Daughter with PDG code " << child->getPDGCode() << " does not have a valid p-value: p=" << child->getPValue() << ", E=" <<
312 child->getEnergy() << " GeV");
313 return false; // error matrix not valid
314 }
315 }
316 mother->setMomentumVertexErrorMatrix(MomentumVertexErrorMatrix);
317 return true;
318}
319
321{
322 B2DEBUG(17, "ParticleKinematicFitterModule: adding a particle to the fitter!");
323
324 if (m_add3CPhoton && index == 0) {
325 if (particle -> getPDGCode() != Const::photon.getPDGCode()) {
326 B2ERROR("In 3C Kinematic fit, the first daughter should be the Unmeasured Photon!");
327 }
328
329 const ECLCluster* cluster = particle->getECLCluster();
330 double startingE = cluster->getEnergy(particle->getECLClusterEHypothesisBit());
331 double startingPhi = cluster->getPhi();
332 double startingTheta = cluster->getTheta();
333
334 ClusterUtils clutls;
335 const auto EPhiThetaCov = clutls.GetCovarianceMatrix3x3FromCluster(cluster);
336 double startingeE = sqrt(fabs(EPhiThetaCov[0][0]));
337 double startingePhi = sqrt(fabs(EPhiThetaCov[1][1]));
338 double startingeTheta = sqrt(fabs(EPhiThetaCov[2][2]));
339
340 B2DEBUG(17, startingPhi << " " << startingTheta << " " << startingePhi << " " << startingeTheta);
341 // create a fit object
342 ParticleFitObject* pfitobject;
343 // memory allocated: it will be deallocated via "delete fo" in doOrcaKinFitFit
344 pfitobject = new JetFitObject(startingE, startingTheta, startingPhi, startingeE, startingeTheta, startingePhi, 0.);
345 pfitobject->setParam(0, startingE, false, false);
347 pfitobject->setParam(1, startingTheta, false, false);
348 else
349 pfitobject->setParam(1, startingTheta, true, false);
350 pfitobject->setParam(2, startingPhi, true, false);
351
352 std::string fitObjectName = "Unmeasured3C";
353 pfitobject->setName(fitObjectName.c_str());
354 ParticleFitObject& fitobject = *pfitobject;
355
356 // add this fit object (=particle) to the constraints
357 addFitObjectToConstraints(fitobject);
358
359 // add fit particle to the fitter
360 fitter.addFitObject(fitobject);
361
362 } else {
363 // four vector
364 CLHEP::HepLorentzVector clheplorentzvector = getCLHEPLorentzVector(particle);
365
366 // error matrix
367 CLHEP::HepSymMatrix clhepmomentumerrormatrix = getCLHEPMomentumErrorMatrix(particle);
368
369 // create the fit object (ParticleFitObject is the base class)
370 ParticleFitObject* pfitobject;
371 // memory allocated: it will be deallocated via "delete fo" in doOrcaKinFitFit
372 pfitobject = new PxPyPzMFitObject(clheplorentzvector, clhepmomentumerrormatrix);
373 std::string fitObjectName = "particle_" + SSTR(index);
374 pfitobject->setName(fitObjectName.c_str());
375 ParticleFitObject& fitobject = *pfitobject;
376
377 // add this fit object (=particle) to the constraints
378 addFitObjectToConstraints(fitobject);
379
380 // add fit particle to the fitter
381 fitter.addFitObject(fitobject);
382 }
383
384 return;
385}
386
388{
389 CLHEP::HepSymMatrix covMatrix(4);
390 TMatrixFSym errMatrix = particle->getMomentumErrorMatrix();
391
392 for (int i = 0; i < 4; i++) {
393 for (int j = 0; j < 4; j++) {
394 covMatrix[i][j] = errMatrix[i][j];
395 }
396 }
397
398 return covMatrix;
399}
400
402{
403 CLHEP::HepSymMatrix covMatrix(7);
404 TMatrixFSym errMatrix = particle->getMomentumVertexErrorMatrix();
405
406 for (int i = 0; i < 7; i++) {
407 for (int j = 0; j < 7; j++) {
408 covMatrix[i][j] = errMatrix[i][j];
409 }
410 }
411
412 return covMatrix;
413}
414
416{
417 CLHEP::HepLorentzVector mom(particle->getPx(), particle->getPy(), particle->getPz(), particle->get4Vector().E());
418 return mom;
419}
420
422{
423 ROOT::Math::PxPyPzEVector mom(fitobject->getPx(), fitobject->getPy(), fitobject->getPz(), fitobject->getE());
424 return mom;
425}
426
427// TMatrixFSym ParticleKinematicFitterModule::getTMatrixFSymMomentumErrorMatrix(ParticleFitObject* fitobject) //fix the warning
429{
430 TMatrixFSym errMatrix(4);
431
432 for (int i = 0; i < 4; i++) {
433 for (int j = i; j < 4; j++) {
434 errMatrix[i][j] = 0.0;
435 }
436 }
437
438 return errMatrix;
439}
440
441// TMatrixFSym ParticleKinematicFitterModule::getTMatrixFSymMomentumVertexErrorMatrix(ParticleFitObject* fitobject) //fix the warning
443{
444 TMatrixFSym errMatrix(7);
445
446 for (int i = 0; i < 7; i++) {
447 for (int j = i; j < 7; j++) {
448 errMatrix[i][j] = 0.0;
449 }
450 }
451
452 return errMatrix;
453}
454
456{
457 if (m_orcaConstraint == "HardBeam") {
458 ROOT::Math::PxPyPzEVector constraints4vector(m_hardConstraintPx.getValue(),
462 return constraints4vector;
463 } else {
464 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
465 }
466
467 // should not reach this point...
468 return ROOT::Math::PxPyPzEVector(0., 0., 0., 0.);
469}
470
472{
473 if (m_orcaConstraint == "HardBeam") {
475 const ROOT::Math::PxPyPzEVector boost = T.getBeamFourMomentum();
476
477 m_hardConstraintPx = MomentumConstraint(0, 1, 0, 0, boost.Px());
478 m_hardConstraintPy = MomentumConstraint(0, 0, 1, 0, boost.Py());
479 m_hardConstraintPz = MomentumConstraint(0, 0, 0, 1, boost.Pz());
480 m_hardConstraintE = MomentumConstraint(1, 0, 0, 0, boost.E());
481
486
487 m_hardConstraintPx.setName("Sum(p_x) [hard]");
488 m_hardConstraintPy.setName("Sum(p_y) [hard]");
489 m_hardConstraintPz.setName("Sum(p_z) [hard]");
490 m_hardConstraintE.setName("Sum(E) [hard]");
491
492 } else if (m_orcaConstraint == "RecoilMass") {
494 const ROOT::Math::PxPyPzEVector boost = T.getBeamFourMomentum();
495
496 m_hardConstraintRecoilMass = RecoilMassConstraint(m_recoilMass, boost.Px(), boost.Py(), boost.Pz(), boost.E());
497
499 m_hardConstraintRecoilMass.setName("Recoil Mass [hard]");
500
501 } else if (m_orcaConstraint == "Mass") {
503
505 m_hardConstraintMass.setName("Mass [hard]");
506 } else {
507 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
508 }
509}
510
512{
513 B2DEBUG(17, "ParticleKinematicFitterModule: Resetting the fitter");
514 fitter.reset();
515}
516
518{
519 if (m_orcaConstraint == "HardBeam") {
524 } else if (m_orcaConstraint == "RecoilMass") {
526 } else if (m_orcaConstraint == "Mass") {
528 } else {
529 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
530 }
531}
532
534{
535 if (m_orcaConstraint == "HardBeam") {
536 fitter.addConstraint(m_hardConstraintPx);
537 fitter.addConstraint(m_hardConstraintPy);
538 fitter.addConstraint(m_hardConstraintPz);
539 fitter.addConstraint(m_hardConstraintE);
540 } else if (m_orcaConstraint == "RecoilMass") {
541 fitter.addConstraint(m_hardConstraintRecoilMass);
542 } else if (m_orcaConstraint == "Mass") {
543 fitter.addConstraint(m_hardConstraintMass);
544 }
545
546 else {
547 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
548 }
549}
550
552{
553 if (m_orcaTracer == "Text") {
554 m_textTracer = new TextTracer(std::cout);
555 fitter.setTracer(m_textTracer);
556 } else if (m_orcaTracer != "None") {
557 B2FATAL("ParticleKinematicFitterModule: " << m_orcaTracer << " is an invalid OrcaKinFit tracer!");
558 }
559}
560
562{
563 B2DEBUG(17, "ParticleKinematicFitterModule::addUnmeasuredGammaToOrcaKinFit: adding an unmeasured photon to the fitter!");
564 // Initialize photon using the existing constraints
565 ROOT::Math::PxPyPzEVector tlv = getLorentzVectorConstraints();
566 double startingE = tlv.E();
567 double startingPhi = tlv.Phi();
568 double startingTheta = tlv.Theta();
569 bool paramFlag = false;
570
571 // create a fit object
572 ParticleFitObject* pfitobject;
573
574 std::string fitObjectName = "UnmeasuredAlongBeam";
575
577 startingTheta = 41.5e-3; // TODO: Read beam crossing angle from db if it's available
578 startingPhi = 0.0;
579 paramFlag = true;
580 } else if (m_fixUnmeasuredPhotonToLER) {
581 startingTheta = TMath::Pi() - 41.5e-3;
582 startingPhi = 0.0;
583 paramFlag = true;
584 } else {
585 fitObjectName = "Unmeasured";
586 }
587
588 // memory allocated: it will be deallocated via "delete fo" in doOrcaKinFitFit
589 pfitobject = new JetFitObject(startingE, startingTheta, startingPhi, 0.0, 0.0, 0.0, 0.);
590 pfitobject->setParam(0, startingE, false, false);
591 pfitobject->setParam(1, startingTheta, paramFlag, paramFlag);
592 pfitobject->setParam(2, startingPhi, paramFlag, paramFlag);
593
594 pfitobject->setName(fitObjectName.c_str());
595 ParticleFitObject& fitobject = *pfitobject;
596
597 // add this fit object (=particle) to the constraints
598 addFitObjectToConstraints(fitobject);
599
600 // add fit particle to the fitter
601 fitter.addFitObject(fitobject);
602}
603
605{
606 std::vector <Belle2::Particle*> bDau = mother->getDaughters();
607 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
608
609 const unsigned nd = bDau.size();
610 unsigned l = 0;
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]);
620 } else {
621 pard.push_back(l);
622 pars.push_back(pard);
623 allparticles.push_back(bDau[ichild]);
624 l++;
625 }
626 }
627
628 if (l == fitObjectContainer->size() - m_addUnmeasuredPhoton) {
629
630 if (fitter.getError() == 0) {
631 for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
632 ROOT::Math::PxPyPzEVector tlv ;
633 TMatrixFSym errMatrixU(7);
634 if (pars[iDaug].size() > 0) {
635 for (unsigned int iChild : pars[iDaug]) {
636 BaseFitObject* fo = fitObjectContainer->at(iChild);
637 auto* fitobject = static_cast<ParticleFitObject*>(fo);
638 ROOT::Math::PxPyPzEVector tlv_sub = getLorentzVector(fitobject);
639 TMatrixFSym errMatrixU_sub = getCovMat7(fitobject);
640 tlv = tlv + tlv_sub;
641 errMatrixU = errMatrixU + errMatrixU_sub;
642 }
643 } else {
644 B2FATAL("ParticleKinematicFitterModule: no fitObject could be used to update the daughter!");
645 }
646 ROOT::Math::XYZVector pos = allparticles[iDaug]->getVertex(); // we don't update the vertex yet
647 TMatrixFSym errMatrix = allparticles[iDaug]->getMomentumVertexErrorMatrix();
648 TMatrixFSym errMatrixMom = allparticles[iDaug]->getMomentumErrorMatrix();
649 TMatrixFSym errMatrixVer = allparticles[iDaug]->getVertexErrorMatrix();
650
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];
654 }
655 }
656
657 allparticles[iDaug]->set4Vector(tlv);
658 allparticles[iDaug]->setVertex(pos);
659 allparticles[iDaug]->setMomentumVertexErrorMatrix(errMatrixU);
660 }
661 }
662
663 return true;
664 } else {
665 B2ERROR("updateOrcaKinFitDaughters: Cannot update daughters, mismatch between number of daughters and number of fitobjects!");
666 return false;
667 }
668}
669
670void ParticleKinematicFitterModule::updateMapOfTrackAndDaughter(unsigned& l, std::vector<std::vector<unsigned>>& pars,
671 std::vector<unsigned>& parm, std::vector<Particle*>& allparticles, const Particle* daughter)
672{
673 std::vector <Belle2::Particle*> dDau = daughter->getDaughters();
674 for (unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
675 const Particle* child = daughter->getDaughter(ichild);
676 std::vector<unsigned> pard;
677 if (child->getNDaughters() > 0) {
678 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
679 parm.insert(parm.end(), pard.begin(), pard.end());
680 pars.push_back(pard);
681 allparticles.push_back(dDau[ichild]);
682 } else {
683 pard.push_back(l);
684 parm.push_back(l);
685 pars.push_back(pard);
686 allparticles.push_back(dDau[ichild]);
687 l++;
688 }
689 }
690}
691
692void ParticleKinematicFitterModule::updateOrcaKinFitMother(BaseFitter& fitter, std::vector<Particle*>& particleChildren,
693 Particle* mother)
694{
695 // get old values
696 ROOT::Math::XYZVector pos = mother->getVertex();
697 TMatrixFSym errMatrix = mother->getMomentumVertexErrorMatrix();
698 float pvalue = mother->getPValue();
699
700 // update momentum vector
701 ROOT::Math::PxPyPzEVector momnew(0., 0., 0., 0.);
702
703 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
704 for (unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
705 BaseFitObject* fo = fitObjectContainer->at(iChild);
706 auto* fitobject = static_cast<ParticleFitObject*>(fo);
707 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
708 momnew += tlv;
709 }
710
711 // update
712 // TODO: use pvalue of the fit or the old one of the mother? use fit covariance matrix?
713 // Maybe here should use the pvalue and errmatrix of the fit ----Yu Hu
714 mother->updateMomentum(momnew, pos, errMatrix, pvalue);
715}
716
718 std::vector<Particle*>& particleChildren, Particle* mother)
719{
720 bool updated = false;
721 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
722
723 for (unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
724 BaseFitObject* fo = fitObjectContainer->at(iChild);
725 auto* fitobject = static_cast<ParticleFitObject*>(fo);
726 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
727
728 // name of extra variables
729 std::string extraVariableParticlePx = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_Px";
730 std::string extraVariableParticlePy = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_Py";
731 std::string extraVariableParticlePz = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_Pz";
732 std::string extraVariableParticleE = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_E";
733 std::string extraVariableParticlePxErr = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_PxErr";
734 std::string extraVariableParticlePyErr = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_PyErr";
735 std::string extraVariableParticlePzErr = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_PzErr";
736 std::string extraVariableParticleEErr = m_prefix + "OrcaKinFit" + fitSuffix + "_" + SSTR(iChild) + "_EErr";
737
738 mother->addExtraInfo(extraVariableParticlePx, tlv.Px());
739 mother->addExtraInfo(extraVariableParticlePy, tlv.Py());
740 mother->addExtraInfo(extraVariableParticlePz, tlv.Pz());
741 mother->addExtraInfo(extraVariableParticleE, tlv.E());
742 mother->addExtraInfo(extraVariableParticlePxErr, getFitObjectError(fitobject, 0));
743 mother->addExtraInfo(extraVariableParticlePyErr, getFitObjectError(fitobject, 1));
744 mother->addExtraInfo(extraVariableParticlePzErr, getFitObjectError(fitobject, 2));
745 mother->addExtraInfo(extraVariableParticleEErr, -1.0);
746
747 }
748
749 return updated;
750}
751
753{
754 //check if it is a PxPyPzMFitObject
755 auto* pxpypzmfitobject = static_cast<PxPyPzMFitObject*>(fitobject);
756 if (pxpypzmfitobject) {
757 return fitobject->getError(ilocal);
758 } else {
759 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
760 }
761}
762
763// this returns the local covariance matrix
765{
766
767 //check if it is a PxPyPzMFitObject
768 auto* pxpypzmfitobject = static_cast<PxPyPzMFitObject*>(fitobject);
769 if (pxpypzmfitobject) {
770
771 TMatrixFSym errMatrix(3);
772
773 //loop over the i-j local variables.
774 for (int i = 0; i < 3; i++) {
775 for (int j = 0; j < 3; j++) {
776 errMatrix[i][j] = pxpypzmfitobject->getCov(i, j);
777 }
778 }
779
780 return errMatrix;
781 } else {
782 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
783 }
784}
785
787{
788 TMatrixFSym fitCovMatrix(3);
789
790 if (strcmp(fitobject->getParamName(0), "E") == 0) {
791 //check if it is a JetFitObject
792 auto* jetfitObject = static_cast<JetFitObject*>(fitobject);
793 if (jetfitObject) {
794
795 fitCovMatrix = getFitObjectCovMat(fitobject);
796 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
797
798 const double energy = tlv.E();
799 const double theta = tlv.Theta();
800 const double phi = tlv.Phi();
801
802 const double st = sin(theta);
803 const double ct = cos(theta);
804 const double sp = sin(phi);
805 const double cp = cos(phi);
806
807 // updated covariance matrix is: A * cov * A^T where A is the Jacobi matrix (use Similarity)
808 TMatrixF A(7, 3);
809 A(0, 0) = cp * st ; // dpx/dE
810 A(0, 1) = energy * cp * ct ; // dpx/dtheta
811 A(0, 2) = -energy * sp * st ; // dpx/dphi
812 A(1, 0) = sp * st ; // dpy/dE
813 A(1, 1) = energy * sp * ct ; // dpz/dtheta
814 A(1, 2) = energy * cp * st ; // dpy/dphi
815 A(2, 0) = ct ; // dpz/dE
816 A(2, 1) = -energy * st ; // dpz/dtheta
817 A(2, 2) = 0 ; // dpz/dphi
818 A(3, 0) = 1.0; // dE/dE
819 A(3, 1) = 0.0; // dE/dphi
820 A(3, 2) = 0.0; // dE/dtheta
821
822 TMatrixFSym D = fitCovMatrix.Similarity(A);
823 return D;
824
825 } else {
826 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
827 }
828 } else {
829 //check if it is a PxPyPzMFitObject
830 auto* pxpypzmfitobject = static_cast<PxPyPzMFitObject*>(fitobject);
831 if (pxpypzmfitobject) {
832
833 fitCovMatrix = getFitObjectCovMat(fitobject);
834
835 // updated covariance matrix is: A * cov * A^T where A is the Jacobi matrix (use Similarity)
836 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
837 TMatrixF A(7, 3);
838 A[0][0] = 1.; // px/dpx
839 A[0][1] = 0.; // px/dpy
840 A[0][2] = 0.; // px/dpz
841 A[1][0] = 0.; // py/dpx
842 A[1][1] = 1.; // py/dpy
843 A[1][2] = 0.; // py/dpz
844 A[2][0] = 0.; // pz/dpx
845 A[2][1] = 0.; // pz/dpy
846 A[2][2] = 1.; // pz/dpz
847 if (tlv.E() > 0.0) {
848 A[3][0] = tlv.Px() / tlv.E(); // E/dpx, E=sqrt(px^2 + py^2 + pz^2 + m^2)
849 A[3][1] = tlv.Py() / tlv.E(); // E/dpy
850 A[3][2] = tlv.Pz() / tlv.E(); // E/dpz
851 }
852
853 TMatrixFSym D = fitCovMatrix.Similarity(A);
854
855 return D;
856 } else {
857 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
858 }
859 }
860}
Class to provide momentum-related information from ECLClusters.
Definition: ClusterUtils.h:36
const TMatrixDSym GetCovarianceMatrix3x3FromCluster(const ECLCluster *cluster)
Returns 3x3 covariance matrix (E, theta, phi)
static const ParticleType photon
photon particle
Definition: Const.h:673
In the store you can park objects that have to be accessed by various modules.
Definition: DataStore.h:51
bool init(const std::string &str)
Initialise the DecayDescriptor from given string.
ECL cluster data.
Definition: ECLCluster.h:27
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
void setName(const char *name_)
Set object's name.
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.
Definition: BaseFitter.h:47
Class for jets with (E, eta, phi) in kinematic fits.
Definition: JetFitObject.h:43
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.
virtual double getValue() const override
Returns the value of the constraint.
A kinematic fitter using the Newton-Raphson method to solve the equations.
Definition: NewFitterGSL.h:51
Description of the fit algorithm and interface:
Definition: OPALFitterGSL.h:88
virtual void addToFOList(ParticleFitObject &fitobject, int flag=1)
Adds one ParticleFitObject objects to the list.
virtual void resetFOList()
Resests ParticleFitObject list.
virtual double getPx() const
Return px.
virtual double getPz() const
Return pz.
virtual double getPy() const
Return py.
virtual double getE() const
Return E.
bool fillFitParticles(Particle *mother, std::vector< Particle * > &particleChildren)
Fills valid particle's children (with valid error matrix) in the vector of Particles that will enter ...
void addFitObjectToConstraints(ParticleFitObject &fitobject)
Adds Orca fit object to the constraints.
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 m_liftPhotonTheta
lift theta constraint of the 3CPhoton.
bool AddFour(Particle *mother)
Added four vectors and calculated a covariance matrix for a combined particles.
std::string m_orcaConstraint
Constraint (softBeam, hardBeam (default))
int m_debugFitterLevel
internal debugging level (for New and Newton fitter only)
virtual void initialize() override
Initialize the Module.
CLHEP::HepSymMatrix getCLHEPMomentumVertexErrorMatrix(Particle *particle)
Returns particle's 7x7 momentum-vertex-error matrix as a HepSymMatrix.
MomentumConstraint m_hardConstraintE
hard beam constraint E
std::string m_prefix
prefix attached to extrainfo names
bool m_debugFitter
activate internal debugging (for New and Newton fitter only)
void addUnmeasuredGammaToOrcaKinFit(BaseFitter &fitter)
Adds an unmeasured gamma (E, phi, theta) to the fit (-3C) stored as EventExtraInfo TODO.
MomentumConstraint m_hardConstraintPz
hard beam constraint pz
float getFitObjectError(ParticleFitObject *fitobject, int ilocal)
Returns fit object error on the parameter ilocal.
ROOT::Math::PxPyPzEVector getLorentzVectorConstraints()
Get constraints (at whatever stage before/after fitting)
bool doOrcaKinFitFit(Particle *p)
Kinematic fit using OrcaKinFit.
bool m_add3CPhoton
add one photon with unmeasured energy to the fit (costs 1 constraints)
TMatrixFSym getTMatrixFSymMomentumErrorMatrix()
Returns particle's 7x7 momentum-error matrix as a TMatrixFSym.
void addConstraintsToFitter(BaseFitter &fitter)
Adds Orca fit object to the constraints.
void addTracerToFitter(BaseFitter &fitter)
Adds tracer to the fitter.
RecoilMassConstraint m_hardConstraintRecoilMass
hard recoil mass constraint
TMatrixFSym getTMatrixFSymMomentumVertexErrorMatrix()
Returns particle's 7x7 momentum-vertex-error matrix as a TMatrixFSym.
StoreObjPtr< EventExtraInfo > m_eventextrainfo
StoreObjPtr for the EventExtraInfo in this mode.
MomentumConstraint m_hardConstraintPy
hard beam constraint py
double m_recoilMass
Recoil mass for RecoilMass constraint.
CLHEP::HepLorentzVector getCLHEPLorentzVector(Particle *particle)
Returns particle's 4-momentum as a HepLorentzVector.
TextTracer * m_textTracer
internal text output variable
DecayDescriptor m_decaydescriptor
Decay descriptor of decays to look for.
TMatrixFSym getCovMat7(ParticleFitObject *fitobject)
Returns covariance matrix.
ROOT::Math::PxPyPzEVector getLorentzVector(ParticleFitObject *fitobject)
bool updateOrcaKinFitDaughters(BaseFitter &fitter, Particle *mother)
Update the daughters: momentum is sum of daughters TODO update covariance matrix.
void setConstraints()
Sets constraints, this is not connect to particles or a fitter at this stage.
bool m_addUnmeasuredPhoton
add one unmeasured photon to the fit (costs 3 constraints)
TMatrixFSym getFitObjectCovMat(ParticleFitObject *fitobject)
Returns covariance matrix.
bool storeOrcaKinFitParticles(const std::string &fitSuffix, BaseFitter &fitter, std::vector< Particle * > &particleChildren, Particle *mother)
store fit object information as ExtraInfo
CLHEP::HepSymMatrix getCLHEPMomentumErrorMatrix(Particle *particle)
Returns particle's 4x4 momentum-error matrix as a HepSymMatrix.
bool m_fixUnmeasuredPhotonToHER
fix the momentum of the unmeasured photon to HER
StoreObjPtr< ParticleList > m_plist
StoreObjPtr for the particle list.
void addParticleToOrcaKinFit(BaseFitter &fitter, Particle *particle, const int index)
Adds given particle to the OrcaKinFit.
void updateOrcaKinFitMother(BaseFitter &fitter, std::vector< Particle * > &particleChildren, Particle *mother)
Update the mother: momentum is sum of daughters TODO update covariance matrix.
void resetFitter(BaseFitter &fitter)
Resets all objects associated with the OrcaKinFit fitter.
bool m_fixUnmeasuredPhotonToLER
fix the momentum of the unmeasured photon to LER
bool doKinematicFit(Particle *p)
Main steering routine for any kinematic fitter.
MomentumConstraint m_hardConstraintPx
hard beam constraint px
Class to hold Lorentz transformations from/to CMS and boost vector.
ROOT::Math::PxPyPzEVector getBeamFourMomentum() const
Returns LAB four-momentum of e+e-, i.e.
Class to store reconstructed particles.
Definition: Particle.h:75
double getPx() const
Returns x component of momentum.
Definition: Particle.h:587
const ECLCluster * getECLCluster() const
Returns the pointer to the ECLCluster object that was used to create this Particle (if ParticleType =...
Definition: Particle.cc:891
double getPz() const
Returns z component of momentum.
Definition: Particle.h:605
double getPValue() const
Returns chi^2 probability of fit if done or -1.
Definition: Particle.h:667
ROOT::Math::XYZVector getVertex() const
Returns vertex position (POCA for charged, IP for neutral FS particles)
Definition: Particle.h:631
double getPy() const
Returns y component of momentum.
Definition: Particle.h:596
unsigned getNDaughters(void) const
Returns number of daughter particles.
Definition: Particle.h:727
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
Definition: Particle.cc:637
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:547
TMatrixFSym getMomentumErrorMatrix() const
Returns the 4x4 momentum error matrix.
Definition: Particle.cc:435
void addExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1336
void setMomentumVertexErrorMatrix(const TMatrixFSym &errMatrix)
Sets 7x7 error matrix.
Definition: Particle.cc:393
void setPValue(double pValue)
Sets chi^2 probability of fit.
Definition: Particle.h:366
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.
Definition: Particle.h:386
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
Definition: Particle.cc:420
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:631
ECLCluster::EHypothesisBit getECLClusterEHypothesisBit() const
Returns the ECLCluster EHypothesisBit for this Particle.
Definition: Particle.h:1001
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
void copyDaughters(Particle *mother)
Function copies all (grand-)^n-daughter particles of the argument mother Particle.
Definition: ParticleCopy.cc:56
Abstract base class for different kinds of events.
Definition: ClusterUtils.h:24
STL namespace.