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