CMS 3D CMS Logo

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

#include <MultShiftMETcorrInputProducer.h>

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

Public Member Functions

 MultShiftMETcorrInputProducer (const edm::ParameterSet &)
 
 ~MultShiftMETcorrInputProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

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_
 
std::vector< int > counts_
 
std::vector< double > etaMax_
 
std::vector< double > etaMin_
 
std::vector< std::unique_ptr< TF1 > > formula_x_
 
std::vector< std::unique_ptr< TF1 > > formula_y_
 
std::string moduleLabel_
 
edm::EDGetTokenT< edm::View< reco::Candidate > > pflow_
 
std::vector< double > sumPt_
 
std::vector< int > type_
 
std::vector< int > varType_
 
edm::EDGetTokenT< edm::View< reco::Vertex > > vertices_
 
edm::EDGetTokenT< edm::ValueMap< float > > weightsToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::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
Robert Schoefbeck, Vienna

Definition at line 32 of file MultShiftMETcorrInputProducer.h.

Constructor & Destructor Documentation

◆ MultShiftMETcorrInputProducer()

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

Definition at line 37 of file MultShiftMETcorrInputProducer.cc.

References mps_setup::append, looper::cfg, cfgCorrParameters_, counts_, etaMax_, etaMin_, formula_x_, formula_y_, mps_fire::i, moduleLabel_, pflow_, HLT_2022v12_cff::srcWeights, AlCaHLTBitMon_QueryRunRegistry::string, sumPt_, type_, findQualityFiles::v, varType_, vertices_, and weightsToken_.

38  : moduleLabel_(cfg.getParameter<std::string>("@module_label")) {
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  edm::InputTag srcWeights = cfg.getParameter<edm::InputTag>("srcWeights");
43  if (!srcWeights.label().empty())
45 
46  cfgCorrParameters_ = cfg.getParameter<std::vector<edm::ParameterSet>>("parameters");
47  etaMin_.clear();
48  etaMax_.clear();
49  type_.clear();
50  varType_.clear();
51 
52  produces<CorrMETData>();
53 
54  for (std::vector<edm::ParameterSet>::const_iterator v = cfgCorrParameters_.begin(); v != cfgCorrParameters_.end();
55  v++) {
56  TString corrPxFormula = v->getParameter<std::string>("fx");
57  TString corrPyFormula = v->getParameter<std::string>("fy");
58  std::vector<double> corrPxParams = v->getParameter<std::vector<double>>("px");
59  std::vector<double> corrPyParams = v->getParameter<std::vector<double>>("py");
60 
61  formula_x_.push_back(std::make_unique<TF1>(
62  std::string(moduleLabel_).append("_").append(v->getParameter<std::string>("name")).append("_corrPx").c_str(),
63  v->getParameter<std::string>("fx").c_str()));
64  formula_y_.push_back(std::make_unique<TF1>(
65  std::string(moduleLabel_).append("_").append(v->getParameter<std::string>("name")).append("_corrPy").c_str(),
66  v->getParameter<std::string>("fy").c_str()));
67 
68  for (unsigned i = 0; i < corrPxParams.size(); i++)
69  formula_x_.back()->SetParameter(i, corrPxParams[i]);
70  for (unsigned i = 0; i < corrPyParams.size(); i++)
71  formula_y_.back()->SetParameter(i, corrPyParams[i]);
72 
73  counts_.push_back(0);
74  sumPt_.push_back(0.);
75  etaMin_.push_back(v->getParameter<double>("etaMin"));
76  etaMax_.push_back(v->getParameter<double>("etaMax"));
77  type_.push_back(v->getParameter<int>("type"));
78  varType_.push_back(v->getParameter<int>("varType"));
79  }
80 }
edm::EDGetTokenT< edm::View< reco::Vertex > > vertices_
edm::EDGetTokenT< edm::ValueMap< float > > weightsToken_
std::vector< edm::ParameterSet > cfgCorrParameters_
std::vector< std::unique_ptr< TF1 > > formula_x_
std::vector< std::unique_ptr< TF1 > > formula_y_
edm::EDGetTokenT< edm::View< reco::Candidate > > pflow_

◆ ~MultShiftMETcorrInputProducer()

MultShiftMETcorrInputProducer::~MultShiftMETcorrInputProducer ( )
override

Definition at line 82 of file MultShiftMETcorrInputProducer.cc.

82 {}

Member Function Documentation

◆ produce()

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

Definition at line 84 of file MultShiftMETcorrInputProducer.cc.

References funct::abs(), c, cfgCorrParameters_, counts_, etaMax_, etaMin_, formula_x_, formula_y_, edm::Event::getByToken(), metFilters_cff::goodVertices, mps_fire::i, edm::EDGetTokenT< T >::isUninitialized(), edm::HandleBase::isValid(), dqmiolumiharvest::j, eostools::move(), pfLinker_cff::particleFlow, pflow_, position, edm::Event::put(), rho, sumPt_, translateTypeToAbsPdgId(), type_, findQualityFiles::v, heppy_batch::val, varType_, vertices_, mps_merge::weight, HLT_2022v12_cff::weights, weightsToken_, and z.

84  {
85  //get primary vertices
87  evt.getByToken(vertices_, hpv);
88  if (!hpv.isValid()) {
89  edm::LogError("MultShiftMETcorrInputProducer::produce") << "could not find vertex collection ";
90  }
91  std::vector<reco::Vertex> goodVertices;
92  for (unsigned i = 0; i < hpv->size(); i++) {
93  if ((*hpv)[i].ndof() > 4 && (fabs((*hpv)[i].z()) <= 24.) && (fabs((*hpv)[i].position().rho()) <= 2.0))
94  goodVertices.push_back((*hpv)[i]);
95  }
96  int ngoodVertices = goodVertices.size();
97 
98  for (unsigned i = 0; i < counts_.size(); i++)
99  counts_[i] = 0;
100  for (unsigned i = 0; i < sumPt_.size(); i++)
101  sumPt_[i] = 0.;
102 
105 
109  for (unsigned i = 0; i < particleFlow->size(); ++i) {
110  const reco::Candidate& c = particleFlow->at(i);
111  for (unsigned j = 0; j < type_.size(); j++) {
113  if ((c.eta() > etaMin_[j]) and (c.eta() < etaMax_[j])) {
114  float weight = (!weightsToken_.isUninitialized()) ? (*weights)[particleFlow->ptrAt(i)] : 1.0;
115  counts_[j] += (weight > 0);
116  sumPt_[j] += c.pt() * weight;
117  continue;
118  }
119  }
120  }
121  }
122 
123  //MM: loop over all constituent types and sum each correction
124  std::unique_ptr<CorrMETData> metCorr(new CorrMETData());
125 
126  double corx = 0.;
127  double cory = 0.;
128 
129  for (std::vector<edm::ParameterSet>::const_iterator v = cfgCorrParameters_.begin(); v != cfgCorrParameters_.end();
130  v++) {
131  unsigned j = v - cfgCorrParameters_.begin();
132 
133  double val(0.);
134  if (varType_[j] == 0) {
135  val = counts_[j];
136  }
137  if (varType_[j] == 1) {
138  val = ngoodVertices;
139  }
140  if (varType_[j] == 2) {
141  val = sumPt_[j];
142  }
143 
144  corx -= formula_x_[j]->Eval(val);
145  cory -= formula_y_[j]->Eval(val);
146 
147  } //end loop over corrections
148 
149  metCorr->mex = corx;
150  metCorr->mey = cory;
151  evt.put(std::move(metCorr), "");
152 }
edm::EDGetTokenT< edm::View< reco::Vertex > > vertices_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
ParticleType
particle types
Definition: PFCandidate.h:44
edm::EDGetTokenT< edm::ValueMap< float > > weightsToken_
std::vector< edm::ParameterSet > cfgCorrParameters_
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:99
goodVertices
The Good vertices collection needed by the tracking failure filter ________||.
Definition: weight.py:1
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
Log< level::Error, false > LogError
static int translateTypeToAbsPdgId(reco::PFCandidate::ParticleType type)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:70
a MET correction term
Definition: CorrMETData.h:14
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::vector< std::unique_ptr< TF1 > > formula_x_
std::vector< std::unique_ptr< TF1 > > formula_y_
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< edm::View< reco::Candidate > > pflow_

◆ translateTypeToAbsPdgId()

int MultShiftMETcorrInputProducer::translateTypeToAbsPdgId ( reco::PFCandidate::ParticleType  type)
staticprivate

Definition at line 15 of file MultShiftMETcorrInputProducer.cc.

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

Referenced by produce().

15  {
16  switch (type) {
18  return 211; // pi+
20  return 11;
22  return 13;
24  return 22;
25  case reco::PFCandidate::ParticleType::h0:
26  return 130; // K_L0
27  case reco::PFCandidate::ParticleType::h_HF:
28  return 1; // dummy pdg code
29  case reco::PFCandidate::ParticleType::egamma_HF:
30  return 2; // dummy pdg code
32  default:
33  return 0;
34  }
35 }
#define X(str)
Definition: MuonsGrabber.cc:38
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4

Member Data Documentation

◆ cfgCorrParameters_

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

Definition at line 46 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ counts_

std::vector<int> MultShiftMETcorrInputProducer::counts_
private

Definition at line 49 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ etaMax_

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

Definition at line 48 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ etaMin_

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

Definition at line 48 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ formula_x_

std::vector<std::unique_ptr<TF1> > MultShiftMETcorrInputProducer::formula_x_
private

Definition at line 51 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ formula_y_

std::vector<std::unique_ptr<TF1> > MultShiftMETcorrInputProducer::formula_y_
private

Definition at line 52 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ moduleLabel_

std::string MultShiftMETcorrInputProducer::moduleLabel_
private

◆ pflow_

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

Definition at line 41 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ sumPt_

std::vector<double> MultShiftMETcorrInputProducer::sumPt_
private

Definition at line 50 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ type_

std::vector<int> MultShiftMETcorrInputProducer::type_
private

◆ varType_

std::vector<int> MultShiftMETcorrInputProducer::varType_
private

Definition at line 49 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ vertices_

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

Definition at line 42 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().

◆ weightsToken_

edm::EDGetTokenT<edm::ValueMap<float> > MultShiftMETcorrInputProducer::weightsToken_
private

Definition at line 44 of file MultShiftMETcorrInputProducer.h.

Referenced by MultShiftMETcorrInputProducer(), and produce().