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 _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 23 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 
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")) {}

Member Function Documentation

◆ analyze()

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

Definition at line 139 of file RivetAnalysis.h.

139  {
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  };

References _fatJetMaxEta, _fatJetMinPt, _fatjets, _jetMaxEta, _jetMinPt, _jets, _leptons, _met, _neutrinos, _phoIsoConeSize, _phoMaxEta, _phoMaxRelIso, _phoMinPt, _photons, funct::abs(), GetRecoTauVFromDQM_MC_cff::cl, PbPb_ZMuSkimMuonDPG_cff::deltaR, MillePedeFileConverter_cfg::e, GeV, AlCaHLTBitMon_ParallelJobs::p, muons2muons_cfi::photon, and PVValHelper::pT.

◆ fatjets()

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

Definition at line 29 of file RivetAnalysis.h.

29 { return _fatjets; }

References _fatjets.

Referenced by ParticleLevelProducer::produce().

◆ finalize()

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

Definition at line 204 of file RivetAnalysis.h.

204 {};

◆ init()

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

Definition at line 76 of file RivetAnalysis.h.

76  {
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  };

References _excludeNeutrinosFromJetClustering, _excludePromptLeptonsFromJetClustering, _fatJetConeSize, _jetConeSize, _lepConeSize, _lepMaxEta, _lepMinPt, _particleMaxEta, _particleMinPt, _usePromptFinalStates, diffTwoXMLs::ALL(), GeV, neutrinos(), NONE, pat::PHOTON, photons(), and PVValHelper::pT.

◆ jets()

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

Definition at line 28 of file RivetAnalysis.h.

28 { return _jets; }

References _jets.

Referenced by ParticleLevelProducer::produce().

◆ leptons()

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

Definition at line 25 of file RivetAnalysis.h.

25 { return _leptons; }

References _leptons.

Referenced by ParticleLevelProducer::produce().

◆ met()

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

Definition at line 30 of file RivetAnalysis.h.

30 { return _met; }

References _met.

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

◆ neutrinos()

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

Definition at line 27 of file RivetAnalysis.h.

27 { return _neutrinos; }

References _neutrinos.

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

◆ photons()

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

Definition at line 26 of file RivetAnalysis.h.

26 { return _photons; }

References _photons.

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

◆ status()

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

Definition at line 206 of file RivetAnalysis.h.

206 { return "VALIDATED"; }

Member Data Documentation

◆ _excludeNeutrinosFromJetClustering

bool Rivet::RivetAnalysis::_excludeNeutrinosFromJetClustering
private

Definition at line 35 of file RivetAnalysis.h.

Referenced by init().

◆ _excludePromptLeptonsFromJetClustering

bool Rivet::RivetAnalysis::_excludePromptLeptonsFromJetClustering
private

Definition at line 34 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.

Referenced by analyze().

◆ _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.

Referenced by analyze().

◆ _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.

Referenced by init().

◆ _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.

Referenced by init().

◆ _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 33 of file RivetAnalysis.h.

Referenced by init().

pat::PHOTON
4:
Definition: ParticleCode.h:22
muons2muons_cfi.photon
photon
Definition: muons2muons_cfi.py:28
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::_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
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
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::_excludeNeutrinosFromJetClustering
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:35
Rivet::RivetAnalysis::_particleMaxEta
double _particleMaxEta
Definition: RivetAnalysis.h:37
PVValHelper::pT
Definition: PVValidationHelpers.h:70
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
Rivet::RivetAnalysis::_fatjets
Jets _fatjets
Definition: RivetAnalysis.h:45
Rivet::RivetAnalysis::neutrinos
Particles neutrinos() const
Definition: RivetAnalysis.h:27
GeV
const double GeV
Definition: MathUtil.h:16
diffTwoXMLs.ALL
def ALL(dt, wheel, station, sector)
Definition: diffTwoXMLs.py:28
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
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
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