Belle II Software  release-06-00-14
ParticleVertexFitterModule.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 
9 // Own include
10 #include <analysis/modules/ParticleVertexFitter/ParticleVertexFitterModule.h>
11 
12 // framework aux
13 #include <framework/gearbox/Unit.h>
14 #include <framework/gearbox/Const.h>
15 #include <framework/logging/Logger.h>
16 #include <framework/particledb/EvtGenDatabasePDG.h>
17 
18 // dataobjects
19 #include <analysis/dataobjects/Particle.h>
20 #include <analysis/dataobjects/Btube.h>
21 // utilities
22 #include <analysis/utility/CLHEPToROOT.h>
23 #include <analysis/utility/PCmsLabTransform.h>
24 #include <analysis/utility/ParticleCopy.h>
25 #include <analysis/utility/ROOTToCLHEP.h>
26 
27 // Magnetic field
28 #include <framework/geometry/BFieldManager.h>
29 
30 #include <TMath.h>
31 
32 using namespace std;
33 
34 namespace Belle2 {
41  //-----------------------------------------------------------------
42  // Register module
43  //-----------------------------------------------------------------
44 
45  REG_MODULE(ParticleVertexFitter)
46 
47  //-----------------------------------------------------------------
48  // Implementation
49  //-----------------------------------------------------------------
50 
52  m_Bfield(0)
53  {
54  // set module description (e.g. insert text)
55  setDescription("Vertex fitter for modular analysis");
56  setPropertyFlags(c_ParallelProcessingCertified);
57 
58  // Add parameters
59  addParam("listName", m_listName, "name of particle list", string(""));
60  addParam("confidenceLevel", m_confidenceLevel,
61  "Confidence level to accept the fit. Particle candidates with "
62  "p-value less than confidenceLevel are removed from the particle "
63  "list. If set to -1, all candidates are kept; if set to 0, "
64  "the candidates failing the fit are removed.",
65  0.001);
66  addParam("vertexFitter", m_vertexFitter, "KFit or Rave", string("KFit"));
67  addParam("fitType", m_fitType, "type of the kinematic fit (vertex, massvertex, mass)", string("vertex"));
68  addParam("withConstraint", m_withConstraint,
69  "additional constraint on vertex: ipprofile, iptube, mother, iptubecut, pointing, btube",
70  string(""));
71  addParam("decayString", m_decayString, "specifies which daughter particles are included in the kinematic fit", string(""));
72  addParam("updateDaughters", m_updateDaughters, "true: update the daughters after the vertex fit", false);
73  addParam("smearing", m_smearing, "smear IP tube width by given length", 0.002);
74  addParam("massConstraintList", m_massConstraintList,
75  "Type::[int]. List of daughter particles to mass constrain with int = pdg code. (only for MassFourCKFit)", {});
76  addParam("massConstraintListParticlename", m_massConstraintListParticlename,
77  "Type::[string]. List of daughter particles to mass constrain with string = particle name. (only for MassFourCKFit)", {});
78  }
79 
80  void ParticleVertexFitterModule::initialize()
81  {
82  // Particle list with name m_listName has to exist
83  m_plist.isRequired(m_listName);
84 
85  // magnetic field
86  m_Bfield = BFieldManager::getField(TVector3(0, 0, 0)).Z() / Unit::T;
87 
88  // RAVE setup
89  if (m_vertexFitter == "Rave")
90  analysis::RaveSetup::initialize(1, m_Bfield);
91 
92  B2DEBUG(1, "ParticleVertexFitterModule : magnetic field = " << m_Bfield);
93 
94 
95  if (m_decayString != "")
96  m_decaydescriptor.init(m_decayString);
97 
98  if ((m_massConstraintList.size()) == 0 && (m_massConstraintListParticlename.size()) > 0) {
99  for (auto& containedParticle : m_massConstraintListParticlename) {
100  TParticlePDG* particletemp = TDatabasePDG::Instance()->GetParticle((containedParticle).c_str());
101  m_massConstraintList.push_back(particletemp->PdgCode());
102  }
103  }
104 
105  B2INFO("ParticleVertexFitter: Performing " << m_fitType << " fit on " << m_listName << " using " << m_vertexFitter);
106  if (m_decayString != "")
107  B2INFO("ParticleVertexFitter: Using specified decay string: " << m_decayString);
108  if (m_withConstraint != "")
109  B2INFO("ParticleVertexFitter: Additional " << m_withConstraint << " will be applied");
110 
111  }
112 
113  void ParticleVertexFitterModule::beginRun()
114  {
115  //TODO: set magnetic field for each run
116  //m_Bfield = BFieldMap::Instance().getBField(TVector3(0,0,0)).Z();
117  //TODO: set IP spot size for each run
118  }
119 
120  void ParticleVertexFitterModule::event()
121  {
122  if (m_vertexFitter == "Rave")
123  analysis::RaveSetup::initialize(1, m_Bfield);
124 
125  m_BeamSpotCenter = m_beamSpotDB->getIPPosition();
126  m_beamSpotCov.ResizeTo(3, 3);
127  TMatrixDSym beamSpotCov(3);
128  if (m_withConstraint == "ipprofile") m_beamSpotCov = m_beamSpotDB->getCovVertex();
129  if (m_withConstraint == "iptube") {
130  if (m_smearing > 0 && m_vertexFitter == "KFit") {
131  ParticleVertexFitterModule::smearBeamSpot(m_smearing);
132  } else {
133  ParticleVertexFitterModule::findConstraintBoost(2.);
134  }
135  }
136  if (m_withConstraint == "iptubecut") { // for development purpose only
137  m_BeamSpotCenter = TVector3(0.001, 0., .013);
138  findConstraintBoost(0.03);
139  }
140  if ((m_vertexFitter == "Rave") && (m_withConstraint == "ipprofile" || m_withConstraint == "iptube"
141  || m_withConstraint == "mother" || m_withConstraint == "iptubecut" || m_withConstraint == "btube"))
142  analysis::RaveSetup::getInstance()->setBeamSpot(m_BeamSpotCenter, m_beamSpotCov);
143  std::vector<unsigned int> toRemove;
144  unsigned int n = m_plist->getListSize();
145  for (unsigned i = 0; i < n; i++) {
146  Particle* particle = m_plist->getParticle(i);
147  m_hasCovMatrix = false;
148  if (m_updateDaughters == true) {
149  if (m_decayString.empty()) ParticleCopy::copyDaughters(particle);
150  else B2ERROR("Daughters update works only when all daughters are selected. Daughters will not be updated");
151  }
152 
153  if (m_withConstraint == "mother") {
154  m_BeamSpotCenter = particle->getVertex();
155  m_beamSpotCov = particle->getVertexErrorMatrix();
156  }
157 
158  TMatrixFSym mother_errMatrix(7);
159  mother_errMatrix = particle->getMomentumVertexErrorMatrix();
160  for (int k = 0; k < 7; k++) {
161  for (int j = 0; j < 7; j++) {
162  if (mother_errMatrix[k][j] > 0) {
163  m_hasCovMatrix = true;
164  }
165  }
166  }
167  bool hasTube = true;
168  if (m_withConstraint == "btube") {
169  Btube* Ver = particle->getRelatedTo<Btube>();
170  if (!Ver) {
171  hasTube = false;
172  toRemove.push_back(particle->getArrayIndex());
173  } else {
174  m_BeamSpotCenter.SetXYZ(Ver->getTubeCenter()(0, 0), Ver->getTubeCenter()(1, 0), Ver->getTubeCenter()(2, 0));
175  m_beamSpotCov = Ver->getTubeMatrix();
176  }
177  }
178  bool ok = false;
179  if (hasTube) {
180  ok = doVertexFit(particle);
181  }
182  if (!ok)
183  particle->setPValue(-1);
184  if (particle->getPValue() < m_confidenceLevel)
185  toRemove.push_back(particle->getArrayIndex());
186  }
187  m_plist->removeParticles(toRemove);
188 
189  //free memory allocated by rave. initialize() would be enough, except that we must clean things up before program end...
190  if (m_vertexFitter == "Rave")
191  analysis::RaveSetup::getInstance()->reset();
192  }
193 
194  bool ParticleVertexFitterModule::doVertexFit(Particle* mother)
195  {
196  // steering starts here
197 
198  if (m_Bfield == 0) {
199  B2FATAL("ParticleVertexFitter: No magnetic field");
200  }
201 
202  if (m_withConstraint != "ipprofile" &&
203  m_withConstraint != "iptube" &&
204  m_withConstraint != "mother" &&
205  m_withConstraint != "iptubecut" &&
206  m_withConstraint != "pointing" &&
207  m_withConstraint != "btube" &&
208  m_withConstraint != "")
209  B2FATAL("ParticleVertexFitter: " << m_withConstraint << " ***invalid Constraint ");
210 
211  bool ok = false;
212  // fits with KFit
213  if (m_vertexFitter == "KFit") {
214 
215  // vertex fit
216  if (m_fitType == "vertex") {
217  if (m_withConstraint == "ipprofile") {
218  ok = doKVertexFit(mother, true, false);
219  } else if (m_withConstraint == "iptube") {
220  ok = doKVertexFit(mother, false, true);
221  } else {
222  ok = doKVertexFit(mother, false, false);
223  }
224  }
225 
226  // mass-constrained vertex fit
227  if (m_fitType == "massvertex") {
228  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
229  B2FATAL("ParticleVertexFitter: Invalid options - mass-constrained fit using KFit does not work with iptube or ipprofile constraint.");
230  } else if (m_withConstraint == "pointing") {
231  ok = doKMassPointingVertexFit(mother);
232  } else {
233  ok = doKMassVertexFit(mother);
234  }
235  }
236 
237  // mass fit
238  if (m_fitType == "mass") {
239  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
240  B2FATAL("ParticleVertexFitter: Invalid options - mass fit using KFit does not work with iptube or ipprofile constraint.");
241  } else {
242  ok = doKMassFit(mother);
243  }
244  }
245 
246  // four C fit
247  if (m_fitType == "fourC") {
248  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
249  B2FATAL("ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
250  } else {
251  ok = doKFourCFit(mother);
252  }
253  }
254 
255  // four mass C fit
256  if (m_fitType == "massfourC") {
257  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
258  B2FATAL("ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
259  } else {
260  ok = doKMassFourCFit(mother);
261  }
262  }
263 
264  // invalid KFit fit type
265  if (m_fitType != "vertex"
266  && m_fitType != "massvertex"
267  && m_fitType != "mass"
268  && m_fitType != "fourC"
269  && m_fitType != "massfourC")
270  B2FATAL("ParticleVertexFitter: " << m_fitType << " ***invalid fit type for the vertex fitter ");
271  }
272 
273  // fits using Rave
274  if (m_vertexFitter == "Rave") {
275  try {
276  ok = doRaveFit(mother);
277  } catch (const rave::CheckedFloatException&) {
278  B2ERROR("Invalid inputs (nan/inf)?");
279  ok = false;
280  }
281  }
282 
283  // invalid fitter
284  if (m_vertexFitter != "KFit" && m_vertexFitter != "Rave")
285  B2FATAL("ParticleVertexFitter: " << m_vertexFitter << " ***invalid vertex fitter ");
286 
287  if (!ok) return false;
288 
289  // steering ends here
290 
291  //if (mother->getPValue() < m_confidenceLevel) return false;
292  return true;
293 
294  }
295 
296  bool ParticleVertexFitterModule::fillFitParticles(const Particle* mother, std::vector<const Particle*>& fitChildren,
297  std::vector<const Particle*>& twoPhotonChildren)
298  {
299  if (m_decayString.empty()) {
300  // if decayString is empty, just use all primary daughters
301  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
302  const Particle* child = mother->getDaughter(ichild);
303  // This if allows to skip the daughters, which cannot be used in the fits, particularly K_L0 from KLM.
304  // Useful for fully-inclusive particles.
305  if (mother->getProperty() == Particle::PropertyFlags::c_IsUnspecified and child->getPValue() < 0) {
306  continue;
307  }
308  fitChildren.push_back(child);
309  }
310  } else {
311  fitChildren = m_decaydescriptor.getSelectionParticles(mother);
312  }
313 
314  auto itr = fitChildren.begin();
315  while (itr != fitChildren.end()) {
316  const Particle* child = *itr;
317 
318  if (child->getPValue() < 0) {
319  B2WARNING("Daughter with PDG code " << child->getPDGCode() << " does not have a valid error matrix.");
320  return false; // error matrix not valid
321  }
322  bool isTwoPhotonParticle = false;
323  if (m_hasCovMatrix == false) {
324  if (child->getPDGCode() == Const::pi0.getPDGCode() or child->getPDGCode() == 221) { // pi0 or eta
325  if (child->getNDaughters() == 2) {
326  if (child->getDaughter(0)->getPDGCode() == Const::photon.getPDGCode()
327  && child->getDaughter(1)->getPDGCode() == Const::photon.getPDGCode()) {
328  isTwoPhotonParticle = true;
329  }
330  }
331  }
332  }
333  if (isTwoPhotonParticle) {
334  // move children from fitChildren to twoPhotonChildren
335  twoPhotonChildren.push_back(child);
336  itr = fitChildren.erase(itr);
337  } else {
338  itr++;
339  }
340  }
341 
342  return true;
343  }
344 
345  bool ParticleVertexFitterModule::redoTwoPhotonDaughterMassFit(Particle* postFit, const Particle* preFit,
346  const analysis::VertexFitKFit& kv)
347  {
348  // TODO: something like setGammaError is necessary
349  // this is just workaround for the moment
350 
351  const Particle* g1Orig = preFit->getDaughter(0);
352  const Particle* g2Orig = preFit->getDaughter(1);
353  Particle g1Temp(g1Orig->get4Vector(), 22);
354  Particle g2Temp(g2Orig->get4Vector(), 22);
355 
356  TMatrixFSym g1ErrMatrix = g1Orig->getMomentumVertexErrorMatrix();
357  TMatrixFSym g2ErrMatrix = g2Orig->getMomentumVertexErrorMatrix();
358 
359  TVector3 pos(kv.getVertex().x(), kv.getVertex().y(), kv.getVertex().z());
360  CLHEP::HepSymMatrix posErrorMatrix = kv.getVertexError();
361 
362  TMatrixFSym errMatrix(3);
363  for (int i = 0; i < 3; i++)
364  for (int j = 0; j < 3; j++)
365  errMatrix(i, j) = posErrorMatrix[i][j];
366 
367  g1ErrMatrix.SetSub(4, errMatrix);
368  g2ErrMatrix.SetSub(4, errMatrix);
369 
370  g1Temp.updateMomentum(g1Orig->get4Vector(), pos, g1ErrMatrix, 1.0);
371  g2Temp.updateMomentum(g2Orig->get4Vector(), pos, g2ErrMatrix, 1.0);
372 
373  // perform the mass fit for the two-photon particle
375  km.setMagneticField(m_Bfield);
376 
377  km.addParticle(&g1Temp);
378  km.addParticle(&g2Temp);
379 
380  km.setVertex(kv.getVertex());
382  km.setInvariantMass(preFit->getPDGMass());
383 
384  int err = km.doFit();
385  if (err != 0) {
386  return false;
387  }
388 
389  // The update of the daughters is disabled for this mass fit.
390  bool updateDaughters = m_updateDaughters;
391  m_updateDaughters = false;
392  bool ok = makeKMassMother(km, postFit);
393  m_updateDaughters = updateDaughters;
394 
395  return ok;
396  }
397 
398  bool ParticleVertexFitterModule::doKVertexFit(Particle* mother, bool ipProfileConstraint, bool ipTubeConstraint)
399  {
400  if ((mother->getNDaughters() < 2 && !ipTubeConstraint) || mother->getNDaughters() < 1) return false;
401 
402  std::vector<const Particle*> fitChildren;
403  std::vector<const Particle*> twoPhotonChildren;
404  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
405 
406  if (!validChildren)
407  return false;
408 
409  if (twoPhotonChildren.size() > 1) {
410  B2FATAL("[ParticleVertexFitterModule::doKVertexFit] Vertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
411  }
412 
413  if ((fitChildren.size() < 2 && !ipTubeConstraint) || fitChildren.size() < 1) {
414  B2WARNING("[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is too low.");
415  return false;
416  }
417 
418  // Initialise the Fitter
420  kv.setMagneticField(m_Bfield);
421 
422  for (auto& child : fitChildren)
423  kv.addParticle(child);
424 
425  if (ipProfileConstraint)
426  addIPProfileToKFit(kv);
427 
428  if (ipTubeConstraint)
429  addIPTubeToKFit(kv);
430 
431  // Perform vertex fit using only the particles with valid error matrices
432  int err = kv.doFit();
433  if (err != 0)
434  return false;
435 
436  bool ok = false;
437  if (twoPhotonChildren.size() == 0)
438  // in the case daughters do not include pi0 - this is it (fit done)
439  ok = makeKVertexMother(kv, mother);
440  else if (twoPhotonChildren.size() == 1) {
441  // there is a daughter reconstructed from two photons so without position information
442  // 1. determine vertex based on all other valid daughters
443  // 2. set position and error matrix of two-photon daughter to previously determined vertex
444  // 3. redo the fit using all particles (including two-photon particle this time)
445 
446  const Particle* twoPhotonDaughter = twoPhotonChildren[0];
447  Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
448  ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
449  if (!ok)
450  return false;
451 
452  // finally perform the fit using all daughter particles
454  kv2.setMagneticField(m_Bfield);
455 
456  for (auto& child : fitChildren)
457  kv2.addParticle(child);
458 
459  kv2.addParticle(&fixedTwoPhotonDaughter);
460 
461  if (ipProfileConstraint)
462  addIPProfileToKFit(kv2);
463 
464  err = kv2.doFit();
465 
466  if (err != 0)
467  return false;
468 
469  ok = makeKVertexMother(kv2, mother);
470  }
471 
472  return ok;
473  }
474 
475  bool ParticleVertexFitterModule::doKMassVertexFit(Particle* mother)
476  {
477  if (mother->getNDaughters() < 2) return false;
478 
479  std::vector<const Particle*> fitChildren;
480  std::vector<const Particle*> twoPhotonChildren;
481  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
482 
483  if (!validChildren)
484  return false;
485 
486  if (twoPhotonChildren.size() > 1) {
487  B2FATAL("[ParticleVertexFitterModule::doKVertexFit] MassVertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
488  }
489 
490  if (fitChildren.size() < 2) {
491  B2WARNING("[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
492  return false;
493  }
494 
495  bool ok = false;
496  if (twoPhotonChildren.size() == 0) {
497  // Initialise the Fitter
499  kmv.setMagneticField(m_Bfield);
500 
501  for (auto child : fitChildren)
502  kmv.addParticle(child);
503 
504  kmv.setInvariantMass(mother->getPDGMass());
505  int err = kmv.doFit();
506  if (err != 0)
507  return false;
508 
509  // in the case daughters do not include particles with two photon daughters like pi0 - this is it (fit done)
510  ok = makeKMassVertexMother(kmv, mother);
511  } else if (twoPhotonChildren.size() == 1) {
512  // there is a daughter reconstructed from two photons so without position information
513  // 1. determine vertex based on all other valid daughters
514  // 2. set position and error matrix of two-photon daughter to previously determined vertex
515  // 3. redo the fit using all particles (including two-photon particle this time)
516 
518  kv.setMagneticField(m_Bfield);
519 
520  for (auto child : fitChildren)
521  kv.addParticle(child);
522 
523  // Perform vertex fit using only the particles with valid error matrices
524  int err = kv.doFit();
525  if (err != 0)
526  return false;
527 
528  const Particle* twoPhotonDaughter = twoPhotonChildren[0];
529  Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
530  ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
531  if (!ok)
532  return false;
533 
534  // finally perform the fit using all daughter particles
536  kmv2.setMagneticField(m_Bfield);
537 
538  for (auto child : fitChildren)
539  kmv2.addParticle(child);
540  kmv2.addParticle(&fixedTwoPhotonDaughter);
541 
542  kmv2.setInvariantMass(mother->getPDGMass());
543  err = kmv2.doFit();
544 
545  if (err != 0)
546  return false;
547 
548  ok = makeKMassVertexMother(kmv2, mother);
549  }
550 
551  return ok;
552 
553  }
554 
555  bool ParticleVertexFitterModule::doKMassPointingVertexFit(Particle* mother)
556  {
557  if (!(mother->hasExtraInfo("prodVertX") && mother->hasExtraInfo("prodVertY") && mother->hasExtraInfo("prodVertZ"))) {
558  return false;
559  }
560 
561  if (mother->getNDaughters() < 2) return false;
562 
563  std::vector<const Particle*> fitChildren;
564  std::vector<const Particle*> twoPhotonChildren;
565  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
566 
567  if (!validChildren)
568  return false;
569 
570  if (twoPhotonChildren.size() > 0) {
571  B2FATAL("[ParticleVertexFitterModule::doKMassPointingVertexFit] MassPointingVertex fit using KFit does not support fit with two-photon daughters (yet).");
572  }
573 
574  if (fitChildren.size() < 2) {
575  B2WARNING("[ParticleVertexFitterModule::doKMassPointingVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
576  return false;
577  }
578 
579  bool ok = false;
580  // Initialise the Fitter
582  kmpv.setMagneticField(m_Bfield);
583 
584  for (auto child : fitChildren)
585  kmpv.addParticle(child);
586 
587  kmpv.setInvariantMass(mother->getPDGMass());
588  HepPoint3D productionVertex(mother->getExtraInfo("prodVertX"),
589  mother->getExtraInfo("prodVertY"),
590  mother->getExtraInfo("prodVertZ"));
591  kmpv.setProductionVertex(productionVertex);
592  int err = kmpv.doFit();
593  if (err != 0) return false;
594 
595  ok = makeKMassPointingVertexMother(kmpv, mother);
596 
597  return ok;
598  }
599 
600  bool ParticleVertexFitterModule::doKMassFit(Particle* mother)
601  {
602  if (mother->getNDaughters() < 2) return false;
603 
605  km.setMagneticField(m_Bfield);
606 
607  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
608  const Particle* child = mother->getDaughter(ichild);
609 
610  if (child->getPValue() < 0) return false; // error matrix not valid
611 
612  km.addParticle(child);
613  }
614 
615  // apply mass constraint
616  km.setInvariantMass(mother->getPDGMass());
617 
618  int err = km.doFit();
619 
620  if (err != 0) return false;
621 
622  bool ok = makeKMassMother(km, mother);
623 
624  return ok;
625  }
626 
627  bool ParticleVertexFitterModule::doKFourCFit(Particle* mother)
628  {
629  if (mother->getNDaughters() < 2) return false;
630 
632  kf.setMagneticField(m_Bfield);
633 
634  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
635  const Particle* child = mother->getDaughter(ichild);
636 
637  if (child->getNDaughters() > 0) {
638  bool err = addChildofParticletoKFit(kf, child);
639  if (!err) return false;
640  } else {
641  if (child->getPValue() < 0) return false; // error matrix not valid
642 
643  kf.addParticle(child);
644  }
645  }
646 
647  // apply four momentum constraint
650 
651  int err = kf.doFit();
652 
653  if (err != 0) return false;
654 
655  bool ok = makeKFourCMother(kf, mother);
656 
657  return ok;
658  }
659 
660  bool ParticleVertexFitterModule::doKMassFourCFit(Particle* mother)
661  {
662  if (mother->getNDaughters() < 2) return false;
663 
665  kf.setMagneticField(m_Bfield);
666 
667  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
668  const Particle* child = mother->getDaughter(ichild);
669 
670  if (child->getNDaughters() > 0) {
671  bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
672  std::abs(child->getPDGCode())) != m_massConstraintList.end();
673  std::vector<unsigned> childId;
674  bool err = addChildofParticletoMassKFit(kf, child, childId);
675  if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
676  if (!err) return false;
677  } else {
678  if (child->getPValue() < 0) return false; // error matrix not valid
679  kf.addParticle(child);
680  }
681  }
682 
683  // apply four momentum constraint
686 
687  int err = kf.doFit();
688 
689  if (err != 0) return false;
690 
691  bool ok = makeMassKFourCMother(kf, mother);
692 
693  return ok;
694  }
695 
696  bool ParticleVertexFitterModule::makeKVertexMother(analysis::VertexFitKFit& kv,
697  Particle* mother)
698  {
699  enum analysis::KFitError::ECode fitError;
700  fitError = kv.updateMother(mother);
701  if (fitError != analysis::KFitError::kNoError)
702  return false;
703  if (m_decayString.empty() && m_updateDaughters == true) {
704  // update daughter momenta as well
705  // the order of daughters in the *fitter is the same as in the mother Particle
706 
707  std::vector<Particle*> daughters = mother->getDaughters();
708 
709  unsigned track_count = kv.getTrackCount();
710  if (daughters.size() != track_count)
711  return false;
712 
713  for (unsigned iChild = 0; iChild < track_count; iChild++) {
714  daughters[iChild]->set4Vector(
715  CLHEPToROOT::getTLorentzVector(kv.getTrackMomentum(iChild)));
716  daughters[iChild]->setVertex(
717  CLHEPToROOT::getTVector3(kv.getTrackPosition(iChild)));
718  daughters[iChild]->setMomentumVertexErrorMatrix(
719  CLHEPToROOT::getTMatrixFSym(kv.getTrackError(iChild)));
720  }
721  }
722 
723  return true;
724  }
725 
726  bool ParticleVertexFitterModule::makeKMassVertexMother(analysis::MassVertexFitKFit& kmv,
727  Particle* mother)
728  {
729  enum analysis::KFitError::ECode fitError;
730  fitError = kmv.updateMother(mother);
731  if (fitError != analysis::KFitError::kNoError)
732  return false;
733  if (m_decayString.empty() && m_updateDaughters == true) {
734  // update daughter momenta as well
735  // the order of daughters in the *fitter is the same as in the mother Particle
736 
737  std::vector<Particle*> daughters = mother->getDaughters();
738 
739  unsigned track_count = kmv.getTrackCount();
740  if (daughters.size() != track_count)
741  return false;
742 
743  for (unsigned iChild = 0; iChild < track_count; iChild++) {
744  daughters[iChild]->set4Vector(
745  CLHEPToROOT::getTLorentzVector(kmv.getTrackMomentum(iChild)));
746  daughters[iChild]->setVertex(
747  CLHEPToROOT::getTVector3(kmv.getTrackPosition(iChild)));
748  daughters[iChild]->setMomentumVertexErrorMatrix(
749  CLHEPToROOT::getTMatrixFSym(kmv.getTrackError(iChild)));
750  }
751  }
752 
753  return true;
754  }
755 
756  bool ParticleVertexFitterModule::makeKMassPointingVertexMother(analysis::MassPointingVertexFitKFit& kmpv,
757  Particle* mother)
758  {
759  enum analysis::KFitError::ECode fitError;
760  fitError = kmpv.updateMother(mother);
761  if (fitError != analysis::KFitError::kNoError) {
762  return false;
763  }
764 
765  if (m_decayString.empty() && m_updateDaughters == true) {
766  // update daughter momenta as well
767  // the order of daughters in the *fitter is the same as in the mother Particle
768 
769  std::vector<Particle*> daughters = mother->getDaughters();
770 
771  unsigned track_count = kmpv.getTrackCount();
772  if (daughters.size() != track_count)
773  return false;
774 
775  for (unsigned iChild = 0; iChild < track_count; iChild++) {
776  daughters[iChild]->set4Vector(
777  CLHEPToROOT::getTLorentzVector(kmpv.getTrackMomentum(iChild)));
778  daughters[iChild]->setVertex(
779  CLHEPToROOT::getTVector3(kmpv.getTrackPosition(iChild)));
780  daughters[iChild]->setMomentumVertexErrorMatrix(
781  CLHEPToROOT::getTMatrixFSym(kmpv.getTrackError(iChild)));
782  }
783  }
784 
785  return true;
786  }
787 
788 
789  bool ParticleVertexFitterModule::makeKMassMother(analysis::MassFitKFit& km,
790  Particle* mother)
791  {
792  enum analysis::KFitError::ECode fitError;
793  fitError = km.updateMother(mother);
794  if (fitError != analysis::KFitError::kNoError)
795  return false;
796  if (m_decayString.empty() && m_updateDaughters == true) {
797  // update daughter momenta as well
798  // the order of daughters in the *fitter is the same as in the mother Particle
799 
800  std::vector<Particle*> daughters = mother->getDaughters();
801 
802  unsigned track_count = km.getTrackCount();
803  if (daughters.size() != track_count)
804  return false;
805 
806  for (unsigned iChild = 0; iChild < track_count; iChild++) {
807  daughters[iChild]->set4Vector(
808  CLHEPToROOT::getTLorentzVector(km.getTrackMomentum(iChild)));
809  daughters[iChild]->setVertex(
810  CLHEPToROOT::getTVector3(km.getTrackPosition(iChild)));
811  daughters[iChild]->setMomentumVertexErrorMatrix(
812  CLHEPToROOT::getTMatrixFSym(km.getTrackError(iChild)));
813  }
814  }
815 
816  return true;
817  }
818 
819 
820 
821  bool ParticleVertexFitterModule::makeKFourCMother(analysis::FourCFitKFit& kf, Particle* mother)
822  {
823  enum analysis::KFitError::ECode fitError;
824  fitError = kf.updateMother(mother);
825  if (fitError != analysis::KFitError::kNoError)
826  return false;
827  mother->addExtraInfo("FourCFitProb", kf.getCHIsq());
828  mother->addExtraInfo("FourCFitChi2", kf.getNDF());
829  if (m_decayString.empty() && m_updateDaughters == true) {
830  // update daughter momenta as well
831  // the order of daughters in the *fitter is the same as in the mother Particle
832 
833  std::vector<Particle*> daughters = mother->getDaughters();
834 
835  const unsigned nd = daughters.size();
836  unsigned l = 0;
837  std::vector<std::vector<unsigned>> pars;
838  std::vector<Particle*> allparticles;
839  for (unsigned ichild = 0; ichild < nd; ichild++) {
840  const Particle* daughter = mother->getDaughter(ichild);
841  std::vector<unsigned> pard;
842  if (daughter->getNDaughters() > 0) {
843  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
844  pars.push_back(pard);
845  allparticles.push_back(daughters[ichild]);
846  } else {
847  pard.push_back(l);
848  pars.push_back(pard);
849  allparticles.push_back(daughters[ichild]);
850  l++;
851  }
852  }
853 
854  unsigned track_count = kf.getTrackCount();
855  if (l != track_count)
856  return false;
857 
858  for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
859  TLorentzVector childMoms;
860  TVector3 childPoss;
861  TMatrixFSym childErrMatrixs(7);
862  for (unsigned int iChild : pars[iDaug]) {
863  childMoms = childMoms +
864  CLHEPToROOT::getTLorentzVector(
865  kf.getTrackMomentum(iChild));
866  childPoss = childPoss +
867  CLHEPToROOT::getTVector3(
868  kf.getTrackPosition(iChild));
869  TMatrixFSym childErrMatrix =
870  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
871  childErrMatrixs = childErrMatrixs + childErrMatrix;
872  }
873  allparticles[iDaug]->set4Vector(childMoms);
874  allparticles[iDaug]->setVertex(childPoss);
875  allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
876  }
877  }
878 
879  return true;
880  }
881 
882  bool ParticleVertexFitterModule::makeMassKFourCMother(analysis::MassFourCFitKFit& kf, Particle* mother)
883  {
884  enum analysis::KFitError::ECode fitError;
885  fitError = kf.updateMother(mother);
886  if (fitError != analysis::KFitError::kNoError)
887  return false;
888  mother->addExtraInfo("MassFourCFitProb", TMath::Prob(kf.getCHIsq(), kf.getNDF()));
889  mother->addExtraInfo("MassFourCFitChi2", kf.getCHIsq());
890  mother->addExtraInfo("MassFourCFitNDF", kf.getNDF());
891  if (m_decayString.empty() && m_updateDaughters == true) {
892  // update daughter momenta as well
893  // the order of daughters in the *fitter is the same as in the mother Particle
894 
895  std::vector<Particle*> daughters = mother->getDaughters();
896 
897  const unsigned nd = daughters.size();
898  unsigned l = 0;
899  std::vector<std::vector<unsigned>> pars;
900  std::vector<Particle*> allparticles;
901  for (unsigned ichild = 0; ichild < nd; ichild++) {
902  const Particle* daughter = mother->getDaughter(ichild);
903  std::vector<unsigned> pard;
904  if (daughter->getNDaughters() > 0) {
905  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
906  pars.push_back(pard);
907  allparticles.push_back(daughters[ichild]);
908  } else {
909  pard.push_back(l);
910  pars.push_back(pard);
911  allparticles.push_back(daughters[ichild]);
912  l++;
913  }
914  }
915 
916  unsigned track_count = kf.getTrackCount();
917  if (l != track_count)
918  return false;
919 
920  for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
921  TLorentzVector childMoms;
922  TVector3 childPoss;
923  TMatrixFSym childErrMatrixs(7);
924  for (unsigned int iChild : pars[iDaug]) {
925  childMoms = childMoms +
926  CLHEPToROOT::getTLorentzVector(
927  kf.getTrackMomentum(iChild));
928  childPoss = childPoss +
929  CLHEPToROOT::getTVector3(
930  kf.getTrackPosition(iChild));
931  TMatrixFSym childErrMatrix =
932  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
933  childErrMatrixs = childErrMatrixs + childErrMatrix;
934  }
935  allparticles[iDaug]->set4Vector(childMoms);
936  allparticles[iDaug]->setVertex(childPoss);
937  allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
938  }
939  }
940 
941  return true;
942  }
943 
944 
945  void ParticleVertexFitterModule::updateMapOfTrackAndDaughter(unsigned& l, std::vector<std::vector<unsigned>>& pars,
946  std::vector<unsigned>& parm, std::vector<Particle*>& allparticles, const Particle* daughter)
947  {
948  std::vector <Belle2::Particle*> childs = daughter->getDaughters();
949  for (unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
950  const Particle* child = daughter->getDaughter(ichild);
951  std::vector<unsigned> pard;
952  if (child->getNDaughters() > 0) {
953  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
954  parm.insert(parm.end(), pard.begin(), pard.end());
955  pars.push_back(pard);
956  allparticles.push_back(childs[ichild]);
957  } else {
958  pard.push_back(l);
959  parm.push_back(l);
960  pars.push_back(pard);
961  allparticles.push_back(childs[ichild]);
962  l++;
963  }
964  }
965  }
966 
967 
968  bool ParticleVertexFitterModule::doRaveFit(Particle* mother)
969  {
970  if ((m_decayString.empty() ||
971  (m_withConstraint == "" && m_fitType != "mass")) && mother->getNDaughters() < 2) return false;
972  if (m_withConstraint == "") analysis::RaveSetup::getInstance()->unsetBeamSpot();
973  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "mother"
974  || m_withConstraint == "iptubecut" || m_withConstraint == "btube")
975  analysis::RaveSetup::getInstance()->setBeamSpot(m_BeamSpotCenter, m_beamSpotCov);
976 
978  if (m_fitType == "mass") rf.setVertFit(false);
979 
980  if (m_decayString.empty()) {
981  rf.addMother(mother);
982  } else {
983  std::vector<const Particle*> tracksVertex = m_decaydescriptor.getSelectionParticles(mother);
984  std::vector<std::string> tracksName = m_decaydescriptor.getSelectionNames();
985 
986  if (allSelectedDaughters(mother, tracksVertex)) {
987  for (auto& itrack : tracksVertex) {
988  if (itrack != mother) rf.addTrack(itrack);
989  }
990  rf.setMother(mother);
991  } else {
992 
994  bool mothSel = false;
995  int nTrk = 0;
996  for (unsigned itrack = 0; itrack < tracksVertex.size(); itrack++) {
997  if (tracksVertex[itrack] != mother) {
998  rsf.addTrack(tracksVertex[itrack]);
999  B2DEBUG(1, "ParticleVertexFitterModule: Adding particle " << tracksName[itrack] << " to vertex fit ");
1000  nTrk++;
1001  }
1002  if (tracksVertex[itrack] == mother) mothSel = true;
1003  }
1004 
1005 
1006  // Fit one particle constrained to originate from the beam spot
1007  bool mothIPfit = false;
1008  if (tracksVertex.size() == 1 && mothSel == true && m_withConstraint != "" && nTrk == 0) {
1009  rsf.addTrack(tracksVertex[0]);
1010  if (tracksVertex[0] != mother)
1011  B2FATAL("ParticleVertexFitterModule: FATAL Error in IP constrained mother fit");
1012  nTrk++;
1013  mothIPfit = true;
1014  }
1015 
1016 
1017  TVector3 pos;
1018  TMatrixDSym RerrMatrix(3);
1019  int nvert = 0;
1020 
1021  // one track fit is not kinematic
1022  if (nTrk == 1) {
1024  for (auto& itrack : tracksVertex) {
1025  rsg.addTrack(itrack);
1026  nvert = rsg.fit("kalman");
1027  if (nvert > 0) {
1028  pos = rsg.getPos(0);
1029  RerrMatrix = rsg.getCov(0);
1030  double prob = rsg.getPValue(0);
1031  TLorentzVector mom(mother->getMomentum(), mother->getEnergy());
1032  TMatrixDSym errMatrix(7);
1033  for (int i = 0; i < 7; i++) {
1034  for (int j = 0; j < 7; j++) {
1035  if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1036  else {errMatrix[i][j] = 0;}
1037  }
1038  }
1039  if (mothIPfit) {
1040  mother->writeExtraInfo("prodVertX", pos.X());
1041  mother->writeExtraInfo("prodVertY", pos.Y());
1042  mother->writeExtraInfo("prodVertZ", pos.Z());
1043  mother->writeExtraInfo("prodVertSxx", RerrMatrix[0][0]);
1044  mother->writeExtraInfo("prodVertSxy", RerrMatrix[0][1]);
1045  mother->writeExtraInfo("prodVertSxz", RerrMatrix[0][2]);
1046  mother->writeExtraInfo("prodVertSyx", RerrMatrix[1][0]);
1047  mother->writeExtraInfo("prodVertSyy", RerrMatrix[1][1]);
1048  mother->writeExtraInfo("prodVertSyz", RerrMatrix[1][2]);
1049  mother->writeExtraInfo("prodVertSzx", RerrMatrix[2][0]);
1050  mother->writeExtraInfo("prodVertSzy", RerrMatrix[2][1]);
1051  mother->writeExtraInfo("prodVertSzz", RerrMatrix[2][2]);
1052  } else {
1053  mother->updateMomentum(mom, pos, errMatrix, prob);
1054  }
1055  return true;
1056  } else {return false;}
1057  }
1058  } else {
1059  nvert = rsf.fit();
1060  }
1061 
1062  if (nvert > 0) {
1063  pos = rsf.getPos();
1064  RerrMatrix = rsf.getVertexErrorMatrix();
1065  double prob = rsf.getPValue();
1066  TLorentzVector mom(mother->getMomentum(), mother->getEnergy());
1067  TMatrixDSym errMatrix(7);
1068  for (int i = 0; i < 7; i++) {
1069  for (int j = 0; j < 7; j++) {
1070  if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1071  else {errMatrix[i][j] = 0;}
1072  }
1073  }
1074  mother->updateMomentum(mom, pos, errMatrix, prob);
1075  } else {return false;}
1076 
1077 
1078  if (mothSel && nTrk > 1) {
1079  analysis::RaveSetup::getInstance()->setBeamSpot(pos, RerrMatrix);
1080  rf.addMother(mother);
1081  int nKfit = rf.fit();
1082  rf.updateMother();
1083  analysis::RaveSetup::getInstance()->unsetBeamSpot();
1084 
1085  if (nKfit > 0) {return true;}
1086  else return false;
1087  } else return true;
1088  }
1089  }
1090 
1091  bool okFT = false;
1092  if (m_fitType == "vertex") {
1093  okFT = true;
1094  int nVert = rf.fit();
1095  rf.updateMother();
1096  if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1097  if (nVert != 1) return false;
1098  }
1099  if (m_fitType == "mass") {
1100  // add protection
1101  okFT = true;
1102  rf.setMassConstFit(true);
1103  rf.setVertFit(false);
1104  int nVert = rf.fit();
1105  rf.updateMother();
1106  if (nVert != 1) return false;
1107  };
1108  if (m_fitType == "massvertex") {
1109  okFT = true;
1110  rf.setMassConstFit(true);
1111  int nVert = rf.fit();
1112  rf.updateMother();
1113  if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1114  if (nVert != 1) return false;
1115  };
1116  if (!okFT) {
1117  B2FATAL("fitType : " << m_fitType << " ***invalid fit type ");
1118  }
1119 
1120  return true;
1121  }
1122 
1123  bool ParticleVertexFitterModule::allSelectedDaughters(const Particle* mother,
1124  const std::vector<const Particle*>& tracksVertex)
1125  {
1126  bool isAll = false;
1127  if (mother->getNDaughters() == 0) return false;
1128 
1129  int nNotIncluded = mother->getNDaughters();
1130 
1131  for (unsigned i = 0; i < mother->getNDaughters(); i++) {
1132  bool dauOk = false;
1133  for (auto& vi : tracksVertex) {
1134  if (vi == mother->getDaughter(i)) {
1135  nNotIncluded = nNotIncluded - 1;
1136  dauOk = true;
1137  }
1138  }
1139  if (!dauOk) {
1140  if (allSelectedDaughters(mother->getDaughter(i), tracksVertex)) nNotIncluded--;
1141  }
1142  }
1143  if (nNotIncluded == 0) isAll = true;
1144  return isAll;
1145  }
1146 
1147  bool ParticleVertexFitterModule::addChildofParticletoKFit(analysis::FourCFitKFit& kf, const Particle* particle)
1148  {
1149  for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1150  const Particle* child = particle->getDaughter(ichild);
1151  if (child->getNDaughters() > 0) addChildofParticletoKFit(kf, child);
1152  else {
1153  if (child->getPValue() < 0) return false; // error matrix not valid
1154 
1155  kf.addParticle(child);
1156  }
1157  }
1158  return true;
1159  }
1160 
1161  bool ParticleVertexFitterModule::addChildofParticletoMassKFit(analysis::MassFourCFitKFit& kf, const Particle* particle,
1162  std::vector<unsigned>& particleId)
1163  {
1164  for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1165  const Particle* child = particle->getDaughter(ichild);
1166  if (child->getNDaughters() > 0) {
1167  bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
1168  std::abs(child->getPDGCode())) != m_massConstraintList.end();
1169  std::vector<unsigned> childId;
1170  addChildofParticletoMassKFit(kf, child, childId);
1171  if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
1172  particleId.insert(particleId.end(), childId.begin(), childId.end());
1173  } else {
1174  if (child->getPValue() < 0) return false; // error matrix not valid
1175  kf.addParticle(child);
1176  particleId.push_back(kf.getTrackCount() - 1);
1177  }
1178  }
1179  return true;
1180  }
1181 
1182  void ParticleVertexFitterModule::addIPProfileToKFit(analysis::VertexFitKFit& kv)
1183  {
1184  HepPoint3D pos(0.0, 0.0, 0.0);
1185  CLHEP::HepSymMatrix covMatrix(3, 0);
1186 
1187  for (int i = 0; i < 3; i++) {
1188  pos[i] = m_BeamSpotCenter(i);
1189  for (int j = 0; j < 3; j++) {
1190  covMatrix[i][j] = m_beamSpotCov(i, j);
1191  }
1192  }
1193 
1194  kv.setIpProfile(pos, covMatrix);
1195  }
1196 
1197  void ParticleVertexFitterModule::addIPTubeToKFit(analysis::VertexFitKFit& kv)
1198  {
1199  CLHEP::HepSymMatrix err(7, 0);
1200 
1201  for (int i = 0; i < 3; i++) {
1202  for (int j = 0; j < 3; j++) {
1203  err[i + 4][j + 4] = m_beamSpotCov(i, j);
1204  }
1205  }
1206 
1207  PCmsLabTransform T;
1208  TLorentzVector iptube_mom = T.getBeamFourMomentum();
1209 
1210  kv.setIpTubeProfile(
1211  ROOTToCLHEP::getHepLorentzVector(iptube_mom),
1212  ROOTToCLHEP::getPoint3D(m_BeamSpotCenter),
1213  err,
1214  0.);
1215  }
1216 
1217  void ParticleVertexFitterModule::findConstraintBoost(double cut)
1218  {
1219  PCmsLabTransform T;
1220 
1221  TVector3 boost = T.getBoostVector();
1222  TVector3 boostDir = boost.Unit();
1223 
1224  TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1225  beamSpotCov(2, 2) = cut * cut;
1226  double thetab = boostDir.Theta();
1227  double phib = boostDir.Phi();
1228 
1229  double stb = TMath::Sin(thetab);
1230  double ctb = TMath::Cos(thetab);
1231  double spb = TMath::Sin(phib);
1232  double cpb = TMath::Cos(phib);
1233 
1234 
1235  TMatrix rz(3, 3); rz(2, 2) = 1;
1236  rz(0, 0) = cpb; rz(0, 1) = spb;
1237  rz(1, 0) = -1 * spb; rz(1, 1) = cpb;
1238 
1239  TMatrix ry(3, 3); ry(1, 1) = 1;
1240  ry(0, 0) = ctb; ry(0, 2) = -1 * stb;
1241  ry(2, 0) = stb; ry(2, 2) = ctb;
1242 
1243  TMatrix r(3, 3); r.Mult(rz, ry);
1244  TMatrix rt(3, 3); rt.Transpose(r);
1245 
1246  TMatrix TubePart(3, 3); TubePart.Mult(rt, beamSpotCov);
1247  TMatrix Tube(3, 3); Tube.Mult(TubePart, r);
1248 
1249  m_beamSpotCov(0, 0) = Tube(0, 0); m_beamSpotCov(0, 1) = Tube(0, 1); m_beamSpotCov(0, 2) = Tube(0, 2);
1250  m_beamSpotCov(1, 0) = Tube(1, 0); m_beamSpotCov(1, 1) = Tube(1, 1); m_beamSpotCov(1, 2) = Tube(1, 2);
1251  m_beamSpotCov(2, 0) = Tube(2, 0); m_beamSpotCov(2, 1) = Tube(2, 1); m_beamSpotCov(2, 2) = Tube(2, 2);
1252  }
1253 
1254  void ParticleVertexFitterModule::smearBeamSpot(double width)
1255  {
1256  TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1257  for (int i = 0; i < 3; i++)
1258  beamSpotCov(i, i) += width * width;
1259 
1260  m_beamSpotCov = beamSpotCov;
1261  }
1263 } // end Belle2 namespace
For each MCParticle with hits in the CDC, this class stores some summarising information on those hit...
Definition: Btube.h:28
TMatrixFSym getTubeMatrix() const
Returns Btube matrix.
Definition: Btube.h:79
Eigen::Matrix< double, 3, 1 > getTubeCenter() const
Returns Btube center.
Definition: Btube.h:65
Base class for Modules.
Definition: Module.h:72
Class to hold Lorentz transformations from/to CMS and boost vector.
TVector3 getBoostVector() const
Returns boost vector (beta=p/E)
TLorentzVector getBeamFourMomentum() const
Returns LAB four-momentum of e+e-.
Class to store reconstructed particles.
Definition: Particle.h:74
float getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
Definition: Particle.cc:621
void updateMomentum(const TLorentzVector &p4, const TVector3 &vertex, const TMatrixFSym &errMatrix, float pValue)
Sets Lorentz vector, position, 7x7 error matrix and p-value.
Definition: Particle.h:340
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:392
TLorentzVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:477
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
Definition: Particle.cc:430
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:639
FourCFitKFit is a derived class from KFitBase to perform 4 momentum-constraint kinematical fit.
Definition: FourCFitKFit.h:30
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode setFourMomentum(const TLorentzVector &m)
Set an 4 Momentum for the FourC-constraint fit.
Definition: FourCFitKFit.cc:77
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a four momentum-constraint fit.
const CLHEP::HepSymMatrix getTrackError(const int id) const
Get an error matrix of the track.
Definition: KFitBase.cc:168
const CLHEP::HepLorentzVector getTrackMomentum(const int id) const
Get a Lorentz vector of the track.
Definition: KFitBase.cc:154
const HepPoint3D getTrackPosition(const int id) const
Get a position of the track.
Definition: KFitBase.cc:161
virtual int getNDF(void) const
Get an NDF of the fit.
Definition: KFitBase.cc:114
enum KFitError::ECode setMagneticField(const double mf)
Change a magnetic field from the default value KFitConst::kDefaultMagneticField.
Definition: KFitBase.cc:93
enum KFitError::ECode addParticle(const Particle *particle)
Add a particle to the fitter.
Definition: KFitBase.cc:59
int getTrackCount(void) const
Get the number of added tracks.
Definition: KFitBase.cc:107
ECode
ECode is a error code enumerate.
Definition: KFitError.h:34
MassFitKFit is a derived class from KFitBase to perform mass-constraint kinematical fit.
Definition: MassFitKFit.h:33
enum KFitError::ECode setVertex(const HepPoint3D &v)
Set an initial vertex position for the mass-constraint fit.
Definition: MassFitKFit.cc:43
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
Definition: MassFitKFit.cc:695
enum KFitError::ECode doFit(void)
Perform a mass-constraint fit.
Definition: MassFitKFit.cc:279
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-constraint fit.
Definition: MassFitKFit.cc:67
enum KFitError::ECode setVertexError(const CLHEP::HepSymMatrix &e)
Set an initial vertex error matrix for the mass-constraint fit.
Definition: MassFitKFit.cc:51
MassFourCFitKFit is a derived class from KFitBase to perform mass and 4 momentum-constraint kinematic...
enum KFitError::ECode addMassConstraint(const double m, std::vector< unsigned > &childTrackId)
Set an invariant mass of daughter particle for the mass-four-momentum-constraint fit.
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode setFourMomentum(const TLorentzVector &m)
Set an 4 Momentum for the mass-four-constraint fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a mass-four-momentum-constraint fit.
MassPointingVertexFitKFit is a derived class from KFitBase It performs a kinematical fit with three c...
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a mass-vertex-pointing constraint fit.
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-vertex-pointing constraint fit.
enum KFitError::ECode setProductionVertex(const HepPoint3D &v)
Set the production vertex of the particle.
MassVertexFitKFit is a derived class from KFitBase to perform mass-vertex-constraint kinematical fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode doFit(void)
Perform a mass-vertex-constraint fit.
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-vertex constraint fit.
The RaveKinematicVertexFitter class is part of the RaveInterface together with RaveSetup.
void addMother(const Particle *aMotherParticlePtr)
All daughters of the argument of this function will be used as input for the vertex fit.
void addTrack(const Particle *aParticlePtr)
add a track (in the format of a Belle2::Particle) to set of tracks that should be fitted to a vertex
void setMother(const Particle *aMotherParticlePtr)
Set Mother particle for Vertex/momentum update.
int fit()
do the kinematic vertex fit with all tracks previously added with the addTrack or addMother function.
void setVertFit(bool isVertFit=true)
Set vertex fit: set false in case of mass fit only.
void setMassConstFit(bool isConstFit=true)
Set mass constrained fit
TVector3 getPos()
get the position of the fitted vertex.
double getPValue()
get the p value of the fitted vertex.
TMatrixDSym getVertexErrorMatrix()
get the covariance matrix (3x3) of the of the fitted vertex position.
void updateDaughters()
update the Daughters particles
The RaveVertexFitter class is part of the RaveInterface together with RaveSetup.
TMatrixDSym getCov(VecSize vertexId=0) const
get the covariance matrix (3x3) of the of the fitted vertex position.
int fit(std::string options="default")
do the vertex fit with all tracks previously added with the addTrack or addMother function.
void addTrack(const Particle *const aParticlePtr)
add a track (in the format of a Belle2::Particle) to set of tracks that should be fitted to a vertex
TVector3 getPos(VecSize vertexId=0) const
get the position of the fitted vertex.
double getPValue(VecSize vertexId=0) const
get the p value of the fitted vertex.
VertexFitKFit is a derived class from KFitBase to perform vertex-constraint kinematical fit.
Definition: VertexFitKFit.h:32
enum KFitError::ECode setIpTubeProfile(const CLHEP::HepLorentzVector &p, const HepPoint3D &x, const CLHEP::HepSymMatrix &e, const double q)
Set a virtual IP-tube track for the vertex constraint fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
const CLHEP::HepSymMatrix getVertexError(void) const
Get a fitted vertex error matrix.
enum KFitError::ECode doFit(void)
Perform a vertex-constraint fit.
enum KFitError::ECode setIpProfile(const HepPoint3D &ip, const CLHEP::HepSymMatrix &ipe)
Set an IP-ellipsoid shape for the vertex constraint fit.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
#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.