Imperial Analysis
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ICSingleMetProducer.cc
Go to the documentation of this file.
2 #include <memory>
3 #include <string>
4 #include <vector>
5 #include "FWCore/Framework/interface/Event.h"
6 #include "FWCore/Framework/interface/EventSetup.h"
7 #include "FWCore/Framework/interface/MakerMacros.h"
8 #include "FWCore/ParameterSet/interface/ParameterSet.h"
9 #include "FWCore/Utilities/interface/InputTag.h"
10 #include "DataFormats/Common/interface/Handle.h"
11 #include "DataFormats/Common/interface/View.h"
12 #include "DataFormats/METReco/interface/MET.h"
18 
19 ICSingleMetProducer::ICSingleMetProducer(const edm::ParameterSet& config)
20  : input_(config.getParameter<edm::InputTag>("input")),
21  branch_(config.getParameter<std::string>("branch")) {
22  consumes<edm::View<reco::MET>>(input_);
23  met_ = new ic::Met();
24  PrintHeaderWithProduces(config, input_, branch_);
25 }
26 
28  delete met_;
29 }
30 
31 void ICSingleMetProducer::produce(edm::Event& event,
32  const edm::EventSetup& setup) {
33  edm::Handle<edm::View<reco::MET> > mets_handle;
34  event.getByLabel(input_, mets_handle);
35 
36  *met_ = ic::Met();
37 
38  if (mets_handle->size() >= 1) {
39  reco::MET const& src = mets_handle->at(0);
40  ic::Met& dest = *met_;
41  dest.set_id(met_hasher_(&src));
42  dest.set_pt(src.pt());
43  dest.set_eta(src.eta());
44  dest.set_phi(src.phi());
45  dest.set_energy(src.energy());
46  dest.set_sum_et(src.sumEt());
47  dest.set_et_sig(src.mEtSig());
48  dest.set_xx_sig(src.getSignificanceMatrix()(0, 0));
49  dest.set_xy_sig(src.getSignificanceMatrix()(0, 1));
50  dest.set_yx_sig(src.getSignificanceMatrix()(1, 0));
51  dest.set_yy_sig(src.getSignificanceMatrix()(1, 1));
52  }
53 }
54 
55 void ICSingleMetProducer::beginJob() {
56  ic::StaticTree::tree_->Branch(branch_.c_str(), &met_);
57 }
58 
59 void ICSingleMetProducer::endJob() {}
60 
61 // define this as a plug-in
Stores a missing transverse energy object and the corresponding significance and corrections.
Definition: Met.hh:15
See documentation here.
void set_xx_sig(double const &xx_sig)
The component of the significance matrix.
Definition: Met.hh:86
void set_yy_sig(double const &yy_sig)
The component of the significance matrix.
Definition: Met.hh:95
void set_et_sig(double const &et_sig)
Signifiance of the missing transverse energy.
Definition: Met.hh:83
void set_phi(double const &phi)
Direct access to .
Definition: Candidate.hh:75
void set_id(std::size_t const &id)
Unique identifier.
Definition: Candidate.hh:66
static TTree * tree_
Definition: StaticTree.hh:13
ICSingleMetProducer(const edm::ParameterSet &)
void set_energy(double const &energy)
Direct access to the energy.
Definition: Candidate.hh:78
void set_yx_sig(double const &yx_sig)
The component of the significance matrix.
Definition: Met.hh:92
Definition: Consumes.h:19
DEFINE_FWK_MODULE(ICSingleMetProducer)
void set_pt(double const &pt)
Direct access to the .
Definition: Candidate.hh:69
void set_eta(double const &eta)
Direct access to .
Definition: Candidate.hh:72
void set_sum_et(double const &sum_et)
Scalar sum of transverse energies for all input objects.
Definition: Met.hh:80
void PrintHeaderWithProduces(edm::ParameterSet const &config, edm::InputTag const &in, std::string branch)
void set_xy_sig(double const &xy_sig)
The component of the significance matrix.
Definition: Met.hh:89