CMS 3D CMS Logo

DQMLumiMonitor.cc
Go to the documentation of this file.
1 /*
2  * \file DQMLumiMonitor.cc
3  * \author S. Dutta
4  * Last Update:
5  *
6  * Description: Pixel Luminosity Monitoring
7  *
8 */
18 #include "TPRegexp.h"
19 
20 
21 // -----------------------------
22 // constructors and destructor
23 // -----------------------------
24 
25 DQMLumiMonitor::DQMLumiMonitor( const edm::ParameterSet& ps ) : parameters_(ps) {
26 
27 
30  pixelClusterInputTag_= consumes<edmNew::DetSetVector<SiPixelCluster> >(parameters_.getParameter<edm::InputTag>("PixelClusterInputTag"));
31  lumiRecordName_ = consumes<LumiSummary,edm::InLumi>(parameters_.getParameter<std::string>("LumiRecordName"));
32 
33  nClusME_ = nullptr;
34  nClusVsLSME_ = nullptr;
35  intLumiVsLSME_= nullptr;
36  corrIntLumiAndClusVsLSME_ = nullptr;
37 }
38 
40 
42 
43  edm::ParameterSet ClusHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1ClusPar");
44  edm::ParameterSet LumiHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1LumiPar");
45  edm::ParameterSet LumiSecHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1LSPar");
46 
47  std::string currentFolder = moduleName_ + "/" + folderName_ ;
48  dbe_->setCurrentFolder(currentFolder);
49 
50  if (nClusME_ == nullptr) nClusME_ = dbe_->book1D("nPlxClus", " Number of Pixel Clusters ",
51  ClusHistoPar.getParameter<int32_t>("Xbins"),
52  ClusHistoPar.getParameter<double>("Xmin"),
53  ClusHistoPar.getParameter<double>("Xmax"));
54  else nClusME_->Reset();
55  if (nClusVsLSME_ == nullptr) nClusVsLSME_ = dbe_->bookProfile("nClusVsLS", " Number of Pixel Cluster Vs LS number",
56  LumiSecHistoPar.getParameter<int32_t>("Xbins"),
57  LumiSecHistoPar.getParameter<double>("Xmin"),
58  LumiSecHistoPar.getParameter<double>("Xmax"),
59  0.0, 0.0, "");
60  else nClusVsLSME_->Reset();
61  if (intLumiVsLSME_ == nullptr) intLumiVsLSME_ = dbe_->bookProfile("intLumiVsLS", " Integrated Luminosity Vs LS number",
62  LumiSecHistoPar.getParameter<int32_t>("Xbins"),
63  LumiSecHistoPar.getParameter<double>("Xmin"),
64  LumiSecHistoPar.getParameter<double>("Xmax"),
65  0.0, 0.0, "");
66  else intLumiVsLSME_->Reset();
67 
68  if (corrIntLumiAndClusVsLSME_== nullptr) corrIntLumiAndClusVsLSME_ = dbe_->bookProfile2D("corrIntLumiAndClusVsLS", " Correlation of nCluster and Integrated Luminosity Vs LS number",
69  LumiSecHistoPar.getParameter<int32_t>("Xbins"),
70  LumiSecHistoPar.getParameter<double>("Xmin"),
71  LumiSecHistoPar.getParameter<double>("Xmax"),
72  LumiHistoPar.getParameter<int32_t>("Xbins"),
73  LumiHistoPar.getParameter<double>("Xmin"),
74  LumiHistoPar.getParameter<double>("Xmax"),
75  0.0, 0.0);
77 }
80  intLumi_ = -1.0;
81  nLumi_ = -1;
82 }
83 
84 void DQMLumiMonitor::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
86 
87 }
89 
90  //Access Pixel Clusters
92  iEvent.getByToken(pixelClusterInputTag_, siPixelClusters);
93 
94  if(!siPixelClusters.isValid()) {
95  edm::LogError("PixelLumiMonotor") << "Could not find Cluster Collection ";
96  return;
97  }
98  unsigned int nClusterPix = (*siPixelClusters).dataSize();
99  nClusME_->Fill(nClusterPix);
100  if (nLumi_ != -1) nClusVsLSME_->Fill(nLumi_, nClusterPix);
101  if (intLumi_ != -1 || nLumi_ != -1) corrIntLumiAndClusVsLSME_->Fill(nLumi_, intLumi_, nClusterPix);
102 }
103 
105  edm::LogInfo("PixelLumiMonotor") <<" Run Number "<<lumiBlock.run() <<" Lumi Section Numnber "<< lumiBlock.luminosityBlock();
106 
107  nLumi_ = lumiBlock.luminosityBlock();
108 
109  // Access Lumi Summary
110  edm::Handle<LumiSummary> lumiSummary_;
111  lumiBlock.getByToken(lumiRecordName_, lumiSummary_);
112  if(lumiSummary_->isValid()){
113  intLumi_ = lumiSummary_->intgDelLumi();
114  edm::LogInfo("PixelLumiMonotor") <<" Luminosity in this Lumi Section " << intLumi_ ;
116  } else {
117  edm::LogError("PixelLumiMonotor") << "No valid data found!";
118  }
119  /*
120  // Access Lumi Details
121  Handle<LumiDetails> lumiDetails;
122  lumiBlock.getByLabel("expressLumiProducer", lumiDetails);
123  if(lumiDetails->isValid()){
124  std::cout<<"valid detail"<<std::endl;
125  }else{
126  std::cout << "no valid lumi detail data" <<std::endl;
127  } */
128 
129 }
130 
131 void DQMLumiMonitor::endRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
132 
133 }
134 
135 
137 
138 }
139 // Define this as a plug-in
T getParameter(std::string const &) const
edm::EDGetTokenT< LumiSummary > lumiRecordName_
MonitorElement * corrIntLumiAndClusVsLSME_
void endRun(edm::Run const &iRun, edm::EventSetup const &iSetup) override
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:1040
edm::ParameterSet parameters_
void endJob() override
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
DQMStore * dbe_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusterInputTag_
void Fill(long long x)
LuminosityBlockNumber_t luminosityBlock() const
void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &eSetup) override
int iEvent
Definition: GenABIO.cc:230
float intgDelLumi() const
Definition: LumiSummary.cc:17
RunNumber_t run() const
MonitorElement * bookProfile(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const char *option="s")
Definition: DQMStore.cc:1372
void beginRun(edm::Run const &iRun, edm::EventSetup const &iSetup) override
void Reset()
reset ME (ie. contents, errors, etc)
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * nClusVsLSME_
DQMLumiMonitor(const edm::ParameterSet &)
std::string moduleName_
~DQMLumiMonitor() override
void beginJob() override
MonitorElement * nClusME_
std::string folderName_
bool isValid() const
Definition: LumiSummary.cc:78
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:746
MonitorElement * intLumiVsLSME_
Definition: Run.h:43
MonitorElement * bookProfile2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, int nchZ, double lowZ, double highZ, const char *option="s")
Definition: DQMStore.cc:1516