CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuScleFitMuonProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MuScleFitMuonProducer
4 // Class: MuScleFitMuonProducer
5 //
13 //
14 // Original Author: Marco De Mattia,40 3-B32,+41227671551,
15 // Created: Tue Jun 22 13:50:22 CEST 2010
16 // $Id: MuScleFitMuonProducer.cc,v 1.1 2010/06/22 17:25:39 demattia Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
27 
30 
33 
37 
41 
43 
45  public:
48 
49  private:
50  virtual void beginJob() ;
51  virtual void produce(edm::Event&, const edm::EventSetup&);
52  virtual void endJob() ;
53 
55  // Contains the numbers taken from the J/Psi
56  // smearFunctionType7 smearFunction;
57 };
58 
60  theMuonLabel_( iConfig.getParameter<edm::InputTag>( "MuonLabel" ) )
61 {
62  produces<reco::MuonCollection>();
63 }
64 
65 
67 {
68 }
69 
70 // ------------ method called to produce the data ------------
72 {
74  iEvent.getByLabel (theMuonLabel_, allMuons);
75 
76  std::auto_ptr<reco::MuonCollection> pOut(new reco::MuonCollection);
77 
78  // Apply any bias and/or smearing to the events
79  for( std::vector<reco::Muon>::const_iterator muon = allMuons->begin(); muon != allMuons->end(); ++muon ) {
80 
81  double pt = muon->pt();
82 
83  // Biasing the MC such that it reproduces the data + the residual error
84  double a_0 = 1.0039;
85  // double a_1 = 0.;
86  // // The miscalibration can be applied only for pt < 626 GeV/c
87  // if( pt < -a_0*a_0/(4*a_1) ) {
88  // pt = (-a_0 + sqrt(a_0*a_0 + 4*a_1*pt))/(2*a_1);
89  // }
90  pt = a_0*pt;
91 
92  // Smearing
93  // std::cout << "smearing muon" << std::endl;
94  // std::vector<double> par;
95  // double * y = 0;
96  TF1 G("G", "[0]*exp(-0.5*pow(x,2)/[1])", -5., 5.);
97  double sigma = 0.00014; // converted from TeV^-1 to GeV^-1
98  double norm = 1/(sqrt(2*TMath::Pi()));
99  G.SetParameter(0,norm);
100  G.SetParameter(1,1);
101  // std::cout << "old pt = " << pt;
102  pt = 1/(1/pt + sigma*G.GetRandom());
103  // pt' = pt + sigma pt^2
104  // pt = pt*(1-G.GetRandom());
105  double eta = muon->eta();
106  double phi = muon->phi();
107  // smearFunction.smear(pt, eta, phi, y, par);
108  // std::cout << " new pt = " << pt << std::endl;
109 
110  reco::Muon * newMuon = muon->clone();
111  newMuon->setP4( reco::Particle::PolarLorentzVector( pt, eta, phi, muon->mass() ) );
112 
113  pOut->push_back(*newMuon);
114  }
115 
116  // put into the Event
117  // std::auto_ptr<reco::MuonCollection> pOut(new reco::MuonCollection(*allMuons));
118  iEvent.put(pOut);
119 }
120 
121 // ------------ method called once each job just before starting event loop ------------
122 void
124 {
125 }
126 
127 // ------------ method called once each job just after ending the event loop ------------
128 void
130 {
131 }
132 
133 //define this as a plug-in
const double Pi
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void produce(edm::Event &, const edm::EventSetup &)
virtual void setP4(const LorentzVector &p4)
set 4-momentum
T eta() const
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
T sqrt(T t)
Definition: SSEVec.h:28
tuple norm
Definition: lumiNorm.py:78
tuple allMuons
Definition: allMuons_cfi.py:3
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Particle.h:28
MuScleFitMuonProducer(const edm::ParameterSet &)
Definition: DDAxes.h:10