rivet is hosted by Hepforge, IPPP Durham
ATLAS_2012_I1125575.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/ChargedFinalState.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 #include "Rivet/Tools/BinnedHistogram.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   class ATLAS_2012_I1125575 : public Analysis {
00011     public:
00012 
00013       /// @name Constructors etc.
00014       //@{
00015 
00016       /// Constructor
00017       ATLAS_2012_I1125575()
00018         : Analysis("ATLAS_2012_I1125575")
00019       {    }
00020 
00021       //@}
00022 
00023 
00024     public:
00025 
00026       /// @name Analysis methods
00027       //@{
00028 
00029       /// Book histograms and initialise projections before the run
00030       void init() {
00031 
00032         /// @todo Initialise and register projections here
00033         const ChargedFinalState jet_input(-2.5, 2.5, 0.5*GeV);
00034         addProjection(jet_input, "JET_INPUT");
00035 
00036         const ChargedFinalState track_input(-1.5, 1.5, 0.5*GeV);
00037         addProjection(track_input, "TRACK_INPUT");
00038 
00039         const FastJets jets02(jet_input, FastJets::ANTIKT, 0.2);
00040         addProjection(jets02, "JETS_02");
00041 
00042         const FastJets jets04(jet_input, FastJets::ANTIKT, 0.4);
00043         addProjection(jets04, "JETS_04");
00044 
00045         const FastJets jets06(jet_input, FastJets::ANTIKT, 0.6);
00046         addProjection(jets06, "JETS_06");
00047 
00048         const FastJets jets08(jet_input, FastJets::ANTIKT, 0.8);
00049         addProjection(jets08, "JETS_08");
00050 
00051         const FastJets jets10(jet_input, FastJets::ANTIKT, 1.0);
00052         addProjection(jets10, "JETS_10");
00053 
00054         // Mean number of tracks
00055         InitializeProfiles(_h_meanNch, 1);
00056 
00057         // Mean of the average track pT in each region
00058         InitializeProfiles(_h_meanPtAvg, 2);
00059 
00060         // Mean of the scalar sum of track pT in each region
00061         InitializeProfiles(_h_meanPtSum, 3);
00062 
00063         // Distribution of Nch, in bins of leading track-jet pT
00064         InitializeHistograms(_h_Nch, 4);
00065 
00066         // Distribution of average track-jet pT, in bins of leading track-jet pT
00067         InitializeHistograms(_h_PtAvg, 5);
00068 
00069         // Distribution of sum of track-jet pT, in bins of leading track-jet pT
00070         InitializeHistograms(_h_PtSum, 6);
00071 
00072         for (int i=0; i<5; ++i) {
00073           _nEvents[i] = 0.0;
00074         }
00075 
00076       }
00077 
00078       void InitializeProfiles(Profile1DPtr plots[5][2], int distribution) {
00079         for (int i=0; i<5; ++i) {
00080           for (int j=0; j<2; ++j) {
00081             plots[i][j] = bookProfile1D(distribution, i+1, j+1);
00082           }
00083         }
00084       }
00085 
00086       void InitializeHistograms(BinnedHistogram<double> plots[5][2], int distribution) {
00087         Scatter2D refscatter = refData(1, 1, 1);
00088         for (int i=0; i<5; ++i) {
00089           for (int y=0; y<2; ++y) {
00090             for (size_t j=0; j<refscatter.numPoints(); ++j) {
00091               int histogram_number = ((j+1)*2)-((y+1)%2);
00092               double low_edge = refscatter.point(j).xMin();
00093               double high_edge = refscatter.point(j).xMax();
00094               plots[i][y].addHistogram(low_edge, high_edge, bookHisto1D(distribution, i+1, histogram_number));
00095             }
00096           }
00097         }
00098       }
00099 
00100 
00101       /// Perform the per-event analysis
00102       void analyze(const Event& event) {
00103         const double weight = event.weight();
00104 
00105         std::vector<Jets*> all_jets;
00106         Jets jets_02 = applyProjection<FastJets>(event, "JETS_02").jetsByPt(4.0*GeV, MAXDOUBLE, -1.5, 1.5, PSEUDORAPIDITY);
00107         all_jets.push_back(&jets_02);
00108         Jets jets_04 = applyProjection<FastJets>(event, "JETS_04").jetsByPt(4.0*GeV, MAXDOUBLE, -1.5, 1.5, PSEUDORAPIDITY);
00109         all_jets.push_back(&jets_04);
00110         Jets jets_06 = applyProjection<FastJets>(event, "JETS_06").jetsByPt(4.0*GeV, MAXDOUBLE, -1.5, 1.5, PSEUDORAPIDITY);
00111         all_jets.push_back(&jets_06);
00112         Jets jets_08 = applyProjection<FastJets>(event, "JETS_08").jetsByPt(4.0*GeV, MAXDOUBLE, -1.5, 1.5, PSEUDORAPIDITY);
00113         all_jets.push_back(&jets_08);
00114         Jets jets_10 = applyProjection<FastJets>(event, "JETS_10").jetsByPt(4.0*GeV, MAXDOUBLE, -1.5, 1.5, PSEUDORAPIDITY);
00115         all_jets.push_back(&jets_10);
00116 
00117         // Count the number of tracks in the away and transverse regions, for each set of jets
00118         double n_ch[5][2] = { {0,0}, {0,0}, {0,0}, {0,0}, {0,0} };
00119         // Also add up the sum pT
00120         double sumpt[5][2] = { {0,0}, {0,0}, {0,0}, {0,0}, {0,0} };
00121         // ptmean = sumpt / n_ch
00122         double ptavg[5][2] = { {0,0}, {0,0}, {0,0}, {0,0}, {0,0} };
00123         // lead jet pT defines which bin we want to fill
00124         double lead_jet_pts[5] = {0.0};
00125 
00126         // Loop over each of the jet radii:
00127         for (int i=0; i<5; ++i) {
00128           if (all_jets[i]->size() < 1) continue;
00129 
00130           // Find the lead jet pT
00131           lead_jet_pts[i] = all_jets[i]->at(0).pT();
00132 
00133           // Loop over each of the charged particles
00134           const Particles& tracks = applyProjection<ChargedFinalState>(event, "TRACK_INPUT").particlesByPt();
00135           foreach(const Particle& t, tracks) {
00136 
00137             // Get the delta-phi between the track and the leading jet
00138             double dphi = deltaPhi(all_jets[i]->at(0), t);
00139 
00140             // Find out which region this puts it in.
00141             // 0 = away region, 1 = transverse region, 2 = toward region
00142             int region = region_index(dphi);
00143 
00144             // If the track is in the toward region, ignore it.
00145             if (region == 2) continue;
00146 
00147             // Otherwise, increment the relevant counters
00148             ++n_ch[i][region];
00149             sumpt[i][region] += t.pT();
00150 
00151           }
00152           // Calculate the pT_avg for the away and transverse regions.
00153           // (And make sure we don't try to divide by zero.)
00154           ptavg[i][0] = (n_ch[i][0] == 0 ? 0.0 : sumpt[i][0] / n_ch[i][0]);
00155           ptavg[i][1] = (n_ch[i][1] == 0 ? 0.0 : sumpt[i][1] / n_ch[i][1]);
00156 
00157           _nEvents[i] += weight;
00158         }
00159 
00160         FillProfiles(_h_meanNch,    n_ch, lead_jet_pts, weight, 1.0 / (2*PI));
00161         FillProfiles(_h_meanPtAvg, ptavg, lead_jet_pts, weight, 1.0);
00162         FillProfiles(_h_meanPtSum, sumpt, lead_jet_pts, weight, 1.0 / (2*PI));
00163 
00164         FillHistograms(_h_Nch,    n_ch, lead_jet_pts, weight);
00165         FillHistograms(_h_PtAvg, ptavg, lead_jet_pts, weight);
00166         FillHistograms(_h_PtSum, sumpt, lead_jet_pts, weight);
00167       }
00168 
00169       void FillProfiles(Profile1DPtr plots[5][2], double var[5][2], double lead_pt[5], double weight, double scale) {
00170         for (int i=0; i<5; ++i) {
00171           double pt = lead_pt[i];
00172           for (int j=0; j<2; ++j) {
00173             double v = var[i][j];
00174             plots[i][j]->fill(pt, v*scale, weight);
00175           }
00176         }
00177       }
00178 
00179       void FillHistograms(BinnedHistogram<double> plots[5][2], double var[5][2], double lead_pt[5], double weight) {
00180         for (int i=0; i<5; ++i) {
00181           double pt = lead_pt[i];
00182           for (int j=0; j<2; ++j) {
00183             double v = var[i][j];
00184             plots[i][j].fill(pt, v, weight);
00185           }
00186         }
00187       }
00188 
00189       inline int region_index(double dphi) {
00190         assert(inRange(dphi, 0.0, PI, CLOSED, CLOSED));
00191         if (dphi < PI/3.0) return 2;
00192         if (dphi < 2*PI/3.0) return 1;
00193         return 0;
00194       }
00195 
00196 
00197       /// Normalise histograms etc., after the run
00198       void finalize() {
00199         FinalizeHistograms(_h_Nch);
00200         FinalizeHistograms(_h_PtAvg);
00201         FinalizeHistograms(_h_PtSum);
00202       }
00203 
00204       void FinalizeHistograms(BinnedHistogram<double> plots[5][2]) {
00205         for (int i=0; i<5; ++i) {
00206           for (int j=0; j<2; ++j) {
00207             vector<Histo1DPtr> histos = plots[i][j].getHistograms();
00208             foreach(Histo1DPtr h, histos) {
00209               scale(h, 1.0/_nEvents[i]);
00210             }
00211           }
00212         }
00213       }
00214 
00215       //@}
00216 
00217 
00218     private:
00219 
00220       // Data members like post-cuts event weight counters go here
00221       double _nEvents[5];
00222 
00223     private:
00224 
00225       Profile1DPtr _h_meanNch[5][2];
00226       Profile1DPtr _h_meanPtAvg[5][2];
00227       Profile1DPtr _h_meanPtSum[5][2];
00228 
00229       BinnedHistogram<double> _h_Nch[5][2];
00230       BinnedHistogram<double> _h_PtAvg[5][2];
00231       BinnedHistogram<double> _h_PtSum[5][2];
00232 
00233   };
00234 
00235 
00236 
00237   // The hook for the plugin system
00238   DECLARE_RIVET_PLUGIN(ATLAS_2012_I1125575);
00239 
00240 }