CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DistortedMETProducer.cc
Go to the documentation of this file.
1 #include <memory>
6 
9 
10 //
11 // class declaration
12 //
14  public:
15  explicit DistortedMETProducer(const edm::ParameterSet&);
17 
18  private:
19  virtual void beginJob() override ;
20  virtual void produce(edm::Event&, const edm::EventSetup&) override;
21  virtual void endJob() override ;
22 
24  double metScaleShift_; // relative shift (0. => no shift)
25 };
26 
30 
33 
34  // What is being produced
35  produces<std::vector<reco::MET> >();
36 
37  // Input products
38  metToken_ = consumes<edm::View<reco::MET> >(pset.getUntrackedParameter<edm::InputTag> ("MetTag", edm::InputTag("met")));
39  // Distortions in MET in Gev**{-1/2}
40  metScaleShift_ = pset.getUntrackedParameter<double> ("MetScaleShift",1.e-3);
41 
42 }
43 
46 }
47 
50 }
51 
54 
57 
58  if (ev.isRealData()) return;
59 
60  // MET collection
61  edm::Handle<edm::View<reco::MET> > metCollection;
62  if (!ev.getByToken(metToken_, metCollection)) {
63  edm::LogError("") << ">>> MET collection does not exist !!!";
64  return;
65  }
66  edm::RefToBase<reco::MET> met = metCollection->refAt(0);
67 
68  std::auto_ptr<reco::METCollection> newmetCollection (new reco::METCollection);
69 
70  double met_et = met->et() * (1. + metScaleShift_);
71  double sum_et = met->sumEt() * (1. + metScaleShift_);
72  double met_phi = met->phi();
73  double met_ex = met_et*cos(met_phi);
74  double met_ey = met_et*sin(met_phi);
75  reco::Particle::LorentzVector met_p4(met_ex, met_ey, 0., met_et);
76  reco::Particle::Point met_vtx(0.,0.,0.);
77  reco::MET* newmet = new reco::MET(sum_et, met_p4, met_vtx);
78 
79  newmetCollection->push_back(*newmet);
80 
81  ev.put(newmetCollection);
82 }
83 
T getUntrackedParameter(std::string const &, T const &) const
tuple met
____________________________________________________________________________||
Definition: CaloMET_cfi.py:4
virtual double et() const
transverse energy
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual float phi() const
momentum azimuthal angle
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool ev
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:23
bool isRealData() const
Definition: EventBase.h:60
virtual void beginJob() override
double sumEt() const
Definition: MET.h:53
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
Definition: MET.h:39
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
math::XYZPoint Point
point in the space
Definition: Particle.h:31
edm::EDGetTokenT< edm::View< reco::MET > > metToken_
virtual void produce(edm::Event &, const edm::EventSetup &) override
virtual void endJob() override
DistortedMETProducer(const edm::ParameterSet &)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27