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