CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Private Attributes
DQMHcalPhiSymAlCaReco Class Reference

#include <DQMHcalPhiSymAlCaReco.h>

Inheritance diagram for DQMHcalPhiSymAlCaReco:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 DQMHcalPhiSymAlCaReco (const edm::ParameterSet &)
 
 ~DQMHcalPhiSymAlCaReco ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c)
 
void beginJob ()
 
void beginLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context)
 
void beginRun (const edm::Run &r, const edm::EventSetup &c)
 
void endJob ()
 
void endLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &c)
 
void endRun (const edm::Run &r, const edm::EventSetup &c)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Private Attributes

DQMStoredbe_
 
int eventCounter_
 
std::string fileName_
 Output file name if required. More...
 
std::string folderName_
 DQM folder name. More...
 
edm::EDGetTokenT
< HBHERecHitCollection
hbherecoMB
 object to monitor More...
 
edm::EDGetTokenT
< HBHERecHitCollection
hbherecoNoise
 
MonitorElementhFEDsize
 
edm::EDGetTokenT
< HFRecHitCollection
hfrecoMB
 
edm::EDGetTokenT
< HFRecHitCollection
hfrecoNoise
 
MonitorElementhHcalIsZS
 
int hiDistr_r_nbin_
 
double hiDistr_x_max_
 
double hiDistr_x_min_
 
int hiDistr_x_nbin_
 
double hiDistr_y_max_
 
double hiDistr_y_min_
 
int hiDistr_y_nbin_
 
MonitorElementhiDistrHBHEsize1D_
 
MonitorElementhiDistrHFsize1D_
 
MonitorElementhiDistrMB2Min2D_
 
MonitorElementhiDistrMB2Pl2D_
 
MonitorElementhiDistrMBMin2D_
 
MonitorElementhiDistrMBPl2D_
 
MonitorElementhiDistrNoise2Min2D_
 
MonitorElementhiDistrNoise2Pl2D_
 
MonitorElementhiDistrNoiseMin2D_
 
MonitorElementhiDistrNoisePl2D_
 
MonitorElementhiDistrVarMBMin2D_
 
MonitorElementhiDistrVarMBPl2D_
 
MonitorElementhiDistrVarNoiseMin2D_
 
MonitorElementhiDistrVarNoisePl2D_
 
MonitorElementhL1Id
 
edm::InputTag horecoMB
 
edm::InputTag horecoNoise
 
double ihbhe_size_
 
double ihf_size_
 
unsigned int period_
 
edm::EDGetTokenT
< FEDRawDataCollection
rawInLabel_
 
bool saveToFile_
 Write to file. More...
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 

Detailed Description

Definition at line 24 of file DQMHcalPhiSymAlCaReco.h.

Constructor & Destructor Documentation

DQMHcalPhiSymAlCaReco::DQMHcalPhiSymAlCaReco ( const edm::ParameterSet ps)

Definition at line 45 of file DQMHcalPhiSymAlCaReco.cc.

References dbe_, fileName_, folderName_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hbherecoMB, hbherecoNoise, hfrecoMB, hfrecoNoise, hiDistr_r_nbin_, hiDistr_x_max_, hiDistr_x_min_, hiDistr_x_nbin_, hiDistr_y_max_, hiDistr_y_min_, hiDistr_y_nbin_, horecoMB, horecoNoise, ihbhe_size_, ihf_size_, cppFunctionSkipper::operator, period_, rawInLabel_, and saveToFile_.

45  :
47 {
49  //
50  // Input from configurator file
51  //
52  folderName_ = ps.getUntrackedParameter<string>("FolderName","ALCAStreamHcalPhiSym");
53 
54  hbherecoMB = consumes<HBHERecHitCollection>(ps.getParameter<edm::InputTag>("hbheInputMB"));
55  horecoMB = ps.getParameter<edm::InputTag>("hoInputMB");
56  hfrecoMB = consumes<HFRecHitCollection>(ps.getParameter<edm::InputTag>("hfInputMB"));
57 
58  hbherecoNoise = consumes<HBHERecHitCollection>(ps.getParameter<edm::InputTag>("hbheInputNoise"));
59  horecoNoise = ps.getParameter<edm::InputTag>("hoInputNoise");
60  hfrecoNoise = consumes<HFRecHitCollection>(ps.getParameter<edm::InputTag>("hfInputNoise"));
61 
62  rawInLabel_ = consumes<FEDRawDataCollection>(ps.getParameter<edm::InputTag>("rawInputLabel"));
63 
64  period_ = ps.getParameter<unsigned int>("period") ;
65 
66  saveToFile_ = ps.getUntrackedParameter<bool>("SaveToFile",false);
67  fileName_ = ps.getUntrackedParameter<string>("FileName","MonitorAlCaHcalPhiSym.root");
68 
69  // histogram parameters
70 
71  // Distribution of rechits in iPhi, iEta
72  hiDistr_y_nbin_ = ps.getUntrackedParameter<int>("hiDistr_y_nbin",72);
73  hiDistr_y_min_ = ps.getUntrackedParameter<double>("hiDistr_y_min",0.5);
74  hiDistr_y_max_ = ps.getUntrackedParameter<double>("hiDistr_y_max",72.5);
75  hiDistr_x_nbin_ = ps.getUntrackedParameter<int>("hiDistr_x_nbin",41);
76  hiDistr_x_min_ = ps.getUntrackedParameter<double>("hiDistr_x_min",0.5);
77  hiDistr_x_max_ = ps.getUntrackedParameter<double>("hiDistr_x_max",41.5);
78  // Check for NZS
79  hiDistr_r_nbin_ = ps.getUntrackedParameter<int>("hiDistr_r_nbin",100);
80  ihbhe_size_ = ps.getUntrackedParameter<double>("ihbhe_size_",5184.);
81  ihf_size_ = ps.getUntrackedParameter<double>("ihf_size_",1728.);
82 
83 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool saveToFile_
Write to file.
edm::EDGetTokenT< HBHERecHitCollection > hbherecoNoise
edm::EDGetTokenT< FEDRawDataCollection > rawInLabel_
edm::EDGetTokenT< HFRecHitCollection > hfrecoNoise
std::string folderName_
DQM folder name.
edm::EDGetTokenT< HFRecHitCollection > hfrecoMB
std::string fileName_
Output file name if required.
edm::EDGetTokenT< HBHERecHitCollection > hbherecoMB
object to monitor
DQMHcalPhiSymAlCaReco::~DQMHcalPhiSymAlCaReco ( )

Definition at line 85 of file DQMHcalPhiSymAlCaReco.cc.

86 {}

Member Function Documentation

void DQMHcalPhiSymAlCaReco::analyze ( const edm::Event e,
const edm::EventSetup c 
)
protectedvirtual

Level-1 event number generated by the TTC system

Implements edm::EDAnalyzer.

Definition at line 309 of file DQMHcalPhiSymAlCaReco.cc.

References edm::SortedCollection< T, SORT >::begin(), FEDRawData::data(), HcalDetId::depth(), edm::SortedCollection< T, SORT >::end(), eventCounter_, edm::false, FEDRawDataCollection::FEDData(), MonitorElement::Fill(), edm::Event::getByToken(), HcalDCCHeader::getSpigotData(), HcalDCCHeader::getSpigotPresent(), hbherecoMB, hbherecoNoise, hFEDsize, hfrecoMB, hfrecoNoise, hHcalIsZS, hiDistrHBHEsize1D_, hiDistrHFsize1D_, hiDistrMB2Min2D_, hiDistrMB2Pl2D_, hiDistrMBMin2D_, hiDistrMBPl2D_, hiDistrNoise2Min2D_, hiDistrNoise2Pl2D_, hiDistrNoiseMin2D_, hiDistrNoisePl2D_, hL1Id, i, HcalDetId::ieta(), ihbhe_size_, ihf_size_, HcalDetId::iphi(), HcalHTRData::isUnsuppressed(), edm::HandleBase::isValid(), gen::k, LogDebug, FEDNumbering::MAXHCALFEDID, FEDNumbering::MINHCALFEDID, FEDNumbering::MINTriggerGTPFEDID, period_, funct::pow(), edm::Handle< T >::product(), rawInLabel_, reco::return(), FEDRawData::size(), edm::SortedCollection< T, SORT >::size(), and HcalDCCHeader::SPIGOT_COUNT.

310  {
311 
312 
313  eventCounter_++;
314 
316  iEvent.getByToken(rawInLabel_,rawIn);
317 
318  if(!rawIn.isValid()){
319  LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
320  return ;
321  }
322 
323  //get HCAL FEDs:
324  std::vector<int> selFEDs;
326  {
327  selFEDs.push_back(i);
328  }
329 
330  // std::cout<<" Size of FED "<<selFEDs.size()<<std::endl;
331 
332  const FEDRawDataCollection *rdc=rawIn.product();
333 
334  bool hcalIsZS = false ;
335  for (unsigned int k=0; k<selFEDs.size(); k++)
336  {
337  const FEDRawData & fedData = rdc->FEDData(selFEDs[k]);
338  // std::cout<<fedData.size()*std::pow(1024.,-1)<<std::endl;
339  hFEDsize->Fill(fedData.size()*std::pow(1024.,-1),1);
340 
341  // get HCAL DCC Header for each FEDRawData
342  const HcalDCCHeader* dccHeader=(const HcalDCCHeader*)(fedData.data());
343 
344  if (!dccHeader) continue;
345 
346  // walk through the HTR data...
347  HcalHTRData htr;
348 
349  int nspigot =0;
350  for (int spigot=0; spigot<HcalDCCHeader::SPIGOT_COUNT; spigot++) {
351  nspigot++;
352 
353  if (!dccHeader->getSpigotPresent(spigot)) continue;
354 
355  // Load the given decoder with the pointer and length from this spigot.
356  dccHeader->getSpigotData(spigot,htr, fedData.size());
357 
358  if(k != 20 && nspigot !=14 ) {
359  if ( !htr.isUnsuppressed() ) { hcalIsZS = true; }
360  }
361  }
362 
363  } // loop over HcalFEDs
364 
365  hHcalIsZS->Fill( hcalIsZS );
366 
367  // get Trigger FED-Id
368  const FEDRawData& fedData = rdc->FEDData(FEDNumbering::MINTriggerGTPFEDID) ;
369  FEDHeader header(fedData.data()) ;
370 
372  hL1Id->Fill( (header.lvl1ID())%period_ );
373 
375  iEvent.getByToken(hbherecoNoise, hbheNS);
376 
377  if(!hbheNS.isValid()){
378  LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
379  return ;
380  }
381 
383  iEvent.getByToken(hbherecoMB, hbheMB);
384 
385  if(!hbheMB.isValid()){
386  LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
387  return ;
388  }
389 
391  iEvent.getByToken(hfrecoNoise, hfNS);
392 
393  if(!hfNS.isValid()){
394  LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
395  return ;
396  }
397 
399  iEvent.getByToken(hfrecoMB, hfMB);
400 
401  if(!hfMB.isValid()){
402  LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
403  return ;
404  }
405 
406  const HBHERecHitCollection HithbheNS = *(hbheNS.product());
407 
408  hiDistrHBHEsize1D_->Fill(HithbheNS.size()/ihbhe_size_);
409 
410 
411  for(HBHERecHitCollection::const_iterator hbheItr=HithbheNS.begin(); hbheItr!=HithbheNS.end(); hbheItr++)
412  {
413  DetId id = (*hbheItr).detid();
414  HcalDetId hid=HcalDetId(id);
415 
416  if(hid.depth() == 1) {
417  if( hid.ieta() > 0 ) {
418  hiDistrNoisePl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
419  hiDistrNoise2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
420  } else {
421  hiDistrNoiseMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
422  hiDistrNoise2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
423  }
424  }
425  }
426 
427  const HBHERecHitCollection HithbheMB = *(hbheMB.product());
428 
429  for(HBHERecHitCollection::const_iterator hbheItr=HithbheMB.begin(); hbheItr!=HithbheMB.end(); hbheItr++)
430  {
431  DetId id = (*hbheItr).detid();
432  HcalDetId hid=HcalDetId(id);
433 
434  if(hid.depth() == 1) {
435  if( hid.ieta() > 0 ) {
436  hiDistrMBPl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
437  hiDistrMB2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
438  } else {
439  hiDistrMBMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
440  hiDistrMB2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
441  }
442  }
443 
444  }
445 
446  const HFRecHitCollection HithfNS = *(hfNS.product());
447 
448  hiDistrHFsize1D_->Fill(HithfNS.size()/ihf_size_);
449 
450  for(HFRecHitCollection::const_iterator hbheItr=HithfNS.begin(); hbheItr!=HithfNS.end(); hbheItr++)
451  {
452 
453  DetId id = (*hbheItr).detid();
454  HcalDetId hid=HcalDetId(id);
455 
456  if(hid.depth() == 1) {
457  if( hid.ieta() > 0 ) {
458  hiDistrNoisePl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
459  hiDistrNoise2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
460  } else {
461  hiDistrNoiseMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
462  hiDistrNoise2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
463  }
464  }
465 
466  }
467 
468  const HFRecHitCollection HithfMB = *(hfMB.product());
469 
470  for(HFRecHitCollection::const_iterator hbheItr=HithfMB.begin(); hbheItr!=HithfMB.end(); hbheItr++)
471  {
472  DetId id = (*hbheItr).detid();
473  HcalDetId hid=HcalDetId(id);
474 
475  if(hid.depth() == 1) {
476  if( hid.ieta() > 0 ) {
477  hiDistrMBPl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
478  hiDistrMB2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
479  } else {
480  hiDistrMBMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
481  hiDistrMB2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
482  }
483  }
484  }
485 
486 
487 } //analyze
#define LogDebug(id)
MonitorElement * hiDistrNoisePl2D_
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< HBHERecHitCollection > hbherecoNoise
std::vector< HBHERecHit >::const_iterator const_iterator
edm::EDGetTokenT< FEDRawDataCollection > rawInLabel_
int getSpigotData(int nspigot, HcalHTRData &decodeTool, int validSize) const
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
void Fill(long long x)
bool isUnsuppressed() const
Is this event an unsuppresed event?
Definition: HcalHTRData.cc:353
MonitorElement * hiDistrMBPl2D_
MonitorElement * hiDistrHFsize1D_
int depth() const
get the tower depth
Definition: HcalDetId.h:40
edm::EDGetTokenT< HFRecHitCollection > hfrecoNoise
int iEvent
Definition: GenABIO.cc:230
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
bool getSpigotPresent(unsigned int nspigot) const
Read the &quot;PRESENT&quot; bit for this spigot.
MonitorElement * hiDistrMB2Pl2D_
bool isValid() const
Definition: HandleBase.h:76
MonitorElement * hiDistrNoise2Min2D_
int k[5][pyjets_maxn]
const_iterator end() const
int iphi() const
get the cell iphi
Definition: HcalDetId.h:38
Definition: DetId.h:18
MonitorElement * hiDistrMBMin2D_
T const * product() const
Definition: Handle.h:81
return(e1-e2)*(e1-e2)+dp *dp
edm::EDGetTokenT< HFRecHitCollection > hfrecoMB
static const int SPIGOT_COUNT
Definition: HcalDCCHeader.h:19
size_type size() const
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
MonitorElement * hiDistrNoise2Pl2D_
MonitorElement * hiDistrNoiseMin2D_
volatile std::atomic< bool > shutdown_flag false
edm::EDGetTokenT< HBHERecHitCollection > hbherecoMB
object to monitor
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const_iterator begin() const
MonitorElement * hiDistrHBHEsize1D_
MonitorElement * hiDistrMB2Min2D_
void DQMHcalPhiSymAlCaReco::beginJob ( void  )
protectedvirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 89 of file DQMHcalPhiSymAlCaReco.cc.

References DQMStore::book1D(), DQMStore::book2D(), dbe_, eventCounter_, folderName_, hFEDsize, hHcalIsZS, hiDistr_r_nbin_, hiDistr_x_max_, hiDistr_x_min_, hiDistr_x_nbin_, hiDistr_y_max_, hiDistr_y_min_, hiDistr_y_nbin_, hiDistrHBHEsize1D_, hiDistrHFsize1D_, hiDistrMB2Min2D_, hiDistrMB2Pl2D_, hiDistrMBMin2D_, hiDistrMBPl2D_, hiDistrNoise2Min2D_, hiDistrNoise2Pl2D_, hiDistrNoiseMin2D_, hiDistrNoisePl2D_, hiDistrVarMBMin2D_, hiDistrVarMBPl2D_, hiDistrVarNoiseMin2D_, hiDistrVarNoisePl2D_, hL1Id, period_, MonitorElement::setAxisTitle(), MonitorElement::setBinLabel(), DQMStore::setCurrentFolder(), SiStripMonitorClusterAlca_cfi::xmax, and SiStripMonitorClusterAlca_cfi::xmin.

89  {
90 
91  // create and cd into new folder
93 
94  eventCounter_ = 0;
95 
96  hFEDsize = dbe_->book1D("hFEDsize","HCAL FED size (kB)",200,-0.5,20.5);
97  hFEDsize->setAxisTitle("kB",1);
98 
99  hHcalIsZS = dbe_->book1D("hHcalIsZS", "Hcal Is ZS", 4, -1.5, 2.5);
100  hHcalIsZS->setBinLabel(2, "NZS");
101  hHcalIsZS->setBinLabel(3, "ZS");
102 
103  char hname[50];
104  sprintf(hname, "L1 Event Number %% %i", period_);
105  hL1Id = dbe_->book1D("hL1Id", hname,4200,-99.5,4099.5);
106  hL1Id->setAxisTitle(hname);
107 
108 
109  // book some histograms 1D
110  double xmin = 0.1;
111  double xmax = 1.1;
113  dbe_->book1D("DistrHBHEsize","Size of HBHE Collection",
115  xmin,
116  xmax
117  );
119  dbe_->book1D("DistrHFsize","Size of HF Collection",
121  xmin,
122  xmax
123  );
124 
125 
126  // First moment
127  hiDistrMBPl2D_ =
128  dbe_->book2D("MBdepthPl1", "iphi- +ieta signal distribution at depth1",
135  );
136 
137  hiDistrMBPl2D_->setAxisTitle("i#phi ", 2);
138  hiDistrMBPl2D_->setAxisTitle("i#eta ", 1);
139 
140 
142  dbe_->book2D("NoisedepthPl1", "iphi-ieta noise distribution at depth1",
143  hiDistr_x_nbin_+1,
144  hiDistr_x_min_-1.,
146  hiDistr_y_nbin_+1,
147  hiDistr_y_min_-1.,
149  );
150 
151  hiDistrNoisePl2D_->setAxisTitle("i#phi ", 2);
152  hiDistrNoisePl2D_->setAxisTitle("i#eta ", 1);
153 // Second moment
155  dbe_->book2D("MB2depthPl1", "iphi- +ieta signal distribution at depth1",
162  );
163 
164  hiDistrMB2Pl2D_->setAxisTitle("i#phi ", 2);
165  hiDistrMB2Pl2D_->setAxisTitle("i#eta ", 1);
166 
167 
169  dbe_->book2D("Noise2depthPl1", "iphi-ieta noise distribution at depth1",
176  );
177 
178  hiDistrNoise2Pl2D_->setAxisTitle("i#phi ", 2);
179  hiDistrNoise2Pl2D_->setAxisTitle("i#eta ", 1);
180 
181 // Variance
183  dbe_->book2D("VarMBdepthPl1", "iphi- +ieta signal distribution at depth1",
190  );
191 
192  hiDistrVarMBPl2D_->setAxisTitle("i#phi ", 2);
193  hiDistrVarMBPl2D_->setAxisTitle("i#eta ", 1);
194 
195 
197  dbe_->book2D("VarNoisedepthPl1", "iphi-ieta noise distribution at depth1",
204  );
205 
206  hiDistrVarNoisePl2D_->setAxisTitle("i#phi ", 2);
207  hiDistrVarNoisePl2D_->setAxisTitle("i#eta ", 1);
208 
209 //==================================================================================
210 // First moment
211  hiDistrMBMin2D_ =
212  dbe_->book2D("MBdepthMin1", "iphi- +ieta signal distribution at depth1",
219  );
220 
221  hiDistrMBMin2D_->setAxisTitle("i#phi ", 2);
222  hiDistrMBMin2D_->setAxisTitle("i#eta ", 1);
223 
224 
226  dbe_->book2D("NoisedepthMin1", "iphi-ieta noise distribution at depth1",
233  );
234 
235  hiDistrNoiseMin2D_->setAxisTitle("i#phi ", 2);
236  hiDistrNoiseMin2D_->setAxisTitle("i#eta ", 1);
237 // Second moment
239  dbe_->book2D("MB2depthMin1", "iphi- +ieta signal distribution at depth1",
246  );
247 
248  hiDistrMB2Min2D_->setAxisTitle("i#phi ", 2);
249  hiDistrMB2Min2D_->setAxisTitle("i#eta ", 1);
250 
251 
253  dbe_->book2D("Noise2depthMin1", "iphi-ieta noise distribution at depth1",
260  );
261 
262  hiDistrNoise2Min2D_->setAxisTitle("i#phi ", 2);
263  hiDistrNoise2Min2D_->setAxisTitle("i#eta ", 1);
264 
265 // Variance
267  dbe_->book2D("VarMBdepthMin1", "iphi- +ieta signal distribution at depth1",
274  );
275 
276  hiDistrVarMBMin2D_->setAxisTitle("i#phi ", 2);
277  hiDistrVarMBMin2D_->setAxisTitle("i#eta ", 1);
278 
279 
281  dbe_->book2D("VarNoisedepthMin1", "iphi-ieta noise distribution at depth1",
288  );
289 
290  hiDistrVarNoiseMin2D_->setAxisTitle("i#phi ", 2);
291  hiDistrVarNoiseMin2D_->setAxisTitle("i#eta ", 1);
292 
293 
294 }
MonitorElement * hiDistrNoisePl2D_
MonitorElement * hiDistrVarMBPl2D_
MonitorElement * hiDistrVarNoiseMin2D_
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:954
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * hiDistrMBPl2D_
MonitorElement * hiDistrHFsize1D_
MonitorElement * hiDistrVarMBMin2D_
MonitorElement * hiDistrVarNoisePl2D_
std::string folderName_
DQM folder name.
MonitorElement * hiDistrMB2Pl2D_
MonitorElement * hiDistrNoise2Min2D_
MonitorElement * hiDistrMBMin2D_
MonitorElement * hiDistrNoise2Pl2D_
MonitorElement * hiDistrNoiseMin2D_
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:1082
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:667
MonitorElement * hiDistrHBHEsize1D_
MonitorElement * hiDistrMB2Min2D_
void DQMHcalPhiSymAlCaReco::beginLuminosityBlock ( const edm::LuminosityBlock lumiSeg,
const edm::EventSetup context 
)
protectedvirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 302 of file DQMHcalPhiSymAlCaReco.cc.

303  {
304 
305 }
void DQMHcalPhiSymAlCaReco::beginRun ( const edm::Run r,
const edm::EventSetup c 
)
protectedvirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 297 of file DQMHcalPhiSymAlCaReco.cc.

297  {
298 // eventCounter_ = 0;
299 }
void DQMHcalPhiSymAlCaReco::endJob ( void  )
protectedvirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 531 of file DQMHcalPhiSymAlCaReco.cc.

References dbe_, eventCounter_, fileName_, MonitorElement::getBinContent(), hiDistr_x_nbin_, hiDistr_y_nbin_, hiDistrMB2Min2D_, hiDistrMB2Pl2D_, hiDistrMBMin2D_, hiDistrMBPl2D_, hiDistrNoise2Min2D_, hiDistrNoise2Pl2D_, hiDistrNoiseMin2D_, hiDistrNoisePl2D_, hiDistrVarMBMin2D_, hiDistrVarMBPl2D_, hiDistrVarNoiseMin2D_, hiDistrVarNoisePl2D_, j, gen::k, DQMStore::save(), saveToFile_, and MonitorElement::setBinContent().

531  {
532  if (saveToFile_) {
533 
534  for(int k=0; k<=hiDistr_x_nbin_;k++)
535  {
536  for(int j=0; j<=hiDistr_y_nbin_;j++)
537  {
538 // First moment
539  float cc1=hiDistrMBPl2D_->getBinContent(k,j);
540  cc1 = cc1 * 1./eventCounter_;
542  float cc2=hiDistrNoisePl2D_->getBinContent(k,j);
543  cc2 = cc2 * 1./eventCounter_;
545  float cc3=hiDistrMBMin2D_->getBinContent(k,j);
546  cc3 = cc3 * 1./eventCounter_;
548  float cc4=hiDistrNoiseMin2D_->getBinContent(k,j);
549  cc4 = cc4 * 1./eventCounter_;
551 // Second moment
552  float cc11=hiDistrMB2Pl2D_->getBinContent(k,j);
553  cc11 = cc11 * 1./eventCounter_;
555  hiDistrVarMBPl2D_->setBinContent(k,j,cc11-cc1*cc1);
556  float cc22=hiDistrNoise2Pl2D_->getBinContent(k,j);
557  cc22 = cc22 * 1./eventCounter_;
559  hiDistrVarNoisePl2D_->setBinContent(k,j,cc22-cc2*cc2);
560  float cc33=hiDistrMB2Min2D_->getBinContent(k,j);
561  cc33 = cc33 * 1./eventCounter_;
563  hiDistrVarMBMin2D_->setBinContent(k,j,cc33-cc3*cc3);
564  float cc44=hiDistrNoise2Min2D_->getBinContent(k,j);
565  cc44 = cc44 * 1./eventCounter_;
567  hiDistrVarNoiseMin2D_->setBinContent(k,j,cc44-cc4*cc4);
568  }
569  }
570  dbe_->save(fileName_);
571  }
572 }
MonitorElement * hiDistrNoisePl2D_
bool saveToFile_
Write to file.
void setBinContent(int binx, double content)
set content of bin (1-D)
MonitorElement * hiDistrVarMBPl2D_
MonitorElement * hiDistrVarNoiseMin2D_
MonitorElement * hiDistrMBPl2D_
MonitorElement * hiDistrVarMBMin2D_
int j
Definition: DBlmapReader.cc:9
MonitorElement * hiDistrVarNoisePl2D_
MonitorElement * hiDistrMB2Pl2D_
MonitorElement * hiDistrNoise2Min2D_
int k[5][pyjets_maxn]
MonitorElement * hiDistrMBMin2D_
double getBinContent(int binx) const
get content of bin (1-D)
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, const uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE", const bool resetMEsAfterWriting=false)
Definition: DQMStore.cc:2540
MonitorElement * hiDistrNoise2Pl2D_
MonitorElement * hiDistrNoiseMin2D_
std::string fileName_
Output file name if required.
MonitorElement * hiDistrMB2Min2D_
void DQMHcalPhiSymAlCaReco::endLuminosityBlock ( const edm::LuminosityBlock lumiSeg,
const edm::EventSetup c 
)
protectedvirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 493 of file DQMHcalPhiSymAlCaReco.cc.

494  {
495 }
void DQMHcalPhiSymAlCaReco::endRun ( const edm::Run r,
const edm::EventSetup c 
)
protectedvirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 497 of file DQMHcalPhiSymAlCaReco.cc.

References eventCounter_, MonitorElement::getBinContent(), hiDistr_x_nbin_, hiDistr_y_nbin_, hiDistrMB2Min2D_, hiDistrMB2Pl2D_, hiDistrMBMin2D_, hiDistrMBPl2D_, hiDistrNoise2Min2D_, hiDistrNoise2Pl2D_, hiDistrNoiseMin2D_, hiDistrNoisePl2D_, hiDistrVarMBMin2D_, hiDistrVarMBPl2D_, hiDistrVarNoiseMin2D_, hiDistrVarNoisePl2D_, j, gen::k, and MonitorElement::setBinContent().

497  {
498 // Keep Variances
499  if(eventCounter_ > 0) {
500  for(int k=0; k<=hiDistr_x_nbin_;k++)
501  {
502  for(int j=0; j<=hiDistr_y_nbin_;j++)
503  {
504 // First moment
505  float cc1=hiDistrMBPl2D_->getBinContent(k,j);
506  cc1 = cc1 * 1./eventCounter_;
507  float cc2=hiDistrNoisePl2D_->getBinContent(k,j);
508  cc2 = cc2 * 1./eventCounter_;
509  float cc3=hiDistrMBMin2D_->getBinContent(k,j);
510  cc3 = cc3 * 1./eventCounter_;
511  float cc4=hiDistrNoiseMin2D_->getBinContent(k,j);
512  cc4 = cc4 * 1./eventCounter_;
513 // Second moment
514  float cc11=hiDistrMB2Pl2D_->getBinContent(k,j);
515  cc11 = cc11 * 1./eventCounter_;
516  hiDistrVarMBPl2D_->setBinContent(k,j,cc11-cc1*cc1);
517  float cc22=hiDistrNoise2Pl2D_->getBinContent(k,j);
518  cc22 = cc22 * 1./eventCounter_;
519  hiDistrVarNoisePl2D_->setBinContent(k,j,cc22-cc2*cc2);
520  float cc33=hiDistrMB2Min2D_->getBinContent(k,j);
521  cc33 = cc33 * 1./eventCounter_;
522  hiDistrVarMBMin2D_->setBinContent(k,j,cc33-cc3*cc3);
523  float cc44=hiDistrNoise2Min2D_->getBinContent(k,j);
524  cc44 = cc44 * 1./eventCounter_;
525  hiDistrVarNoiseMin2D_->setBinContent(k,j,cc44-cc4*cc4);
526  }
527  }
528  }
529 }
MonitorElement * hiDistrNoisePl2D_
void setBinContent(int binx, double content)
set content of bin (1-D)
MonitorElement * hiDistrVarMBPl2D_
MonitorElement * hiDistrVarNoiseMin2D_
MonitorElement * hiDistrMBPl2D_
MonitorElement * hiDistrVarMBMin2D_
int j
Definition: DBlmapReader.cc:9
MonitorElement * hiDistrVarNoisePl2D_
MonitorElement * hiDistrMB2Pl2D_
MonitorElement * hiDistrNoise2Min2D_
int k[5][pyjets_maxn]
MonitorElement * hiDistrMBMin2D_
double getBinContent(int binx) const
get content of bin (1-D)
MonitorElement * hiDistrNoise2Pl2D_
MonitorElement * hiDistrNoiseMin2D_
MonitorElement * hiDistrMB2Min2D_

Member Data Documentation

DQMStore* DQMHcalPhiSymAlCaReco::dbe_
private

Definition at line 52 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), DQMHcalPhiSymAlCaReco(), and endJob().

int DQMHcalPhiSymAlCaReco::eventCounter_
private

Definition at line 53 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

std::string DQMHcalPhiSymAlCaReco::fileName_
private

Output file name if required.

Definition at line 114 of file DQMHcalPhiSymAlCaReco.h.

Referenced by DQMHcalPhiSymAlCaReco(), and endJob().

std::string DQMHcalPhiSymAlCaReco::folderName_
private

DQM folder name.

Definition at line 105 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), and DQMHcalPhiSymAlCaReco().

edm::EDGetTokenT<HBHERecHitCollection> DQMHcalPhiSymAlCaReco::hbherecoMB
private

object to monitor

Definition at line 94 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and DQMHcalPhiSymAlCaReco().

edm::EDGetTokenT<HBHERecHitCollection> DQMHcalPhiSymAlCaReco::hbherecoNoise
private

Definition at line 98 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and DQMHcalPhiSymAlCaReco().

MonitorElement* DQMHcalPhiSymAlCaReco::hFEDsize
private

Definition at line 76 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and beginJob().

edm::EDGetTokenT<HFRecHitCollection> DQMHcalPhiSymAlCaReco::hfrecoMB
private

Definition at line 96 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and DQMHcalPhiSymAlCaReco().

edm::EDGetTokenT<HFRecHitCollection> DQMHcalPhiSymAlCaReco::hfrecoNoise
private

Definition at line 100 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and DQMHcalPhiSymAlCaReco().

MonitorElement* DQMHcalPhiSymAlCaReco::hHcalIsZS
private

Definition at line 77 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and beginJob().

int DQMHcalPhiSymAlCaReco::hiDistr_r_nbin_
private

Definition at line 87 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), and DQMHcalPhiSymAlCaReco().

double DQMHcalPhiSymAlCaReco::hiDistr_x_max_
private

Definition at line 85 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), and DQMHcalPhiSymAlCaReco().

double DQMHcalPhiSymAlCaReco::hiDistr_x_min_
private

Definition at line 84 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), and DQMHcalPhiSymAlCaReco().

int DQMHcalPhiSymAlCaReco::hiDistr_x_nbin_
private

Definition at line 81 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), DQMHcalPhiSymAlCaReco(), endJob(), and endRun().

double DQMHcalPhiSymAlCaReco::hiDistr_y_max_
private

Definition at line 83 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), and DQMHcalPhiSymAlCaReco().

double DQMHcalPhiSymAlCaReco::hiDistr_y_min_
private

Definition at line 82 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), and DQMHcalPhiSymAlCaReco().

int DQMHcalPhiSymAlCaReco::hiDistr_y_nbin_
private

Definition at line 80 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), DQMHcalPhiSymAlCaReco(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrHBHEsize1D_
private

Definition at line 73 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrHFsize1D_
private

Definition at line 74 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrMB2Min2D_
private

Definition at line 65 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrMB2Pl2D_
private

Definition at line 63 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrMBMin2D_
private

Definition at line 60 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrMBPl2D_
private

Definition at line 58 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrNoise2Min2D_
private

Definition at line 66 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrNoise2Pl2D_
private

Definition at line 64 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrNoiseMin2D_
private

Definition at line 61 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrNoisePl2D_
private

Definition at line 59 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrVarMBMin2D_
private

Definition at line 70 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrVarMBPl2D_
private

Definition at line 68 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrVarNoiseMin2D_
private

Definition at line 71 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hiDistrVarNoisePl2D_
private

Definition at line 69 of file DQMHcalPhiSymAlCaReco.h.

Referenced by beginJob(), endJob(), and endRun().

MonitorElement* DQMHcalPhiSymAlCaReco::hL1Id
private

Definition at line 78 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and beginJob().

edm::InputTag DQMHcalPhiSymAlCaReco::horecoMB
private

Definition at line 95 of file DQMHcalPhiSymAlCaReco.h.

Referenced by DQMHcalPhiSymAlCaReco().

edm::InputTag DQMHcalPhiSymAlCaReco::horecoNoise
private

Definition at line 99 of file DQMHcalPhiSymAlCaReco.h.

Referenced by DQMHcalPhiSymAlCaReco().

double DQMHcalPhiSymAlCaReco::ihbhe_size_
private

Definition at line 88 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and DQMHcalPhiSymAlCaReco().

double DQMHcalPhiSymAlCaReco::ihf_size_
private

Definition at line 89 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and DQMHcalPhiSymAlCaReco().

unsigned int DQMHcalPhiSymAlCaReco::period_
private

Definition at line 111 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), beginJob(), and DQMHcalPhiSymAlCaReco().

edm::EDGetTokenT<FEDRawDataCollection> DQMHcalPhiSymAlCaReco::rawInLabel_
private

Definition at line 102 of file DQMHcalPhiSymAlCaReco.h.

Referenced by analyze(), and DQMHcalPhiSymAlCaReco().

bool DQMHcalPhiSymAlCaReco::saveToFile_
private

Write to file.

Definition at line 108 of file DQMHcalPhiSymAlCaReco.h.

Referenced by DQMHcalPhiSymAlCaReco(), and endJob().