9#include <analysis/modules/TagVertex/TagVertexModule.h>
15#include <framework/gearbox/Unit.h>
16#include <framework/gearbox/Const.h>
17#include <framework/logging/Logger.h>
20#include <analysis/dataobjects/RestOfEvent.h>
21#include <analysis/dataobjects/FlavorTaggerInfo.h>
24#include <analysis/utility/PCmsLabTransform.h>
25#include <analysis/variables/TrackVariables.h>
26#include <analysis/utility/ParticleCopy.h>
27#include <analysis/utility/CLHEPToROOT.h>
28#include <analysis/utility/ROOTToCLHEP.h>
29#include <analysis/utility/DistanceTools.h>
30#include <analysis/utility/RotationTools.h>
33#include <analysis/VertexFitting/KFit/VertexFitKFit.h>
36#include <mdst/dataobjects/HitPatternVXD.h>
39#include <framework/geometry/BFieldManager.h>
43#include <Math/Vector4D.h>
49static const double realNaN = std::numeric_limits<double>::quiet_NaN();
52static const ROOT::Math::XYZVector vecNaN(realNaN, realNaN, realNaN);
55static const double arrayNaN[] = {
56 realNaN, realNaN, realNaN,
57 realNaN, realNaN, realNaN,
58 realNaN, realNaN, realNaN,
62static const TMatrixDSym matNaN(3, arrayNaN);
65using RotationTools::rotateTensor;
66using RotationTools::rotateTensorInv;
67using RotationTools::toSymMatrix;
68using RotationTools::toVec;
69using RotationTools::getUnitOrthogonal;
81 m_Bfield(0), m_fitTruthStatus(0), m_rollbackStatus(0), m_fitPval(0), m_mcTagLifeTime(-1), m_mcPDG(0), m_mcLifeTimeReco(-1),
82 m_deltaT(0), m_deltaTErr(0), m_mcDeltaTau(0), m_mcDeltaT(0),
83 m_FitType(0), m_tagVl(0),
84 m_truthTagVl(0), m_tagVlErr(0), m_tagVol(0), m_truthTagVol(0), m_tagVolErr(0), m_tagVNDF(0), m_tagVChi2(0), m_tagVChi2IP(0),
85 m_kFitReqReducedChi2(5),
95 "required confidence level of fit to keep particles in the list. Note that even with confidenceLevel == 0.0, errors during the fit might discard Particles in the list. confidenceLevel = -1 if an error occurs during the fit",
98 "'': no MC association. breco: use standard Breco MC association. internal: use internal MC association",
string(
"breco"));
100 "Choose the type of the constraint: noConstraint, IP (tag tracks constrained to be within the beam spot), tube (long tube along the BTag line of flight, only for fully reconstruced B rec), boost (long tube along the Upsilon(4S) boost direction), (breco)",
103 "Choose how to reconstruct the tracks on the tag side: standard, standard_PXD",
104 string(
"standard_PXD"));
108 "TRUE when requesting MC Information from the tracks performing the vertex fit",
false);
110 "Minimum number of PXD hits for a track to be used in the vertex fit", 0);
112 "Fitter used for the tag vertex fit: Rave or KFit",
string(
"KFit"));
114 "The required chi2/ndf to accept the kFit result, if it is higher, iteration procedure is applied", 5.0);
116 "Use the true track parameters in the vertex fit",
false);
118 "Use rolled back non-primary tracks",
false);
127 B2INFO(
"TagVertexModule : magnetic field = " <<
m_Bfield);
141 B2FATAL(
"TagVertexModule: invalid fitting algorithm (must be set to either Rave or KFit).");
143 B2FATAL(
"TagVertexModule: invalid fitting option (useRollBack and useTruthInFit cannot be simultaneously set to true).");
146 B2FATAL(
"TagVertexModule : the singleTrack option is temporarily broken.");
158 B2ERROR(
"TagVertexModule: ParticleList " <<
m_listName <<
" not found");
165 std::vector<unsigned int> toRemove;
167 for (
unsigned i = 0; i <
m_plist->getListSize(); ++i) {
177 toRemove.push_back(particle->getArrayIndex());
182 particle->addRelationTo(ver);
244 m_plist->removeParticles(toRemove);
271 B2ERROR(
"TagVertex: No magnetic field");
284 double bg = beta /
sqrt(1 - beta * beta);
289 double lB0 = tauB * bg * c;
300 B2ERROR(
"TagVertex: Invalid constraintType selected");
305 B2ERROR(
"TagVertex: No correct fit constraint");
314 double minPVal = (
m_fitAlgo !=
"KFit") ? 0.001 : 0.;
352 if (Breco->
getPValue() < 0.)
return make_pair(vecNaN, matNaN);
354 TMatrixDSym beamSpotCov(3);
364 TMatrixDSym PerrMatrix(7);
366 for (
int i = 0; i < 3; ++i) {
367 for (
int j = 0; j < 3; ++j) {
369 PerrMatrix(i, j) = (beamSpotCov(i, j) + TerrMatrix(i, j)) * pmag / xmag;
371 PerrMatrix(i, j) = TerrMatrix(i, j);
373 PerrMatrix(i + 4, j + 4) = TerrMatrix(i + 4, j + 4);
377 PerrMatrix(3, 3) = 0.;
385 if (BRecoRes->
getPValue() < 0)
return make_pair(vecNaN, matNaN);
394 ROOT::Math::PxPyPzMVector pBrecEstimate(BvertDiff.X(), BvertDiff.Y(), BvertDiff.Z(), Breco->
getPDGMass());
396 pBtagEstimate.SetPxPyPzE(-pBtagEstimate.px(), -pBtagEstimate.py(), -pBtagEstimate.pz(), pBtagEstimate.E());
400 TMatrixD TubeZ = rotateTensorInv(pBrecEstimate.Vect(), errFinal);
402 TubeZ(2, 2) = cut * cut;
403 TubeZ(2, 0) = 0; TubeZ(0, 2) = 0;
404 TubeZ(2, 1) = 0; TubeZ(1, 2) = 0;
408 TMatrixD Tube = rotateTensor(
B2Vector3D(pBtagEstimate.Vect()), TubeZ);
418 B2WARNING(
"In TagVertexModule::findConstraintBTube: cannot get a proper vertex for BReco. BTube constraint replaced by Boost.");
426 if (tubecreatorBCopy->
getPValue() < 0)
return make_pair(vecNaN, matNaN);
430 ROOT::Math::PxPyPzEVector pBrec = tubecreatorBCopy->
get4Vector();
441 pBtag.SetPxPyPzE(-pBtag.px(), -pBtag.py(), -pBtag.pz(), pBtag.E());
453 B2DEBUG(10,
"Brec direction after fit: " <<
printVector(
float(1. / tubecreatorBCopy->
getP()) * tubecreatorBCopy->
getMomentum()));
458 B2DEBUG(10,
"Brec PV covariance: " <<
printMatrix(pv));
459 B2DEBUG(10,
"BTag direction: " <<
printVector((1. / pBtag.P())*pBtag.Vect()));
463 TMatrixD longerror(3, 3); longerror(2, 2) = cut * cut;
467 TMatrixD longerrorRotated = rotateTensor(
B2Vector3D(pBtag.Vect()), longerror);
470 TMatrixD pvNew = TMatrixD(pv) + longerrorRotated;
473 ROOT::Math::XYZVector constraintCenter = tubecreatorBCopy->
getVertex();
484 B2DEBUG(10,
"IPTube covariance: " <<
printMatrix(pvNew));
495 return make_pair(constraintCenter, toSymMatrix(pvNew));
503 TMatrixD longerror(3, 3); longerror(2, 2) = cut * cut;
504 longerror(0, 0) = longerror(1, 1) = d * d;
508 TMatrixD longerrorRotated = rotateTensor(boostDir, longerror);
512 TMatrixD Tube = TMatrixD(beamSpotCov) + longerrorRotated;
517 return make_pair(constraintCenter, toSymMatrix(Tube));
521static double getProperLifeTime(
const MCParticle* mc)
523 double beta = mc->getMomentum().R() / mc->getEnergy();
524 return 1e3 * mc->getLifetime() *
sqrt(1 - pow(beta, 2));
530 vector<const MCParticle*> mcBs;
532 if (abs(mc.getPDG()) == abs(Breco->
getPDGCode()))
536 if (mcBs.size() < 2)
return;
538 if (mcBs.size() > 2) {
539 B2WARNING(
"TagVertexModule:: Too many Bs found in MC");
548 if (!isReco(mcBs[0]) && !isReco(mcBs[1])) {
553 if (!isReco(mcBs[0]) && isReco(mcBs[1]))
554 swap(mcBs[0], mcBs[1]);
557 if (isReco(mcBs[0]) && isReco(mcBs[1])) {
558 double dist0 = (mcBs[0]->getDecayVertex() - Breco->
getVertex()).Mag2();
559 double dist1 = (mcBs[1]->getDecayVertex() - Breco->
getVertex()).Mag2();
561 swap(mcBs[0], mcBs[1]);
566 m_mcTagV = mcBs[1]->getDecayVertex();
575 bool isDecMode =
true;
577 const std::vector<Belle2::Particle*> recDau = Breco->
getDaughters();
578 const std::vector<Belle2::MCParticle*> genDau = Bgen->
getDaughters();
580 if (recDau.size() > 0 && genDau.size() > 0) {
581 for (
auto dauRec : recDau) {
583 for (
auto dauGen : genDau) {
584 if (dauGen->getPDG() == dauRec->getPDGCode())
587 if (!isDau) isDecMode =
false;
590 if (recDau.size() == 0) {
592 }
else {isDecMode =
false;}
604 std::vector<const Particle*> fitParticles;
606 if (!roe)
return fitParticles;
609 if (ROEParticles.size() == 0)
return fitParticles;
611 for (
auto& ROEParticle : ROEParticles) {
612 HitPatternVXD roeTrackPattern = ROEParticle->getTrackFitResult()->getHitPatternVXD();
615 fitParticles.push_back(ROEParticle);
624 vector<ParticleAndWeight> particleAndWeights;
626 for (
const Particle* particle : tagParticles) {
627 ROOT::Math::PxPyPzEVector mom = particle->get4Vector();
628 if (!isfinite(mom.mag2()))
continue;
632 particleAndWeight.
weight = -1111.;
633 particleAndWeight.
particle = particle;
638 particleAndWeights.push_back(particleAndWeight);
641 return particleAndWeights;
653 unsigned n = particleAndWeights.size();
654 sort(particleAndWeights.begin(), particleAndWeights.end(),
661 for (
unsigned i = 0; i < n; ++i) {
694 for (
const auto& pw : particleAndWeights) {
709 rFit.
addTrack(pw.particle->getTrackFitResult());
711 }
catch (
const rave::CheckedFloatException&) {
712 B2ERROR(
"Exception caught in TagVertexModule::makeGeneralFitRave(): Invalid inputs (nan/inf)?");
721 isGoodFit = rFit.
fit(
"avf");
723 if (isGoodFit < 1)
return false;
724 }
catch (
const rave::CheckedFloatException&) {
725 B2ERROR(
"Exception caught in TagVertexModule::makeGeneralFitRave(): Invalid inputs (nan/inf)?");
731 for (
unsigned int i(0); i < particleAndWeights.size() && isGoodFit >= 1; ++i)
732 particleAndWeights.at(i).weight = rFit.
getWeight(i);
759 CLHEP::HepSymMatrix err(7, 0);
761 err.sub(5, ROOTToCLHEP::getHepSymMatrix(
m_pvCov));
762 kFit.setIpTubeProfile(
774 for (
auto& pawi : particleAndWeights) {
777 if (pawi.mcParticle) {
778 addedOK = kFit.addTrack(
779 ROOTToCLHEP::getHepLorentzVector(pawi.mcParticle->get4Vector()),
780 ROOTToCLHEP::getPoint3DFromB2Vector(
getTruePoca(pawi)),
781 ROOTToCLHEP::getHepSymMatrix(pawi.particle->getMomentumVertexErrorMatrix()),
782 pawi.particle->getCharge());
787 if (pawi.mcParticle) {
788 addedOK = kFit.addTrack(
789 ROOTToCLHEP::getHepLorentzVector(pawi.mcParticle->get4Vector()),
791 ROOTToCLHEP::getHepSymMatrix(pawi.particle->getMomentumVertexErrorMatrix()),
792 pawi.particle->getCharge());
797 addedOK = kFit.addParticle(pawi.particle);
803 B2WARNING(
"TagVertexModule::makeGeneralFitKFit: failed to add a track");
809 int nTracksAdded = kFit.getTrackCount();
812 if ((nTracksAdded < 2 &&
m_constraintType ==
"noConstraint") || nTracksAdded < 1)
815 int isGoodFit = kFit.doFit();
824 int largest_chi2_trackid = -1;
825 double largest_track_chi2 = -1;
826 for (
int i = 0; i < kFit.getTrackCount(); ++i) {
827 double track_chi2 = kFit.getTrackCHIsq(i);
828 if (track_chi2 > largest_track_chi2) {
829 largest_track_chi2 = track_chi2;
830 largest_chi2_trackid = i;
833 return largest_chi2_trackid;
848 for (
int iteration_counter = 0; iteration_counter < 100; ++iteration_counter) {
856 if (nTracks !=
int(particleAndWeights.size()))
857 B2ERROR(
"TagVertexModule: Different number of tracks in kFit and particles");
863 int badTrackID = getLargestChi2ID(kFitTemp);
864 if (0 <= badTrackID && badTrackID <
int(particleAndWeights.size()))
865 particleAndWeights.erase(particleAndWeights.begin() + badTrackID);
867 B2ERROR(
"TagVertexModule: Obtained badTrackID is not within limits");
877 fillTagVinfo(CLHEPToROOT::getXYZVector(kFit.getVertex()),
878 CLHEPToROOT::getTMatrixDSym(kFit.getVertexError()));
891 ROOT::Math::XYZVector boostDir = ROOT::Math::XYZVector(boost.
Unit());
897 double dl = dVert.Dot(boostDir);
902 double MCdl = MCdVert.Dot(boostDir);
926 TVectorD oVec = toVec(oboost);
946 int nvert = rsg.
fit(fitType);
961 B2ERROR(
"In TagVertexModule::getTrackWithTrueCoordinate: no MC particle set");
980 B2ERROR(
"In TagVertexModule::getTruePoca: no MC particle set");
981 return ROOT::Math::XYZVector(0., 0., 0.);
1005 B2ERROR(
"In TagVertexModule::getTruePoca: no MC particle set");
1006 return ROOT::Math::XYZVector(0., 0., 0.);
1050 std::ostringstream oss;
1052 oss <<
"(" << std::setw(w) << vec.X() <<
", " << std::setw(w) << vec.Y() <<
", " << std::setw(w) << vec.Z() <<
")" << std::endl;
1059 std::ostringstream oss;
1061 for (
int i = 0; i < mat.GetNrows(); ++i) {
1062 for (
int j = 0; j < mat.GetNcols(); ++j) {
1063 oss << std::setw(w) << mat(i, j) <<
" ";
1073 std::ostringstream oss;
1075 for (
int i = 0; i < mat.GetNrows(); ++i) {
1076 for (
int j = 0; j < mat.GetNcols(); ++j) {
1077 oss << std::setw(w) << mat(i, j) <<
" ";
DataType Mag() const
The magnitude (rho in spherical coordinate system).
DataType Mag2() const
The magnitude squared (rho^2 in spherical coordinate system).
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
static const double speedOfLight
[cm/ns]
Hit pattern of the VXD within a track.
unsigned short getNPXDHits() const
Get total number of hits in the PXD.
A Class to store the Monte Carlo particle information.
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
ROOT::Math::XYZVector getProductionVertex() const
Return production vertex position.
ROOT::Math::PxPyPzEVector get4Vector() const
Return 4Vector of particle.
int getPDG() const
Return PDG code of particle.
ROOT::Math::XYZVector getMomentum() const
Return momentum.
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Class to store reconstructed particles.
TMatrixFSym getVertexErrorMatrix() const
Returns the 3x3 position error sub-matrix.
double getPValue() const
Returns chi^2 probability of fit if done or -1.
ROOT::Math::XYZVector getVertex() const
Returns vertex position (POCA for charged, IP for neutral FS particles)
int getPDGCode(void) const
Returns PDG code.
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
void setMomentumVertexErrorMatrix(const TMatrixFSym &errMatrix)
Sets 7x7 error matrix.
ROOT::Math::XYZVector getMomentum() const
Returns momentum vector.
void setPValue(double pValue)
Sets chi^2 probability of fit.
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
const TrackFitResult * getTrackFitResult() const
Returns the pointer to the TrackFitResult that was used to create this Particle (ParticleType == c_Tr...
double getMomentumMagnitude() const
Returns momentum magnitude.
double getP() const
Returns momentum magnitude (same as getMomentumMagnitude but with shorter name)
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
This is a general purpose class for collecting reconstructed MDST data objects that are not used in r...
std::vector< const Particle * > getChargedParticles(const std::string &maskName=c_defaultMaskName, unsigned int pdg=0, bool unpackComposite=true) const
Get charged particles from ROE mask.
static constexpr const char * c_defaultMaskName
Default mask name.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Accessor to arrays stored in the data store.
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
int m_fitTruthStatus
Store info about whether the fit was performed with the truth info 0 fit performed with measured para...
TMatrixDSym m_constraintCov
constraint to be used in the tag vertex fit
double m_truthTagVol
MC tagV component in the direction orthogonal to the boost.
std::vector< const Particle * > m_tagParticles
tracks of the rest of the event
double m_tagVl
tagV component in the boost direction
bool doVertexFit(const Particle *Breco)
central method for the tag side vertex fit
std::vector< const Particle * > getTagTracks_standardAlgorithm(const Particle *Breco, int nPXDHits) const
performs the fit using the standard algorithm - using all tracks in RoE The user can specify a reques...
double m_truthTagVl
MC tagV component in the boost direction
std::pair< ROOT::Math::XYZVector, TMatrixDSym > findConstraintBoost(double cut) const
calculate the standard constraint for the vertex fit on the tag side
bool m_useTruthInFit
Set to true if the tag fit is to be made with the TRUE tag track momentum and position.
std::vector< double > m_raveWeights
Store the weights used by Rave in the vtx fit so that they can be accessed later.
TrackFitResult getTrackWithRollBackCoordinates(ParticleAndWeight const &paw)
If the fit has to be done with the rolled back tracks, Rave or KFit is fed with a track where the pos...
virtual void initialize() override
Initialize the Module.
void fillTagVinfo(const ROOT::Math::XYZVector &tagVpos, const TMatrixDSym &tagVposErr)
Fill tagV vertex info.
ROOT::Math::XYZVector m_mcVertReco
generated Breco decay vertex
TMatrixDSym m_pvCov
covariance matrix of the PV (useful with tube and KFit)
double m_mcDeltaT
generated DeltaT with boost-direction approximation
static std::string printVector(const ROOT::Math::XYZVector &vec)
Print a XYZVector (useful for debugging)
ROOT::Math::XYZVector m_tagV
tag side fit result
virtual void event() override
Event processor.
std::pair< ROOT::Math::XYZVector, TMatrixDSym > findConstraintBTube(const Particle *Breco, double cut)
calculate constraint for the vertex fit on the tag side using the B tube (cylinder along the expected...
Particle * doVertexFitForBTube(const Particle *mother, std::string fitType) const
it returns an intersection between B rec and beam spot (= origin of BTube)
std::string m_listName
Breco particle list name.
std::vector< ParticleAndWeight > getParticlesAndWeights(const std::vector< const Particle * > &tagParticles) const
Get a list of particles with attached weight and associated MC particle.
bool m_useRollBack
Set to true if the tag fit is to be made with the tag track position rolled back to mother B.
double m_tagVlErr
Error of the tagV component in the boost direction
std::string m_roeMaskName
ROE particles from this mask will be used for vertex fitting.
double m_tagVChi2
chi^2 value of the tag vertex fit result
static ROOT::Math::XYZVector getTruePoca(ParticleAndWeight const &paw)
This finds the point on the true particle trajectory closest to the measured track position.
void BtagMCVertex(const Particle *Breco)
get the vertex of the MC B particle associated to Btag.
std::pair< ROOT::Math::XYZVector, TMatrixDSym > findConstraint(const Particle *Breco, double cut) const
calculate the constraint for the vertex fit on the tag side using Breco information
bool m_mcInfo
true if user wants to retrieve MC information out from the tracks used in the fit
double m_kFitReqReducedChi2
The required chi2/ndf to accept the kFit result, if it is higher, iteration procedure is applied.
std::string m_useMCassociation
No MC association or standard Breco particle or internal MCparticle association.
ROOT::Math::XYZVector m_constraintCenter
centre position of the constraint for the tag Vertex fit
double m_tagVolErr
Error of the tagV component in the direction orthogonal to the boost.
double m_mcTagLifeTime
generated tag side life time of B-decay
ROOT::Math::XYZVector m_BeamSpotCenter
Beam spot position.
double m_tagVNDF
Number of degrees of freedom in the tag vertex fit.
double m_deltaTErr
reconstructed DeltaT error
std::string m_fitAlgo
Algorithm used for the tag fit (Rave or KFit)
bool makeGeneralFit()
TO DO: tag side vertex fit in the case of semileptonic tag side decay.
ROOT::Math::XYZVector getRollBackPoca(ParticleAndWeight const &paw)
This shifts the position of tracks by the vector difference of mother B and production point of track...
virtual void beginRun() override
Called when entering a new run.
double m_mcDeltaTau
generated DeltaT
double m_fitPval
P value of the tag side fit result.
StoreArray< TagVertex > m_verArray
StoreArray of TagVertexes.
DBObjPtr< BeamSpot > m_beamSpotDB
Beam spot database object.
void deltaT(const Particle *Breco)
calculate DeltaT and MC-DeltaT (rec - tag) in ps from Breco and Btag vertices DT = Dl / gamma beta c ...
std::vector< const Particle * > m_raveParticles
tracks given to rave for the track fit (after removing Kshorts
void fillParticles(std::vector< ParticleAndWeight > &particleAndWeights)
Fill sorted list of particles into external variable.
TagVertexModule()
Constructor.
int m_reqPXDHits
N of PXD hits for a track to be used.
double m_confidenceLevel
required fit confidence level
double m_tagVChi2IP
IP component of the chi^2 of the tag vertex fit result.
double m_tagVol
tagV component in the direction orthogonal to the boost
analysis::VertexFitKFit doSingleKfit(std::vector< ParticleAndWeight > &particleAndWeights)
performs single KFit on particles stored in particleAndWeights this function can be iterated several ...
std::string m_constraintType
Choose constraint: noConstraint, IP, tube, boost, (breco)
void resetReturnParams()
Reset all parameters that are computed in each event and then used to compute tuple variables.
ROOT::Math::PxPyPzEVector m_tagMomentum
B tag momentum computed from fully reconstructed B sig.
bool makeGeneralFitRave()
make the vertex fit on the tag side: RAVE AVF tracks coming from Ks removed all other tracks used
TrackFitResult getTrackWithTrueCoordinates(ParticleAndWeight const &paw) const
If the fit has to be done with the truth info, Rave is fed with a track where the momentum is replace...
StoreArray< MCParticle > m_mcParticles
StoreArray of MCParticles.
double m_mcLifeTimeReco
generated Breco life time
static bool compBrecoBgen(const Particle *Breco, const MCParticle *Bgen)
compare Breco with the two MC B particles
double m_Bfield
magnetic field from data base
TMatrixDSym m_tagVErrMatrix
Error matrix of the tag side fit result.
int m_mcPDG
generated tag side B flavor
StoreObjPtr< ParticleList > m_plist
input particle list
static std::string printMatrix(const TMatrixD &mat)
Print a TMatrix (useful for debugging)
ROOT::Math::XYZVector m_mcTagV
generated tag side vertex
int m_rollbackStatus
Store info about whether the fit was performed with the rolled back tracks 0 fit performed with measu...
std::string m_trackFindingType
Choose how to find the tag tracks: standard, standard_PXD.
bool m_verbose
choose if you want to print extra infos
int m_FitType
fit algo used
bool makeGeneralFitKFit()
make the vertex fit on the tag side: KFit tracks coming from Ks removed all other tracks used
std::vector< const MCParticle * > m_raveMCParticles
Store the MC particles corresponding to each track used by Rave in the vtx fit.
TMatrixDSym m_BeamSpotCov
size of the beam spot == covariance matrix on the beam spot position
double m_deltaT
reconstructed DeltaT
TagVertex data object: contains Btag Vertex and DeltaT.
void setConstraintType(const std::string &constraintType)
Set the type of the constraint for the tag fit.
void setTagVlErr(float TagVlErr)
Set the error of the tagV component in the boost direction.
void setTruthTagVl(float TruthTagVl)
Set the MC tagV component in the boost direction.
void setTagVertex(const B2Vector3D &TagVertex)
Set BTag Vertex.
void setConstraintCenter(const B2Vector3D &constraintCenter)
Set the centre of the constraint for the tag fit.
void setTruthTagVol(float TruthTagVol)
Set the tagV component in the direction orthogonal to the boost.
void setMCTagBFlavor(int mcTagBFlavor)
Set generated Btag PDG code.
void setTagVolErr(float TagVolErr)
Set the error of the tagV component in the direction orthogonal to the boost.
void setMCTagVertex(const B2Vector3D &mcTagVertex)
Set generated BTag Vertex.
void setTagVNDF(float TagVNDF)
Set the number of degrees of freedom in the tag vertex fit.
void setDeltaTErr(float DeltaTErr)
Set DeltaTErr.
void setNTracks(int nTracks)
Set number of tracks used in the fit.
void setTagVChi2(float TagVChi2)
Set the chi^2 value of the tag vertex fit result.
void setMCDeltaT(float mcDeltaT)
Set generated DeltaT (in kin.
void setRollBackStatus(int backStatus)
Set the status of the fit performed with the rolled back tracks.
void setVertexFitMCParticles(const std::vector< const MCParticle * > &vtxFitMCParticles)
Set a vector of pointers to the MC p'cles corresponding to the tracks in the tag vtx fit.
void setTagVol(float TagVol)
Set the tagV component in the direction orthogonal to the boost.
void setDeltaT(float DeltaT)
Set DeltaT.
void setRaveWeights(const std::vector< double > &raveWeights)
Set the weights used by Rave in the tag vtx fit.
void setTagVertexPval(float TagVertexPval)
Set BTag Vertex P value.
void setMCDeltaTau(float mcDeltaTau)
Set generated DeltaT.
void setTagVl(float TagVl)
Set the tagV component in the boost direction.
void setTagVertexErrMatrix(const TMatrixDSym &TagVertexErrMatrix)
Set BTag Vertex (3x3) error matrix.
void setConstraintCov(const TMatrixDSym &constraintCov)
Set the covariance matrix of the constraint for the tag fit.
void setFitType(float FitType)
Set fit algo type.
void setVertexFitParticles(const std::vector< const Particle * > &vtxFitParticles)
Set a vector of pointers to the tracks used in the tag vtx fit.
void setTagVChi2IP(float TagVChi2IP)
Set the IP component of the chi^2 value of the tag vertex fit result.
void setFitTruthStatus(int truthStatus)
Set the status of the fit performed with the truth info of the tracks.
Values of the result of a track fit with a given particle hypothesis.
float getNDF() const
Getter for number of degrees of freedom of the track fit.
short getChargeSign() const
Return track charge (1 or -1).
double getPValue() const
Getter for Chi2 Probability of the track fit.
TMatrixDSym getCovariance6() const
Position and Momentum Covariance Matrix.
Const::ParticleType getParticleType() const
Getter for ParticleType of the mass hypothesis of the track fit.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
ROOT::Math::XYZVector getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
virtual int getNDF(void) const
Get an NDF of the fit.
bool isFitted(void) const
Return false if fit is not performed yet or performed fit is failed; otherwise true.
int getTrackCount(void) const
Get the number of added tracks.
void unsetBeamSpot()
unset beam spot constraint
static void initialize(int verbosity=1, double MagneticField=1.5)
Set everything up so everything needed for vertex fitting is there.
static RaveSetup * getInstance()
get the pointer to the instance to get/set any of options stored in RaveSetup
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...
void reset()
frees memory allocated by initialize().
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
double getNdf(VecSize vertexId=0) const
get the number of degrees of freedom (NDF) of the fitted vertex.
double getChi2(VecSize vertexId=0) const
get the χ² of the fitted vertex.
B2Vector3D getPos(VecSize vertexId=0) const
get the position of the fitted vertex.
void updateDaughters()
update the Daughters particles
double getWeight(int trackId, VecSize vertexId=0) const
get the weight Rave assigned to a specific input track.
double getPValue(VecSize vertexId=0) const
get the p value of the fitted vertex.
VertexFitKFit is a derived class from KFitBase to perform vertex-constraint kinematical fit.
double getCHIsq(void) const override
Get a chi-square of the fit.
static const double realNaN
This collects the B-meson properties in the hadronic B-decays It is used for the Ecms calibration in ...
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
B2Vector3< double > B2Vector3D
typedef for common usage with double
double sqrt(double a)
sqrt for double
Particle * copyParticle(const Particle *original)
Function takes argument Particle and creates a copy of it and copies of all its (grand-)^n-daughters.
Abstract base class for different kinds of events.
this struct is used to store and sort the tag tracks
const Particle * particle
tag track fit result with pion mass hypo, for sorting purposes
const MCParticle * mcParticle
mc particle matched to the tag track, for sorting purposes
double weight
rave weight associated to the track, for sorting purposes