CMS 3D CMS Logo

JetCorrectorDBReader.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: JetCorrectorDBReader
4 // Class:
5 //
13 //
14 // Original Author: Benedikt Hegner
15 // Created: Tue Mar 09 01:32:51 CET 2010
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
34 //
35 // class declaration
36 //
37 
39 public:
40  explicit JetCorrectorDBReader(const edm::ParameterSet&);
41  ~JetCorrectorDBReader() override;
42 
43 
44 private:
45  void beginJob() override ;
46  void analyze(const edm::Event&, const edm::EventSetup&) override;
47  void endJob() override ;
48 
51 };
52 
53 
55 {
56  mPayloadName = iConfig.getUntrackedParameter<std::string>("payloadName");
57  mGlobalTag = iConfig.getUntrackedParameter<std::string>("globalTag");
58  mPrintScreen = iConfig.getUntrackedParameter<bool>("printScreen");
59  mCreateTextFile = iConfig.getUntrackedParameter<bool>("createTextFile");
60 }
61 
62 
64 {
65 
66 }
67 
69 {
71  std::cout <<"Inspecting JEC payload with label: "<< mPayloadName <<std::endl;
72  iSetup.get<JetCorrectionsRecord>().get(mPayloadName,JetCorParamsColl);
73  std::vector<JetCorrectorParametersCollection::key_type> keys;
74  JetCorParamsColl->validKeys( keys );
75  for ( std::vector<JetCorrectorParametersCollection::key_type>::const_iterator ibegin = keys.begin(),
76  iend = keys.end(), ikey = ibegin; ikey != iend; ++ikey ) {
77  std::cout<<"-------------------------------------------------" << std::endl;
78  std::cout<<"Processing key = " << *ikey << std::endl;
79  std::cout<<"object label: "<<JetCorParamsColl->findLabel(*ikey)<<std::endl;
80  JetCorrectorParameters const & JetCorParams = (*JetCorParamsColl)[*ikey];
81 
82  if (mCreateTextFile)
83  {
84  std::cout<<"Creating txt file: "<<mGlobalTag+"_"+mPayloadName+"_"+JetCorParamsColl->findLabel(*ikey)+".txt"<<std::endl;
85  JetCorParams.printFile(mGlobalTag+"_"+JetCorParamsColl->findLabel(*ikey)+"_"+mPayloadName+".txt");
86  }
87  if (mPrintScreen)
88  JetCorParams.printScreen();
89  }
90 }
91 
92 void
94 {
95 }
96 
97 void
99 {
100 }
101 
102 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
JetCorrectorDBReader(const edm::ParameterSet &)
static std::string findLabel(key_type k)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void printFile(const std::string &fFileName) const
void analyze(const edm::Event &, const edm::EventSetup &) override
T get() const
Definition: EventSetup.h:71
void validKeys(std::vector< key_type > &keys) const