CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATDiObjectProxy.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatUtils_interface_PATDiObjectProxy_h
2 #define PhysicsTools_PatUtils_interface_PATDiObjectProxy_h
3 
6 
15 
16 namespace pat {
17 
18 /* Now we implement PATDiObjectProxy with typeid & static_casts on the fly. */
20 
21  public:
27  cand1_(&c1), cand2_(&c2), type1_(&typeid(c1)), type2_(&typeid(c2)), totalP4ok_(false), totalP4_() {}
28 
30  const reco::Candidate & cand1() const { return *cand1_; }
32  const reco::Candidate & cand2() const { return *cand2_; }
33 
35  double deltaR() const { return ::deltaR(*cand1_, *cand2_); }
37  double deltaPhi() const { return ::deltaPhi(cand1_->phi(), cand2_->phi()); }
38 
40  // Implementation notice: return by const reference, not by value,
41  // as it's easier for Reflex.
43  if (!totalP4ok_) {
44  totalP4_ = cand1_->p4() + cand2_->p4();
45  totalP4ok_ = true;
46  }
47  return totalP4_;
48  }
49 
51  const Electron & ele() const { return tryGetOne_<Electron>(); }
53  const Muon & mu() const { return tryGetOne_<Muon>(); }
55  const Tau & tau() const { return tryGetOne_<Tau>(); }
57  const Photon & gam() const { return tryGetOne_<Photon>(); }
59  const Jet & jet() const { return tryGetOne_<Jet>(); }
61  const MET & met() const { return tryGetOne_<MET>(); }
63  const GenericParticle & part() const { return tryGetOne_<GenericParticle>(); }
65  const PFParticle & pf() const { return tryGetOne_<PFParticle>(); }
66 
68  const Electron & ele1() const { return tryGet_<Electron>(cand1_, type1_); }
70  const Muon & mu1() const { return tryGet_<Muon>(cand1_, type1_); }
72  const Tau & tau1() const { return tryGet_<Tau>(cand1_, type1_); }
74  const Photon & gam1() const { return tryGet_<Photon>(cand1_, type1_); }
76  const Jet & jet1() const { return tryGet_<Jet>(cand1_, type1_); }
78  const MET & met1() const { return tryGet_<MET>(cand1_, type1_); }
80  const GenericParticle & part1() const { return tryGet_<GenericParticle>(cand1_, type1_); }
82  const PFParticle & pf1() const { return tryGet_<PFParticle>(cand1_, type1_); }
83 
85  const Electron & ele2() const { return tryGet_<Electron>(cand2_, type2_); }
87  const Muon & mu2() const { return tryGet_<Muon>(cand2_, type2_); }
89  const Tau & tau2() const { return tryGet_<Tau>(cand2_, type2_); }
91  const Photon & gam2() const { return tryGet_<Photon>(cand2_, type2_); }
93  const Jet & jet2() const { return tryGet_<Jet>(cand2_, type2_); }
95  const MET & met2() const { return tryGet_<MET>(cand2_, type2_); }
97  const GenericParticle & part2() const { return tryGet_<GenericParticle>(cand2_, type2_); }
99  const PFParticle & pf2() const { return tryGet_<PFParticle>(cand2_, type2_); }
100 
101  private:
102 
103  template<typename T>
104  const T & tryGet_(const reco::Candidate *ptr, const std::type_info *type) const {
105  if (typeid(T) != *type) {
106  throw cms::Exception("Type Error") << "pat::DiObjectProxy: the object of the pair is not of the type you request.\n"
107  << " Item Index in pair: " << (ptr == cand1_ ? "first" : "second") << "\n"
108  << " Requested TypeID : " << ClassName<T>::name() << "\n"
109  << " Found TypeID : " << className(*ptr) << "\n";
110  }
111  return static_cast<const T &>(*ptr);
112  }
113 
114  template<typename T>
115  const T & tryGetOne_() const {
116  if (typeid(T) == *type1_) {
117  if (typeid(T) == *type2_) {
118  throw cms::Exception("Type Error") << "pat::DiObjectProxy: " <<
119  "you can't get use methods that get a particle by type if the two are of the same type!\n" <<
120  " Requested Type:" << ClassName<T>::name() << "\n";
121  }
122  return static_cast<const T &>(*cand1_);
123  } else {
124  if (typeid(T) != *type2_) {
125  throw cms::Exception("Type Error") << "pat::DiObjectProxy: " <<
126  "you can't get use methods that get a particle by type if neither of the two is of that type!\n" <<
127  " Requested Type:" << ClassName<T>::name() << "\n" <<
128  " Type of first :" << className(*cand1_) << "\n" <<
129  " Type of second:" << className(*cand2_) << "\n";
130  }
131  return static_cast<const T &>(*cand2_);
132  }
133  }
134 
136  const std::type_info *type1_, *type2_;
137 
138  mutable bool totalP4ok_;
140 
141 
142 };
143 
144 }
145 #endif
type
Definition: HCALResponse.h:22
Analysis-level MET class.
Definition: MET.h:42
Analysis-level Photon class.
Definition: Photon.h:46
double deltaR() const
Get the angular separation.
const GenericParticle & part1() const
Get the first item, if it&#39;s a PAT GenericParticle (throw exception otherwise)
const reco::Candidate & cand2() const
Gets the second Candidate.
const PFParticle & pf1() const
Get the first item, if it&#39;s a PAT PFParticle (throw exception otherwise)
const reco::Candidate & cand1() const
Gets the first Candidate.
const reco::Candidate * cand2_
const PFParticle & pf2() const
Get the second item, if it&#39;s a PAT PFParticle (throw exception otherwise)
const reco::Candidate * cand1_
const Electron & ele() const
Get the PAT Electron, if the pair contains one and only one PAT Electron (throw exception otherwise) ...
const GenericParticle & part2() const
Get the second item, if it&#39;s a PAT GenericParticle (throw exception otherwise)
const std::type_info * type1_
const Muon & mu2() const
Get the second item, if it&#39;s a PAT Muon (throw exception otherwise)
const GenericParticle & part() const
Get the PAT GenericParticle, if the pair contains one and only one PAT GenericParticle (throw excepti...
const Photon & gam1() const
Get the first item, if it&#39;s a PAT Photon (throw exception otherwise)
Analysis-level Generic Particle class (e.g. for hadron or muon not fully reconstructed) ...
const Electron & ele1() const
Get the first item, if it&#39;s a PAT Electron (throw exception otherwise)
const Tau & tau2() const
Get the second item, if it&#39;s a PAT Tau (throw exception otherwise)
const Muon & mu() const
Get the PAT Muon, if the pair contains one and only one PAT Muon (throw exception otherwise) ...
const T & tryGetOne_() const
static const std::string & name()
Definition: ClassName.h:37
const Electron & ele2() const
Get the second item, if it&#39;s a PAT Electron (throw exception otherwise)
DiObjectProxy(const reco::Candidate &c1, const reco::Candidate &c2)
const MET & met2() const
Get the second item, if it&#39;s a PAT MET (throw exception otherwise)
const MET & met1() const
Get the first item, if it&#39;s a PAT MET (throw exception otherwise)
Analysis-level tau class.
Definition: Tau.h:50
const Photon & gam2() const
Get the second item, if it&#39;s a PAT Photon (throw exception otherwise)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
const Jet & jet2() const
Get the second item, if it&#39;s a PAT Jet (throw exception otherwise)
const std::type_info * type2_
const T & tryGet_(const reco::Candidate *ptr, const std::type_info *type) const
const Jet & jet1() const
Get the first item, if it&#39;s a PAT Jet (throw exception otherwise)
const reco::Candidate::LorentzVector & totalP4() const
Get the total four momentum.
const Jet & jet() const
Get the PAT Jet, if the pair contains one and only one PAT Jet (throw exception otherwise) ...
const MET & met() const
Get the PAT MET, if the pair contains one and only one PAT MET (throw exception otherwise) ...
reco::Candidate::LorentzVector totalP4_
Analysis-level electron class.
Definition: Electron.h:52
Analysis-level calorimeter jet class.
Definition: Jet.h:71
const Tau & tau() const
Get the PAT Tau, if the pair contains one and only one PAT Tau (throw exception otherwise) ...
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
Analysis-level class for reconstructed particles.
Definition: PFParticle.h:37
double deltaPhi() const
Get the phi separation.
const Photon & gam() const
Get the PAT Photon, if the pair contains one and only one PAT Photon (throw exception otherwise) ...
const Muon & mu1() const
Get the first item, if it&#39;s a PAT Muon (throw exception otherwise)
long double T
Analysis-level muon class.
Definition: Muon.h:51
const Tau & tau1() const
Get the first item, if it&#39;s a PAT Tau (throw exception otherwise)
std::string className(const T &t)
Definition: ClassName.h:30
DiObjectProxy()
Default constructor, requested by ROOT. NEVER use a default constructed item!
virtual double phi() const =0
momentum azimuthal angle
const PFParticle & pf() const
Get the PAT PFParticle, if the pair contains one and only one PAT PFParticle (throw exception otherwi...
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector