CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
BPHKinematicFit Class Reference

#include <BPHKinematicFit.h>

Inheritance diagram for BPHKinematicFit:
BPHDecayVertex BPHDecayMomentum BPHRecoCandidate BPHPlusMinusCandidate

Public Member Functions

double constrMass () const
 retrieve the constraint More...
 
double constrSigma () const
 
virtual const RefCountedKinematicVertex currentDecayVertex () const
 
virtual const RefCountedKinematicParticle currentParticle () const
 
virtual bool isEmpty () const
 
virtual bool isValidFit () const
 
virtual const RefCountedKinematicTreekinematicTree () const
 perform the kinematic fit and get the result More...
 
virtual const RefCountedKinematicTreekinematicTree (const std::string &name, double mass, double sigma) const
 
virtual const RefCountedKinematicTreekinematicTree (const std::string &name, double mass) const
 
virtual const RefCountedKinematicTreekinematicTree (const std::string &name, KinematicConstraint *kc) const
 
virtual const RefCountedKinematicTreekinematicTree (const std::string &name, MultiTrackKinematicConstraint *kc) const
 
virtual const std::vector< RefCountedKinematicParticle > & kinParticles () const
 get kinematic particles More...
 
virtual std::vector< RefCountedKinematicParticlekinParticles (const std::vector< std::string > &names) const
 
virtual ParticleMass mass () const
 
virtual const math::XYZTLorentzVectorp4 () const
 compute total momentum after the fit More...
 
virtual void resetKinematicFit () const
 reset the kinematic fit More...
 
void setConstraint (double mass, double sigma)
 apply a mass constraint More...
 
 ~BPHKinematicFit () override
 
- Public Member Functions inherited from BPHDecayVertex
const reco::TrackgetTrack (const reco::Candidate *cand) const
 get Track for a daughter More...
 
const std::string & getTrackSearchList (const reco::Candidate *cand) const
 retrieve track search list More...
 
reco::TransientTrackgetTransientTrack (const reco::Candidate *cand) const
 get TransientTrack for a daughter More...
 
const std::vector< const reco::Track * > & tracks () const
 get list of Tracks More...
 
const std::vector< reco::TransientTrack > & transientTracks () const
 get list of TransientTracks More...
 
virtual bool validTracks () const
 check for valid reconstructed vertex More...
 
virtual bool validVertex () const
 
virtual const reco::Vertexvertex () const
 get reconstructed vertex More...
 
 ~BPHDecayVertex () override
 
- Public Member Functions inherited from BPHDecayMomentum
virtual const std::vector< std::string > & compNames () const
 
virtual const pat::CompositeCandidatecomposite () const
 get a composite by the simple sum of simple particles More...
 
virtual const std::vector< BPHRecoConstCandPtr > & daughComp () const
 
virtual const std::vector< const reco::Candidate * > & daughFull () const
 
virtual const std::vector< const reco::Candidate * > & daughters () const
 
virtual const std::vector< std::string > & daugNames () const
 
virtual BPHRecoConstCandPtr getComp (const std::string &name) const
 
virtual const reco::CandidategetDaug (const std::string &name) const
 
virtual const reco::CandidateoriginalReco (const reco::Candidate *daug) const
 get the original particle from the clone More...
 
virtual ~BPHDecayMomentum ()
 

Protected Member Functions

virtual void addK (const std::string &name, const reco::Candidate *daug, double mass=-1.0, double sigma=-1.0)
 
virtual void addK (const std::string &name, const reco::Candidate *daug, const std::string &searchList, double mass=-1.0, double sigma=-1.0)
 
virtual void addK (const std::string &name, const BPHRecoConstCandPtr &comp)
 add a previously reconstructed particle giving it a name More...
 
 BPHKinematicFit ()
 
 BPHKinematicFit (const BPHKinematicFit *ptr)
 
void setNotUpdated () const override
 
- Protected Member Functions inherited from BPHDecayVertex
virtual void addV (const std::string &name, const reco::Candidate *daug, const std::string &searchList, double mass)
 
virtual void addV (const std::string &name, const BPHRecoConstCandPtr &comp)
 add a previously reconstructed particle giving it a name More...
 
 BPHDecayVertex (const edm::EventSetup *es)
 
 BPHDecayVertex (const BPHDecayVertex *ptr, const edm::EventSetup *es)
 
- Protected Member Functions inherited from BPHDecayMomentum
virtual void addP (const std::string &name, const reco::Candidate *daug, double mass=-1.0)
 
virtual void addP (const std::string &name, const BPHRecoConstCandPtr &comp)
 add a previously reconstructed particle giving it a name More...
 
 BPHDecayMomentum ()
 
 BPHDecayMomentum (const std::map< std::string, Component > &daugMap)
 
 BPHDecayMomentum (const std::map< std::string, Component > &daugMap, const std::map< std::string, BPHRecoConstCandPtr > compMap)
 
const std::vector< Component > & componentList () const
 

Private Member Functions

virtual void buildParticles () const
 
virtual void fitMomentum () const
 

Private Attributes

std::vector< RefCountedKinematicParticleallParticles
 
std::map< const reco::Candidate *, double > dMSig
 
std::map< const reco::Candidate *, RefCountedKinematicParticlekinMap
 
RefCountedKinematicTree kinTree
 
double massConst
 
double massSigma
 
bool oldFit
 
bool oldKPs
 
bool oldMom
 
math::XYZTLorentzVector totalMomentum
 

Detailed Description

Description: Highest-level base class to encapsulate kinematic fit operations

Author
Paolo Ronchese INFN Padova high-level base class to perform a kinematic fit

Definition at line 35 of file BPHKinematicFit.h.

Constructor & Destructor Documentation

BPHKinematicFit::~BPHKinematicFit ( )
override

Constructor is protected this object can exist only as part of a derived classDestructor

Definition at line 83 of file BPHKinematicFit.cc.

83 {}
BPHKinematicFit::BPHKinematicFit ( )
protected

Definition at line 41 of file BPHKinematicFit.cc.

42  : BPHDecayVertex(nullptr),
43  massConst(-1.0),
44  massSigma(-1.0),
45  oldKPs(true),
46  oldFit(true),
47  oldMom(true),
48  kinTree(nullptr) {}
BPHDecayVertex(const edm::EventSetup *es)
RefCountedKinematicTree kinTree
BPHKinematicFit::BPHKinematicFit ( const BPHKinematicFit ptr)
protected

Definition at line 50 of file BPHKinematicFit.cc.

References HltBtagPostValidation_cff::c, BPHDecayMomentum::Component::cand, BPHDecayMomentum::componentList(), BPHDecayMomentum::daughComp(), BPHDecayMomentum::daughters(), BPHDecayMomentum::dMap, dMSig, mps_fire::i, dqmiolumiharvest::j, list(), visualization-live-secondInstance_cfg::m, BPHDecayMomentum::Component::msig, dqmiodumpmetadata::n, and BPHDecayMomentum::originalReco().

51  : BPHDecayVertex(ptr, nullptr),
52  massConst(-1.0),
53  massSigma(-1.0),
54  oldKPs(true),
55  oldFit(true),
56  oldMom(true),
57  kinTree(nullptr) {
58  map<const reco::Candidate*, const reco::Candidate*> iMap;
59  const vector<const reco::Candidate*>& daug = daughters();
60  const vector<Component>& list = ptr->componentList();
61  int i;
62  int n = daug.size();
63  for (i = 0; i < n; ++i) {
64  const reco::Candidate* cand = daug[i];
65  iMap[originalReco(cand)] = cand;
66  }
67  for (i = 0; i < n; ++i) {
68  const Component& c = list[i];
69  dMSig[iMap[c.cand]] = c.msig;
70  }
71  const vector<BPHRecoConstCandPtr>& dComp = daughComp();
72  int j;
73  int m = dComp.size();
74  for (j = 0; j < m; ++j) {
75  const map<const reco::Candidate*, double>& dMap = dComp[j]->dMSig;
76  dMSig.insert(dMap.begin(), dMap.end());
77  }
78 }
std::map< std::string, const reco::Candidate * > dMap
virtual const std::vector< const reco::Candidate * > & daughters() const
std::map< const reco::Candidate *, double > dMSig
const std::vector< Component > & componentList() const
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
BPHDecayVertex(const edm::EventSetup *es)
RefCountedKinematicTree kinTree
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger list("!*","!HLTx*"if it matches 2 triggers or more) will accept the event if all the matching triggers are FAIL.It will reject the event if any of the triggers are PASS or EXCEPTION(this matches the behavior of"!*"before the partial wildcard feature was incorporated).Triggers which are in the READY state are completely ignored.(READY should never be returned since the trigger paths have been run

Member Function Documentation

virtual void BPHKinematicFit::addK ( const std::string &  name,
const reco::Candidate daug,
double  mass = -1.0,
double  sigma = -1.0 
)
protectedvirtual

add a simple particle giving it a name particles are cloned, eventually specifying a different mass and a sigma

Referenced by p4(), and BPHPlusMinusCandidate::~BPHPlusMinusCandidate().

virtual void BPHKinematicFit::addK ( const std::string &  name,
const reco::Candidate daug,
const std::string &  searchList,
double  mass = -1.0,
double  sigma = -1.0 
)
protectedvirtual

add a simple particle and specify a criterion to search for the associated track

virtual void BPHKinematicFit::addK ( const std::string &  name,
const BPHRecoConstCandPtr comp 
)
protectedvirtual

add a previously reconstructed particle giving it a name

void BPHKinematicFit::buildParticles ( ) const
privatevirtual

Definition at line 339 of file BPHKinematicFit.cc.

References allParticles, BPHDecayMomentum::daughFull(), dMSig, BPHDecayVertex::getTransientTrack(), kinMap, mass(), reco::Candidate::mass(), dqmiodumpmetadata::n, oldKPs, KinematicParticleFactoryFromTransientTrack::particle(), and groupFilesInBlocks::tt.

Referenced by kinParticles().

339  {
340  kinMap.clear();
341  allParticles.clear();
342  const vector<const reco::Candidate*>& daug = daughFull();
344  int n = daug.size();
345  allParticles.reserve(n);
346  float chi = 0.0;
347  float ndf = 0.0;
348  while (n--) {
349  const reco::Candidate* cand = daug[n];
350  ParticleMass mass = cand->mass();
351  float sigma = dMSig.find(cand)->second;
352  if (sigma < 0)
353  sigma = 1.0e-7;
355  if (tt != nullptr)
356  allParticles.push_back(kinMap[cand] = pFactory.particle(*tt, mass, chi, ndf, sigma));
357  }
358  oldKPs = false;
359  return;
360 }
reco::TransientTrack * getTransientTrack(const reco::Candidate *cand) const
get TransientTrack for a daughter
double ParticleMass
Definition: ParticleMass.h:4
virtual ParticleMass mass() const
std::map< const reco::Candidate *, RefCountedKinematicParticle > kinMap
std::map< const reco::Candidate *, double > dMSig
std::vector< RefCountedKinematicParticle > allParticles
virtual double mass() const =0
mass
virtual const std::vector< const reco::Candidate * > & daughFull() const
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
double BPHKinematicFit::constrMass ( ) const

retrieve the constraint

Definition at line 95 of file BPHKinematicFit.cc.

References massConst.

95 { return massConst; }
double BPHKinematicFit::constrSigma ( ) const

Definition at line 97 of file BPHKinematicFit.cc.

References massSigma.

97 { return massSigma; }
const RefCountedKinematicVertex BPHKinematicFit::currentDecayVertex ( ) const
virtual

Definition at line 292 of file BPHKinematicFit.cc.

References isEmpty(), and kinTree.

292  {
293  if (isEmpty())
294  return RefCountedKinematicVertex(nullptr);
295  return kinTree->currentDecayVertex();
296 }
virtual bool isEmpty() const
ReferenceCountingPointer< KinematicVertex > RefCountedKinematicVertex
RefCountedKinematicTree kinTree
const RefCountedKinematicParticle BPHKinematicFit::currentParticle ( ) const
virtual

Definition at line 286 of file BPHKinematicFit.cc.

References isEmpty(), and kinTree.

Referenced by fitMomentum(), isValidFit(), and mass().

286  {
287  if (isEmpty())
288  return RefCountedKinematicParticle(nullptr);
289  return kinTree->currentParticle();
290 }
virtual bool isEmpty() const
RefCountedKinematicTree kinTree
ReferenceCountingPointer< KinematicParticle > RefCountedKinematicParticle
void BPHKinematicFit::fitMomentum ( ) const
privatevirtual

Definition at line 362 of file BPHKinematicFit.cc.

References AlCaHLTBitMon_QueryRunRegistry::comp, currentParticle(), BPHDecayMomentum::daughComp(), BPHDecayMomentum::daughters(), MillePedeFileConverter_cfg::e, KinematicState::globalMomentum(), isValidFit(), visualization-live-secondInstance_cfg::m, KinematicState::mass(), dqmiodumpmetadata::n, oldMom, mathSSE::sqrt(), totalMomentum, x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by p4().

362  {
363  if (isValidFit()) {
364  const KinematicState& ks = currentParticle()->currentState();
365  GlobalVector tm = ks.globalMomentum();
366  double x = tm.x();
367  double y = tm.y();
368  double z = tm.z();
369  double m = ks.mass();
370  double e = sqrt((x * x) + (y * y) + (z * z) + (m * m));
371  totalMomentum.SetPxPyPzE(x, y, z, e);
372  } else {
373  edm::LogPrint("FitNotFound") << "BPHKinematicFit::fitMomentum: "
374  << "simple momentum sum computed";
376  const vector<const reco::Candidate*>& daug = daughters();
377  int n = daug.size();
378  while (n--)
379  tm += daug[n]->p4();
380  const vector<BPHRecoConstCandPtr>& comp = daughComp();
381  int m = comp.size();
382  while (m--)
383  tm += comp[m]->p4();
384  totalMomentum = tm;
385  }
386  oldMom = false;
387  return;
388 }
T y() const
Definition: PV3DBase.h:60
GlobalVector globalMomentum() const
virtual const std::vector< const reco::Candidate * > & daughters() const
virtual const RefCountedKinematicParticle currentParticle() const
ParticleMass mass() const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
virtual bool isValidFit() const
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
math::XYZTLorentzVector totalMomentum
T x() const
Definition: PV3DBase.h:59
bool BPHKinematicFit::isEmpty ( void  ) const
virtual

Definition at line 272 of file BPHKinematicFit.cc.

References kinematicTree(), and kinTree.

Referenced by plotting.Plot::clone(), currentDecayVertex(), and currentParticle().

272  {
273  kinematicTree();
274  if (kinTree.get() == nullptr)
275  return true;
276  return kinTree->isEmpty();
277 }
virtual const RefCountedKinematicTree & kinematicTree() const
perform the kinematic fit and get the result
RefCountedKinematicTree kinTree
bool BPHKinematicFit::isValidFit ( ) const
virtual

Definition at line 279 of file BPHKinematicFit.cc.

References currentParticle().

Referenced by fitMomentum().

279  {
281  if (kPart.get() == nullptr)
282  return false;
283  return kPart->currentState().isValid();
284 }
virtual const RefCountedKinematicParticle currentParticle() const
const RefCountedKinematicTree & BPHKinematicFit::kinematicTree ( ) const
virtual
virtual const RefCountedKinematicTree& BPHKinematicFit::kinematicTree ( const std::string &  name,
double  mass,
double  sigma 
) const
virtual
virtual const RefCountedKinematicTree& BPHKinematicFit::kinematicTree ( const std::string &  name,
double  mass 
) const
virtual
virtual const RefCountedKinematicTree& BPHKinematicFit::kinematicTree ( const std::string &  name,
KinematicConstraint kc 
) const
virtual
virtual const RefCountedKinematicTree& BPHKinematicFit::kinematicTree ( const std::string &  name,
MultiTrackKinematicConstraint kc 
) const
virtual
const vector< RefCountedKinematicParticle > & BPHKinematicFit::kinParticles ( ) const
virtual

get kinematic particles

Definition at line 99 of file BPHKinematicFit.cc.

References allParticles, buildParticles(), BPHDecayMomentum::daughFull(), BPHDecayMomentum::getDaug(), mps_fire::i, dqmiolumiharvest::j, kinMap, kp, visualization-live-secondInstance_cfg::m, dqmiodumpmetadata::n, names, oldKPs, unpackData-CaloStage2::pname, and muonDTDigis_cfi::pset.

Referenced by kinematicTree().

99  {
100  if (oldKPs)
101  buildParticles();
102  return allParticles;
103 }
virtual void buildParticles() const
std::vector< RefCountedKinematicParticle > allParticles
virtual std::vector<RefCountedKinematicParticle> BPHKinematicFit::kinParticles ( const std::vector< std::string > &  names) const
virtual
ParticleMass BPHKinematicFit::mass ( ) const
virtual

Definition at line 298 of file BPHKinematicFit.cc.

References currentParticle(), KinematicState::isValid(), and KinematicState::mass().

Referenced by Particle.Particle::__str__(), DiObject.DiMuon::__str__(), BPHPlusMinusCandidate::build(), BPHRecoCandidate::build(), buildParticles(), kinematicTree(), p4(), setConstraint(), and BPHPlusMinusCandidate::~BPHPlusMinusCandidate().

298  {
300  if (kPart.get() == nullptr)
301  return -1.0;
302  const KinematicState kStat = kPart->currentState();
303  if (kStat.isValid())
304  return kStat.mass();
305  return -1.0;
306 }
bool isValid() const
virtual const RefCountedKinematicParticle currentParticle() const
ParticleMass mass() const
const math::XYZTLorentzVector & BPHKinematicFit::p4 ( ) const
virtual
void BPHKinematicFit::resetKinematicFit ( ) const
virtual

reset the kinematic fit

Definition at line 267 of file BPHKinematicFit.cc.

References oldFit, oldKPs, and oldMom.

Referenced by setNotUpdated().

267  {
268  oldKPs = oldFit = oldMom = true;
269  return;
270 }
void BPHKinematicFit::setConstraint ( double  mass,
double  sigma 
)

apply a mass constraint

Operations

Definition at line 88 of file BPHKinematicFit.cc.

References mass(), massConst, massSigma, oldFit, and oldMom.

Referenced by BPHOniaToMuMuBuilder::extractList().

88  {
89  oldFit = oldMom = true;
90  massConst = mass;
91  massSigma = sigma;
92  return;
93 }
virtual ParticleMass mass() const
void BPHKinematicFit::setNotUpdated ( ) const
overrideprotectedvirtual

Reimplemented from BPHDecayVertex.

Reimplemented in BPHPlusMinusCandidate.

Definition at line 333 of file BPHKinematicFit.cc.

References resetKinematicFit(), and BPHDecayVertex::setNotUpdated().

Referenced by BPHPlusMinusCandidate::setNotUpdated().

333  {
336  return;
337 }
void setNotUpdated() const override
virtual void resetKinematicFit() const
reset the kinematic fit

Member Data Documentation

std::vector<RefCountedKinematicParticle> BPHKinematicFit::allParticles
mutableprivate

Definition at line 115 of file BPHKinematicFit.h.

Referenced by buildParticles(), kinematicTree(), and kinParticles().

std::map<const reco::Candidate*, double> BPHKinematicFit::dMSig
private

Definition at line 108 of file BPHKinematicFit.h.

Referenced by BPHKinematicFit(), buildParticles(), and p4().

std::map<const reco::Candidate*, RefCountedKinematicParticle> BPHKinematicFit::kinMap
mutableprivate

Definition at line 114 of file BPHKinematicFit.h.

Referenced by buildParticles(), and kinParticles().

RefCountedKinematicTree BPHKinematicFit::kinTree
mutableprivate

Definition at line 116 of file BPHKinematicFit.h.

Referenced by currentDecayVertex(), currentParticle(), isEmpty(), and kinematicTree().

double BPHKinematicFit::massConst
private

Definition at line 104 of file BPHKinematicFit.h.

Referenced by constrMass(), kinematicTree(), and setConstraint().

double BPHKinematicFit::massSigma
private

Definition at line 105 of file BPHKinematicFit.h.

Referenced by constrSigma(), kinematicTree(), and setConstraint().

bool BPHKinematicFit::oldFit
mutableprivate

Definition at line 112 of file BPHKinematicFit.h.

Referenced by kinematicTree(), resetKinematicFit(), and setConstraint().

bool BPHKinematicFit::oldKPs
mutableprivate

Definition at line 111 of file BPHKinematicFit.h.

Referenced by buildParticles(), kinParticles(), and resetKinematicFit().

bool BPHKinematicFit::oldMom
mutableprivate

Definition at line 113 of file BPHKinematicFit.h.

Referenced by fitMomentum(), p4(), resetKinematicFit(), and setConstraint().

math::XYZTLorentzVector BPHKinematicFit::totalMomentum
mutableprivate

Definition at line 117 of file BPHKinematicFit.h.

Referenced by fitMomentum(), and p4().