CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BPHKinematicFit.h
Go to the documentation of this file.
1 #ifndef HeavyFlavorAnalysis_RecoDecay_BPHKinematicFit_h
2 #define HeavyFlavorAnalysis_RecoDecay_BPHKinematicFit_h
3 
12 //----------------------
13 // Base Class Headers --
14 //----------------------
16 
17 //------------------------------------
18 // Collaborating Class Declarations --
19 //------------------------------------
22 
25 
26 //---------------
27 // C++ Headers --
28 //---------------
29 #include <string>
30 #include <vector>
31 #include <map>
32 #include <set>
33 
34 // ---------------------
35 // -- Class Interface --
36 // ---------------------
37 
38 class BPHKinematicFit : public virtual BPHDecayVertex {
39 public:
43  // deleted copy constructor and assignment operator
44  BPHKinematicFit(const BPHKinematicFit& x) = delete;
45  BPHKinematicFit& operator=(const BPHKinematicFit& x) = delete;
46 
49  ~BPHKinematicFit() override;
50 
54  void setConstraint(double mass, double sigma);
57  double constrMass() const;
58  double constrSigma() const;
60  void setIndependentFit(const std::string& name, bool flag = true, double mass = -1.0, double sigma = -1.0);
61 
63  virtual const std::vector<RefCountedKinematicParticle>& kinParticles() const;
64  virtual std::vector<RefCountedKinematicParticle> kinParticles(const std::vector<std::string>& names) const;
65 
67  virtual const RefCountedKinematicTree& kinematicTree() const;
68  virtual const RefCountedKinematicTree& kinematicTree(const std::string& name, double mass, double sigma) const;
69  virtual const RefCountedKinematicTree& kinematicTree(const std::string& name, double mass) const;
70  virtual const RefCountedKinematicTree& kinematicTree(const std::string& name) const;
71  virtual const RefCountedKinematicTree& kinematicTree(const std::string& name, KinematicConstraint* kc) const;
72  virtual const RefCountedKinematicTree& kinematicTree(const std::string& name,
74 
76  virtual void resetKinematicFit() const;
77 
79  virtual bool isEmpty() const;
80  virtual bool isValidFit() const;
81 
83  virtual const RefCountedKinematicParticle currentParticle() const;
84  virtual const RefCountedKinematicVertex currentDecayVertex() const;
85 
87  virtual const RefCountedKinematicParticle topParticle() const;
88  virtual const RefCountedKinematicVertex topDecayVertex() const;
89  virtual ParticleMass mass() const;
90 
92  virtual const math::XYZTLorentzVector& p4() const;
93 
95  double getMassSigma(const reco::Candidate* cand) const;
96 
98  bool getIndependentFit(const std::string& name) const;
99 
100 protected:
101  // constructors
102  BPHKinematicFit();
103  // pointer used to retrieve informations from other bases
104  BPHKinematicFit(const BPHKinematicFit* ptr);
105 
109  virtual void addK(const std::string& name, const reco::Candidate* daug, double mass = -1.0, double sigma = -1.0);
112  virtual void addK(const std::string& name,
113  const reco::Candidate* daug,
114  const std::string& searchList,
115  double mass = -1.0,
116  double sigma = -1.0);
118  virtual void addK(const std::string& name, const BPHRecoConstCandPtr& comp);
119 
120  // utility function used to cash reconstruction results
121  void setNotUpdated() const override;
122 
123 private:
124  // mass constraint
125  double massConst;
126  double massSigma;
127 
128  // map linking daughters to mass sigma
129  std::map<const reco::Candidate*, double> dMSig;
130 
131  // map to handle composite daughters as single particles
132  struct FlyingParticle {
133  bool flag = false;
134  double mass = -1.0;
135  double sigma = -1.0;
136  };
137  std::map<const BPHRecoCandidate*, FlyingParticle> cKinP;
138 
139  // temporary particle set
140  mutable std::vector<BPHRecoConstCandPtr> tmpList;
141 
142  // reconstruction results cache
143  mutable bool oldKPs;
144  mutable bool oldFit;
145  mutable bool oldMom;
146  mutable std::map<const reco::Candidate*, RefCountedKinematicParticle> kinMap;
147  mutable std::map<const BPHRecoCandidate*, RefCountedKinematicParticle> kCDMap;
148  mutable std::vector<RefCountedKinematicParticle> allParticles;
151 
152  // build kin particles, perform the fit and compute the total momentum
153  virtual void buildParticles() const;
154  virtual void addParticles(std::vector<RefCountedKinematicParticle>& kl,
155  std::map<const reco::Candidate*, RefCountedKinematicParticle>& km,
156  std::map<const BPHRecoCandidate*, RefCountedKinematicParticle>& cm) const;
157  virtual void getParticles(const std::string& moth,
158  const std::string& daug,
159  std::vector<RefCountedKinematicParticle>& kl,
160  std::set<RefCountedKinematicParticle>& ks) const;
161  virtual void getParticles(const std::string& moth,
162  const std::vector<std::string>& daug,
163  std::vector<RefCountedKinematicParticle>& kl,
164  std::set<RefCountedKinematicParticle>& ks) const;
165  virtual unsigned int numParticles(const BPHKinematicFit* cand = nullptr) const;
167  std::vector<RefCountedKinematicParticle>& kl,
168  std::set<RefCountedKinematicParticle>& ks);
169  virtual const BPHKinematicFit* splitKP(const std::string& name,
170  std::vector<RefCountedKinematicParticle>* kComp,
171  std::vector<RefCountedKinematicParticle>* kTail = nullptr) const;
172  virtual const RefCountedKinematicTree& kinematicTree(const std::vector<RefCountedKinematicParticle>& kPart,
174  virtual void fitMomentum() const;
175 };
176 
177 #endif
virtual const math::XYZTLorentzVector & p4() const
compute total momentum after the fit
virtual unsigned int numParticles(const BPHKinematicFit *cand=nullptr) const
virtual void buildParticles() const
virtual void addK(const std::string &name, const reco::Candidate *daug, double mass=-1.0, double sigma=-1.0)
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
std::map< const BPHRecoCandidate *, RefCountedKinematicParticle > kCDMap
int kp
double ParticleMass
Definition: ParticleMass.h:4
virtual bool isEmpty() const
get fit status
virtual ParticleMass mass() const
virtual const BPHKinematicFit * splitKP(const std::string &name, std::vector< RefCountedKinematicParticle > *kComp, std::vector< RefCountedKinematicParticle > *kTail=nullptr) const
virtual void getParticles(const std::string &moth, const std::string &daug, std::vector< RefCountedKinematicParticle > &kl, std::set< RefCountedKinematicParticle > &ks) const
const std::string names[nVars_]
virtual const RefCountedKinematicParticle currentParticle() const
get current particle
std::vector< BPHRecoConstCandPtr > tmpList
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double constrSigma() const
std::map< const reco::Candidate *, RefCountedKinematicParticle > kinMap
double constrMass() const
retrieve the constraint
void setNotUpdated() const override
virtual void fitMomentum() const
virtual const RefCountedKinematicTree & kinematicTree() const
perform the kinematic fit and get the result
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
std::map< const reco::Candidate *, double > dMSig
virtual const RefCountedKinematicVertex currentDecayVertex() const
BPHKinematicFit & operator=(const BPHKinematicFit &x)=delete
virtual const std::vector< RefCountedKinematicParticle > & kinParticles() const
get kinematic particles
static void insertParticle(RefCountedKinematicParticle &kp, std::vector< RefCountedKinematicParticle > &kl, std::set< RefCountedKinematicParticle > &ks)
virtual void resetKinematicFit() const
reset the kinematic fit
std::vector< RefCountedKinematicParticle > allParticles
virtual bool isValidFit() const
virtual void addParticles(std::vector< RefCountedKinematicParticle > &kl, std::map< const reco::Candidate *, RefCountedKinematicParticle > &km, std::map< const BPHRecoCandidate *, RefCountedKinematicParticle > &cm) const
RefCountedKinematicTree kinTree
std::map< const BPHRecoCandidate *, FlyingParticle > cKinP
virtual const RefCountedKinematicVertex topDecayVertex() const
void setConstraint(double mass, double sigma)
apply a mass constraint
math::XYZTLorentzVector totalMomentum
virtual const RefCountedKinematicParticle topParticle() const
get top particle
~BPHKinematicFit() override
bool getIndependentFit(const std::string &name) const
retrieve independent fit flag
double getMassSigma(const reco::Candidate *cand) const
retrieve particle mass sigma