Belle II Software development
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
77}
78
80{
81 m_eventextrainfo.registerInDataStore();
82
83 if (m_decayString != "") {
85 B2INFO("ParticleKinematicFitter: Using specified decay string: " << m_decayString);
86 }
87
88 m_plist.isRequired(m_listName);
89}
90
92{
93 B2INFO("ParticleKinematicFitterModule::terminate");
94}
95
97{
98 B2DEBUG(17, "ParticleKinematicFitterModule::event");
99
100 unsigned int n = m_plist->getListSize();
101
102 for (unsigned i = 0; i < n; i++) {
103 Particle* particle = m_plist->getParticle(i);
104
105 if (m_updateDaughters == true) {
106 if (m_decayString.empty()) ParticleCopy::copyDaughters(particle);
107 else B2ERROR("Daughters update works only when all daughters are selected. Daughters will not be updated");
108 }
109
110 bool ok = doKinematicFit(particle);
111
112 if (!ok) particle->setPValue(-1.);
113 }
114}
115
117{
118 B2DEBUG(17, "ParticleKinematicFitterModule::doKinematicFit");
119
120 bool ok = false;
121
122 // fitting with OrcaKinFit
123 if (m_kinematicFitter == "OrcaKinFit") {
124
125 // select subset of particles for the fit
126 if (m_decayString != "") {
127 B2FATAL("ParticleKinematicFitterModule: OrcaKinFit does not support yet selection of daughters via decay string!") ;
128 }
129
130 // check requested fit engine
131 if (!(m_orcaFitterEngine == "OPALFitterGSL" or
132 m_orcaFitterEngine == "NewtonFitterGSL" or
133 m_orcaFitterEngine == "NewFitterGSL")) {
134 B2FATAL("ParticleKinematicFitterModule: " << m_orcaFitterEngine << " is an invalid OrcaKinFit fitter engine!");
135 }
136
137 // check requested constraint
138 if (!(m_orcaConstraint == "HardBeam" or
139 m_orcaConstraint == "RecoilMass" or
140 m_orcaConstraint == "Mass")) {
141 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
142 }
143
144 // basic check is good, go to fitting routine
145 ok = doOrcaKinFitFit(mother);
146 } else { // invalid fitter
147 B2FATAL("ParticleKinematicFitter: " << m_kinematicFitter << " is an invalid kinematic fitter!");
148 }
149
150 if (!ok) return false;
151
152 return true;
153
154}
155
157{
158 if (mother->getNDaughters() <= 1) {
159 B2WARNING("ParticleKinematicFitterModule: Cannot fit with " << mother->getNDaughters() << " daughters.");
160 return false;
161 }
162
163 // fill particles
164 std::vector<Particle*> particleChildren;
165 bool validChildren = fillFitParticles(mother, particleChildren);
166
167 if (!validChildren) {
168 B2WARNING("ParticleKinematicFitterModule: Cannot find valid children for the fit.");
169 return false;
170 }
171
172 // set fit engine
173 BaseFitter* pfitter;
174
175 // internal debugger
176 int debugfitter = 0;
177 if (m_debugFitter) debugfitter = m_debugFitterLevel;
178
179 // choose minimization
180 if (m_orcaFitterEngine == "OPALFitterGSL") {
181 pfitter = new OPALFitterGSL(); // OPAL fitter has no debugger
182 } else if (m_orcaFitterEngine == "NewtonFitterGSL") {
183 pfitter = new NewtonFitterGSL();
184 (dynamic_cast<NewtonFitterGSL*>(pfitter))->setDebug(debugfitter);
185 } else if (m_orcaFitterEngine == "NewFitterGSL") {
186 pfitter = new NewFitterGSL();
187 (dynamic_cast<NewFitterGSL*>(pfitter))->setDebug(debugfitter);
188 } else {
189 B2FATAL("ParticleKinematicFitterModule: " << m_orcaFitterEngine << " is an invalid OrcaKinFit fitter engine!");
190 return false;
191 }
192
193 if (!pfitter) return false;
194 BaseFitter& fitter(*pfitter);
195
196 // reset fitter
197 resetFitter(fitter);
198
199 // set constraints (not connected to a fitter or particles at this point!)
201
202 // add fit particles from particle list to the fitter and to all constraints
203 for (unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
204 addParticleToOrcaKinFit(fitter, particleChildren[iChild], iChild);
205 }
206
207 // add unmeasured photon to the fitter and to all constraints
209
210 // add constraints to the fitter
212
213 // add tracers to the fitter
214 addTracerToFitter(fitter);
215
216 //store information before the fit
217 storeOrcaKinFitParticles("Measured", fitter, particleChildren, mother);
218
219 double prob = fitter.fit();
220 double chi2 = fitter.getChi2();
221 int niter = fitter.getIterations();
222 int ndof = fitter.getDoF();
223 int errorcode = fitter.getError();
224
225 B2DEBUG(17, "ParticleKinematicFitterModule: -------------------------------------------");
226 B2DEBUG(17, "ParticleKinematicFitterModule: Fit result of OrcaKinFit using " << m_orcaFitterEngine);
227 B2DEBUG(17, "ParticleKinematicFitterModule: prob " << prob);
228 B2DEBUG(17, "ParticleKinematicFitterModule: chi2 " << chi2);
229 B2DEBUG(17, "ParticleKinematicFitterModule: iterations " << niter);
230 B2DEBUG(17, "ParticleKinematicFitterModule: ndf " << ndof);
231 B2DEBUG(17, "ParticleKinematicFitterModule: errorcode " << errorcode);
232 B2DEBUG(17, "ParticleKinematicFitterModule: -------------------------------------------");
233
234 // default update mother information
235 if (m_updateMother) updateOrcaKinFitMother(fitter, particleChildren, mother);
236
237 // update daughter information if that is requested
238 if (m_updateDaughters) updateOrcaKinFitDaughters(fitter, mother);
239
240 // store information after the fit
241 storeOrcaKinFitParticles("Fitted", fitter, particleChildren, mother);
242
243 //store general fit results
244 mother->addExtraInfo("OrcaKinFitProb", prob);
245 mother->setPValue(prob);
246 mother->addExtraInfo("OrcaKinFitChi2", chi2);
247 mother->addExtraInfo("OrcaKinFitErrorCode", errorcode);
248
249 // 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?
250 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
251 for (auto fo : *fitObjectContainer) {
253 const std::string name = fo->getName();
254 if (name.find("Unmeasured") != std::string::npos) {
255 auto* fitobject = static_cast<ParticleFitObject*>(fo);
256 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
257 mother->addExtraInfo("OrcaKinFit" + name + "Theta", tlv.Theta());
258 mother->addExtraInfo("OrcaKinFit" + name + "Phi", tlv.Phi());
259 mother->addExtraInfo("OrcaKinFit" + name + "E", tlv.E());
260
261 // Uncertainty
262 // const double err0 = getFitObjectError(fitobject, 0);
263 mother->addExtraInfo("OrcaKinFit" + name + "ErrorTheta", getFitObjectError(fitobject, 1));
264 mother->addExtraInfo("OrcaKinFit" + name + "ErrorPhi", getFitObjectError(fitobject, 2));
265 mother->addExtraInfo("OrcaKinFit" + name + "ErrorE", getFitObjectError(fitobject, 0));
266 }
267 }
268 delete fo;
269 }
270
271 delete pfitter;
272 delete m_textTracer;
273 return true;
274}
275
276bool ParticleKinematicFitterModule::fillFitParticles(Particle* mother, std::vector<Particle*>& particleChildren)
277{
278 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
279 auto* child = const_cast<Particle*>(mother->getDaughter(ichild));
280
281 if (child->getNDaughters() > 0) {
282 bool err = fillFitParticles(child, particleChildren);
283 if (!err) {
284 B2WARNING("ParticleKinematicFitterModule: Cannot find valid children for the fit.");
285 return false;
286 }
287 } else if (child->getPValue() > 0) {
288 particleChildren.push_back(child);
289 } else {
290 B2ERROR("Daughter with PDG code " << child->getPDGCode() << " does not have a valid p-value: p=" << child->getPValue() << ", E=" <<
291 child->getEnergy() << " GeV");
292 return false; // error matrix not valid
293 }
294 }
295 return true;
296}
297
299{
300 TMatrixFSym MomentumVertexErrorMatrix(7);
301 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
302 auto* child = const_cast<Particle*>(mother->getDaughter(ichild));
303
304 if (child->getPValue() > 0) {
305 MomentumVertexErrorMatrix += child->getMomentumVertexErrorMatrix();
306 } else if (child->getNDaughters() > 0) {
307 AddFour(child);
308 MomentumVertexErrorMatrix += child->getMomentumVertexErrorMatrix();
309 } else {
310 B2ERROR("Daughter with PDG code " << child->getPDGCode() << " does not have a valid p-value: p=" << child->getPValue() << ", E=" <<
311 child->getEnergy() << " GeV");
312 return false; // error matrix not valid
313 }
314 }
315 mother->setMomentumVertexErrorMatrix(MomentumVertexErrorMatrix);
316 return true;
317}
318
320{
321 B2DEBUG(17, "ParticleKinematicFitterModule: adding a particle to the fitter!");
322
323 if (m_add3CPhoton && index == 0) {
324 if (particle -> getPDGCode() != Const::photon.getPDGCode()) {
325 B2ERROR("In 3C Kinematic fit, the first daughter should be the Unmeasured Photon!");
326 }
327
328 const ECLCluster* cluster = particle->getECLCluster();
329 double startingE = cluster->getEnergy(particle->getECLClusterEHypothesisBit());
330 double startingPhi = cluster->getPhi();
331 double startingTheta = cluster->getTheta();
332
333 ClusterUtils clutls;
334 const auto EPhiThetaCov = clutls.GetCovarianceMatrix3x3FromCluster(cluster);
335 double startingeE = sqrt(fabs(EPhiThetaCov[0][0]));
336 double startingePhi = sqrt(fabs(EPhiThetaCov[1][1]));
337 double startingeTheta = sqrt(fabs(EPhiThetaCov[2][2]));
338
339 B2DEBUG(17, startingPhi << " " << startingTheta << " " << startingePhi << " " << startingeTheta);
340 // create a fit object
341 ParticleFitObject* pfitobject;
342 // memory allocated: it will be deallocated via "delete fo" in doOrcaKinFitFit
343 pfitobject = new JetFitObject(startingE, startingTheta, startingPhi, startingeE, startingeTheta, startingePhi, 0.);
344 pfitobject->setParam(0, startingE, false, false);
346 pfitobject->setParam(1, startingTheta, false, false);
347 else
348 pfitobject->setParam(1, startingTheta, true, false);
349 pfitobject->setParam(2, startingPhi, true, false);
350
351 std::string fitObjectName = "Unmeasured3C";
352 pfitobject->setName(fitObjectName.c_str());
353 ParticleFitObject& fitobject = *pfitobject;
354
355 // add this fit object (=particle) to the constraints
356 addFitObjectToConstraints(fitobject);
357
358 // add fit particle to the fitter
359 fitter.addFitObject(fitobject);
360
361 } else {
362 // four vector
363 CLHEP::HepLorentzVector clheplorentzvector = getCLHEPLorentzVector(particle);
364
365 // error matrix
366 CLHEP::HepSymMatrix clhepmomentumerrormatrix = getCLHEPMomentumErrorMatrix(particle);
367
368 // create the fit object (ParticleFitObject is the base class)
369 ParticleFitObject* pfitobject;
370 // memory allocated: it will be deallocated via "delete fo" in doOrcaKinFitFit
371 pfitobject = new PxPyPzMFitObject(clheplorentzvector, clhepmomentumerrormatrix);
372 std::string fitObjectName = "particle_" + SSTR(index);
373 pfitobject->setName(fitObjectName.c_str());
374 ParticleFitObject& fitobject = *pfitobject;
375
376 // add this fit object (=particle) to the constraints
377 addFitObjectToConstraints(fitobject);
378
379 // add fit particle to the fitter
380 fitter.addFitObject(fitobject);
381 }
382
383 return;
384}
385
387{
388 CLHEP::HepSymMatrix covMatrix(4);
389 TMatrixFSym errMatrix = particle->getMomentumErrorMatrix();
390
391 for (int i = 0; i < 4; i++) {
392 for (int j = 0; j < 4; j++) {
393 covMatrix[i][j] = errMatrix[i][j];
394 }
395 }
396
397 return covMatrix;
398}
399
401{
402 CLHEP::HepSymMatrix covMatrix(7);
403 TMatrixFSym errMatrix = particle->getMomentumVertexErrorMatrix();
404
405 for (int i = 0; i < 7; i++) {
406 for (int j = 0; j < 7; j++) {
407 covMatrix[i][j] = errMatrix[i][j];
408 }
409 }
410
411 return covMatrix;
412}
413
415{
416 CLHEP::HepLorentzVector mom(particle->getPx(), particle->getPy(), particle->getPz(), particle->get4Vector().E());
417 return mom;
418}
419
421{
422 ROOT::Math::PxPyPzEVector mom(fitobject->getPx(), fitobject->getPy(), fitobject->getPz(), fitobject->getE());
423 return mom;
424}
425
426// TMatrixFSym ParticleKinematicFitterModule::getTMatrixFSymMomentumErrorMatrix(ParticleFitObject* fitobject) //fix the warning
428{
429 TMatrixFSym errMatrix(4);
430
431 for (int i = 0; i < 4; i++) {
432 for (int j = i; j < 4; j++) {
433 errMatrix[i][j] = 0.0;
434 }
435 }
436
437 return errMatrix;
438}
439
440// TMatrixFSym ParticleKinematicFitterModule::getTMatrixFSymMomentumVertexErrorMatrix(ParticleFitObject* fitobject) //fix the warning
442{
443 TMatrixFSym errMatrix(7);
444
445 for (int i = 0; i < 7; i++) {
446 for (int j = i; j < 7; j++) {
447 errMatrix[i][j] = 0.0;
448 }
449 }
450
451 return errMatrix;
452}
453
455{
456 if (m_orcaConstraint == "HardBeam") {
457 ROOT::Math::PxPyPzEVector constraints4vector(m_hardConstraintPx.getValue(),
461 return constraints4vector;
462 } else {
463 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
464 }
465
466 // should not reach this point...
467 return ROOT::Math::PxPyPzEVector(0., 0., 0., 0.);
468}
469
471{
472 if (m_orcaConstraint == "HardBeam") {
474 const ROOT::Math::PxPyPzEVector boost = T.getBeamFourMomentum();
475
476 m_hardConstraintPx = MomentumConstraint(0, 1, 0, 0, boost.Px());
477 m_hardConstraintPy = MomentumConstraint(0, 0, 1, 0, boost.Py());
478 m_hardConstraintPz = MomentumConstraint(0, 0, 0, 1, boost.Pz());
479 m_hardConstraintE = MomentumConstraint(1, 0, 0, 0, boost.E());
480
485
486 m_hardConstraintPx.setName("Sum(p_x) [hard]");
487 m_hardConstraintPy.setName("Sum(p_y) [hard]");
488 m_hardConstraintPz.setName("Sum(p_z) [hard]");
489 m_hardConstraintE.setName("Sum(E) [hard]");
490
491 } else if (m_orcaConstraint == "RecoilMass") {
493 const ROOT::Math::PxPyPzEVector boost = T.getBeamFourMomentum();
494
495 m_hardConstraintRecoilMass = RecoilMassConstraint(m_recoilMass, boost.Px(), boost.Py(), boost.Pz(), boost.E());
496
498 m_hardConstraintRecoilMass.setName("Recoil Mass [hard]");
499
500 } else if (m_orcaConstraint == "Mass") {
502
504 m_hardConstraintMass.setName("Mass [hard]");
505 } else {
506 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
507 }
508}
509
511{
512 B2DEBUG(17, "ParticleKinematicFitterModule: Resetting the fitter");
513 fitter.reset();
514}
515
517{
518 if (m_orcaConstraint == "HardBeam") {
523 } else if (m_orcaConstraint == "RecoilMass") {
525 } else if (m_orcaConstraint == "Mass") {
527 } else {
528 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
529 }
530}
531
533{
534 if (m_orcaConstraint == "HardBeam") {
535 fitter.addConstraint(m_hardConstraintPx);
536 fitter.addConstraint(m_hardConstraintPy);
537 fitter.addConstraint(m_hardConstraintPz);
538 fitter.addConstraint(m_hardConstraintE);
539 } else if (m_orcaConstraint == "RecoilMass") {
540 fitter.addConstraint(m_hardConstraintRecoilMass);
541 } else if (m_orcaConstraint == "Mass") {
542 fitter.addConstraint(m_hardConstraintMass);
543 }
544
545 else {
546 B2FATAL("ParticleKinematicFitterModule: " << m_orcaConstraint << " is an invalid OrcaKinFit constraint!");
547 }
548}
549
551{
552 if (m_orcaTracer == "Text") {
553 m_textTracer = new TextTracer(std::cout);
554 fitter.setTracer(m_textTracer);
555 } else if (m_orcaTracer != "None") {
556 B2FATAL("ParticleKinematicFitterModule: " << m_orcaTracer << " is an invalid OrcaKinFit tracer!");
557 }
558}
559
561{
562 B2DEBUG(17, "ParticleKinematicFitterModule::addUnmeasuredGammaToOrcaKinFit: adding an unmeasured photon to the fitter!");
563 // Initialize photon using the existing constraints
564 ROOT::Math::PxPyPzEVector tlv = getLorentzVectorConstraints();
565 double startingE = tlv.E();
566 double startingPhi = tlv.Phi();
567 double startingTheta = tlv.Theta();
568 bool paramFlag = false;
569
570 // create a fit object
571 ParticleFitObject* pfitobject;
572
573 std::string fitObjectName = "UnmeasuredAlongBeam";
574
576 startingTheta = 41.5e-3; // TODO: Read beam crossing angle from db if it's available
577 startingPhi = 0.0;
578 paramFlag = true;
579 } else if (m_fixUnmeasuredPhotonToLER) {
580 startingTheta = TMath::Pi() - 41.5e-3;
581 startingPhi = 0.0;
582 paramFlag = true;
583 } else {
584 fitObjectName = "Unmeasured";
585 }
586
587 // memory allocated: it will be deallocated via "delete fo" in doOrcaKinFitFit
588 pfitobject = new JetFitObject(startingE, startingTheta, startingPhi, 0.0, 0.0, 0.0, 0.);
589 pfitobject->setParam(0, startingE, false, false);
590 pfitobject->setParam(1, startingTheta, paramFlag, paramFlag);
591 pfitobject->setParam(2, startingPhi, paramFlag, paramFlag);
592
593 pfitobject->setName(fitObjectName.c_str());
594 ParticleFitObject& fitobject = *pfitobject;
595
596 // add this fit object (=particle) to the constraints
597 addFitObjectToConstraints(fitobject);
598
599 // add fit particle to the fitter
600 fitter.addFitObject(fitobject);
601}
602
604{
605 std::vector <Belle2::Particle*> bDau = mother->getDaughters();
606 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
607
608 const unsigned nd = bDau.size();
609 unsigned l = 0;
610 std::vector<std::vector<unsigned>> pars;
611 std::vector<Particle*> allparticles;
612 for (unsigned ichild = 0; ichild < nd; ichild++) {
613 const Particle* daughter = mother->getDaughter(ichild);
614 std::vector<unsigned> pard;
615 if (daughter->getNDaughters() > 0) {
616 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
617 pars.push_back(pard);
618 allparticles.push_back(bDau[ichild]);
619 } else {
620 pard.push_back(l);
621 pars.push_back(pard);
622 allparticles.push_back(bDau[ichild]);
623 l++;
624 }
625 }
626
627 if (l == fitObjectContainer->size() - m_addUnmeasuredPhoton) {
628
629 if (fitter.getError() == 0) {
630 for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
631 ROOT::Math::PxPyPzEVector tlv ;
632 TMatrixFSym errMatrixU(7);
633 if (pars[iDaug].size() > 0) {
634 for (unsigned int iChild : pars[iDaug]) {
635 BaseFitObject* fo = fitObjectContainer->at(iChild);
636 auto* fitobject = static_cast<ParticleFitObject*>(fo);
637 ROOT::Math::PxPyPzEVector tlv_sub = getLorentzVector(fitobject);
638 TMatrixFSym errMatrixU_sub = getCovMat7(fitobject);
639 tlv = tlv + tlv_sub;
640 errMatrixU = errMatrixU + errMatrixU_sub;
641 }
642 } else {
643 B2FATAL("ParticleKinematicFitterModule: no fitObject could be used to update the daughter!");
644 }
645 ROOT::Math::XYZVector pos = allparticles[iDaug]->getVertex(); // we don't update the vertex yet
646 TMatrixFSym errMatrix = allparticles[iDaug]->getMomentumVertexErrorMatrix();
647 TMatrixFSym errMatrixMom = allparticles[iDaug]->getMomentumErrorMatrix();
648 TMatrixFSym errMatrixVer = allparticles[iDaug]->getVertexErrorMatrix();
649
650 for (int i = 0; i < 3; i++) {
651 for (int j = i; j < 3; j++) {
652 errMatrixU[i + 4][j + 4] = errMatrixVer[i][j];
653 }
654 }
655
656 allparticles[iDaug]->set4Vector(tlv);
657 allparticles[iDaug]->setVertex(pos);
658 allparticles[iDaug]->setMomentumVertexErrorMatrix(errMatrixU);
659 }
660 }
661
662 return true;
663 } else {
664 B2ERROR("updateOrcaKinFitDaughters: Cannot update daughters, mismatch between number of daughters and number of fitobjects!");
665 return false;
666 }
667}
668
669void ParticleKinematicFitterModule::updateMapOfTrackAndDaughter(unsigned& l, std::vector<std::vector<unsigned>>& pars,
670 std::vector<unsigned>& parm, std::vector<Particle*>& allparticles, const Particle* daughter)
671{
672 std::vector <Belle2::Particle*> dDau = daughter->getDaughters();
673 for (unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
674 const Particle* child = daughter->getDaughter(ichild);
675 std::vector<unsigned> pard;
676 if (child->getNDaughters() > 0) {
677 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
678 parm.insert(parm.end(), pard.begin(), pard.end());
679 pars.push_back(pard);
680 allparticles.push_back(dDau[ichild]);
681 } else {
682 pard.push_back(l);
683 parm.push_back(l);
684 pars.push_back(pard);
685 allparticles.push_back(dDau[ichild]);
686 l++;
687 }
688 }
689}
690
691void ParticleKinematicFitterModule::updateOrcaKinFitMother(BaseFitter& fitter, std::vector<Particle*>& particleChildren,
692 Particle* mother)
693{
694 // get old values
695 ROOT::Math::XYZVector pos = mother->getVertex();
696 TMatrixFSym errMatrix = mother->getMomentumVertexErrorMatrix();
697 float pvalue = mother->getPValue();
698
699 // update momentum vector
700 ROOT::Math::PxPyPzEVector momnew(0., 0., 0., 0.);
701
702 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
703 for (unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
704 BaseFitObject* fo = fitObjectContainer->at(iChild);
705 auto* fitobject = static_cast<ParticleFitObject*>(fo);
706 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
707 momnew += tlv;
708 }
709
710 // update
711 // TODO: use pvalue of the fit or the old one of the mother? use fit covariance matrix?
712 // Maybe here should use the pvalue and errmatrix of the fit ----Yu Hu
713 mother->updateMomentum(momnew, pos, errMatrix, pvalue);
714}
715
717 std::vector<Particle*>& particleChildren, Particle* mother)
718{
719 bool updated = false;
720 std::vector <BaseFitObject*>* fitObjectContainer = fitter.getFitObjects();
721
722 for (unsigned iChild = 0; iChild < particleChildren.size(); iChild++) {
723 BaseFitObject* fo = fitObjectContainer->at(iChild);
724 auto* fitobject = static_cast<ParticleFitObject*>(fo);
725 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
726
727 // name of extra variables
728 std::string extraVariableParticlePx = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_Px";
729 std::string extraVariableParticlePy = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_Py";
730 std::string extraVariableParticlePz = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_Pz";
731 std::string extraVariableParticleE = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_E";
732 std::string extraVariableParticlePxErr = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_PxErr";
733 std::string extraVariableParticlePyErr = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_PyErr";
734 std::string extraVariableParticlePzErr = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_PzErr";
735 std::string extraVariableParticleEErr = "OrcaKinFit" + prefix + "_" + SSTR(iChild) + "_EErr";
736
737 mother->addExtraInfo(extraVariableParticlePx, tlv.Px());
738 mother->addExtraInfo(extraVariableParticlePy, tlv.Py());
739 mother->addExtraInfo(extraVariableParticlePz, tlv.Pz());
740 mother->addExtraInfo(extraVariableParticleE, tlv.E());
741 mother->addExtraInfo(extraVariableParticlePxErr, getFitObjectError(fitobject, 0));
742 mother->addExtraInfo(extraVariableParticlePyErr, getFitObjectError(fitobject, 1));
743 mother->addExtraInfo(extraVariableParticlePzErr, getFitObjectError(fitobject, 2));
744 mother->addExtraInfo(extraVariableParticleEErr, -1.0);
745
746 }
747
748 return updated;
749}
750
752{
753 //check if it is a PxPyPzMFitObject
754 auto* pxpypzmfitobject = static_cast<PxPyPzMFitObject*>(fitobject);
755 if (pxpypzmfitobject) {
756 return fitobject->getError(ilocal);
757 } else {
758 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
759 }
760}
761
762// this returns the local covariance matrix
764{
765
766 //check if it is a PxPyPzMFitObject
767 auto* pxpypzmfitobject = static_cast<PxPyPzMFitObject*>(fitobject);
768 if (pxpypzmfitobject) {
769
770 TMatrixFSym errMatrix(3);
771
772 //loop over the i-j local variables.
773 for (int i = 0; i < 3; i++) {
774 for (int j = 0; j < 3; j++) {
775 errMatrix[i][j] = pxpypzmfitobject->getCov(i, j);
776 }
777 }
778
779 return errMatrix;
780 } else {
781 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
782 }
783}
784
786{
787 TMatrixFSym fitCovMatrix(3);
788
789 if (strcmp(fitobject->getParamName(0), "E") == 0) {
790 //check if it is a JetFitObject
791 auto* jetfitObject = static_cast<JetFitObject*>(fitobject);
792 if (jetfitObject) {
793
794 fitCovMatrix = getFitObjectCovMat(fitobject);
795 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
796
797 const double energy = tlv.E();
798 const double theta = tlv.Theta();
799 const double phi = tlv.Phi();
800
801 const double st = sin(theta);
802 const double ct = cos(theta);
803 const double sp = sin(phi);
804 const double cp = cos(phi);
805
806 // updated covariance matrix is: A * cov * A^T where A is the Jacobi matrix (use Similarity)
807 TMatrixF A(7, 3);
808 A(0, 0) = cp * st ; // dpx/dE
809 A(0, 1) = energy * cp * ct ; // dpx/dtheta
810 A(0, 2) = -energy * sp * st ; // dpx/dphi
811 A(1, 0) = sp * st ; // dpy/dE
812 A(1, 1) = energy * sp * ct ; // dpz/dtheta
813 A(1, 2) = energy * cp * st ; // dpy/dphi
814 A(2, 0) = ct ; // dpz/dE
815 A(2, 1) = -energy * st ; // dpz/dtheta
816 A(2, 2) = 0 ; // dpz/dphi
817 A(3, 0) = 1.0; // dE/dE
818 A(3, 1) = 0.0; // dE/dphi
819 A(3, 2) = 0.0; // dE/dtheta
820
821 TMatrixFSym D = fitCovMatrix.Similarity(A);
822 return D;
823
824 } else {
825 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
826 }
827 } else {
828 //check if it is a PxPyPzMFitObject
829 auto* pxpypzmfitobject = static_cast<PxPyPzMFitObject*>(fitobject);
830 if (pxpypzmfitobject) {
831
832 fitCovMatrix = getFitObjectCovMat(fitobject);
833
834 // updated covariance matrix is: A * cov * A^T where A is the Jacobi matrix (use Similarity)
835 ROOT::Math::PxPyPzEVector tlv = getLorentzVector(fitobject);
836 TMatrixF A(7, 3);
837 A[0][0] = 1.; // px/dpx
838 A[0][1] = 0.; // px/dpy
839 A[0][2] = 0.; // px/dpz
840 A[1][0] = 0.; // py/dpx
841 A[1][1] = 1.; // py/dpy
842 A[1][2] = 0.; // py/dpz
843 A[2][0] = 0.; // pz/dpx
844 A[2][1] = 0.; // pz/dpy
845 A[2][2] = 1.; // pz/dpz
846 if (tlv.E() > 0.0) {
847 A[3][0] = tlv.Px() / tlv.E(); // E/dpx, E=sqrt(px^2 + py^2 + pz^2 + m^2)
848 A[3][1] = tlv.Py() / tlv.E(); // E/dpy
849 A[3][2] = tlv.Pz() / tlv.E(); // E/dpz
850 }
851
852 TMatrixFSym D = fitCovMatrix.Similarity(A);
853
854 return D;
855 } else {
856 B2FATAL("ParticleKinematicFitterModule: not implemented yet");
857 }
858 }
859}
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
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.
CLHEP::HepSymMatrix getCLHEPMomentumErrorMatrix(Particle *particle)
Returns particle's 4x4 momentum-error matrix as a HepSymMatrix.
double m_invMass
Invariant mass for Mass constraint.
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 storeOrcaKinFitParticles(const std::string &prefix, BaseFitter &fitter, std::vector< Particle * > &particleChildren, Particle *mother)
store fit object information as ExtraInfo
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 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
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
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
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
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.
STL namespace.