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