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
 
std::vector< DressedLepton > leptons () const
 
Vector3 met () const
 
ParticleVector neutrinos () const
 
ParticleVector photons () const
 
 RivetAnalysis (const edm::ParameterSet &pset)
 

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
 
std::vector< DressedLepton > _leptons
 
Vector3 _met
 
ParticleVector _neutrinos
 
double _particleMaxEta
 
double _particleMinPt
 
double _phoIsoConeSize
 
double _phoMaxEta
 
double _phoMaxRelIso
 
double _phoMinPt
 
ParticleVector _photons
 
bool _usePromptFinalStates
 

Detailed Description

Definition at line 23 of file RivetAnalysis.h.

Constructor & Destructor Documentation

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

Definition at line 51 of file RivetAnalysis.h.

51  : Analysis("RivetAnalysis"),
52  _usePromptFinalStates(pset.getParameter<bool>("usePromptFinalStates")),
53  _excludePromptLeptonsFromJetClustering(pset.getParameter<bool>("excludePromptLeptonsFromJetClustering")),
54  _excludeNeutrinosFromJetClustering(pset.getParameter<bool>("excludeNeutrinosFromJetClustering")),
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"))
75  {
76  }
T getParameter(std::string const &) const
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:36
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:35

Member Function Documentation

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

Definition at line 145 of file RivetAnalysis.h.

References _phoIsoConeSize, _phoMaxRelIso, _phoMinPt, funct::abs(), GetRecoTauVFromDQM_MC_cff::cl, deltaR(), MillePedeFileConverter_cfg::e, event(), GeV, AlCaHLTBitMon_ParallelJobs::p, HadronAndPartonSelector_cfi::particles, and muons2muons_cfi::photon.

145  {
146  _jets.clear();
147  _fatjets.clear();
148  _leptons.clear();
149  _photons.clear();
150  _neutrinos.clear();
151 
152  // Get analysis objects from projections
153  Cut jet_cut = (Cuts::abseta < _jetMaxEta) and (Cuts::pT > _jetMinPt*GeV);
154  Cut fatjet_cut = (Cuts::abseta < _fatJetMaxEta) and (Cuts::pT > _fatJetMinPt*GeV);
155 
156  _leptons = applyProjection<DressedLeptons>(event, "DressedLeptons").dressedLeptons();
157  // search tau ancestors
158  for ( auto & lepton : _leptons ) {
159  const auto & cl = lepton.constituentLepton();
160  for ( auto & p : cl.ancestors()) {
161  if (p.abspid() == 15) {
162  p.setMomentum(p.momentum()*10e-20);
163  lepton.addPhoton(p, false);
164  }
165  }
166  }
167 
168  Particles fsparticles = applyProjection<FinalState>(event,"FS").particles();
169 
170  for ( auto & photon : applyProjection<FinalState>(event, "Photons").particlesByPt() ) {
171 
172  if (photon.pt() < _phoMinPt)
173  continue;
174 
175  if (abs(photon.eta()) > _phoMaxEta)
176  continue;
177 
178  double photonptsum = 0;
179 
180  for (auto &fsparticle : fsparticles) {
181 
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 
198  _jets = applyProjection<FastJets>(event, "Jets").jetsByPt(jet_cut);
199  _fatjets = applyProjection<FastJets>(event, "FatJets").jetsByPt(fatjet_cut);
200  _neutrinos = applyProjection<FinalState>(event, "Neutrinos").particlesByPt();
201  _met = applyProjection<MissingMomentum>(event, "MET").missingMomentum().p3();
202  };
const double GeV
Definition: MathUtil.h:16
std::vector< DressedLepton > _leptons
Definition: RivetAnalysis.h:45
ParticleVector _photons
Definition: RivetAnalysis.h:46
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
ParticleVector _neutrinos
Definition: RivetAnalysis.h:46
Jets Rivet::RivetAnalysis::fatjets ( ) const
inline

Definition at line 30 of file RivetAnalysis.h.

References _fatjets.

Referenced by ParticleLevelProducer::produce().

30 {return _fatjets;}
void Rivet::RivetAnalysis::finalize ( void  )
inlineoverride

Definition at line 205 of file RivetAnalysis.h.

205 {};
void Rivet::RivetAnalysis::init ( void  )
inlineoverride

Definition at line 79 of file RivetAnalysis.h.

References GeV, neutrinos(), pat::PHOTON, and photons().

79  {
80  // Cuts
81  Cut particle_cut = (Cuts::abseta < _particleMaxEta) and (Cuts::pT > _particleMinPt*GeV);
82  Cut lepton_cut = (Cuts::abseta < _lepMaxEta) and (Cuts::pT > _lepMinPt*GeV);
83 
84  // Generic final state
85  FinalState fs(particle_cut);
86 
87  addProjection(fs, "FS");
88 
89  // Dressed leptons
90  ChargedLeptons charged_leptons(fs);
91  IdentifiedFinalState photons(fs);
92  photons.acceptIdPair(PID::PHOTON);
93 
94  PromptFinalState prompt_leptons(charged_leptons);
95  prompt_leptons.acceptMuonDecays(true);
96  prompt_leptons.acceptTauDecays(true);
97 
98  PromptFinalState prompt_photons(photons);
99  prompt_photons.acceptMuonDecays(true);
100  prompt_photons.acceptTauDecays(true);
101 
102  // useDecayPhotons=true allows for photons with tau ancestor,
103  // photons from hadrons are vetoed by the PromptFinalState;
104  // will be default DressedLeptons behaviour for Rivet >= 2.5.4
105  DressedLeptons dressed_leptons(prompt_photons, prompt_leptons, _lepConeSize,
106  lepton_cut, /*cluster*/ true, /*useDecayPhotons*/ true);
107  if (not _usePromptFinalStates)
108  dressed_leptons = DressedLeptons(photons, charged_leptons, _lepConeSize,
109  lepton_cut, /*cluster*/ true, /*useDecayPhotons*/ true);
110  addProjection(dressed_leptons, "DressedLeptons");
111 
112  // Photons
113  addProjection(photons, "Photons");
114 
115  // Jets
116  VetoedFinalState fsForJets(fs);
117  if (_usePromptFinalStates and _excludePromptLeptonsFromJetClustering)
118  fsForJets.addVetoOnThisFinalState(dressed_leptons);
119  JetAlg::InvisiblesStrategy invisiblesStrategy = JetAlg::DECAY_INVISIBLES;
121  invisiblesStrategy = JetAlg::NO_INVISIBLES;
122  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _jetConeSize,
123  JetAlg::ALL_MUONS, invisiblesStrategy), "Jets");
124 
125  // FatJets
126  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _fatJetConeSize), "FatJets");
127 
128  // Neutrinos
129  IdentifiedFinalState neutrinos(fs);
130  neutrinos.acceptNeutrinos();
131  if (_usePromptFinalStates) {
132  PromptFinalState prompt_neutrinos(neutrinos);
133  prompt_neutrinos.acceptMuonDecays(true);
134  prompt_neutrinos.acceptTauDecays(true);
135  addProjection(prompt_neutrinos, "Neutrinos");
136  }
137  else
138  addProjection(neutrinos, "Neutrinos");
139 
140  // MET
141  addProjection(MissingMomentum(fs), "MET");
142  };
const double GeV
Definition: MathUtil.h:16
ParticleVector neutrinos() const
Definition: RivetAnalysis.h:28
ParticleVector photons() const
Definition: RivetAnalysis.h:27
bool _excludeNeutrinosFromJetClustering
Definition: RivetAnalysis.h:36
bool _excludePromptLeptonsFromJetClustering
Definition: RivetAnalysis.h:35
Jets Rivet::RivetAnalysis::jets ( ) const
inline

Definition at line 29 of file RivetAnalysis.h.

References _jets.

Referenced by ParticleLevelProducer::produce().

29 {return _jets;}
std::vector<DressedLepton> Rivet::RivetAnalysis::leptons ( void  ) const
inline

Definition at line 26 of file RivetAnalysis.h.

References _leptons.

Referenced by ParticleLevelProducer::produce().

26 {return _leptons;}
std::vector< DressedLepton > _leptons
Definition: RivetAnalysis.h:45
Vector3 Rivet::RivetAnalysis::met ( ) const
inline

Definition at line 31 of file RivetAnalysis.h.

References _met.

Referenced by ParticleLevelProducer::produce().

31 {return _met;}
ParticleVector Rivet::RivetAnalysis::neutrinos ( ) const
inline

Definition at line 28 of file RivetAnalysis.h.

References _neutrinos.

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

28 {return _neutrinos;}
ParticleVector _neutrinos
Definition: RivetAnalysis.h:46
ParticleVector Rivet::RivetAnalysis::photons ( ) const
inline

Definition at line 27 of file RivetAnalysis.h.

References _photons.

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

27 {return _photons;}
ParticleVector _photons
Definition: RivetAnalysis.h:46

Member Data Documentation

bool Rivet::RivetAnalysis::_excludeNeutrinosFromJetClustering
private

Definition at line 36 of file RivetAnalysis.h.

bool Rivet::RivetAnalysis::_excludePromptLeptonsFromJetClustering
private

Definition at line 35 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_fatJetConeSize
private

Definition at line 41 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_fatJetMaxEta
private

Definition at line 41 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_fatJetMinPt
private

Definition at line 41 of file RivetAnalysis.h.

Jets Rivet::RivetAnalysis::_fatjets
private

Definition at line 47 of file RivetAnalysis.h.

Referenced by fatjets().

double Rivet::RivetAnalysis::_jetConeSize
private

Definition at line 40 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_jetMaxEta
private

Definition at line 40 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_jetMinPt
private

Definition at line 40 of file RivetAnalysis.h.

Jets Rivet::RivetAnalysis::_jets
private

Definition at line 47 of file RivetAnalysis.h.

Referenced by jets().

double Rivet::RivetAnalysis::_lepConeSize
private

Definition at line 39 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_lepMaxEta
private

Definition at line 39 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_lepMinPt
private

Definition at line 39 of file RivetAnalysis.h.

std::vector<DressedLepton> Rivet::RivetAnalysis::_leptons
private

Definition at line 45 of file RivetAnalysis.h.

Referenced by leptons().

Vector3 Rivet::RivetAnalysis::_met
private

Definition at line 48 of file RivetAnalysis.h.

Referenced by met().

ParticleVector Rivet::RivetAnalysis::_neutrinos
private

Definition at line 46 of file RivetAnalysis.h.

Referenced by neutrinos().

double Rivet::RivetAnalysis::_particleMaxEta
private

Definition at line 38 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_particleMinPt
private

Definition at line 38 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_phoIsoConeSize
private

Definition at line 43 of file RivetAnalysis.h.

Referenced by analyze().

double Rivet::RivetAnalysis::_phoMaxEta
private

Definition at line 43 of file RivetAnalysis.h.

double Rivet::RivetAnalysis::_phoMaxRelIso
private

Definition at line 43 of file RivetAnalysis.h.

Referenced by analyze().

double Rivet::RivetAnalysis::_phoMinPt
private

Definition at line 43 of file RivetAnalysis.h.

Referenced by analyze().

ParticleVector Rivet::RivetAnalysis::_photons
private

Definition at line 46 of file RivetAnalysis.h.

Referenced by photons().

bool Rivet::RivetAnalysis::_usePromptFinalStates
private

Definition at line 34 of file RivetAnalysis.h.