CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/TopQuarkAnalysis/TopHitFit/src/Lepjets_Event.cc

Go to the documentation of this file.
00001 //
00002 // $Id: Lepjets_Event.cc,v 1.1 2011/05/26 09:47:00 mseidel Exp $
00003 //
00004 // File: src/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 // CMSSW File      : src/Lepjets_Event.cc
00009 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
00010 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
00011 //
00012 
00013 
00036 #include "TopQuarkAnalysis/TopHitFit/interface/Lepjets_Event.h"
00037 #include <algorithm>
00038 #include <functional>
00039 #include <cmath>
00040 #include <cassert>
00041 
00042 
00043 using std::stable_sort;
00044 using std::less;
00045 using std::not2;
00046 using std::remove_if;
00047 using std::abs;
00048 
00049 
00050 namespace hitfit {
00051 
00052 
00053 Lepjets_Event::Lepjets_Event (int runnum, int evnum)
00054 //
00055 // Purpose: Constructor.
00056 //
00057 // Inputs:
00058 //   runnum -      The run number.
00059 //   evnum -       The event number.
00060 //
00061   : _zvertex (0),
00062     _isMC(false),
00063     _runnum (runnum),
00064     _evnum (evnum),
00065     _dlb (-1),
00066     _dnn (-1)
00067 {
00068 }
00069 
00070 
00071 const int& Lepjets_Event::runnum () const
00072 //
00073 // Purpose: Return the run number.
00074 //
00075 // Returns:
00076 //   The run number.
00077 //
00078 {
00079   return _runnum;
00080 }
00081 
00082 
00083 int& Lepjets_Event::runnum ()
00084 //
00085 // Purpose: Return the run number.
00086 //
00087 // Returns:
00088 //   The run number.
00089 //
00090 {
00091   return _runnum;
00092 }
00093 
00094 
00095 const int& Lepjets_Event::evnum () const
00096 //
00097 // Purpose: Return the event number.
00098 //
00099 // Returns:
00100 //   The event number.
00101 //
00102 {
00103   return _evnum;
00104 }
00105 
00106 
00107 int& Lepjets_Event::evnum ()
00108 //
00109 // Purpose: Return the event number.
00110 //
00111 // Returns:
00112 //   The event number.
00113 //
00114 {
00115   return _evnum;
00116 }
00117 
00118 
00119 std::vector<Lepjets_Event_Lep>::size_type Lepjets_Event::nleps () const
00120 //
00121 // Purpose: Return the length of the lepton list.
00122 //
00123 // Returns:
00124 //   The length of the lepton list.
00125 //
00126 {
00127   return _leps.size ();
00128 }
00129 
00130 
00131 std::vector<Lepjets_Event_Jet>::size_type Lepjets_Event::njets () const
00132 //
00133 // Purpose: Return the length of the jet list.
00134 //
00135 // Returns:
00136 //   The length of the jet list.
00137 //
00138 {
00139   return _jets.size ();
00140 }
00141 
00142 
00143 Lepjets_Event_Lep& Lepjets_Event::lep (std::vector<Lepjets_Event_Lep>::size_type  i)
00144 //
00145 // Purpose: Return the Ith lepton.
00146 //
00147 // Inputs:
00148 //   i -           The lepton index (counting from 0).
00149 //
00150 // Returns:
00151 //   The Ith lepton.
00152 //
00153 {
00154   assert (i < _leps.size());
00155   return _leps[i];
00156 }
00157 
00158 
00159 Lepjets_Event_Jet& Lepjets_Event::jet (std::vector<Lepjets_Event_Jet>::size_type i)
00160 //
00161 // Purpose: Return the Ith jet.
00162 //
00163 // Inputs:
00164 //   i -           The jet index (counting from 0).
00165 //
00166 // Returns:
00167 //   The Ith jet.
00168 //
00169 {
00170   assert (i < _jets.size());
00171   return _jets[i];
00172 }
00173 
00174 
00175 const Lepjets_Event_Lep& Lepjets_Event::lep (std::vector<Lepjets_Event_Lep>::size_type i) const
00176 //
00177 // Purpose: Return the Ith lepton.
00178 //
00179 // Inputs:
00180 //   i -           The lepton index (counting from 0).
00181 //
00182 // Returns:
00183 //   The Ith lepton.
00184 //
00185 {
00186   assert (i < _leps.size());
00187   return _leps[i];
00188 }
00189 
00190 
00191 const Lepjets_Event_Jet& Lepjets_Event::jet (std::vector<Lepjets_Event_Jet>::size_type i) const
00192 //
00193 // Purpose: Return the Ith jet.
00194 //
00195 // Inputs:
00196 //   i -           The jet index (counting from 0).
00197 //
00198 // Returns:
00199 //   The Ith jet.
00200 //
00201 {
00202   assert (i < _jets.size());
00203   return _jets[i];
00204 }
00205 
00206 
00207 Fourvec& Lepjets_Event::met ()
00208 //
00209 // Purpose: Return the missing Et.
00210 //
00211 // Returns:
00212 //   The missing Et.
00213 //
00214 {
00215   return _met;
00216 }
00217 
00218 
00219 const Fourvec& Lepjets_Event::met () const
00220 //
00221 // Purpose: Return the missing Et.
00222 //
00223 // Returns:
00224 //   The missing Et.
00225 //
00226 {
00227   return _met;
00228 }
00229 
00230 
00231 Resolution& Lepjets_Event::kt_res ()
00232 //
00233 // Purpose: Return the kt resolution.
00234 //
00235 // Returns:
00236 //   The kt resolution.
00237 //
00238 {
00239   return _kt_res;
00240 }
00241 
00242 
00243 const Resolution& Lepjets_Event::kt_res () const
00244 //
00245 // Purpose: Return the kt resolution.
00246 //
00247 // Returns:
00248 //   The kt resolution.
00249 //
00250 {
00251   return _kt_res;
00252 }
00253 
00254 
00255 double Lepjets_Event::zvertex () const
00256 //
00257 // Purpose: Return the z-vertex.
00258 //
00259 // Returns:
00260 //   The z-vertex.
00261 //
00262 {
00263   return _zvertex;
00264 }
00265 
00266 
00267 double& Lepjets_Event::zvertex ()
00268 //
00269 // Purpose: Return the z-vertex.
00270 //
00271 // Returns:
00272 //   The z-vertex.
00273 //
00274 {
00275   return _zvertex;
00276 }
00277 
00278 
00279 bool Lepjets_Event::isMC () const
00280 //
00281 // Purpose: Return the isMC flag.
00282 //
00283 // Returns:
00284 //   The isMC flag.
00285 //
00286 {
00287   return _isMC;
00288 }
00289 
00290 
00291 void Lepjets_Event::setMC (bool isMC)
00292 //
00293 // Purpose: set isMC flag.
00294 //
00295 // Returns:
00296 //   nothing
00297 //
00298 {
00299   _isMC = isMC;
00300 }
00301 
00302 double Lepjets_Event::dlb () const
00303 //
00304 // Purpose: Return the LB discriminant.
00305 //
00306 // Returns:
00307 //   The LB discriminant.
00308 //
00309 {
00310   return _dlb;
00311 }
00312 
00313 
00314 double& Lepjets_Event::dlb ()
00315 //
00316 // Purpose: Return the LB discriminant.
00317 //
00318 // Returns:
00319 //   The LB discriminant.
00320 //
00321 {
00322   return _dlb;
00323 }
00324 
00325 
00326 double Lepjets_Event::dnn () const
00327 //
00328 // Purpose: Return the NN discriminant.
00329 //
00330 // Returns:
00331 //   The NN discriminant.
00332 //
00333 {
00334   return _dnn;
00335 }
00336 
00337 
00338 double& Lepjets_Event::dnn ()
00339 //
00340 // Purpose: Return the NN discriminant.
00341 //
00342 // Returns:
00343 //   The NN discriminant.
00344 //
00345 {
00346   return _dnn;
00347 }
00348 
00349 
00350 Fourvec Lepjets_Event::sum (int type) const
00351 //
00352 // Purpose: Sum all objects with type code TYPE.
00353 //
00354 // Inputs:
00355 //   type -        The type code to match.
00356 //
00357 // Returns:
00358 //   The sum of all objects with type code TYPE.
00359 //
00360 {
00361   Fourvec out;
00362   for (std::vector<Lepjets_Event_Lep>::size_type  i=0; i < _leps.size(); i++)
00363     if (_leps[i].type() == type)
00364       out += _leps[i].p();
00365   for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++)
00366     if (_jets[i].type() == type)
00367       out += _jets[i].p();
00368   return out;
00369 }
00370 
00371 
00372 Fourvec Lepjets_Event::kt () const
00373 //
00374 // Purpose: Calculate kt --- sum of all objects plus missing Et.
00375 //
00376 // Returns:
00377 //   The event kt.
00378 {
00379   Fourvec v = _met;
00380   for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++)
00381     v += _leps[i].p();
00382   for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++)
00383     v += _jets[i].p();
00384   return v;
00385 }
00386 
00387 
00388 void Lepjets_Event::add_lep (const Lepjets_Event_Lep& lep)
00389 //
00390 // Purpose: Add a lepton to the event.
00391 //
00392 // Inputs:
00393 //   lep -         The lepton to add.
00394 //
00395 {
00396   _leps.push_back (lep);
00397 }
00398 
00399 
00400 void Lepjets_Event::add_jet (const Lepjets_Event_Jet& jet)
00401 //
00402 // Purpose: Add a jet to the event.
00403 //
00404 // Inputs:
00405 //   jet -         The jet to add.
00406 //
00407 {
00408   _jets.push_back (jet);
00409 }
00410 
00411 
00412 void Lepjets_Event::smear (CLHEP::HepRandomEngine& engine, bool smear_dir /*= false*/)
00413 //
00414 // Purpose: Smear the objects in the event according to their resolutions.
00415 //
00416 // Inputs:
00417 //   engine -      The underlying RNG.
00418 //   smear_dir -   If false, smear the momentum only.
00419 //
00420 {
00421   Fourvec before, after;
00422   for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
00423     before += _leps[i].p();
00424     _leps[i].smear (engine, smear_dir);
00425     after += _leps[i].p();
00426   }
00427   for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
00428     before += _jets[i].p();
00429     _jets[i].smear (engine, smear_dir);
00430     after += _jets[i].p();
00431   }
00432 
00433   Fourvec kt = _met + before;
00434   kt(Fourvec::X) = _kt_res.pick (kt(Fourvec::X), kt(Fourvec::X), engine);
00435   kt(Fourvec::Y) = _kt_res.pick (kt(Fourvec::Y), kt(Fourvec::Y), engine);
00436   _met = kt - after;
00437 }
00438 
00439 
00440 void Lepjets_Event::sort ()
00441 //
00442 // Purpose: Sort the objects in the event in order of descending pt.
00443 //
00444 {
00445   std::stable_sort (_leps.begin(), _leps.end(), not2 (less<Lepjets_Event_Lep> ()));
00446   std::stable_sort (_jets.begin(), _jets.end(), not2 (less<Lepjets_Event_Lep> ()));
00447 }
00448 
00449 
00450 std::vector<int> Lepjets_Event::jet_types() const
00451 //
00452 // Purpose: Return the jet types of the event
00453 //
00454 {
00455   std::vector<int> ret;
00456   for (std::vector<Lepjets_Event_Jet>::size_type ijet =  0 ;
00457        ijet != _jets.size() ;
00458        ijet++) {
00459     ret.push_back(jet(ijet).type());
00460   }
00461   return ret;
00462 }
00463 
00464 
00465 bool Lepjets_Event::set_jet_types(const std::vector<int>& _jet_types)
00466 //
00467 // Purpose: Set the jet types of the event
00468 // Return false if it fails, trus if it succeeds
00469 //
00470 {
00471   if (_jets.size() != _jet_types.size()) {
00472     return false;
00473   }
00474   bool saw_hadw1 = false;
00475   for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < njets(); i++) {
00476     int t = _jet_types[i];
00477     if (t == hadw1_label) {
00478       if (saw_hadw1)
00479         t = hadw2_label;
00480       saw_hadw1 = true;
00481     }
00482     jet (i).type() = t;
00483   }
00484   return true;
00485 }
00486 
00487 
00488 namespace {
00489 
00490 
00491 struct Lepjets_Event_Cutter
00492 //
00493 // Purpose: Helper for cutting on objects.
00494 //
00495 {
00496   Lepjets_Event_Cutter (double pt_cut, double eta_cut)
00497     : _pt_cut (pt_cut), _eta_cut (eta_cut)
00498   {}
00499   bool operator() (const Lepjets_Event_Lep& l) const;
00500   double _pt_cut;
00501   double _eta_cut;
00502 };
00503 
00504 
00505 bool Lepjets_Event_Cutter::operator () (const Lepjets_Event_Lep& l) const
00506 //
00507 // Purpose: Object cut evaluator.
00508 //
00509 {
00510   return ! (l.p().perp() > _pt_cut && abs (l.p().pseudoRapidity()) < _eta_cut);
00511 }
00512 
00513 
00514 } // unnamed namespace
00515 
00516 
00517 int Lepjets_Event::cut_leps (double pt_cut, double eta_cut)
00518 //
00519 // Purpose: Remove all leptons failing the pt and eta cuts.
00520 //
00521 // Inputs:
00522 //   pt_cut -      Pt cut.  Remove objects with pt less than this.
00523 //   eta_cut -     Eta cut.  Remove objects with abs(eta) larger than this.
00524 //
00525 // Returns:
00526 //   The number of leptons remaining after the cuts.
00527 //
00528 {
00529   _leps.erase (remove_if (_leps.begin(), _leps.end(),
00530                           Lepjets_Event_Cutter (pt_cut, eta_cut)),
00531                _leps.end ());
00532   return _leps.size ();
00533 }
00534 
00535 
00536 int Lepjets_Event::cut_jets (double pt_cut, double eta_cut)
00537 //
00538 // Purpose: Remove all jets failing the pt and eta cuts.
00539 //
00540 // Inputs:
00541 //   pt_cut -      Pt cut.  Remove objects with pt less than this.
00542 //   eta_cut -     Eta cut.  Remove objects with abs(eta) larger than this.
00543 //
00544 // Returns:
00545 //   The number of jets remaining after the cuts.
00546 //
00547 {
00548   _jets.erase (remove_if (_jets.begin(), _jets.end(),
00549                           Lepjets_Event_Cutter (pt_cut, eta_cut)),
00550                _jets.end ());
00551   return _jets.size ();
00552 }
00553 
00554 
00555 void Lepjets_Event::trimjets (std::vector<Lepjets_Event_Jet>::size_type n)
00556 //
00557 // Purpose: Remove all but the first N jets.
00558 //
00559 // Inputs:
00560 //   n -           The number of jets to keep.
00561 //
00562 {
00563   if (n >= _jets.size())
00564     return;
00565   _jets.erase (_jets.begin() + n, _jets.end());
00566 }
00567 
00568 
00569 std::ostream& Lepjets_Event::dump (std::ostream& s, bool full /*=false*/) const
00570 //
00571 // Purpose: Dump out this object.
00572 //
00573 // Inputs:
00574 //   s -           The stream to which to write.
00575 //   full -        If true, dump all information for this object.
00576 //
00577 // Returns:
00578 //   The stream S.
00579 //
00580 {
00581   s << "Run: " << _runnum << "  Event: " << _evnum << "\n";
00582   s << "Leptons:\n";
00583   for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
00584     s << "  ";
00585     _leps[i].dump (s, full);
00586     s << "\n";
00587   }
00588   s << "Jets:\n";
00589   for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
00590     s << "  ";
00591     _jets[i].dump (s, full);
00592     s << "\n";
00593   }
00594   s << "Missing Et: " << _met << "\n";
00595   if (_zvertex != 0)
00596     s << "z-vertex: " << _zvertex << "\n";
00597   return s;
00598 }
00599 
00600 
00601 std::string Lepjets_Event::jet_permutation() const
00602 //
00603 // Purpose: Return a string representation of the jet permutation
00604 //          g - isr/gluon
00605 //          b - leptonic b
00606 //          B - hadronic b
00607 //          w - hadronic W
00608 //          h - Higgs to b-bbar
00609 //          ? - Unknown
00610 {
00611     std::string permutation;
00612     for (size_t jet = 0 ; jet != _jets.size() ; ++jet) {
00613         permutation = permutation + hitfit::jetTypeString(_jets[jet].type());
00614     }
00615     return permutation;
00616 }
00617 
00626 std::ostream& operator<< (std::ostream& s, const Lepjets_Event& ev)
00627 //
00628 // Purpose: Dump out this object.
00629 //
00630 // Inputs:
00631 //   s -           The stream to which to write.
00632 //   ev -          The object to dump.
00633 //
00634 // Returns:
00635 //   The stream S.
00636 //
00637 {
00638   return ev.dump (s);
00639 }
00640 
00641 
00642 } // namespace hitfit
00643 
00644 
00645 
00646 
00647