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  std::vector<DressedLepton> 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:
37 
42 
43  std::vector<DressedLepton> _leptons;
44  ParticleVector _photons, _neutrinos;
45  Jets _jets, _fatjets;
46  Vector3 _met;
47 
48  public:
49  RivetAnalysis(const edm::ParameterSet& pset) : Analysis("RivetAnalysis"),
50  _usePromptFinalStates(pset.getParameter<bool>("usePromptFinalStates")),
51  _excludePromptLeptonsFromJetClustering(pset.getParameter<bool>("excludePromptLeptonsFromJetClustering")),
52  _excludeNeutrinosFromJetClustering(pset.getParameter<bool>("excludeNeutrinosFromJetClustering")),
53 
54  _particleMinPt (pset.getParameter<double>("particleMinPt")),
55  _particleMaxEta (pset.getParameter<double>("particleMaxEta")),
56 
57  _lepConeSize (pset.getParameter<double>("lepConeSize")),
58  _lepMinPt (pset.getParameter<double>("lepMinPt")),
59  _lepMaxEta (pset.getParameter<double>("lepMaxEta")),
60 
61  _jetConeSize (pset.getParameter<double>("jetConeSize")),
62  _jetMinPt (pset.getParameter<double>("jetMinPt")),
63  _jetMaxEta (pset.getParameter<double>("jetMaxEta")),
64 
65  _fatJetConeSize (pset.getParameter<double>("fatJetConeSize")),
66  _fatJetMinPt (pset.getParameter<double>("fatJetMinPt")),
67  _fatJetMaxEta (pset.getParameter<double>("fatJetMaxEta"))
68  {
69  }
70 
71  // Initialize Rivet projections
72  void init() {
73  // Cuts
74  Cut particle_cut = (Cuts::abseta < _particleMaxEta) and (Cuts::pT > _particleMinPt*GeV);
75  Cut lepton_cut = (Cuts::abseta < _lepMaxEta) and (Cuts::pT > _lepMinPt*GeV);
76 
77  // Generic final state
78  FinalState fs(particle_cut);
79 
80  // Dressed leptons
81  ChargedLeptons charged_leptons(fs);
82  IdentifiedFinalState photons(fs);
83  photons.acceptIdPair(PID::PHOTON);
84 
85  PromptFinalState prompt_leptons(charged_leptons);
86  prompt_leptons.acceptMuonDecays(true);
87  prompt_leptons.acceptTauDecays(true);
88 
89  PromptFinalState prompt_photons(photons);
90  prompt_photons.acceptMuonDecays(true);
91  prompt_photons.acceptTauDecays(true);
92 
93  // useDecayPhotons=true allows for photons with tau ancestor,
94  // photons from hadrons are vetoed by the PromptFinalState;
95  // will be default DressedLeptons behaviour for Rivet >= 2.5.4
96  DressedLeptons dressed_leptons(prompt_photons, prompt_leptons, _lepConeSize,
97  lepton_cut, /*cluster*/ true, /*useDecayPhotons*/ true);
98  if (not _usePromptFinalStates)
99  dressed_leptons = DressedLeptons(photons, charged_leptons, _lepConeSize,
100  lepton_cut, /*cluster*/ true, /*useDecayPhotons*/ true);
101  addProjection(dressed_leptons, "DressedLeptons");
102 
103  // Photons
104  if (_usePromptFinalStates) {
105  // We remove the photons used up for lepton dressing in this case
106  VetoedFinalState vetoed_prompt_photons(prompt_photons);
107  vetoed_prompt_photons.addVetoOnThisFinalState(dressed_leptons);
108  addProjection(vetoed_prompt_photons, "Photons");
109  }
110  else
111  addProjection(photons, "Photons");
112 
113  // Jets
114  VetoedFinalState fsForJets(fs);
115  if (_usePromptFinalStates and _excludePromptLeptonsFromJetClustering)
116  fsForJets.addVetoOnThisFinalState(dressed_leptons);
117  JetAlg::InvisiblesStrategy invisiblesStrategy = JetAlg::DECAY_INVISIBLES;
118  if (_excludeNeutrinosFromJetClustering)
119  invisiblesStrategy = JetAlg::NO_INVISIBLES;
120  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _jetConeSize,
121  JetAlg::ALL_MUONS, invisiblesStrategy), "Jets");
122 
123  // FatJets
124  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _fatJetConeSize), "FatJets");
125 
126  // Neutrinos
127  IdentifiedFinalState neutrinos(fs);
128  neutrinos.acceptNeutrinos();
129  if (_usePromptFinalStates) {
130  PromptFinalState prompt_neutrinos(neutrinos);
131  prompt_neutrinos.acceptMuonDecays(true);
132  prompt_neutrinos.acceptTauDecays(true);
133  addProjection(prompt_neutrinos, "Neutrinos");
134  }
135  else
136  addProjection(neutrinos, "Neutrinos");
137 
138  // MET
139  addProjection(MissingMomentum(fs), "MET");
140  };
141 
142  // Apply Rivet projections
143  void analyze(const Event& event) {
144  _jets.clear();
145  _fatjets.clear();
146  _leptons.clear();
147  _photons.clear();
148  _neutrinos.clear();
149 
150  // Get analysis objects from projections
151  Cut jet_cut = (Cuts::abseta < _jetMaxEta) and (Cuts::pT > _jetMinPt*GeV);
152  Cut fatjet_cut = (Cuts::abseta < _fatJetMaxEta) and (Cuts::pT > _fatJetMinPt*GeV);
153 
154  _leptons = applyProjection<DressedLeptons>(event, "DressedLeptons").dressedLeptons();
155  _jets = applyProjection<FastJets>(event, "Jets").jetsByPt(jet_cut);
156  _fatjets = applyProjection<FastJets>(event, "FatJets").jetsByPt(fatjet_cut);
157  _photons = applyProjection<FinalState>(event, "Photons").particlesByPt();
158  _neutrinos = applyProjection<FinalState>(event, "Neutrinos").particlesByPt();
159  _met = applyProjection<MissingMomentum>(event, "MET").missingMomentum().p3();
160  };
161 
162  // Do nothing here
163  void finalize() {};
164 
165  };
166 
167 }
168 
169 #endif
const double GeV
Definition: MathUtil.h:16
Vector3 met() const
Definition: RivetAnalysis.h:31
std::vector< DressedLepton > _leptons
Definition: RivetAnalysis.h:43
Jets fatjets() const
Definition: RivetAnalysis.h:30
ParticleVector _photons
Definition: RivetAnalysis.h:44
ParticleVector neutrinos() const
Definition: RivetAnalysis.h:28
RivetAnalysis(const edm::ParameterSet &pset)
Definition: RivetAnalysis.h:49
Jets jets() const
Definition: RivetAnalysis.h:29
ParticleVector photons() const
Definition: RivetAnalysis.h:27
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:36
std::vector< DressedLepton > leptons() const
Definition: RivetAnalysis.h:26
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
void analyze(const Event &event)
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:35
ParticleVector _neutrinos
Definition: RivetAnalysis.h:44
Definition: event.py:1