CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
DQMSourceEleCalib Class Reference

#include <DQMSourceEleCalib.h>

Inheritance diagram for DQMSourceEleCalib:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

 DQMSourceEleCalib (const edm::ParameterSet &)
 
 ~DQMSourceEleCalib () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 

Private Member Functions

void fillAroundBarrel (const EcalRecHitCollection *, int, int)
 fills local occupancy graphs More...
 
void fillAroundEndcap (const EcalRecHitCollection *, int, int)
 
DetId findMaxHit (const std::vector< std::pair< DetId, float >> &, const EcalRecHitCollection *, const EcalRecHitCollection *)
 find the MOX More...
 

Private Attributes

MonitorElementElectronsNumber_
 Number of electrons. More...
 
MonitorElementESCoP_
 ESCoP. More...
 
int eventCounter_
 
std::string fileName_
 Output file name if required. More...
 
std::string folderName_
 DQM folder name. More...
 
MonitorElementHitsVsAssociatedHits_
 recHits over associated recHits More...
 
MonitorElementLocalOccupancyEB_
 
MonitorElementLocalOccupancyEE_
 
MonitorElementOccupancyEB_
 Occupancy. More...
 
MonitorElementOccupancyEEM_
 
MonitorElementOccupancyEEP_
 
unsigned int prescaleFactor_
 Monitor every prescaleFactor_ events. More...
 
edm::EDGetTokenT< EcalRecHitCollectionproductMonitoredEB_
 object to monitor More...
 
edm::EDGetTokenT< EcalRecHitCollectionproductMonitoredEE_
 object to monitor More...
 
edm::EDGetTokenT< reco::GsfElectronCollectionproductMonitoredElectrons_
 electrons to monitor More...
 
MonitorElementrecHitsPerElectron_
 Number of recHits per electron. More...
 
bool saveToFile_
 Write to file. More...
 

Detailed Description

Definition at line 27 of file DQMSourceEleCalib.h.

Constructor & Destructor Documentation

DQMSourceEleCalib::DQMSourceEleCalib ( const edm::ParameterSet ps)

Definition at line 37 of file DQMSourceEleCalib.cc.

References fileName_, folderName_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), prescaleFactor_, productMonitoredEB_, productMonitoredEE_, productMonitoredElectrons_, and saveToFile_.

37  : eventCounter_(0) {
38  folderName_ = ps.getUntrackedParameter<string>("FolderName", "ALCAStreamEcalSingleEle");
39  productMonitoredEB_ = consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("AlCaStreamEBTag"));
40  productMonitoredEE_ = consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("AlCaStreamEETag"));
41 
42  saveToFile_ = ps.getUntrackedParameter<bool>("SaveToFile", false);
43  fileName_ = ps.getUntrackedParameter<string>("FileName", "MonitorAlCaEcalSingleEle.root");
44  productMonitoredElectrons_ = consumes<reco::GsfElectronCollection>(ps.getParameter<InputTag>("electronCollection"));
45  prescaleFactor_ = ps.getUntrackedParameter<unsigned int>("prescaleFactor", 1);
46 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::string fileName_
Output file name if required.
edm::EDGetTokenT< EcalRecHitCollection > productMonitoredEE_
object to monitor
bool saveToFile_
Write to file.
edm::EDGetTokenT< EcalRecHitCollection > productMonitoredEB_
object to monitor
std::string folderName_
DQM folder name.
unsigned int prescaleFactor_
Monitor every prescaleFactor_ events.
edm::EDGetTokenT< reco::GsfElectronCollection > productMonitoredElectrons_
electrons to monitor
DQMSourceEleCalib::~DQMSourceEleCalib ( )
override

Definition at line 48 of file DQMSourceEleCalib.cc.

48 {}

Member Function Documentation

void DQMSourceEleCalib::analyze ( const edm::Event e,
const edm::EventSetup c 
)
overrideprotected

Definition at line 73 of file DQMSourceEleCalib.cc.

References edm::SortedCollection< T, SORT >::begin(), DetId::det(), EcalBarrel, EcalEndcap, ElectronsNumber_, edm::SortedCollection< T, SORT >::end(), ESCoP_, eventCounter_, MonitorElement::Fill(), fillAroundBarrel(), fillAroundEndcap(), findMaxHit(), edm::Event::getByToken(), HitsVsAssociatedHits_, triggerObjects_cff::id, EBDetId::ieta(), EBDetId::iphi(), edm::HandleBase::isValid(), EEDetId::ix(), EEDetId::iy(), Max(), OccupancyEB_, OccupancyEEM_, OccupancyEEP_, edm::Handle< T >::product(), productMonitoredEB_, productMonitoredEE_, productMonitoredElectrons_, recHitsPerElectron_, edm::SortedCollection< T, SORT >::size(), DetId::subdetId(), and ecaldqm::zside().

73  {
74  // if (eventCounter_% prescaleFactor_ ) return; //FIXME
75  eventCounter_++;
76  int numberOfHits = 0;
77  int numberOfElectrons = 0;
78  int numberOfAssociatedHits = 0;
79  // reads the recHits
82 
83  iEvent.getByToken(productMonitoredEB_, rhEB);
84  iEvent.getByToken(productMonitoredEE_, rhEE);
85 
87 
88  // reads the electrons
90  iEvent.getByToken(productMonitoredElectrons_, pElectrons);
91 
92  if (pElectrons.isValid()) {
93  ElectronsNumber_->Fill(pElectrons->size() + 0.1);
94  numberOfElectrons = pElectrons->size();
95  for (reco::GsfElectronCollection::const_iterator eleIt = pElectrons->begin(); eleIt != pElectrons->end(); ++eleIt) {
96  ESCoP_->Fill(eleIt->eSuperClusterOverP());
97  numberOfAssociatedHits += eleIt->superCluster()->size();
98  DetId Max = findMaxHit(eleIt->superCluster()->hitsAndFractions(), rhEB.product(), rhEE.product());
99  if (!Max.det())
100  continue;
101  if (Max.subdetId() == EcalBarrel) {
102  EBDetId EBMax(Max);
103  fillAroundBarrel(rhEB.product(), EBMax.ieta(), EBMax.iphi());
104  }
105  if (Max.subdetId() == EcalEndcap) {
106  EEDetId EEMax(Max);
107  fillAroundEndcap(rhEE.product(), EEMax.ix(), EEMax.iy());
108  }
109  }
110  } // is valid electron
111 
112  // fill EB histos
113  if (rhEB.isValid()) {
114  numberOfHits += rhEB->size();
115  for (itb = rhEB->begin(); itb != rhEB->end(); ++itb) {
116  EBDetId id(itb->id());
117  OccupancyEB_->Fill(id.iphi(), id.ieta());
118  } // Eb rechits
119  } // is Valid
120  if (rhEE.isValid()) {
121  numberOfHits += rhEE->size();
122  for (itb = rhEE->begin(); itb != rhEE->end(); ++itb) {
123  EEDetId id(itb->id());
124  if (id.zside() > 0) {
125  OccupancyEEP_->Fill(id.ix(), id.iy());
126  } // zside>0
127  else if (id.zside() < 0) {
128  OccupancyEEM_->Fill(id.ix(), id.iy());
129  } // zside<0
130 
131  } // EE reChit
132  } // is Valid
133  if (numberOfElectrons)
134  recHitsPerElectron_->Fill((double)numberOfHits / ((double)numberOfElectrons));
135  if (numberOfHits)
136  HitsVsAssociatedHits_->Fill((double)numberOfAssociatedHits / ((double)numberOfHits));
137 } // end of the analyzer
MonitorElement * ESCoP_
ESCoP.
DetId findMaxHit(const std::vector< std::pair< DetId, float >> &, const EcalRecHitCollection *, const EcalRecHitCollection *)
find the MOX
std::vector< EcalRecHit >::const_iterator const_iterator
int zside(DetId const &)
edm::EDGetTokenT< EcalRecHitCollection > productMonitoredEE_
object to monitor
void Fill(long long x)
MonitorElement * OccupancyEEM_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< EcalRecHitCollection > productMonitoredEB_
object to monitor
MonitorElement * HitsVsAssociatedHits_
recHits over associated recHits
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
MonitorElement * OccupancyEEP_
bool isValid() const
Definition: HandleBase.h:74
void fillAroundEndcap(const EcalRecHitCollection *, int, int)
MonitorElement * recHitsPerElectron_
Number of recHits per electron.
const_iterator end() const
T Max(T a, T b)
Definition: MathUtil.h:44
void fillAroundBarrel(const EcalRecHitCollection *, int, int)
fills local occupancy graphs
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:74
MonitorElement * OccupancyEB_
Occupancy.
size_type size() const
MonitorElement * ElectronsNumber_
Number of electrons.
const_iterator begin() const
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39
edm::EDGetTokenT< reco::GsfElectronCollection > productMonitoredElectrons_
electrons to monitor
void DQMSourceEleCalib::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  irun,
edm::EventSetup const &  isetup 
)
overrideprotected

Definition at line 51 of file DQMSourceEleCalib.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), ElectronsNumber_, ESCoP_, folderName_, HitsVsAssociatedHits_, LocalOccupancyEB_, LocalOccupancyEE_, OccupancyEB_, OccupancyEEM_, OccupancyEEP_, recHitsPerElectron_, and DQMStore::IBooker::setCurrentFolder().

53  {
54  // create and cd into new folder
56 
57  recHitsPerElectron_ = ibooker.book1D("recHitsPerElectron_", "recHitPerElectron", 200, 0, 200);
58  ElectronsNumber_ = ibooker.book1D("ElectronsNumber_", "electrons in the event", 40, 0, 40);
59  ESCoP_ = ibooker.book1D("ESCoP", "ESCoP", 50, 0, 5);
60 
61  OccupancyEB_ = ibooker.book2D("OccupancyEB_", "OccupancyEB", 360, 1, 361, 171, -85, 86);
62  OccupancyEEP_ = ibooker.book2D("OccupancyEEP_", "Occupancy EE Plus", 100, 1, 101, 100, 1, 101);
63  OccupancyEEM_ = ibooker.book2D("OccupancyEEM_", "Occupancy EE Minus", 100, 1, 101, 100, 1, 101);
64  HitsVsAssociatedHits_ = ibooker.book1D("HitsVsAssociatedHits_", "HitsVsAssociatedHits", 100, 0, 5);
65  LocalOccupancyEB_ = ibooker.book2D("LocalOccupancyEB_", "Local occupancy Barrel", 9, -4, 5, 9, -4, 5);
66  LocalOccupancyEE_ = ibooker.book2D("LocalOccupancyEE_", "Local occupancy Endcap", 9, -4, 5, 9, -4, 5);
67 }
MonitorElement * ESCoP_
ESCoP.
MonitorElement * LocalOccupancyEE_
MonitorElement * OccupancyEEM_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * HitsVsAssociatedHits_
recHits over associated recHits
std::string folderName_
DQM folder name.
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * OccupancyEEP_
MonitorElement * recHitsPerElectron_
Number of recHits per electron.
MonitorElement * OccupancyEB_
Occupancy.
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
MonitorElement * ElectronsNumber_
Number of electrons.
MonitorElement * LocalOccupancyEB_
void DQMSourceEleCalib::fillAroundBarrel ( const EcalRecHitCollection recHits,
int  eta,
int  phi 
)
private

fills local occupancy graphs

Definition at line 179 of file DQMSourceEleCalib.cc.

References edm::SortedCollection< T, SORT >::begin(), HTMLExport::elem(), edm::SortedCollection< T, SORT >::end(), PVValHelper::eta, MonitorElement::Fill(), EBDetId::ieta(), EBDetId::iphi(), LocalOccupancyEB_, and phi.

Referenced by analyze().

179  {
180  for (EcalRecHitCollection::const_iterator elem = recHits->begin(); elem != recHits->end(); ++elem) {
181  EBDetId elementId = elem->id();
182  LocalOccupancyEB_->Fill(elementId.ieta() - eta, elementId.iphi() - phi, elem->energy());
183  }
184  return;
185 }
std::vector< EcalRecHit >::const_iterator const_iterator
void Fill(long long x)
int iphi() const
get the crystal iphi
Definition: EBDetId.h:51
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
def elem(elemtype, innerHTML='', html_class='', kwargs)
Definition: HTMLExport.py:19
const_iterator end() const
const_iterator begin() const
MonitorElement * LocalOccupancyEB_
void DQMSourceEleCalib::fillAroundEndcap ( const EcalRecHitCollection recHits,
int  ics,
int  ips 
)
private

Definition at line 189 of file DQMSourceEleCalib.cc.

References edm::SortedCollection< T, SORT >::begin(), HTMLExport::elem(), edm::SortedCollection< T, SORT >::end(), MonitorElement::Fill(), EEDetId::ix(), EEDetId::iy(), and LocalOccupancyEE_.

Referenced by analyze().

189  {
190  for (EcalRecHitCollection::const_iterator elem = recHits->begin(); elem != recHits->end(); ++elem) {
191  EEDetId elementId = elem->id();
192  LocalOccupancyEE_->Fill(elementId.ix() - ics, elementId.iy() - ips, elem->energy());
193  }
194  return;
195 }
int ix() const
Definition: EEDetId.h:77
MonitorElement * LocalOccupancyEE_
std::vector< EcalRecHit >::const_iterator const_iterator
void Fill(long long x)
int iy() const
Definition: EEDetId.h:83
def elem(elemtype, innerHTML='', html_class='', kwargs)
Definition: HTMLExport.py:19
const_iterator end() const
const_iterator begin() const
DetId DQMSourceEleCalib::findMaxHit ( const std::vector< std::pair< DetId, float >> &  v1,
const EcalRecHitCollection EBhits,
const EcalRecHitCollection EEhits 
)
private

find the MOX

Definition at line 143 of file DQMSourceEleCalib.cc.

References EcalBarrel, edm::SortedCollection< T, SORT >::end(), edm::SortedCollection< T, SORT >::find(), and MTVHistoProducerAlgoForTrackerBlock_cfi::maxHit.

Referenced by analyze().

145  {
146  double currEnergy = 0.;
147  DetId maxHit;
148  for (std::vector<std::pair<DetId, float>>::const_iterator idsIt = v1.begin(); idsIt != v1.end(); ++idsIt) {
149  if (idsIt->first.subdetId() == EcalBarrel) {
151  itrechit = EBhits->find((*idsIt).first);
152  if (itrechit == EBhits->end()) {
153  edm::LogInfo("reading") << "[findMaxHit] rechit not found! ";
154  continue;
155  }
156  // FIXME: wnat to use the fraction i.e. .second??
157  if (itrechit->energy() > currEnergy) {
158  currEnergy = itrechit->energy();
159  maxHit = (*idsIt).first;
160  }
161  } else {
163  itrechit = EEhits->find((*idsIt).first);
164  if (itrechit == EEhits->end()) {
165  edm::LogInfo("reading") << "[findMaxHit] rechit not found! ";
166  continue;
167  }
168 
169  // FIXME: wnat to use the fraction i.e. .second??
170  if (itrechit->energy() > currEnergy) {
171  currEnergy = itrechit->energy();
172  maxHit = (*idsIt).first;
173  }
174  }
175  }
176  return maxHit;
177 }
std::vector< EcalRecHit >::const_iterator const_iterator
const_iterator end() const
Definition: DetId.h:18
iterator find(key_type k)

Member Data Documentation

MonitorElement* DQMSourceEleCalib::ElectronsNumber_
private

Number of electrons.

Definition at line 51 of file DQMSourceEleCalib.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DQMSourceEleCalib::ESCoP_
private

ESCoP.

Definition at line 53 of file DQMSourceEleCalib.h.

Referenced by analyze(), and bookHistograms().

int DQMSourceEleCalib::eventCounter_
private

Definition at line 46 of file DQMSourceEleCalib.h.

Referenced by analyze().

std::string DQMSourceEleCalib::fileName_
private

Output file name if required.

Definition at line 82 of file DQMSourceEleCalib.h.

Referenced by DQMSourceEleCalib().

std::string DQMSourceEleCalib::folderName_
private

DQM folder name.

Definition at line 76 of file DQMSourceEleCalib.h.

Referenced by bookHistograms(), and DQMSourceEleCalib().

MonitorElement* DQMSourceEleCalib::HitsVsAssociatedHits_
private

recHits over associated recHits

Definition at line 62 of file DQMSourceEleCalib.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DQMSourceEleCalib::LocalOccupancyEB_
private

Definition at line 58 of file DQMSourceEleCalib.h.

Referenced by bookHistograms(), and fillAroundBarrel().

MonitorElement* DQMSourceEleCalib::LocalOccupancyEE_
private

Definition at line 59 of file DQMSourceEleCalib.h.

Referenced by bookHistograms(), and fillAroundEndcap().

MonitorElement* DQMSourceEleCalib::OccupancyEB_
private

Occupancy.

Definition at line 55 of file DQMSourceEleCalib.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DQMSourceEleCalib::OccupancyEEM_
private

Definition at line 57 of file DQMSourceEleCalib.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DQMSourceEleCalib::OccupancyEEP_
private

Definition at line 56 of file DQMSourceEleCalib.h.

Referenced by analyze(), and bookHistograms().

unsigned int DQMSourceEleCalib::prescaleFactor_
private

Monitor every prescaleFactor_ events.

Definition at line 73 of file DQMSourceEleCalib.h.

Referenced by DQMSourceEleCalib().

edm::EDGetTokenT<EcalRecHitCollection> DQMSourceEleCalib::productMonitoredEB_
private

object to monitor

Definition at line 65 of file DQMSourceEleCalib.h.

Referenced by analyze(), and DQMSourceEleCalib().

edm::EDGetTokenT<EcalRecHitCollection> DQMSourceEleCalib::productMonitoredEE_
private

object to monitor

Definition at line 68 of file DQMSourceEleCalib.h.

Referenced by analyze(), and DQMSourceEleCalib().

edm::EDGetTokenT<reco::GsfElectronCollection> DQMSourceEleCalib::productMonitoredElectrons_
private

electrons to monitor

Definition at line 70 of file DQMSourceEleCalib.h.

Referenced by analyze(), and DQMSourceEleCalib().

MonitorElement* DQMSourceEleCalib::recHitsPerElectron_
private

Number of recHits per electron.

Definition at line 49 of file DQMSourceEleCalib.h.

Referenced by analyze(), and bookHistograms().

bool DQMSourceEleCalib::saveToFile_
private

Write to file.

Definition at line 79 of file DQMSourceEleCalib.h.

Referenced by DQMSourceEleCalib().