CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EmbeddingKineReweightProducer.cc
Go to the documentation of this file.
2 
4 
6 
8 
10 
11 #include <vector>
12 #include <algorithm>
13 #include <math.h>
14 
16 {
17  srcGenDiTaus_ = cfg.getParameter<edm::InputTag>("srcGenDiTaus");
18 
20  if ( inputFileName.location() == edm::FileInPath::Unknown)
21  throw cms::Exception("EmbeddingReweightProducer")
22  << " Failed to find File = " << inputFileName << " !!\n";
23  std::auto_ptr<TFile> inputFile(new TFile(inputFileName.fullPath().data()));
24 
25  edm::ParameterSet cfgLUTs = cfg.getParameter<edm::ParameterSet>("lutNames");
26  std::vector<std::string> variables = cfgLUTs.getParameterNamesForType<std::string>();
27  for ( std::vector<std::string>::const_iterator variable = variables.begin();
28  variable != variables.end(); ++variable ) {
29  std::string lutName = cfgLUTs.getParameter<std::string>(*variable);
30  lutEntryType* lutEntry =
31  new lutEntryType(*inputFile, *variable, lutName);
32  lutEntries_.push_back(lutEntry);
33  }
34 
35  minWeight_ = cfg.getParameter<double>("minWeight");
36  maxWeight_ = cfg.getParameter<double>("maxWeight");
37 
38  verbosity_ = ( cfg.exists("verbosity") ) ?
39  cfg.getParameter<int>("verbosity") : 0;
40 
41  for ( std::vector<lutEntryType*>::iterator lutEntry = lutEntries_.begin();
42  lutEntry != lutEntries_.end(); ++lutEntry ) {
43  produces<double>((*lutEntry)->variableName_);
44  }
45 }
46 
48 {
49  for ( std::vector<lutEntryType*>::iterator it = lutEntries_.begin();
50  it != lutEntries_.end(); ++it ) {
51  delete (*it);
52  }
53 }
54 
56 {
57  if ( verbosity_ ) {
58  std::cout << "<EmbeddingKineReweightProducer::produce>:" << std::endl;
59  }
60 
63  evt.getByLabel(srcGenDiTaus_, genDiTaus);
64  if ( genDiTaus->size() != 1 )
65  throw cms::Exception("EmbeddingKineReweightProducer")
66  << "Failed to find unique genDiTau object !!\n";
67  const reco::Candidate& genDiTau = genDiTaus->front();
68  if ( verbosity_ ) {
69  std::cout << "diTau: Pt = " << genDiTau.pt() << ", eta = " << genDiTau.eta() << ", phi = " << genDiTau.phi() << ", mass = " << genDiTau.mass() << std::endl;
70  }
71 
72  for ( std::vector<lutEntryType*>::const_iterator lutEntry = lutEntries_.begin();
73  lutEntry != lutEntries_.end(); ++lutEntry ) {
74  double weight = (**lutEntry)(genDiTau);
75  if ( weight < minWeight_ ) weight = minWeight_;
76  if ( weight > maxWeight_ ) weight = maxWeight_;
77  if ( verbosity_ ) {
78  std::cout << " " << (*lutEntry)->variableName_ << " = " << weight << std::endl;
79  }
80  std::auto_ptr<double> weightPtr(new double(weight));
81  evt.put(weightPtr, (*lutEntry)->variableName_);
82  }
83 }
84 
86 
88 
89 
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:259
virtual double pt() const =0
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual double mass() const =0
mass
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< std::string > getParameterNamesForType(bool trackiness=true) const
Definition: ParameterSet.h:194
edm::View< reco::Candidate > CandidateView
std::vector< lutEntryType * > lutEntries_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:413
void produce(edm::Event &, const edm::EventSetup &)
tuple cout
Definition: gather_cfg.py:121
int weight
Definition: histoStyle.py:50
EmbeddingKineReweightProducer(const edm::ParameterSet &)
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity