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