CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/HLTrigger/JetMET/src/HLTMhtProducer.cc

Go to the documentation of this file.
00001 
00008 #include "HLTrigger/JetMET/interface/HLTMhtProducer.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 #include "DataFormats/Common/interface/View.h"
00011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "DataFormats/Math/interface/deltaPhi.h"
00016 #include "DataFormats/Math/interface/LorentzVector.h"
00017 #include "DataFormats/Math/interface/Point3D.h"
00018 //#include "RecoMET/METProducers/interface/METProducer.h"
00019 #include "DataFormats/METReco/interface/MET.h"
00020 #include "DataFormats/METReco/interface/METFwd.h"
00021 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00022 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00023 #include "FWCore/Utilities/interface/InputTag.h"
00024 #include <vector>
00025 
00026 
00027 //
00028 // constructors and destructor
00029 //
00030 HLTMhtProducer::HLTMhtProducer(const edm::ParameterSet& iConfig)
00031 {
00032   inputJetTag_ = iConfig.getParameter< edm::InputTag > ("inputJetTag");
00033   minPtJet_= iConfig.getParameter<double> ("minPtJet");
00034   etaJet_= iConfig.getParameter<double> ("etaJet");
00035   usePt_= iConfig.getParameter<bool>("usePt");
00036 
00037   //register your products
00038   produces<reco::METCollection>();
00039 }
00040 
00041 HLTMhtProducer::~HLTMhtProducer(){}
00042 
00043 void HLTMhtProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00044   edm::ParameterSetDescription desc;
00045   desc.add<edm::InputTag>("inputJetTag",edm::InputTag("hltMCJetCorJetIcone5HF07"));
00046   desc.add<double>("minPtJet",20.0);
00047   desc.add<double>("etaJet",9999.0);
00048   desc.add<bool>("usePt",true);
00049   descriptions.add("hltMhtProducer",desc);
00050 }
00051 
00052 // ------------ method called to produce the data  ------------
00053 void
00054   HLTMhtProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00055 {
00056   using namespace std;
00057   using namespace edm;
00058   using namespace reco;
00059 
00060   auto_ptr<reco::METCollection> result (new reco::METCollection); 
00061 
00062   math::XYZPoint vtx(0,0,0);
00063   
00064   //Handle<CaloJetCollection> recocalojets;
00065   edm::Handle<edm::View<reco::Jet> > jets;
00066   iEvent.getByLabel(inputJetTag_,jets);
00067 
00068   // look at all candidates,  check cuts and add to result object
00069   double mhtx=0., mhty=0., mht;
00070   double jetVar;
00071   
00072   if(jets->size() > 0){
00073     // events with at least one jet
00074     //for (CaloJetCollection::const_iterator jet = jets->begin(); jet != jets->end(); jet++) {
00075     for(edm::View<reco::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); jet++ ) {
00076       jetVar = jet->pt();
00077       if (!usePt_) jetVar = jet->et();
00078 
00079       //---get MHT
00080       if (jetVar > minPtJet_ && std::abs(jet->eta()) < etaJet_) {
00081         mhtx -= jetVar*cos(jet->phi());
00082         mhty -= jetVar*sin(jet->phi());
00083       }
00084     }
00085     mht = sqrt(mhtx*mhtx + mhty*mhty);
00086 
00087     math::XYZTLorentzVector mhtVec(mhtx,mhty,0,mht);
00088     reco::MET mhtobj(mhtVec,vtx);
00089     result->push_back( mhtobj );
00090     
00091   } // events with at least one jet
00092   
00093     
00094   // put object into the Event
00095   iEvent.put(result);
00096 
00097 }