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