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
 
ParticleVector leptons () const
 
Vector3 met () const
 
ParticleVector neutrinos () const
 
ParticleVector 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
 
ParticleVector _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 50 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"))
74  {
75  }
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, boostedElectronIsolation_cff::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").particlesByPt();
157  // search tau ancestors
158  Particles promptleptons = applyProjection<PromptFinalState>(event, "PromptLeptons").particles();
159  for ( auto & lepton : _leptons ) {
160  const auto & cl = lepton.constituents().front(); // this has no ancestors anymore :(
161  for ( auto const & pl : promptleptons ) {
162  if (cl.momentum() == pl.momentum()) {
163  for ( auto & p : pl.ancestors()) {
164  if (p.abspid() == 15) {
165  p.setMomentum(p.momentum()*10e-20);
166  lepton.addConstituent(p, false);
167  }
168  }
169  break;
170  }
171  }
172  }
173 
174  Particles fsparticles = applyProjection<FinalState>(event,"FS").particles();
175 
176  for ( auto & photon : applyProjection<FinalState>(event, "Photons").particlesByPt() ) {
177 
178  if (photon.pt() < _phoMinPt)
179  continue;
180 
181  if (abs(photon.eta()) > _phoMaxEta)
182  continue;
183 
184  double photonptsum = 0;
185 
186  for (auto &fsparticle : fsparticles) {
187 
188  if (deltaR(fsparticle, photon) == 0)
189  continue;
190 
191  if (deltaR(fsparticle, photon) > _phoIsoConeSize)
192  continue;
193 
194  photonptsum += fsparticle.pt();
195  }
196 
197  if (photonptsum/photon.pt() > _phoMaxRelIso)
198  continue;
199 
200  _photons.push_back(photon);
201 
202  }
203 
204  _jets = applyProjection<FastJets>(event, "Jets").jetsByPt(jet_cut);
205  _fatjets = applyProjection<FastJets>(event, "FatJets").jetsByPt(fatjet_cut);
206  _neutrinos = applyProjection<FinalState>(event, "Neutrinos").particlesByPt();
207  _met = applyProjection<MissingMomentum>(event, "MET").missingMomentum().p3();
208  };
const double GeV
Definition: MathUtil.h:16
ParticleVector _photons
Definition: RivetAnalysis.h:45
ParticleVector _leptons
Definition: RivetAnalysis.h:45
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
ParticleVector _neutrinos
Definition: RivetAnalysis.h:45
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 211 of file RivetAnalysis.h.

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

Definition at line 78 of file RivetAnalysis.h.

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

78  {
79  // Cuts
80  Cut particle_cut = (Cuts::abseta < _particleMaxEta) and (Cuts::pT > _particleMinPt*GeV);
81  Cut lepton_cut = (Cuts::abseta < _lepMaxEta) and (Cuts::pT > _lepMinPt*GeV);
82 
83  // Generic final state
84  FinalState fs(particle_cut);
85 
86  addProjection(fs, "FS");
87 
88  // Dressed leptons
89  ChargedLeptons charged_leptons(fs);
90  IdentifiedFinalState photons(fs);
91  photons.acceptIdPair(PID::PHOTON);
92 
93  PromptFinalState prompt_leptons(charged_leptons);
94  prompt_leptons.acceptMuonDecays(true);
95  prompt_leptons.acceptTauDecays(true);
96  addProjection(prompt_leptons, "PromptLeptons");
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, /*useDecayPhotons*/ true);
107  if (not _usePromptFinalStates)
108  dressed_leptons = DressedLeptons(photons, charged_leptons, _lepConeSize,
109  lepton_cut, /*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;}
ParticleVector Rivet::RivetAnalysis::leptons ( void  ) const
inline

Definition at line 26 of file RivetAnalysis.h.

References _leptons.

Referenced by ParticleLevelProducer::produce().

26 {return _leptons;}
ParticleVector _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:45
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:45
std::string Rivet::RivetAnalysis::status ( void  ) const
inlineoverride

Definition at line 213 of file RivetAnalysis.h.

213 { return "VALIDATED"; }

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 46 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 46 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.

ParticleVector Rivet::RivetAnalysis::_leptons
private

Definition at line 45 of file RivetAnalysis.h.

Referenced by leptons().

Vector3 Rivet::RivetAnalysis::_met
private

Definition at line 47 of file RivetAnalysis.h.

Referenced by met().

ParticleVector Rivet::RivetAnalysis::_neutrinos
private

Definition at line 45 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 45 of file RivetAnalysis.h.

Referenced by photons().

bool Rivet::RivetAnalysis::_usePromptFinalStates
private

Definition at line 34 of file RivetAnalysis.h.