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  public:
25  Particles leptons() const { return _leptons; }
26  Particles photons() const { return _photons; }
27  Particles neutrinos() const { return _neutrinos; }
28  Jets jets() const { return _jets; }
29  Jets fatjets() const { return _fatjets; }
30  Vector3 met() const { return _met; }
31 
32  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 
55  _particleMinPt(pset.getParameter<double>("particleMinPt")),
56  _particleMaxEta(pset.getParameter<double>("particleMaxEta")),
57 
58  _lepConeSize(pset.getParameter<double>("lepConeSize")),
59  _lepMinPt(pset.getParameter<double>("lepMinPt")),
60  _lepMaxEta(pset.getParameter<double>("lepMaxEta")),
61 
62  _jetConeSize(pset.getParameter<double>("jetConeSize")),
63  _jetMinPt(pset.getParameter<double>("jetMinPt")),
64  _jetMaxEta(pset.getParameter<double>("jetMaxEta")),
65 
66  _fatJetConeSize(pset.getParameter<double>("fatJetConeSize")),
67  _fatJetMinPt(pset.getParameter<double>("fatJetMinPt")),
68  _fatJetMaxEta(pset.getParameter<double>("fatJetMaxEta")),
69 
70  _phoMinPt(pset.getParameter<double>("phoMinPt")),
71  _phoMaxEta(pset.getParameter<double>("phoMaxEta")),
72  _phoIsoConeSize(pset.getParameter<double>("phoIsoConeSize")),
73  _phoMaxRelIso(pset.getParameter<double>("phoMaxRelIso")) {}
74 
75  // Initialize Rivet projections
76  void init() override {
77  // Cuts
78  Cut particle_cut = (Cuts::abseta < _particleMaxEta) and (Cuts::pT > _particleMinPt * GeV);
79  Cut lepton_cut = (Cuts::abseta < _lepMaxEta) and (Cuts::pT > _lepMinPt * GeV);
80 
81  // Generic final state
82  FinalState fs(particle_cut);
83 
84  declare(fs, "FS");
85 
86  // Dressed leptons
87  ChargedLeptons charged_leptons(fs);
88  IdentifiedFinalState photons(fs);
89  photons.acceptIdPair(PID::PHOTON);
90 
91  PromptFinalState prompt_leptons(charged_leptons);
92  prompt_leptons.acceptMuonDecays(true);
93  prompt_leptons.acceptTauDecays(true);
94  declare(prompt_leptons, "PromptLeptons");
95 
96  PromptFinalState prompt_photons(photons);
97  prompt_photons.acceptMuonDecays(true);
98  prompt_photons.acceptTauDecays(true);
99 
100  // useDecayPhotons=true allows for photons with tau ancestor,
101  // photons from hadrons are vetoed by the PromptFinalState;
102  // will be default DressedLeptons behaviour for Rivet >= 2.5.4
103  DressedLeptons dressed_leptons(
104  prompt_photons, prompt_leptons, _lepConeSize, lepton_cut, /*useDecayPhotons*/ true);
105  if (not _usePromptFinalStates)
106  dressed_leptons = DressedLeptons(photons, charged_leptons, _lepConeSize, lepton_cut, /*useDecayPhotons*/ true);
107  declare(dressed_leptons, "DressedLeptons");
108 
109  declare(photons, "Photons");
110 
111  // Jets
112  VetoedFinalState fsForJets(fs);
114  fsForJets.addVetoOnThisFinalState(dressed_leptons);
115  JetAlg::Invisibles invisiblesStrategy = JetAlg::Invisibles::DECAY;
117  invisiblesStrategy = JetAlg::Invisibles::NONE;
118  declare(FastJets(fsForJets, FastJets::ANTIKT, _jetConeSize, JetAlg::Muons::ALL, invisiblesStrategy), "Jets");
119 
120  // FatJets
121  declare(FastJets(fsForJets, FastJets::ANTIKT, _fatJetConeSize), "FatJets");
122 
123  // Neutrinos
124  IdentifiedFinalState neutrinos(fs);
125  neutrinos.acceptNeutrinos();
126  if (_usePromptFinalStates) {
127  PromptFinalState prompt_neutrinos(neutrinos);
128  prompt_neutrinos.acceptMuonDecays(true);
129  prompt_neutrinos.acceptTauDecays(true);
130  declare(prompt_neutrinos, "Neutrinos");
131  } else
132  declare(neutrinos, "Neutrinos");
133 
134  // MET
135  declare(MissingMomentum(fs), "MET");
136  };
137 
138  // Apply Rivet projections
139  void analyze(const Event& event) override {
140  _jets.clear();
141  _fatjets.clear();
142  _leptons.clear();
143  _photons.clear();
144  _neutrinos.clear();
145 
146  // Get analysis objects from projections
147  Cut jet_cut = (Cuts::abseta < _jetMaxEta) and (Cuts::pT > _jetMinPt * GeV);
148  Cut fatjet_cut = (Cuts::abseta < _fatJetMaxEta) and (Cuts::pT > _fatJetMinPt * GeV);
149 
150  _leptons = apply<DressedLeptons>(event, "DressedLeptons").particlesByPt();
151 
152  // search tau ancestors
153  Particles promptleptons = apply<PromptFinalState>(event, "PromptLeptons").particles();
154  for (auto& lepton : _leptons) {
155  const auto& cl = lepton.constituents().front(); // this has no ancestors anymore :(
156  for (auto const& pl : promptleptons) {
157  if (cl.momentum() == pl.momentum()) {
158  for (auto& p : pl.ancestors()) {
159  if (p.abspid() == 15) {
160  p.setMomentum(p.momentum() * 10e-20);
161  lepton.addConstituent(p, false);
162  }
163  }
164  break;
165  }
166  }
167  }
168 
169  // Photons
170  Particles fsparticles = apply<FinalState>(event, "FS").particles();
171 
172  for (auto& photon : apply<FinalState>(event, "Photons").particlesByPt()) {
173  if (photon.pt() < _phoMinPt)
174  continue;
175 
176  if (abs(photon.eta()) > _phoMaxEta)
177  continue;
178 
179  double photonptsum = 0;
180 
181  for (auto& fsparticle : fsparticles) {
182  if (deltaR(fsparticle, photon) == 0)
183  continue;
184 
185  if (deltaR(fsparticle, photon) > _phoIsoConeSize)
186  continue;
187 
188  photonptsum += fsparticle.pt();
189  }
190 
191  if (photonptsum / photon.pt() > _phoMaxRelIso)
192  continue;
193 
194  _photons.push_back(photon);
195  }
196 
197  _jets = apply<FastJets>(event, "Jets").jetsByPt(jet_cut);
198  _fatjets = apply<FastJets>(event, "FatJets").jetsByPt(fatjet_cut);
199  _neutrinos = apply<FinalState>(event, "Neutrinos").particlesByPt();
200  _met = apply<MissingMomentum>(event, "MET").missingMomentum().p3();
201  };
202 
203  // Do nothing here
204  void finalize() override{};
205 
206  std::string status() const override { return "VALIDATED"; }
207  };
208 
209 } // namespace Rivet
210 
211 #endif
pat::PHOTON
4:
Definition: ParticleCode.h:22
muons2muons_cfi.photon
photon
Definition: muons2muons_cfi.py:28
electrons_cff.bool
bool
Definition: electrons_cff.py:393
Rivet::RivetAnalysis::_phoMaxEta
double _phoMaxEta
Definition: RivetAnalysis.h:42
Rivet::RivetAnalysis::_lepMaxEta
double _lepMaxEta
Definition: RivetAnalysis.h:38
Rivet::RivetAnalysis::_fatJetConeSize
double _fatJetConeSize
Definition: RivetAnalysis.h:40
Rivet::RivetAnalysis
Definition: RivetAnalysis.h:23
Rivet::RivetAnalysis::_usePromptFinalStates
bool _usePromptFinalStates
Definition: RivetAnalysis.h:33
Rivet::RivetAnalysis::_fatJetMinPt
double _fatJetMinPt
Definition: RivetAnalysis.h:40
Rivet::RivetAnalysis::_neutrinos
Particles _neutrinos
Definition: RivetAnalysis.h:44
Rivet
Definition: RivetAnalysis.h:21
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
Rivet::RivetAnalysis::leptons
Particles leptons() const
Definition: RivetAnalysis.h:25
Rivet::RivetAnalysis::_excludePromptLeptonsFromJetClustering
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:34
Rivet::RivetAnalysis::_lepMinPt
double _lepMinPt
Definition: RivetAnalysis.h:38
Rivet::RivetAnalysis::_phoIsoConeSize
double _phoIsoConeSize
Definition: RivetAnalysis.h:42
Rivet::RivetAnalysis::_jetMaxEta
double _jetMaxEta
Definition: RivetAnalysis.h:39
Rivet::RivetAnalysis::_jetConeSize
double _jetConeSize
Definition: RivetAnalysis.h:39
Rivet::RivetAnalysis::_phoMaxRelIso
double _phoMaxRelIso
Definition: RivetAnalysis.h:42
Rivet::RivetAnalysis::photons
Particles photons() const
Definition: RivetAnalysis.h:26
GetRecoTauVFromDQM_MC_cff.cl
cl
Definition: GetRecoTauVFromDQM_MC_cff.py:38
Rivet::RivetAnalysis::fatjets
Jets fatjets() const
Definition: RivetAnalysis.h:29
Rivet::RivetAnalysis::_excludeNeutrinosFromJetClustering
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:35
Rivet::RivetAnalysis::init
void init() override
Definition: RivetAnalysis.h:76
Rivet::RivetAnalysis::_particleMaxEta
double _particleMaxEta
Definition: RivetAnalysis.h:37
PVValHelper::pT
Definition: PVValidationHelpers.h:70
Event
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Rivet::RivetAnalysis::_fatjets
Jets _fatjets
Definition: RivetAnalysis.h:45
Rivet::RivetAnalysis::neutrinos
Particles neutrinos() const
Definition: RivetAnalysis.h:27
Rivet::RivetAnalysis::jets
Jets jets() const
Definition: RivetAnalysis.h:28
edm::ParameterSet
Definition: ParameterSet.h:47
GeV
const double GeV
Definition: MathUtil.h:16
diffTwoXMLs.ALL
def ALL(dt, wheel, station, sector)
Definition: diffTwoXMLs.py:28
Rivet::RivetAnalysis::status
std::string status() const override
Definition: RivetAnalysis.h:206
Rivet::RivetAnalysis::RivetAnalysis
RivetAnalysis(const edm::ParameterSet &pset)
Definition: RivetAnalysis.h:49
Rivet::RivetAnalysis::finalize
void finalize() override
Definition: RivetAnalysis.h:204
METSkim_cff.Jets
Jets
Definition: METSkim_cff.py:17
Rivet::RivetAnalysis::_photons
Particles _photons
Definition: RivetAnalysis.h:44
Rivet::RivetAnalysis::_lepConeSize
double _lepConeSize
Definition: RivetAnalysis.h:38
Rivet::RivetAnalysis::_met
Vector3 _met
Definition: RivetAnalysis.h:46
Rivet::RivetAnalysis::_phoMinPt
double _phoMinPt
Definition: RivetAnalysis.h:42
Rivet::RivetAnalysis::_fatJetMaxEta
double _fatJetMaxEta
Definition: RivetAnalysis.h:40
Rivet::RivetAnalysis::analyze
void analyze(const Event &event) override
Definition: RivetAnalysis.h:139
NONE
Definition: TkAlStyle.cc:47
Rivet::RivetAnalysis::_leptons
Particles _leptons
Definition: RivetAnalysis.h:44
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
Rivet::RivetAnalysis::met
Vector3 met() const
Definition: RivetAnalysis.h:30
event
Definition: event.py:1
Rivet::RivetAnalysis::_jets
Jets _jets
Definition: RivetAnalysis.h:45
Rivet::RivetAnalysis::_jetMinPt
double _jetMinPt
Definition: RivetAnalysis.h:39
Rivet::RivetAnalysis::_particleMinPt
double _particleMinPt
Definition: RivetAnalysis.h:37
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37