rivet is hosted by Hepforge, IPPP Durham
ATLAS_2011_S9126244.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetYODA.hh"
00004 #include "Rivet/Tools/Logging.hh"
00005 #include "Rivet/Projections/FinalState.hh"
00006 #include "Rivet/Projections/FastJets.hh"
00007 
00008 #include "Rivet/Tools/BinnedHistogram.hh"
00009 
00010 namespace Rivet {
00011 
00012   struct ATLAS_2011_S9126244_Plots {
00013 
00014     int selectionType; // The third value in the axis code d**-**-**
00015     std::string intermediateHistName;
00016 
00017     // Gap Fraction Vs Delta Y Plot Setup
00018     int m_gapFractionDeltaYHistIndex;
00019     std::vector<double> m_gapFractionDeltaYSlices;
00020     BinnedHistogram<double> _h_gapVsDeltaYVeto;
00021     BinnedHistogram<double> _h_gapVsDeltaYInc;
00022 
00023     // Gap Fraction Vs PtBar Plot Setup
00024     int m_gapFractionPtBarHistIndex;
00025     std::vector<double> m_gapFractionPtBarSlices;
00026     BinnedHistogram<double> _h_gapVsPtBarVeto;
00027     BinnedHistogram<double> _h_gapVsPtBarInc;
00028 
00029     // Gap Fraction Vs Q0 Plot Setup
00030     int m_gapFractionQ0HistIndex;
00031     std::vector<double> m_gapFractionQ0SlicesPtBar;
00032     std::vector<double> m_gapFractionQ0SlicesDeltaY;
00033     std::vector<Histo1DPtr> _h_vetoPt;
00034     std::vector<Scatter2DPtr> _d_vetoPtGapFraction;
00035     std::vector<double> _h_vetoPtTotalSum;
00036 
00037     // Average NJet Vs DeltaY Setup
00038     int m_avgNJetDeltaYHistIndex;
00039     std::vector<double> m_avgNJetDeltaYSlices;
00040     std::vector<Profile1DPtr> _p_avgJetVsDeltaY;
00041 
00042     // Average NJet Vs PtBar Setup
00043     int m_avgNJetPtBarHistIndex;
00044     std::vector<double> m_avgNJetPtBarSlices;
00045     std::vector<Profile1DPtr> _p_avgJetVsPtBar;
00046   };
00047 
00048 
00049   class ATLAS_2011_S9126244 : public Analysis {
00050   public:
00051 
00052     /// Constructor
00053     ATLAS_2011_S9126244()
00054       : Analysis("ATLAS_2011_S9126244")
00055     {    }
00056 
00057 
00058   public:
00059 
00060     /// Book histograms and initialise projections before the run
00061     void init() {
00062 
00063       // Initialize the lone projection required
00064       addProjection(FastJets(FinalState(), FastJets::ANTIKT, 0.6), "AntiKtJets06");
00065 
00066       // Make Q0 bins 0->20 then 20 to 195.0 in steps of 5
00067       m_q0BinEdges += 0.0;
00068       for (unsigned int x=0; x<36; x++){
00069         m_q0BinEdges += 20.0 + x*5.0 ;
00070       }
00071 
00072       // Initialize plots for each selection type
00073       m_selectionPlots[0].intermediateHistName = "highestPt";
00074       m_selectionPlots[0].selectionType = 1;
00075       m_selectionPlots[0].m_gapFractionDeltaYHistIndex = 6;
00076       m_selectionPlots[0].m_gapFractionPtBarHistIndex = 1;
00077       m_selectionPlots[0].m_gapFractionQ0HistIndex = 13;
00078       m_selectionPlots[0].m_avgNJetDeltaYHistIndex = 37;
00079       m_selectionPlots[0].m_avgNJetPtBarHistIndex = 26;
00080       m_selectionPlots[0].m_gapFractionDeltaYSlices += 70.0, 90.0, 120.0, 150.0, 180.0, 210.0, 240.0, 270.0;
00081       m_selectionPlots[0].m_gapFractionPtBarSlices += 1.0, 2.0, 3.0, 4.0, 5.0, 6.0;
00082       m_selectionPlots[0].m_gapFractionQ0SlicesPtBar += 70.0, 90.0, 120.0, 150.0, 210.0, 240.0;
00083       m_selectionPlots[0].m_gapFractionQ0SlicesDeltaY += 2.0, 3.0, 4.0, 5.0;
00084       m_selectionPlots[0].m_avgNJetPtBarSlices += 1.0, 2.0, 3.0, 4.0, 5.0;
00085       m_selectionPlots[0].m_avgNJetDeltaYSlices += 70.0, 90.0, 120.0, 150.0, 180.0, 210.0, 240.0, 270.0;
00086       initializePlots(m_selectionPlots[0]);
00087 
00088       m_selectionPlots[1].intermediateHistName = "forwardBackward";
00089       m_selectionPlots[1].selectionType = 2;
00090       m_selectionPlots[1].m_gapFractionDeltaYHistIndex = 6;
00091       m_selectionPlots[1].m_gapFractionPtBarHistIndex = 1;
00092       m_selectionPlots[1].m_gapFractionQ0HistIndex = 13;
00093       m_selectionPlots[1].m_avgNJetDeltaYHistIndex = 37;
00094       m_selectionPlots[1].m_avgNJetPtBarHistIndex = 26;
00095       m_selectionPlots[1].m_gapFractionDeltaYSlices += 70.0, 90.0, 120.0, 150.0, 180.0, 210.0, 240.0, 270.0;
00096       m_selectionPlots[1].m_gapFractionPtBarSlices += 1.0, 2.0, 3.0, 4.0, 5.0, 6.0;
00097       m_selectionPlots[1].m_gapFractionQ0SlicesPtBar += 70.0, 90.0, 120.0, 150.0, 210.0, 240.0;
00098       m_selectionPlots[1].m_gapFractionQ0SlicesDeltaY += 2.0, 3.0, 4.0, 5.0;
00099       m_selectionPlots[1].m_avgNJetPtBarSlices += 1.0, 2.0, 3.0, 4.0, 5.0;
00100       m_selectionPlots[1].m_avgNJetDeltaYSlices += 70.0, 90.0, 120.0, 150.0, 180.0, 210.0, 240.0, 270.0;
00101       initializePlots(m_selectionPlots[1]);
00102 
00103       m_selectionPlots[2].intermediateHistName = "forwardBackward_PtBarVeto";
00104       m_selectionPlots[2].selectionType = 1;
00105       m_selectionPlots[2].m_gapFractionDeltaYHistIndex = 19;
00106       m_selectionPlots[2].m_avgNJetDeltaYHistIndex = 30;
00107       m_selectionPlots[2].m_gapFractionDeltaYSlices += 70.0, 90.0, 120.0, 150.0, 180.0, 210.0, 240.0, 270.0;
00108       m_selectionPlots[2].m_avgNJetDeltaYSlices += 70.0, 90.0, 120.0, 150.0, 180.0, 210.0, 240.0, 270.0;
00109       initializePlots(m_selectionPlots[2]);
00110     }
00111 
00112 
00113     void initializePlots(ATLAS_2011_S9126244_Plots& plots) {
00114 
00115       // Gap Fraction Vs DeltaY
00116       for (int x=0; x<((int)plots.m_gapFractionDeltaYSlices.size()-1); x++) {
00117         std::stringstream vetoHistName;
00118         std::stringstream inclusiveHistName;
00119 
00120         vetoHistName << "gapDeltaYVeto_" << plots.intermediateHistName << "_" << x;
00121         inclusiveHistName << "gapDeltaYInclusive_" << plots.intermediateHistName << "_" << x;
00122 
00123         plots._h_gapVsDeltaYVeto.addHistogram(plots.m_gapFractionDeltaYSlices[x], plots.m_gapFractionDeltaYSlices[x+1],
00124             bookHisto1D(plots.m_gapFractionDeltaYHistIndex+x, 1, plots.selectionType, vetoHistName.str()));
00125         plots._h_gapVsDeltaYInc.addHistogram(plots.m_gapFractionDeltaYSlices[x], plots.m_gapFractionDeltaYSlices[x+1],
00126             bookHisto1D(plots.m_gapFractionDeltaYHistIndex+x, 1, plots.selectionType, inclusiveHistName.str()));
00127       }
00128 
00129       // Average NJet Vs DeltaY
00130       for (int x=0; x<((int)plots.m_avgNJetDeltaYSlices.size()-1); x++) {
00131         plots._p_avgJetVsDeltaY += bookProfile1D(plots.m_avgNJetDeltaYHistIndex+x, 1, plots.selectionType);
00132       }
00133 
00134       // Gap Fraction Vs PtBar
00135       for (int x=0; x<((int)plots.m_gapFractionPtBarSlices.size()-1); x++) {
00136 
00137         std::stringstream vetoHistName;
00138         std::stringstream inclusiveHistName;
00139 
00140         vetoHistName << "gapPtBarVeto_" << plots.intermediateHistName << "_" << x;
00141         inclusiveHistName << "gapPtBarInclusive_" << plots.intermediateHistName << "_" << x;
00142 
00143         plots._h_gapVsPtBarVeto.addHistogram(plots.m_gapFractionPtBarSlices[x], plots.m_gapFractionPtBarSlices[x+1],
00144             bookHisto1D(plots.m_gapFractionPtBarHistIndex+x, 1, plots.selectionType, vetoHistName.str()));
00145         plots._h_gapVsPtBarInc.addHistogram(plots.m_gapFractionPtBarSlices[x], plots.m_gapFractionPtBarSlices[x+1],
00146             bookHisto1D(plots.m_gapFractionPtBarHistIndex+x, 1, plots.selectionType, inclusiveHistName.str()));
00147       }
00148 
00149       // Average NJet Vs PtBar
00150       for (int x=0; x<((int)plots.m_avgNJetPtBarSlices.size()-1); x++) {
00151         plots._p_avgJetVsPtBar += bookProfile1D(plots.m_avgNJetPtBarHistIndex+x, 1, plots.selectionType);
00152       }
00153 
00154       // Gap fraction Vs Q0
00155       int q0PlotCount = 0;
00156       for (int x=0; x<((int)plots.m_gapFractionQ0SlicesPtBar.size()/2); x++) {
00157         for (int y=0; y<((int)plots.m_gapFractionQ0SlicesDeltaY.size()/2); y++) {
00158           std::stringstream vetoPtHistName;
00159           std::stringstream vetoPtGapDataPointName;
00160 
00161           vetoPtHistName << "vetoPt_" << plots.intermediateHistName << "_" << q0PlotCount;
00162           vetoPtGapDataPointName << "gapQ0GapFractionDataPoints_" << plots.intermediateHistName << "_" << q0PlotCount;
00163 
00164           plots._h_vetoPt += bookHisto1D(vetoPtHistName.str(),
00165                                              m_q0BinEdges);
00166           plots._d_vetoPtGapFraction += bookScatter2D(plots.m_gapFractionQ0HistIndex+q0PlotCount, 1, plots.selectionType);
00167           plots._h_vetoPtTotalSum += 0.0;
00168           q0PlotCount++;
00169         }
00170       }
00171     }
00172 
00173 
00174     /// Perform the per-event analysis
00175     void analyze(const Event& event) {
00176       const double weight = event.weight();
00177 
00178       // Get minimal list of jets needed to be considered
00179       double minimumJetPtBar = 50.0*GeV; //of interval defining jets
00180 
00181       vector<FourMomentum> acceptJets;
00182       foreach (const Jet& jet, applyProjection<FastJets>(event, "AntiKtJets06").jetsByPt(20.0*GeV)) {
00183         if (fabs(jet.momentum().rapidity()) < 4.4) {
00184           acceptJets.push_back(jet.momentum());
00185         }
00186       }
00187 
00188       // If can't form an interval drop out of the analysis early
00189       if (acceptJets.size() < 2) {
00190         return;
00191       }
00192 
00193       // Analyze leading jet case
00194       if ((acceptJets[0].pT() + acceptJets[1].pT())/2.0 > minimumJetPtBar) {
00195         analyzeJets(acceptJets, m_selectionPlots[0], weight, 20.0*GeV);
00196       }
00197 
00198       // Re-order jets to have forward backward selection
00199       unsigned int minRapidityJet = 0;
00200       unsigned int maxRapidityJet = 0;
00201 
00202       for (size_t j=1; j<acceptJets.size(); j++) {
00203         if (acceptJets[j].rapidity() >
00204             acceptJets[maxRapidityJet].rapidity()) {
00205           maxRapidityJet=j;
00206         }
00207 
00208         if (acceptJets[j].rapidity() <
00209             acceptJets[minRapidityJet].rapidity()) {
00210           minRapidityJet=j;
00211         }
00212       }
00213 
00214       vector<FourMomentum> fwdBkwdJets;
00215       fwdBkwdJets.push_back(acceptJets[maxRapidityJet]);
00216       fwdBkwdJets.push_back(acceptJets[minRapidityJet]);
00217 
00218       for (size_t j=0; j<acceptJets.size(); j++) {
00219         if (j==minRapidityJet or j==maxRapidityJet){
00220           continue;
00221         }
00222         fwdBkwdJets.push_back(acceptJets[j]);
00223       }
00224 
00225       if ((fwdBkwdJets[0].pT() + fwdBkwdJets[1].pT())/2.0 > minimumJetPtBar) {
00226         //Use most forward/backward jets in rapidity to define the interval
00227         analyzeJets(fwdBkwdJets, m_selectionPlots[1], weight,
00228                     20.0*GeV);
00229 
00230         //As before but now using PtBar of interval to define veto threshold
00231         analyzeJets(fwdBkwdJets, m_selectionPlots[2], weight,
00232                     (fwdBkwdJets[0].pT()+fwdBkwdJets[1].pT())/2.0);
00233       }
00234     }
00235 
00236 
00237     // Fill plots!
00238     void analyzeJets(vector<FourMomentum>& jets, ATLAS_2011_S9126244_Plots& plots,
00239                      const double weight, double vetoPtThreshold){
00240       // Calculate the interval size, ptBar and veto Pt (if any)
00241       double intervalSize;
00242       double ptBar;
00243       double vetoPt = 0.0*GeV;
00244 
00245       intervalSize = fabs(jets[0].rapidity()-jets[1].rapidity());
00246       ptBar = (jets[0].pT()+jets[1].pT())/2.0;
00247 
00248       double minY;
00249       double maxY;
00250       if (jets[0].rapidity() > jets[1].rapidity()) {
00251         minY = jets[1].rapidity();
00252         maxY = jets[0].rapidity();
00253       } else {
00254         minY = jets[0].rapidity();
00255         maxY = jets[1].rapidity();
00256       }
00257 
00258       for (size_t j=2; j<jets.size(); j++){
00259         if (jets[j].rapidity() > minY &&
00260             jets[j].rapidity() < maxY &&
00261             jets[j].pT() > vetoPt){
00262           vetoPt = jets[j].pT();
00263         }
00264       }
00265 
00266       // Fill the gap fraction vs delta Y histograms
00267       plots._h_gapVsDeltaYInc.fill(ptBar/GeV, intervalSize, weight);
00268       if (vetoPt < vetoPtThreshold) {
00269         plots._h_gapVsDeltaYVeto.fill(ptBar/GeV, intervalSize, weight);
00270       }
00271 
00272       // Fill the gap fraction vs pt Bar histograms
00273       plots._h_gapVsPtBarInc.fill(intervalSize, ptBar/GeV,  weight);
00274       if (vetoPt < vetoPtThreshold) {
00275         plots._h_gapVsPtBarVeto.fill(intervalSize, ptBar/GeV, weight);
00276       }
00277 
00278       // Count the number of veto jets present
00279       int vetoJetsCount=0;
00280       for (size_t j=2; j<jets.size(); j++){
00281         if (jets[j].rapidity() > minY &&
00282             jets[j].rapidity() < maxY &&
00283             jets[j].pT() > vetoPtThreshold){
00284           vetoJetsCount += 1;
00285         }
00286       }
00287 
00288       // Fill the avg NJet, deltaY slices
00289       for (int i=0; i<(int)plots.m_avgNJetPtBarSlices.size()-1; i++) {
00290         if ( intervalSize >= plots.m_avgNJetPtBarSlices[i] &&
00291              intervalSize < plots.m_avgNJetPtBarSlices[i+1]) {
00292          plots._p_avgJetVsPtBar[i]->fill(ptBar/GeV, vetoJetsCount, weight);
00293         }
00294       }
00295 
00296       // Fill the avg NJet, ptBar slices
00297       for (int i=0; i<(int)plots.m_avgNJetDeltaYSlices.size()-1; i++) {
00298         if ( ptBar/GeV >= plots.m_avgNJetDeltaYSlices[i] &&
00299              ptBar/GeV < plots.m_avgNJetDeltaYSlices[i+1] ) {
00300           plots._p_avgJetVsDeltaY[i]->fill(intervalSize, vetoJetsCount, weight);
00301         }
00302       }
00303 
00304       // Fill the veto pt plots
00305       int q0PlotCount = 0;
00306       for (int x=0; x<((int)plots.m_gapFractionQ0SlicesPtBar.size()/2); x++) {
00307         for (int y=0; y<((int)plots.m_gapFractionQ0SlicesDeltaY.size()/2); y++) {
00308           // Check if it should be filled
00309           if ( ptBar/GeV < plots.m_gapFractionQ0SlicesPtBar[x*2] ||
00310                ptBar/GeV >= plots.m_gapFractionQ0SlicesPtBar[x*2+1] ) {
00311             q0PlotCount++;
00312             continue;
00313           }
00314 
00315           if ( intervalSize < plots.m_gapFractionQ0SlicesDeltaY[y*2] ||
00316                intervalSize >= plots.m_gapFractionQ0SlicesDeltaY[y*2+1] ) {
00317             q0PlotCount++;
00318             continue;
00319           }
00320 
00321           plots._h_vetoPt[q0PlotCount]->fill(vetoPt, weight);
00322           plots._h_vetoPtTotalSum[q0PlotCount] += weight;
00323 
00324           q0PlotCount++;
00325         }
00326       }
00327     }
00328 
00329 
00330     /// Derive final distributions for each selection
00331     void finalize() {
00332       /// \todo YODA divide and binEdges
00333       /// foreach (const ATLAS_2011_S9126244_Plots& plots, m_selectionPlots) {
00334       ///    // Calculate the gap fraction for each slice
00335       ///    for (size_t x=0; x<plots._h_gapVsDeltaYVeto.getHistograms().size(); x++) {
00336       ///      histogramFactory().divide(histoPath(makeAxisCode(plots.m_gapFractionDeltaYHistIndex+x, 1, plots.selectionType)),
00337       ///                                *(plots._h_gapVsDeltaYVeto.getHistograms()[x]),
00338       ///                                *(plots._h_gapVsDeltaYInc.getHistograms()[x]));
00339       ///      histogramFactory().destroy(plots._h_gapVsDeltaYVeto.getHistograms()[x]);
00340       ///      histogramFactory().destroy(plots._h_gapVsDeltaYInc.getHistograms()[x]);
00341       ///    }
00342 
00343       ///    for (size_t x=0; x<plots._h_gapVsPtBarVeto.getHistograms().size(); x++) {
00344       ///      histogramFactory().divide(histoPath(makeAxisCode(plots.m_gapFractionPtBarHistIndex+x, 1, plots.selectionType)),
00345       ///                                *(plots._h_gapVsPtBarVeto.getHistograms()[x]),
00346       ///                                *(plots._h_gapVsPtBarInc.getHistograms()[x]));
00347       ///      histogramFactory().destroy(plots._h_gapVsPtBarVeto.getHistograms()[x]);
00348       ///      histogramFactory().destroy(plots._h_gapVsPtBarInc.getHistograms()[x]);
00349       ///    }
00350 
00351       ///    for (size_t h=0; h<plots._d_vetoPtGapFraction.size(); h++) {
00352       ///      //Get the number of bins needed for this slice
00353       ///      const BinEdges q0Edges = binEdges(plots.m_gapFractionQ0HistIndex+h, 1, plots.selectionType);
00354       ///      finalizeQ0GapFraction(plots._h_vetoPtTotalSum[h],
00355       ///                            plots._d_vetoPtGapFraction[h],
00356       ///                            plots._h_vetoPt[h],
00357       ///                            q0Edges.size());
00358       ///    }
00359       /// }
00360     }
00361 
00362 
00363     /// \todo YODA
00364     /// void finalizeQ0GapFraction(double totalWeightSum,
00365     ///                            Scatter2DPtr gapFractionDP,
00366     ///                            Histo1DPtr vetoPtHist,
00367     ///                            int binNumber) {
00368     ///   double vetoPtWeightSum = 0.0;
00369     ///   for (int x=0; x<binNumber-1; x++) {
00370     ///    vetoPtWeightSum += vetoPtHist->bin(x).area();
00371 
00372     ///     //Alternatively try saving as data points
00373     ///     IDataPoint* currentPoint = gapFractionDP->point(x);
00374     ///     IMeasurement* xCoord = currentPoint->coordinate(0);
00375     ///     IMeasurement* yCoord = currentPoint->coordinate(1);
00376 
00377     ///     //Calculate the efficiency uncertainty
00378     ///     double efficiency = vetoPtWeightSum/totalWeightSum;
00379     ///     double efficiencyError = std::sqrt(efficiency*(1.0-efficiency)/totalWeightSum);
00380     ///     if (totalWeightSum==0.) efficiency = efficiencyError = 0.;
00381 
00382     ///     xCoord->setValue(m_q0BinEdges[x+1]);
00383     ///     xCoord->setErrorPlus(2.5);
00384     ///     xCoord->setErrorMinus(2.5);
00385     ///     yCoord->setValue(efficiency);
00386     ///     yCoord->setErrorPlus(efficiencyError);
00387     ///     yCoord->setErrorMinus(efficiencyError);
00388     ///   }
00389     ///   histogramFactory().destroy(vetoPtHist);
00390     /// }
00391 
00392 
00393   private:
00394 
00395     // Only need to define the q0 binning here
00396     std::vector<double> m_q0BinEdges;
00397 
00398 
00399   private:
00400 
00401     // Structure containing complete set of plots
00402     ATLAS_2011_S9126244_Plots m_selectionPlots[3];
00403 
00404   };
00405 
00406 
00407   // The hook for the plugin system
00408   DECLARE_RIVET_PLUGIN(ATLAS_2011_S9126244);
00409 
00410 }