rivet is hosted by Hepforge, IPPP Durham
ATLAS_2012_I1183818.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/ChargedFinalState.hh"
00005 #include "Rivet/Projections/FastJets.hh"
00006 
00007 
00008 /// @author Peter Wijeratne <paw@hep.ucl.ac.uk>
00009 /// @author Robindra Prabhu <prabhu@cern.ch>
00010 namespace Rivet {
00011 
00012   // A very basic analysis sensitive to ET flow in minbias and dijet events
00013   class ATLAS_2012_I1183818 : public Analysis {
00014 
00015   public:
00016 
00017     ATLAS_2012_I1183818()
00018       : Analysis("ATLAS_2012_I1183818")
00019     {}
00020 
00021 
00022   public:
00023 
00024     void init() {
00025 
00026       const FinalState cnfs(-4.8, 4.8, 0*MeV);
00027       const ChargedFinalState cfs(-2.5, 2.5, 250*MeV);
00028       addProjection(cnfs, "FS");
00029       addProjection(cfs, "CFS");
00030 
00031       const FastJets jetsAntiKt4(cnfs, FastJets::ANTIKT, 0.4);
00032       addProjection(jetsAntiKt4, "AntiKt4Jets");
00033 
00034       // ------- MINBIAS HISTOGRAMS --------
00035       //
00036       // MB event counter
00037       m_chargedEvents = 0.0;
00038 
00039       _h_ETflowEta = bookHisto1D(1, 1, 1);
00040       _h_SumETbin1 = bookHisto1D(3, 1, 1);
00041       _h_SumETbin2 = bookHisto1D(4, 1, 1);
00042       _h_SumETbin3 = bookHisto1D(5, 1, 1);
00043       _h_SumETbin4 = bookHisto1D(6, 1, 1);
00044       _h_SumETbin5 = bookHisto1D(7, 1, 1);
00045       _h_SumETbin6 = bookHisto1D(8, 1, 1);
00046 
00047       // ------- DIJET HISTOGRAMS --------
00048       //
00049       // Dijet event counter
00050       m_events_dijets = 0.0;
00051 
00052       // sumET
00053       _h_transETflowEta = bookHisto1D( 2, 1, 1);
00054       _h_transSumETbin1 = bookHisto1D( 9, 1, 1);
00055       _h_transSumETbin2 = bookHisto1D(10, 1, 1);
00056       _h_transSumETbin3 = bookHisto1D(11, 1, 1);
00057       _h_transSumETbin4 = bookHisto1D(12, 1, 1);
00058       _h_transSumETbin5 = bookHisto1D(13, 1, 1);
00059       _h_transSumETbin6 = bookHisto1D(14, 1, 1);
00060 
00061 
00062     }
00063 
00064 
00065     void analyze(const Event& event) {
00066 
00067       const double weight = event.weight();
00068 
00069       const FinalState& cfs = applyProjection<FinalState>(event, "CFS");
00070 
00071       bool isCharged = false;
00072       if (cfs.size() >= 2) {  // event selection: > 2 charged particles with pT > 250.MeV and |eta| < 2.5
00073         isCharged = true;
00074         m_chargedEvents += weight;
00075       }
00076 
00077       const FinalState& cnfs = applyProjection<FinalState>(event, "FS");
00078 
00079       Particles particles;
00080       foreach( const Particle& p, cnfs.particles() ) {
00081         // enforce truth selection representing detected particle sensitivity
00082         double pp = p.momentum().p().mod();
00083         if (PID::threeCharge(p.pdgId()) != 0 && pp < 0.5*GeV) continue;
00084         if (PID::threeCharge(p.pdgId()) == 0 && pp < 0.2*GeV) continue;
00085 
00086         particles.push_back(p);
00087       }
00088 
00089 
00090       // get jets
00091       const FastJets& jetsAntiKt4 = applyProjection<FastJets>(event, "AntiKt4Jets");
00092       const Jets& jets = jetsAntiKt4.jetsByPt(20.0*GeV);
00093 
00094       // initialise sumET variables
00095       double sumETbin1 = 0;
00096       double sumETbin2 = 0;
00097       double sumETbin3 = 0;
00098       double sumETbin4 = 0;
00099       double sumETbin5 = 0;
00100       double sumETbin6 = 0;
00101 
00102       // if (passes event selection)
00103       if (isCharged) {
00104 
00105         foreach( const Particle& p, particles ) {
00106 
00107           ///calculate variables
00108           double ET = p.momentum().Et()/GeV;
00109           double eta = fabs(p.eta());
00110 
00111           // fill histograms
00112           _h_ETflowEta->fill(eta, weight*ET);
00113 
00114           if      (eta <  0.8) sumETbin1 += ET;
00115           else if (eta <  1.6) sumETbin2 += ET;
00116           else if (eta <  2.4) sumETbin3 += ET;
00117           else if (eta <  3.2) sumETbin4 += ET;
00118           else if (eta <  4.0) sumETbin5 += ET;
00119           else if (eta <= 4.8) sumETbin6 += ET;
00120 
00121         } // end of foreach
00122 
00123         _h_SumETbin1->fill(sumETbin1, weight);
00124         _h_SumETbin2->fill(sumETbin2, weight);
00125         _h_SumETbin3->fill(sumETbin3, weight);
00126         _h_SumETbin4->fill(sumETbin4, weight);
00127         _h_SumETbin5->fill(sumETbin5, weight);
00128         _h_SumETbin6->fill(sumETbin6, weight);
00129       }
00130 
00131       // --- do dijet analysis ---
00132 
00133       if ( jets.size() >= 2                       && // require at least two jets
00134            jets[0].momentum().Et() >= 20.*GeV     && // require two leading jets to pass ET cuts
00135            jets[1].momentum().Et() >= 20.*GeV     &&
00136            fabs(jets[0].eta()) < 2.5   && // require leading jets to be central
00137            fabs(jets[1].eta()) < 2.5   &&
00138            deltaPhi(jets[0], jets[1]) > 2.5       && // require back-to-back topology
00139            jets[1].momentum().Et()/jets[0].momentum().Et() >= 0.5) { //require ET-balance
00140 
00141         // found an event that satisfies dijet selection, now fill histograms...
00142         // initialise dijet sumET variables
00143         double trans_sumET_bin1 = 0.;
00144         double trans_sumET_bin2 = 0.;
00145         double trans_sumET_bin3 = 0.;
00146         double trans_sumET_bin4 = 0.;
00147         double trans_sumET_bin5 = 0.;
00148         double trans_sumET_bin6 = 0.;
00149 
00150         m_events_dijets += weight;
00151 
00152         // loop over all particles and check their relation to leading jet
00153         foreach( const Particle& particle, particles ) {
00154 
00155           // calculate variables
00156           double dPhi = deltaPhi( jets[0], particle.momentum() );
00157           double ET   = particle.momentum().Et()/GeV;
00158           double eta  = fabs(particle.eta());
00159 
00160           // Transverse region
00161           if ( dPhi > 1./3.*M_PI && dPhi < 2./3.*M_PI ) {
00162             _h_transETflowEta->fill( eta, weight*ET );
00163             if      (eta <  0.8) { trans_sumET_bin1 += ET; }
00164             else if (eta <  1.6) { trans_sumET_bin2 += ET; }
00165             else if (eta <  2.4) { trans_sumET_bin3 += ET; }
00166             else if (eta <  3.2) { trans_sumET_bin4 += ET; }
00167             else if (eta <  4.0) { trans_sumET_bin5 += ET; }
00168             else if (eta <= 4.8) { trans_sumET_bin6 += ET; }
00169           }
00170 
00171         } // end loop over particles
00172 
00173         _h_transSumETbin1->fill( trans_sumET_bin1, weight);
00174         _h_transSumETbin2->fill( trans_sumET_bin2, weight);
00175         _h_transSumETbin3->fill( trans_sumET_bin3, weight);
00176         _h_transSumETbin4->fill( trans_sumET_bin4, weight);
00177         _h_transSumETbin5->fill( trans_sumET_bin5, weight);
00178         _h_transSumETbin6->fill( trans_sumET_bin6, weight);
00179       } // end of dijet selection cuts
00180 
00181     }
00182 
00183 
00184     void finalize() {
00185       /// several scale factors here:
00186       /// 1. nEvents (m_chargedEvents)
00187       /// 2. phase-space (2*M_PI)
00188       /// 3. double binning due to symmetrisation (2)
00189       scale( _h_ETflowEta, 1./m_chargedEvents/(4.*M_PI) );
00190       scale( _h_SumETbin1, 1./m_chargedEvents );
00191       scale( _h_SumETbin2, 1./m_chargedEvents );
00192       scale( _h_SumETbin3, 1./m_chargedEvents );
00193       scale( _h_SumETbin4, 1./m_chargedEvents );
00194       scale( _h_SumETbin5, 1./m_chargedEvents );
00195       scale( _h_SumETbin6, 1./m_chargedEvents );
00196 
00197       //Dijet analysis
00198 
00199       // Dijet scale factors:
00200       //1. number of events passing dijet selection
00201       //2. phase-space: 1. / 2/3*M_PI
00202       //3. double binning due to symmetrisation in |eta| plots : 1/2
00203       scale( _h_transETflowEta, 1./m_events_dijets * 1./(4./3.*M_PI) );
00204       scale( _h_transSumETbin1, 1./m_events_dijets );
00205       scale( _h_transSumETbin2, 1./m_events_dijets );
00206       scale( _h_transSumETbin3, 1./m_events_dijets );
00207       scale( _h_transSumETbin4, 1./m_events_dijets );
00208       scale( _h_transSumETbin5, 1./m_events_dijets );
00209       scale( _h_transSumETbin6, 1./m_events_dijets );
00210     }
00211 
00212   private:
00213 
00214     // Event counts
00215     double m_chargedEvents;
00216     double m_events_dijets;
00217 
00218     // Minbias-analysis: variable + histograms
00219     Histo1DPtr _h_ETflowEta;
00220     Histo1DPtr _h_SumETbin1;
00221     Histo1DPtr _h_SumETbin2;
00222     Histo1DPtr _h_SumETbin3;
00223     Histo1DPtr _h_SumETbin4;
00224     Histo1DPtr _h_SumETbin5;
00225     Histo1DPtr _h_SumETbin6;
00226 
00227     // Transverse region
00228     Histo1DPtr _h_transETflowEta;
00229     Histo1DPtr _h_transSumETbin1;
00230     Histo1DPtr _h_transSumETbin2;
00231     Histo1DPtr _h_transSumETbin3;
00232     Histo1DPtr _h_transSumETbin4;
00233     Histo1DPtr _h_transSumETbin5;
00234     Histo1DPtr _h_transSumETbin6;
00235 
00236   };
00237 
00238 
00239   DECLARE_RIVET_PLUGIN(ATLAS_2012_I1183818);
00240 
00241 }