Belle II Software  light-2205-abys
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/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 
23 // utilities
24 #include <analysis/utility/CLHEPToROOT.h>
25 #include <analysis/utility/PCmsLabTransform.h>
26 #include <analysis/utility/ParticleCopy.h>
27 #include <analysis/utility/ROOTToCLHEP.h>
28 
29 // Magnetic field
30 #include <framework/geometry/BFieldManager.h>
31 
32 #include <TVector.h>
33 #include <TRotation.h>
34 #include <TMath.h>
35 
36 using namespace std;
37 using namespace Belle2;
38 
39 
40 //-----------------------------------------------------------------
41 // Register module
42 //-----------------------------------------------------------------
43 
44 REG_MODULE(ParticleVertexFitter);
45 
46 //-----------------------------------------------------------------
47 // Implementation
48 //-----------------------------------------------------------------
49 
50 ParticleVertexFitterModule::ParticleVertexFitterModule() : Module(),
51  m_Bfield(0)
52 {
53  // set module description (e.g. insert text)
54  setDescription("Vertex fitter for modular analysis");
56 
57  // Add parameters
58  addParam("listName", m_listName, "name of particle list", string(""));
59  addParam("confidenceLevel", m_confidenceLevel,
60  "Confidence level to accept the fit. Particle candidates with "
61  "p-value less than confidenceLevel are removed from the particle "
62  "list. If set to -1, all candidates are kept; if set to 0, "
63  "the candidates failing the fit are removed.",
64  0.001);
65  addParam("vertexFitter", m_vertexFitter, "KFit or Rave", string("KFit"));
66  addParam("fitType", m_fitType, "type of the kinematic fit (vertex, massvertex, mass)", string("vertex"));
67  addParam("withConstraint", m_withConstraint,
68  "additional constraint on vertex: ipprofile, iptube, mother, iptubecut, pointing, btube",
69  string(""));
70  addParam("decayString", m_decayString, "specifies which daughter particles are included in the kinematic fit", string(""));
71  addParam("updateDaughters", m_updateDaughters, "true: update the daughters after the vertex fit", false);
72  addParam("smearing", m_smearing, "smear IP tube width by given length", 0.002);
73  addParam("recoilMass", m_recoilMass, "recoil invariant mass (GeV)", 0.);
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 
81 {
82  // Particle list with name m_listName has to exist
83  m_plist.isRequired(m_listName);
84 
85  // magnetic field
87 
88  // RAVE setup
89  if (m_vertexFitter == "Rave")
91 
92  B2DEBUG(1, "ParticleVertexFitterModule : magnetic field = " << m_Bfield);
93 
94 
95  if (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 
114 {
115  //TODO: set magnetic field for each run
116  //m_Bfield = BFieldMap::Instance().getBField(B2Vector3D(0,0,0)).Z();
117  //TODO: set IP spot size for each run
118 }
119 
121 {
122  if (m_vertexFitter == "Rave")
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") {
132  } else {
134  }
135  }
136  if (m_withConstraint == "iptubecut") { // for development purpose only
137  m_BeamSpotCenter = B2Vector3D(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"))
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 = B2Vector3D(particle->getVertex().x(), particle->getVertex().y(), particle->getVertex().z());
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")
192 }
193 
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  if (m_decayString != "" and m_fitType != "vertex")
216  B2FATAL("ParticleVertexFitter: KFit does not support yet selection of daughters via decay string except for vertex fit!");
217 
218  // vertex fit
219  if (m_fitType == "vertex") {
220  if (m_withConstraint == "ipprofile") {
221  ok = doKVertexFit(mother, true, false);
222  } else if (m_withConstraint == "iptube") {
223  ok = doKVertexFit(mother, false, true);
224  } else {
225  ok = doKVertexFit(mother, false, false);
226  }
227  }
228 
229  // mass-constrained vertex fit
230  if (m_fitType == "massvertex") {
231  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
232  B2FATAL("ParticleVertexFitter: Invalid options - mass-constrained fit using KFit does not work with iptube or ipprofile constraint.");
233  } else if (m_withConstraint == "pointing") {
234  ok = doKMassPointingVertexFit(mother);
235  } else {
236  ok = doKMassVertexFit(mother);
237  }
238  }
239 
240  // mass fit
241  if (m_fitType == "mass") {
242  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
243  B2FATAL("ParticleVertexFitter: Invalid options - mass fit using KFit does not work with iptube or ipprofile constraint.");
244  } else {
245  ok = doKMassFit(mother);
246  }
247  }
248 
249  // four C fit
250  if (m_fitType == "fourC") {
251  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
252  B2FATAL("ParticleVertexFitter: Invalid options - four C fit using KFit does not work with iptube or ipprofile constraint.");
253  } else {
254  ok = doKFourCFit(mother);
255  }
256  }
257 
258  // four mass C fit
259  if (m_fitType == "massfourC") {
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 = doKMassFourCFit(mother);
264  }
265  }
266 
267  // recoil mass C fit
268  if (m_fitType == "recoilmass") {
269  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "iptubecut") {
270  B2FATAL("ParticleVertexFitter: Invalid options - recoil mass fit using KFit does not work with iptube or ipprofile constraint.");
271  } else {
272  ok = doKRecoilMassFit(mother);
273  }
274  }
275 
276  // invalid KFit fit type
277  if (m_fitType != "vertex"
278  && m_fitType != "massvertex"
279  && m_fitType != "mass"
280  && m_fitType != "fourC"
281  && m_fitType != "massfourC"
282  && m_fitType != "recoilmass")
283  B2FATAL("ParticleVertexFitter: " << m_fitType << " ***invalid fit type for the vertex fitter ");
284  }
285 
286  // fits using Rave
287  if (m_vertexFitter == "Rave") {
288  try {
289  ok = doRaveFit(mother);
290  } catch (const rave::CheckedFloatException&) {
291  B2ERROR("Invalid inputs (nan/inf)?");
292  ok = false;
293  }
294  }
295 
296  // invalid fitter
297  if (m_vertexFitter != "KFit" && m_vertexFitter != "Rave")
298  B2FATAL("ParticleVertexFitter: " << m_vertexFitter << " ***invalid vertex fitter ");
299 
300  if (!ok) return false;
301 
302  // steering ends here
303 
304  //if (mother->getPValue() < m_confidenceLevel) return false;
305  return true;
306 
307 }
308 
309 bool ParticleVertexFitterModule::fillFitParticles(const Particle* mother, std::vector<const Particle*>& fitChildren,
310  std::vector<const Particle*>& twoPhotonChildren)
311 {
312  if (m_decayString.empty()) {
313  // if decayString is empty, just use all primary daughters
314  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
315  const Particle* child = mother->getDaughter(ichild);
316  // This if allows to skip the daughters, which cannot be used in the fits, particularly K_L0 from KLM.
317  // Useful for fully-inclusive particles.
318  if (mother->getProperty() == Particle::PropertyFlags::c_IsUnspecified and child->getPValue() < 0) {
319  continue;
320  }
321  fitChildren.push_back(child);
322  }
323  } else {
324  fitChildren = m_decaydescriptor.getSelectionParticles(mother);
325  }
326 
327  auto itr = fitChildren.begin();
328  while (itr != fitChildren.end()) {
329  const Particle* child = *itr;
330 
331  if (child->getPValue() < 0) {
332  B2WARNING("Daughter with PDG code " << child->getPDGCode() << " does not have a valid error matrix.");
333  return false; // error matrix not valid
334  }
335  bool isTwoPhotonParticle = false;
336  if (m_hasCovMatrix == false) {
337  if (child->getPDGCode() == Const::pi0.getPDGCode() or child->getPDGCode() == 221) { // pi0 or eta
338  if (child->getNDaughters() == 2) {
339  if (child->getDaughter(0)->getPDGCode() == Const::photon.getPDGCode()
340  && child->getDaughter(1)->getPDGCode() == Const::photon.getPDGCode()) {
341  isTwoPhotonParticle = true;
342  }
343  }
344  }
345  }
346  if (isTwoPhotonParticle) {
347  // move children from fitChildren to twoPhotonChildren
348  twoPhotonChildren.push_back(child);
349  itr = fitChildren.erase(itr);
350  } else {
351  itr++;
352  }
353  }
354 
355  return true;
356 }
357 
358 bool ParticleVertexFitterModule::fillNotFitParticles(const Particle* mother, std::vector<const Particle*>& notFitChildren,
359  const std::vector<const Particle*>& fitChildren)
360 {
361  if (fitChildren.empty())
362  B2WARNING("[ParticleVertexFitterModule::fillNotFitParticles] fitChildren is empty! Please call fillFitParticles firstly");
363  if (!notFitChildren.empty())
364  B2WARNING("[ParticleVertexFitterModule::fillNotFitParticles] notFitChildren is NOT empty!"
365  << " The function should be called only once");
366 
367  if (m_decayString.empty())
368  // if decayString is empty, just use all primary daughters
369  return true;
370 
371  std::function<bool(const Particle*)> funcCheckInFit =
372  [&funcCheckInFit, &notFitChildren, fitChildren](const Particle * part) {
373 
374  // check if the given particle in fitChildren
375  // if it is included, return true
376  if (std::find(fitChildren.begin(), fitChildren.end(), part) != fitChildren.end())
377  return true;
378 
379  // if not, firstly check if particle has children
380  if (part->getNDaughters() == 0)
381  // if it has no children (=final-state-particle), return false
382  return false;
383 
384  // here, the given particle is not in fitChildren and has children
385  bool isAnyChildrenInFit = false;
386  vector<const Particle*> notFitChildren_tmp;
387  for (unsigned ichild = 0; ichild < part->getNDaughters(); ichild++) {
388  // call funcCheckInFit recursively for all children
389  const Particle* child = part->getDaughter(ichild);
390  bool isChildrenInFit = funcCheckInFit(child);
391  isAnyChildrenInFit = isChildrenInFit or isAnyChildrenInFit;
392 
393  // if the child is not in fitChildren, fill the child in a temporary vector
394  if (!isChildrenInFit)
395  notFitChildren_tmp.push_back(child);
396  }
397 
398  // if there are a sister in fitChildren, the children in the temporary vector will be filled in notFitChildren
399  if (isAnyChildrenInFit)
400  notFitChildren.insert(notFitChildren.end(), notFitChildren_tmp.begin(), notFitChildren_tmp.end());
401 
402  // if no children in fitChildren, the given particle should be filled instead of all children.
403 
404  return isAnyChildrenInFit;
405  };
406 
407 
408  // call funcCheckInFit for all primary children
409  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
410  const Particle* child = mother->getDaughter(ichild);
411  bool isGivenParticleOrAnyChildrenInFit = funcCheckInFit(child);
412  if (!isGivenParticleOrAnyChildrenInFit)
413  notFitChildren.push_back(child);
414  }
415 
416  return true;
417 }
418 
420  const analysis::VertexFitKFit& kv)
421 {
422  // TODO: something like setGammaError is necessary
423  // this is just workaround for the moment
424 
425  const Particle* g1Orig = preFit->getDaughter(0);
426  const Particle* g2Orig = preFit->getDaughter(1);
427  Particle g1Temp(g1Orig->get4Vector(), 22);
428  Particle g2Temp(g2Orig->get4Vector(), 22);
429 
430  TMatrixFSym g1ErrMatrix = g1Orig->getMomentumVertexErrorMatrix();
431  TMatrixFSym g2ErrMatrix = g2Orig->getMomentumVertexErrorMatrix();
432 
433  ROOT::Math::XYZVector pos(kv.getVertex().x(), kv.getVertex().y(), kv.getVertex().z());
434  CLHEP::HepSymMatrix posErrorMatrix = kv.getVertexError();
435 
436  TMatrixFSym errMatrix(3);
437  for (int i = 0; i < 3; i++)
438  for (int j = 0; j < 3; j++)
439  errMatrix(i, j) = posErrorMatrix[i][j];
440 
441  g1ErrMatrix.SetSub(4, errMatrix);
442  g2ErrMatrix.SetSub(4, errMatrix);
443 
444  g1Temp.updateMomentum(g1Orig->get4Vector(), pos, g1ErrMatrix, 1.0);
445  g2Temp.updateMomentum(g2Orig->get4Vector(), pos, g2ErrMatrix, 1.0);
446 
447  // perform the mass fit for the two-photon particle
450 
451  km.addParticle(&g1Temp);
452  km.addParticle(&g2Temp);
453 
454  km.setVertex(kv.getVertex());
456  km.setInvariantMass(preFit->getPDGMass());
457 
458  int err = km.doFit();
459  if (err != 0) {
460  return false;
461  }
462 
463  // The update of the daughters is disabled for this mass fit.
464  bool updateDaughters = m_updateDaughters;
465  m_updateDaughters = false;
466  bool ok = makeKMassMother(km, postFit);
467  m_updateDaughters = updateDaughters;
468 
469  return ok;
470 }
471 
472 bool ParticleVertexFitterModule::doKVertexFit(Particle* mother, bool ipProfileConstraint, bool ipTubeConstraint)
473 {
474  if ((mother->getNDaughters() < 2 && !ipTubeConstraint) || mother->getNDaughters() < 1) return false;
475 
476  std::vector<const Particle*> fitChildren;
477  std::vector<const Particle*> twoPhotonChildren;
478  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
479 
480  if (!validChildren)
481  return false;
482 
483  std::vector<const Particle*> notFitChildren;
484  fillNotFitParticles(mother, notFitChildren, fitChildren);
485 
486 
487  if (twoPhotonChildren.size() > 1) {
488  B2FATAL("[ParticleVertexFitterModule::doKVertexFit] Vertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
489  }
490 
491  if ((fitChildren.size() < 2 && !ipTubeConstraint) || fitChildren.size() < 1) {
492  B2WARNING("[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is too low.");
493  return false;
494  }
495 
496  // Initialise the Fitter
499 
500  for (auto& child : fitChildren)
501  kv.addParticle(child);
502 
503  if (ipProfileConstraint)
504  addIPProfileToKFit(kv);
505 
506  if (ipTubeConstraint)
507  addIPTubeToKFit(kv);
508 
509  // Perform vertex fit using only the particles with valid error matrices
510  int err = kv.doFit();
511  if (err != 0)
512  return false;
513 
514  double chi2_track = getChi2TracksLBoost(kv);
515  unsigned track_count = kv.getTrackCount();
516  mother->writeExtraInfo("chiSquared_trackL", chi2_track);
517  mother->writeExtraInfo("kFit_nTracks", track_count);
518 
519  bool ok = false;
520  if (twoPhotonChildren.size() == 0)
521  // in the case daughters do not include pi0 - this is it (fit done)
522  ok = makeKVertexMother(kv, mother);
523  else if (twoPhotonChildren.size() == 1) {
524  // there is a daughter reconstructed from two photons so without position information
525  // 1. determine vertex based on all other valid daughters
526  // 2. set position and error matrix of two-photon daughter to previously determined vertex
527  // 3. redo the fit using all particles (including two-photon particle this time)
528 
529  const Particle* twoPhotonDaughter = twoPhotonChildren[0];
530  Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
531  ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
532  if (!ok)
533  return false;
534 
535  // finally perform the fit using all daughter particles
538 
539  for (auto& child : fitChildren)
540  kv2.addParticle(child);
541 
542  kv2.addParticle(&fixedTwoPhotonDaughter);
543 
544  if (ipProfileConstraint)
545  addIPProfileToKFit(kv2);
546 
547  err = kv2.doFit();
548 
549  if (err != 0)
550  return false;
551 
552  ok = makeKVertexMother(kv2, mother);
553  }
554 
555  // update 4-vector using not-fit-particles
556  ROOT::Math::PxPyPzEVector total4Vector(mother->get4Vector());
557  for (auto& child : notFitChildren)
558  total4Vector += child->get4Vector();
559  mother->set4Vector(total4Vector);
560 
561  return ok;
562 }
563 
565 {
566  if (mother->getNDaughters() < 2) return false;
567 
568  std::vector<const Particle*> fitChildren;
569  std::vector<const Particle*> twoPhotonChildren;
570  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
571 
572  if (!validChildren)
573  return false;
574 
575  if (twoPhotonChildren.size() > 1) {
576  B2FATAL("[ParticleVertexFitterModule::doKVertexFit] MassVertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
577  }
578 
579  if (fitChildren.size() < 2) {
580  B2WARNING("[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
581  return false;
582  }
583 
584  bool ok = false;
585  if (twoPhotonChildren.size() == 0) {
586  // Initialise the Fitter
589 
590  for (auto child : fitChildren)
591  kmv.addParticle(child);
592 
593  kmv.setInvariantMass(mother->getPDGMass());
594  int err = kmv.doFit();
595  if (err != 0)
596  return false;
597 
598  // in the case daughters do not include particles with two photon daughters like pi0 - this is it (fit done)
599  ok = makeKMassVertexMother(kmv, mother);
600  } else if (twoPhotonChildren.size() == 1) {
601  // there is a daughter reconstructed from two photons so without position information
602  // 1. determine vertex based on all other valid daughters
603  // 2. set position and error matrix of two-photon daughter to previously determined vertex
604  // 3. redo the fit using all particles (including two-photon particle this time)
605 
608 
609  for (auto child : fitChildren)
610  kv.addParticle(child);
611 
612  // Perform vertex fit using only the particles with valid error matrices
613  int err = kv.doFit();
614  if (err != 0)
615  return false;
616 
617  const Particle* twoPhotonDaughter = twoPhotonChildren[0];
618  Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
619  ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
620  if (!ok)
621  return false;
622 
623  // finally perform the fit using all daughter particles
626 
627  for (auto child : fitChildren)
628  kmv2.addParticle(child);
629  kmv2.addParticle(&fixedTwoPhotonDaughter);
630 
631  kmv2.setInvariantMass(mother->getPDGMass());
632  err = kmv2.doFit();
633 
634  if (err != 0)
635  return false;
636 
637  ok = makeKMassVertexMother(kmv2, mother);
638  }
639 
640  return ok;
641 
642 }
643 
645 {
646  if (!(mother->hasExtraInfo("prodVertX") && mother->hasExtraInfo("prodVertY") && mother->hasExtraInfo("prodVertZ"))) {
647  return false;
648  }
649 
650  if (mother->getNDaughters() < 2) return false;
651 
652  std::vector<const Particle*> fitChildren;
653  std::vector<const Particle*> twoPhotonChildren;
654  bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
655 
656  if (!validChildren)
657  return false;
658 
659  if (twoPhotonChildren.size() > 0) {
660  B2FATAL("[ParticleVertexFitterModule::doKMassPointingVertexFit] MassPointingVertex fit using KFit does not support fit with two-photon daughters (yet).");
661  }
662 
663  if (fitChildren.size() < 2) {
664  B2WARNING("[ParticleVertexFitterModule::doKMassPointingVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
665  return false;
666  }
667 
668  bool ok = false;
669  // Initialise the Fitter
672 
673  for (auto child : fitChildren)
674  kmpv.addParticle(child);
675 
676  kmpv.setInvariantMass(mother->getPDGMass());
677  HepPoint3D productionVertex(mother->getExtraInfo("prodVertX"),
678  mother->getExtraInfo("prodVertY"),
679  mother->getExtraInfo("prodVertZ"));
680  kmpv.setProductionVertex(productionVertex);
681  int err = kmpv.doFit();
682  if (err != 0) return false;
683 
684  ok = makeKMassPointingVertexMother(kmpv, mother);
685 
686  return ok;
687 }
688 
690 {
691  if (mother->getNDaughters() < 2) return false;
692 
695 
696  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
697  const Particle* child = mother->getDaughter(ichild);
698 
699  if (child->getPValue() < 0) return false; // error matrix not valid
700 
701  km.addParticle(child);
702  }
703 
704  // apply mass constraint
705  km.setInvariantMass(mother->getPDGMass());
706 
707  int err = km.doFit();
708 
709  if (err != 0) return false;
710 
711  bool ok = makeKMassMother(km, mother);
712 
713  return ok;
714 }
715 
717 {
718  if (mother->getNDaughters() < 2) return false;
719 
722 
723  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
724  const Particle* child = mother->getDaughter(ichild);
725 
726  if (child->getNDaughters() > 0) {
727  bool err = addChildofParticletoKFit(kf, child);
728  if (!err) return false;
729  } else {
730  if (child->getPValue() < 0) return false; // error matrix not valid
731 
732  kf.addParticle(child);
733  }
734  }
735 
736  // apply four momentum constraint
739 
740  int err = kf.doFit();
741 
742  if (err != 0) return false;
743 
744  bool ok = makeKFourCMother(kf, mother);
745 
746  return ok;
747 }
748 
750 {
751  if (mother->getNDaughters() < 2) return false;
752 
755 
756  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
757  const Particle* child = mother->getDaughter(ichild);
758 
759  if (child->getNDaughters() > 0) {
760  bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
761  std::abs(child->getPDGCode())) != m_massConstraintList.end();
762  std::vector<unsigned> childId;
763  bool err = addChildofParticletoMassKFit(kf, child, childId);
764  if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
765  if (!err) return false;
766  } else {
767  if (child->getPValue() < 0) return false; // error matrix not valid
768  kf.addParticle(child);
769  }
770  }
771 
772  // apply four momentum constraint
775 
776  int err = kf.doFit();
777 
778  if (err != 0) return false;
779 
780  bool ok = makeMassKFourCMother(kf, mother);
781 
782  return ok;
783 }
784 
786 {
789 
790  for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
791  const Particle* child = mother->getDaughter(ichild);
792 
793  if (child->getPValue() < 0) return false; // error matrix not valid
794 
795  kf.addParticle(child);
796  }
797 
798  // apply four momentum constraint
801 
802  // apply recoil mass constraint
804 
805  int err = kf.doFit();
806 
807  if (err != 0) return false;
808 
809  bool ok = makeKRecoilMassMother(kf, mother);
810 
811  return ok;
812 }
813 
815  Particle* mother)
816 {
817  enum analysis::KFitError::ECode fitError;
818  fitError = kv.updateMother(mother);
819  if (fitError != analysis::KFitError::kNoError)
820  return false;
821  if (m_decayString.empty() && m_updateDaughters == true) {
822  // update daughter momenta as well
823  // the order of daughters in the *fitter is the same as in the mother Particle
824 
825  std::vector<Particle*> daughters = mother->getDaughters();
826 
827  unsigned track_count = kv.getTrackCount();
828  if (daughters.size() != track_count)
829  return false;
830 
831  for (unsigned iChild = 0; iChild < track_count; iChild++) {
832  daughters[iChild]->set4Vector(
833  CLHEPToROOT::getLorentzVector(kv.getTrackMomentum(iChild)));
834  daughters[iChild]->setVertex(
835  CLHEPToROOT::getXYZVector(kv.getTrackPosition(iChild)));
836  daughters[iChild]->setMomentumVertexErrorMatrix(
837  CLHEPToROOT::getTMatrixFSym(kv.getTrackError(iChild)));
838  }
839  }
840 
841  return true;
842 }
843 
845  Particle* mother)
846 {
847  enum analysis::KFitError::ECode fitError;
848  fitError = kmv.updateMother(mother);
849  if (fitError != analysis::KFitError::kNoError)
850  return false;
851  if (m_decayString.empty() && m_updateDaughters == true) {
852  // update daughter momenta as well
853  // the order of daughters in the *fitter is the same as in the mother Particle
854 
855  std::vector<Particle*> daughters = mother->getDaughters();
856 
857  unsigned track_count = kmv.getTrackCount();
858  if (daughters.size() != track_count)
859  return false;
860 
861  for (unsigned iChild = 0; iChild < track_count; iChild++) {
862  daughters[iChild]->set4Vector(
863  CLHEPToROOT::getLorentzVector(kmv.getTrackMomentum(iChild)));
864  daughters[iChild]->setVertex(
865  CLHEPToROOT::getXYZVector(kmv.getTrackPosition(iChild)));
866  daughters[iChild]->setMomentumVertexErrorMatrix(
867  CLHEPToROOT::getTMatrixFSym(kmv.getTrackError(iChild)));
868  }
869  }
870 
871  return true;
872 }
873 
875  Particle* mother)
876 {
877  enum analysis::KFitError::ECode fitError;
878  fitError = kmpv.updateMother(mother);
879  if (fitError != analysis::KFitError::kNoError) {
880  return false;
881  }
882 
883  if (m_decayString.empty() && m_updateDaughters == true) {
884  // update daughter momenta as well
885  // the order of daughters in the *fitter is the same as in the mother Particle
886 
887  std::vector<Particle*> daughters = mother->getDaughters();
888 
889  unsigned track_count = kmpv.getTrackCount();
890  if (daughters.size() != track_count)
891  return false;
892 
893  for (unsigned iChild = 0; iChild < track_count; iChild++) {
894  daughters[iChild]->set4Vector(
895  CLHEPToROOT::getLorentzVector(kmpv.getTrackMomentum(iChild)));
896  daughters[iChild]->setVertex(
897  CLHEPToROOT::getXYZVector(kmpv.getTrackPosition(iChild)));
898  daughters[iChild]->setMomentumVertexErrorMatrix(
899  CLHEPToROOT::getTMatrixFSym(kmpv.getTrackError(iChild)));
900  }
901  }
902 
903  return true;
904 }
905 
906 
908  Particle* mother)
909 {
910  enum analysis::KFitError::ECode fitError;
911  fitError = km.updateMother(mother);
912  if (fitError != analysis::KFitError::kNoError)
913  return false;
914  if (m_decayString.empty() && m_updateDaughters == true) {
915  // update daughter momenta as well
916  // the order of daughters in the *fitter is the same as in the mother Particle
917 
918  std::vector<Particle*> daughters = mother->getDaughters();
919 
920  unsigned track_count = km.getTrackCount();
921  if (daughters.size() != track_count)
922  return false;
923 
924  for (unsigned iChild = 0; iChild < track_count; iChild++) {
925  daughters[iChild]->set4Vector(
926  CLHEPToROOT::getLorentzVector(km.getTrackMomentum(iChild)));
927  daughters[iChild]->setVertex(
928  CLHEPToROOT::getXYZVector(km.getTrackPosition(iChild)));
929  daughters[iChild]->setMomentumVertexErrorMatrix(
930  CLHEPToROOT::getTMatrixFSym(km.getTrackError(iChild)));
931  }
932  }
933 
934  return true;
935 }
936 
937 
938 
940 {
941  enum analysis::KFitError::ECode fitError;
942  fitError = kf.updateMother(mother);
943  if (fitError != analysis::KFitError::kNoError)
944  return false;
945  mother->addExtraInfo("FourCFitProb", kf.getCHIsq());
946  mother->addExtraInfo("FourCFitChi2", kf.getNDF());
947  if (m_decayString.empty() && m_updateDaughters == true) {
948  // update daughter momenta as well
949  // the order of daughters in the *fitter is the same as in the mother Particle
950 
951  std::vector<Particle*> daughters = mother->getDaughters();
952 
953  const unsigned nd = daughters.size();
954  unsigned l = 0;
955  std::vector<std::vector<unsigned>> pars;
956  std::vector<Particle*> allparticles;
957  for (unsigned ichild = 0; ichild < nd; ichild++) {
958  const Particle* daughter = mother->getDaughter(ichild);
959  std::vector<unsigned> pard;
960  if (daughter->getNDaughters() > 0) {
961  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
962  pars.push_back(pard);
963  allparticles.push_back(daughters[ichild]);
964  } else {
965  pard.push_back(l);
966  pars.push_back(pard);
967  allparticles.push_back(daughters[ichild]);
968  l++;
969  }
970  }
971 
972  unsigned track_count = kf.getTrackCount();
973  if (l != track_count)
974  return false;
975 
976  for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
977  ROOT::Math::PxPyPzEVector childMoms;
978  ROOT::Math::XYZVector childPoss;
979  TMatrixFSym childErrMatrixs(7);
980  for (unsigned int iChild : pars[iDaug]) {
981  childMoms = childMoms +
982  CLHEPToROOT::getLorentzVector(
983  kf.getTrackMomentum(iChild));
984  childPoss = childPoss +
985  CLHEPToROOT::getXYZVector(
986  kf.getTrackPosition(iChild));
987  TMatrixFSym childErrMatrix =
988  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
989  childErrMatrixs = childErrMatrixs + childErrMatrix;
990  }
991  allparticles[iDaug]->set4Vector(childMoms);
992  allparticles[iDaug]->setVertex(childPoss);
993  allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
994  }
995  }
996 
997  return true;
998 }
999 
1001 {
1002  enum analysis::KFitError::ECode fitError;
1003  fitError = kf.updateMother(mother);
1004  if (fitError != analysis::KFitError::kNoError)
1005  return false;
1006  mother->addExtraInfo("MassFourCFitProb", TMath::Prob(kf.getCHIsq(), kf.getNDF()));
1007  mother->addExtraInfo("MassFourCFitChi2", kf.getCHIsq());
1008  mother->addExtraInfo("MassFourCFitNDF", kf.getNDF());
1009  if (m_decayString.empty() && m_updateDaughters == true) {
1010  // update daughter momenta as well
1011  // the order of daughters in the *fitter is the same as in the mother Particle
1012 
1013  std::vector<Particle*> daughters = mother->getDaughters();
1014 
1015  const unsigned nd = daughters.size();
1016  unsigned l = 0;
1017  std::vector<std::vector<unsigned>> pars;
1018  std::vector<Particle*> allparticles;
1019  for (unsigned ichild = 0; ichild < nd; ichild++) {
1020  const Particle* daughter = mother->getDaughter(ichild);
1021  std::vector<unsigned> pard;
1022  if (daughter->getNDaughters() > 0) {
1023  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
1024  pars.push_back(pard);
1025  allparticles.push_back(daughters[ichild]);
1026  } else {
1027  pard.push_back(l);
1028  pars.push_back(pard);
1029  allparticles.push_back(daughters[ichild]);
1030  l++;
1031  }
1032  }
1033 
1034  unsigned track_count = kf.getTrackCount();
1035  if (l != track_count)
1036  return false;
1037 
1038  for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
1039  ROOT::Math::PxPyPzEVector childMoms;
1040  ROOT::Math::XYZVector childPoss;
1041  TMatrixFSym childErrMatrixs(7);
1042  for (unsigned int iChild : pars[iDaug]) {
1043  childMoms = childMoms +
1044  CLHEPToROOT::getLorentzVector(
1045  kf.getTrackMomentum(iChild));
1046  childPoss = childPoss +
1047  CLHEPToROOT::getXYZVector(
1048  kf.getTrackPosition(iChild));
1049  TMatrixFSym childErrMatrix =
1050  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
1051  childErrMatrixs = childErrMatrixs + childErrMatrix;
1052  }
1053  allparticles[iDaug]->set4Vector(childMoms);
1054  allparticles[iDaug]->setVertex(childPoss);
1055  allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
1056  }
1057  }
1058 
1059  return true;
1060 }
1061 
1063 {
1064  enum analysis::KFitError::ECode fitError;
1065  fitError = kf.updateMother(mother);
1066  if (fitError != analysis::KFitError::kNoError)
1067  return false;
1068  mother->addExtraInfo("RecoilMassFitProb", TMath::Prob(kf.getCHIsq(), kf.getNDF()));
1069  mother->addExtraInfo("RecoilMassFitChi2", kf.getCHIsq());
1070  mother->addExtraInfo("RecoilMassFitNDF", kf.getNDF());
1071  if (m_decayString.empty() && m_updateDaughters == true) {
1072  // update daughter momenta as well
1073  // the order of daughters in the *fitter is the same as in the mother Particle
1074 
1075  std::vector<Particle*> daughters = mother->getDaughters();
1076 
1077  unsigned track_count = kf.getTrackCount();
1078  if (daughters.size() != track_count)
1079  return false;
1080 
1081  for (unsigned iChild = 0; iChild < track_count; iChild++) {
1082  daughters[iChild]->set4Vector(
1083  CLHEPToROOT::getLorentzVector(kf.getTrackMomentum(iChild)));
1084  daughters[iChild]->setVertex(
1085  CLHEPToROOT::getXYZVector(kf.getTrackPosition(iChild)));
1086  daughters[iChild]->setMomentumVertexErrorMatrix(
1087  CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild)));
1088  }
1089  }
1090 
1091  return true;
1092 }
1093 
1094 
1095 void ParticleVertexFitterModule::updateMapOfTrackAndDaughter(unsigned& l, std::vector<std::vector<unsigned>>& pars,
1096  std::vector<unsigned>& parm, std::vector<Particle*>& allparticles, const Particle* daughter)
1097 {
1098  std::vector <Belle2::Particle*> childs = daughter->getDaughters();
1099  for (unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
1100  const Particle* child = daughter->getDaughter(ichild);
1101  std::vector<unsigned> pard;
1102  if (child->getNDaughters() > 0) {
1103  updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
1104  parm.insert(parm.end(), pard.begin(), pard.end());
1105  pars.push_back(pard);
1106  allparticles.push_back(childs[ichild]);
1107  } else {
1108  pard.push_back(l);
1109  parm.push_back(l);
1110  pars.push_back(pard);
1111  allparticles.push_back(childs[ichild]);
1112  l++;
1113  }
1114  }
1115 }
1116 
1117 
1119 {
1120  if ((m_decayString.empty() ||
1121  (m_withConstraint == "" && m_fitType != "mass")) && mother->getNDaughters() < 2) return false;
1123  if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "mother"
1124  || m_withConstraint == "iptubecut" || m_withConstraint == "btube")
1126 
1128  if (m_fitType == "mass") rf.setVertFit(false);
1129 
1130  if (m_decayString.empty()) {
1131  rf.addMother(mother);
1132  } else {
1133  std::vector<const Particle*> tracksVertex = m_decaydescriptor.getSelectionParticles(mother);
1134  std::vector<std::string> tracksName = m_decaydescriptor.getSelectionNames();
1135 
1136  if (allSelectedDaughters(mother, tracksVertex)) {
1137  for (auto& itrack : tracksVertex) {
1138  if (itrack != mother) rf.addTrack(itrack);
1139  }
1140  rf.setMother(mother);
1141  } else {
1142 
1144  bool mothSel = false;
1145  int nTrk = 0;
1146  for (unsigned itrack = 0; itrack < tracksVertex.size(); itrack++) {
1147  if (tracksVertex[itrack] != mother) {
1148  rsf.addTrack(tracksVertex[itrack]);
1149  B2DEBUG(1, "ParticleVertexFitterModule: Adding particle " << tracksName[itrack] << " to vertex fit ");
1150  nTrk++;
1151  }
1152  if (tracksVertex[itrack] == mother) mothSel = true;
1153  }
1154 
1155 
1156  // Fit one particle constrained to originate from the beam spot
1157  bool mothIPfit = false;
1158  if (tracksVertex.size() == 1 && mothSel == true && m_withConstraint != "" && nTrk == 0) {
1159  rsf.addTrack(tracksVertex[0]);
1160  if (tracksVertex[0] != mother)
1161  B2FATAL("ParticleVertexFitterModule: FATAL Error in IP constrained mother fit");
1162  nTrk++;
1163  mothIPfit = true;
1164  }
1165 
1166 
1167  ROOT::Math::XYZVector pos;
1168  TMatrixDSym RerrMatrix(3);
1169  int nvert = 0;
1170 
1171  // one track fit is not kinematic
1172  if (nTrk == 1) {
1174  for (auto& itrack : tracksVertex) {
1175  rsg.addTrack(itrack);
1176  nvert = rsg.fit("kalman");
1177  if (nvert > 0) {
1178  pos = rsg.getPos(0);
1179  RerrMatrix = rsg.getCov(0);
1180  double prob = rsg.getPValue(0);
1181  ROOT::Math::PxPyPzEVector mom(mother->get4Vector());
1182  TMatrixDSym errMatrix(7);
1183  for (int i = 0; i < 7; i++) {
1184  for (int j = 0; j < 7; j++) {
1185  if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1186  else {errMatrix[i][j] = 0;}
1187  }
1188  }
1189  if (mothIPfit) {
1190  mother->writeExtraInfo("prodVertX", pos.X());
1191  mother->writeExtraInfo("prodVertY", pos.Y());
1192  mother->writeExtraInfo("prodVertZ", pos.Z());
1193  mother->writeExtraInfo("prodVertSxx", RerrMatrix[0][0]);
1194  mother->writeExtraInfo("prodVertSxy", RerrMatrix[0][1]);
1195  mother->writeExtraInfo("prodVertSxz", RerrMatrix[0][2]);
1196  mother->writeExtraInfo("prodVertSyx", RerrMatrix[1][0]);
1197  mother->writeExtraInfo("prodVertSyy", RerrMatrix[1][1]);
1198  mother->writeExtraInfo("prodVertSyz", RerrMatrix[1][2]);
1199  mother->writeExtraInfo("prodVertSzx", RerrMatrix[2][0]);
1200  mother->writeExtraInfo("prodVertSzy", RerrMatrix[2][1]);
1201  mother->writeExtraInfo("prodVertSzz", RerrMatrix[2][2]);
1202  } else {
1203  mother->updateMomentum(mom, pos, errMatrix, prob);
1204  }
1205  return true;
1206  } else {return false;}
1207  }
1208  } else {
1209  nvert = rsf.fit();
1210  }
1211 
1212  if (nvert > 0) {
1213  pos = rsf.getPos();
1214  RerrMatrix = rsf.getVertexErrorMatrix();
1215  double prob = rsf.getPValue();
1216  ROOT::Math::PxPyPzEVector mom(mother->get4Vector());
1217  TMatrixDSym errMatrix(7);
1218  for (int i = 0; i < 7; i++) {
1219  for (int j = 0; j < 7; j++) {
1220  if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1221  else {errMatrix[i][j] = 0;}
1222  }
1223  }
1224  mother->updateMomentum(mom, pos, errMatrix, prob);
1225  } else {return false;}
1226 
1227 
1228  if (mothSel && nTrk > 1) {
1229  analysis::RaveSetup::getInstance()->setBeamSpot(B2Vector3D(pos.x(), pos.y(), pos.z()), RerrMatrix);
1230  rf.addMother(mother);
1231  int nKfit = rf.fit();
1232  rf.updateMother();
1234 
1235  if (nKfit > 0) {return true;}
1236  else return false;
1237  } else return true;
1238  }
1239  }
1240 
1241  bool okFT = false;
1242  if (m_fitType == "vertex") {
1243  okFT = true;
1244  int nVert = rf.fit();
1245  rf.updateMother();
1246  if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1247  if (nVert != 1) return false;
1248  }
1249  if (m_fitType == "mass") {
1250  // add protection
1251  okFT = true;
1252  rf.setMassConstFit(true);
1253  rf.setVertFit(false);
1254  int nVert = rf.fit();
1255  rf.updateMother();
1256  if (nVert != 1) return false;
1257  };
1258  if (m_fitType == "massvertex") {
1259  okFT = true;
1260  rf.setMassConstFit(true);
1261  int nVert = rf.fit();
1262  rf.updateMother();
1263  if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1264  if (nVert != 1) return false;
1265  };
1266  if (!okFT) {
1267  B2FATAL("fitType : " << m_fitType << " ***invalid fit type ");
1268  }
1269 
1270  return true;
1271 }
1272 
1274  const std::vector<const Particle*>& tracksVertex)
1275 {
1276  bool isAll = false;
1277  if (mother->getNDaughters() == 0) return false;
1278 
1279  int nNotIncluded = mother->getNDaughters();
1280 
1281  for (unsigned i = 0; i < mother->getNDaughters(); i++) {
1282  bool dauOk = false;
1283  for (auto& vi : tracksVertex) {
1284  if (vi == mother->getDaughter(i)) {
1285  nNotIncluded = nNotIncluded - 1;
1286  dauOk = true;
1287  }
1288  }
1289  if (!dauOk) {
1290  if (allSelectedDaughters(mother->getDaughter(i), tracksVertex)) nNotIncluded--;
1291  }
1292  }
1293  if (nNotIncluded == 0) isAll = true;
1294  return isAll;
1295 }
1296 
1298 {
1299  for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1300  const Particle* child = particle->getDaughter(ichild);
1301  if (child->getNDaughters() > 0) addChildofParticletoKFit(kf, child);
1302  else {
1303  if (child->getPValue() < 0) return false; // error matrix not valid
1304 
1305  kf.addParticle(child);
1306  }
1307  }
1308  return true;
1309 }
1310 
1312  std::vector<unsigned>& particleId)
1313 {
1314  for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1315  const Particle* child = particle->getDaughter(ichild);
1316  if (child->getNDaughters() > 0) {
1317  bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
1318  std::abs(child->getPDGCode())) != m_massConstraintList.end();
1319  std::vector<unsigned> childId;
1320  addChildofParticletoMassKFit(kf, child, childId);
1321  if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
1322  particleId.insert(particleId.end(), childId.begin(), childId.end());
1323  } else {
1324  if (child->getPValue() < 0) return false; // error matrix not valid
1325  kf.addParticle(child);
1326  particleId.push_back(kf.getTrackCount() - 1);
1327  }
1328  }
1329  return true;
1330 }
1331 
1333 {
1334  HepPoint3D pos(0.0, 0.0, 0.0);
1335  CLHEP::HepSymMatrix covMatrix(3, 0);
1336 
1337  for (int i = 0; i < 3; i++) {
1338  pos[i] = m_BeamSpotCenter(i);
1339  for (int j = 0; j < 3; j++) {
1340  covMatrix[i][j] = m_beamSpotCov(i, j);
1341  }
1342  }
1343 
1344  kv.setIpProfile(pos, covMatrix);
1345 }
1346 
1348 {
1349  CLHEP::HepSymMatrix err(7, 0);
1350 
1351  for (int i = 0; i < 3; i++) {
1352  for (int j = 0; j < 3; j++) {
1353  err[i + 4][j + 4] = m_beamSpotCov(i, j);
1354  }
1355  }
1356 
1357  PCmsLabTransform T;
1358  ROOT::Math::PxPyPzEVector iptube_mom = T.getBeamFourMomentum();
1359 
1360  kv.setIpTubeProfile(
1361  ROOTToCLHEP::getHepLorentzVector(iptube_mom),
1362  ROOTToCLHEP::getPoint3DFromB2Vector(m_BeamSpotCenter),
1363  err,
1364  0.);
1365 }
1366 
1368 {
1369  PCmsLabTransform T;
1370 
1371  B2Vector3D boost = T.getBoostVector();
1372  B2Vector3D boostDir = boost.Unit();
1373 
1374  TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1375  beamSpotCov(2, 2) = cut * cut;
1376  double thetab = boostDir.Theta();
1377  double phib = boostDir.Phi();
1378 
1379  double stb = TMath::Sin(thetab);
1380  double ctb = TMath::Cos(thetab);
1381  double spb = TMath::Sin(phib);
1382  double cpb = TMath::Cos(phib);
1383 
1384 
1385  TMatrix rz(3, 3); rz(2, 2) = 1;
1386  rz(0, 0) = cpb; rz(0, 1) = spb;
1387  rz(1, 0) = -1 * spb; rz(1, 1) = cpb;
1388 
1389  TMatrix ry(3, 3); ry(1, 1) = 1;
1390  ry(0, 0) = ctb; ry(0, 2) = -1 * stb;
1391  ry(2, 0) = stb; ry(2, 2) = ctb;
1392 
1393  TMatrix r(3, 3); r.Mult(rz, ry);
1394  TMatrix rt(3, 3); rt.Transpose(r);
1395 
1396  TMatrix TubePart(3, 3); TubePart.Mult(rt, beamSpotCov);
1397  TMatrix Tube(3, 3); Tube.Mult(TubePart, r);
1398 
1399  m_beamSpotCov(0, 0) = Tube(0, 0); m_beamSpotCov(0, 1) = Tube(0, 1); m_beamSpotCov(0, 2) = Tube(0, 2);
1400  m_beamSpotCov(1, 0) = Tube(1, 0); m_beamSpotCov(1, 1) = Tube(1, 1); m_beamSpotCov(1, 2) = Tube(1, 2);
1401  m_beamSpotCov(2, 0) = Tube(2, 0); m_beamSpotCov(2, 1) = Tube(2, 1); m_beamSpotCov(2, 2) = Tube(2, 2);
1402 }
1403 
1405 {
1406  TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1407  for (int i = 0; i < 3; i++)
1408  beamSpotCov(i, i) += width * width;
1409 
1410  m_beamSpotCov = beamSpotCov;
1411 }
1412 
1414 {
1415  double chi2TrackL = 0;
1416 
1417  for (int iTrack = 0; iTrack < kv.getTrackCount(); iTrack++) {
1418 
1419  analysis::KFitTrack trk_i = kv.getTrack(iTrack); // KFitTrack contains parameters before/after fit.
1420 
1421  TMatrixFSym err = CLHEPToROOT::getTMatrixFSym(trk_i.getError(analysis::KFitConst::kBeforeFit)); // px, py, pz, E, x, y, z
1422 
1423  B2Vector3D x_before = CLHEPToROOT::getXYZVector(trk_i.getPosition(analysis::KFitConst::kBeforeFit));
1424  B2Vector3D x_after = CLHEPToROOT::getXYZVector(trk_i.getPosition());
1425  B2Vector3D dPos = x_after - x_before;
1426 
1427  PCmsLabTransform T;
1428  B2Vector3D boost3 = T.getBoostVector().Unit();
1429  TVectorD boostD(0, 6, 0., 0., 0., 0., boost3.X(), boost3.Y(), boost3.Z(), "END");
1430 
1431  double dLBoost = dPos.Dot(boost3);
1432 
1433  chi2TrackL += TMath::Power(dLBoost, 2) / err.Similarity(boostD);
1434  }
1435  return chi2TrackL;
1436 }
DataType Phi() const
The azimuth angle.
Definition: B2Vector3.h:145
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:429
DataType Theta() const
The polar angle.
Definition: B2Vector3.h:147
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:425
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:427
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
Definition: B2Vector3.h:263
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
Definition: B2Vector3.h:284
void SetXYZ(DataType x, DataType y, DataType z)
set all coordinates using data type
Definition: B2Vector3.h:454
static B2Vector3D getFieldInTesla(const B2Vector3D &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:70
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:354
static const ParticleType pi0
neutral pion particle
Definition: Const.h:555
static const ParticleType photon
photon particle
Definition: Const.h:554
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-.
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:74
TMatrixFSym getVertexErrorMatrix() const
Returns the 3x3 position error sub-matrix.
Definition: Particle.cc:479
void writeExtraInfo(const std::string &name, const double value)
Sets the user defined extraInfo.
Definition: Particle.cc:1297
double getPValue() const
Returns chi^2 probability of fit if done or -1.
Definition: Particle.h:596
ROOT::Math::XYZVector getVertex() const
Returns vertex position (POCA for charged, IP for neutral FS particles)
Definition: Particle.h:560
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
Definition: Particle.cc:1255
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:403
int getProperty() const
Returns particle property as a bit pattern The values are defined in the PropertyFlags enum and descr...
Definition: Particle.h:447
unsigned getNDaughters(void) const
Returns number of daughter particles.
Definition: Particle.h:656
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
Definition: Particle.cc:660
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
Definition: Particle.cc:636
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:488
void addExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1325
void set4Vector(const ROOT::Math::PxPyPzEVector &p4)
Sets Lorentz vector.
Definition: Particle.h:276
void setPValue(double pValue)
Sets chi^2 probability of fit.
Definition: Particle.h:331
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:351
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
Definition: Particle.cc:452
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:654
double getExtraInfo(const std::string &name) const
Return given value if set.
Definition: Particle.cc:1278
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
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:76
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:43
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
Definition: MassFitKFit.cc:693
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 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.
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.
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.
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 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:502
void copyDaughters(Particle *mother)
Function copies all (grand-)^n-daughter particles of the argument mother Particle.
Definition: ParticleCopy.cc:51
Abstract base class for different kinds of events.
Definition: ClusterUtils.h:23
static const int kBeforeFit
Input parameter to specify before-fit when setting/getting a track attribute.
Definition: KFitConst.h:34