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