CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
DistortedMuonProducer Class Reference
Inheritance diagram for DistortedMuonProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 DistortedMuonProducer (const edm::ParameterSet &)
 
 ~DistortedMuonProducer () override
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducer () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void beginJob () override
 
void endJob () override
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

std::vector< double > efficiencyRatioOverMC_
 
std::vector< double > etaBinEdges_
 
edm::EDGetTokenT< reco::GenParticleMatchgenMatchMapToken_
 
edm::EDGetTokenT< edm::View< reco::Muon > > muonToken_
 
std::vector< double > relativeShiftOnPt_
 
std::vector< double > relativeUncertaintyOnPt_
 
std::vector< double > shiftOnOneOverPt_
 
std::vector< double > uncertaintyOnOneOverPt_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 17 of file DistortedMuonProducer.cc.

Constructor & Destructor Documentation

DistortedMuonProducer::DistortedMuonProducer ( const edm::ParameterSet pset)
explicit

Definition at line 49 of file DistortedMuonProducer.cc.

References efficiencyRatioOverMC_, etaBinEdges_, genMatchMapToken_, edm::ParameterSet::getUntrackedParameter(), mps_fire::i, muonToken_, relativeShiftOnPt_, relativeUncertaintyOnPt_, shiftOnOneOverPt_, and uncertaintyOnOneOverPt_.

49  {
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 }
T getUntrackedParameter(std::string const &, T const &) const
std::vector< double > etaBinEdges_
edm::EDGetTokenT< reco::GenParticleMatch > genMatchMapToken_
std::vector< double > efficiencyRatioOverMC_
std::vector< double > shiftOnOneOverPt_
edm::EDGetTokenT< edm::View< reco::Muon > > muonToken_
std::vector< double > relativeUncertaintyOnPt_
std::vector< double > relativeShiftOnPt_
std::vector< double > uncertaintyOnOneOverPt_
DistortedMuonProducer::~DistortedMuonProducer ( )
override

Definition at line 110 of file DistortedMuonProducer.cc.

110  {
111 }

Member Function Documentation

void DistortedMuonProducer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 114 of file DistortedMuonProducer.cc.

114  {
115 }
void DistortedMuonProducer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 118 of file DistortedMuonProducer.cc.

118  {
119 }
void DistortedMuonProducer::produce ( edm::Event ev,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 122 of file DistortedMuonProducer.cc.

References reco::Muon::clone(), DEFINE_FWK_MODULE, efficiencyRatioOverMC_, reco::LeafCandidate::eta(), conversionPostprocessing_cfi::etaBin, etaBinEdges_, genMatchMapToken_, edm::Event::getByToken(), mps_fire::i, edm::Ref< C, T, F >::isNull(), edm::EventBase::isRealData(), LogTrace, reco::LeafCandidate::mass(), eostools::move(), RPCpg::mu, HiRecoMuon_cff::muonCollection, muonToken_, pileupCalc::nbins, MuonAlignmentFromReference_cff::newmuons, reco::LeafCandidate::phi(), reco::LeafCandidate::pt(), edm::Event::put(), relativeShiftOnPt_, relativeUncertaintyOnPt_, reco::LeafCandidate::setP4(), shiftOnOneOverPt_, and uncertaintyOnOneOverPt_.

Referenced by JSONExport.JsonExport::export(), HTMLExport.HTMLExport::export(), and HTMLExport.HTMLExportStatic::export().

122  {
123 
124  if (ev.isRealData()) return;
125 
126  // Muon collection
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::unique_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
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(std::move(newmuons));
220 }
std::vector< double > etaBinEdges_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
double eta() const final
momentum pseudorapidity
edm::EDGetTokenT< reco::GenParticleMatch > genMatchMapToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
double pt() const final
transverse momentum
std::vector< double > efficiencyRatioOverMC_
bool isRealData() const
Definition: EventBase.h:64
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
std::vector< double > shiftOnOneOverPt_
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:250
#define LogTrace(id)
std::vector< double > relativeShiftOnPt_
std::vector< double > uncertaintyOnOneOverPt_
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Particle.h:23
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
Muon * clone() const override
create a clone
def move(src, dest)
Definition: eostools.py:510
double mass() const final
mass

Member Data Documentation

std::vector<double> DistortedMuonProducer::efficiencyRatioOverMC_
private

Definition at line 36 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().

std::vector<double> DistortedMuonProducer::etaBinEdges_
private

Definition at line 29 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().

edm::EDGetTokenT<reco::GenParticleMatch> DistortedMuonProducer::genMatchMapToken_
private

Definition at line 28 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().

edm::EDGetTokenT<edm::View<reco::Muon> > DistortedMuonProducer::muonToken_
private

Definition at line 27 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().

std::vector<double> DistortedMuonProducer::relativeShiftOnPt_
private

Definition at line 32 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().

std::vector<double> DistortedMuonProducer::relativeUncertaintyOnPt_
private

Definition at line 34 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().

std::vector<double> DistortedMuonProducer::shiftOnOneOverPt_
private

Definition at line 31 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().

std::vector<double> DistortedMuonProducer::uncertaintyOnOneOverPt_
private

Definition at line 33 of file DistortedMuonProducer.cc.

Referenced by DistortedMuonProducer(), and produce().