Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 //
0002 //
0003 // File: src/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 // CMSSW File      : src/Lepjets_Event.cc
0008 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
0009 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
0010 //
0011 
0012 /**
0013     @file Lepjets_Event.cc
0014 
0015     @brief Represent a simple event consisting of lepton(s) and jet(s).
0016     See the documentation of header file Lepjets_Event.h for details.
0017 
0018     @author Scott Stuart Snyder <snyder@bnl.gov>
0019 
0020     @par Creation date:
0021     Jul 2000.
0022 
0023     @par Modification History:
0024     Apr 2009: Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>:
0025     Imported to CMSSW.<br>
0026     Nov 2009: Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>:
0027     Added doxygen tags for automatic generation of documentation.
0028 
0029     @par Terms of Usage:
0030     With consent for the original author (Scott Snyder).
0031 
0032  */
0033 
0034 #include "TopQuarkAnalysis/TopHitFit/interface/Lepjets_Event.h"
0035 #include <algorithm>
0036 #include <functional>
0037 #include <cmath>
0038 #include <cassert>
0039 
0040 using std::abs;
0041 using std::remove_if;
0042 using std::stable_sort;
0043 
0044 namespace hitfit {
0045 
0046   Lepjets_Event::Lepjets_Event(int runnum, int evnum)
0047       //
0048       // Purpose: Constructor.
0049       //
0050       // Inputs:
0051       //   runnum -      The run number.
0052       //   evnum -       The event number.
0053       //
0054       : _zvertex(0), _isMC(false), _runnum(runnum), _evnum(evnum), _dlb(-1), _dnn(-1) {}
0055 
0056   const int& Lepjets_Event::runnum() const
0057   //
0058   // Purpose: Return the run number.
0059   //
0060   // Returns:
0061   //   The run number.
0062   //
0063   {
0064     return _runnum;
0065   }
0066 
0067   int& Lepjets_Event::runnum()
0068   //
0069   // Purpose: Return the run number.
0070   //
0071   // Returns:
0072   //   The run number.
0073   //
0074   {
0075     return _runnum;
0076   }
0077 
0078   const int& Lepjets_Event::evnum() const
0079   //
0080   // Purpose: Return the event number.
0081   //
0082   // Returns:
0083   //   The event number.
0084   //
0085   {
0086     return _evnum;
0087   }
0088 
0089   int& Lepjets_Event::evnum()
0090   //
0091   // Purpose: Return the event number.
0092   //
0093   // Returns:
0094   //   The event number.
0095   //
0096   {
0097     return _evnum;
0098   }
0099 
0100   std::vector<Lepjets_Event_Lep>::size_type Lepjets_Event::nleps() const
0101   //
0102   // Purpose: Return the length of the lepton list.
0103   //
0104   // Returns:
0105   //   The length of the lepton list.
0106   //
0107   {
0108     return _leps.size();
0109   }
0110 
0111   std::vector<Lepjets_Event_Jet>::size_type Lepjets_Event::njets() const
0112   //
0113   // Purpose: Return the length of the jet list.
0114   //
0115   // Returns:
0116   //   The length of the jet list.
0117   //
0118   {
0119     return _jets.size();
0120   }
0121 
0122   Lepjets_Event_Lep& Lepjets_Event::lep(std::vector<Lepjets_Event_Lep>::size_type i)
0123   //
0124   // Purpose: Return the Ith lepton.
0125   //
0126   // Inputs:
0127   //   i -           The lepton index (counting from 0).
0128   //
0129   // Returns:
0130   //   The Ith lepton.
0131   //
0132   {
0133     assert(i < _leps.size());
0134     return _leps[i];
0135   }
0136 
0137   Lepjets_Event_Jet& Lepjets_Event::jet(std::vector<Lepjets_Event_Jet>::size_type i)
0138   //
0139   // Purpose: Return the Ith jet.
0140   //
0141   // Inputs:
0142   //   i -           The jet index (counting from 0).
0143   //
0144   // Returns:
0145   //   The Ith jet.
0146   //
0147   {
0148     assert(i < _jets.size());
0149     return _jets[i];
0150   }
0151 
0152   const Lepjets_Event_Lep& Lepjets_Event::lep(std::vector<Lepjets_Event_Lep>::size_type i) const
0153   //
0154   // Purpose: Return the Ith lepton.
0155   //
0156   // Inputs:
0157   //   i -           The lepton index (counting from 0).
0158   //
0159   // Returns:
0160   //   The Ith lepton.
0161   //
0162   {
0163     assert(i < _leps.size());
0164     return _leps[i];
0165   }
0166 
0167   const Lepjets_Event_Jet& Lepjets_Event::jet(std::vector<Lepjets_Event_Jet>::size_type i) const
0168   //
0169   // Purpose: Return the Ith jet.
0170   //
0171   // Inputs:
0172   //   i -           The jet index (counting from 0).
0173   //
0174   // Returns:
0175   //   The Ith jet.
0176   //
0177   {
0178     assert(i < _jets.size());
0179     return _jets[i];
0180   }
0181 
0182   Fourvec& Lepjets_Event::met()
0183   //
0184   // Purpose: Return the missing Et.
0185   //
0186   // Returns:
0187   //   The missing Et.
0188   //
0189   {
0190     return _met;
0191   }
0192 
0193   const Fourvec& Lepjets_Event::met() const
0194   //
0195   // Purpose: Return the missing Et.
0196   //
0197   // Returns:
0198   //   The missing Et.
0199   //
0200   {
0201     return _met;
0202   }
0203 
0204   Resolution& Lepjets_Event::kt_res()
0205   //
0206   // Purpose: Return the kt resolution.
0207   //
0208   // Returns:
0209   //   The kt resolution.
0210   //
0211   {
0212     return _kt_res;
0213   }
0214 
0215   const Resolution& Lepjets_Event::kt_res() const
0216   //
0217   // Purpose: Return the kt resolution.
0218   //
0219   // Returns:
0220   //   The kt resolution.
0221   //
0222   {
0223     return _kt_res;
0224   }
0225 
0226   double Lepjets_Event::zvertex() const
0227   //
0228   // Purpose: Return the z-vertex.
0229   //
0230   // Returns:
0231   //   The z-vertex.
0232   //
0233   {
0234     return _zvertex;
0235   }
0236 
0237   double& Lepjets_Event::zvertex()
0238   //
0239   // Purpose: Return the z-vertex.
0240   //
0241   // Returns:
0242   //   The z-vertex.
0243   //
0244   {
0245     return _zvertex;
0246   }
0247 
0248   bool Lepjets_Event::isMC() const
0249   //
0250   // Purpose: Return the isMC flag.
0251   //
0252   // Returns:
0253   //   The isMC flag.
0254   //
0255   {
0256     return _isMC;
0257   }
0258 
0259   void Lepjets_Event::setMC(bool isMC)
0260   //
0261   // Purpose: set isMC flag.
0262   //
0263   // Returns:
0264   //   nothing
0265   //
0266   {
0267     _isMC = isMC;
0268   }
0269 
0270   double Lepjets_Event::dlb() const
0271   //
0272   // Purpose: Return the LB discriminant.
0273   //
0274   // Returns:
0275   //   The LB discriminant.
0276   //
0277   {
0278     return _dlb;
0279   }
0280 
0281   double& Lepjets_Event::dlb()
0282   //
0283   // Purpose: Return the LB discriminant.
0284   //
0285   // Returns:
0286   //   The LB discriminant.
0287   //
0288   {
0289     return _dlb;
0290   }
0291 
0292   double Lepjets_Event::dnn() const
0293   //
0294   // Purpose: Return the NN discriminant.
0295   //
0296   // Returns:
0297   //   The NN discriminant.
0298   //
0299   {
0300     return _dnn;
0301   }
0302 
0303   double& Lepjets_Event::dnn()
0304   //
0305   // Purpose: Return the NN discriminant.
0306   //
0307   // Returns:
0308   //   The NN discriminant.
0309   //
0310   {
0311     return _dnn;
0312   }
0313 
0314   Fourvec Lepjets_Event::sum(int type) const
0315   //
0316   // Purpose: Sum all objects with type code TYPE.
0317   //
0318   // Inputs:
0319   //   type -        The type code to match.
0320   //
0321   // Returns:
0322   //   The sum of all objects with type code TYPE.
0323   //
0324   {
0325     Fourvec out;
0326     for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++)
0327       if (_leps[i].type() == type)
0328         out += _leps[i].p();
0329     for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++)
0330       if (_jets[i].type() == type)
0331         out += _jets[i].p();
0332     return out;
0333   }
0334 
0335   Fourvec Lepjets_Event::kt() const
0336   //
0337   // Purpose: Calculate kt --- sum of all objects plus missing Et.
0338   //
0339   // Returns:
0340   //   The event kt.
0341   {
0342     Fourvec v = _met;
0343     for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++)
0344       v += _leps[i].p();
0345     for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++)
0346       v += _jets[i].p();
0347     return v;
0348   }
0349 
0350   void Lepjets_Event::add_lep(const Lepjets_Event_Lep& lep)
0351   //
0352   // Purpose: Add a lepton to the event.
0353   //
0354   // Inputs:
0355   //   lep -         The lepton to add.
0356   //
0357   {
0358     _leps.push_back(lep);
0359   }
0360 
0361   void Lepjets_Event::add_jet(const Lepjets_Event_Jet& jet)
0362   //
0363   // Purpose: Add a jet to the event.
0364   //
0365   // Inputs:
0366   //   jet -         The jet to add.
0367   //
0368   {
0369     _jets.push_back(jet);
0370   }
0371 
0372   void Lepjets_Event::smear(CLHEP::HepRandomEngine& engine, bool smear_dir /*= false*/)
0373   //
0374   // Purpose: Smear the objects in the event according to their resolutions.
0375   //
0376   // Inputs:
0377   //   engine -      The underlying RNG.
0378   //   smear_dir -   If false, smear the momentum only.
0379   //
0380   {
0381     Fourvec before, after;
0382     for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++) {
0383       before += _leps[i].p();
0384       _leps[i].smear(engine, smear_dir);
0385       after += _leps[i].p();
0386     }
0387     for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++) {
0388       before += _jets[i].p();
0389       _jets[i].smear(engine, smear_dir);
0390       after += _jets[i].p();
0391     }
0392 
0393     Fourvec kt = _met + before;
0394     kt(Fourvec::X) = _kt_res.pick(kt(Fourvec::X), kt(Fourvec::X), engine);
0395     kt(Fourvec::Y) = _kt_res.pick(kt(Fourvec::Y), kt(Fourvec::Y), engine);
0396     _met = kt - after;
0397   }
0398 
0399   void Lepjets_Event::sort()
0400   //
0401   // Purpose: Sort the objects in the event in order of descending pt.
0402   //
0403   {
0404     std::stable_sort(_leps.begin(), _leps.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
0405     std::stable_sort(_jets.begin(), _jets.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
0406   }
0407 
0408   std::vector<int> Lepjets_Event::jet_types() const
0409   //
0410   // Purpose: Return the jet types of the event
0411   //
0412   {
0413     std::vector<int> ret;
0414     for (std::vector<Lepjets_Event_Jet>::size_type ijet = 0; ijet != _jets.size(); ijet++) {
0415       ret.push_back(jet(ijet).type());
0416     }
0417     return ret;
0418   }
0419 
0420   bool Lepjets_Event::set_jet_types(const std::vector<int>& _jet_types)
0421   //
0422   // Purpose: Set the jet types of the event
0423   // Return false if it fails, trus if it succeeds
0424   //
0425   {
0426     if (_jets.size() != _jet_types.size()) {
0427       return false;
0428     }
0429     bool saw_hadw1 = false;
0430     for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < njets(); i++) {
0431       int t = _jet_types[i];
0432       if (t == hadw1_label) {
0433         if (saw_hadw1)
0434           t = hadw2_label;
0435         saw_hadw1 = true;
0436       }
0437       jet(i).type() = t;
0438     }
0439     return true;
0440   }
0441 
0442   namespace {
0443 
0444     struct Lepjets_Event_Cutter
0445     //
0446     // Purpose: Helper for cutting on objects.
0447     //
0448     {
0449       Lepjets_Event_Cutter(double pt_cut, double eta_cut) : _pt_cut(pt_cut), _eta_cut(eta_cut) {}
0450       bool operator()(const Lepjets_Event_Lep& l) const;
0451       double _pt_cut;
0452       double _eta_cut;
0453     };
0454 
0455     bool Lepjets_Event_Cutter::operator()(const Lepjets_Event_Lep& l) const
0456     //
0457     // Purpose: Object cut evaluator.
0458     //
0459     {
0460       return !(l.p().perp() > _pt_cut && abs(l.p().pseudoRapidity()) < _eta_cut);
0461     }
0462 
0463   }  // unnamed namespace
0464 
0465   int Lepjets_Event::cut_leps(double pt_cut, double eta_cut)
0466   //
0467   // Purpose: Remove all leptons failing the pt and eta cuts.
0468   //
0469   // Inputs:
0470   //   pt_cut -      Pt cut.  Remove objects with pt less than this.
0471   //   eta_cut -     Eta cut.  Remove objects with abs(eta) larger than this.
0472   //
0473   // Returns:
0474   //   The number of leptons remaining after the cuts.
0475   //
0476   {
0477     _leps.erase(remove_if(_leps.begin(), _leps.end(), Lepjets_Event_Cutter(pt_cut, eta_cut)), _leps.end());
0478     return _leps.size();
0479   }
0480 
0481   int Lepjets_Event::cut_jets(double pt_cut, double eta_cut)
0482   //
0483   // Purpose: Remove all jets failing the pt and eta cuts.
0484   //
0485   // Inputs:
0486   //   pt_cut -      Pt cut.  Remove objects with pt less than this.
0487   //   eta_cut -     Eta cut.  Remove objects with abs(eta) larger than this.
0488   //
0489   // Returns:
0490   //   The number of jets remaining after the cuts.
0491   //
0492   {
0493     _jets.erase(remove_if(_jets.begin(), _jets.end(), Lepjets_Event_Cutter(pt_cut, eta_cut)), _jets.end());
0494     return _jets.size();
0495   }
0496 
0497   void Lepjets_Event::trimjets(std::vector<Lepjets_Event_Jet>::size_type n)
0498   //
0499   // Purpose: Remove all but the first N jets.
0500   //
0501   // Inputs:
0502   //   n -           The number of jets to keep.
0503   //
0504   {
0505     if (n >= _jets.size())
0506       return;
0507     _jets.erase(_jets.begin() + n, _jets.end());
0508   }
0509 
0510   std::ostream& Lepjets_Event::dump(std::ostream& s, bool full /*=false*/) const
0511   //
0512   // Purpose: Dump out this object.
0513   //
0514   // Inputs:
0515   //   s -           The stream to which to write.
0516   //   full -        If true, dump all information for this object.
0517   //
0518   // Returns:
0519   //   The stream S.
0520   //
0521   {
0522     s << "Run: " << _runnum << "  Event: " << _evnum << "\n";
0523     s << "Leptons:\n";
0524     for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++) {
0525       s << "  ";
0526       _leps[i].dump(s, full);
0527       s << "\n";
0528     }
0529     s << "Jets:\n";
0530     for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++) {
0531       s << "  ";
0532       _jets[i].dump(s, full);
0533       s << "\n";
0534     }
0535     s << "Missing Et: " << _met << "\n";
0536     if (_zvertex != 0)
0537       s << "z-vertex: " << _zvertex << "\n";
0538     return s;
0539   }
0540 
0541   std::string Lepjets_Event::jet_permutation() const
0542   //
0543   // Purpose: Return a string representation of the jet permutation
0544   //          g - isr/gluon
0545   //          b - leptonic b
0546   //          B - hadronic b
0547   //          w - hadronic W
0548   //          h - Higgs to b-bbar
0549   //          ? - Unknown
0550   {
0551     std::string permutation;
0552     for (size_t jet = 0; jet != _jets.size(); ++jet) {
0553       permutation = permutation + hitfit::jetTypeString(_jets[jet].type());
0554     }
0555     return permutation;
0556   }
0557 
0558   /**
0559     @brief Output stream operator, print the content of this Lepjets_Event
0560     to an output stream.
0561 
0562     @param s The output stream to which to write.
0563 
0564     @param ev The instance of Lepjets_Event to be printed.
0565  */
0566   std::ostream& operator<<(std::ostream& s, const Lepjets_Event& ev)
0567   //
0568   // Purpose: Dump out this object.
0569   //
0570   // Inputs:
0571   //   s -           The stream to which to write.
0572   //   ev -          The object to dump.
0573   //
0574   // Returns:
0575   //   The stream S.
0576   //
0577   {
0578     return ev.dump(s);
0579   }
0580 
0581 }  // namespace hitfit