00001 // 00002 // $Id: Lepjets_Event.h,v 1.1 2011/05/26 09:46:53 mseidel Exp $ 00003 // 00004 // File: hitfit/Lepjets_Event.h 00005 // Purpose: Represent a simple `event' consisting of leptons and jets. 00006 // Created: Jul, 2000, sss, based on run 1 mass analysis code. 00007 // 00008 // An instance of this class holds a list of `leptons' (as represented 00009 // by Lepjets_Event_Lep) and `jets' (as represented by Lepjets_Event_Jet). 00010 // We also record: 00011 // 00012 // - Missing Et 00013 // - z-vertex 00014 // - run and event number 00015 // 00016 // CMSSW File : interface/Lepjets_Event.h 00017 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0 00018 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch> 00019 // 00020 00021 00043 #ifndef HITFIT_LEPJETS_EVENT_H 00044 #define HITFIT_LEPJETS_EVENT_H 00045 00046 00047 #include "TopQuarkAnalysis/TopHitFit/interface/Lepjets_Event_Jet.h" 00048 #include "TopQuarkAnalysis/TopHitFit/interface/Lepjets_Event_Lep.h" 00049 #include <vector> 00050 #include <iosfwd> 00051 00052 00053 namespace hitfit { 00054 00055 00067 class Lepjets_Event 00068 // 00069 // Purpose: Represent a simple `event' consisting of leptons and jets. 00070 // 00071 { 00072 public: 00073 00074 // Constructor. 00082 Lepjets_Event (int runnum, int evnum); 00083 00084 // Get the run and event number. 00085 00089 const int& runnum () const; 00090 00094 int& runnum(); 00095 00099 const int& evnum () const; 00100 00104 int& evnum(); 00105 00106 // Get the length of the lepton and jet lists. 00110 std::vector<Lepjets_Event_Lep>::size_type nleps () const; 00111 00115 std::vector<Lepjets_Event_Jet>::size_type njets () const; 00116 00117 // Access leptons and jets. 00118 00124 Lepjets_Event_Lep& lep (std::vector<Lepjets_Event_Lep>::size_type i); 00125 00131 Lepjets_Event_Jet& jet (std::vector<Lepjets_Event_Jet>::size_type i); 00132 00138 const Lepjets_Event_Lep& lep (std::vector<Lepjets_Event_Lep>::size_type i) const; 00139 00145 const Lepjets_Event_Jet& jet (std::vector<Lepjets_Event_Jet>::size_type i) const; 00146 00147 // Access missing Et. 00151 Fourvec& met (); 00152 00156 const Fourvec& met () const; 00157 00158 // Access kt resolution. 00159 00163 Resolution& kt_res (); 00164 00168 const Resolution& kt_res () const; 00169 00170 // Access the z-vertex. 00171 00175 double zvertex () const; 00176 00180 double& zvertex (); 00181 00182 // Access the isMC flag. 00186 bool isMC () const; 00187 00191 void setMC (bool isMC); 00192 00193 // Access the discriminants. 00198 double& dlb (); 00199 00204 double dlb () const; 00205 00210 double& dnn (); 00211 00216 double dnn () const; 00217 00218 // Sum all objects (leptons or jets) with type TYPE. 00225 Fourvec sum (int type) const; 00226 00227 // Calculate kt --- sum of all objects plus missing Et. 00232 Fourvec kt () const; 00233 00234 // Add new objects to the event. 00240 void add_lep (const Lepjets_Event_Lep& lep); 00241 00247 void add_jet (const Lepjets_Event_Jet& jet); 00248 00249 // Smear the objects in the event according to their resolutions. 00258 void smear (CLHEP::HepRandomEngine& engine, bool smear_dir = false); 00259 00260 // Sort according to pt. 00265 void sort (); 00266 00267 // Get jet types 00271 std::vector<int> jet_types() const; 00272 00273 // Set jet types 00277 bool set_jet_types(const std::vector<int>&); 00278 00279 // Remove objects failing pt and eta cuts. 00290 int cut_leps (double pt_cut, double eta_cut); 00291 00302 int cut_jets (double pt_cut, double eta_cut); 00303 00304 // Remove all but the first N jets. 00310 void trimjets (std::vector<Lepjets_Event_Jet>::size_type n); 00311 00312 // Dump this object. 00323 std::ostream& dump (std::ostream& s, bool full = false) const; 00324 00335 std::string jet_permutation() const; 00336 00337 private: 00338 // The lepton and jet lists. 00339 00343 std::vector<Lepjets_Event_Lep> _leps; 00344 00348 std::vector<Lepjets_Event_Jet> _jets; 00349 00350 // Other event state. 00354 Fourvec _met; 00355 00359 Resolution _kt_res; 00360 00364 double _zvertex; 00365 00369 bool _isMC; 00370 00374 int _runnum; 00375 00379 int _evnum; 00380 00384 double _dlb; 00385 00389 double _dnn; 00390 }; 00391 00392 00393 // Print the object. 00394 std::ostream& operator<< (std::ostream& s, const Lepjets_Event& ev); 00395 00396 00397 00398 } // namespace hitfit 00399 00400 00401 #endif // not HITFIT_LEPJETS_EVENT_H 00402