CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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)
 
Jets fatjets () const
 
void finalize ()
 
void init ()
 
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
 
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 49 of file RivetAnalysis.h.

49  : Analysis("RivetAnalysis"),
50  _usePromptFinalStates(pset.getParameter<bool>("usePromptFinalStates")),
51  _excludePromptLeptonsFromJetClustering(pset.getParameter<bool>("excludePromptLeptonsFromJetClustering")),
52  _excludeNeutrinosFromJetClustering(pset.getParameter<bool>("excludeNeutrinosFromJetClustering")),
53 
54  _particleMinPt (pset.getParameter<double>("particleMinPt")),
55  _particleMaxEta (pset.getParameter<double>("particleMaxEta")),
56 
57  _lepConeSize (pset.getParameter<double>("lepConeSize")),
58  _lepMinPt (pset.getParameter<double>("lepMinPt")),
59  _lepMaxEta (pset.getParameter<double>("lepMaxEta")),
60 
61  _jetConeSize (pset.getParameter<double>("jetConeSize")),
62  _jetMinPt (pset.getParameter<double>("jetMinPt")),
63  _jetMaxEta (pset.getParameter<double>("jetMaxEta")),
64 
65  _fatJetConeSize (pset.getParameter<double>("fatJetConeSize")),
66  _fatJetMinPt (pset.getParameter<double>("fatJetMinPt")),
67  _fatJetMaxEta (pset.getParameter<double>("fatJetMaxEta"))
68  {
69  }
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)
inline

Definition at line 143 of file RivetAnalysis.h.

References _fatJetMinPt, _fatjets, _jetMinPt, _jets, _leptons, _met, _neutrinos, _photons, event(), and GeV.

143  {
144  _jets.clear();
145  _fatjets.clear();
146  _leptons.clear();
147  _photons.clear();
148  _neutrinos.clear();
149 
150  // Get analysis objects from projections
151  Cut jet_cut = (Cuts::abseta < _jetMaxEta) and (Cuts::pT > _jetMinPt*GeV);
152  Cut fatjet_cut = (Cuts::abseta < _fatJetMaxEta) and (Cuts::pT > _fatJetMinPt*GeV);
153 
154  _leptons = applyProjection<DressedLeptons>(event, "DressedLeptons").dressedLeptons();
155  _jets = applyProjection<FastJets>(event, "Jets").jetsByPt(jet_cut);
156  _fatjets = applyProjection<FastJets>(event, "FatJets").jetsByPt(fatjet_cut);
157  _photons = applyProjection<FinalState>(event, "Photons").particlesByPt();
158  _neutrinos = applyProjection<FinalState>(event, "Neutrinos").particlesByPt();
159  _met = applyProjection<MissingMomentum>(event, "MET").missingMomentum().p3();
160  };
const double GeV
Definition: MathUtil.h:16
std::vector< DressedLepton > _leptons
Definition: RivetAnalysis.h:43
ParticleVector _photons
Definition: RivetAnalysis.h:44
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:44
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  )
inline

Definition at line 163 of file RivetAnalysis.h.

163 {};
void Rivet::RivetAnalysis::init ( void  )
inline

Definition at line 72 of file RivetAnalysis.h.

References _excludeNeutrinosFromJetClustering, _excludePromptLeptonsFromJetClustering, _fatJetConeSize, _jetConeSize, _lepConeSize, _lepMinPt, _particleMinPt, _usePromptFinalStates, GeV, neutrinos(), pat::PHOTON, and photons().

72  {
73  // Cuts
74  Cut particle_cut = (Cuts::abseta < _particleMaxEta) and (Cuts::pT > _particleMinPt*GeV);
75  Cut lepton_cut = (Cuts::abseta < _lepMaxEta) and (Cuts::pT > _lepMinPt*GeV);
76 
77  // Generic final state
78  FinalState fs(particle_cut);
79 
80  // Dressed leptons
81  ChargedLeptons charged_leptons(fs);
82  IdentifiedFinalState photons(fs);
83  photons.acceptIdPair(PID::PHOTON);
84 
85  PromptFinalState prompt_leptons(charged_leptons);
86  prompt_leptons.acceptMuonDecays(true);
87  prompt_leptons.acceptTauDecays(true);
88 
89  PromptFinalState prompt_photons(photons);
90  prompt_photons.acceptMuonDecays(true);
91  prompt_photons.acceptTauDecays(true);
92 
93  // useDecayPhotons=true allows for photons with tau ancestor,
94  // photons from hadrons are vetoed by the PromptFinalState;
95  // will be default DressedLeptons behaviour for Rivet >= 2.5.4
96  DressedLeptons dressed_leptons(prompt_photons, prompt_leptons, _lepConeSize,
97  lepton_cut, /*cluster*/ true, /*useDecayPhotons*/ true);
98  if (not _usePromptFinalStates)
99  dressed_leptons = DressedLeptons(photons, charged_leptons, _lepConeSize,
100  lepton_cut, /*cluster*/ true, /*useDecayPhotons*/ true);
101  addProjection(dressed_leptons, "DressedLeptons");
102 
103  // Photons
104  if (_usePromptFinalStates) {
105  // We remove the photons used up for lepton dressing in this case
106  VetoedFinalState vetoed_prompt_photons(prompt_photons);
107  vetoed_prompt_photons.addVetoOnThisFinalState(dressed_leptons);
108  addProjection(vetoed_prompt_photons, "Photons");
109  }
110  else
111  addProjection(photons, "Photons");
112 
113  // Jets
114  VetoedFinalState fsForJets(fs);
115  if (_usePromptFinalStates and _excludePromptLeptonsFromJetClustering)
116  fsForJets.addVetoOnThisFinalState(dressed_leptons);
117  JetAlg::InvisiblesStrategy invisiblesStrategy = JetAlg::DECAY_INVISIBLES;
119  invisiblesStrategy = JetAlg::NO_INVISIBLES;
120  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _jetConeSize,
121  JetAlg::ALL_MUONS, invisiblesStrategy), "Jets");
122 
123  // FatJets
124  addProjection(FastJets(fsForJets, FastJets::ANTIKT, _fatJetConeSize), "FatJets");
125 
126  // Neutrinos
127  IdentifiedFinalState neutrinos(fs);
128  neutrinos.acceptNeutrinos();
129  if (_usePromptFinalStates) {
130  PromptFinalState prompt_neutrinos(neutrinos);
131  prompt_neutrinos.acceptMuonDecays(true);
132  prompt_neutrinos.acceptTauDecays(true);
133  addProjection(prompt_neutrinos, "Neutrinos");
134  }
135  else
136  addProjection(neutrinos, "Neutrinos");
137 
138  // MET
139  addProjection(MissingMomentum(fs), "MET");
140  };
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:43
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:44
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:44

Member Data Documentation

bool Rivet::RivetAnalysis::_excludeNeutrinosFromJetClustering
private

Definition at line 36 of file RivetAnalysis.h.

Referenced by init().

bool Rivet::RivetAnalysis::_excludePromptLeptonsFromJetClustering
private

Definition at line 35 of file RivetAnalysis.h.

Referenced by init().

double Rivet::RivetAnalysis::_fatJetConeSize
private

Definition at line 41 of file RivetAnalysis.h.

Referenced by init().

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.

Referenced by analyze().

Jets Rivet::RivetAnalysis::_fatjets
private

Definition at line 45 of file RivetAnalysis.h.

Referenced by analyze(), and fatjets().

double Rivet::RivetAnalysis::_jetConeSize
private

Definition at line 40 of file RivetAnalysis.h.

Referenced by init().

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.

Referenced by analyze().

Jets Rivet::RivetAnalysis::_jets
private

Definition at line 45 of file RivetAnalysis.h.

Referenced by analyze(), and jets().

double Rivet::RivetAnalysis::_lepConeSize
private

Definition at line 39 of file RivetAnalysis.h.

Referenced by init().

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.

Referenced by init().

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

Definition at line 43 of file RivetAnalysis.h.

Referenced by analyze(), and leptons().

Vector3 Rivet::RivetAnalysis::_met
private

Definition at line 46 of file RivetAnalysis.h.

Referenced by analyze(), and met().

ParticleVector Rivet::RivetAnalysis::_neutrinos
private

Definition at line 44 of file RivetAnalysis.h.

Referenced by analyze(), and 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.

Referenced by init().

ParticleVector Rivet::RivetAnalysis::_photons
private

Definition at line 44 of file RivetAnalysis.h.

Referenced by analyze(), and photons().

bool Rivet::RivetAnalysis::_usePromptFinalStates
private

Definition at line 34 of file RivetAnalysis.h.

Referenced by init().