CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DQMHcalPhiSymHLT.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 /*
4  Description: <one line class summary>
5 
6  Implementation:
7  <Notes on implementation>
8 */
9 //
10 // Original Author: Grigory Safronov
11 // Created: Thu Sep 10 08:28:14 CEST 2009
12 // $Id$
13 //
14 //
15 
16 
17 #include <memory>
18 
21 
24 
26 
30 
34 
36 public:
37  explicit DQMHcalPhiSymHLT(const edm::ParameterSet&);
39 
44 
45 
46  std::string folderName_;
47  std::string outRootFileName_;
50 
52 
53 private:
54  virtual void beginJob() ;
55  virtual void analyze(const edm::Event&, const edm::EventSetup&);
56  virtual void endJob() ;
57 
59  int iEvt;
60 };
61 
62 
64 {
65  folderName_ = iConfig.getParameter<std::string>("folderName");
66  outRootFileName_=iConfig.getParameter<std::string>("outputRootFileName");
67  rawInLabel_=iConfig.getParameter<edm::InputTag>("rawInputLabel");
68  saveToRootFile_=iConfig.getParameter<bool>("SaveToRootFile");
69 
70  iEvt=0;
71 }
72 
74 {
75 }
76 
77 
78 void
80 {
81  if (iEvt==0) firstLumiSec=iEvent.luminosityBlock();
82  iEvt++;
83 
84  std::auto_ptr<FEDRawDataCollection> producedData(new FEDRawDataCollection);
85 
87  iEvent.getByLabel(rawInLabel_,rawIn);
88 
89  std::vector<int> selFEDs;
90 
91  //get HCAL FEDs:
93  {
94  selFEDs.push_back(i);
95  }
96 
97  const FEDRawDataCollection *rdc=rawIn.product();
98 
99  //calculate full HCAL data size:
100  size_t hcalSize=0;
101  for (unsigned int k=0; k<selFEDs.size(); k++)
102  {
103  const FEDRawData & fedData = rdc->FEDData(selFEDs[k]);
104  hcalSize+=fedData.size();
105  hFEDsize->Fill(fedData.size()*pow(1024,-1),1);
106  }
107  hHCALsize->Fill(hcalSize*pow(1024,-1),1);
108  hHCALvsLumiSec->Fill(iEvent.luminosityBlock()-firstLumiSec,hcalSize*pow(1024,-1),1);
109 
110  //calculate full data size:
111  size_t fullSize=0;
112  for (int j=0; j<FEDNumbering::MAXFEDID; ++j )
113  {
114  const FEDRawData & fedData = rdc->FEDData(j);
115  fullSize+=fedData.size();
116  }
117  hFULLsize->Fill(fullSize*pow(1024,-1),1);
118 }
119 
120 
121 // ------------ method called once each job just before starting event loop ------------
122 void
124 {
127 
128  hFEDsize=dbe_->book1D("hFEDsize","HCAL FED size (kB)",1000,0,100);
129  hFEDsize->setAxisTitle("kB",1);
130 
131  hHCALsize=dbe_->book1D("hHCALsize","HCAL data size (kB)",1000,0,1000);
132  hHCALsize->setAxisTitle("kB",1);
133 
134  hFULLsize=dbe_->book1D("hFULLsize","Full data size (kB)",1000,0,2000);
135  hFULLsize->setAxisTitle("kB",1);
136 
137  hHCALvsLumiSec=dbe_->book2D("hHCALvsLumiSec","HCAL data size (kB) vs. internal lumi block number",10000,0,10000,1000,0,1000);
138  hHCALvsLumiSec->setAxisTitle("kB",2);
139  hHCALvsLumiSec->setAxisTitle("internal luminosity block number",1);
140 }
141 
142 void
144 
145 {
146  if(dbe_&&saveToRootFile_)
147  {
149  }
150 }
151 
152 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual void beginJob()
MonitorElement * hFULLsize
virtual void endJob()
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:514
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::string outRootFileName_
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:1898
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
MonitorElement * hFEDsize
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
void Fill(long long x)
DQMHcalPhiSymHLT(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:243
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
int j
Definition: DBlmapReader.cc:9
edm::InputTag rawInLabel_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
int k[5][pyjets_maxn]
MonitorElement * hHCALsize
T const * product() const
Definition: Handle.h:74
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::string folderName_
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:642
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
MonitorElement * hHCALvsLumiSec