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