CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTMhtProducer.cc
Go to the documentation of this file.
1 
16 //#include "RecoMET/METProducers/interface/METProducer.h"
22 #include <vector>
23 
24 
25 //
26 // constructors and destructor
27 //
29 {
30  inputJetTag_ = iConfig.getParameter< edm::InputTag > ("inputJetTag");
31  minPtJet_= iConfig.getParameter<double> ("minPtJet");
32  etaJet_= iConfig.getParameter<double> ("etaJet");
33  usePt_= iConfig.getParameter<bool>("usePt");
34 
35  m_theJetToken = consumes<edm::View<reco::Jet>>(inputJetTag_);
36  //register your products
37  produces<reco::METCollection>();
38 }
39 
41 
44  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("hltMCJetCorJetIcone5HF07"));
45  desc.add<double>("minPtJet",20.0);
46  desc.add<double>("etaJet",9999.0);
47  desc.add<bool>("usePt",true);
48  descriptions.add("hltMhtProducer",desc);
49 }
50 
51 // ------------ method called to produce the data ------------
52 void
54 {
55  using namespace std;
56  using namespace edm;
57  using namespace reco;
58 
59  auto_ptr<reco::METCollection> result (new reco::METCollection);
60 
61  math::XYZPoint vtx(0,0,0);
62 
63  //Handle<CaloJetCollection> recocalojets;
66 
67  // look at all candidates, check cuts and add to result object
68  double mhtx=0., mhty=0., mht;
69  double jetVar;
70 
71  if(jets->size() > 0){
72  // events with at least one jet
73  //for (CaloJetCollection::const_iterator jet = jets->begin(); jet != jets->end(); jet++) {
74  for(edm::View<reco::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); jet++ ) {
75  jetVar = jet->pt();
76  if (!usePt_) jetVar = jet->et();
77 
78  //---get MHT
79  if (jetVar > minPtJet_ && std::abs(jet->eta()) < etaJet_) {
80  mhtx -= jetVar*cos(jet->phi());
81  mhty -= jetVar*sin(jet->phi());
82  }
83  }
84  mht = sqrt(mhtx*mhtx + mhty*mhty);
85 
86  math::XYZTLorentzVector mhtVec(mhtx,mhty,0,mht);
87  reco::MET mhtobj(mhtVec,vtx);
88  result->push_back( mhtobj );
89 
90  } // events with at least one jet
91 
92 
93  // put object into the Event
94  iEvent.put(result);
95 
96 }
T getParameter(std::string const &) const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:23
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:243
edm::InputTag inputJetTag_
edm::EDGetTokenT< edm::View< reco::Jet > > m_theJetToken
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
Definition: MET.h:32
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
tuple result
Definition: query.py:137
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
HLTMhtProducer(const edm::ParameterSet &)
virtual void produce(edm::Event &, const edm::EventSetup &)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void add(std::string const &label, ParameterSetDescription const &psetDescription)