CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DistortedMuonProducer.cc
Go to the documentation of this file.
1 #include <memory>
6 
9 
10 //
11 // class declaration
12 //
14  public:
17 
18  private:
19  virtual void beginJob() ;
20  virtual void produce(edm::Event&, const edm::EventSetup&);
21  virtual void endJob() ;
22 
25  std::vector<double> etaBinEdges_;
26 
27  std::vector<double> shiftOnOneOverPt_; // in [1/GeV]
28  std::vector<double> relativeShiftOnPt_;
29  std::vector<double> uncertaintyOnOneOverPt_; // in [1/GeV]
30  std::vector<double> relativeUncertaintyOnPt_;
31 
32  std::vector<double> efficiencyRatioOverMC_;
33 };
34 
42 
43 #include <CLHEP/Random/RandFlat.h>
44 #include <CLHEP/Random/RandGauss.h>
45 
48 
49  // What is being produced
50  produces<std::vector<reco::Muon> >();
51 
52  // Input products
53  muonTag_ = pset.getUntrackedParameter<edm::InputTag> ("MuonTag", edm::InputTag("muons"));
54  genMatchMapTag_ = pset.getUntrackedParameter<edm::InputTag> ("GenMatchMapTag", edm::InputTag("genMatchMap"));
55 
56  // Eta edges
57  std::vector<double> defEtaEdges;
58  defEtaEdges.push_back(-999999.);
59  defEtaEdges.push_back(999999.);
60  etaBinEdges_ = pset.getUntrackedParameter<std::vector<double> > ("EtaBinEdges",defEtaEdges);
61  unsigned int ninputs_expected = etaBinEdges_.size()-1;
62 
63  // Distortions in muon momentum
64  std::vector<double> defDistortion;
65  defDistortion.push_back(0.);
66 
67  shiftOnOneOverPt_ = pset.getUntrackedParameter<std::vector<double> > ("ShiftOnOneOverPt",defDistortion); // in [1/GeV]
68  if (shiftOnOneOverPt_.size()==1 && ninputs_expected>1) {
69  for (unsigned int i=1; i<ninputs_expected; i++){ shiftOnOneOverPt_.push_back(shiftOnOneOverPt_[0]);}
70  }
71 
72  relativeShiftOnPt_ = pset.getUntrackedParameter<std::vector<double> > ("RelativeShiftOnPt",defDistortion);
73  if (relativeShiftOnPt_.size()==1 && ninputs_expected>1) {
74  for (unsigned int i=1; i<ninputs_expected; i++){ relativeShiftOnPt_.push_back(relativeShiftOnPt_[0]);}
75  }
76 
77  uncertaintyOnOneOverPt_ = pset.getUntrackedParameter<std::vector<double> > ("UncertaintyOnOneOverPt",defDistortion); // in [1/GeV]
78  if (uncertaintyOnOneOverPt_.size()==1 && ninputs_expected>1) {
79  for (unsigned int i=1; i<ninputs_expected; i++){ uncertaintyOnOneOverPt_.push_back(uncertaintyOnOneOverPt_[0]);}
80  }
81 
82  relativeUncertaintyOnPt_ = pset.getUntrackedParameter<std::vector<double> > ("RelativeUncertaintyOnPt",defDistortion);
83  if (relativeUncertaintyOnPt_.size()==1 && ninputs_expected>1) {
84  for (unsigned int i=1; i<ninputs_expected; i++){ relativeUncertaintyOnPt_.push_back(relativeUncertaintyOnPt_[0]);}
85  }
86 
87  // Data/MC efficiency ratios
88  std::vector<double> defEfficiencyRatio;
89  defEfficiencyRatio.push_back(1.);
90  efficiencyRatioOverMC_ = pset.getUntrackedParameter<std::vector<double> > ("EfficiencyRatioOverMC",defEfficiencyRatio);
91  if (efficiencyRatioOverMC_.size()==1 && ninputs_expected>1) {
92  for (unsigned int i=1; i<ninputs_expected; i++){ efficiencyRatioOverMC_.push_back(efficiencyRatioOverMC_[0]);}
93  }
94 
95  // Send a warning if there are inconsistencies in vector sizes !!
96  bool effWrong = efficiencyRatioOverMC_.size()!=ninputs_expected;
97  bool momWrong = shiftOnOneOverPt_.size()!=ninputs_expected
98  || relativeShiftOnPt_.size()!=ninputs_expected
99  || uncertaintyOnOneOverPt_.size()!=ninputs_expected
100  || relativeUncertaintyOnPt_.size()!=ninputs_expected;
101  if ( effWrong and momWrong) {
102  edm::LogError("") << "WARNING: DistortedMuonProducer : Size of some parameters do not match the EtaBinEdges vector!!";
103  }
104 
105 }
106 
109 }
110 
113 }
114 
117 }
118 
121 
122  if (ev.isRealData()) return;
123 
124  // Muon collection
125  edm::Handle<edm::View<reco::Muon> > muonCollection;
126  if (!ev.getByLabel(muonTag_, muonCollection)) {
127  edm::LogError("") << ">>> Muon collection does not exist !!!";
128  return;
129  }
130 
132  if (!ev.getByLabel(genMatchMapTag_, genMatchMap)) {
133  edm::LogError("") << ">>> Muon-GenParticle match map does not exist !!!";
134  return;
135  }
136 
137  unsigned int muonCollectionSize = muonCollection->size();
138 
139  std::auto_ptr<reco::MuonCollection> newmuons (new reco::MuonCollection);
140 
141  for (unsigned int i=0; i<muonCollectionSize; i++) {
142  edm::RefToBase<reco::Muon> mu = muonCollection->refAt(i);
143 
144  double ptgen = mu->pt();
145  double etagen = mu->eta();
146  reco::GenParticleRef gen = (*genMatchMap)[mu];
147  if( !gen.isNull()) {
148  ptgen = gen->pt();
149  etagen = gen->eta();
150  LogTrace("") << ">>> Muon-GenParticle match found; ptmu= " << mu->pt() << ", ptgen= " << ptgen;
151  } else {
152  LogTrace("") << ">>> MUON-GENPARTICLE MATCH NOT FOUND!!!";
153  }
154 
155  // Initialize parameters
156  double effRatio = 0.;
157  double shift1 = 0.;
158  double shift2 = 0.;
159  double sigma1 = 0.;
160  double sigma2 = 0.;
161 
162  // Find out which eta bin should be used
163  unsigned int nbins = etaBinEdges_.size()-1;
164  unsigned int etaBin = nbins;
165  if (etagen>etaBinEdges_[0] && etagen<etaBinEdges_[nbins]) {
166  for (unsigned int j=1; j<=nbins; ++j) {
167  if (etagen>etaBinEdges_[j]) continue;
168  etaBin = j-1;
169  break;
170  }
171  }
172  if (etaBin<nbins) {
173  LogTrace("") << ">>> etaBin: " << etaBin << ", for etagen =" << etagen;
174  } else {
175  // Muon is rejected if outside the considered eta range
176  LogTrace("") << ">>> Muon outside eta range: reject it; etagen = " << etagen;
177  continue;
178  }
179 
180  // Set shifts
181  shift1 = shiftOnOneOverPt_[etaBin];
182  shift2 = relativeShiftOnPt_[etaBin];
183  LogTrace("") << "\tshiftOnOneOverPt= " << shift1*100 << " [%]";
184  LogTrace("") << "\trelativeShiftOnPt= " << shift2*100 << " [%]";
185 
186  // Set resolutions
187  sigma1 = uncertaintyOnOneOverPt_[etaBin];
188  sigma2 = relativeUncertaintyOnPt_[etaBin];
189  LogTrace("") << "\tuncertaintyOnOneOverPt= " << sigma1 << " [1/GeV]";
190  LogTrace("") << "\trelativeUncertaintyOnPt= " << sigma2*100 << " [%]";
191 
192  // Set efficiency ratio
193  effRatio = efficiencyRatioOverMC_[etaBin];
194  LogTrace("") << "\tefficiencyRatioOverMC= " << effRatio;
195 
196  // Reject muons according to efficiency ratio
197  double rndf = CLHEP::RandFlat::shoot();
198  if (rndf>effRatio) continue;
199 
200  // Gaussian Random numbers for smearing
201  double rndg1 = CLHEP::RandGauss::shoot();
202  double rndg2 = CLHEP::RandGauss::shoot();
203 
204  // New muon
205  double ptmu = mu->pt();
206  ptmu += ptgen * ( shift1*ptgen + shift2 + sigma1*rndg1*ptgen + sigma2*rndg2);
207  reco::Muon* newmu = mu->clone();
208  newmu->setP4 (
210  ptmu, mu->eta(), mu->phi(), mu->mass()
211  )
212  );
213  newmuons->push_back(*newmu);
214 
215  }
216 
217  ev.put(newmuons);
218 }
219 
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< double > etaBinEdges_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
DistortedMuonProducer(const edm::ParameterSet &)
virtual void setP4(const LorentzVector &p4)
set 4-momentum
std::vector< double > efficiencyRatioOverMC_
bool isRealData() const
Definition: EventBase.h:60
virtual double eta() const
momentum pseudorapidity
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
virtual double mass() const
mass
bool isNull() const
Checks for null.
Definition: Ref.h:247
std::vector< double > shiftOnOneOverPt_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
int j
Definition: DBlmapReader.cc:9
const int mu
Definition: Constants.h:23
std::vector< double > relativeUncertaintyOnPt_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#define LogTrace(id)
Muon * clone() const
create a clone
Definition: Muon.cc:40
virtual double pt() const
transverse momentum
std::vector< double > relativeShiftOnPt_
std::vector< double > uncertaintyOnOneOverPt_
virtual void produce(edm::Event &, const edm::EventSetup &)
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Particle.h:27
virtual double phi() const
momentum azimuthal angle