CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
MultShiftMETcorrDBInputProducer Class Reference

#include <MultShiftMETcorrDBInputProducer.h>

Inheritance diagram for MultShiftMETcorrDBInputProducer:
edm::stream::EDProducer<>

Public Member Functions

 MultShiftMETcorrDBInputProducer (const edm::ParameterSet &)
 
 ~MultShiftMETcorrDBInputProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Static Private Member Functions

static int translateTypeToAbsPdgId (reco::PFCandidate::ParticleType type)
 

Private Attributes

std::vector< edm::ParameterSetcfgCorrParameters_
 
int counts_
 
std::vector< double > etaMax_
 
std::vector< double > etaMin_
 
std::unique_ptr< TF1 > formula_x_
 
std::unique_ptr< TF1 > formula_y_
 
bool mIsData
 
std::string moduleLabel_
 
std::string mPayloadName
 
std::string mSampleType
 
edm::EDGetTokenT< edm::View< reco::Candidate > > pflow_
 
double sumPt_
 
edm::EDGetTokenT< edm::View< reco::Vertex > > vertices_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Compute MET correction to compensate systematic shift of MET in x/y-direction (cf. https://indico.cern.ch/getFile.py/access?contribId=1&resId=0&materialId=slides&confId=174318 )

Authors
SangEun Lee,
Robert Schoefbeck, Vienna

Definition at line 32 of file MultShiftMETcorrDBInputProducer.h.

Constructor & Destructor Documentation

MultShiftMETcorrDBInputProducer::MultShiftMETcorrDBInputProducer ( const edm::ParameterSet cfg)
explicit

Definition at line 31 of file MultShiftMETcorrDBInputProducer.cc.

References etaMax_, etaMin_, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), mIsData, mPayloadName, mSampleType, pflow_, AlCaHLTBitMon_QueryRunRegistry::string, and vertices_.

31  :
32  moduleLabel_(cfg.getParameter<std::string>("@module_label"))
33 {
34 
35  mPayloadName = cfg.getParameter<std::string>("payloadName");
36  mSampleType = (cfg.exists("sampleType")) ? cfg.getParameter< std::string >("sampleType"): "MC";
37  mIsData = cfg.getParameter< bool >("isData");
38 
39  pflow_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter< edm::InputTag >("srcPFlow") );
40  vertices_ = consumes<edm::View<reco::Vertex> >( cfg.getParameter< edm::InputTag >("vertexCollection") );
41 
42  etaMin_.clear();
43  etaMax_.clear();
44 
45  produces<CorrMETData>();
46 
47 
48 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::View< reco::Candidate > > pflow_
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< edm::View< reco::Vertex > > vertices_
MultShiftMETcorrDBInputProducer::~MultShiftMETcorrDBInputProducer ( )
override

Definition at line 50 of file MultShiftMETcorrDBInputProducer.cc.

51 {
52 }

Member Function Documentation

void MultShiftMETcorrDBInputProducer::produce ( edm::Event evt,
const edm::EventSetup es 
)
overrideprivate

Definition at line 54 of file MultShiftMETcorrDBInputProducer.cc.

References funct::abs(), EnergyCorrector::c, counts_, DEFINE_FWK_MODULE, MEtXYcorrectParameters::definitions(), reco::Candidate::eta(), Exception, MEtXYcorrectParametersCollection::findLabel(), MEtXYcorrectParameters::Definitions::formula(), formula_x_, formula_y_, edm::EventSetup::get(), edm::Event::getByToken(), metFilters_cff::goodVertices, mps_fire::i, MEtXYcorrectParametersCollection::isShiftData(), MEtXYcorrectParametersCollection::isShiftDY(), MEtXYcorrectParametersCollection::isShiftMC(), MEtXYcorrectParametersCollection::isShiftTTJets(), MEtXYcorrectParametersCollection::isShiftWJets(), edm::HandleBase::isValid(), relativeConstraints::keys, mIsData, eostools::move(), mPayloadName, mSampleType, MEtXYcorrectParameters::Record::nParameters(), MEtXYcorrectParameters::Record::parameter(), pfLinker_cff::particleFlow, MEtXYcorrectParameters::Definitions::parVar(), reco::Candidate::pdgId(), pflow_, position, reco::Candidate::pt(), MEtXYcorrectParameters::Definitions::PtclType(), edm::Event::put(), MEtXYcorrectParameters::record(), rho, AlCaHLTBitMon_QueryRunRegistry::string, sumPt_, translateTypeToAbsPdgId(), heppy_batch::val, MEtXYcorrectParametersCollection::validKeys(), vertices_, MEtXYcorrectParameters::Record::xMax(), MEtXYcorrectParameters::Record::xMin(), and z.

55 {
56  // Get para.s from DB
58  es.get<MEtXYcorrectRecord>().get(mPayloadName,MEtXYcorParaColl);
59 
60  // get the sections from Collection (pair of section and METCorr.Par class)
61  std::vector<MEtXYcorrectParametersCollection::key_type> keys;
62  // save level to keys for each METParameter in METParameter collection
63  MEtXYcorParaColl->validKeys( keys );
64 
65 
66  //get primary vertices
68  evt.getByToken( vertices_, hpv );
69  if(!hpv.isValid()) {
70  edm::LogError("MultShiftMETcorrDBInputProducer::produce") << "could not find vertex collection ";
71  }
72  std::vector<reco::Vertex> goodVertices;
73  for (unsigned i = 0; i < hpv->size(); i++) {
74  if ( (*hpv)[i].ndof() > 4 &&
75  ( fabs((*hpv)[i].z()) <= 24. ) &&
76  ( fabs((*hpv)[i].position().rho()) <= 2.0 ) )
77  goodVertices.push_back((*hpv)[i]);
78  }
79  int ngoodVertices = goodVertices.size();
80 
82  evt.getByToken(pflow_, particleFlow);
83 
84 
85  //loop over all constituent types and sum each correction
86  //std::auto_ptr<CorrMETData> metCorr(new CorrMETData());
87  std::unique_ptr<CorrMETData> metCorr(new CorrMETData());
88 
89  double corx=0.;
90  double cory=0.;
91 
92 
93  // check DB
94  for ( std::vector<MEtXYcorrectParametersCollection::key_type>::const_iterator ikey = keys.begin();
95  ikey != keys.end(); ++ikey ) {
96  if(mIsData)
97  {
98  if(!MEtXYcorParaColl->isShiftData(*ikey) )
99  throw cms::Exception("MultShiftMETcorrDBInputProducer::produce")
100  << "DB is not for Data. Set proper option: \"corrPfMetXYMultDB.isData\" !!\n";
101  }else{
102  if( MEtXYcorParaColl->isShiftData(*ikey) )
103  throw cms::Exception("MultShiftMETcorrDBInputProducer::produce")
104  << "DB is for Data. Set proper option: \"corrPfMetXYMultDB.isData\" !!\n";
105  }
106  }
107 
108  for ( std::vector<MEtXYcorrectParametersCollection::key_type>::const_iterator ikey = keys.begin();
109  ikey != keys.end(); ++ikey ) {
110 
111  if( !mIsData){
112 
113  if(mSampleType == "MC"){
114  if(!MEtXYcorParaColl->isShiftMC(*ikey)) continue;
115  }else if(mSampleType == "DY"){
116  if(!MEtXYcorParaColl->isShiftDY(*ikey)) continue;
117  }else if(mSampleType == "TTJets"){
118  if(!MEtXYcorParaColl->isShiftTTJets(*ikey)) continue;
119  }else if(mSampleType == "WJets"){
120  if(!MEtXYcorParaColl->isShiftWJets(*ikey)) continue;
121  }else throw cms::Exception("MultShiftMETcorrDBInputProducer::produce")
122  << "SampleType: "<<mSampleType<<" is not reserved !!!\n";
123  }
124 
125  std::string sectionName= MEtXYcorParaColl->findLabel(*ikey);
126  MEtXYcorrectParameters const & MEtXYcorParams = (*MEtXYcorParaColl)[*ikey];
127 
128  counts_ = 0;
129  sumPt_ = 0;
130 
131  for (unsigned i = 0; i < particleFlow->size(); ++i) {
132  const reco::Candidate& c = particleFlow->at(i);
134  if ((c.eta()>MEtXYcorParams.record(0).xMin(0)) and (c.eta()<MEtXYcorParams.record(0).xMax(0))) {
135  counts_ +=1;
136  sumPt_ +=c.pt();
137  continue;
138  }
139  }
140  }
141  double val(0.);
142  unsigned parVar = MEtXYcorParams.definitions().parVar(0);
143 
144  if ( parVar ==0) {
145  val = counts_;
146 
147  }else if ( parVar ==1) {
148  val = ngoodVertices;
149 
150  }else if ( parVar ==2) {
151  val = sumPt_;
152 
153  }else{
154  throw cms::Exception("MultShiftMETcorrDBInputProducer::produce")
155  << "parVar: "<<parVar<<" is not reserved !!!\n";
156  }
157 
158  formula_x_.reset( new TF1("corrPx", MEtXYcorParams.definitions().formula().c_str()));
159  formula_y_.reset( new TF1("corrPy", MEtXYcorParams.definitions().formula().c_str()));
160 
161  for( unsigned i(0); i<MEtXYcorParams.record(0).nParameters(); i++)
162  {
163  formula_x_->SetParameter(i,MEtXYcorParams.record(0).parameter(i));
164  }
165  for( unsigned i(0); i<MEtXYcorParams.record(1).nParameters(); i++)
166  {
167  formula_y_->SetParameter(i,MEtXYcorParams.record(1).parameter(i));
168  }
169 
170  corx -= formula_x_->Eval(val);
171  cory -= formula_y_->Eval(val);
172 
173  } //end loop over corrections
174 
175 
176  metCorr->mex = corx;
177  metCorr->mey = cory;
178  evt.put(std::move(metCorr), "");
179 
180 }
float parameter(unsigned fIndex) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
ParticleType
particle types
Definition: PFCandidate.h:45
void validKeys(std::vector< key_type > &keys) const
edm::EDGetTokenT< edm::View< reco::Candidate > > pflow_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
const Definitions & definitions() const
float xMax(unsigned fVar) const
goodVertices
The Good vertices collection needed by the tracking failure filter ________||.
const Record & record(unsigned fBin) const
static int translateTypeToAbsPdgId(reco::PFCandidate::ParticleType type)
virtual int pdgId() const =0
PDG identifier.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:74
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
static std::string findLabel(key_type k)
a MET correction term
Definition: CorrMETData.h:14
float xMin(unsigned fVar) const
static int position[264][3]
Definition: ReadPGInfo.cc:509
T get() const
Definition: EventSetup.h:63
edm::EDGetTokenT< edm::View< reco::Vertex > > vertices_
std::vector< unsigned > parVar() const
def move(src, dest)
Definition: eostools.py:510
int MultShiftMETcorrDBInputProducer::translateTypeToAbsPdgId ( reco::PFCandidate::ParticleType  type)
staticprivate

Definition at line 16 of file MultShiftMETcorrDBInputProducer.cc.

References MillePedeFileConverter_cfg::e, CustomPhysics_cfi::gamma, h, RPCpg::mu, and X.

Referenced by produce().

16  {
17  switch( type ) {
18  case reco::PFCandidate::ParticleType::h: return 211; // pi+
19  case reco::PFCandidate::ParticleType::e: return 11;
22  case reco::PFCandidate::ParticleType::h0: return 130; // K_L0
23  case reco::PFCandidate::ParticleType::h_HF: return 1; // dummy pdg code
24  case reco::PFCandidate::ParticleType::egamma_HF: return 2; // dummy pdg code
26  default: return 0;
27  }
28 }
type
Definition: HCALResponse.h:21
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
#define X(str)
Definition: MuonsGrabber.cc:48
const int mu
Definition: Constants.h:22

Member Data Documentation

std::vector<edm::ParameterSet> MultShiftMETcorrDBInputProducer::cfgCorrParameters_
private

Definition at line 52 of file MultShiftMETcorrDBInputProducer.h.

int MultShiftMETcorrDBInputProducer::counts_
private

Definition at line 55 of file MultShiftMETcorrDBInputProducer.h.

Referenced by produce().

std::vector<double> MultShiftMETcorrDBInputProducer::etaMax_
private

Definition at line 54 of file MultShiftMETcorrDBInputProducer.h.

Referenced by MultShiftMETcorrDBInputProducer().

std::vector<double> MultShiftMETcorrDBInputProducer::etaMin_
private

Definition at line 54 of file MultShiftMETcorrDBInputProducer.h.

Referenced by MultShiftMETcorrDBInputProducer().

std::unique_ptr< TF1 > MultShiftMETcorrDBInputProducer::formula_x_
private

Definition at line 57 of file MultShiftMETcorrDBInputProducer.h.

Referenced by produce().

std::unique_ptr< TF1 > MultShiftMETcorrDBInputProducer::formula_y_
private

Definition at line 58 of file MultShiftMETcorrDBInputProducer.h.

Referenced by produce().

bool MultShiftMETcorrDBInputProducer::mIsData
private

Definition at line 50 of file MultShiftMETcorrDBInputProducer.h.

Referenced by MultShiftMETcorrDBInputProducer(), and produce().

std::string MultShiftMETcorrDBInputProducer::moduleLabel_
private
std::string MultShiftMETcorrDBInputProducer::mPayloadName
private

Definition at line 48 of file MultShiftMETcorrDBInputProducer.h.

Referenced by MultShiftMETcorrDBInputProducer(), and produce().

std::string MultShiftMETcorrDBInputProducer::mSampleType
private

Definition at line 49 of file MultShiftMETcorrDBInputProducer.h.

Referenced by MultShiftMETcorrDBInputProducer(), and produce().

edm::EDGetTokenT<edm::View<reco::Candidate> > MultShiftMETcorrDBInputProducer::pflow_
private

Definition at line 45 of file MultShiftMETcorrDBInputProducer.h.

Referenced by MultShiftMETcorrDBInputProducer(), and produce().

double MultShiftMETcorrDBInputProducer::sumPt_
private

Definition at line 56 of file MultShiftMETcorrDBInputProducer.h.

Referenced by produce().

edm::EDGetTokenT<edm::View<reco::Vertex> > MultShiftMETcorrDBInputProducer::vertices_
private

Definition at line 46 of file MultShiftMETcorrDBInputProducer.h.

Referenced by MultShiftMETcorrDBInputProducer(), and produce().