CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
Rivet::RivetAnalysis Class Reference

#include <RivetAnalysis.h>

Inheritance diagram for Rivet::RivetAnalysis:

Public Member Functions

void analyze (const Event &event) override
 
Jets fatjets () const
 
void finalize () override
 
void init () override
 
Jets jets () const
 
Particles leptons () const
 
Vector3 met () const
 
Particles neutrinos () const
 
Particles photons () const
 
 RivetAnalysis (const edm::ParameterSet &pset)
 
std::string status () const override
 

Private Attributes

bool _doJetClustering
 
bool _excludeNeutrinosFromJetClustering
 
bool _excludePromptLeptonsFromJetClustering
 
double _fatJetConeSize
 
double _fatJetMaxEta
 
double _fatJetMinPt
 
Jets _fatjets
 
double _jetConeSize
 
double _jetMaxEta
 
double _jetMinPt
 
Jets _jets
 
double _lepConeSize
 
double _lepMaxEta
 
double _lepMinPt
 
Particles _leptons
 
Vector3 _met
 
Particles _neutrinos
 
double _particleMaxEta
 
double _particleMinPt
 
double _phoIsoConeSize
 
double _phoMaxEta
 
double _phoMaxRelIso
 
double _phoMinPt
 
Particles _photons
 
bool _usePromptFinalStates
 

Detailed Description

Definition at line 22 of file RivetAnalysis.h.

Constructor & Destructor Documentation

◆ RivetAnalysis()

Rivet::RivetAnalysis::RivetAnalysis ( const edm::ParameterSet pset)
inline

Definition at line 49 of file RivetAnalysis.h.

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")) {}
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:34
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:33

Member Function Documentation

◆ analyze()

void Rivet::RivetAnalysis::analyze ( const Event event)
inlineoverride

Definition at line 135 of file RivetAnalysis.h.

References _doJetClustering, _fatJetMinPt, _fatjets, _jetMinPt, _jets, _leptons, _met, _neutrinos, _phoIsoConeSize, _phoMaxEta, _phoMaxRelIso, _phoMinPt, _photons, funct::abs(), haddnano::cl, eleIsoSequence_cff::deltaR, MillePedeFileConverter_cfg::e, edmPickEvents::event, reco::isLepton(), metsig::jet, AlCaHLTBitMon_ParallelJobs::p, primaryVertexAssociation_cfi::particles, displacedMuons_cfi::photon, and DiDispStaMuonMonitor_cfi::pt.

135  {
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  };
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:13
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: event.py:1

◆ fatjets()

Jets Rivet::RivetAnalysis::fatjets ( ) const
inline

Definition at line 28 of file RivetAnalysis.h.

References _fatjets.

Referenced by ParticleLevelProducer::produce().

28 { return _fatjets; }

◆ finalize()

void Rivet::RivetAnalysis::finalize ( void  )
inlineoverride

Definition at line 218 of file RivetAnalysis.h.

218 {}

◆ init()

void Rivet::RivetAnalysis::init ( void  )
inlineoverride

Definition at line 77 of file RivetAnalysis.h.

References _excludeNeutrinosFromJetClustering, _excludePromptLeptonsFromJetClustering, _fatJetConeSize, _jetConeSize, _lepConeSize, _lepMinPt, _particleMinPt, _usePromptFinalStates, diffTwoXMLs::ALL(), make_classfiles::fs, neutrinos(), NONE, pat::PHOTON, and photons().

77  {
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  };
def ALL(dt, wheel, station, sector)
Definition: diffTwoXMLs.py:28
Particles photons() const
Definition: RivetAnalysis.h:25
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:34
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:33
Particles neutrinos() const
Definition: RivetAnalysis.h:26
Definition: TkAlStyle.h:43

◆ jets()

Jets Rivet::RivetAnalysis::jets ( ) const
inline

Definition at line 27 of file RivetAnalysis.h.

References _jets.

Referenced by ParticleLevelProducer::produce().

27 { return _jets; }

◆ leptons()

Particles Rivet::RivetAnalysis::leptons ( void  ) const
inline

Definition at line 24 of file RivetAnalysis.h.

References _leptons.

Referenced by ParticleLevelProducer::produce().

24 { return _leptons; }

◆ met()

Vector3 Rivet::RivetAnalysis::met ( ) const
inline

Definition at line 29 of file RivetAnalysis.h.

References _met.

Referenced by objects.METAnalyzer.METAnalyzer::applyDeltaMet(), and ParticleLevelProducer::produce().

29 { return _met; }

◆ neutrinos()

Particles Rivet::RivetAnalysis::neutrinos ( ) const
inline

Definition at line 26 of file RivetAnalysis.h.

References _neutrinos.

Referenced by init(), and ParticleLevelProducer::produce().

26 { return _neutrinos; }

◆ photons()

Particles Rivet::RivetAnalysis::photons ( ) const
inline

Definition at line 25 of file RivetAnalysis.h.

References _photons.

Referenced by init(), and ParticleLevelProducer::produce().

25 { return _photons; }

◆ status()

std::string Rivet::RivetAnalysis::status ( void  ) const
inlineoverride

Definition at line 220 of file RivetAnalysis.h.

220 { return "VALIDATED"; }

Member Data Documentation

◆ _doJetClustering

bool Rivet::RivetAnalysis::_doJetClustering
private

Definition at line 35 of file RivetAnalysis.h.

Referenced by analyze().

◆ _excludeNeutrinosFromJetClustering

bool Rivet::RivetAnalysis::_excludeNeutrinosFromJetClustering
private

Definition at line 34 of file RivetAnalysis.h.

Referenced by init().

◆ _excludePromptLeptonsFromJetClustering

bool Rivet::RivetAnalysis::_excludePromptLeptonsFromJetClustering
private

Definition at line 33 of file RivetAnalysis.h.

Referenced by init().

◆ _fatJetConeSize

double Rivet::RivetAnalysis::_fatJetConeSize
private

Definition at line 40 of file RivetAnalysis.h.

Referenced by init().

◆ _fatJetMaxEta

double Rivet::RivetAnalysis::_fatJetMaxEta
private

Definition at line 40 of file RivetAnalysis.h.

◆ _fatJetMinPt

double Rivet::RivetAnalysis::_fatJetMinPt
private

Definition at line 40 of file RivetAnalysis.h.

Referenced by analyze().

◆ _fatjets

Jets Rivet::RivetAnalysis::_fatjets
private

Definition at line 45 of file RivetAnalysis.h.

Referenced by analyze(), and fatjets().

◆ _jetConeSize

double Rivet::RivetAnalysis::_jetConeSize
private

Definition at line 39 of file RivetAnalysis.h.

Referenced by init().

◆ _jetMaxEta

double Rivet::RivetAnalysis::_jetMaxEta
private

Definition at line 39 of file RivetAnalysis.h.

◆ _jetMinPt

double Rivet::RivetAnalysis::_jetMinPt
private

Definition at line 39 of file RivetAnalysis.h.

Referenced by analyze().

◆ _jets

Jets Rivet::RivetAnalysis::_jets
private

Definition at line 45 of file RivetAnalysis.h.

Referenced by analyze(), and jets().

◆ _lepConeSize

double Rivet::RivetAnalysis::_lepConeSize
private

Definition at line 38 of file RivetAnalysis.h.

Referenced by init().

◆ _lepMaxEta

double Rivet::RivetAnalysis::_lepMaxEta
private

Definition at line 38 of file RivetAnalysis.h.

◆ _lepMinPt

double Rivet::RivetAnalysis::_lepMinPt
private

Definition at line 38 of file RivetAnalysis.h.

Referenced by init().

◆ _leptons

Particles Rivet::RivetAnalysis::_leptons
private

Definition at line 44 of file RivetAnalysis.h.

Referenced by analyze(), and leptons().

◆ _met

Vector3 Rivet::RivetAnalysis::_met
private

Definition at line 46 of file RivetAnalysis.h.

Referenced by analyze(), and met().

◆ _neutrinos

Particles Rivet::RivetAnalysis::_neutrinos
private

Definition at line 44 of file RivetAnalysis.h.

Referenced by analyze(), and neutrinos().

◆ _particleMaxEta

double Rivet::RivetAnalysis::_particleMaxEta
private

Definition at line 37 of file RivetAnalysis.h.

◆ _particleMinPt

double Rivet::RivetAnalysis::_particleMinPt
private

Definition at line 37 of file RivetAnalysis.h.

Referenced by init().

◆ _phoIsoConeSize

double Rivet::RivetAnalysis::_phoIsoConeSize
private

Definition at line 42 of file RivetAnalysis.h.

Referenced by analyze().

◆ _phoMaxEta

double Rivet::RivetAnalysis::_phoMaxEta
private

Definition at line 42 of file RivetAnalysis.h.

Referenced by analyze().

◆ _phoMaxRelIso

double Rivet::RivetAnalysis::_phoMaxRelIso
private

Definition at line 42 of file RivetAnalysis.h.

Referenced by analyze().

◆ _phoMinPt

double Rivet::RivetAnalysis::_phoMinPt
private

Definition at line 42 of file RivetAnalysis.h.

Referenced by analyze().

◆ _photons

Particles Rivet::RivetAnalysis::_photons
private

Definition at line 44 of file RivetAnalysis.h.

Referenced by analyze(), and photons().

◆ _usePromptFinalStates

bool Rivet::RivetAnalysis::_usePromptFinalStates
private

Definition at line 32 of file RivetAnalysis.h.

Referenced by init().