rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2011_I889807

$B/\bar{B}$ angular correlations based on secondary vertex reconstruction in $pp$ collisions
Experiment: CMS (LHC)
Inspire ID: 889807
Status: VALIDATED
Authors:
  • Lukas Wehrli
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Inclusive QCD at 7 TeV. A $\hat{p_\perp}$ cut (or similar) is recommended since a leading jet $p_\perp > 56$ GeV is required.

The differential $B\bar{B}$ cross-section is measured as a function of the opening angle $\Delta{R}$ and $\Delta\phi$ using data collected with the CMS detector during 2010 and corresponding to an integrated luminosity of 3.1 pb$^{-1}$. The measurement is performed for three different event energy scales, characterized by the transverse momentum of the leading jet in the event (above 56 GeV, above 84 GeV and above 120 GeV). Simulated events are normalised in the region $\Delta{R} > 2.4$ and $\Delta\phi > 3/4\pi$ respectively.

Source code: CMS_2011_I889807.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/UnstableParticles.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6
  7namespace Rivet {
  8
  9
 10  /// B-Bbar angular correlations based on secondary vertex reconstruction
 11  class CMS_2011_I889807 : public Analysis {
 12  public:
 13
 14    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2011_I889807);
 15
 16
 17    /// @name Analysis methods
 18    /// @{
 19
 20    void init() {
 21      FinalState fs;
 22      FastJets jetproj(fs, JetAlg::ANTIKT, 0.5);
 23      jetproj.useInvisibles();
 24      declare(jetproj, "Jets");
 25
 26      UnstableParticles ufs;
 27      declare(ufs, "UFS");
 28
 29      // Book histograms
 30      book(_h_dsigma_dR_56GeV ,1,1,1);
 31      book(_h_dsigma_dR_84GeV ,2,1,1);
 32      book(_h_dsigma_dR_120GeV ,3,1,1);
 33      book(_h_dsigma_dPhi_56GeV ,4,1,1);
 34      book(_h_dsigma_dPhi_84GeV ,5,1,1);
 35      book(_h_dsigma_dPhi_120GeV ,6,1,1);
 36
 37      book(_c["MCDR56"],     "_MCDR56");
 38      book(_c["MCDR84"],     "_MCDR84");
 39      book(_c["MCDR120"],    "_MCDR120");
 40      book(_c["MCDPhi56"],   "_MCDPhi56");
 41      book(_c["MCDPhi84"],   "_MCDPhi84");
 42      book(_c["MCDPhi120"], "_MCDPhi120");
 43    }
 44
 45
 46    /// Perform the per-event analysis
 47    void analyze(const Event& event) {
 48
 49      const Jets& jets = apply<FastJets>(event,"Jets").jetsByPt();
 50      const UnstableParticles& ufs = apply<UnstableParticles>(event, "UFS");
 51
 52      // Find the leading jet pT and eta
 53      if (jets.size() == 0) vetoEvent;
 54      const double ljpT = jets[0].pT();
 55      const double ljeta = jets[0].eta();
 56      MSG_DEBUG("Leading jet pT / eta: " << ljpT << " / " << ljeta);
 57
 58      // Minimum requirement for event
 59      if (ljpT > 56*GeV && fabs(ljeta) < 3.0) {
 60        // Find B hadrons in event
 61        int nb = 0; //counters for all B and independent B hadrons
 62        double etaB1 = 7.7, etaB2 = 7.7;
 63        double phiB1 = 7.7, phiB2 = 7.7;
 64        double pTB1 = 7.7, pTB2 = 7.7;
 65
 66        for (const Particle& p : ufs.particles()) {
 67          int aid = p.abspid();
 68          if (aid/100 == 5 || aid/1000==5) {
 69            // 2J+1 == 1 (mesons) or 2 (baryons)
 70            if (aid%10 == 1 || aid%10 == 2) {
 71              // No B decaying to B
 72              if (aid != 5222 && aid != 5112 && aid != 5212 && aid != 5322) {
 73                if (nb==0) {
 74                  etaB1 = p.eta();
 75                  phiB1 = p.phi();
 76                  pTB1 = p.pT();
 77                } else if (nb==1) {
 78                  etaB2 = p.eta();
 79                  phiB2 = p.phi();
 80                  pTB2 = p.pT();
 81                }
 82                nb++;
 83              }
 84            }
 85            MSG_DEBUG("ID " << aid <<  " B hadron");
 86          }
 87        }
 88
 89        if (nb==2 && pTB1 > 15*GeV && pTB2 > 15*GeV && fabs(etaB1) < 2.0 && fabs(etaB2) < 2.0) {
 90          double dPhi = deltaPhi(phiB1, phiB2);
 91          double dR = deltaR(etaB1, phiB1, etaB2, phiB2);
 92          MSG_DEBUG("DR/DPhi " << dR << " " << dPhi);
 93
 94          // MC counters
 95          if (dR > 2.4) _c["MCDR56"]->fill();
 96          if (dR > 2.4 && ljpT > 84*GeV) _c["MCDR84"]->fill();
 97          if (dR > 2.4 && ljpT > 120*GeV) _c["MCDR120"]->fill();
 98          if (dPhi > 3.*PI/4.) _c["MCDPhi56"]->fill();
 99          if (dPhi > 3.*PI/4. && ljpT > 84*GeV) _c["MCDPhi84"]->fill();
100          if (dPhi > 3.*PI/4. && ljpT > 120*GeV) _c["MCDPhi120"]->fill();
101
102          _h_dsigma_dR_56GeV->fill(dR);
103          if (ljpT > 84*GeV) _h_dsigma_dR_84GeV->fill(dR);
104          if (ljpT > 120*GeV) _h_dsigma_dR_120GeV->fill(discEdge(dR));
105          _h_dsigma_dPhi_56GeV->fill(dPhi);
106          if (ljpT > 84*GeV) _h_dsigma_dPhi_84GeV->fill(dPhi);
107          if (ljpT > 120*GeV) _h_dsigma_dPhi_120GeV->fill(dPhi);
108        }
109      }
110    }
111
112
113    /// Normalise histograms etc., after the run
114    void finalize() {
115      MSG_DEBUG("crossSection " << crossSection()/picobarn << " sumOfWeights " << sumOfWeights());
116
117      // Hardcoded bin widths
118      double DRbin = 0.4;
119      double DPhibin = PI/8.0;
120      // Find out the correct numbers
121      double nDataDR56 = 25862.20;
122      double nDataDR84 = 5675.55;
123      double nDataDR120 = 1042.72;
124      double nDataDPhi56 = 24220.00;
125      double nDataDPhi84 = 4964.00;
126      double nDataDPhi120 = 919.10;
127      double normDR56 = safediv(nDataDR56, dbl(*_c["MCDR56"]), crossSection()/picobarn/sumOfWeights());
128      double normDR84 = safediv(nDataDR84, dbl(*_c["MCDR84"]), crossSection()/picobarn/sumOfWeights());
129      double normDR120 = safediv(nDataDR120, dbl(*_c["MCDR120"]), crossSection()/picobarn/sumOfWeights());
130      double normDPhi56 = safediv(nDataDPhi56, dbl(*_c["MCDPhi56"]), crossSection()/picobarn/sumOfWeights());
131      double normDPhi84 = safediv(nDataDPhi84, dbl(*_c["MCDPhi84"]), crossSection()/picobarn/sumOfWeights());
132      double normDPhi120 = safediv(nDataDPhi120, dbl(*_c["MCDPhi120"]), crossSection()/picobarn/sumOfWeights());
133      scale(_h_dsigma_dR_56GeV, normDR56*DRbin);
134      scale(_h_dsigma_dR_84GeV, normDR84*DRbin);
135      scale(_h_dsigma_dR_120GeV, normDR120*DRbin);
136      scale(_h_dsigma_dPhi_56GeV, normDPhi56*DPhibin);
137      scale(_h_dsigma_dPhi_84GeV, normDPhi84*DPhibin);
138      scale(_h_dsigma_dPhi_120GeV, normDPhi120*DPhibin);
139      for (auto& b : _h_dsigma_dR_120GeV->bins()) {
140        b.scaleW(1.0/_rapaxis.width(b.index()));
141      }
142    }
143
144    string discEdge(const double value) const {
145      const size_t idx = _rapaxis.index(value);
146      return _h_dsigma_dR_120GeV->bin(idx).xEdge();
147    }
148
149    /// @}
150
151
152  private:
153
154    /// Counters
155    map<string, CounterPtr> _c;
156
157    /// @name Histograms
158    /// @{
159    Histo1DPtr _h_dsigma_dR_56GeV, _h_dsigma_dR_84GeV;
160    BinnedHistoPtr<string> _h_dsigma_dR_120GeV; // Why does this one have a different type?
161    Histo1DPtr _h_dsigma_dPhi_56GeV, _h_dsigma_dPhi_84GeV, _h_dsigma_dPhi_120GeV;
162    YODA::Axis<double> _rapaxis{ 0.0, 0.4, 0.8, 1.2, 1.6, 2.0, 2.4, 2.8, 3.2, 3.6, 4.0 };
163    /// @}
164
165  };
166
167
168
169  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2011_I889807, CMS_2011_S8973270);
170
171}