CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetCorrFactorsProducer.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatAlgos_JetCorrFactorsProducer_h
2 #define PhysicsTools_PatAlgos_JetCorrFactorsProducer_h
3 
39 #include <map>
40 #include <string>
41 #include <boost/shared_ptr.hpp>
42 
45 
50 
55 
56 
57 namespace pat {
58 
60  public:
64  typedef std::map<JetCorrFactors::Flavor, std::vector<std::string> > FlavorCorrLevelMap;
65 
66  public:
68  explicit JetCorrFactorsProducer(const edm::ParameterSet& cfg);
72  virtual void produce(edm::Event& event, const edm::EventSetup& setup);
74  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
75 
76  private:
78  bool flavorDependent() const { return (levels_.size()>1); };
80  std::vector<JetCorrectorParameters> params(const JetCorrectorParametersCollection& parameters, const std::vector<std::string>& levels) const;
84  std::vector<std::string> expand(const std::vector<std::string>& levels, const JetCorrFactors::Flavor& flavor);
86  float evaluate(edm::View<reco::Jet>::const_iterator& jet, boost::shared_ptr<FactorizedJetCorrector>& corrector, int level);
88  int numberOf(const edm::Handle<std::vector<reco::Vertex> >& primaryVertices);
90  std::string payload();
91 
92  private:
94  bool emf_;
98  std::string type_;
100  std::string label_;
102  std::string payload_;
108  bool useNPV_;
109  bool useRho_;
121  };
122 
123  inline int
124  JetCorrFactorsProducer::numberOf(const edm::Handle<std::vector<reco::Vertex> >& primaryVertices)
125  {
126  int npv=0;
127  for(std::vector<reco::Vertex>::const_iterator pv=primaryVertices->begin(); pv!=primaryVertices->end(); ++pv){
128  if(pv->ndof()>=4) ++npv;
129  }
130  return npv;
131  }
132 }
133 
134 #endif
dictionary parameters
Definition: Parameters.py:2
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
std::string payload()
map jet algorithm to payload in DB
std::string label_
label of jec factors
bool emf_
use electromagnetic fraction for jet energy corrections or not (will only have an effect for jets Cal...
std::string payload_
label of payload
JetCorrFactorsProducer(const edm::ParameterSet &cfg)
default constructor
bool flavorDependent() const
return true if the jec levels contain at least one flavor dependent correction level ...
int numberOf(const edm::Handle< std::vector< reco::Vertex > > &primaryVertices)
determines the number of valid primary vertices for the standard L1Offset correction of JetMET ...
virtual void produce(edm::Event &event, const edm::EventSetup &setup)
everything that needs to be done per event
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
Produces a ValueMap between JetCorrFactors and the to the originating reco jets.
edm::ValueMap< pat::JetCorrFactors > JetCorrFactorsMap
value map for JetCorrFactors (to be written into the event)
std::vector< std::string > expand(const std::vector< std::string > &levels, const JetCorrFactors::Flavor &flavor)
std::vector< JetCorrectorParameters > params(const JetCorrectorParametersCollection &parameters, const std::vector< std::string > &levels) const
return the jec parameters as input to the FactorizedJetCorrector for different flavors ...
~JetCorrFactorsProducer()
default destructor
edm::InputTag src_
input jet collection
float evaluate(edm::View< reco::Jet >::const_iterator &jet, boost::shared_ptr< FactorizedJetCorrector > &corrector, int level)
evaluate jet correction factor up to a given level
tuple level
Definition: testEve_cfg.py:34
edm::InputTag primaryVertices_
label for L1Offset primaryVertex collection
edm::InputTag rho_
label for L1FastJet energy density parameter rho
bool useNPV_
use the NPV and rho with the JEC? (used for L1Offset/L1FastJet and L1FastJet, resp.)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
description of configuration file parameters
std::string type_
type of flavor dependent JEC factors (only &#39;J&#39; and &#39;T&#39; are allowed)
std::map< JetCorrFactors::Flavor, std::vector< std::string > > FlavorCorrLevelMap
map of correction levels to different flavors
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")