ATLAS_2015_I1364361.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/FinalState.hh" 00004 #include "Rivet/Projections/FastJets.hh" 00005 00006 namespace Rivet { 00007 00008 00009 /// @brief Higgs differential cross-section combination between the ATLAS measurements in the yy and 4l channels 00010 /// 00011 /// Computes Higgs transverse momentum, rapidity, jet multiplicity and leading jet pT. 00012 /// 00013 /// @author Michaela Queitsch-Maitland <michaela.queitsch-maitland@cern.ch> 00014 /// @author Dag Gillberg <dag.gillberg@cern.ch> 00015 /// @author Florian Bernlochner <florian.bernlochner@cern.ch> 00016 /// @author Sarah Heim <sarah.heim@cern.ch> 00017 class ATLAS_2015_I1364361 : public Analysis { 00018 public: 00019 00020 /// Constructor 00021 DEFAULT_RIVET_ANALYSIS_CTOR(ATLAS_2015_I1364361); 00022 00023 00024 /// Book histograms and initialise projections before the run 00025 void init() { 00026 00027 // All final state particles 00028 const FinalState fs; 00029 declare(fs, "FS"); 00030 00031 // Histograms with data bins 00032 _h_pTH_incl = bookHisto1D(1,1,1); 00033 _h_yH_incl = bookHisto1D(2,1,1); 00034 _h_Njets_incl = bookHisto1D(3,1,1); 00035 _h_pTj1_incl = bookHisto1D(4,1,1); 00036 } 00037 00038 00039 /// Perform the per-event analysis 00040 void analyze(const Event& event) { 00041 00042 // Get the final state particles ordered by pT 00043 const Particles& fs = apply<FinalState>(event, "FS").particlesByPt(); 00044 00045 // Find a stable Higgs (mandatory) 00046 const auto higgsiter = std::find_if(fs.begin(), fs.end(), [](const Particle& p){ return p.pid() == PID::HIGGSBOSON; }); 00047 if (higgsiter == fs.end()) vetoEvent; 00048 const Particle& higgs = *higgsiter; 00049 00050 // bool stable_higgs = false; 00051 // const Particle* higgs = 0; 00052 // for (const Particle& p : fs) { 00053 // if (p.pid() == PID::HIGGSBOSON) { 00054 // stable_higgs = true; 00055 // higgs = &p; 00056 // break; 00057 // } 00058 // } 00059 // if (!stable_higgs) { 00060 // MSG_WARNING("FATAL: No stable Higgs found in event record.\n"); 00061 // vetoEvent; 00062 // } 00063 00064 00065 // Loop over final state particles and fill various particle vectors 00066 Particles leptons, photons, jet_ptcls; 00067 foreach ( const Particle& ptcl, fs ) { 00068 // Do not include the Higgs in jet finding! 00069 if ( ptcl.pid() == PID::HIGGSBOSON ) continue; 00070 // Neutrinos not from hadronisation 00071 if ( ptcl.isNeutrino() && !ptcl.fromHadron() ) continue; 00072 // Electrons and muons not from hadronisation 00073 if ( ( ptcl.abspid() == PID::ELECTRON || ptcl.abspid() == PID::MUON ) && !ptcl.fromHadron() ) { 00074 leptons.push_back(ptcl); 00075 continue; 00076 } 00077 // Photons not from hadronisation 00078 if ( ptcl.abspid() == PID::PHOTON && !ptcl.fromHadron() ) { 00079 photons.push_back(ptcl); 00080 continue; 00081 } 00082 // Add particle to jet inputs 00083 jet_ptcls.push_back(ptcl); 00084 } 00085 00086 // Match FS photons to leptons within cone R=0.1 00087 // If they are not 'dressing' photons, add to jet particle vector 00088 foreach ( const Particle& ph, photons ) { 00089 bool fsr_photon = false; 00090 foreach ( const Particle& lep, leptons ) { 00091 if ( deltaR(ph.momentum(),lep.momentum()) < 0.1 ){ 00092 fsr_photon=true; 00093 continue; 00094 } 00095 } 00096 if ( !fsr_photon ) jet_ptcls.push_back(ph); 00097 } 00098 00099 // Let's build the jets! By hand... 00100 const PseudoJets pjs_in = mkPseudoJets(jet_ptcls); 00101 const fastjet::JetDefinition jdef(fastjet::antikt_algorithm, 0.4); 00102 const Jets alljets = mkJets(fastjet::ClusterSequence(pjs_in, jdef).inclusive_jets()); 00103 const Jets jets = sortByPt(filterBy(alljets, Cuts::pT > 30*GeV && Cuts::absrap < 4.4)); 00104 // FastJets jet_pro(FastJets::ANTIKT, 0.4); 00105 // jet_pro.calc(jet_ptcls); 00106 // Jets jets = jet_pro.jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 4.4); 00107 00108 const double weight = event.weight(); 00109 _h_pTH_incl->fill(higgs.pT(), weight); 00110 _h_yH_incl->fill(higgs.absrap(), weight); 00111 _h_Njets_incl->fill(jets.size() > 3 ? 3 : jets.size(), weight); 00112 _h_pTj1_incl->fill(jets.empty() ? 0 : jets[0].pT(), weight); 00113 } 00114 00115 00116 /// Normalise histograms etc., after the run 00117 void finalize() { 00118 const double xs = crossSectionPerEvent(); 00119 scale(_h_pTH_incl, xs); 00120 scale(_h_yH_incl, xs); 00121 scale(_h_Njets_incl, xs); 00122 scale(_h_pTj1_incl, xs); 00123 } 00124 00125 00126 private: 00127 00128 Histo1DPtr _h_pTH_incl, _h_yH_incl, _h_Njets_incl, _h_pTj1_incl; 00129 00130 }; 00131 00132 00133 // The hook for the plugin system 00134 DECLARE_RIVET_PLUGIN(ATLAS_2015_I1364361); 00135 00136 00137 } Generated on Tue Dec 13 2016 16:32:35 for The Rivet MC analysis system by ![]() |