rivet is hosted by Hepforge, IPPP Durham
ATLAS_2014_I1307243.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Tools/Logging.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 #include "Rivet/Tools/BinnedHistogram.hh"
00006 // // STL
00007 // #include <algorithm>
00008 // #include <map>
00009 
00010 namespace Rivet {
00011 
00012   /// @brief ATLAS azimuthal decorrelation with jet veto analysis
00013   /// @author James Robinson <james.robinson@cern.ch>
00014   class ATLAS_2014_I1307243 : public Analysis {
00015   public:
00016 
00017     /// Constructor
00018     ATLAS_2014_I1307243()
00019       : Analysis( "ATLAS_2014_I1307243" )
00020       , _dy_max(8)
00021       , _nEventsInAcceptance(0)
00022       , _sumOfAcceptedWeights(0.)
00023     {
00024       // Cannot do this in initialiser list without C++11
00025       _fiducialRegions += 2010, 2011;
00026       _vetoScale += 20*GeV, 30*GeV;
00027       _yFiducial += 4.4, 2.4;
00028       _gapCategories += "inclusive", "gap";
00029     }
00030 
00031     /// Book histograms and initialise projections before the run
00032     void init() {
00033       /// Initialise and register projections here
00034       FinalState fs;
00035       FastJets fastJets(fs, FastJets::ANTIKT, 0.6);
00036       fastJets.useInvisibles(true);
00037       addProjection(fastJets, "AntiKt6JetsWithInvisibles");
00038 
00039 
00040       /// Book histograms
00041       foreach( const std::string &gapCategory, _gapCategories ) {
00042         int gapCategoryOffset( gapCategory == "inclusive" ? 0 : 1 );
00043 
00044         // Temporary inclusive and gap histograms
00045         _h_tmp_events_dy[gapCategory] = bookHisto1D(1, 1, 1);
00046         _h_tmp_events_dy[gapCategory]->setPath("/TMP/" + toString(gapCategory) + "_events_dy");
00047         _h_tmp_events_pTbar[gapCategory] = bookHisto1D(2, 1, 1);
00048         _h_tmp_events_pTbar[gapCategory]->setPath("/TMP/" + toString(gapCategory) + "_events_pTbar");
00049 
00050 
00051         // Azimuthal moment histograms
00052         _h_profiled_cosDeltaPhi_dy[gapCategory]       = bookProfile1D( 5+4*gapCategoryOffset, 1, 1);
00053         _h_profiled_cosDeltaPhi_pTbar[gapCategory]    = bookProfile1D( 6+4*gapCategoryOffset, 1, 1);
00054         _h_C2C1_dy[gapCategory]                       = bookScatter2D( 7+4*gapCategoryOffset, 1, 1, false);
00055         _h_C2C1_pTbar[gapCategory]                    = bookScatter2D( 8+4*gapCategoryOffset, 1, 1, false);
00056         _h_profiled_cosTwoDeltaPhi_dy[gapCategory]    = bookProfile1D(37+2*gapCategoryOffset, 1, 1);
00057         _h_profiled_cosTwoDeltaPhi_pTbar[gapCategory] = bookProfile1D(38+2*gapCategoryOffset, 1, 1);
00058 
00059         // Gap fraction vs. Q0 and cross-section in dy slices
00060         for (size_t dyLow = 0; dyLow < _dy_max; ++dyLow ) {
00061           Histo1DPtr _h_tmp_events_Q0_single_dySlice = bookHisto1D( 29+dyLow, 1, 1);
00062           _h_tmp_events_Q0_single_dySlice->setPath("/TMP/" + toString(gapCategory) + "_events_dySlice_" + toString(dyLow) + "_" + toString(dyLow+1) + "_Q0");
00063           _h_tmp_events_Q0_dySlices[gapCategory].addHistogram( dyLow, dyLow+1, _h_tmp_events_Q0_single_dySlice );
00064           _h_crossSection_dphi_dySlices[gapCategory].addHistogram( dyLow, dyLow+1, bookHisto1D( 13+(_dy_max*gapCategoryOffset)+dyLow, 1, 1));
00065         }
00066 
00067       }
00068 
00069       // Number of jets in rapidity interval
00070       _h_profiled_nJets_rapidity_interval_dy    = bookProfile1D( 3, 1, 1);
00071       _h_profiled_nJets_rapidity_interval_pTbar = bookProfile1D( 4, 1, 1);
00072     }
00073 
00074 
00075     /// Perform the per-event analysis
00076     void analyze(const Event& event) {
00077 
00078       // Get the event weight
00079       const double weight( event.weight() );
00080       bool eventAccepted( false );
00081 
00082       for (int iFiducialRegion = 0; iFiducialRegion < 2; ++iFiducialRegion ) {
00083 
00084         // Retrieve all anti-kt R=0.6 jets above _pTMin and inside |_yFiducial|
00085         const Jets akt6Jets = applyProjection<JetAlg>(event, "AntiKt6JetsWithInvisibles").jetsByPt( Cuts::absrap < _yFiducial.at(iFiducialRegion) );
00086         // If there are fewer than 2 jets then bail
00087         if ( akt6Jets.size() < 2 ) { vetoEvent; }
00088 
00089         // Require jets to be above {60, 50} GeV
00090         if ( akt6Jets.at(0).momentum().pT() < 60.*GeV || akt6Jets.at(1).momentum().pT() < 50.*GeV ) { vetoEvent; }
00091 
00092         // Identify gap boundaries
00093         double yMin( std::min( akt6Jets.at(0).momentum().rapidity(), akt6Jets.at(1).momentum().rapidity() ) );
00094         double yMax( std::max( akt6Jets.at(0).momentum().rapidity(), akt6Jets.at(1).momentum().rapidity() ) );
00095 
00096         // Determine azimuthal decorrelation quantities
00097         const double dy( yMax - yMin );
00098         const double dphi( mapAngle0ToPi( akt6Jets.at(0).momentum().phi() - akt6Jets.at(1).momentum().phi() ) );
00099         const double pTbar( (akt6Jets.at(0).momentum().pT() + akt6Jets.at(1).momentum().pT())/2.0 );
00100 
00101         // Impose minimum dy for the 2011 phase space
00102         if ( _fiducialRegions.at(iFiducialRegion) == 2011 && dy < 1.0 ) { vetoEvent; }
00103 
00104         // Construct gap candidates sample
00105         Jets gapCandidates;
00106         foreach( const Jet &j, akt6Jets ) {
00107           if ( inRange( j.momentum().rapidity(), yMin, yMax, OPEN, OPEN ) ) {
00108             gapCandidates.push_back( j );
00109           }
00110         }
00111 
00112         // Determine gap quantities
00113         unsigned int nJets_rapidity_interval( 0 );
00114         double maximumGapQ0( 0. );
00115         foreach( const Jet &jet, gapCandidates ) {
00116           const double pT( jet.momentum().pT() );
00117           if ( pT > _vetoScale.at(iFiducialRegion) ) { ++nJets_rapidity_interval; }
00118           if ( pT > maximumGapQ0 ) { maximumGapQ0 = pT; }
00119         }
00120 
00121         // Fill histograms
00122         if ( weight < 0.0 ) {
00123           MSG_DEBUG( "Negative weight " << weight << "found!" );
00124         }
00125         fillHistograms( _fiducialRegions.at(iFiducialRegion), dy, pTbar, dphi, nJets_rapidity_interval, maximumGapQ0, weight );
00126         eventAccepted = true;
00127       }
00128 
00129       // Count number of accepted events
00130       if ( eventAccepted ) {
00131         _nEventsInAcceptance++;
00132         _sumOfAcceptedWeights += weight;
00133       }
00134       return;
00135     }
00136 
00137     void fillHistograms( const unsigned int &fiducialRegion, const double &dy, const double &pTbar, const double &dphi, const double &nJets_rapidity_interval, const double &maximumGapQ0, const double &weight) {
00138       // Determine gap category
00139       std::vector<std::string> eventGapCategories = boost::assign::list_of("inclusive");
00140       if ( nJets_rapidity_interval == 0 ) { eventGapCategories.push_back("gap"); }
00141 
00142       // Fill histograms relevant for comparison with 2010 data
00143       if ( fiducialRegion == _fiducialRegions.at(0) ) {
00144         // Fill inclusive and gap histograms
00145         foreach( const std::string &gapCategory, eventGapCategories ) {
00146           _h_tmp_events_dy[gapCategory]->fill( dy, weight );
00147           _h_crossSection_dphi_dySlices[gapCategory].fill( dy, dphi / M_PI, weight );
00148           _h_profiled_cosDeltaPhi_dy[gapCategory]->fill( dy, cos(M_PI - dphi), weight );
00149           _h_profiled_cosTwoDeltaPhi_dy[gapCategory]->fill( dy, cos(2.0 * dphi), weight );
00150         }
00151         // Fill profiled nJets_rapidity_interval
00152         _h_profiled_nJets_rapidity_interval_dy->fill( dy, nJets_rapidity_interval, weight );
00153         // Fill Q0 histograms - can fill multiple points per event
00154         foreach( const YODA::HistoBin1D Q0_bin, _h_tmp_events_Q0_dySlices["inclusive"].getHistograms().at(0)->bins() ) {
00155           const double Q0( Q0_bin.xMid() );
00156           _h_tmp_events_Q0_dySlices["inclusive"].fill(dy, Q0, weight);
00157           if ( maximumGapQ0 <= Q0 ) { _h_tmp_events_Q0_dySlices["gap"].fill(dy, Q0, weight); }
00158         }
00159 
00160       // Fill histograms relevant for comparison with 2011 data
00161       } else if ( fiducialRegion == _fiducialRegions.at(1) ) {
00162         // Fill inclusive and gap histograms
00163         foreach( const std::string &gapCategory, eventGapCategories ) {
00164           _h_tmp_events_pTbar[gapCategory]->fill( pTbar, weight );
00165           _h_profiled_cosDeltaPhi_pTbar[gapCategory]->fill( pTbar, cos(M_PI - dphi), weight );
00166           _h_profiled_cosTwoDeltaPhi_pTbar[gapCategory]->fill( pTbar, cos(2.0 * dphi), weight );
00167         }
00168         // Fill profiled nJets_rapidity_interval
00169         _h_profiled_nJets_rapidity_interval_pTbar->fill( pTbar, nJets_rapidity_interval, weight );
00170       }
00171       return;
00172     }
00173 
00174     /// Normalise histograms etc., after the run
00175     void finalize() {
00176       /// Normalise cross-section plots to correct cross-section
00177       const double xs_pb( crossSection()/picobarn );
00178       const double nEventsGenerated( sumOfWeights() );
00179       double xs_norm_factor( xs_pb/nEventsGenerated );
00180       const double ySpan(1); // all dy spans are 1
00181       foreach( const string &gapCategory, _gapCategories ) {
00182         _h_crossSection_dphi_dySlices[gapCategory].scale(xs_norm_factor/ySpan/M_PI, this);
00183       }
00184 
00185       /// Create C2/C1 scatter from profiles
00186       foreach( const std::string &gapCategory, _gapCategories ) {
00187         divide( _h_profiled_cosTwoDeltaPhi_dy[gapCategory], _h_profiled_cosDeltaPhi_dy[gapCategory], _h_C2C1_dy[gapCategory] );
00188         divide( _h_profiled_cosTwoDeltaPhi_pTbar[gapCategory], _h_profiled_cosDeltaPhi_pTbar[gapCategory], _h_C2C1_pTbar[gapCategory] );
00189       }
00190 
00191       /// Fill simple gap fractions
00192       Scatter2DPtr h_gap_fraction_dy    = bookScatter2D( 1, 1, 1);
00193       Scatter2DPtr h_gap_fraction_pTbar = bookScatter2D( 2, 1, 1, false);
00194       Rivet::Analysis::efficiency(    *_h_tmp_events_dy["gap"].get(),    *_h_tmp_events_dy["inclusive"].get(),    h_gap_fraction_dy );
00195       Rivet::Analysis::efficiency( *_h_tmp_events_pTbar["gap"].get(), *_h_tmp_events_pTbar["inclusive"].get(), h_gap_fraction_pTbar );
00196 
00197       /// Register and fill Q0 gap fractions
00198       for( unsigned int dyLow = 0; dyLow < _dy_max; ++dyLow ) {
00199         Scatter2DPtr h_gap_fraction_Q0 = bookScatter2D( 29+dyLow, 1, 1, false);
00200         Rivet::Analysis::efficiency( *_h_tmp_events_Q0_dySlices["gap"].getHistograms().at(dyLow).get(), *_h_tmp_events_Q0_dySlices["inclusive"].getHistograms().at(dyLow).get(), h_gap_fraction_Q0 );
00201       }
00202 
00203       /// Print summary information
00204       MSG_INFO( "Cross-Section/pb     : " << xs_pb );
00205       MSG_INFO( "Sum of weights       : " << sumOfWeights() << "  (" << _sumOfAcceptedWeights << " accepted)" );
00206       MSG_INFO( "nEvents              : " << numEvents() << "  (" << _nEventsInAcceptance << " accepted)" );
00207       MSG_INFO( "Normalisation factor : " << xs_norm_factor );
00208       return;
00209     }
00210 
00211   private:
00212     /// Member variables
00213     std::vector<unsigned int> _fiducialRegions;
00214     std::vector<double> _vetoScale, _yFiducial;
00215     std::vector<std::string> _gapCategories;
00216 
00217     unsigned int _dy_max;
00218     unsigned int _nEventsInAcceptance;
00219     double _sumOfAcceptedWeights;
00220 
00221     /// Histograms
00222     // Number of events : gap and non-gap : necessary input for gap fractions but not saved as output
00223     std::map< std::string, Histo1DPtr > _h_tmp_events_dy;
00224     std::map< std::string, Histo1DPtr > _h_tmp_events_pTbar;
00225     std::map< std::string, BinnedHistogram<double> > _h_tmp_events_Q0_dySlices;
00226 
00227     // Number of jets in rapidity interval
00228     Profile1DPtr _h_profiled_nJets_rapidity_interval_dy;
00229     Profile1DPtr _h_profiled_nJets_rapidity_interval_pTbar;
00230 
00231     // Azimuthal moment histograms
00232     std::map< std::string, Profile1DPtr > _h_profiled_cosDeltaPhi_dy;
00233     std::map< std::string, Profile1DPtr > _h_profiled_cosDeltaPhi_pTbar;
00234     std::map< std::string, Profile1DPtr > _h_profiled_cosTwoDeltaPhi_dy;
00235     std::map< std::string, Profile1DPtr > _h_profiled_cosTwoDeltaPhi_pTbar;
00236     std::map< std::string, Scatter2DPtr > _h_C2C1_dy;
00237     std::map< std::string, Scatter2DPtr > _h_C2C1_pTbar;
00238 
00239     // Cross-section vs. #Delta#phi
00240     std::map< std::string, BinnedHistogram<double> > _h_crossSection_dphi_dySlices;
00241   };
00242 
00243   // The hook for the plugin system
00244   DECLARE_RIVET_PLUGIN(ATLAS_2014_I1307243);
00245 
00246 }