1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
| // -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ZFinder.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
//#include "Rivet/ParticleName.hh"
namespace Rivet {
/// Underlying event activity in the Drell-Yan process at 13 TeV
class CMS_2017_I1635889 : public Analysis {
public:
/// Constructor
CMS_2017_I1635889()
: Analysis("CMS_2017_I1635889")
{ }
/// Initialization
void init() {
/// @note Using a bare muon Z (but with a clustering radius!?)
Cut cut = Cuts::abseta < 2.4 && Cuts::pT > 10*GeV;
ZFinder zfinder(FinalState(), cut, PID::MUON, 81*GeV, 101*GeV, 0.2, ZFinder::ClusterPhotons::NONE);
declare(zfinder, "ZFinder");
ChargedFinalState cfs(zfinder.remainingFinalState());
declare(cfs, "cfs");
book(_h_Nchg_towards_pTmumu , 1, 1, 1);
book(_h_Nchg_transverse_pTmumu , 2, 1, 1);
book(_h_Nchg_away_pTmumu , 3, 1, 1);
book(_h_pTsum_towards_pTmumu , 4, 1, 1);
book(_h_pTsum_transverse_pTmumu , 5, 1, 1);
book(_h_pTsum_away_pTmumu , 6, 1, 1);
}
/// Perform the per-event analysis
void analyze(const Event& event) {
const ZFinder& zfinder = applyProjection<ZFinder>(event, "ZFinder");
if (zfinder.bosons().size() != 1) vetoEvent;
if (zfinder.constituents()[0].pT()<20 && zfinder.constituents()[1].pT()<20)vetoEvent;
//std::cout<<"pt[0] = "<<zfinder.constituents()[0].pT()<<"pt[1] = "<<zfinder.constituents()[1].pT()<<std::endl;
double Zpt = zfinder.bosons()[0].pT()/GeV;
double Zphi = zfinder.bosons()[0].phi();
//double Zmass = zfinder.bosons()[0].mass()/GeV;
Particles particles = applyProjection<ChargedFinalState>(event, "cfs").particlesByPt(Cuts::pT > 0.5*GeV && Cuts::abseta <2.0);
int nTowards = 0;
int nTransverse = 0;
int nAway = 0;
double ptSumTowards = 0;
double ptSumTransverse = 0;
double ptSumAway = 0;
for (const Particle& p : particles) {
double dphi = fabs(deltaPhi(Zphi, p.phi()));
double pT = p.pT();
if ( dphi < M_PI/3 ) {
nTowards++;
ptSumTowards += pT;
} else if ( dphi < 2.*M_PI/3 ) {
nTransverse++;
ptSumTransverse += pT;
} else {
nAway++;
ptSumAway += pT;
}
} // Loop over particles
const double area = 8./3.*M_PI;
_h_Nchg_towards_pTmumu-> fill(Zpt, 1./area * nTowards);
_h_Nchg_transverse_pTmumu-> fill(Zpt, 1./area * nTransverse);
_h_Nchg_away_pTmumu-> fill(Zpt, 1./area * nAway);
_h_pTsum_towards_pTmumu-> fill(Zpt, 1./area * ptSumTowards);
_h_pTsum_transverse_pTmumu-> fill(Zpt, 1./area * ptSumTransverse);
_h_pTsum_away_pTmumu-> fill(Zpt, 1./area * ptSumAway);
}
/// Normalise histograms etc., after the run
void finalize() {
}
private:
/// @name Histogram objects
//@{
Profile1DPtr _h_Nchg_towards_pTmumu;
Profile1DPtr _h_Nchg_transverse_pTmumu;
Profile1DPtr _h_Nchg_away_pTmumu;
Profile1DPtr _h_pTsum_towards_pTmumu;
Profile1DPtr _h_pTsum_transverse_pTmumu;
Profile1DPtr _h_pTsum_away_pTmumu;
//@}
};
// Hook for the plugin system
RIVET_DECLARE_PLUGIN(CMS_2017_I1635889);
}
|