CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Type1PFMET.cc
Go to the documentation of this file.
1 
2 // user include files
4 
8 
9 //using namespace std;
10 
11 using namespace reco;
12 
13 // PRODUCER CONSTRUCTORS ------------------------------------------
15 {
16  tokenUncorMet = consumes<METCollection>(iConfig.getParameter<edm::InputTag>("inputUncorMetLabel"));
17  tokenUncorJets = consumes<PFJetCollection>(iConfig.getParameter<edm::InputTag>("inputUncorJetsTag"));
18  correctorToken = consumes<JetCorrector>(iConfig.getParameter<edm::InputTag>("corrector"));
19  jetPTthreshold = iConfig.getParameter<double>("jetPTthreshold");
20  jetEMfracLimit = iConfig.getParameter<double>("jetEMfracLimit");
21  jetMufracLimit = iConfig.getParameter<double>("jetMufracLimit");
22  produces<METCollection>();
23 }
25 
26 // PRODUCER DESTRUCTORS -------------------------------------------
28 
29 // PRODUCER METHODS -----------------------------------------------
31 {
32  using namespace edm;
33  Handle<PFJetCollection> inputUncorJets;
34  iEvent.getByToken( tokenUncorJets, inputUncorJets );
36  iEvent.getByToken( correctorToken, corrector );
37  Handle<METCollection> inputUncorMet; //Define Inputs
38  iEvent.getByToken( tokenUncorMet, inputUncorMet ); //Get Inputs
39  std::auto_ptr<METCollection> output( new METCollection() ); //Create empty output
40  run( *(inputUncorMet.product()), *(corrector.product()), *(inputUncorJets.product()),
41  jetPTthreshold, jetEMfracLimit, jetMufracLimit,
42  &*output ); //Invoke the algorithm
43  iEvent.put( output ); //Put output into Event
44 }
45 
46 void Type1PFMET::run(const METCollection& uncorMET,
48  const PFJetCollection& uncorJet,
49  double jetPTthreshold,
50  double jetEMfracLimit,
51  double jetMufracLimit,
52  METCollection* corMET)
53 {
54  if (!corMET) {
55  std::cerr << "Type1METAlgo_run-> undefined output MET collection. Stop. " << std::endl;
56  return;
57  }
58 
59  double DeltaPx = 0.0;
60  double DeltaPy = 0.0;
61  double DeltaSumET = 0.0;
62  // ---------------- Calculate jet corrections, but only for those uncorrected jets
63  // ---------------- which are above the given threshold. This requires that the
64  // ---------------- uncorrected jets be matched with the corrected jets.
65  for( PFJetCollection::const_iterator jet = uncorJet.begin(); jet != uncorJet.end(); ++jet) {
66  if( jet->pt() > jetPTthreshold ) {
67  double emEFrac =
68  jet->chargedEmEnergyFraction() + jet->neutralEmEnergyFraction();
69  double muEFrac = jet->chargedMuEnergyFraction();
70  if( emEFrac < jetEMfracLimit
71  && muEFrac < jetMufracLimit ) {
72  double corr = corrector.correction (*jet) - 1.; // correction itself
73  DeltaPx += jet->px() * corr;
74  DeltaPy += jet->py() * corr;
75  DeltaSumET += jet->et() * corr;
76  }
77  }
78  }
79  //----------------- Calculate and set deltas for new MET correction
81  delta.mex = - DeltaPx; //correction to MET (from Jets) is negative,
82  delta.mey = - DeltaPy; //since MET points in direction opposite of jets
83  delta.sumet = DeltaSumET;
84  //----------------- Fill holder with corrected MET (= uncorrected + delta) values
85  const MET* u = &(uncorMET.front());
86  double corrMetPx = u->px()+delta.mex;
87  double corrMetPy = u->py()+delta.mey;
88  MET::LorentzVector correctedMET4vector( corrMetPx, corrMetPy, 0.,
89  sqrt (corrMetPx*corrMetPx + corrMetPy*corrMetPy)
90  );
91  //----------------- get previous corrections and push into new corrections
92  std::vector<CorrMETData> corrections = u->mEtCorr();
93  corrections.push_back( delta );
94  //----------------- Push onto MET Collection
95  MET result = MET(u->sumEt()+delta.sumet,
96  corrections,
97  correctedMET4vector,
98  u->vertex() );
99  corMET->push_back(result);
100 
101  return;
102 }
103 
104 // DEFINE_FWK_MODULE(Type1PFMET); //define this as a plug-in
105 
106 
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual const Point & vertex() const
vertex position (overwritten by PF...)
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:23
void run(const reco::METCollection &uncorMET, const reco::JetCorrector &corrector, const reco::PFJetCollection &uncorJet, double jetPTthreshold, double jetEMfracLimit, double jetMufracLimit, reco::METCollection *corMET)
Definition: Type1PFMET.cc:46
int iEvent
Definition: GenABIO.cc:230
double sumEt() const
Definition: MET.h:56
tuple corrector
Definition: mvaPFMET_cff.py:86
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
Definition: MET.h:42
T sqrt(T t)
Definition: SSEVec.h:48
tuple result
Definition: query.py:137
Collection of MET.
double sumet
Definition: CorrMETData.h:20
JetCorrectorParameters corr
Definition: classes.h:5
std::vector< CorrMETData > mEtCorr() const
Definition: MET.h:70
virtual ~Type1PFMET()
Definition: Type1PFMET.cc:27
virtual double px() const
x coordinate of momentum vector
a MET correction term
Definition: CorrMETData.h:14
std::vector< PFJet > PFJetCollection
collection of PFJet objects
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: Type1PFMET.cc:30
double mey
Definition: CorrMETData.h:18
double mex
Definition: CorrMETData.h:17
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
virtual double py() const
y coordinate of momentum vector