rivet is hosted by Hepforge, IPPP Durham
ATLAS_2012_I1199269.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/IdentifiedFinalState.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 
00006 namespace Rivet {
00007 
00008 
00009   /// @brief Measurement of isolated diphoton + X differential cross-sections
00010   ///
00011   /// Inclusive isolated gamma gamma cross-sections, differential in M(gg), pT(gg),
00012   /// dphi(gg), cos(theta*)_CS
00013   ///
00014   /// @author Giovanni Marchiori
00015   ///
00016   class ATLAS_2012_I1199269 : public Analysis {
00017   public:
00018 
00019     /// Constructor
00020     ATLAS_2012_I1199269()
00021       : Analysis("ATLAS_2012_I1199269")
00022     {    }
00023 
00024 
00025     /// Book histograms and initialise projections before the run
00026     void init() {
00027 
00028       FinalState fs;
00029       declare(fs, "FS");
00030 
00031       FastJets fj(fs, FastJets::KT, 0.5);
00032       fj.useJetArea(new fastjet::AreaDefinition(fastjet::VoronoiAreaSpec()));
00033       declare(fj, "KtJetsD05");
00034 
00035       IdentifiedFinalState photonfs(Cuts::abseta < 2.37 && Cuts::pT > 22*GeV);
00036       photonfs.acceptId(PID::PHOTON);
00037       declare(photonfs, "Photon");
00038 
00039       _h_M            = bookHisto1D(1, 1, 1);
00040       _h_pT           = bookHisto1D(2, 1, 1);
00041       _h_dPhi         = bookHisto1D(3, 1, 1);
00042       _h_cosThetaStar = bookHisto1D(4, 1, 1);
00043     }
00044 
00045 
00046     /// Perform the per-event analysis
00047     void analyze(const Event& event) {
00048 
00049       // Require at least 2 photons in final state
00050       const Particles photons = apply<IdentifiedFinalState>(event, "Photon").particlesByPt();
00051       if (photons.size() < 2) vetoEvent;
00052 
00053       // Get jets, and corresponding jet areas
00054       vector<vector<double> > ptDensities(_eta_bins_areaoffset.size()-1);
00055       const auto clust_seq_area = apply<FastJets>(event, "KtJetsD05").clusterSeqArea();
00056       for (const Jet& jet : apply<FastJets>(event, "KtJetsD05").jets()) {
00057         const double area = clust_seq_area->area(jet); // implicit .pseudojet()
00058         if (area < 1e-3) continue;
00059         const int ieta = binIndex(jet.abseta(), _eta_bins_areaoffset);
00060         if (ieta != -1) ptDensities[ieta].push_back(jet.pT()/area);
00061       }
00062 
00063       // Compute median jet properties over the jets in the event
00064       vector<double> vptDensity; //, vsigma, vNjets;
00065       for (size_t b = 0; b < _eta_bins_areaoffset.size()-1; ++b) {
00066         vptDensity += ptDensities[b].empty() ? 0 : median(ptDensities[b]);
00067       }
00068 
00069 
00070       // Loop over photons and fill vector of isolated ones
00071       Particles isolated_photons;
00072       for (const Particle& photon : photons) {
00073         /// Remove photons in ECAL crack region
00074         if (inRange(photon.abseta(), 1.37, 1.52)) continue;
00075         // Compute isolation via particles within an R=0.4 cone of the photon
00076         const Particles& fs = apply<FinalState>(event, "FS").particles();
00077         FourMomentum mom_in_EtCone;
00078         for (const Particle& p : fs) {
00079           // Reject if not in cone
00080           if (deltaR(photon, p) > 0.4) continue;
00081           // Reject if in the 5x7 cell central core
00082           if (fabs(deltaEta(photon, p)) < 0.025 * 5 * 0.5 &&
00083               fabs(deltaPhi(photon, p)) < PI/128. * 7 * 0.5) continue;
00084           // Sum momentum
00085           mom_in_EtCone += p.momentum();
00086         }
00087         // Now figure out the correction (area*density)
00088         const double ETCONE_AREA = PI*sqr(0.4) - (7*.025)*(5*PI/128.); // cone area - central core rectangle
00089         const double correction = vptDensity[binIndex(photon.abseta(), _eta_bins_areaoffset)] * ETCONE_AREA;
00090 
00091         // Discard the photon if there is more than 4 GeV of cone activity
00092         // NOTE: Shouldn't need to subtract photon itself (it's in the central core)
00093         // NOTE: using expected cut at hadron/particle level, not at reco level
00094         if (mom_in_EtCone.Et() - correction > 4*GeV) continue;
00095         // Add isolated photon to list
00096         isolated_photons.push_back(photon);
00097       }
00098 
00099       // Require at least two isolated photons and select leading pT pair
00100       if (isolated_photons.size() < 2) vetoEvent;
00101       sortByPt(isolated_photons);
00102       const FourMomentum& y1 = isolated_photons[0].momentum();
00103       const FourMomentum& y2 = isolated_photons[1].momentum();
00104 
00105       // Leading photon should have pT > 25 GeV
00106       if (y1.pT() < 25*GeV) vetoEvent;
00107 
00108       // Require the two photons to be separated by dR > 0.4
00109       if (deltaR(y1, y2) < 0.4) vetoEvent;
00110 
00111       // Compute diphoton vector and fill histos
00112       const double weight = event.weight();
00113       FourMomentum yy = y1 + y2;
00114       const double costhetayy = 2 * y1.pT() * y2.pT() * sinh(y1.eta() - y2.eta()) / yy.mass() / add_quad(yy.mass(), yy.pT());
00115       _h_M->fill(yy.mass()/GeV, weight);
00116       _h_pT->fill(yy.pT()/GeV, weight);
00117       _h_dPhi->fill(mapAngle0ToPi(y1.phi() - y2.phi()), weight);
00118       _h_cosThetaStar->fill(costhetayy, weight);
00119     }
00120 
00121 
00122     /// Normalise histograms etc., after the run
00123     void finalize() {
00124       scale(_h_M, crossSection()/sumOfWeights());
00125       scale(_h_pT, crossSection()/sumOfWeights());
00126       scale(_h_dPhi, crossSection()/sumOfWeights());
00127       scale(_h_cosThetaStar, crossSection()/sumOfWeights());
00128     }
00129 
00130 
00131   private:
00132 
00133     Histo1DPtr _h_M, _h_pT, _h_dPhi, _h_cosThetaStar;
00134 
00135     const vector<double> _eta_bins_areaoffset = {0.0, 1.5, 3.0};
00136 
00137   };
00138 
00139 
00140   DECLARE_RIVET_PLUGIN(ATLAS_2012_I1199269);
00141 
00142 }