Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:17

0001 //
0002 //
0003 // File: hitfit/Lepjets_Event.h
0004 // Purpose: Represent a simple `event' consisting of leptons and jets.
0005 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
0006 //
0007 // An instance of this class holds a list of `leptons' (as represented
0008 // by Lepjets_Event_Lep) and `jets' (as represented by Lepjets_Event_Jet).
0009 // We also record:
0010 //
0011 //   - Missing Et
0012 //   - z-vertex
0013 //   - run and event number
0014 //
0015 // CMSSW File      : interface/Lepjets_Event.h
0016 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
0017 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
0018 //
0019 
0020 /**
0021     @file Lepjets_Event.h
0022 
0023     @brief Represent a simple event consisting of lepton(s) and jet(s).
0024 
0025     @author Scott Stuart Snyder <snyder@bnl.gov>
0026 
0027     @par Creation date:
0028     Jul 2000.
0029 
0030     @par Modification History:
0031     Apr 2009: Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>:
0032     Imported to CMSSW.<br>
0033     Nov 2009: Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>:
0034     Added doxygen tags for automatic generation of documentation.
0035 
0036     @par Terms of Usage:
0037     With consent for the original author (Scott Snyder).
0038 
0039  */
0040 
0041 #ifndef HITFIT_LEPJETS_EVENT_H
0042 #define HITFIT_LEPJETS_EVENT_H
0043 
0044 #include "TopQuarkAnalysis/TopHitFit/interface/Lepjets_Event_Jet.h"
0045 #include "TopQuarkAnalysis/TopHitFit/interface/Lepjets_Event_Lep.h"
0046 #include <vector>
0047 #include <iosfwd>
0048 
0049 namespace hitfit {
0050 
0051   /**
0052     @class Lepjets_Event
0053 
0054     @brief Represent a simple event consisting of lepton(s) and jet(s).
0055     An instance of this class holds a list of leptons (as represented by
0056     the Lepjets_Event_Lep class) and jets (as represented by Lepjets_Event_Jet
0057     class).  Also recorded are:
0058     - Missing transverse energy.
0059     -  \f$ z- \f$ vertex (Irrelevant for non-D0 experiment)
0060     - Run and event number (Irrelevant for non-D0 experiment)
0061  */
0062   class Lepjets_Event
0063   //
0064   // Purpose: Represent a simple `event' consisting of leptons and jets.
0065   //
0066   {
0067   public:
0068     // Constructor.
0069     /**
0070       @brief Constructor.
0071 
0072       @param runnum The run number.
0073 
0074       @param evnum The event number.
0075    */
0076     Lepjets_Event(int runnum, int evnum);
0077 
0078     // Get the run and event number.
0079 
0080     /**
0081      @brief Return a constant reference to the run number.
0082    */
0083     const int& runnum() const;
0084 
0085     /**
0086      @brief Return a reference to the run number.
0087    */
0088     int& runnum();
0089 
0090     /**
0091      @brief Return a constant reference to the event number.
0092    */
0093     const int& evnum() const;
0094 
0095     /**
0096      @brief Return a reference to the event number.
0097    */
0098     int& evnum();
0099 
0100     // Get the length of the lepton and jet lists.
0101     /**
0102      @brief Return the number of leptons in the event.
0103    */
0104     std::vector<Lepjets_Event_Lep>::size_type nleps() const;
0105 
0106     /**
0107      @brief Return the number of jets in the event.
0108    */
0109     std::vector<Lepjets_Event_Jet>::size_type njets() const;
0110 
0111     // Access leptons and jets.
0112 
0113     /**
0114      @brief Return a reference to lepton at index position <i>i</i>.
0115 
0116      @param i The lepton index position.
0117    */
0118     Lepjets_Event_Lep& lep(std::vector<Lepjets_Event_Lep>::size_type i);
0119 
0120     /**
0121      @brief Return a reference to jet at index position <i>i</i>.
0122 
0123      @param i The jet index position.
0124    */
0125     Lepjets_Event_Jet& jet(std::vector<Lepjets_Event_Jet>::size_type i);
0126 
0127     /**
0128      @brief Return a constant reference to lepton at index position <i>i</i>.
0129 
0130      @param i The lepton index position.
0131    */
0132     const Lepjets_Event_Lep& lep(std::vector<Lepjets_Event_Lep>::size_type i) const;
0133 
0134     /**
0135      @brief Return a constant reference to jet at index position <i>i</i>.
0136 
0137      @param i The jet index position.
0138    */
0139     const Lepjets_Event_Jet& jet(std::vector<Lepjets_Event_Jet>::size_type i) const;
0140 
0141     // Access missing Et.
0142     /**
0143      @brief Return a reference to the missing transverse energy.
0144    */
0145     Fourvec& met();
0146 
0147     /**
0148      @brief Return a constant reference to the missing transverse energy.
0149    */
0150     const Fourvec& met() const;
0151 
0152     // Access kt resolution.
0153 
0154     /**
0155      @brief Return a reference to the  \f$ k_{T} \f$  resolution.
0156    */
0157     Resolution& kt_res();
0158 
0159     /**
0160      @brief Return a const reference to the  \f$ k_{T} \f$  resolution.
0161    */
0162     const Resolution& kt_res() const;
0163 
0164     // Access the z-vertex.
0165 
0166     /**
0167      @brief Return the value of z-vertex.
0168    */
0169     double zvertex() const;
0170 
0171     /**
0172      @brief Return a reference to the value of z-vertex.
0173    */
0174     double& zvertex();
0175 
0176     // Access the isMC flag.
0177     /**
0178      @brief Return the Monte Carlo flag.
0179    */
0180     bool isMC() const;
0181 
0182     /**
0183      @brief Set the Monte Carlo flag.
0184    */
0185     void setMC(bool isMC);
0186 
0187     // Access the discriminants.
0188     /**
0189      @brief Return a reference to the value of low-bias (LB) discriminant
0190      (Irrelevant for non-D0 experiment).
0191    */
0192     double& dlb();
0193 
0194     /**
0195      @brief Return the value of low-bias (LB) discriminant
0196      (Irrelevant for non-D0 experiment).
0197    */
0198     double dlb() const;
0199 
0200     /**
0201      @brief Return a reference to the value of neural network (NN) discriminant
0202      (Irrelevant for non-D0 experiment).
0203    */
0204     double& dnn();
0205 
0206     /**
0207      @brief Return a the value of neural network (NN) discriminant
0208      (Irrelevant for non-D0 experiment).
0209    */
0210     double dnn() const;
0211 
0212     // Sum all objects (leptons or jets) with type TYPE.
0213     /**
0214      @brief Return the sum of all objects' four-momentum which have
0215      a particular type.
0216 
0217      @param type The type code of the objects to be summed up.
0218    */
0219     Fourvec sum(int type) const;
0220 
0221     // Calculate kt --- sum of all objects plus missing Et.
0222     /**
0223      @brief Return the sum of all objects' four-momentum and
0224      missing transverse energy.
0225    */
0226     Fourvec kt() const;
0227 
0228     // Add new objects to the event.
0229     /**
0230      @brief Add a new lepton to the event.
0231 
0232      @param lep The lepton to be added.
0233    */
0234     void add_lep(const Lepjets_Event_Lep& lep);
0235 
0236     /**
0237      @brief Add a new jet to the event.
0238 
0239      @param jet The jet to be added.
0240    */
0241     void add_jet(const Lepjets_Event_Jet& jet);
0242 
0243     // Smear the objects in the event according to their resolutions.
0244     /**
0245      @brief Smear the objects in the event according to their resolutions.
0246 
0247      @param engine The underlying random number generator.
0248 
0249      @param smear_dir If <b>TRUE</b>, also smear the object's direction.<br>
0250      If <b>FALSE</b>, then only smear the magnitude of three-momentum.
0251    */
0252     void smear(CLHEP::HepRandomEngine& engine, bool smear_dir = false);
0253 
0254     // Sort according to pt.
0255     /**
0256      @brief Sort objects in the event according to their transverse momentum
0257       \f$ p_{T} \f$ .
0258    */
0259     void sort();
0260 
0261     // Get jet types
0262     /**
0263      @brief Return the jet types in the event.
0264    */
0265     std::vector<int> jet_types() const;
0266 
0267     // Set jet types
0268     /**
0269      @brief Set the jet types in the event.
0270    */
0271     bool set_jet_types(const std::vector<int>&);
0272 
0273     // Remove objects failing pt and eta cuts.
0274     /**
0275      @brief Remove leptons which fail transverse momentum  \f$ p_{T} \f$ 
0276      and pseudorapidity  \f$ \eta \f$  cut.
0277 
0278      @param pt_cut Remove leptons which have transverse momentum  \f$ p_{T} \f$ 
0279      less than this value, in GeV.
0280 
0281      @param eta_cut Remove leptons which have absolute
0282      pseudorapidity  \f$ |\eta| \f$  more than this value.
0283    */
0284     int cut_leps(double pt_cut, double eta_cut);
0285 
0286     /**
0287      @brief Remove jets which fail transverse momentum  \f$ p_{T} \f$
0288      and pseudorapidity  \f$ \eta \f$  cut.
0289 
0290      @param pt_cut Remove jets which have transverse momentum  \f$ p_{T} \f$
0291      less than this value, in GeV.
0292 
0293      @param eta_cut Remove jetss which have absolute
0294      pseudorapidity  \f$ |\eta| \f$  more than this value.
0295    */
0296     int cut_jets(double pt_cut, double eta_cut);
0297 
0298     // Remove all but the first N jets.
0299     /**
0300      @brief Remove all but the first <i>n</i> jets.
0301 
0302      @param n The number of jets to keep.
0303    */
0304     void trimjets(std::vector<Lepjets_Event_Jet>::size_type n);
0305 
0306     // Dump this object.
0307     /**
0308      @brief Print the content of this object.
0309 
0310      @param s The output stream to which to write
0311 
0312      @param full If <b>TRUE</b>, print all information about this instance
0313      of Lepjets_Event.<br>
0314      If <b>FALSE</b>, print partial information about this instance
0315      of Lepjets_Event.
0316    */
0317     std::ostream& dump(std::ostream& s, bool full = false) const;
0318 
0319     /**
0320      @brief Return a string representing the jet permutation. The following
0321      notation is used for each type of jet:
0322      - g ISR/gluon.
0323      - b leptonic  \f$ b- \f$ quark.
0324      - B hadronic  \f$ b- \f$ quark.
0325      - w hadronic jet from  \f$ W- \f$ boson.
0326      - H  \f$ b- \f$ jet from Higgs boson.
0327      - ? Unknown.
0328    */
0329     std::string jet_permutation() const;
0330 
0331   private:
0332     // The lepton and jet lists.
0333 
0334     /**
0335      The list of leptons in the event.
0336    */
0337     std::vector<Lepjets_Event_Lep> _leps;
0338 
0339     /**
0340      The list of jets in the event.
0341    */
0342     std::vector<Lepjets_Event_Jet> _jets;
0343 
0344     // Other event state.
0345     /**
0346      Missing transverse energy.
0347    */
0348     Fourvec _met;
0349 
0350     /**
0351      The  \f$ k_{T} \f$  resolution.
0352    */
0353     Resolution _kt_res;
0354 
0355     /**
0356      The  \f$ z- \f$ vertex of the event.
0357    */
0358     double _zvertex;
0359 
0360     /**
0361      The Monte Calro flag.
0362    */
0363     bool _isMC;
0364 
0365     /**
0366      The run number.
0367    */
0368     int _runnum;
0369 
0370     /**
0371      The event number.
0372    */
0373     int _evnum;
0374 
0375     /**
0376      The low-bias (LB) discriminant.
0377    */
0378     double _dlb;
0379 
0380     /**
0381      The neural network (NN) discriminant.
0382    */
0383     double _dnn;
0384   };
0385 
0386   // Print the object.
0387   std::ostream& operator<<(std::ostream& s, const Lepjets_Event& ev);
0388 
0389 }  // namespace hitfit
0390 
0391 #endif  // not HITFIT_LEPJETS_EVENT_H