Belle II Software development
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
40using namespace std;
41using namespace Belle2;
42
43
44//-----------------------------------------------------------------
45// Register module
46//-----------------------------------------------------------------
47
48REG_MODULE(ParticleVertexFitter);
49
50//-----------------------------------------------------------------
51// Implementation
52//-----------------------------------------------------------------
53
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}
81
83{
84 // Particle list with name m_listName has to exist
85 m_plist.isRequired(m_listName);
86
87 // magnetic field
88 m_Bfield = BFieldManager::getFieldInTesla(ROOT::Math::XYZVector(0, 0, 0)).Z();
89
90 // RAVE setup
91 if (m_vertexFitter == "Rave")
93
94 B2DEBUG(1, "ParticleVertexFitterModule : magnetic field = " << m_Bfield);
95
96
97 if (m_decayString != "")
99
100 B2INFO("ParticleVertexFitter: Performing " << m_fitType << " fit on " << m_listName << " using " << m_vertexFitter);
101 if (m_decayString != "")
102 B2INFO("ParticleVertexFitter: Using specified decay string: " << m_decayString);
103 if (m_withConstraint != "")
104 B2INFO("ParticleVertexFitter: Additional " << m_withConstraint << " will be applied");
105
106}
107
109{
110 //TODO: set magnetic field for each run
111 //m_Bfield = BFieldMap::Instance().getBField(B2Vector3D(0,0,0)).Z();
112 //TODO: set IP spot size for each run
113}
114
116{
117 if (m_vertexFitter == "Rave")
119
120 m_BeamSpotCenter = m_beamSpotDB->getIPPosition();
121 m_beamSpotCov.ResizeTo(3, 3);
122 TMatrixDSym beamSpotCov(3);
123 if (m_withConstraint == "ipprofile") m_beamSpotCov = m_beamSpotDB->getCovVertex();
124 if (m_withConstraint == "iptube") {
125 if (m_smearing > 0 && m_vertexFitter == "KFit") {
127 } else {
129 }
130 }
131 if (m_withConstraint == "iptubecut") { // for development purpose only
132 m_BeamSpotCenter = B2Vector3D(0.001, 0., .013);
134 }
135 if ((m_vertexFitter == "Rave") && (m_withConstraint == "ipprofile" || m_withConstraint == "iptube"
136 || m_withConstraint == "mother" || m_withConstraint == "iptubecut" || m_withConstraint == "btube"))
138
139 std::vector<unsigned int> toRemove;
140 unsigned int nParticles = m_plist->getListSize();
141
142 for (unsigned iPart = 0; iPart < nParticles; iPart++) {
143 Particle* particle = m_plist->getParticle(iPart);
144 m_hasCovMatrix = false;
145 if (m_updateDaughters == true) {
146 if (m_decayString.empty() || m_vertexFitter == "KFit")
148 else B2ERROR("Daughters update works only when all daughters are selected. Daughters will not be updated");
149 }
150
151 if (m_withConstraint == "mother") {
152 m_BeamSpotCenter = B2Vector3D(particle->getVertex().x(), particle->getVertex().y(), particle->getVertex().z());
153 m_beamSpotCov = particle->getVertexErrorMatrix();
154 }
155
156 TMatrixFSym mother_errMatrix(7);
157 mother_errMatrix = particle->getMomentumVertexErrorMatrix();
158 for (int k = 0; k < 7; k++) {
159 for (int j = 0; j < 7; j++) {
160 if (mother_errMatrix[k][j] > 0) {
161 m_hasCovMatrix = true;
162 }
163 }
164 }
165
166 bool hasTube = true;
167 if (m_withConstraint == "btube") {
168 Btube* Ver = particle->getRelatedTo<Btube>();
169 if (!Ver) {
170 hasTube = false;
171 toRemove.push_back(particle->getArrayIndex());
172 } else {
173 m_BeamSpotCenter.SetXYZ(Ver->getTubeCenter()(0, 0), Ver->getTubeCenter()(1, 0), Ver->getTubeCenter()(2, 0));
175 }
176 }
177 bool ok = false;
178 if (hasTube) {
179 ok = doVertexFit(particle);
180 }
181 if (!ok)
182 particle->setPValue(-1);
183 if (particle->getPValue() < m_confidenceLevel)
184 toRemove.push_back(particle->getArrayIndex());
185
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
309bool 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
358bool 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
472bool 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 if (mother->getV0()) {
501 HepPoint3D V0vertex_heppoint(mother->getV0()->getFittedVertexX(),
502 mother->getV0()->getFittedVertexY(),
503 mother->getV0()->getFittedVertexZ());
504 kv.setInitialVertex(V0vertex_heppoint);
505 }
506
507 for (auto& child : fitChildren)
508 kv.addParticle(child);
509
510 if (ipProfileConstraint)
512
513 if (ipTubeConstraint)
514 addIPTubeToKFit(kv);
515
516 // Perform vertex fit using only the particles with valid error matrices
517 int err = kv.doFit();
518 if (err != 0)
519 return false;
520
521 double chi2_track = getChi2TracksLBoost(kv);
522 unsigned track_count = kv.getTrackCount();
523 mother->writeExtraInfo("chiSquared_trackL", chi2_track);
524 mother->writeExtraInfo("kFit_nTracks", track_count);
525
526 bool ok = false;
527 if (twoPhotonChildren.size() == 0)
528 // in the case daughters do not include pi0 - this is it (fit done)
529 ok = makeKVertexMother(kv, mother);
530 else if (twoPhotonChildren.size() == 1) {
531 // there is a daughter reconstructed from two photons so without position information
532 // 1. determine vertex based on all other valid daughters
533 // 2. set position and error matrix of two-photon daughter to previously determined vertex
534 // 3. redo the fit using all particles (including two-photon particle this time)
535
536 const Particle* twoPhotonDaughter = twoPhotonChildren[0];
537 Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
538 ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
539 if (!ok)
540 return false;
541
542 // finally perform the fit using all daughter particles
545
546 for (auto& child : fitChildren)
547 kv2.addParticle(child);
548
549 kv2.addParticle(&fixedTwoPhotonDaughter);
550
551 if (ipProfileConstraint)
553
554 err = kv2.doFit();
555
556 if (err != 0)
557 return false;
558
559 ok = makeKVertexMother(kv2, mother);
560 }
561
562 // update 4-vector using not-fit-particles
563 ROOT::Math::PxPyPzEVector total4Vector(mother->get4Vector());
564 for (auto& child : notFitChildren)
565 total4Vector += child->get4Vector();
566 mother->set4Vector(total4Vector);
567
568 return ok;
569}
570
572{
573 if (mother->getNDaughters() < 2) return false;
574
575 std::vector<const Particle*> fitChildren;
576 std::vector<const Particle*> twoPhotonChildren;
577 bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
578
579 if (!validChildren)
580 return false;
581
582 if (twoPhotonChildren.size() > 1) {
583 B2FATAL("[ParticleVertexFitterModule::doKVertexFit] MassVertex fit using KFit does not support fit with multiple particles decaying to two photons like pi0 (yet).");
584 }
585
586 if (fitChildren.size() < 2) {
587 B2WARNING("[ParticleVertexFitterModule::doKVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
588 return false;
589 }
590
591 bool ok = false;
592 if (twoPhotonChildren.size() == 0) {
593 // Initialise the Fitter
596
597 if (mother->getV0()) {
598 HepPoint3D V0vertex_heppoint(mother->getV0()->getFittedVertexX(),
599 mother->getV0()->getFittedVertexY(),
600 mother->getV0()->getFittedVertexZ());
601 kmv.setInitialVertex(V0vertex_heppoint);
602 }
603
604 for (auto child : fitChildren)
605 kmv.addParticle(child);
606
607 kmv.setInvariantMass(mother->getPDGMass());
608 int err = kmv.doFit();
609 if (err != 0)
610 return false;
611
612 // in the case daughters do not include particles with two photon daughters like pi0 - this is it (fit done)
613 ok = makeKMassVertexMother(kmv, mother);
614 } else if (twoPhotonChildren.size() == 1) {
615 // there is a daughter reconstructed from two photons so without position information
616 // 1. determine vertex based on all other valid daughters
617 // 2. set position and error matrix of two-photon daughter to previously determined vertex
618 // 3. redo the fit using all particles (including two-photon particle this time)
619
622
623 for (auto child : fitChildren)
624 kv.addParticle(child);
625
626 // Perform vertex fit using only the particles with valid error matrices
627 int err = kv.doFit();
628 if (err != 0)
629 return false;
630
631 const Particle* twoPhotonDaughter = twoPhotonChildren[0];
632 Particle fixedTwoPhotonDaughter(twoPhotonDaughter->get4Vector(), twoPhotonDaughter->getPDGCode());
633 ok = redoTwoPhotonDaughterMassFit(&fixedTwoPhotonDaughter, twoPhotonDaughter, kv);
634 if (!ok)
635 return false;
636
637 // finally perform the fit using all daughter particles
640
641 for (auto child : fitChildren)
642 kmv2.addParticle(child);
643 kmv2.addParticle(&fixedTwoPhotonDaughter);
644
645 kmv2.setInvariantMass(mother->getPDGMass());
646 err = kmv2.doFit();
647
648 if (err != 0)
649 return false;
650
651 ok = makeKMassVertexMother(kmv2, mother);
652 }
653
654 return ok;
655
656}
657
659{
660 if (!(mother->hasExtraInfo("prodVertX") && mother->hasExtraInfo("prodVertY") && mother->hasExtraInfo("prodVertZ"))) {
661 return false;
662 }
663
664 if (mother->getNDaughters() < 2) return false;
665
666 std::vector<const Particle*> fitChildren;
667 std::vector<const Particle*> twoPhotonChildren;
668 bool validChildren = fillFitParticles(mother, fitChildren, twoPhotonChildren);
669
670 if (!validChildren)
671 return false;
672
673 if (twoPhotonChildren.size() > 0) {
674 B2FATAL("[ParticleVertexFitterModule::doKMassPointingVertexFit] MassPointingVertex fit using KFit does not support fit with two-photon daughters (yet).");
675 }
676
677 if (fitChildren.size() < 2) {
678 B2WARNING("[ParticleVertexFitterModule::doKMassPointingVertexFit] Number of particles with valid error matrix entering the vertex fit using KFit is less than 2.");
679 return false;
680 }
681
682 bool ok = false;
683 // Initialise the Fitter
686
687 for (auto child : fitChildren)
688 kmpv.addParticle(child);
689
690 kmpv.setInvariantMass(mother->getPDGMass());
691 HepPoint3D productionVertex(mother->getExtraInfo("prodVertX"),
692 mother->getExtraInfo("prodVertY"),
693 mother->getExtraInfo("prodVertZ"));
694 kmpv.setProductionVertex(productionVertex);
695 int err = kmpv.doFit();
696 if (err != 0) return false;
697
698 ok = makeKMassPointingVertexMother(kmpv, mother);
699
700 return ok;
701}
702
704{
705 if (mother->getNDaughters() < 2) return false;
706
709
710 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
711 const Particle* child = mother->getDaughter(ichild);
712
713 if (child->getPValue() < 0) return false; // error matrix not valid
714
715 km.addParticle(child);
716 }
717
718 // apply mass constraint
719 km.setInvariantMass(mother->getPDGMass());
720
721 int err = km.doFit();
722
723 if (err != 0) return false;
724
725 bool ok = makeKMassMother(km, mother);
726
727 return ok;
728}
729
731{
732 if (mother->getNDaughters() < 2) return false;
733
736
737 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
738 const Particle* child = mother->getDaughter(ichild);
739
740 if (child->getNDaughters() > 0) {
741 bool err = addChildofParticletoKFit(kf, child);
742 if (!err) return false;
743 } else {
744 if (child->getPValue() < 0) return false; // error matrix not valid
745
746 kf.addParticle(child);
747 }
748 }
749
750 // apply four momentum constraint
753
754 int err = kf.doFit();
755
756 if (err != 0) return false;
757
758 bool ok = makeKFourCMother(kf, mother);
759
760 return ok;
761}
762
764{
765 if (mother->getNDaughters() < 2) return false;
766
769
770 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
771 const Particle* child = mother->getDaughter(ichild);
772
773 if (child->getNDaughters() > 0) {
774 bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
775 std::abs(child->getPDGCode())) != m_massConstraintList.end();
776 std::vector<unsigned> childId;
777 bool err = addChildofParticletoMassKFit(kf, child, childId);
778 if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
779 if (!err) return false;
780 } else {
781 if (child->getPValue() < 0) return false; // error matrix not valid
782 kf.addParticle(child);
783 }
784 }
785
786 // apply four momentum constraint
789
790 int err = kf.doFit();
791
792 if (err != 0) return false;
793
794 bool ok = makeMassKFourCMother(kf, mother);
795
796 return ok;
797}
798
800{
803
804 for (unsigned ichild = 0; ichild < mother->getNDaughters(); ichild++) {
805 const Particle* child = mother->getDaughter(ichild);
806
807 if (child->getPValue() < 0) return false; // error matrix not valid
808
809 kf.addParticle(child);
810 }
811
812 // apply four momentum constraint
815
816 // apply recoil mass constraint
818
819 int err = kf.doFit();
820
821 if (err != 0) return false;
822
823 bool ok = makeKRecoilMassMother(kf, mother);
824
825 return ok;
826}
827
829 Particle* mother)
830{
831 enum analysis::KFitError::ECode fitError;
832 fitError = kv.updateMother(mother);
833 if (fitError != analysis::KFitError::kNoError)
834 return false;
835 if (m_decayString.empty() && m_updateDaughters == true) {
836 // update daughter momenta as well
837 // the order of daughters in the *fitter is the same as in the mother Particle
838
839 std::vector<Particle*> daughters = mother->getDaughters();
840
841 unsigned track_count = kv.getTrackCount();
842 if (daughters.size() != track_count)
843 return false;
844
845 for (unsigned iChild = 0; iChild < track_count; iChild++) {
846 double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
847 double dx = kv.getVertex().x() - kv.getTrackPosition(iChild).x();
848 double dy = kv.getVertex().y() - kv.getTrackPosition(iChild).y();
849
850 ROOT::Math::PxPyPzEVector i4Vector(kv.getTrackMomentum(iChild).x() - a * dy,
851 kv.getTrackMomentum(iChild).y() + a * dx,
852 kv.getTrackMomentum(iChild).z(),
853 kv.getTrackMomentum(iChild).t());
854 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
855
856 daughters[iChild]->setVertex(
857 CLHEPToROOT::getXYZVector(kv.getTrackPosition(iChild)));
858 daughters[iChild]->setMomentumVertexErrorMatrix(
859 CLHEPToROOT::getTMatrixFSym(kv.getTrackError(iChild)));
860 }
861
862 } else if (m_updateDaughters == true) { // if decayString is not empty
863 // first, update only the fit children
864 std::vector<const Particle*> fitChildren = m_decaydescriptor.getSelectionParticles(mother);
865
866 unsigned track_count = kv.getTrackCount();
867 if (fitChildren.size() != track_count)
868 return false;
869
870 for (unsigned iChild = 0; iChild < track_count; iChild++) {
871 auto daughter = const_cast<Particle*>(fitChildren[iChild]);
872
873 double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughter->getCharge();
874 double dx = kv.getVertex().x() - kv.getTrackPosition(iChild).x();
875 double dy = kv.getVertex().y() - kv.getTrackPosition(iChild).y();
876
877 ROOT::Math::PxPyPzEVector i4Vector(kv.getTrackMomentum(iChild).x() - a * dy,
878 kv.getTrackMomentum(iChild).y() + a * dx,
879 kv.getTrackMomentum(iChild).z(),
880 kv.getTrackMomentum(iChild).t());
881 daughter->set4VectorDividingByMomentumScaling(i4Vector);
882
883 daughter->setVertex(CLHEPToROOT::getXYZVector(kv.getTrackPosition(iChild)));
884 daughter->setMomentumVertexErrorMatrix(CLHEPToROOT::getTMatrixFSym(kv.getTrackError(iChild)));
885 }
886
887 // then, update other particles that have a fit-child in decay
888 std::function<bool(Particle*)> funcUpdateMomentum =
889 [&funcUpdateMomentum, fitChildren](Particle * part) {
890
891 if (part->getNDaughters() == 0) {
892 // check if part is included in fitChildren
893 if (std::find(fitChildren.begin(), fitChildren.end(), part) != fitChildren.end())
894 return true;
895 else
896 return false;
897 }
898
899 bool includeFitChildren = false;
900
901 // Update daughters' momentum
902 for (auto daughter : part->getDaughters())
903 includeFitChildren = funcUpdateMomentum(daughter) || includeFitChildren;
904
905 if (includeFitChildren) {
906 // Using updated daughters, update part's momentum
907 ROOT::Math::PxPyPzEVector sum4Vector;
908 for (auto daughter : part->getDaughters())
909 sum4Vector += daughter->get4Vector();
910
911 part->set4VectorDividingByMomentumScaling(sum4Vector);
912 }
913
914 return includeFitChildren;
915 };
916
917 // Update all daughters
918 for (auto daughter : mother->getDaughters())
919 funcUpdateMomentum(daughter);
920
921 }
922
923
924 return true;
925}
926
928 Particle* mother)
929{
930 enum analysis::KFitError::ECode fitError;
931 fitError = kmv.updateMother(mother);
932 if (fitError != analysis::KFitError::kNoError)
933 return false;
934 if (m_decayString.empty() && m_updateDaughters == true) {
935 // update daughter momenta as well
936 // the order of daughters in the *fitter is the same as in the mother Particle
937
938 std::vector<Particle*> daughters = mother->getDaughters();
939
940 unsigned track_count = kmv.getTrackCount();
941 if (daughters.size() != track_count)
942 return false;
943
944 for (unsigned iChild = 0; iChild < track_count; iChild++) {
945 double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
946 double dx = kmv.getVertex().x() - kmv.getTrackPosition(iChild).x();
947 double dy = kmv.getVertex().y() - kmv.getTrackPosition(iChild).y();
948
949 ROOT::Math::PxPyPzEVector i4Vector(kmv.getTrackMomentum(iChild).x() - a * dy,
950 kmv.getTrackMomentum(iChild).y() + a * dx,
951 kmv.getTrackMomentum(iChild).z(),
952 kmv.getTrackMomentum(iChild).t());
953 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
954
955 daughters[iChild]->setVertex(
956 CLHEPToROOT::getXYZVector(kmv.getTrackPosition(iChild)));
957 daughters[iChild]->setMomentumVertexErrorMatrix(
958 CLHEPToROOT::getTMatrixFSym(kmv.getTrackError(iChild)));
959 }
960 }
961
962 return true;
963}
964
966 Particle* mother)
967{
968 enum analysis::KFitError::ECode fitError;
969 fitError = kmpv.updateMother(mother);
970 if (fitError != analysis::KFitError::kNoError) {
971 return false;
972 }
973
974 if (m_decayString.empty() && m_updateDaughters == true) {
975 // update daughter momenta as well
976 // the order of daughters in the *fitter is the same as in the mother Particle
977
978 std::vector<Particle*> daughters = mother->getDaughters();
979
980 unsigned track_count = kmpv.getTrackCount();
981 if (daughters.size() != track_count)
982 return false;
983
984 for (unsigned iChild = 0; iChild < track_count; iChild++) {
985 double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
986 double dx = kmpv.getVertex().x() - kmpv.getTrackPosition(iChild).x();
987 double dy = kmpv.getVertex().y() - kmpv.getTrackPosition(iChild).y();
988
989 ROOT::Math::PxPyPzEVector i4Vector(kmpv.getTrackMomentum(iChild).x() - a * dy,
990 kmpv.getTrackMomentum(iChild).y() + a * dx,
991 kmpv.getTrackMomentum(iChild).z(),
992 kmpv.getTrackMomentum(iChild).t());
993 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
994
995 daughters[iChild]->setVertex(
996 CLHEPToROOT::getXYZVector(kmpv.getTrackPosition(iChild)));
997 daughters[iChild]->setMomentumVertexErrorMatrix(
998 CLHEPToROOT::getTMatrixFSym(kmpv.getTrackError(iChild)));
999 }
1000 }
1001
1002 return true;
1003}
1004
1005
1007 Particle* mother)
1008{
1009 enum analysis::KFitError::ECode fitError;
1010 fitError = km.updateMother(mother);
1011 if (fitError != analysis::KFitError::kNoError)
1012 return false;
1013 if (m_decayString.empty() && m_updateDaughters == true) {
1014 // update daughter momenta as well
1015 // the order of daughters in the *fitter is the same as in the mother Particle
1016
1017 std::vector<Particle*> daughters = mother->getDaughters();
1018
1019 unsigned track_count = km.getTrackCount();
1020 if (daughters.size() != track_count)
1021 return false;
1022
1023 for (unsigned iChild = 0; iChild < track_count; iChild++) {
1024 double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
1025 double dx = km.getVertex().x() - km.getTrackPosition(iChild).x();
1026 double dy = km.getVertex().y() - km.getTrackPosition(iChild).y();
1027
1028 ROOT::Math::PxPyPzEVector i4Vector(km.getTrackMomentum(iChild).x() - a * dy,
1029 km.getTrackMomentum(iChild).y() + a * dx,
1030 km.getTrackMomentum(iChild).z(),
1031 km.getTrackMomentum(iChild).t());
1032 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1033
1034 daughters[iChild]->setVertex(
1035 CLHEPToROOT::getXYZVector(km.getTrackPosition(iChild)));
1036 daughters[iChild]->setMomentumVertexErrorMatrix(
1037 CLHEPToROOT::getTMatrixFSym(km.getTrackError(iChild)));
1038 }
1039 }
1040
1041 return true;
1042}
1043
1044
1045
1047{
1048 enum analysis::KFitError::ECode fitError;
1049 fitError = kf.updateMother(mother);
1050 if (fitError != analysis::KFitError::kNoError)
1051 return false;
1052 mother->addExtraInfo("FourCFitProb", kf.getCHIsq());
1053 mother->addExtraInfo("FourCFitChi2", kf.getNDF());
1054 if (m_decayString.empty() && m_updateDaughters == true) {
1055 // update daughter momenta as well
1056 // the order of daughters in the *fitter is the same as in the mother Particle
1057
1058 std::vector<Particle*> daughters = mother->getDaughters();
1059
1060 const unsigned nd = daughters.size();
1061 unsigned l = 0;
1062 std::vector<std::vector<unsigned>> pars;
1063 std::vector<Particle*> allparticles;
1064 for (unsigned ichild = 0; ichild < nd; ichild++) {
1065 const Particle* daughter = mother->getDaughter(ichild);
1066 std::vector<unsigned> pard;
1067 if (daughter->getNDaughters() > 0) {
1068 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
1069 pars.push_back(pard);
1070 allparticles.push_back(daughters[ichild]);
1071 } else {
1072 pard.push_back(l);
1073 pars.push_back(pard);
1074 allparticles.push_back(daughters[ichild]);
1075 l++;
1076 }
1077 }
1078
1079 unsigned track_count = kf.getTrackCount();
1080 if (l != track_count)
1081 return false;
1082
1083 for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
1084 ROOT::Math::PxPyPzEVector childMoms;
1085 ROOT::Math::XYZVector childPoss;
1086 TMatrixFSym childErrMatrixs(7);
1087 for (unsigned int iChild : pars[iDaug]) {
1088 childMoms = childMoms +
1089 CLHEPToROOT::getLorentzVector(
1090 kf.getTrackMomentum(iChild));
1091 childPoss = childPoss +
1092 CLHEPToROOT::getXYZVector(
1093 kf.getTrackPosition(iChild));
1094 TMatrixFSym childErrMatrix =
1095 CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
1096 childErrMatrixs = childErrMatrixs + childErrMatrix;
1097 }
1098 allparticles[iDaug]->set4Vector(childMoms);
1099 allparticles[iDaug]->setVertex(childPoss);
1100 allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
1101 }
1102 }
1103
1104 return true;
1105}
1106
1108{
1109 enum analysis::KFitError::ECode fitError;
1110 fitError = kf.updateMother(mother);
1111 if (fitError != analysis::KFitError::kNoError)
1112 return false;
1113 mother->addExtraInfo("MassFourCFitProb", TMath::Prob(kf.getCHIsq(), kf.getNDF()));
1114 mother->addExtraInfo("MassFourCFitChi2", kf.getCHIsq());
1115 mother->addExtraInfo("MassFourCFitNDF", kf.getNDF());
1116 if (m_decayString.empty() && m_updateDaughters == true) {
1117 // update daughter momenta as well
1118 // the order of daughters in the *fitter is the same as in the mother Particle
1119
1120 std::vector<Particle*> daughters = mother->getDaughters();
1121
1122 const unsigned nd = daughters.size();
1123 unsigned l = 0;
1124 std::vector<std::vector<unsigned>> pars;
1125 std::vector<Particle*> allparticles;
1126 for (unsigned ichild = 0; ichild < nd; ichild++) {
1127 const Particle* daughter = mother->getDaughter(ichild);
1128 std::vector<unsigned> pard;
1129 if (daughter->getNDaughters() > 0) {
1130 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, daughter);
1131 pars.push_back(pard);
1132 allparticles.push_back(daughters[ichild]);
1133 } else {
1134 pard.push_back(l);
1135 pars.push_back(pard);
1136 allparticles.push_back(daughters[ichild]);
1137 l++;
1138 }
1139 }
1140
1141 unsigned track_count = kf.getTrackCount();
1142 if (l != track_count)
1143 return false;
1144
1145 for (unsigned iDaug = 0; iDaug < allparticles.size(); iDaug++) {
1146 ROOT::Math::PxPyPzEVector childMoms;
1147 ROOT::Math::XYZVector childPoss;
1148 TMatrixFSym childErrMatrixs(7);
1149 for (unsigned int iChild : pars[iDaug]) {
1150 childMoms = childMoms +
1151 CLHEPToROOT::getLorentzVector(
1152 kf.getTrackMomentum(iChild));
1153 childPoss = childPoss +
1154 CLHEPToROOT::getXYZVector(
1155 kf.getTrackPosition(iChild));
1156 TMatrixFSym childErrMatrix =
1157 CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild));
1158 childErrMatrixs = childErrMatrixs + childErrMatrix;
1159 }
1160 allparticles[iDaug]->set4Vector(childMoms);
1161 allparticles[iDaug]->setVertex(childPoss);
1162 allparticles[iDaug]->setMomentumVertexErrorMatrix(childErrMatrixs);
1163 }
1164 }
1165
1166 return true;
1167}
1168
1170{
1171 enum analysis::KFitError::ECode fitError;
1172 fitError = kf.updateMother(mother);
1173 if (fitError != analysis::KFitError::kNoError)
1174 return false;
1175 mother->addExtraInfo("RecoilMassFitProb", TMath::Prob(kf.getCHIsq(), kf.getNDF()));
1176 mother->addExtraInfo("RecoilMassFitChi2", kf.getCHIsq());
1177 mother->addExtraInfo("RecoilMassFitNDF", kf.getNDF());
1178 if (m_decayString.empty() && m_updateDaughters == true) {
1179 // update daughter momenta as well
1180 // the order of daughters in the *fitter is the same as in the mother Particle
1181
1182 std::vector<Particle*> daughters = mother->getDaughters();
1183
1184 unsigned track_count = kf.getTrackCount();
1185 if (daughters.size() != track_count)
1186 return false;
1187
1188 for (unsigned iChild = 0; iChild < track_count; iChild++) {
1189 double a = -1 * Belle2::Const::speedOfLight * 1e-4 * m_Bfield * daughters[iChild]->getCharge();
1190 double dx = kf.getVertex().x() - kf.getTrackPosition(iChild).x();
1191 double dy = kf.getVertex().y() - kf.getTrackPosition(iChild).y();
1192
1193 ROOT::Math::PxPyPzEVector i4Vector(kf.getTrackMomentum(iChild).x() - a * dy,
1194 kf.getTrackMomentum(iChild).y() + a * dx,
1195 kf.getTrackMomentum(iChild).z(),
1196 kf.getTrackMomentum(iChild).t());
1197 daughters[iChild]->set4VectorDividingByMomentumScaling(i4Vector);
1198
1199 daughters[iChild]->setVertex(
1200 CLHEPToROOT::getXYZVector(kf.getTrackPosition(iChild)));
1201 daughters[iChild]->setMomentumVertexErrorMatrix(
1202 CLHEPToROOT::getTMatrixFSym(kf.getTrackError(iChild)));
1203 }
1204 }
1205
1206 return true;
1207}
1208
1209
1210void ParticleVertexFitterModule::updateMapOfTrackAndDaughter(unsigned& l, std::vector<std::vector<unsigned>>& pars,
1211 std::vector<unsigned>& parm, std::vector<Particle*>& allparticles, const Particle* daughter)
1212{
1213 std::vector <Belle2::Particle*> daughters = daughter->getDaughters();
1214 for (unsigned ichild = 0; ichild < daughter->getNDaughters(); ichild++) {
1215 const Particle* child = daughter->getDaughter(ichild);
1216 std::vector<unsigned> pard;
1217 if (child->getNDaughters() > 0) {
1218 updateMapOfTrackAndDaughter(l, pars, pard, allparticles, child);
1219 parm.insert(parm.end(), pard.begin(), pard.end());
1220 pars.push_back(pard);
1221 allparticles.push_back(daughters[ichild]);
1222 } else {
1223 pard.push_back(l);
1224 parm.push_back(l);
1225 pars.push_back(pard);
1226 allparticles.push_back(daughters[ichild]);
1227 l++;
1228 }
1229 }
1230}
1231
1232
1234{
1235 if ((m_decayString.empty() ||
1236 (m_withConstraint == "" && m_fitType != "mass")) && mother->getNDaughters() < 2) return false;
1238 if (m_withConstraint == "ipprofile" || m_withConstraint == "iptube" || m_withConstraint == "mother"
1239 || m_withConstraint == "iptubecut" || m_withConstraint == "btube")
1241
1243 if (m_fitType == "mass") rf.setVertFit(false);
1244
1245 if (m_decayString.empty()) {
1246 rf.addMother(mother);
1247 } else {
1248 std::vector<const Particle*> tracksVertex = m_decaydescriptor.getSelectionParticles(mother);
1249 std::vector<std::string> tracksName = m_decaydescriptor.getSelectionNames();
1250
1251 if (allSelectedDaughters(mother, tracksVertex)) {
1252 for (auto& itrack : tracksVertex) {
1253 if (itrack != mother) rf.addTrack(itrack);
1254 }
1255 rf.setMother(mother);
1256 } else {
1257
1259 bool mothSel = false;
1260 int nTrk = 0;
1261 for (unsigned itrack = 0; itrack < tracksVertex.size(); itrack++) {
1262 if (tracksVertex[itrack] != mother) {
1263 rsf.addTrack(tracksVertex[itrack]);
1264 B2DEBUG(1, "ParticleVertexFitterModule: Adding particle " << tracksName[itrack] << " to vertex fit ");
1265 nTrk++;
1266 }
1267 if (tracksVertex[itrack] == mother) mothSel = true;
1268 }
1269
1270
1271 // Fit one particle constrained to originate from the beam spot
1272 bool mothIPfit = false;
1273 if (tracksVertex.size() == 1 && mothSel == true && m_withConstraint != "" && nTrk == 0) {
1274 rsf.addTrack(tracksVertex[0]);
1275 if (tracksVertex[0] != mother)
1276 B2FATAL("ParticleVertexFitterModule: FATAL Error in IP constrained mother fit");
1277 nTrk++;
1278 mothIPfit = true;
1279 }
1280
1281
1282 ROOT::Math::XYZVector pos;
1283 TMatrixDSym RerrMatrix(3);
1284 int nvert = 0;
1285
1286 // one track fit is not kinematic
1287 if (nTrk == 1) {
1289 for (auto& itrack : tracksVertex) {
1290 rsg.addTrack(itrack);
1291 nvert = rsg.fit("kalman");
1292 if (nvert > 0) {
1293 pos = rsg.getPos(0);
1294 RerrMatrix = rsg.getCov(0);
1295 double prob = rsg.getPValue(0);
1296 ROOT::Math::PxPyPzEVector mom(mother->get4Vector());
1297 TMatrixDSym errMatrix(7);
1298 for (int i = 0; i < 7; i++) {
1299 for (int j = 0; j < 7; j++) {
1300 if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1301 else {errMatrix[i][j] = 0;}
1302 }
1303 }
1304 if (mothIPfit) {
1305 mother->writeExtraInfo("prodVertX", pos.X());
1306 mother->writeExtraInfo("prodVertY", pos.Y());
1307 mother->writeExtraInfo("prodVertZ", pos.Z());
1308 mother->writeExtraInfo("prodVertSxx", RerrMatrix[0][0]);
1309 mother->writeExtraInfo("prodVertSxy", RerrMatrix[0][1]);
1310 mother->writeExtraInfo("prodVertSxz", RerrMatrix[0][2]);
1311 mother->writeExtraInfo("prodVertSyx", RerrMatrix[1][0]);
1312 mother->writeExtraInfo("prodVertSyy", RerrMatrix[1][1]);
1313 mother->writeExtraInfo("prodVertSyz", RerrMatrix[1][2]);
1314 mother->writeExtraInfo("prodVertSzx", RerrMatrix[2][0]);
1315 mother->writeExtraInfo("prodVertSzy", RerrMatrix[2][1]);
1316 mother->writeExtraInfo("prodVertSzz", RerrMatrix[2][2]);
1317 } else {
1318 mother->updateMomentum(mom, pos, errMatrix, prob);
1319 }
1320 return true;
1321 } else {return false;}
1322 }
1323 } else {
1324 nvert = rsf.fit();
1325 }
1326
1327 if (nvert > 0) {
1328 pos = rsf.getPos();
1329 RerrMatrix = rsf.getVertexErrorMatrix();
1330 double prob = rsf.getPValue();
1331 ROOT::Math::PxPyPzEVector mom(mother->get4Vector());
1332 TMatrixDSym errMatrix(7);
1333 for (int i = 0; i < 7; i++) {
1334 for (int j = 0; j < 7; j++) {
1335 if (i > 3 && j > 3) {errMatrix[i][j] = RerrMatrix[i - 4][j - 4];}
1336 else {errMatrix[i][j] = 0;}
1337 }
1338 }
1339 mother->updateMomentum(mom, pos, errMatrix, prob);
1340 } else {return false;}
1341
1342
1343 if (mothSel && nTrk > 1) {
1344 analysis::RaveSetup::getInstance()->setBeamSpot(B2Vector3D(pos.x(), pos.y(), pos.z()), RerrMatrix);
1345 rf.addMother(mother);
1346 int nKfit = rf.fit();
1347 rf.updateMother();
1349
1350 if (nKfit > 0) {return true;}
1351 else return false;
1352 } else return true;
1353 }
1354 }
1355
1356 bool okFT = false;
1357 if (m_fitType == "vertex") {
1358 okFT = true;
1359 int nVert = rf.fit();
1360 rf.updateMother();
1361 if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1362 if (nVert != 1) return false;
1363 }
1364 if (m_fitType == "mass") {
1365 // add protection
1366 okFT = true;
1367 rf.setMassConstFit(true);
1368 rf.setVertFit(false);
1369 int nVert = rf.fit();
1370 rf.updateMother();
1371 if (nVert != 1) return false;
1372 };
1373 if (m_fitType == "massvertex") {
1374 okFT = true;
1375 rf.setMassConstFit(true);
1376 int nVert = rf.fit();
1377 rf.updateMother();
1378 if (m_decayString.empty() && m_updateDaughters == true) rf.updateDaughters();
1379 if (nVert != 1) return false;
1380 };
1381 if (!okFT) {
1382 B2FATAL("fitType : " << m_fitType << " ***invalid fit type ");
1383 }
1384
1385 return true;
1386}
1387
1389 const std::vector<const Particle*>& tracksVertex)
1390{
1391 bool isAll = false;
1392 if (mother->getNDaughters() == 0) return false;
1393
1394 int nNotIncluded = mother->getNDaughters();
1395
1396 for (unsigned i = 0; i < mother->getNDaughters(); i++) {
1397 bool dauOk = false;
1398 for (auto& vi : tracksVertex) {
1399 if (vi == mother->getDaughter(i)) {
1400 nNotIncluded = nNotIncluded - 1;
1401 dauOk = true;
1402 }
1403 }
1404 if (!dauOk) {
1405 if (allSelectedDaughters(mother->getDaughter(i), tracksVertex)) nNotIncluded--;
1406 }
1407 }
1408 if (nNotIncluded == 0) isAll = true;
1409 return isAll;
1410}
1411
1413{
1414 for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1415 const Particle* child = particle->getDaughter(ichild);
1416 if (child->getNDaughters() > 0) addChildofParticletoKFit(kf, child);
1417 else {
1418 if (child->getPValue() < 0) return false; // error matrix not valid
1419
1420 kf.addParticle(child);
1421 }
1422 }
1423 return true;
1424}
1425
1427 std::vector<unsigned>& particleId)
1428{
1429 for (unsigned ichild = 0; ichild < particle->getNDaughters(); ichild++) {
1430 const Particle* child = particle->getDaughter(ichild);
1431 if (child->getNDaughters() > 0) {
1432 bool massconstraint = std::find(m_massConstraintList.begin(), m_massConstraintList.end(),
1433 std::abs(child->getPDGCode())) != m_massConstraintList.end();
1434 std::vector<unsigned> childId;
1435 addChildofParticletoMassKFit(kf, child, childId);
1436 if (massconstraint) kf.addMassConstraint(child->getPDGMass(), childId);
1437 particleId.insert(particleId.end(), childId.begin(), childId.end());
1438 } else {
1439 if (child->getPValue() < 0) return false; // error matrix not valid
1440 kf.addParticle(child);
1441 particleId.push_back(kf.getTrackCount() - 1);
1442 }
1443 }
1444 return true;
1445}
1446
1448{
1449 HepPoint3D pos(0.0, 0.0, 0.0);
1450 CLHEP::HepSymMatrix covMatrix(3, 0);
1451
1452 for (int i = 0; i < 3; i++) {
1453 pos[i] = m_BeamSpotCenter(i);
1454 for (int j = 0; j < 3; j++) {
1455 covMatrix[i][j] = m_beamSpotCov(i, j);
1456 }
1457 }
1458
1459 kv.setIpProfile(pos, covMatrix);
1460}
1461
1463{
1464 CLHEP::HepSymMatrix err(7, 0);
1465
1466 for (int i = 0; i < 3; i++) {
1467 for (int j = 0; j < 3; j++) {
1468 err[i + 4][j + 4] = m_beamSpotCov(i, j);
1469 }
1470 }
1471
1473 ROOT::Math::PxPyPzEVector iptube_mom = T.getBeamFourMomentum();
1474
1476 ROOTToCLHEP::getHepLorentzVector(iptube_mom),
1477 ROOTToCLHEP::getPoint3DFromB2Vector(m_BeamSpotCenter),
1478 err,
1479 0.);
1480}
1481
1483{
1485
1486 B2Vector3D boost = T.getBoostVector();
1487 B2Vector3D boostDir = boost.Unit();
1488
1489 TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1490 beamSpotCov(2, 2) = cut * cut;
1491 double thetab = boostDir.Theta();
1492 double phib = boostDir.Phi();
1493
1494 double stb = TMath::Sin(thetab);
1495 double ctb = TMath::Cos(thetab);
1496 double spb = TMath::Sin(phib);
1497 double cpb = TMath::Cos(phib);
1498
1499
1500 TMatrix rz(3, 3); rz(2, 2) = 1;
1501 rz(0, 0) = cpb; rz(0, 1) = spb;
1502 rz(1, 0) = -1 * spb; rz(1, 1) = cpb;
1503
1504 TMatrix ry(3, 3); ry(1, 1) = 1;
1505 ry(0, 0) = ctb; ry(0, 2) = -1 * stb;
1506 ry(2, 0) = stb; ry(2, 2) = ctb;
1507
1508 TMatrix r(3, 3); r.Mult(rz, ry);
1509 TMatrix rt(3, 3); rt.Transpose(r);
1510
1511 TMatrix TubePart(3, 3); TubePart.Mult(rt, beamSpotCov);
1512 TMatrix Tube(3, 3); Tube.Mult(TubePart, r);
1513
1514 m_beamSpotCov(0, 0) = Tube(0, 0); m_beamSpotCov(0, 1) = Tube(0, 1); m_beamSpotCov(0, 2) = Tube(0, 2);
1515 m_beamSpotCov(1, 0) = Tube(1, 0); m_beamSpotCov(1, 1) = Tube(1, 1); m_beamSpotCov(1, 2) = Tube(1, 2);
1516 m_beamSpotCov(2, 0) = Tube(2, 0); m_beamSpotCov(2, 1) = Tube(2, 1); m_beamSpotCov(2, 2) = Tube(2, 2);
1517}
1518
1520{
1521 TMatrixDSym beamSpotCov = m_beamSpotDB->getCovVertex();
1522 for (int i = 0; i < 3; i++)
1523 beamSpotCov(i, i) += width * width;
1524
1525 m_beamSpotCov = beamSpotCov;
1526}
1527
1529{
1530 double chi2TrackL = 0;
1531
1532 for (int iTrack = 0; iTrack < kv.getTrackCount(); iTrack++) {
1533
1534 analysis::KFitTrack trk_i = kv.getTrack(iTrack); // KFitTrack contains parameters before/after fit.
1535
1536 TMatrixFSym err = CLHEPToROOT::getTMatrixFSym(trk_i.getError(analysis::KFitConst::kBeforeFit)); // px, py, pz, E, x, y, z
1537
1538 B2Vector3D x_before = CLHEPToROOT::getXYZVector(trk_i.getPosition(analysis::KFitConst::kBeforeFit));
1539 B2Vector3D x_after = CLHEPToROOT::getXYZVector(trk_i.getPosition());
1540 B2Vector3D dPos = x_after - x_before;
1541
1543 B2Vector3D boost3 = T.getBoostVector().Unit();
1544 TVectorD boostD(0, 6, 0., 0., 0., 0., boost3.X(), boost3.Y(), boost3.Z(), "END");
1545
1546 double dLBoost = dPos.Dot(boost3);
1547
1548 chi2TrackL += TMath::Power(dLBoost, 2) / err.Similarity(boostD);
1549 }
1550 return chi2TrackL;
1551}
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
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
Definition: B2Vector3.h:290
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
Definition: B2Vector3.h:269
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:60
For each MCParticle with hits in the CDC, this class stores some summarising information on those hit...
Definition: Btube.h:22
TMatrixFSym getTubeMatrix() const
Returns Btube matrix.
Definition: Btube.h:73
Eigen::Matrix< double, 3, 1 > getTubeCenter() const
Returns Btube center.
Definition: Btube.h:59
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType pi0
neutral pion particle
Definition: Const.h:674
static const double speedOfLight
[cm/ns]
Definition: Const.h:695
static const ParticleType photon
photon particle
Definition: Const.h:673
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.
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.
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::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.
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:76
void writeExtraInfo(const std::string &name, const double value)
Sets the user defined extraInfo.
Definition: Particle.cc:1393
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
Definition: Particle.cc:1351
const V0 * getV0() const
Returns the pointer to the V0 object that was used to create this Particle (if ParticleType == c_V0).
Definition: Particle.cc:956
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:465
int getProperty() const
Returns particle property as a bit pattern The values are defined in the PropertyFlags enum and descr...
Definition: Particle.h:518
unsigned getNDaughters(void) const
Returns number of daughter particles.
Definition: Particle.h:747
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
Definition: Particle.cc:668
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
Definition: Particle.cc:635
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:567
void addExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1421
void set4Vector(const ROOT::Math::PxPyPzEVector &p4)
Sets Lorentz vector.
Definition: Particle.h:282
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:397
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
Definition: Particle.cc:451
@ c_IsUnspecified
Ordinary particles.
Definition: Particle.h:120
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:662
double getExtraInfo(const std::string &name) const
Return given value if set.
Definition: Particle.cc:1374
double getFittedVertexX() const
Get the X coordinate of the fitted vertex.
Definition: V0.h:75
double getFittedVertexY() const
Get the Y coordinate of the fitted vertex.
Definition: V0.h:78
double getFittedVertexZ() const
Get the Z coordinate of the fitted vertex.
Definition: V0.h:81
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
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:55
Abstract base class for different kinds of events.
STL namespace.
static const int kBeforeFit
Input parameter to specify before-fit when setting/getting a track attribute.
Definition: KFitConst.h:35