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 
41 }
42 
44 
45  edm::ParameterSet ClusHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1ClusPar");
46  edm::ParameterSet LumiHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1LumiPar");
47  edm::ParameterSet LumiSecHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1LSPar");
48 
49  std::string currentFolder = moduleName_ + "/" + folderName_ ;
50  dbe_->setCurrentFolder(currentFolder);
51 
52  if (nClusME_ == nullptr) nClusME_ = dbe_->book1D("nPlxClus", " Number of Pixel Clusters ",
53  ClusHistoPar.getParameter<int32_t>("Xbins"),
54  ClusHistoPar.getParameter<double>("Xmin"),
55  ClusHistoPar.getParameter<double>("Xmax"));
56  else nClusME_->Reset();
57  if (nClusVsLSME_ == nullptr) nClusVsLSME_ = dbe_->bookProfile("nClusVsLS", " Number of Pixel Cluster Vs LS number",
58  LumiSecHistoPar.getParameter<int32_t>("Xbins"),
59  LumiSecHistoPar.getParameter<double>("Xmin"),
60  LumiSecHistoPar.getParameter<double>("Xmax"),
61  0.0, 0.0, "");
62  else nClusVsLSME_->Reset();
63  if (intLumiVsLSME_ == nullptr) intLumiVsLSME_ = dbe_->bookProfile("intLumiVsLS", " Integrated Luminosity Vs LS number",
64  LumiSecHistoPar.getParameter<int32_t>("Xbins"),
65  LumiSecHistoPar.getParameter<double>("Xmin"),
66  LumiSecHistoPar.getParameter<double>("Xmax"),
67  0.0, 0.0, "");
68  else intLumiVsLSME_->Reset();
69 
70  if (corrIntLumiAndClusVsLSME_== nullptr) corrIntLumiAndClusVsLSME_ = dbe_->bookProfile2D("corrIntLumiAndClusVsLS", " Correlation of nCluster and Integrated Luminosity Vs LS number",
71  LumiSecHistoPar.getParameter<int32_t>("Xbins"),
72  LumiSecHistoPar.getParameter<double>("Xmin"),
73  LumiSecHistoPar.getParameter<double>("Xmax"),
74  LumiHistoPar.getParameter<int32_t>("Xbins"),
75  LumiHistoPar.getParameter<double>("Xmin"),
76  LumiHistoPar.getParameter<double>("Xmax"),
77  0.0, 0.0);
79 }
82  intLumi_ = -1.0;
83  nLumi_ = -1;
84 }
85 
86 void DQMLumiMonitor::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
88 
89 }
91 
92  //Access Pixel Clusters
94  iEvent.getByToken(pixelClusterInputTag_, siPixelClusters);
95 
96  if(!siPixelClusters.isValid()) {
97  edm::LogError("PixelLumiMonotor") << "Could not find Cluster Collection ";
98  return;
99  }
100  unsigned int nClusterPix = (*siPixelClusters).dataSize();
101  nClusME_->Fill(nClusterPix);
102  if (nLumi_ != -1) nClusVsLSME_->Fill(nLumi_, nClusterPix);
103  if (intLumi_ != -1 || nLumi_ != -1) corrIntLumiAndClusVsLSME_->Fill(nLumi_, intLumi_, nClusterPix);
104 }
105 
107  edm::LogInfo("PixelLumiMonotor") <<" Run Number "<<lumiBlock.run() <<" Lumi Section Numnber "<< lumiBlock.luminosityBlock();
108 
109  nLumi_ = lumiBlock.luminosityBlock();
110 
111  // Access Lumi Summary
112  edm::Handle<LumiSummary> lumiSummary_;
113  lumiBlock.getByToken(lumiRecordName_, lumiSummary_);
114  if(lumiSummary_->isValid()){
115  intLumi_ = lumiSummary_->intgDelLumi();
116  edm::LogInfo("PixelLumiMonotor") <<" Luminosity in this Lumi Section " << intLumi_ ;
118  } else {
119  edm::LogError("PixelLumiMonotor") << "No valid data found!";
120  }
121  /*
122  // Access Lumi Details
123  Handle<LumiDetails> lumiDetails;
124  lumiBlock.getByLabel("expressLumiProducer", lumiDetails);
125  if(lumiDetails->isValid()){
126  std::cout<<"valid detail"<<std::endl;
127  }else{
128  std::cout << "no valid lumi detail data" <<std::endl;
129  } */
130 
131 }
132 
133 void DQMLumiMonitor::endRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
134 
135 }
136 
137 
139 
140 }
141 // 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:1035
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:508
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:1367
void beginRun(edm::Run const &iRun, edm::EventSetup const &iSetup) override
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 Reset(void)
reset ME (ie. contents, errors, etc)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:741
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:1511