rivet is hosted by Hepforge, IPPP Durham
CMS_2015_I1310737.cc
Go to the documentation of this file.
00001 #include "Rivet/Analysis.hh"
00002 #include "Rivet/Projections/FastJets.hh"
00003 #include "Rivet/Projections/ChargedFinalState.hh"
00004 #include "Rivet/Projections/InvMassFinalState.hh"
00005 #include "Rivet/Projections/VisibleFinalState.hh"
00006 #include "Rivet/Projections/VetoedFinalState.hh"
00007 #include "Rivet/Projections/ZFinder.hh"
00008 
00009 namespace Rivet {
00010 
00011 
00012   class CMS_2015_I1310737 : public Analysis {
00013   public:
00014 
00015     /// Constructor
00016     CMS_2015_I1310737()
00017       : Analysis("CMS_2015_I1310737")
00018     {  }
00019 
00020 
00021     /// Book histograms and initialise projections before the run
00022     void init() {
00023 
00024       FinalState fs; ///< @todo No cuts?
00025       VisibleFinalState visfs(fs);
00026 
00027       ZFinder zeeFinder(fs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::ELECTRON, 71.0*GeV, 111.0*GeV);
00028       addProjection(zeeFinder, "ZeeFinder");
00029 
00030       ZFinder zmumuFinder(fs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::MUON, 71.0*GeV, 111.0*GeV);
00031       addProjection(zmumuFinder, "ZmumuFinder");
00032 
00033       VetoedFinalState jetConstits(visfs);
00034       jetConstits.addVetoOnThisFinalState(zeeFinder);
00035       jetConstits.addVetoOnThisFinalState(zmumuFinder);
00036 
00037       FastJets akt05Jets(jetConstits, FastJets::ANTIKT, 0.5);
00038       addProjection(akt05Jets, "AntiKt05Jets");
00039 
00040 
00041       _h_excmult_jets_tot = bookHisto1D(1, 1, 1);
00042       _h_incmult_jets_tot = bookHisto1D(2, 1, 1);
00043       _h_leading_jet_pt_tot = bookHisto1D(3, 1, 1);
00044       _h_second_jet_pt_tot = bookHisto1D(4, 1, 1);
00045       _h_third_jet_pt_tot = bookHisto1D(5, 1, 1);
00046       _h_fourth_jet_pt_tot = bookHisto1D(6, 1, 1);
00047       _h_leading_jet_eta_tot = bookHisto1D(7, 1, 1);
00048       _h_second_jet_eta_tot = bookHisto1D(8, 1, 1);
00049       _h_third_jet_eta_tot = bookHisto1D(9, 1, 1);
00050       _h_fourth_jet_eta_tot = bookHisto1D(10, 1, 1);
00051       _h_ht1_tot = bookHisto1D(11, 1, 1);
00052       _h_ht2_tot = bookHisto1D(12, 1, 1);
00053       _h_ht3_tot = bookHisto1D(13, 1, 1);
00054       _h_ht4_tot = bookHisto1D(14, 1, 1);
00055     }
00056 
00057 
00058     /// Perform the per-event analysis
00059     void analyze(const Event& event) {;
00060 
00061       const ZFinder& zeeFS = applyProjection<ZFinder>(event, "ZeeFinder");
00062       const ZFinder& zmumuFS = applyProjection<ZFinder>(event, "ZmumuFinder");
00063 
00064       const Particles& zees = zeeFS.bosons();
00065       const Particles& zmumus = zmumuFS.bosons();
00066 
00067       // We did not find exactly one Z. No good.
00068       if (zees.size() + zmumus.size() != 1) {
00069         MSG_DEBUG("Did not find exactly one good Z candidate");
00070         vetoEvent;
00071       }
00072 
00073       // Find the (dressed!) leptons
00074       const Particles& dressedLeptons = zees.size() ? zeeFS.constituents() : zmumuFS.constituents();
00075 
00076       // Cluster jets
00077       // NB. Veto has already been applied on leptons and photons used for dressing
00078       const FastJets& fj = applyProjection<FastJets>(event, "AntiKt05Jets");
00079       const Jets& jets = fj.jetsByPt(Cuts::abseta < 2.4 && Cuts::pT > 30*GeV);
00080 
00081       // Perform lepton-jet overlap and HT calculation
00082       double ht = 0;
00083       Jets goodjets;
00084       foreach (const Jet& j, jets) {
00085         // Decide if this jet is "good", i.e. isolated from the leptons
00086         /// @todo Nice use-case for any() and a C++11 lambda
00087         bool overlap = false;
00088         foreach (const Particle& l, dressedLeptons) {
00089           if (Rivet::deltaR(j, l) < 0.5) {
00090             overlap = true;
00091             break;
00092           }
00093         }
00094 
00095         // Fill HT and good-jets collection
00096         if (overlap) continue;
00097         goodjets.push_back(j);
00098         ht += j.pT();
00099       }
00100 
00101       // We don't care about events with no isolated jets
00102       if (goodjets.empty()) {
00103         MSG_DEBUG("No jets in event");
00104         vetoEvent;
00105       }
00106 
00107 
00108       /////////////////
00109 
00110 
00111       // Weight to be used for histo filling
00112       const double w = 0.5 * event.weight();
00113 
00114       // Fill jet number integral histograms
00115       _h_excmult_jets_tot->fill(goodjets.size(), w);
00116       /// @todo Could be better computed by toIntegral transform on exclusive histo
00117       for (size_t iJet = 1; iJet <= goodjets.size(); iJet++ )
00118         _h_incmult_jets_tot->fill(iJet, w);
00119 
00120       // Fill leading jet histograms
00121       const Jet& j1 = goodjets[0];
00122       _h_leading_jet_pt_tot->fill(j1.pT()/GeV, w);
00123       _h_leading_jet_eta_tot->fill(j1.abseta(), w);
00124       _h_ht1_tot->fill(ht/GeV, w);
00125 
00126       // Fill 2nd jet histograms
00127       if (goodjets.size() < 2) return;
00128       const Jet& j2 = goodjets[1];
00129       _h_second_jet_pt_tot->fill(j2.pT()/GeV, w);
00130       _h_second_jet_eta_tot->fill(j2.abseta(), w);
00131       _h_ht2_tot->fill(ht/GeV, w);
00132 
00133       // Fill 3rd jet histograms
00134       if (goodjets.size() < 3) return;
00135       const Jet& j3 = goodjets[2];
00136       _h_third_jet_pt_tot->fill(j3.pT()/GeV, w);
00137       _h_third_jet_eta_tot->fill(j3.abseta(), w);
00138       _h_ht3_tot->fill(ht/GeV, w);
00139 
00140       // Fill 4th jet histograms
00141       if (goodjets.size() < 4) return;
00142       const Jet& j4 = goodjets[3];
00143       _h_fourth_jet_pt_tot->fill(j4.pT()/GeV, w);
00144       _h_fourth_jet_eta_tot->fill(j4.abseta(), w);
00145       _h_ht4_tot->fill(ht/GeV, w);
00146     }
00147 
00148 
00149     /// Normalise histograms etc., after the run
00150     void finalize() {
00151 
00152       const double norm = (sumOfWeights() != 0) ? crossSection()/sumOfWeights() : 1.0;
00153 
00154       MSG_INFO("Cross section = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(3) << crossSection() << " pb");
00155       MSG_INFO("# Events      = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(3) << numEvents() );
00156       MSG_INFO("SumW          = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(3) << sumOfWeights());
00157       MSG_INFO("Norm factor   = " << std::setfill(' ') << std::setw(14) << std::fixed << std::setprecision(6) << norm);
00158 
00159       scale(_h_excmult_jets_tot, norm );
00160       scale(_h_incmult_jets_tot, norm );
00161       scale(_h_leading_jet_pt_tot, norm );
00162       scale(_h_second_jet_pt_tot, norm );
00163       scale(_h_third_jet_pt_tot, norm );
00164       scale(_h_fourth_jet_pt_tot, norm );
00165       scale(_h_leading_jet_eta_tot, norm );
00166       scale(_h_second_jet_eta_tot, norm );
00167       scale(_h_third_jet_eta_tot, norm );
00168       scale(_h_fourth_jet_eta_tot, norm );
00169       scale(_h_ht1_tot, norm );
00170       scale(_h_ht2_tot, norm );
00171       scale(_h_ht3_tot, norm );
00172       scale(_h_ht4_tot, norm );
00173     }
00174 
00175 
00176   private:
00177 
00178     /// @name Histograms
00179 
00180     Histo1DPtr _h_excmult_jets_tot,  _h_incmult_jets_tot;
00181     Histo1DPtr _h_leading_jet_pt_tot, _h_second_jet_pt_tot, _h_third_jet_pt_tot, _h_fourth_jet_pt_tot;
00182     Histo1DPtr _h_leading_jet_eta_tot, _h_second_jet_eta_tot, _h_third_jet_eta_tot, _h_fourth_jet_eta_tot;
00183     Histo1DPtr _h_ht1_tot, _h_ht2_tot, _h_ht3_tot, _h_ht4_tot;
00184 
00185   };
00186 
00187 
00188   DECLARE_RIVET_PLUGIN(CMS_2015_I1310737);
00189 
00190 
00191 }