CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Classes | Public Member Functions | Protected Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
BPHKinematicFit Class Reference

#include <BPHKinematicFit.h>

Inheritance diagram for BPHKinematicFit:
BPHDecayVertex BPHDecayMomentum BPHRecoCandidate BPHPlusMinusCandidate

Classes

struct  FlyingParticle
 

Public Member Functions

 BPHKinematicFit (const BPHKinematicFit &x)=delete
 
double constrMass () const
 retrieve the constraint More...
 
double constrSigma () const
 
virtual const
RefCountedKinematicVertex 
currentDecayVertex () const
 
virtual const
RefCountedKinematicParticle 
currentParticle () const
 get current particle More...
 
bool getIndependentFit (const std::string &name) const
 retrieve independent fit flag More...
 
double getMassSigma (const reco::Candidate *cand) const
 retrieve particle mass sigma More...
 
virtual bool isEmpty () const
 get fit status More...
 
virtual bool isValidFit () const
 
virtual const
RefCountedKinematicTree
kinematicTree () const
 perform the kinematic fit and get the result More...
 
virtual const
RefCountedKinematicTree
kinematicTree (const std::string &name, double mass, double sigma) const
 
virtual const
RefCountedKinematicTree
kinematicTree (const std::string &name, double mass) const
 
virtual const
RefCountedKinematicTree
kinematicTree (const std::string &name) const
 
virtual const
RefCountedKinematicTree
kinematicTree (const std::string &name, KinematicConstraint *kc) const
 
virtual const
RefCountedKinematicTree
kinematicTree (const std::string &name, MultiTrackKinematicConstraint *kc) const
 
virtual const std::vector
< RefCountedKinematicParticle > & 
kinParticles () const
 get kinematic particles More...
 
virtual std::vector
< RefCountedKinematicParticle
kinParticles (const std::vector< std::string > &names) const
 
virtual ParticleMass mass () const
 
BPHKinematicFitoperator= (const BPHKinematicFit &x)=delete
 
virtual const
math::XYZTLorentzVector
p4 () 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...
 
void setIndependentFit (const std::string &name, bool flag=true, double mass=-1.0, double sigma=-1.0)
 set a decaying daughter as an unique particle fitted independently More...
 
virtual const
RefCountedKinematicVertex 
topDecayVertex () const
 
virtual const
RefCountedKinematicParticle 
topParticle () const
 get top particle More...
 
 ~BPHKinematicFit () override
 
- Public Member Functions inherited from BPHDecayVertex
 BPHDecayVertex (const BPHDecayVertex &x)=delete
 
const edm::EventSetupgetEventSetup () const
 retrieve EventSetup More...
 
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...
 
BPHDecayVertexoperator= (const BPHDecayVertex &x)=delete
 
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 (VertexFitter< 5 > *fitter=nullptr, const reco::BeamSpot *bs=nullptr, const GlobalPoint *priorPos=nullptr, const GlobalError *priorError=nullptr) const
 get reconstructed vertex More...
 
 ~BPHDecayVertex () override
 
- Public Member Functions inherited from BPHDecayMomentum
 BPHDecayMomentum (const BPHDecayMomentum &x)=delete
 
virtual const std::vector
< std::string > & 
compNames () const
 
virtual const
pat::CompositeCandidate
composite () 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
 
BPHDecayMomentumoperator= (const BPHDecayMomentum &x)=delete
 
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
 
virtual void fill (BPHRecoCandidate *ptr, int level) const =0
 

Private Member Functions

virtual void addParticles (std::vector< RefCountedKinematicParticle > &kl, std::map< const reco::Candidate *, RefCountedKinematicParticle > &km, std::map< const BPHRecoCandidate *, RefCountedKinematicParticle > &cm) const
 
virtual void buildParticles () const
 
virtual void fitMomentum () const
 
virtual void getParticles (const std::string &moth, const std::string &daug, std::vector< RefCountedKinematicParticle > &kl, std::set< RefCountedKinematicParticle > &ks) const
 
virtual void getParticles (const std::string &moth, const std::vector< std::string > &daug, std::vector< RefCountedKinematicParticle > &kl, std::set< RefCountedKinematicParticle > &ks) const
 
virtual const
RefCountedKinematicTree
kinematicTree (const std::vector< RefCountedKinematicParticle > &kPart, MultiTrackKinematicConstraint *kc) const
 
virtual unsigned int numParticles (const BPHKinematicFit *cand=nullptr) const
 
virtual const BPHKinematicFitsplitKP (const std::string &name, std::vector< RefCountedKinematicParticle > *kComp, std::vector< RefCountedKinematicParticle > *kTail=nullptr) const
 

Static Private Member Functions

static void insertParticle (RefCountedKinematicParticle &kp, std::vector< RefCountedKinematicParticle > &kl, std::set< RefCountedKinematicParticle > &ks)
 

Private Attributes

std::vector
< RefCountedKinematicParticle
allParticles
 
std::map< const
BPHRecoCandidate
*, FlyingParticle
cKinP
 
std::map< const
reco::Candidate *, double > 
dMSig
 
std::map< const
BPHRecoCandidate
*, RefCountedKinematicParticle
kCDMap
 
std::map< const
reco::Candidate
*, RefCountedKinematicParticle
kinMap
 
RefCountedKinematicTree kinTree
 
double massConst
 
double massSigma
 
bool oldFit
 
bool oldKPs
 
bool oldMom
 
std::vector< BPHRecoConstCandPtrtmpList
 
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 38 of file BPHKinematicFit.h.

Constructor & Destructor Documentation

BPHKinematicFit::BPHKinematicFit ( const BPHKinematicFit x)
delete

Constructors are protected this object can exist only as part of a derived class

BPHKinematicFit::~BPHKinematicFit ( )
override

Destructor

Definition at line 85 of file BPHKinematicFit.cc.

85 {}
BPHKinematicFit::BPHKinematicFit ( )
protected

Definition at line 39 of file BPHKinematicFit.cc.

40  : BPHDecayVertex(nullptr),
41  massConst(-1.0),
42  massSigma(-1.0),
43  oldKPs(true),
44  oldFit(true),
45  oldMom(true),
46  kinTree(nullptr) {}
BPHDecayVertex(const BPHDecayVertex &x)=delete
RefCountedKinematicTree kinTree
BPHKinematicFit::BPHKinematicFit ( const BPHKinematicFit ptr)
protected

Definition at line 48 of file BPHKinematicFit.cc.

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

49  : BPHDecayVertex(ptr, nullptr),
50  massConst(-1.0),
51  massSigma(-1.0),
52  oldKPs(true),
53  oldFit(true),
54  oldMom(true),
55  kinTree(nullptr) {
56  map<const reco::Candidate*, const reco::Candidate*> iMap;
57  const vector<const reco::Candidate*>& daug = daughters();
58  const vector<Component>& list = ptr->componentList();
59  int i;
60  int n = daug.size();
61  for (i = 0; i < n; ++i) {
62  const reco::Candidate* cand = daug[i];
63  iMap[originalReco(cand)] = cand;
64  }
65  for (i = 0; i < n; ++i) {
66  const Component& c = list[i];
67  dMSig[iMap[c.cand]] = c.msig;
68  }
69  const vector<BPHRecoConstCandPtr>& dComp = daughComp();
70  int j;
71  int m = dComp.size();
72  for (j = 0; j < m; ++j) {
73  const BPHRecoCandidate* rc = dComp[j].get();
74  const map<const reco::Candidate*, double>& dMap = rc->dMSig;
75  const map<const BPHRecoCandidate*, FlyingParticle>& cMap = rc->cKinP;
76  dMSig.insert(dMap.begin(), dMap.end());
77  cKinP.insert(cMap.begin(), cMap.end());
78  cKinP[rc];
79  }
80 }
const edm::EventSetup & c
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
BPHDecayVertex(const BPHDecayVertex &x)=delete
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
RefCountedKinematicTree kinTree
std::map< const BPHRecoCandidate *, FlyingParticle > cKinP
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
std::map< std::string, BPHRecoConstCandPtr > cMap

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 BPHRecoCandidate::add().

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::addParticles ( std::vector< RefCountedKinematicParticle > &  kl,
std::map< const reco::Candidate *, RefCountedKinematicParticle > &  km,
std::map< const BPHRecoCandidate *, RefCountedKinematicParticle > &  cm 
) const
privatevirtual

Definition at line 371 of file BPHKinematicFit.cc.

References cKinP, BPHRecoCandidate::clone(), BPHDecayMomentum::daughComp(), BPHDecayMomentum::daughters(), dMSig, BPHKinematicFit::FlyingParticle::flag, personalPlayback::fp, BPHDecayVertex::getTransientTrack(), isEmpty(), isValidFit(), visualization-live-secondInstance_cfg::m, reco::Candidate::mass(), BPHKinematicFit::FlyingParticle::mass, dqmiodumpmetadata::n, KinematicParticleFactoryFromTransientTrack::particle(), setConstraint(), BPHKinematicFit::FlyingParticle::sigma, tmpList, topParticle(), and groupFilesInBlocks::tt.

Referenced by buildParticles().

373  {
374  const vector<const reco::Candidate*>& daug = daughters();
376  int n = daug.size();
377  float chi = 0.0;
378  float ndf = 0.0;
379  while (n--) {
380  const reco::Candidate* cand = daug[n];
381  ParticleMass mass = cand->mass();
382  float sigma = dMSig.find(cand)->second;
383  if (sigma < 0)
384  sigma = 1.0e-7;
386  if (tt != nullptr)
387  kl.push_back(km[cand] = pFactory.particle(*tt, mass, chi, ndf, sigma));
388  }
389  const vector<BPHRecoConstCandPtr>& comp = daughComp();
390  int m = comp.size();
391  while (m--) {
392  const BPHRecoCandidate* cptr = comp[m].get();
393  const FlyingParticle& fp = cKinP.at(cptr);
394  if (fp.flag) {
395  BPHRecoCandidate* tptr = cptr->clone();
396  double mass = fp.mass;
397  double sigma = fp.sigma;
398  if ((mass > 0.0) && (sigma > 0.0))
399  tptr->setConstraint(mass, sigma);
400  tmpList.push_back(BPHRecoConstCandPtr(tptr));
401  if (tptr->isEmpty())
402  return;
403  if (!tptr->isValidFit())
404  return;
405  kl.push_back(cm[cptr] = tptr->topParticle());
406  } else {
407  cptr->addParticles(kl, km, cm);
408  }
409  }
410  return;
411 }
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
reco::TransientTrack * getTransientTrack(const reco::Candidate *cand) const
get TransientTrack for a daughter
virtual double mass() const =0
mass
double ParticleMass
Definition: ParticleMass.h:4
virtual bool isEmpty() const
get fit status
virtual ParticleMass mass() const
virtual const std::vector< const reco::Candidate * > & daughters() const
std::vector< BPHRecoConstCandPtr > tmpList
std::map< const reco::Candidate *, double > dMSig
virtual bool isValidFit() const
virtual BPHRecoCandidate * clone(int level=-1) const
std::map< const BPHRecoCandidate *, FlyingParticle > cKinP
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
void setConstraint(double mass, double sigma)
apply a mass constraint
virtual const RefCountedKinematicParticle topParticle() const
get top particle
void BPHKinematicFit::buildParticles ( ) const
privatevirtual

Definition at line 361 of file BPHKinematicFit.cc.

References addParticles(), allParticles, BPHDecayMomentum::daughFull(), kCDMap, kinMap, oldKPs, and findQualityFiles::size.

Referenced by kinParticles().

361  {
362  kinMap.clear();
363  kCDMap.clear();
364  allParticles.clear();
365  allParticles.reserve(daughFull().size());
367  oldKPs = false;
368  return;
369 }
std::map< const BPHRecoCandidate *, RefCountedKinematicParticle > kCDMap
std::map< const reco::Candidate *, RefCountedKinematicParticle > kinMap
std::vector< RefCountedKinematicParticle > allParticles
virtual void addParticles(std::vector< RefCountedKinematicParticle > &kl, std::map< const reco::Candidate *, RefCountedKinematicParticle > &km, std::map< const BPHRecoCandidate *, RefCountedKinematicParticle > &cm) const
virtual const std::vector< const reco::Candidate * > & daughFull() const
tuple size
Write out results.
double BPHKinematicFit::constrMass ( ) const

retrieve the constraint

Definition at line 99 of file BPHKinematicFit.cc.

References massConst.

Referenced by BPHRecoCandidate::fill().

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

Definition at line 101 of file BPHKinematicFit.cc.

References massSigma.

Referenced by BPHRecoCandidate::fill().

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

Definition at line 277 of file BPHKinematicFit.cc.

References isEmpty(), and kinTree.

277  {
278  if (isEmpty())
279  return RefCountedKinematicVertex(nullptr);
280  return kinTree->currentDecayVertex();
281 }
virtual bool isEmpty() const
get fit status
ReferenceCountingPointer< KinematicVertex > RefCountedKinematicVertex
RefCountedKinematicTree kinTree
const RefCountedKinematicParticle BPHKinematicFit::currentParticle ( ) const
virtual

get current particle

Definition at line 271 of file BPHKinematicFit.cc.

References isEmpty(), and kinTree.

271  {
272  if (isEmpty())
273  return RefCountedKinematicParticle(nullptr);
274  return kinTree->currentParticle();
275 }
virtual bool isEmpty() const
get fit status
RefCountedKinematicTree kinTree
ReferenceCountingPointer< KinematicParticle > RefCountedKinematicParticle
void BPHKinematicFit::fitMomentum ( ) const
privatevirtual

Definition at line 538 of file BPHKinematicFit.cc.

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

Referenced by p4().

538  {
539  if (isValidFit()) {
540  const KinematicState& ks = topParticle()->currentState();
541  GlobalVector tm = ks.globalMomentum();
542  double x = tm.x();
543  double y = tm.y();
544  double z = tm.z();
545  double m = ks.mass();
546  double e = sqrt((x * x) + (y * y) + (z * z) + (m * m));
547  totalMomentum.SetPxPyPzE(x, y, z, e);
548  } else {
549  edm::LogPrint("FitNotFound") << "BPHKinematicFit::fitMomentum: "
550  << "simple momentum sum computed";
552  const vector<const reco::Candidate*>& daug = daughters();
553  int n = daug.size();
554  while (n--)
555  tm += daug[n]->p4();
556  const vector<BPHRecoConstCandPtr>& comp = daughComp();
557  int m = comp.size();
558  while (m--)
559  tm += comp[m]->p4();
560  totalMomentum = tm;
561  }
562  oldMom = false;
563  return;
564 }
T y() const
Definition: PV3DBase.h:60
GlobalVector globalMomentum() const
virtual const std::vector< const reco::Candidate * > & daughters() 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
Log< level::Warning, true > LogPrint
virtual bool isValidFit() const
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
math::XYZTLorentzVector totalMomentum
virtual const RefCountedKinematicParticle topParticle() const
get top particle
T x() const
Definition: PV3DBase.h:59
bool BPHKinematicFit::getIndependentFit ( const std::string &  name) const

retrieve independent fit flag

Definition at line 321 of file BPHKinematicFit.cc.

References cKinP, AlCaHLTBitMon_QueryRunRegistry::comp, and BPHDecayMomentum::getComp().

Referenced by BPHRecoCandidate::fill().

321  {
322  const BPHRecoCandidate* comp = getComp(name).get();
323  map<const BPHRecoCandidate*, FlyingParticle>::const_iterator iter = cKinP.find(comp);
324  return (iter != cKinP.end() ? iter->second.flag : false);
325 }
std::map< const BPHRecoCandidate *, FlyingParticle > cKinP
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
double BPHKinematicFit::getMassSigma ( const reco::Candidate cand) const

retrieve particle mass sigma

Definition at line 315 of file BPHKinematicFit.cc.

References dMSig.

Referenced by BPHRecoCandidate::fill().

315  {
316  map<const reco::Candidate*, double>::const_iterator iter = dMSig.find(cand);
317  return (iter != dMSig.end() ? iter->second : -1);
318 }
std::map< const reco::Candidate *, double > dMSig
virtual void BPHKinematicFit::getParticles ( const std::string &  moth,
const std::string &  daug,
std::vector< RefCountedKinematicParticle > &  kl,
std::set< RefCountedKinematicParticle > &  ks 
) const
privatevirtual
virtual void BPHKinematicFit::getParticles ( const std::string &  moth,
const std::vector< std::string > &  daug,
std::vector< RefCountedKinematicParticle > &  kl,
std::set< RefCountedKinematicParticle > &  ks 
) const
privatevirtual
void BPHKinematicFit::insertParticle ( RefCountedKinematicParticle kp,
std::vector< RefCountedKinematicParticle > &  kl,
std::set< RefCountedKinematicParticle > &  ks 
)
staticprivate

Definition at line 477 of file BPHKinematicFit.cc.

479  {
480  if (ks.find(kp) != ks.end())
481  return;
482  kl.push_back(kp);
483  ks.insert(kp);
484  return;
485 }
bool BPHKinematicFit::isEmpty ( void  ) const
virtual

get fit status

Definition at line 256 of file BPHKinematicFit.cc.

References kinematicTree(), and kinTree.

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

256  {
257  kinematicTree();
258  if (kinTree.get() == nullptr)
259  return true;
260  return kinTree->isEmpty();
261 }
virtual const RefCountedKinematicTree & kinematicTree() const
perform the kinematic fit and get the result
RefCountedKinematicTree kinTree
bool BPHKinematicFit::isValidFit ( ) const
virtual

Definition at line 263 of file BPHKinematicFit.cc.

References topParticle().

Referenced by addParticles(), and fitMomentum().

263  {
265  if (kPart.get() == nullptr)
266  return false;
267  return kPart->currentState().isValid();
268 }
virtual const RefCountedKinematicParticle topParticle() const
get top particle
const RefCountedKinematicTree & BPHKinematicFit::kinematicTree ( ) const
virtual

perform the kinematic fit and get the result

Definition at line 160 of file BPHKinematicFit.cc.

References kinTree, massConst, massSigma, and oldFit.

Referenced by isEmpty().

160  {
161  if (oldFit)
162  return kinematicTree("", massConst, massSigma);
163  return kinTree;
164 }
virtual const RefCountedKinematicTree & kinematicTree() const
perform the kinematic fit and get the result
RefCountedKinematicTree kinTree
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) 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 RefCountedKinematicTree & BPHKinematicFit::kinematicTree ( const std::vector< RefCountedKinematicParticle > &  kPart,
MultiTrackKinematicConstraint kc 
) const
privatevirtual

Definition at line 525 of file BPHKinematicFit.cc.

References cppFunctionSkipper::exception, KinematicConstrainedVertexFitter::fit(), and kinTree.

526  {
527  try {
529  kinTree = cvf.fit(kPart, kc);
530  } catch (std::exception const&) {
531  edm::LogPrint("FitFailed") << "BPHKinematicFit::kinematicTree: "
532  << "kin fit reset";
533  kinTree = RefCountedKinematicTree(nullptr);
534  }
535  return kinTree;
536 }
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part)
Log< level::Warning, true > LogPrint
ReferenceCountingPointer< KinematicTree > RefCountedKinematicTree
RefCountedKinematicTree kinTree
const vector< RefCountedKinematicParticle > & BPHKinematicFit::kinParticles ( ) const
virtual

get kinematic particles

Definition at line 125 of file BPHKinematicFit.cc.

References allParticles, buildParticles(), and oldKPs.

Referenced by splitKP().

125  {
126  if (oldKPs)
127  buildParticles();
128  return allParticles;
129 }
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 297 of file BPHKinematicFit.cc.

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

Referenced by Particle.Particle::__str__(), DiObject.DiMuon::__str__(), BPHRecoCandidate::add(), BPHPlusMinusCandidate::build(), BPHRecoCandidate::build(), BPHDecayToV0DiffMassBuilder::buildCandidate(), setConstraint(), and setIndependentFit().

297  {
299  if (kPart.get() == nullptr)
300  return -1.0;
301  const KinematicState kStat = kPart->currentState();
302  if (kStat.isValid())
303  return kStat.mass();
304  return -1.0;
305 }
bool isValid() const
ParticleMass mass() const
virtual const RefCountedKinematicParticle topParticle() const
get top particle
unsigned int BPHKinematicFit::numParticles ( const BPHKinematicFit cand = nullptr) const
privatevirtual

Definition at line 461 of file BPHKinematicFit.cc.

References cKinP, BPHDecayMomentum::compNames(), BPHDecayMomentum::daughters(), BPHDecayMomentum::getComp(), and dqmiodumpmetadata::n.

Referenced by splitKP().

461  {
462  if (cand == nullptr)
463  cand = this;
464  unsigned int n = cand->daughters().size();
465  const vector<string>& cnames = cand->compNames();
466  int i = cnames.size();
467  while (i) {
468  const BPHRecoCandidate* comp = cand->getComp(cnames[--i]).get();
469  if (cKinP.at(comp).flag)
470  ++n;
471  else
472  n += numParticles(comp);
473  }
474  return n;
475 }
virtual unsigned int numParticles(const BPHKinematicFit *cand=nullptr) const
virtual const std::vector< const reco::Candidate * > & daughters() const
std::map< const BPHRecoCandidate *, FlyingParticle > cKinP
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
virtual const std::vector< std::string > & compNames() const
BPHKinematicFit& BPHKinematicFit::operator= ( const BPHKinematicFit x)
delete
const math::XYZTLorentzVector & BPHKinematicFit::p4 ( ) const
virtual

compute total momentum after the fit

Definition at line 308 of file BPHKinematicFit.cc.

References fitMomentum(), oldMom, and totalMomentum.

Referenced by Tau.Tau::dxy_approx(), Tau.Tau::dz(), and Lepton.Lepton::p4WithFSR().

308  {
309  if (oldMom)
310  fitMomentum();
311  return totalMomentum;
312 }
virtual void fitMomentum() const
math::XYZTLorentzVector totalMomentum
void BPHKinematicFit::resetKinematicFit ( ) const
virtual

reset the kinematic fit

Definition at line 250 of file BPHKinematicFit.cc.

References oldFit, oldKPs, and oldMom.

Referenced by BPHDecayToResFlyingBuilder::build(), and setNotUpdated().

250  {
251  oldKPs = oldFit = oldMom = true;
252  return;
253 }
void BPHKinematicFit::setConstraint ( double  mass,
double  sigma 
)

apply a mass constraint

Operations

Definition at line 91 of file BPHKinematicFit.cc.

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

Referenced by addParticles(), BPHOniaToMuMuBuilder::extractList(), and BPHRecoCandidate::fill().

91  {
92  oldKPs = oldFit = oldMom = true;
93  massConst = mass;
94  massSigma = sigma;
95  return;
96 }
virtual ParticleMass mass() const
void BPHKinematicFit::setIndependentFit ( const std::string &  name,
bool  flag = true,
double  mass = -1.0,
double  sigma = -1.0 
)

set a decaying daughter as an unique particle fitted independently

Definition at line 104 of file BPHKinematicFit.cc.

References cKinP, AlCaHLTBitMon_QueryRunRegistry::comp, BPHKinematicFit::FlyingParticle::flag, personalPlayback::fp, BPHDecayMomentum::getComp(), mass(), BPHKinematicFit::FlyingParticle::mass, mergeVDriftHistosByStation::name, oldFit, oldKPs, oldMom, and BPHKinematicFit::FlyingParticle::sigma.

Referenced by BPHDecayToResFlyingBuilder::build(), and BPHRecoCandidate::fill().

104  {
105  string::size_type pos = name.find('/');
106  if (pos != string::npos) {
107  edm::LogPrint("WrongRequest") << "BPHKinematicFit::setIndependentFit: "
108  << "cascade decay specification not admitted " << name;
109  return;
110  }
111  const BPHRecoCandidate* comp = getComp(name).get();
112  if (comp == nullptr) {
113  edm::LogPrint("ParticleNotFound") << "BPHKinematicFit::setIndependentFit: " << name << " daughter not found";
114  return;
115  }
116  oldKPs = oldFit = oldMom = true;
117  FlyingParticle& fp = cKinP[comp];
118  fp.flag = flag;
119  fp.mass = mass;
120  fp.sigma = sigma;
121  return;
122 }
virtual ParticleMass mass() const
uint16_t size_type
Log< level::Warning, true > LogPrint
std::map< const BPHRecoCandidate *, FlyingParticle > cKinP
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
void BPHKinematicFit::setNotUpdated ( ) const
overrideprotectedvirtual

Reimplemented from BPHDecayVertex.

Reimplemented in BPHPlusMinusCandidate.

Definition at line 354 of file BPHKinematicFit.cc.

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

Referenced by BPHPlusMinusCandidate::setNotUpdated().

354  {
357  return;
358 }
virtual void resetKinematicFit() const
reset the kinematic fit
void setNotUpdated() const override
const BPHKinematicFit * BPHKinematicFit::splitKP ( const std::string &  name,
std::vector< RefCountedKinematicParticle > *  kComp,
std::vector< RefCountedKinematicParticle > *  kTail = nullptr 
) const
privatevirtual

Definition at line 487 of file BPHKinematicFit.cc.

References allParticles, cKinP, AlCaHLTBitMon_QueryRunRegistry::comp, BPHDecayMomentum::getComp(), kinParticles(), kinTree, mergeVDriftHistosByStation::name, numParticles(), and oldFit.

489  {
490  kinTree = RefCountedKinematicTree(nullptr);
491  oldFit = false;
492  kinParticles();
493  if (allParticles.size() != numParticles())
494  return nullptr;
495  kComp->clear();
496  if (kTail == nullptr)
497  kTail = kComp;
498  else
499  kTail->clear();
500  if (name.empty()) {
501  *kComp = allParticles;
502  return this;
503  }
504  const BPHRecoCandidate* comp = getComp(name).get();
505  int ns;
506  if (comp != nullptr) {
507  ns = (cKinP.at(comp).flag ? 1 : numParticles(comp));
508  } else {
509  edm::LogPrint("ParticleNotFound") << "BPHKinematicFit::splitKP: " << name << " daughter not found";
510  *kTail = allParticles;
511  return nullptr;
512  }
513  vector<string> nfull(2);
514  nfull[0] = name;
515  nfull[1] = "*";
516  vector<RefCountedKinematicParticle> kPart = kinParticles(nfull);
517  vector<RefCountedKinematicParticle>::const_iterator iter = kPart.begin();
518  vector<RefCountedKinematicParticle>::const_iterator imid = iter + ns;
519  vector<RefCountedKinematicParticle>::const_iterator iend = kPart.end();
520  kComp->insert(kComp->end(), iter, imid);
521  kTail->insert(kTail->end(), imid, iend);
522  return comp;
523 }
virtual unsigned int numParticles(const BPHKinematicFit *cand=nullptr) const
virtual const std::vector< RefCountedKinematicParticle > & kinParticles() const
get kinematic particles
Log< level::Warning, true > LogPrint
ReferenceCountingPointer< KinematicTree > RefCountedKinematicTree
std::vector< RefCountedKinematicParticle > allParticles
RefCountedKinematicTree kinTree
std::map< const BPHRecoCandidate *, FlyingParticle > cKinP
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
const RefCountedKinematicVertex BPHKinematicFit::topDecayVertex ( ) const
virtual

Definition at line 290 of file BPHKinematicFit.cc.

References isEmpty(), and kinTree.

Referenced by BPHDecayToResFlyingBuilder::build().

290  {
291  if (isEmpty())
292  return RefCountedKinematicVertex(nullptr);
293  kinTree->movePointerToTheTop();
294  return kinTree->currentDecayVertex();
295 }
virtual bool isEmpty() const
get fit status
ReferenceCountingPointer< KinematicVertex > RefCountedKinematicVertex
RefCountedKinematicTree kinTree
const RefCountedKinematicParticle BPHKinematicFit::topParticle ( ) const
virtual

get top particle

Definition at line 284 of file BPHKinematicFit.cc.

References isEmpty(), and kinTree.

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

284  {
285  if (isEmpty())
286  return RefCountedKinematicParticle(nullptr);
287  return kinTree->topParticle();
288 }
virtual bool isEmpty() const
get fit status
RefCountedKinematicTree kinTree
ReferenceCountingPointer< KinematicParticle > RefCountedKinematicParticle

Member Data Documentation

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

Definition at line 148 of file BPHKinematicFit.h.

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

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

Definition at line 129 of file BPHKinematicFit.h.

Referenced by addParticles(), BPHKinematicFit(), and getMassSigma().

std::map<const BPHRecoCandidate*, RefCountedKinematicParticle> BPHKinematicFit::kCDMap
mutableprivate

Definition at line 147 of file BPHKinematicFit.h.

Referenced by buildParticles().

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

Definition at line 146 of file BPHKinematicFit.h.

Referenced by buildParticles().

RefCountedKinematicTree BPHKinematicFit::kinTree
mutableprivate
double BPHKinematicFit::massConst
private

Definition at line 125 of file BPHKinematicFit.h.

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

double BPHKinematicFit::massSigma
private

Definition at line 126 of file BPHKinematicFit.h.

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

bool BPHKinematicFit::oldFit
mutableprivate
bool BPHKinematicFit::oldKPs
mutableprivate
bool BPHKinematicFit::oldMom
mutableprivate
std::vector<BPHRecoConstCandPtr> BPHKinematicFit::tmpList
mutableprivate

Definition at line 140 of file BPHKinematicFit.h.

Referenced by addParticles().

math::XYZTLorentzVector BPHKinematicFit::totalMomentum
mutableprivate

Definition at line 150 of file BPHKinematicFit.h.

Referenced by fitMomentum(), and p4().