Belle II Software  release-08-01-10
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 header.
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/geometry/B2Vector3.h>
16 #include <framework/logging/Logger.h>
17 #include <framework/particledb/EvtGenDatabasePDG.h>
18 
19 // dataobjects
20 #include <analysis/dataobjects/Particle.h>
21 #include <analysis/dataobjects/Btube.h>
22 #include <mdst/dataobjects/V0.h>
23 
24 // utilities
25 #include <analysis/utility/CLHEPToROOT.h>
26 #include <analysis/utility/PCmsLabTransform.h>
27 #include <analysis/utility/ParticleCopy.h>
28 #include <analysis/utility/ROOTToCLHEP.h>
29 
30 // Magnetic field
31 #include <framework/geometry/BFieldManager.h>
32 
33 // KFit
34 #include <analysis/VertexFitting/KFit/KFitConst.h>
35 
36 #include <TVector.h>
37 #include <TRotation.h>
38 #include <TMath.h>
39 
40 using namespace std;
41 using namespace Belle2;
42 
43 
44 //-----------------------------------------------------------------
45 // Register module
46 //-----------------------------------------------------------------
47 
48 REG_MODULE(ParticleVertexFitter);
49 
50 //-----------------------------------------------------------------
51 // Implementation
52 //-----------------------------------------------------------------
53 
54 ParticleVertexFitterModule::ParticleVertexFitterModule() : Module(),
55  m_Bfield(0)
56 {
57  // set module description (e.g. insert text)
58  setDescription("Vertex fitter for modular analysis");
60 
61  // Add parameters
62  addParam("listName", m_listName, "name of particle list", string(""));
63  addParam("confidenceLevel", m_confidenceLevel,
64  "Confidence level to accept the fit. Particle candidates with "
65  "p-value less than confidenceLevel are removed from the particle "
66  "list. If set to -1, all candidates are kept; if set to 0, "
67  "the candidates failing the fit are removed.",
68  0.001);
69  addParam("vertexFitter", m_vertexFitter, "KFit or Rave", string("KFit"));
70  addParam("fitType", m_fitType, "type of the kinematic fit (vertex, massvertex, mass)", string("vertex"));
71  addParam("withConstraint", m_withConstraint,
72  "additional constraint on vertex: ipprofile, iptube, mother, iptubecut, pointing, btube",
73  string(""));
74  addParam("decayString", m_decayString, "specifies which daughter particles are included in the kinematic fit", string(""));
75  addParam("updateDaughters", m_updateDaughters, "true: update the daughters after the vertex fit", false);
76  addParam("smearing", m_smearing, "smear IP tube width by given length", 0.002);
77  addParam("recoilMass", m_recoilMass, "recoil invariant mass (GeV)", 0.);
78  addParam("massConstraintList", m_massConstraintList,
79  "Type::[int]. List of daughter particles to mass constrain with int = pdg code. (only for MassFourCKFit)", {});
80  addParam("massConstraintListParticlename", m_massConstraintListParticlename,
81  "Type::[string]. List of daughter particles to mass constrain with string = particle name. (only for MassFourCKFit)", {});
82 }
83 
85 {
86  // Particle list with name m_listName has to exist
87  m_plist.isRequired(m_listName);
88 
89  // magnetic field
90  m_Bfield = BFieldManager::getFieldInTesla(ROOT::Math::XYZVector(0, 0, 0)).Z();
91 
92  // RAVE setup
93  if (m_vertexFitter == "Rave")
95 
96  B2DEBUG(1, "ParticleVertexFitterModule : magnetic field = " << m_Bfield);
97 
98 
99  if (m_decayString != "")
101 
102  if ((m_massConstraintList.size()) == 0 && (m_massConstraintListParticlename.size()) > 0) {
103  for (auto& containedParticle : m_massConstraintListParticlename) {
104  TParticlePDG* particletemp = TDatabasePDG::Instance()->GetParticle((containedParticle).c_str());
105  m_massConstraintList.push_back(particletemp->PdgCode());
106  }
107  }
108 
109  B2INFO("ParticleVertexFitter: Performing " << m_fitType << " fit on " << m_listName << " using " << m_vertexFitter);
110  if (m_decayString != "")
111  B2INFO("ParticleVertexFitter: Using specified decay string: " << m_decayString);
112  if (m_withConstraint != "")
113  B2INFO("ParticleVertexFitter: Additional " << m_withConstraint << " will be applied");
114 
115 }
116 
118 {
119  //TODO: set magnetic field for each run
120  //m_Bfield = BFieldMap::Instance().getBField(B2Vector3D(0,0,0)).Z();
121  //TODO: set IP spot size for each run
122 }
123 
125 {
126  if (m_vertexFitter == "Rave")
128 
129  m_BeamSpotCenter = m_beamSpotDB->getIPPosition();
130  m_beamSpotCov.ResizeTo(3, 3);
131  TMatrixDSym beamSpotCov(3);
132  if (m_withConstraint == "ipprofile") m_beamSpotCov = m_beamSpotDB->getCovVertex();
133  if (m_withConstraint == "iptube") {
134  if (m_smearing > 0 && m_vertexFitter == "KFit") {
136  } else {
138  }
139  }
140  if (m_withConstraint == "iptubecut") { // for development purpose only
141  m_BeamSpotCenter = B2Vector3D(0.001, 0., .013);
142  findConstraintBoost(0.03);
143  }
144  if ((m_vertexFitter == "Rave") && (m_withConstraint == "ipprofile" || m_withConstraint == "iptube"
145  || m_withConstraint == "mother" || m_withConstraint == "iptubecut" || m_withConstraint == "btube"))
147 
148  std::vector<unsigned int> toRemove;
149  unsigned int nParticles = m_plist->getListSize();
150 
151  for (unsigned iPart = 0; iPart < nParticles; iPart++) {
152  Particle* particle = m_plist->getParticle(iPart);
153  m_hasCovMatrix = false;
154  if (m_updateDaughters == true) {
155  if (m_decayString.empty() || m_vertexFitter == "KFit")
156  ParticleCopy::copyDaughters(particle);
157  else B2ERROR("Daughters update works only when all daughters are selected. Daughters will not be updated");
158  }
159 
160  if (m_withConstraint == "mother") {
161  m_BeamSpotCenter = B2Vector3D(particle->getVertex().x(), particle->getVertex().y(), particle->getVertex().z());
162  m_beamSpotCov = particle->getVertexErrorMatrix();
163  }
164 
165  TMatrixFSym mother_errMatrix(7);
166  mother_errMatrix = particle->getMomentumVertexErrorMatrix();
167  for (int k = 0; k < 7; k++) {
168  for (int j = 0; j < 7; j++) {
169  if (mother_errMatrix[k][j] > 0) {
170  m_hasCovMatrix = true;
171  }
172  }
173  }
174 
175  bool hasTube = true;
176  if (m_withConstraint == "btube") {
177  Btube* Ver = particle->getRelatedTo<Btube>();
178  if (!Ver) {
179  hasTube = false;
180  toRemove.push_back(particle->getArrayIndex());
181  } else {
182  m_BeamSpotCenter.SetXYZ(Ver->getTubeCenter()(0, 0), Ver->getTubeCenter()(1, 0), Ver->getTubeCenter()(2, 0));
183  m_beamSpotCov = Ver->getTubeMatrix();
184  }
185  }
186  bool ok = false;
187  if (hasTube) {
188  ok = doVertexFit(particle);
189  }
190  if (!ok)
191  particle->setPValue(-1);
192  if (particle->getPValue() < m_confidenceLevel)
193  toRemove.push_back(particle->getArrayIndex());
194 
195  }
196  m_plist->removeParticles(toRemove);
197 
198  //free memory allocated by rave. initialize() would be enough, except that we must clean things up before program end...
199  if (m_vertexFitter == "Rave")
201 }
202 
204 {
205  // steering starts here
206 
207  if (m_Bfield == 0) {
208  B2FATAL("ParticleVertexFitter: No magnetic field");
209  }
210 
211  if (m_withConstraint != "ipprofile" &&
212  m_withConstraint != "iptube" &&
213  m_withConstraint != "mother" &&
214  m_withConstraint != "iptubecut" &&
215  m_withConstraint != "pointing" &&
216  m_withConstraint != "btube" &&
217  m_withConstraint != "")
218  B2FATAL("ParticleVertexFitter: " << m_withConstraint << " ***invalid Constraint ");
219 
220  bool ok = false;
221  // fits with KFit
222  if (m_vertexFitter == "KFit") {
223 
224  if (m_decayString != "" and m_fitType != "vertex")
225  B2FATAL("ParticleVertexFitter: KFit does not support yet selection of daughters via decay string except for vertex fit!");
226 
227  // vertex fit
228  if (m_fitType == "vertex") {
229  if (m_withConstraint == "ipprofile") {
230  ok = doKVertexFit(mother, true, false);
231  } else if (m_withConstraint == "iptube") {
232  ok = doKVertexFit(mother, false, true);
233  } else {
234  ok = doKVertexFit(mother, false, false);
235  }
236  }
237 
238  // mass-constrained vertex fit
239  if (m_fitType == "massvertex") {
240  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
241  B2FATAL("ParticleVertexFitter: Invalid options - mass-constrained fit using KFit does not work with iptube or ipprofile constraint.");
242  } else if (m_withConstraint == "pointing") {
243  ok = doKMassPointingVertexFit(mother);
244  } else {
245  ok = doKMassVertexFit(mother);
246  }
247  }
248 
249  // mass fit
250  if (m_fitType == "mass") {
251  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
252  B2FATAL("ParticleVertexFitter: Invalid options - mass fit using KFit does not work with iptube or ipprofile constraint.");
253  } else {
254  ok = doKMassFit(mother);
255  }
256  }
257 
258  // four C fit
259  if (m_fitType == "fourC") {
260  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
261  B2FATAL("ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
262  } else {
263  ok = doKFourCFit(mother);
264  }
265  }
266 
267  // four mass C fit
268  if (m_fitType == "massfourC") {
269  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
270  B2FATAL("ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
271  } else {
272  ok = doKMassFourCFit(mother);
273  }
274  }
275 
276  // recoil mass C fit
277  if (m_fitType == "recoilmass") {
278  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
279  B2FATAL("ParticleVertexFitter: Invalid options - recoil mass fit using KFit does not work with iptube or ipprofile constraint.");
280  } else {
281  ok = doKRecoilMassFit(mother);
282  }
283  }
284 
285  // invalid KFit fit type
286  if (m_fitType != "vertex"
287  && m_fitType != "massvertex"
288  && m_fitType != "mass"
289  && m_fitType != "fourC"
290  && m_fitType != "massfourC"
291  && m_fitType != "recoilmass")
292  B2FATAL("ParticleVertexFitter: " << m_fitType << " ***invalid fit type for the vertex fitter ");
293  }
294 
295  // fits using Rave
296  if (m_vertexFitter == "Rave") {
297  try {
298  ok = doRaveFit(mother);
299  } catch (const rave::CheckedFloatException&) {
300  B2ERROR("Invalid inputs (nan/inf)?");
301  ok = false;
302  }
303  }
304 
305  // invalid fitter
306  if (m_vertexFitter != "KFit" && m_vertexFitter != "Rave")
307  B2FATAL("ParticleVertexFitter: " << m_vertexFitter << " ***invalid vertex fitter ");
308 
309  if (!ok) return false;
310 
311  // steering ends here
312 
313  //if (mother->getPValue() < m_confidenceLevel) return false;
314  return true;
315 
316 }
317 
318 bool ParticleVertexFitterModule::fillFitParticles(const Particle* mother, std::vector<const Particle*>& fitChildren,
319  std::vector<const Particle*>& twoPhotonChildren)
320 {
321  if (m_decayString.empty()) {
322  // if decayString is empty, just use all primary daughters
323  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
324  const Particle* child = mother->getDaughter(ichild);
325  // This if allows to skip the daughters, which cannot be used in the fits, particularly K_L0 from KLM.
326  // Useful for fully-inclusive particles.
327  if (mother->getProperty() == Particle::PropertyFlags::c_IsUnspecified and child->getPValue() < 0) {
328  continue;
329  }
330  fitChildren.push_back(child);
331  }
332  } else {
333  fitChildren = m_decaydescriptor.getSelectionParticles(mother);
334  }
335 
336  auto itr = fitChildren.begin();
337  while (itr != fitChildren.end()) {
338  const Particle* child = *itr;
339 
340  if (child->getPValue() < 0) {
341  B2WARNING("Daughter with PDG code " << child->getPDGCode() << " does not have a valid error matrix.");
342  return false; // error matrix not valid
343  }
344  bool isTwoPhotonParticle = false;
345  if (m_hasCovMatrix == false) {
346  if (child->getPDGCode() == Const::pi0.getPDGCode() or child->getPDGCode() == 221) { // pi0 or eta
347  if (child->getNDaughters() == 2) {
348  if (child->getDaughter(0)->getPDGCode() == Const::photon.getPDGCode()
349  && child->getDaughter(1)->getPDGCode() == Const::photon.getPDGCode()) {
350  isTwoPhotonParticle = true;
351  }
352  }
353  }
354  }
355  if (isTwoPhotonParticle) {
356  // move children from fitChildren to twoPhotonChildren
357  twoPhotonChildren.push_back(child);
358  itr = fitChildren.erase(itr);
359  } else {
360  itr++;
361  }
362  }
363 
364  return true;
365 }
366 
367 bool ParticleVertexFitterModule::fillNotFitParticles(const Particle* mother, std::vector<const Particle*>& notFitChildren,
368  const std::vector<const Particle*>& fitChildren)
369 {
370  if (fitChildren.empty())
371  B2WARNING("[ParticleVertexFitterModule::fillNotFitParticles] fitChildren is empty! Please call fillFitParticles firstly");
372  if (!notFitChildren.empty())
373  B2WARNING("[ParticleVertexFitterModule::fillNotFitParticles] notFitChildren is NOT empty!"
374  << " The function should be called only once");
375 
376  if (m_decayString.empty())
377  // if decayString is empty, just use all primary daughters
378  return true;
379 
380  std::function<bool(const Particle*)> funcCheckInFit =
381  [&funcCheckInFit, &notFitChildren, fitChildren](const Particle * part) {
382 
383  // check if the given particle in fitChildren
384  // if it is included, return true
385  if (std::find(fitChildren.begin(), fitChildren.end(), part) != fitChildren.end())
386  return true;
387 
388  // if not, firstly check if particle has children
389  if (part->getNDaughters() == 0)
390  // if it has no children (=final-state-particle), return false
391  return false;
392 
393  // here, the given particle is not in fitChildren and has children
394  bool isAnyChildrenInFit = false;
395  vector<const Particle*> notFitChildren_tmp;
396  for (unsigned ichild = 0; ichild < part->getNDaughters(); ichild++) {
397  // call funcCheckInFit recursively for all children
398  const Particle* child = part->getDaughter(ichild);
399  bool isChildrenInFit = funcCheckInFit(child);
400  isAnyChildrenInFit = isChildrenInFit or isAnyChildrenInFit;
401 
402  // if the child is not in fitChildren, fill the child in a temporary vector
403  if (!isChildrenInFit)
404  notFitChildren_tmp.push_back(child);
405  }
406 
407  // if there are a sister in fitChildren, the children in the temporary vector will be filled in notFitChildren
408  if (isAnyChildrenInFit)
409  notFitChildren.insert(notFitChildren.end(), notFitChildren_tmp.begin(), notFitChildren_tmp.end());
410 
411  // if no children in fitChildren, the given particle should be filled instead of all children.
412 
413  return isAnyChildrenInFit;
414  };
415 
416 
417  // call funcCheckInFit for all primary children
418  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
419  const Particle* child = mother->getDaughter(ichild);
420  bool isGivenParticleOrAnyChildrenInFit = funcCheckInFit(child);
421  if (!isGivenParticleOrAnyChildrenInFit)
422  notFitChildren.push_back(child);
423  }
424 
425  return true;
426 }
427 
429  const analysis::VertexFitKFit& kv)
430 {
431  // TODO: something like setGammaError is necessary
432  // this is just workaround for the moment
433 
434  const Particle* g1Orig = preFit->getDaughter(0);
435  const Particle* g2Orig = preFit->getDaughter(1);
436  Particle g1Temp(g1Orig->get4Vector(), 22);
437  Particle g2Temp(g2Orig->get4Vector(), 22);
438 
439  TMatrixFSym g1ErrMatrix = g1Orig->getMomentumVertexErrorMatrix();
440  TMatrixFSym g2ErrMatrix = g2Orig->getMomentumVertexErrorMatrix();
441 
442  ROOT::Math::XYZVector pos(kv.getVertex().x(), kv.getVertex().y(), kv.getVertex().z());
443  CLHEP::HepSymMatrix posErrorMatrix = kv.getVertexError();
444 
445  TMatrixFSym errMatrix(3);
446  for (int i = 0; i < 3; i++)
447  for (int j = 0; j < 3; j++)
448  errMatrix(i, j) = posErrorMatrix[i][j];
449 
450  g1ErrMatrix.SetSub(4, errMatrix);
451  g2ErrMatrix.SetSub(4, errMatrix);
452 
453  g1Temp.updateMomentum(g1Orig->get4Vector(), pos, g1ErrMatrix, 1.0);
454  g2Temp.updateMomentum(g2Orig->get4Vector(), pos, g2ErrMatrix, 1.0);
455 
456  // perform the mass fit for the two-photon particle
459 
460  km.addParticle(&g1Temp);
461  km.addParticle(&g2Temp);
462 
463  km.setVertex(kv.getVertex());
465  km.setInvariantMass(preFit->getPDGMass());
466 
467  int err = km.doFit();
468  if (err != 0) {
469  return false;
470  }
471 
472  // The update of the daughters is disabled for this mass fit.
473  bool updateDaughters = m_updateDaughters;
474  m_updateDaughters = false;
475  bool ok = makeKMassMother(km, postFit);
476  m_updateDaughters = updateDaughters;
477 
478  return ok;
479 }
480 
481 bool ParticleVertexFitterModule::doKVertexFit(Particle* mother, bool ipProfileConstraint, bool ipTubeConstraint)
482 {
483  if ((mother->getNDaughters() < 2 && !ipTubeConstraint) || mother->getNDaughters() < 1) return false;
484 
485  std::vector<const Particle*> fitChildren;
486  std::vector<const Particle*> twoPhotonChildren;
487  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
488 
489  if (!validChildren)
490  return false;
491 
492  std::vector<const Particle*> notFitChildren;
493  fillNotFitParticles(mother, notFitChildren, fitChildren);
494 
495 
496  if (twoPhotonChildren.size() > 1) {
497  B2FATAL("[ParticleVertexFitterModule::doKVertexFit] Vertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
498  }
499 
500  if ((fitChildren.size() < 2 && !ipTubeConstraint) || fitChildren.size() < 1) {
501  B2WARNING("[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is too low.");
502  return false;
503  }
504 
505  // Initialise the Fitter
508 
509  if (mother->getV0()) {
510  HepPoint3D V0vertex_heppoint(mother->getV0()->getFittedVertexX(),
511  mother->getV0()->getFittedVertexY(),
512  mother->getV0()->getFittedVertexZ());
513  kv.setInitialVertex(V0vertex_heppoint);
514  }
515 
516  for (auto& child : fitChildren)
517  kv.addParticle(child);
518 
519  if (ipProfileConstraint)
520  addIPProfileToKFit(kv);
521 
522  if (ipTubeConstraint)
523  addIPTubeToKFit(kv);
524 
525  // Perform vertex fit using only the particles with valid error matrices
526  int err = kv.doFit();
527  if (err != 0)
528  return false;
529 
530  double chi2_track = getChi2TracksLBoost(kv);
531  unsigned track_count = kv.getTrackCount();
532  mother->writeExtraInfo("chiSquared_trackL", chi2_track);
533  mother->writeExtraInfo("kFit_nTracks", track_count);
534 
535  bool ok = false;
536  if (twoPhotonChildren.size() == 0)
537  // in the case daughters do not include pi0 - this is it (fit done)
538  ok = makeKVertexMother(kv, mother);
539  else if (twoPhotonChildren.size() == 1) {
540  // there is a daughter reconstructed from two photons so without position information
541  // 1. determine vertex based on all other valid daughters
542  // 2. set position and error matrix of two-photon daughter to previously determined vertex
543  // 3. redo the fit using all particles (including two-photon particle this time)
544 
545  const Particle* twoPhotonDaughter = twoPhotonChildren[0];
546  Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
547  ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
548  if (!ok)
549  return false;
550 
551  // finally perform the fit using all daughter particles
554 
555  for (auto& child : fitChildren)
556  kv2.addParticle(child);
557 
558  kv2.addParticle(&fixedTwoPhotonDaughter);
559 
560  if (ipProfileConstraint)
561  addIPProfileToKFit(kv2);
562 
563  err = kv2.doFit();
564 
565  if (err != 0)
566  return false;
567 
568  ok = makeKVertexMother(kv2, mother);
569  }
570 
571  // update 4-vector using not-fit-particles
572  ROOT::Math::PxPyPzEVector total4Vector(mother->get4Vector());
573  for (auto& child : notFitChildren)
574  total4Vector += child->get4Vector();
575  mother->set4Vector(total4Vector);
576 
577  return ok;
578 }
579 
581 {
582  if (mother->getNDaughters() < 2) return false;
583 
584  std::vector<const Particle*> fitChildren;
585  std::vector<const Particle*> twoPhotonChildren;
586  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
587 
588  if (!validChildren)
589  return false;
590 
591  if (twoPhotonChildren.size() > 1) {
592  B2FATAL("[ParticleVertexFitterModule::doKVertexFit] MassVertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
593  }
594 
595  if (fitChildren.size() < 2) {
596  B2WARNING("[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
597  return false;
598  }
599 
600  bool ok = false;
601  if (twoPhotonChildren.size() == 0) {
602  // Initialise the Fitter
605 
606  if (mother->getV0()) {
607  HepPoint3D V0vertex_heppoint(mother->getV0()->getFittedVertexX(),
608  mother->getV0()->getFittedVertexY(),
609  mother->getV0()->getFittedVertexZ());
610  kmv.setInitialVertex(V0vertex_heppoint);
611  }
612 
613  for (auto child : fitChildren)
614  kmv.addParticle(child);
615 
616  kmv.setInvariantMass(mother->getPDGMass());
617  int err = kmv.doFit();
618  if (err != 0)
619  return false;
620 
621  // in the case daughters do not include particles with two photon daughters like pi0 - this is it (fit done)
622  ok = makeKMassVertexMother(kmv, mother);
623  } else if (twoPhotonChildren.size() == 1) {
624  // there is a daughter reconstructed from two photons so without position information
625  // 1. determine vertex based on all other valid daughters
626  // 2. set position and error matrix of two-photon daughter to previously determined vertex
627  // 3. redo the fit using all particles (including two-photon particle this time)
628 
631 
632  for (auto child : fitChildren)
633  kv.addParticle(child);
634 
635  // Perform vertex fit using only the particles with valid error matrices
636  int err = kv.doFit();
637  if (err != 0)
638  return false;
639 
640  const Particle* twoPhotonDaughter = twoPhotonChildren[0];
641  Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
642  ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
643  if (!ok)
644  return false;
645 
646  // finally perform the fit using all daughter particles
649 
650  for (auto child : fitChildren)
651  kmv2.addParticle(child);
652  kmv2.addParticle(&fixedTwoPhotonDaughter);
653 
654  kmv2.setInvariantMass(mother->getPDGMass());
655  err = kmv2.doFit();
656 
657  if (err != 0)
658  return false;
659 
660  ok = makeKMassVertexMother(kmv2, mother);
661  }
662 
663  return ok;
664 
665 }
666 
668 {
669  if (!(mother->hasExtraInfo("prodVertX") && mother->hasExtraInfo("prodVertY") && mother->hasExtraInfo("prodVertZ"))) {
670  return false;
671  }
672 
673  if (mother->getNDaughters() < 2) return false;
674 
675  std::vector<const Particle*> fitChildren;
676  std::vector<const Particle*> twoPhotonChildren;
677  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
678 
679  if (!validChildren)
680  return false;
681 
682  if (twoPhotonChildren.size() > 0) {
683  B2FATAL("[ParticleVertexFitterModule::doKMassPointingVertexFit] MassPointingVertex fit using KFit does not support fit with two-photon daughters (yet).");
684  }
685 
686  if (fitChildren.size() < 2) {
687  B2WARNING("[ParticleVertexFitterModule::doKMassPointingVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
688  return false;
689  }
690 
691  bool ok = false;
692  // Initialise the Fitter
695 
696  for (auto child : fitChildren)
697  kmpv.addParticle(child);
698 
699  kmpv.setInvariantMass(mother->getPDGMass());
700  HepPoint3D productionVertex(mother->getExtraInfo("prodVertX"),
701  mother->getExtraInfo("prodVertY"),
702  mother->getExtraInfo("prodVertZ"));
703  kmpv.setProductionVertex(productionVertex);
704  int err = kmpv.doFit();
705  if (err != 0) return false;
706 
707  ok = makeKMassPointingVertexMother(kmpv, mother);
708 
709  return ok;
710 }
711 
713 {
714  if (mother->getNDaughters() < 2) return false;
715 
718 
719  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
720  const Particle* child = mother->getDaughter(ichild);
721 
722  if (child->getPValue() < 0) return false; // error matrix not valid
723 
724  km.addParticle(child);
725  }
726 
727  // apply mass constraint
728  km.setInvariantMass(mother->getPDGMass());
729 
730  int err = km.doFit();
731 
732  if (err != 0) return false;
733 
734  bool ok = makeKMassMother(km, mother);
735 
736  return ok;
737 }
738 
740 {
741  if (mother->getNDaughters() < 2) return false;
742 
745 
746  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
747  const Particle* child = mother->getDaughter(ichild);
748 
749  if (child->getNDaughters() > 0) {
750  bool err = addChildofParticletoKFit(kf, child);
751  if (!err) return false;
752  } else {
753  if (child->getPValue() < 0) return false; // error matrix not valid
754 
755  kf.addParticle(child);
756  }
757  }
758 
759  // apply four momentum constraint
762 
763  int err = kf.doFit();
764 
765  if (err != 0) return false;
766 
767  bool ok = makeKFourCMother(kf, mother);
768 
769  return ok;
770 }
771 
773 {
774  if (mother->getNDaughters() < 2) return false;
775 
778 
779  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
780  const Particle* child = mother->getDaughter(ichild);
781 
782  if (child->getNDaughters() > 0) {
783  bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
784  std::abs(child->getPDGCode())) != m_massConstraintList.end();
785  std::vector<unsigned> childId;
786  bool err = addChildofParticletoMassKFit(kf, child, childId);
787  if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
788  if (!err) return false;
789  } else {
790  if (child->getPValue() < 0) return false; // error matrix not valid
791  kf.addParticle(child);
792  }
793  }
794 
795  // apply four momentum constraint
798 
799  int err = kf.doFit();
800 
801  if (err != 0) return false;
802 
803  bool ok = makeMassKFourCMother(kf, mother);
804 
805  return ok;
806 }
807 
809 {
812 
813  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
814  const Particle* child = mother->getDaughter(ichild);
815 
816  if (child->getPValue() < 0) return false; // error matrix not valid
817 
818  kf.addParticle(child);
819  }
820 
821  // apply four momentum constraint
824 
825  // apply recoil mass constraint
827 
828  int err = kf.doFit();
829 
830  if (err != 0) return false;
831 
832  bool ok = makeKRecoilMassMother(kf, mother);
833 
834  return ok;
835 }
836 
838  Particle* mother)
839 {
840  enum analysis::KFitError::ECode fitError;
841  fitError = kv.updateMother(mother);
842  if (fitError != analysis::KFitError::kNoError)
843  return false;
844  if (m_decayString.empty() && m_updateDaughters == true) {
845  // update daughter momenta as well
846  // the order of daughters in the *fitter is the same as in the mother Particle
847 
848  std::vector<Particle*> daughters = mother->getDaughters();
849 
850  unsigned track_count = kv.getTrackCount();
851  if (daughters.size() != track_count)
852  return false;
853 
854  for (unsigned iChild = 0; iChild < track_count; iChild++) {
855  double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
856  double dx = kv.getVertex().x() - kv.getTrackPosition(iChild).x();
857  double dy = kv.getVertex().y() - kv.getTrackPosition(iChild).y();
858 
859  ROOT::Math::PxPyPzEVector i4Vector(kv.getTrackMomentum(iChild).x() - a * dy,
860  kv.getTrackMomentum(iChild).y() + a * dx,
861  kv.getTrackMomentum(iChild).z(),
862  kv.getTrackMomentum(iChild).t());
863  daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
864 
865  daughters[iChild]->setVertex(
866  CLHEPToROOT::getXYZVector(kv.getTrackPosition(iChild)));
867  daughters[iChild]->setMomentumVertexErrorMatrix(
868  CLHEPToROOT::getTMatrixFSym(kv.getTrackError(iChild)));
869  }
870 
871  } else if (m_updateDaughters == true) { // if decayString is not empty
872  // first, update only the fit children
873  std::vector<const Particle*> fitChildren = m_decaydescriptor.getSelectionParticles(mother);
874 
875  unsigned track_count = kv.getTrackCount();
876  if (fitChildren.size() != track_count)
877  return false;
878 
879  for (unsigned iChild = 0; iChild < track_count; iChild++) {
880  auto daughter = const_cast<Particle*>(fitChildren[iChild]);
881 
882  double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughter->getCharge();
883  double dx = kv.getVertex().x() - kv.getTrackPosition(iChild).x();
884  double dy = kv.getVertex().y() - kv.getTrackPosition(iChild).y();
885 
886  ROOT::Math::PxPyPzEVector i4Vector(kv.getTrackMomentum(iChild).x() - a * dy,
887  kv.getTrackMomentum(iChild).y() + a * dx,
888  kv.getTrackMomentum(iChild).z(),
889  kv.getTrackMomentum(iChild).t());
890  daughter->set4VectorDividingByMomentumScaling(i4Vector);
891 
892  daughter->setVertex(CLHEPToROOT::getXYZVector(kv.getTrackPosition(iChild)));
893  daughter->setMomentumVertexErrorMatrix(CLHEPToROOT::getTMatrixFSym(kv.getTrackError(iChild)));
894  }
895 
896  // then, update other particles that have a fit-child in decay
897  std::function<bool(Particle*)> funcUpdateMomentum =
898  [&funcUpdateMomentum, fitChildren](Particle * part) {
899 
900  if (part->getNDaughters() == 0) {
901  // check if part is included in fitChildren
902  if (std::find(fitChildren.begin(), fitChildren.end(), part) != fitChildren.end())
903  return true;
904  else
905  return false;
906  }
907 
908  bool includeFitChildren = false;
909 
910  // Update daughters' momentum
911  for (auto daughter : part->getDaughters())
912  includeFitChildren = funcUpdateMomentum(daughter) || includeFitChildren;
913 
914  if (includeFitChildren) {
915  // Using updated daughters, update part's momentum
916  ROOT::Math::PxPyPzEVector sum4Vector;
917  for (auto daughter : part->getDaughters())
918  sum4Vector += daughter->get4Vector();
919 
920  part->set4VectorDividingByMomentumScaling(sum4Vector);
921  }
922 
923  return includeFitChildren;
924  };
925 
926  // Update all daughters
927  for (auto daughter : mother->getDaughters())
928  funcUpdateMomentum(daughter);
929 
930  }
931 
932 
933  return true;
934 }
935 
937  Particle* mother)
938 {
939  enum analysis::KFitError::ECode fitError;
940  fitError = kmv.updateMother(mother);
941  if (fitError != analysis::KFitError::kNoError)
942  return false;
943  if (m_decayString.empty() && m_updateDaughters == true) {
944  // update daughter momenta as well
945  // the order of daughters in the *fitter is the same as in the mother Particle
946 
947  std::vector<Particle*> daughters = mother->getDaughters();
948 
949  unsigned track_count = kmv.getTrackCount();
950  if (daughters.size() != track_count)
951  return false;
952 
953  for (unsigned iChild = 0; iChild < track_count; iChild++) {
954  double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
955  double dx = kmv.getVertex().x() - kmv.getTrackPosition(iChild).x();
956  double dy = kmv.getVertex().y() - kmv.getTrackPosition(iChild).y();
957 
958  ROOT::Math::PxPyPzEVector i4Vector(kmv.getTrackMomentum(iChild).x() - a * dy,
959  kmv.getTrackMomentum(iChild).y() + a * dx,
960  kmv.getTrackMomentum(iChild).z(),
961  kmv.getTrackMomentum(iChild).t());
962  daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
963 
964  daughters[iChild]->setVertex(
965  CLHEPToROOT::getXYZVector(kmv.getTrackPosition(iChild)));
966  daughters[iChild]->setMomentumVertexErrorMatrix(
967  CLHEPToROOT::getTMatrixFSym(kmv.getTrackError(iChild)));
968  }
969  }
970 
971  return true;
972 }
973 
975  Particle* mother)
976 {
977  enum analysis::KFitError::ECode fitError;
978  fitError = kmpv.updateMother(mother);
979  if (fitError != analysis::KFitError::kNoError) {
980  return false;
981  }
982 
983  if (m_decayString.empty() && m_updateDaughters == true) {
984  // update daughter momenta as well
985  // the order of daughters in the *fitter is the same as in the mother Particle
986 
987  std::vector<Particle*> daughters = mother->getDaughters();
988 
989  unsigned track_count = kmpv.getTrackCount();
990  if (daughters.size() != track_count)
991  return false;
992 
993  for (unsigned iChild = 0; iChild < track_count; iChild++) {
994  double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
995  double dx = kmpv.getVertex().x() - kmpv.getTrackPosition(iChild).x();
996  double dy = kmpv.getVertex().y() - kmpv.getTrackPosition(iChild).y();
997 
998  ROOT::Math::PxPyPzEVector i4Vector(kmpv.getTrackMomentum(iChild).x() - a * dy,
999  kmpv.getTrackMomentum(iChild).y() + a * dx,
1000  kmpv.getTrackMomentum(iChild).z(),
1001  kmpv.getTrackMomentum(iChild).t());
1002  daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1003 
1004  daughters[iChild]->setVertex(
1005  CLHEPToROOT::getXYZVector(kmpv.getTrackPosition(iChild)));
1006  daughters[iChild]->setMomentumVertexErrorMatrix(
1007  CLHEPToROOT::getTMatrixFSym(kmpv.getTrackError(iChild)));
1008  }
1009  }
1010 
1011  return true;
1012 }
1013 
1014 
1016  Particle* mother)
1017 {
1018  enum analysis::KFitError::ECode fitError;
1019  fitError = km.updateMother(mother);
1020  if (fitError != analysis::KFitError::kNoError)
1021  return false;
1022  if (m_decayString.empty() && m_updateDaughters == true) {
1023  // update daughter momenta as well
1024  // the order of daughters in the *fitter is the same as in the mother Particle
1025 
1026  std::vector<Particle*> daughters = mother->getDaughters();
1027 
1028  unsigned track_count = km.getTrackCount();
1029  if (daughters.size() != track_count)
1030  return false;
1031 
1032  for (unsigned iChild = 0; iChild < track_count; iChild++) {
1033  double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
1034  double dx = km.getVertex().x() - km.getTrackPosition(iChild).x();
1035  double dy = km.getVertex().y() - km.getTrackPosition(iChild).y();
1036 
1037  ROOT::Math::PxPyPzEVector i4Vector(km.getTrackMomentum(iChild).x() - a * dy,
1038  km.getTrackMomentum(iChild).y() + a * dx,
1039  km.getTrackMomentum(iChild).z(),
1040  km.getTrackMomentum(iChild).t());
1041  daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1042 
1043  daughters[iChild]->setVertex(
1044  CLHEPToROOT::getXYZVector(km.getTrackPosition(iChild)));
1045  daughters[iChild]->setMomentumVertexErrorMatrix(
1046  CLHEPToROOT::getTMatrixFSym(km.getTrackError(iChild)));
1047  }
1048  }
1049 
1050  return true;
1051 }
1052 
1053 
1054 
1056 {
1057  enum analysis::KFitError::ECode fitError;
1058  fitError = kf.updateMother(mother);
1059  if (fitError != analysis::KFitError::kNoError)
1060  return false;
1061  mother->addExtraInfo("FourCFitProb", kf.getCHIsq());
1062  mother->addExtraInfo("FourCFitChi2", kf.getNDF());
1063  if (m_decayString.empty() && m_updateDaughters == true) {
1064  // update daughter momenta as well
1065  // the order of daughters in the *fitter is the same as in the mother Particle
1066 
1067  std::vector<Particle*> daughters = mother->getDaughters();
1068 
1069  const unsigned nd = daughters.size();
1070  unsigned l = 0;
1071  std::vector<std::vector<unsigned>> pars;
1072  std::vector<Particle*> allparticles;
1073  for (unsigned ichild = 0; ichild < nd; ichild++) {
1074  const Particle* daughter = mother->getDaughter(ichild);
1075  std::vector<unsigned> pard;
1076  if (daughter->getNDaughters() > 0) {
1077  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
1078  pars.push_back(pard);
1079  allparticles.push_back(daughters[ichild]);
1080  } else {
1081  pard.push_back(l);
1082  pars.push_back(pard);
1083  allparticles.push_back(daughters[ichild]);
1084  l++;
1085  }
1086  }
1087 
1088  unsigned track_count = kf.getTrackCount();
1089  if (l != track_count)
1090  return false;
1091 
1092  for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
1093  ROOT::Math::PxPyPzEVector childMoms;
1094  ROOT::Math::XYZVector childPoss;
1095  TMatrixFSym childErrMatrixs(7);
1096  for (unsigned int iChild : pars[iDaug]) {
1097  childMoms = childMoms +
1098  CLHEPToROOT::getLorentzVector(
1099  kf.getTrackMomentum(iChild));
1100  childPoss = childPoss +
1101  CLHEPToROOT::getXYZVector(
1102  kf.getTrackPosition(iChild));
1103  TMatrixFSym childErrMatrix =
1104  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
1105  childErrMatrixs = childErrMatrixs + childErrMatrix;
1106  }
1107  allparticles[iDaug]->set4Vector(childMoms);
1108  allparticles[iDaug]->setVertex(childPoss);
1109  allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
1110  }
1111  }
1112 
1113  return true;
1114 }
1115 
1117 {
1118  enum analysis::KFitError::ECode fitError;
1119  fitError = kf.updateMother(mother);
1120  if (fitError != analysis::KFitError::kNoError)
1121  return false;
1122  mother->addExtraInfo("MassFourCFitProb", TMath::Prob(kf.getCHIsq(), kf.getNDF()));
1123  mother->addExtraInfo("MassFourCFitChi2", kf.getCHIsq());
1124  mother->addExtraInfo("MassFourCFitNDF", kf.getNDF());
1125  if (m_decayString.empty() && m_updateDaughters == true) {
1126  // update daughter momenta as well
1127  // the order of daughters in the *fitter is the same as in the mother Particle
1128 
1129  std::vector<Particle*> daughters = mother->getDaughters();
1130 
1131  const unsigned nd = daughters.size();
1132  unsigned l = 0;
1133  std::vector<std::vector<unsigned>> pars;
1134  std::vector<Particle*> allparticles;
1135  for (unsigned ichild = 0; ichild < nd; ichild++) {
1136  const Particle* daughter = mother->getDaughter(ichild);
1137  std::vector<unsigned> pard;
1138  if (daughter->getNDaughters() > 0) {
1139  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
1140  pars.push_back(pard);
1141  allparticles.push_back(daughters[ichild]);
1142  } else {
1143  pard.push_back(l);
1144  pars.push_back(pard);
1145  allparticles.push_back(daughters[ichild]);
1146  l++;
1147  }
1148  }
1149 
1150  unsigned track_count = kf.getTrackCount();
1151  if (l != track_count)
1152  return false;
1153 
1154  for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
1155  ROOT::Math::PxPyPzEVector childMoms;
1156  ROOT::Math::XYZVector childPoss;
1157  TMatrixFSym childErrMatrixs(7);
1158  for (unsigned int iChild : pars[iDaug]) {
1159  childMoms = childMoms +
1160  CLHEPToROOT::getLorentzVector(
1161  kf.getTrackMomentum(iChild));
1162  childPoss = childPoss +
1163  CLHEPToROOT::getXYZVector(
1164  kf.getTrackPosition(iChild));
1165  TMatrixFSym childErrMatrix =
1166  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
1167  childErrMatrixs = childErrMatrixs + childErrMatrix;
1168  }
1169  allparticles[iDaug]->set4Vector(childMoms);
1170  allparticles[iDaug]->setVertex(childPoss);
1171  allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
1172  }
1173  }
1174 
1175  return true;
1176 }
1177 
1179 {
1180  enum analysis::KFitError::ECode fitError;
1181  fitError = kf.updateMother(mother);
1182  if (fitError != analysis::KFitError::kNoError)
1183  return false;
1184  mother->addExtraInfo("RecoilMassFitProb", TMath::Prob(kf.getCHIsq(), kf.getNDF()));
1185  mother->addExtraInfo("RecoilMassFitChi2", kf.getCHIsq());
1186  mother->addExtraInfo("RecoilMassFitNDF", kf.getNDF());
1187  if (m_decayString.empty() && m_updateDaughters == true) {
1188  // update daughter momenta as well
1189  // the order of daughters in the *fitter is the same as in the mother Particle
1190 
1191  std::vector<Particle*> daughters = mother->getDaughters();
1192 
1193  unsigned track_count = kf.getTrackCount();
1194  if (daughters.size() != track_count)
1195  return false;
1196 
1197  for (unsigned iChild = 0; iChild < track_count; iChild++) {
1198  double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
1199  double dx = kf.getVertex().x() - kf.getTrackPosition(iChild).x();
1200  double dy = kf.getVertex().y() - kf.getTrackPosition(iChild).y();
1201 
1202  ROOT::Math::PxPyPzEVector i4Vector(kf.getTrackMomentum(iChild).x() - a * dy,
1203  kf.getTrackMomentum(iChild).y() + a * dx,
1204  kf.getTrackMomentum(iChild).z(),
1205  kf.getTrackMomentum(iChild).t());
1206  daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1207 
1208  daughters[iChild]->setVertex(
1209  CLHEPToROOT::getXYZVector(kf.getTrackPosition(iChild)));
1210  daughters[iChild]->setMomentumVertexErrorMatrix(
1211  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild)));
1212  }
1213  }
1214 
1215  return true;
1216 }
1217 
1218 
1219 void ParticleVertexFitterModule::updateMapOfTrackAndDaughter(unsigned& l, std::vector<std::vector<unsigned>>& pars,
1220  std::vector<unsigned>& parm, std::vector<Particle*>& allparticles, const Particle* daughter)
1221 {
1222  std::vector <Belle2::Particle*> childs = daughter->getDaughters();
1223  for (unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
1224  const Particle* child = daughter->getDaughter(ichild);
1225  std::vector<unsigned> pard;
1226  if (child->getNDaughters() > 0) {
1227  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
1228  parm.insert(parm.end(), pard.begin(), pard.end());
1229  pars.push_back(pard);
1230  allparticles.push_back(childs[ichild]);
1231  } else {
1232  pard.push_back(l);
1233  parm.push_back(l);
1234  pars.push_back(pard);
1235  allparticles.push_back(childs[ichild]);
1236  l++;
1237  }
1238  }
1239 }
1240 
1241 
1243 {
1244  if ((m_decayString.empty() ||
1245  (m_withConstraint == "" && m_fitType != "mass")) && mother->getNDaughters() < 2) return false;
1247  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "mother"
1248  || m_withConstraint == "iptubecut" || m_withConstraint == "btube")
1250 
1252  if (m_fitType == "mass") rf.setVertFit(false);
1253 
1254  if (m_decayString.empty()) {
1255  rf.addMother(mother);
1256  } else {
1257  std::vector<const Particle*> tracksVertex = m_decaydescriptor.getSelectionParticles(mother);
1258  std::vector<std::string> tracksName = m_decaydescriptor.getSelectionNames();
1259 
1260  if (allSelectedDaughters(mother, tracksVertex)) {
1261  for (auto& itrack : tracksVertex) {
1262  if (itrack != mother) rf.addTrack(itrack);
1263  }
1264  rf.setMother(mother);
1265  } else {
1266 
1268  bool mothSel = false;
1269  int nTrk = 0;
1270  for (unsigned itrack = 0; itrack < tracksVertex.size(); itrack++) {
1271  if (tracksVertex[itrack] != mother) {
1272  rsf.addTrack(tracksVertex[itrack]);
1273  B2DEBUG(1, "ParticleVertexFitterModule: Adding particle " << tracksName[itrack] << " to vertex fit ");
1274  nTrk++;
1275  }
1276  if (tracksVertex[itrack] == mother) mothSel = true;
1277  }
1278 
1279 
1280  // Fit one particle constrained to originate from the beam spot
1281  bool mothIPfit = false;
1282  if (tracksVertex.size() == 1 && mothSel == true && m_withConstraint != "" && nTrk == 0) {
1283  rsf.addTrack(tracksVertex[0]);
1284  if (tracksVertex[0] != mother)
1285  B2FATAL("ParticleVertexFitterModule: FATAL Error in IP constrained mother fit");
1286  nTrk++;
1287  mothIPfit = true;
1288  }
1289 
1290 
1291  ROOT::Math::XYZVector pos;
1292  TMatrixDSym RerrMatrix(3);
1293  int nvert = 0;
1294 
1295  // one track fit is not kinematic
1296  if (nTrk == 1) {
1298  for (auto& itrack : tracksVertex) {
1299  rsg.addTrack(itrack);
1300  nvert = rsg.fit("kalman");
1301  if (nvert > 0) {
1302  pos = rsg.getPos(0);
1303  RerrMatrix = rsg.getCov(0);
1304  double prob = rsg.getPValue(0);
1305  ROOT::Math::PxPyPzEVector mom(mother->get4Vector());
1306  TMatrixDSym errMatrix(7);
1307  for (int i = 0; i < 7; i++) {
1308  for (int j = 0; j < 7; j++) {
1309  if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1310  else {errMatrix[i][j] = 0;}
1311  }
1312  }
1313  if (mothIPfit) {
1314  mother->writeExtraInfo("prodVertX", pos.X());
1315  mother->writeExtraInfo("prodVertY", pos.Y());
1316  mother->writeExtraInfo("prodVertZ", pos.Z());
1317  mother->writeExtraInfo("prodVertSxx", RerrMatrix[0][0]);
1318  mother->writeExtraInfo("prodVertSxy", RerrMatrix[0][1]);
1319  mother->writeExtraInfo("prodVertSxz", RerrMatrix[0][2]);
1320  mother->writeExtraInfo("prodVertSyx", RerrMatrix[1][0]);
1321  mother->writeExtraInfo("prodVertSyy", RerrMatrix[1][1]);
1322  mother->writeExtraInfo("prodVertSyz", RerrMatrix[1][2]);
1323  mother->writeExtraInfo("prodVertSzx", RerrMatrix[2][0]);
1324  mother->writeExtraInfo("prodVertSzy", RerrMatrix[2][1]);
1325  mother->writeExtraInfo("prodVertSzz", RerrMatrix[2][2]);
1326  } else {
1327  mother->updateMomentum(mom, pos, errMatrix, prob);
1328  }
1329  return true;
1330  } else {return false;}
1331  }
1332  } else {
1333  nvert = rsf.fit();
1334  }
1335 
1336  if (nvert > 0) {
1337  pos = rsf.getPos();
1338  RerrMatrix = rsf.getVertexErrorMatrix();
1339  double prob = rsf.getPValue();
1340  ROOT::Math::PxPyPzEVector mom(mother->get4Vector());
1341  TMatrixDSym errMatrix(7);
1342  for (int i = 0; i < 7; i++) {
1343  for (int j = 0; j < 7; j++) {
1344  if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1345  else {errMatrix[i][j] = 0;}
1346  }
1347  }
1348  mother->updateMomentum(mom, pos, errMatrix, prob);
1349  } else {return false;}
1350 
1351 
1352  if (mothSel && nTrk > 1) {
1353  analysis::RaveSetup::getInstance()->setBeamSpot(B2Vector3D(pos.x(), pos.y(), pos.z()), RerrMatrix);
1354  rf.addMother(mother);
1355  int nKfit = rf.fit();
1356  rf.updateMother();
1358 
1359  if (nKfit > 0) {return true;}
1360  else return false;
1361  } else return true;
1362  }
1363  }
1364 
1365  bool okFT = false;
1366  if (m_fitType == "vertex") {
1367  okFT = true;
1368  int nVert = rf.fit();
1369  rf.updateMother();
1370  if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1371  if (nVert != 1) return false;
1372  }
1373  if (m_fitType == "mass") {
1374  // add protection
1375  okFT = true;
1376  rf.setMassConstFit(true);
1377  rf.setVertFit(false);
1378  int nVert = rf.fit();
1379  rf.updateMother();
1380  if (nVert != 1) return false;
1381  };
1382  if (m_fitType == "massvertex") {
1383  okFT = true;
1384  rf.setMassConstFit(true);
1385  int nVert = rf.fit();
1386  rf.updateMother();
1387  if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1388  if (nVert != 1) return false;
1389  };
1390  if (!okFT) {
1391  B2FATAL("fitType : " << m_fitType << " ***invalid fit type ");
1392  }
1393 
1394  return true;
1395 }
1396 
1398  const std::vector<const Particle*>& tracksVertex)
1399 {
1400  bool isAll = false;
1401  if (mother->getNDaughters() == 0) return false;
1402 
1403  int nNotIncluded = mother->getNDaughters();
1404 
1405  for (unsigned i = 0; i < mother->getNDaughters(); i++) {
1406  bool dauOk = false;
1407  for (auto& vi : tracksVertex) {
1408  if (vi == mother->getDaughter(i)) {
1409  nNotIncluded = nNotIncluded - 1;
1410  dauOk = true;
1411  }
1412  }
1413  if (!dauOk) {
1414  if (allSelectedDaughters(mother->getDaughter(i), tracksVertex)) nNotIncluded--;
1415  }
1416  }
1417  if (nNotIncluded == 0) isAll = true;
1418  return isAll;
1419 }
1420 
1422 {
1423  for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1424  const Particle* child = particle->getDaughter(ichild);
1425  if (child->getNDaughters() > 0) addChildofParticletoKFit(kf, child);
1426  else {
1427  if (child->getPValue() < 0) return false; // error matrix not valid
1428 
1429  kf.addParticle(child);
1430  }
1431  }
1432  return true;
1433 }
1434 
1436  std::vector<unsigned>& particleId)
1437 {
1438  for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1439  const Particle* child = particle->getDaughter(ichild);
1440  if (child->getNDaughters() > 0) {
1441  bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
1442  std::abs(child->getPDGCode())) != m_massConstraintList.end();
1443  std::vector<unsigned> childId;
1444  addChildofParticletoMassKFit(kf, child, childId);
1445  if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
1446  particleId.insert(particleId.end(), childId.begin(), childId.end());
1447  } else {
1448  if (child->getPValue() < 0) return false; // error matrix not valid
1449  kf.addParticle(child);
1450  particleId.push_back(kf.getTrackCount() - 1);
1451  }
1452  }
1453  return true;
1454 }
1455 
1457 {
1458  HepPoint3D pos(0.0, 0.0, 0.0);
1459  CLHEP::HepSymMatrix covMatrix(3, 0);
1460 
1461  for (int i = 0; i < 3; i++) {
1462  pos[i] = m_BeamSpotCenter(i);
1463  for (int j = 0; j < 3; j++) {
1464  covMatrix[i][j] = m_beamSpotCov(i, j);
1465  }
1466  }
1467 
1468  kv.setIpProfile(pos, covMatrix);
1469 }
1470 
1472 {
1473  CLHEP::HepSymMatrix err(7, 0);
1474 
1475  for (int i = 0; i < 3; i++) {
1476  for (int j = 0; j < 3; j++) {
1477  err[i + 4][j + 4] = m_beamSpotCov(i, j);
1478  }
1479  }
1480 
1481  PCmsLabTransform T;
1482  ROOT::Math::PxPyPzEVector iptube_mom = T.getBeamFourMomentum();
1483 
1484  kv.setIpTubeProfile(
1485  ROOTToCLHEP::getHepLorentzVector(iptube_mom),
1486  ROOTToCLHEP::getPoint3DFromB2Vector(m_BeamSpotCenter),
1487  err,
1488  0.);
1489 }
1490 
1492 {
1493  PCmsLabTransform T;
1494 
1495  B2Vector3D boost = T.getBoostVector();
1496  B2Vector3D boostDir = boost.Unit();
1497 
1498  TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1499  beamSpotCov(2, 2) = cut * cut;
1500  double thetab = boostDir.Theta();
1501  double phib = boostDir.Phi();
1502 
1503  double stb = TMath::Sin(thetab);
1504  double ctb = TMath::Cos(thetab);
1505  double spb = TMath::Sin(phib);
1506  double cpb = TMath::Cos(phib);
1507 
1508 
1509  TMatrix rz(3, 3); rz(2, 2) = 1;
1510  rz(0, 0) = cpb; rz(0, 1) = spb;
1511  rz(1, 0) = -1 * spb; rz(1, 1) = cpb;
1512 
1513  TMatrix ry(3, 3); ry(1, 1) = 1;
1514  ry(0, 0) = ctb; ry(0, 2) = -1 * stb;
1515  ry(2, 0) = stb; ry(2, 2) = ctb;
1516 
1517  TMatrix r(3, 3); r.Mult(rz, ry);
1518  TMatrix rt(3, 3); rt.Transpose(r);
1519 
1520  TMatrix TubePart(3, 3); TubePart.Mult(rt, beamSpotCov);
1521  TMatrix Tube(3, 3); Tube.Mult(TubePart, r);
1522 
1523  m_beamSpotCov(0, 0) = Tube(0, 0); m_beamSpotCov(0, 1) = Tube(0, 1); m_beamSpotCov(0, 2) = Tube(0, 2);
1524  m_beamSpotCov(1, 0) = Tube(1, 0); m_beamSpotCov(1, 1) = Tube(1, 1); m_beamSpotCov(1, 2) = Tube(1, 2);
1525  m_beamSpotCov(2, 0) = Tube(2, 0); m_beamSpotCov(2, 1) = Tube(2, 1); m_beamSpotCov(2, 2) = Tube(2, 2);
1526 }
1527 
1529 {
1530  TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1531  for (int i = 0; i < 3; i++)
1532  beamSpotCov(i, i) += width * width;
1533 
1534  m_beamSpotCov = beamSpotCov;
1535 }
1536 
1538 {
1539  double chi2TrackL = 0;
1540 
1541  for (int iTrack = 0; iTrack < kv.getTrackCount(); iTrack++) {
1542 
1543  analysis::KFitTrack trk_i = kv.getTrack(iTrack); // KFitTrack contains parameters before/after fit.
1544 
1545  TMatrixFSym err = CLHEPToROOT::getTMatrixFSym(trk_i.getError(analysis::KFitConst::kBeforeFit)); // px, py, pz, E, x, y, z
1546 
1547  B2Vector3D x_before = CLHEPToROOT::getXYZVector(trk_i.getPosition(analysis::KFitConst::kBeforeFit));
1548  B2Vector3D x_after = CLHEPToROOT::getXYZVector(trk_i.getPosition());
1549  B2Vector3D dPos = x_after - x_before;
1550 
1551  PCmsLabTransform T;
1552  B2Vector3D boost3 = T.getBoostVector().Unit();
1553  TVectorD boostD(0, 6, 0., 0., 0., 0., boost3.X(), boost3.Y(), boost3.Z(), "END");
1554 
1555  double dLBoost = dPos.Dot(boost3);
1556 
1557  chi2TrackL += TMath::Power(dLBoost, 2) / err.Similarity(boostD);
1558  }
1559  return chi2TrackL;
1560 }
DataType Phi() const
The azimuth angle.
Definition: B2Vector3.h:151
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
DataType Theta() const
The polar angle.
Definition: B2Vector3.h:153
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
Definition: B2Vector3.h:269
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
Definition: B2Vector3.h:290
void SetXYZ(DataType x, DataType y, DataType z)
set all coordinates using data type
Definition: B2Vector3.h:464
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:61
For each MCParticle with hits in the CDC, this class stores some summarising information on those hit...
Definition: Btube.h:27
TMatrixFSym getTubeMatrix() const
Returns Btube matrix.
Definition: Btube.h:78
Eigen::Matrix< double, 3, 1 > getTubeCenter() const
Returns Btube center.
Definition: Btube.h:64
int getPDGCode() const
PDG code.
Definition: Const.h:464
static const ParticleType pi0
neutral pion particle
Definition: Const.h:665
static const double speedOfLight
[cm/ns]
Definition: Const.h:686
static const ParticleType photon
photon particle
Definition: Const.h:664
bool init(const std::string &str)
Initialise the DecayDescriptor from given string.
std::vector< std::string > getSelectionNames()
Return list of human readable names of selected particles.
std::vector< const Particle * > getSelectionParticles(const Particle *particle)
Get a vector of pointers with selected daughters in the decay tree.
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
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.
B2Vector3D getBoostVector() const
Returns boost vector (beta=p/E)
TMatrixDSym m_beamSpotCov
Beam spot covariance matrix.
bool makeKMassMother(analysis::MassFitKFit &kv, Particle *p)
Update mother particle after mass fit using KFit.
std::string m_withConstraint
additional constraint on vertex
bool doKMassPointingVertexFit(Particle *p)
Mass-constrained vertex fit with additional pointing constraint using KFit.
void smearBeamSpot(double width)
smear beam spot covariance
bool makeMassKFourCMother(analysis::MassFourCFitKFit &kv, Particle *p)
Update mother particle after MassFourC fit using KFit.
bool doKVertexFit(Particle *p, bool ipProfileConstraint, bool ipTubeConstraint)
Unconstrained vertex fit using KFit.
bool doKMassVertexFit(Particle *p)
Mass-constrained vertex fit using KFit.
virtual void initialize() override
Initialize the Module.
void addIPTubeToKFit(analysis::VertexFitKFit &kv)
Adds IPTube constraint to the vertex fit using KFit.
bool addChildofParticletoMassKFit(analysis::MassFourCFitKFit &kf, const Particle *particle, std::vector< unsigned > &particleId)
Adds given particle's child to the MassFourCFitKFit.
virtual void event() override
Event processor.
bool makeKMassPointingVertexMother(analysis::MassPointingVertexFitKFit &kv, Particle *p)
Update mother particle after mass-constrained vertex fit with additional pointing constraint using KF...
bool m_updateDaughters
flag for daughters update
std::string m_decayString
daughter particles selection
std::vector< std::string > m_massConstraintListParticlename
Name of the particles to be mass constraint (massfourC)
std::string m_listName
particle list name
std::vector< int > m_massConstraintList
PDG codes of the particles to be mass constraint (massfourC)
bool m_hasCovMatrix
flag for mother covariance matrix (PseudoFitter)
bool makeKVertexMother(analysis::VertexFitKFit &kv, Particle *p)
Update mother particle after unconstrained vertex fit using KFit.
bool redoTwoPhotonDaughterMassFit(Particle *postFit, const Particle *preFit, const analysis::VertexFitKFit &kv)
Combines preFit particle and vertex information from vertex fit kv to create new postFit particle.
bool makeKMassVertexMother(analysis::MassVertexFitKFit &kv, Particle *p)
Update mother particle after mass-constrained vertex fit using KFit.
bool makeKFourCMother(analysis::FourCFitKFit &kv, Particle *p)
Update mother particle after FourC fit using KFit.
bool fillFitParticles(const Particle *mother, std::vector< const Particle * > &fitChildren, std::vector< const Particle * > &twoPhotonChildren)
Fills valid particle's children (with valid error matrix) in the vector of Particles that will enter ...
bool doKFourCFit(Particle *p)
FourC fit using KFit.
B2Vector3D m_BeamSpotCenter
Beam spot position.
virtual void beginRun() override
Called when entering a new run.
bool doKRecoilMassFit(Particle *p)
RecoilMass fit using KFit.
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.
DBObjPtr< BeamSpot > m_beamSpotDB
Beam spot database object.
bool doVertexFit(Particle *p)
Main steering routine.
std::string m_vertexFitter
Vertex Fitter name.
bool doRaveFit(Particle *mother)
Fit using Rave.
bool makeKRecoilMassMother(analysis::RecoilMassKFit &kf, Particle *p)
Update mother particle after RecoilMass fit using KFit.
bool doKMassFit(Particle *p)
Mass fit using KFit.
double m_recoilMass
recoil mass for constraint
double m_confidenceLevel
required fit confidence level
bool doKMassFourCFit(Particle *p)
MassFourC fit using KFit.
DecayDescriptor m_decaydescriptor
Decay descriptor of decays to look for.
bool fillNotFitParticles(const Particle *mother, std::vector< const Particle * > &notFitChildren, const std::vector< const Particle * > &fitChildren)
Fills valid particle's children (with valid error matrix) in the vector of Particles that will not en...
void findConstraintBoost(double cut)
calculate iptube constraint (quasi cylinder along boost direction) for RAVE fit
double m_Bfield
magnetic field from data base
bool addChildofParticletoKFit(analysis::FourCFitKFit &kv, const Particle *particle)
Adds given particle's child to the FourCFitKFit.
double getChi2TracksLBoost(const analysis::VertexFitKFit &kv)
calculate the chi2 using only lboost information of tracks
double m_smearing
smearing width applied to IP tube
std::string m_fitType
type of the kinematic fit
void addIPProfileToKFit(analysis::VertexFitKFit &kv)
Adds IPProfile constraint to the vertex fit using KFit.
StoreObjPtr< ParticleList > m_plist
particle list
bool allSelectedDaughters(const Particle *mother, const std::vector< const Particle * > &tracksVertex)
check if all the Daughters (o grand-daughters) are selected for the vertex fit
Class to store reconstructed particles.
Definition: Particle.h:75
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:426
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
Definition: Particle.cc:608
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:517
void updateMomentum(const ROOT::Math::PxPyPzEVector &p4, const ROOT::Math::XYZVector &vertex, const TMatrixFSym &errMatrix, double pValue)
Sets Lorentz vector, position, 7x7 error matrix and p-value.
Definition: Particle.h:358
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
Definition: Particle.cc:424
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:635
FourCFitKFit is a derived class from KFitBase to perform 4 momentum-constraint kinematical fit.
Definition: FourCFitKFit.h:30
enum KFitError::ECode setFourMomentum(const ROOT::Math::PxPyPzEVector &m)
Set an 4 Momentum for the FourC-constraint fit.
Definition: FourCFitKFit.cc:77
double getCHIsq(void) const override
Get a chi-square of the fit.
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
const KFitTrack getTrack(const int id) const
Get a specified track object.
Definition: KFitBase.cc:175
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
KFitTrack is a container of the track information (Lorentz vector, position, and error matrix),...
Definition: KFitTrack.h:38
const CLHEP::HepSymMatrix getError(const int flag=KFitConst::kAfterFit) const
Get an error matrix of the track.
Definition: KFitTrack.cc:172
const HepPoint3D getPosition(const int flag=KFitConst::kAfterFit) const
Get a position of the track.
Definition: KFitTrack.cc:164
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:44
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
Definition: MassFitKFit.cc:691
enum KFitError::ECode doFit(void)
Perform a mass-constraint fit.
Definition: MassFitKFit.cc:280
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-constraint fit.
Definition: MassFitKFit.cc:68
enum KFitError::ECode setVertexError(const CLHEP::HepSymMatrix &e)
Set an initial vertex error matrix for the mass-constraint fit.
Definition: MassFitKFit.cc:52
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
Definition: MassFitKFit.cc:137
MassFourCFitKFit is a derived class from KFitBase to perform mass and 4 momentum-constraint kinematic...
enum KFitError::ECode setFourMomentum(const ROOT::Math::PxPyPzEVector &m)
Set an 4 Momentum for the mass-four-constraint fit.
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 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.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
MassVertexFitKFit is a derived class from KFitBase to perform mass-vertex-constraint kinematical fit.
enum KFitError::ECode setInitialVertex(const HepPoint3D &v)
Set an initial vertex point for the mass-vertex constraint 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.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
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.
ROOT::Math::XYZVector getPos()
get the position of the fitted vertex.
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
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
void unsetBeamSpot()
unset beam spot constraint
Definition: RaveSetup.cc:80
static void initialize(int verbosity=1, double MagneticField=1.5)
Set everything up so everything needed for vertex fitting is there.
Definition: RaveSetup.cc:33
static RaveSetup * getInstance()
get the pointer to the instance to get/set any of options stored in RaveSetup
Definition: RaveSetup.h:43
void setBeamSpot(const B2Vector3D &beamSpot, const TMatrixDSym &beamSpotCov)
The beam spot position and covariance is known you can set it here so that and a vertex in the beam s...
Definition: RaveSetup.cc:72
void reset()
frees memory allocated by initialize().
Definition: RaveSetup.cc:58
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
B2Vector3D 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.
RecoilMassKFit is a derived class from KFitBase to perform a kinematical fit with a recoil mass const...
enum KFitError::ECode setFourMomentum(const ROOT::Math::PxPyPzEVector &m)
Set a recoil mass .
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode setRecoilMass(const double m)
Set an invariant mass for the four momentum-constraint fit.
enum KFitError::ECode doFit(void)
Perform a recoil-mass constraint fit.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
VertexFitKFit is a derived class from KFitBase to perform vertex-constraint kinematical fit.
Definition: VertexFitKFit.h:34
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 setInitialVertex(const HepPoint3D &v)
Set an initial vertex point for the vertex-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.
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
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
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.
static const int kBeforeFit
Input parameter to specify before-fit when setting/getting a track attribute.
Definition: KFitConst.h:35