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/JetAlg.hh"
13 #include "Rivet/Projections/ChargedLeptons.hh"
14 #include "Rivet/Projections/PromptFinalState.hh"
15 #include "Rivet/Projections/DressedLeptons.hh"
16 #include "Rivet/Projections/VetoedFinalState.hh"
17 #include "Rivet/Projections/IdentifiedFinalState.hh"
18 #include "Rivet/Projections/MissingMomentum.hh"
19 #include "Rivet/Tools/RivetHepMC.hh"
20 
21 namespace Rivet {
22 
23  class RivetAnalysis : public Analysis {
24 
25  public:
26  ParticleVector leptons() const {return _leptons;}
27  ParticleVector photons() const {return _photons;}
28  ParticleVector neutrinos() const {return _neutrinos;}
29  Jets jets() const {return _jets;}
30  Jets fatjets() const {return _fatjets;}
31  Vector3 met() const {return _met;}
32 
33  private:
38 
43 
45 
46  ParticleVector _leptons, _photons, _neutrinos;
47  Jets _jets, _fatjets;
48  Vector3 _met;
49 
50  public:
51  RivetAnalysis(const edm::ParameterSet& pset) : Analysis("RivetAnalysis"),
52  _usePromptFinalStates(pset.getParameter<bool>("usePromptFinalStates")),
53  _excludePromptLeptonsFromJetClustering(pset.getParameter<bool>("excludePromptLeptonsFromJetClustering")),
54  _excludeNeutrinosFromJetClustering(pset.getParameter<bool>("excludeNeutrinosFromJetClustering")),
55  _doJetClustering(pset.getParameter<bool>("doJetClustering")),
56 
57  _particleMinPt (pset.getParameter<double>("particleMinPt")),
58  _particleMaxEta (pset.getParameter<double>("particleMaxEta")),
59 
60  _lepConeSize (pset.getParameter<double>("lepConeSize")),
61  _lepMinPt (pset.getParameter<double>("lepMinPt")),
62  _lepMaxEta (pset.getParameter<double>("lepMaxEta")),
63 
64  _jetConeSize (pset.getParameter<double>("jetConeSize")),
65  _jetMinPt (pset.getParameter<double>("jetMinPt")),
66  _jetMaxEta (pset.getParameter<double>("jetMaxEta")),
67 
68  _fatJetConeSize (pset.getParameter<double>("fatJetConeSize")),
69  _fatJetMinPt (pset.getParameter<double>("fatJetMinPt")),
70  _fatJetMaxEta (pset.getParameter<double>("fatJetMaxEta")),
71 
72  _phoMinPt (pset.getParameter<double>("phoMinPt")),
73  _phoMaxEta (pset.getParameter<double>("phoMaxEta")),
74  _phoIsoConeSize (pset.getParameter<double>("phoIsoConeSize")),
75  _phoMaxRelIso (pset.getParameter<double>("phoMaxRelIso"))
76  {
77  }
78 
79  // Initialize Rivet projections
80  void init() override {
81  // Cuts
82  Cut particle_cut = (Cuts::abseta < _particleMaxEta) and (Cuts::pT > _particleMinPt*GeV);
83  Cut lepton_cut = (Cuts::abseta < _lepMaxEta) and (Cuts::pT > _lepMinPt*GeV);
84 
85  // Generic final state
86  FinalState fs(particle_cut);
87 
88  addProjection(fs, "FS");
89 
90  // Dressed leptons
91  ChargedLeptons charged_leptons(fs);
92  IdentifiedFinalState photons(fs);
93  photons.acceptIdPair(PID::PHOTON);
94 
95  PromptFinalState prompt_leptons(charged_leptons);
96  prompt_leptons.acceptMuonDecays(true);
97  prompt_leptons.acceptTauDecays(true);
98  addProjection(prompt_leptons, "PromptLeptons");
99 
100  PromptFinalState prompt_photons(photons);
101  prompt_photons.acceptMuonDecays(true);
102  prompt_photons.acceptTauDecays(true);
103 
104  // useDecayPhotons=true allows for photons with tau ancestor,
105  // photons from hadrons are vetoed by the PromptFinalState;
106  // will be default DressedLeptons behaviour for Rivet >= 2.5.4
107  DressedLeptons dressed_leptons(prompt_photons, prompt_leptons, _lepConeSize,
108  lepton_cut, /*useDecayPhotons*/ true);
109  if (not _usePromptFinalStates)
110  dressed_leptons = DressedLeptons(photons, charged_leptons, _lepConeSize,
111  lepton_cut, /*useDecayPhotons*/ true);
112  addProjection(dressed_leptons, "DressedLeptons");
113 
114  // Photons
115  addProjection(photons, "Photons");
116 
117  // Jets
118  VetoedFinalState fsForJets(fs);
119  if (_usePromptFinalStates and _excludePromptLeptonsFromJetClustering)
120  fsForJets.addVetoOnThisFinalState(dressed_leptons);
121  JetAlg::InvisiblesStrategy invisiblesStrategy = JetAlg::DECAY_INVISIBLES;
122  if (_excludeNeutrinosFromJetClustering)
123  invisiblesStrategy = JetAlg::NO_INVISIBLES;
124  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _jetConeSize,
125  JetAlg::ALL_MUONS, invisiblesStrategy), "Jets");
126 
127  // FatJets
128  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _fatJetConeSize), "FatJets");
129 
130  // Neutrinos
131  IdentifiedFinalState neutrinos(fs);
132  neutrinos.acceptNeutrinos();
133  if (_usePromptFinalStates) {
134  PromptFinalState prompt_neutrinos(neutrinos);
135  prompt_neutrinos.acceptMuonDecays(true);
136  prompt_neutrinos.acceptTauDecays(true);
137  addProjection(prompt_neutrinos, "Neutrinos");
138  }
139  else
140  addProjection(neutrinos, "Neutrinos");
141 
142  // MET
143  addProjection(MissingMomentum(fs), "MET");
144  };
145 
146  // Apply Rivet projections
147  void analyze(const Event& event) override {
148  _jets.clear();
149  _fatjets.clear();
150  _leptons.clear();
151  _photons.clear();
152  _neutrinos.clear();
153 
154  // Get analysis objects from projections
155  Cut jet_cut = (Cuts::abseta < _jetMaxEta) and (Cuts::pT > _jetMinPt*GeV);
156  Cut fatjet_cut = (Cuts::abseta < _fatJetMaxEta) and (Cuts::pT > _fatJetMinPt*GeV);
157 
158  _leptons = applyProjection<DressedLeptons>(event, "DressedLeptons").particlesByPt();
159  // search tau ancestors
160  Particles promptleptons = applyProjection<PromptFinalState>(event, "PromptLeptons").particles();
161  for ( auto & lepton : _leptons ) {
162  const auto & cl = lepton.constituents().front(); // this has no ancestors anymore :(
163  for ( auto const & pl : promptleptons ) {
164  if (cl.momentum() == pl.momentum()) {
165  for ( auto & p : pl.ancestors()) {
166  if (p.abspid() == 15) {
167  p.setMomentum(p.momentum()*10e-20);
168  lepton.addConstituent(p, false);
169  }
170  }
171  break;
172  }
173  }
174  }
175 
176  Particles fsparticles = applyProjection<FinalState>(event,"FS").particles();
177 
178  for ( auto & photon : applyProjection<FinalState>(event, "Photons").particlesByPt() ) {
179 
180  if (photon.pt() < _phoMinPt)
181  continue;
182 
183  if (abs(photon.eta()) > _phoMaxEta)
184  continue;
185 
186  double photonptsum = 0;
187 
188  for (auto &fsparticle : fsparticles) {
189 
190  if (deltaR(fsparticle, photon) == 0)
191  continue;
192 
193  if (deltaR(fsparticle, photon) > _phoIsoConeSize)
194  continue;
195 
196  photonptsum += fsparticle.pt();
197  }
198 
199  if (photonptsum/photon.pt() > _phoMaxRelIso)
200  continue;
201 
202  _photons.push_back(photon);
203 
204  }
205 
206  if (_doJetClustering) {
207  _jets = apply<FastJets>(event, "Jets").jetsByPt(jet_cut);
208  _fatjets = apply<FastJets>(event, "FatJets").jetsByPt(fatjet_cut);
209  }
210  _neutrinos = applyProjection<FinalState>(event, "Neutrinos").particlesByPt();
211  _met = applyProjection<MissingMomentum>(event, "MET").missingMomentum().p3();
212 
213  // check for leptonic decays of tag particles
214  for (auto& jet : _jets) {
215  for (auto& pt : jet.tags()) {
216  for (auto& p : pt.children(isLepton)) {
217  pt.addConstituent(p, false);
218  }
219  }
220  }
221  for (auto& jet : _fatjets) {
222  for (auto& pt : jet.tags()) {
223  for (auto& p : pt.children(isLepton)) {
224  pt.addConstituent(p, false);
225  }
226  }
227  }
228  };
229 
230  // Do nothing here
231  void finalize() override {};
232 
233  std::string status() const override { return "VALIDATED"; }
234 
235  };
236 
237 }
238 
239 #endif
const double GeV
Definition: MathUtil.h:16
Vector3 met() const
Definition: RivetAnalysis.h:31
void finalize() override
void init() override
Definition: RivetAnalysis.h:80
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:19
Jets fatjets() const
Definition: RivetAnalysis.h:30
ParticleVector _photons
Definition: RivetAnalysis.h:46
ParticleVector neutrinos() const
Definition: RivetAnalysis.h:28
RivetAnalysis(const edm::ParameterSet &pset)
Definition: RivetAnalysis.h:51
Jets jets() const
Definition: RivetAnalysis.h:29
ParticleVector _leptons
Definition: RivetAnalysis.h:46
ParticleVector photons() const
Definition: RivetAnalysis.h:27
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:36
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void analyze(const Event &event) override
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
std::string status() const override
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:35
ParticleVector _neutrinos
Definition: RivetAnalysis.h:46
ParticleVector leptons() const
Definition: RivetAnalysis.h:26
Definition: event.py:1