CMS 3D CMS Logo

RivetAnalysis.h
Go to the documentation of this file.
1 #ifndef GeneratorInterface_RivetInterface_RivetAnalysis_H
2 #define GeneratorInterface_RivetInterface_RivetAnalysis_H
3 
5 
6 #include "Rivet/Analysis.hh"
7 #include "Rivet/Projections/FinalState.hh"
8 #include "Rivet/Particle.hh"
9 #include "Rivet/Particle.fhh"
10 #include "Rivet/Event.hh"
11 #include "Rivet/Projections/FastJets.hh"
12 #include "Rivet/Projections/ChargedLeptons.hh"
13 #include "Rivet/Projections/PromptFinalState.hh"
14 #include "Rivet/Projections/LeptonFinder.hh"
15 #include "Rivet/Projections/VetoedFinalState.hh"
16 #include "Rivet/Projections/IdentifiedFinalState.hh"
17 #include "Rivet/Projections/MissingMomentum.hh"
18 #include "Rivet/Tools/RivetHepMC.hh"
19 
20 namespace Rivet {
21 
22  class RivetAnalysis : public Analysis {
23  public:
24  Particles leptons() const { return _leptons; }
25  Particles photons() const { return _photons; }
26  Particles neutrinos() const { return _neutrinos; }
27  Jets jets() const { return _jets; }
28  Jets fatjets() const { return _fatjets; }
29  Vector3 met() const { return _met; }
30 
31  private:
36 
41 
43 
46  Vector3 _met;
47 
48  public:
50  : Analysis("RivetAnalysis"),
51  _usePromptFinalStates(pset.getParameter<bool>("usePromptFinalStates")),
52  _excludePromptLeptonsFromJetClustering(pset.getParameter<bool>("excludePromptLeptonsFromJetClustering")),
53  _excludeNeutrinosFromJetClustering(pset.getParameter<bool>("excludeNeutrinosFromJetClustering")),
54  _doJetClustering(pset.getParameter<bool>("doJetClustering")),
55 
56  _particleMinPt(pset.getParameter<double>("particleMinPt")),
57  _particleMaxEta(pset.getParameter<double>("particleMaxEta")),
58 
59  _lepConeSize(pset.getParameter<double>("lepConeSize")),
60  _lepMinPt(pset.getParameter<double>("lepMinPt")),
61  _lepMaxEta(pset.getParameter<double>("lepMaxEta")),
62 
63  _jetConeSize(pset.getParameter<double>("jetConeSize")),
64  _jetMinPt(pset.getParameter<double>("jetMinPt")),
65  _jetMaxEta(pset.getParameter<double>("jetMaxEta")),
66 
67  _fatJetConeSize(pset.getParameter<double>("fatJetConeSize")),
68  _fatJetMinPt(pset.getParameter<double>("fatJetMinPt")),
69  _fatJetMaxEta(pset.getParameter<double>("fatJetMaxEta")),
70 
71  _phoMinPt(pset.getParameter<double>("phoMinPt")),
72  _phoMaxEta(pset.getParameter<double>("phoMaxEta")),
73  _phoIsoConeSize(pset.getParameter<double>("phoIsoConeSize")),
74  _phoMaxRelIso(pset.getParameter<double>("phoMaxRelIso")) {}
75 
76  // Initialize Rivet projections
77  void init() override {
78  // Cuts
79  Cut particle_cut = (Cuts::abseta < _particleMaxEta) and (Cuts::pT > _particleMinPt * GeV);
80  Cut lepton_cut = (Cuts::abseta < _lepMaxEta) and (Cuts::pT > _lepMinPt * GeV);
81 
82  // Generic final state
83  FinalState fs(particle_cut);
84 
85  declare(fs, "FS");
86 
87  // Dressed leptons
88  ChargedLeptons charged_leptons(fs);
89  IdentifiedFinalState photons(fs);
90  photons.acceptIdPair(PID::PHOTON);
91 
92  PromptFinalState prompt_leptons(charged_leptons);
93  prompt_leptons.acceptMuonDecays(true);
94  prompt_leptons.acceptTauDecays(true);
95  declare(prompt_leptons, "PromptLeptons");
96 
97  PromptFinalState prompt_photons(photons);
98  prompt_photons.acceptMuonDecays(true);
99  prompt_photons.acceptTauDecays(true);
100 
101  LeptonFinder prompt_dressed_leptons(prompt_leptons, prompt_photons, _lepConeSize, lepton_cut);
102  LeptonFinder dressed_leptons(charged_leptons, photons, _lepConeSize, lepton_cut);
103  declare(_usePromptFinalStates ? prompt_dressed_leptons : dressed_leptons, "DressedLeptons");
104 
105  declare(photons, "Photons");
106 
107  // Jets
108  VetoedFinalState fsForJets(fs);
110  fsForJets.addVetoOnThisFinalState(prompt_dressed_leptons);
111  JetInvisibles invisiblesStrategy = JetInvisibles::DECAY;
113  invisiblesStrategy = JetInvisibles::NONE;
114  declare(FastJets(fsForJets, JetAlg::ANTIKT, _jetConeSize, JetMuons::ALL, invisiblesStrategy), "Jets");
115 
116  // FatJets
117  declare(FastJets(fsForJets, JetAlg::ANTIKT, _fatJetConeSize, JetMuons::ALL, invisiblesStrategy), "FatJets");
118 
119  // Neutrinos
120  IdentifiedFinalState neutrinos(fs);
121  neutrinos.acceptNeutrinos();
122  if (_usePromptFinalStates) {
123  PromptFinalState prompt_neutrinos(neutrinos);
124  prompt_neutrinos.acceptMuonDecays(true);
125  prompt_neutrinos.acceptTauDecays(true);
126  declare(prompt_neutrinos, "Neutrinos");
127  } else
128  declare(neutrinos, "Neutrinos");
129 
130  // MET
131  declare(MissingMomentum(fs), "MET");
132  };
133 
134  // Apply Rivet projections
135  void analyze(const Event& event) override {
136  _jets.clear();
137  _fatjets.clear();
138  _leptons.clear();
139  _photons.clear();
140  _neutrinos.clear();
141 
142  // Get analysis objects from projections
143  Cut jet_cut = (Cuts::abseta < _jetMaxEta) and (Cuts::pT > _jetMinPt * GeV);
144  Cut fatjet_cut = (Cuts::abseta < _fatJetMaxEta) and (Cuts::pT > _fatJetMinPt * GeV);
145 
146  _leptons = apply<LeptonFinder>(event, "DressedLeptons").particlesByPt();
147 
148  // search tau ancestors
149  Particles promptleptons = apply<PromptFinalState>(event, "PromptLeptons").particles();
150  for (auto& lepton : _leptons) {
151  const auto& cl = lepton.constituents().front(); // this has no ancestors anymore :(
152  for (auto const& pl : promptleptons) {
153  if (cl.momentum() == pl.momentum()) {
154  for (auto& p : pl.ancestors()) {
155  if (p.abspid() == 15) {
156  p.setMomentum(p.momentum() * 1e-20);
157  lepton.addConstituent(p, false);
158  }
159  }
160  break;
161  }
162  }
163  }
164 
165  // Photons
166  Particles fsparticles = apply<FinalState>(event, "FS").particles();
167 
168  for (auto& photon : apply<FinalState>(event, "Photons").particlesByPt()) {
169  if (photon.pt() < _phoMinPt)
170  continue;
171 
172  if (abs(photon.eta()) > _phoMaxEta)
173  continue;
174 
175  double photonptsum = 0;
176 
177  for (auto& fsparticle : fsparticles) {
178  if (deltaR(fsparticle, photon) == 0)
179  continue;
180 
181  if (deltaR(fsparticle, photon) > _phoIsoConeSize)
182  continue;
183 
184  photonptsum += fsparticle.pt();
185  }
186 
187  if (photonptsum / photon.pt() > _phoMaxRelIso)
188  continue;
189 
190  _photons.push_back(photon);
191  }
192 
193  if (_doJetClustering) {
194  _jets = apply<FastJets>(event, "Jets").jetsByPt(jet_cut);
195  _fatjets = apply<FastJets>(event, "FatJets").jetsByPt(fatjet_cut);
196  }
197  _neutrinos = apply<FinalState>(event, "Neutrinos").particlesByPt();
198  _met = apply<MissingMomentum>(event, "MET").missingMomentum().p3();
199 
200  // check for leptonic decays of tag particles
201  for (auto& jet : _jets) {
202  for (auto& pt : jet.tags()) {
203  for (auto& p : pt.children(isLepton)) {
204  pt.addConstituent(p, false);
205  }
206  }
207  }
208  for (auto& jet : _fatjets) {
209  for (auto& pt : jet.tags()) {
210  for (auto& p : pt.children(isLepton)) {
211  pt.addConstituent(p, false);
212  }
213  }
214  }
215  };
216 
217  // Do nothing here
218  void finalize() override{};
219 
220  std::string status() const override { return "VALIDATED"; }
221  };
222 
223 } // namespace Rivet
224 
225 #endif
void finalize() override
void init() override
Definition: RivetAnalysis.h:77
Particles leptons() const
Definition: RivetAnalysis.h:24
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:13
def ALL(dt, wheel, station, sector)
Definition: diffTwoXMLs.py:28
Particles photons() const
Definition: RivetAnalysis.h:25
RivetAnalysis(const edm::ParameterSet &pset)
Definition: RivetAnalysis.h:49
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:34
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void analyze(const Event &event) override
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:33
Jets fatjets() const
Definition: RivetAnalysis.h:28
std::string status() const override
Particles neutrinos() const
Definition: RivetAnalysis.h:26
Definition: TkAlStyle.h:43
Vector3 met() const
Definition: RivetAnalysis.h:29
Definition: event.py:1