CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/PhysicsTools/PatExamples/plugins/JetEnergyShift.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/Event.h"
00002 #include "FWCore/Framework/interface/EDProducer.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 
00031 class JetEnergyShift : public edm::EDProducer {
00032 
00033  public:
00035   explicit JetEnergyShift(const edm::ParameterSet&);
00037   ~JetEnergyShift(){};
00038   
00039  private:
00041   virtual void produce(edm::Event&, const edm::EventSetup&);
00042 
00043  private:
00045   edm::InputTag inputJets_;
00047   edm::InputTag inputMETs_;
00049   std::string outputJets_;
00051   std::string outputMETs_;
00053   double scaleFactor_;
00055   double jetPTThresholdForMET_;
00057   double jetEMLimitForMET_;
00058 };
00059 
00060 
00061 #include "DataFormats/PatCandidates/interface/Jet.h"
00062 #include "DataFormats/PatCandidates/interface/MET.h"
00063 
00064 JetEnergyShift::JetEnergyShift(const edm::ParameterSet& cfg):
00065   inputJets_           (cfg.getParameter<edm::InputTag>("inputJets"           )),
00066   inputMETs_           (cfg.getParameter<edm::InputTag>("inputMETs"           )),
00067   scaleFactor_         (cfg.getParameter<double>       ("scaleFactor"         )),
00068   jetPTThresholdForMET_(cfg.getParameter<double>       ("jetPTThresholdForMET")),
00069   jetEMLimitForMET_    (cfg.getParameter<double>       ("jetEMLimitForMET"    ))
00070 {
00071   // use label of input to create label for output
00072   outputJets_ = inputJets_.label();
00073   outputMETs_ = inputMETs_.label();
00074   // register products
00075   produces<std::vector<pat::Jet> >(outputJets_);
00076   produces<std::vector<pat::MET> >(outputMETs_);
00077 }
00078 
00079 void
00080 JetEnergyShift::produce(edm::Event& event, const edm::EventSetup& setup)
00081 {
00082   edm::Handle<std::vector<pat::Jet> > jets;
00083   event.getByLabel(inputJets_, jets);
00084 
00085   edm::Handle<std::vector<pat::MET> > mets;
00086   event.getByLabel(inputMETs_, mets);
00087   
00088   std::auto_ptr<std::vector<pat::Jet> > pJets(new std::vector<pat::Jet>);
00089   std::auto_ptr<std::vector<pat::MET> > pMETs(new std::vector<pat::MET>);
00090 
00091   double dPx    = 0.;
00092   double dPy    = 0.;
00093   double dSumEt = 0.;
00094 
00095   for(std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); ++jet) {
00096     pat::Jet scaledJet = *jet;
00097     scaledJet.scaleEnergy( scaleFactor_ );
00098     pJets->push_back( scaledJet );
00099     // consider jet scale shift only if the raw jet pt and emf 
00100     // is above the thresholds given in the module definition
00101     if(jet->correctedJet("raw").pt() > jetPTThresholdForMET_
00102        && jet->emEnergyFraction() < jetEMLimitForMET_) {
00103       dPx    += scaledJet.px() - jet->px();
00104       dPy    += scaledJet.py() - jet->py();
00105       dSumEt += scaledJet.et() - jet->et();
00106     }
00107   }
00108 
00109   // scale MET accordingly
00110   pat::MET met = *(mets->begin());
00111   double scaledMETPx = met.px() - dPx;
00112   double scaledMETPy = met.py() - dPy;
00113   pat::MET scaledMET(reco::MET(met.sumEt()+dSumEt, reco::MET::LorentzVector(scaledMETPx, scaledMETPy, 0, sqrt(scaledMETPx*scaledMETPx+scaledMETPy*scaledMETPy)), reco::MET::Point(0,0,0)));
00114   pMETs->push_back( scaledMET );
00115   event.put(pJets, outputJets_);
00116   event.put(pMETs, outputMETs_);
00117 }
00118 
00119 #include "FWCore/Framework/interface/MakerMacros.h"
00120 DEFINE_FWK_MODULE( JetEnergyShift );