CMS 3D CMS Logo

ShiftedParticleProducer.cc
Go to the documentation of this file.
3 
5 {
6  moduleLabel_=cfg.getParameter<std::string>("@module_label");
7  srcToken_ = consumes<CandidateView>(cfg.getParameter<edm::InputTag>("src"));
8  shiftBy_ = cfg.getParameter<double>("shiftBy");
9 
10  if ( cfg.exists("binning") ) {
11  typedef std::vector<edm::ParameterSet> vParameterSet;
12  vParameterSet cfgBinning = cfg.getParameter<vParameterSet>("binning");
13  for ( vParameterSet::const_iterator cfgBinningEntry = cfgBinning.begin();
14  cfgBinningEntry != cfgBinning.end(); ++cfgBinningEntry ) {
15  binning_.push_back(new binningEntryType(*cfgBinningEntry,moduleLabel_));
16  }
17  } else {
18  std::string uncertainty = cfg.getParameter<std::string>("uncertainty");
19  binning_.push_back(new binningEntryType(uncertainty, moduleLabel_));
20  }
21 
22  produces<reco::CandidateCollection>();
23 }
24 
26 {
27  for(std::vector<binningEntryType*>::const_iterator it = binning_.begin();
28  it != binning_.end(); ++it ) {
29  delete (*it);
30  }
31 }
32 
33 void
35 {
36  edm::Handle<CandidateView> originalParticles;
37  evt.getByToken(srcToken_, originalParticles);
38 
39  auto shiftedParticles = std::make_unique<reco::CandidateCollection>();
40 
41  for(CandidateView::const_iterator originalParticle = originalParticles->begin();
42  originalParticle != originalParticles->end(); ++originalParticle ) {
43 
44  double uncertainty = getUncShift(originalParticle);
45  double shift = shiftBy_*uncertainty;
46 
47  reco::Candidate::LorentzVector shiftedParticleP4 = originalParticle->p4();
48  //leave 0*nan = 0
49  if (! (edm::isNotFinite(shift) && shiftedParticleP4.mag2()==0)) shiftedParticleP4 *= (1. + shift);
50 
51  std::unique_ptr<reco::Candidate> shiftedParticle = std::make_unique<reco::LeafCandidate>(*originalParticle);
52  shiftedParticle->setP4(shiftedParticleP4);
53 
54  shiftedParticles->push_back( shiftedParticle.release() );
55  }
56 
57  evt.put(std::move(shiftedParticles));
58 }
59 
60 double
62  double valx=0;
63  double valy=0;
64  for(std::vector<binningEntryType*>::iterator binningEntry=binning_.begin();
65  binningEntry!=binning_.end(); ++binningEntry ) {
66  if( (!(*binningEntry)->binSelection_) || (*(*binningEntry)->binSelection_)(*originalParticle) ) {
67 
68  if( (*binningEntry)->energyDep_ ) valx=originalParticle->energy();
69  else valx=originalParticle->pt();
70 
71  valy=originalParticle->eta();
72  return (*binningEntry)->binUncFormula_->Eval(valx, valy);
73  }
74  }
75  return 0;
76 }
77 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< CandidateView > srcToken_
bool exists(std::string const &parameterName) const
checks if a parameter exists
void produce(edm::Event &evt, const edm::EventSetup &es) override
ShiftedParticleProducer(const edm::ParameterSet &cfg)
const_iterator begin() const
bool isNotFinite(T x)
Definition: isFinite.h:10
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
double getUncShift(const CandidateView::const_iterator &originalParticle)
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
static unsigned int const shift
const_iterator end() const
def move(src, dest)
Definition: eostools.py:510
std::vector< binningEntryType * > binning_