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
HcalSimHitStudy Class Reference

#include <HcalSimHitStudy.h>

Inheritance diagram for HcalSimHitStudy:
DQMEDAnalyzer edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > > edm::stream::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

virtual void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
 
 HcalSimHitStudy (const edm::ParameterSet &ps)
 
 ~HcalSimHitStudy ()
 
- Public Member Functions inherited from DQMEDAnalyzer
virtual void beginRun (edm::Run const &, edm::EventSetup const &) final
 
virtual void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer (void)
 
virtual void endLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, dqmDetails::NoCache *) const final
 
virtual void endRunSummary (edm::Run const &, edm::EventSetup const &, dqmDetails::NoCache *) const final
 
uint32_t streamId () const
 
- Public Member Functions inherited from edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > >
 EDAnalyzer ()=default
 
- Public Member Functions inherited from edm::stream::EDAnalyzerBase
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDAnalyzerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 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
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) 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 analyzeHits (std::vector< PCaloHit > &)
 
- Protected Member Functions inherited from edm::stream::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- 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

bool checkHit_
 
std::string g4Label
 
std::string hcalHits
 
MonitorElementmeAllNHit_
 
MonitorElementmeBadDetHit_
 
MonitorElementmeBadIdHit_
 
MonitorElementmeBadSubHit_
 
MonitorElementmeDepthHit_
 
MonitorElementmeDetectHit_
 
MonitorElementmeEnergyHit_
 
MonitorElementmeEtaHit_
 
MonitorElementmeHBDepHit_
 
MonitorElementmeHBEneHit2_
 
MonitorElementmeHBEneHit_
 
MonitorElementmeHBEtaHit_
 
MonitorElementmeHBL10Ene_
 
MonitorElementmeHBL10EneP_
 
MonitorElementmeHBNHit_
 
MonitorElementmeHBPhiHit_
 
MonitorElementmeHBTimHit_
 
MonitorElementmeHEDepHit_
 
MonitorElementmeHEEneHit2_
 
MonitorElementmeHEEneHit_
 
MonitorElementmeHEEtaHit_
 
MonitorElementmeHEL10Ene_
 
MonitorElementmeHEL10EneP_
 
MonitorElementmeHENHit_
 
MonitorElementmeHEPhiHit_
 
MonitorElementmeHETimHit_
 
MonitorElementmeHFDepHit_
 
MonitorElementmeHFEneHit2_
 
MonitorElementmeHFEneHit_
 
MonitorElementmeHFEtaHit_
 
MonitorElementmeHFL10Ene_
 
MonitorElementmeHFL10EneP_
 
MonitorElementmeHFNHit_
 
MonitorElementmeHFPhiHit_
 
MonitorElementmeHFTimHit_
 
MonitorElementmeHODepHit_
 
MonitorElementmeHOEneHit2_
 
MonitorElementmeHOEneHit_
 
MonitorElementmeHOEtaHit_
 
MonitorElementmeHOL10Ene_
 
MonitorElementmeHOL10EneP_
 
MonitorElementmeHONHit_
 
MonitorElementmeHOPhiHit_
 
MonitorElementmeHOTimHit_
 
MonitorElementmePhiHit_
 
MonitorElementmeSubdetHit_
 
MonitorElementmeTimeHit_
 
MonitorElementmeTimeWHit_
 
std::string outFile_
 
edm::EDGetTokenT
< edm::PCaloHitContainer
tok_hits_
 
bool verbose_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > >
typedef CacheContexts< T...> CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T...> HasAbility
 
typedef
CacheTypes::LuminosityBlockCache 
LuminosityBlockCache
 
typedef
LuminosityBlockContextT
< LuminosityBlockCache,
RunCache, GlobalCache
LuminosityBlockContext
 
typedef
CacheTypes::LuminosityBlockSummaryCache 
LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache,
GlobalCache
RunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDAnalyzerBase
typedef EDAnalyzerAdaptorBase ModuleType
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static std::shared_ptr
< dqmDetails::NoCache
globalBeginLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, LuminosityBlockContext const *)
 
static std::shared_ptr
< dqmDetails::NoCache
globalBeginRunSummary (edm::Run const &, edm::EventSetup const &, RunContext const *)
 
static void globalEndLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, LuminosityBlockContext const *, dqmDetails::NoCache *)
 
static void globalEndRunSummary (edm::Run const &, edm::EventSetup const &, RunContext const *, dqmDetails::NoCache *)
 
- Static Public Member Functions inherited from edm::stream::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Detailed Description

Definition at line 28 of file HcalSimHitStudy.h.

Constructor & Destructor Documentation

HcalSimHitStudy::HcalSimHitStudy ( const edm::ParameterSet ps)

Definition at line 6 of file HcalSimHitStudy.cc.

References checkHit_, g4Label, edm::ParameterSet::getUntrackedParameter(), hcalHits, HLT_25ns14e33_v1_cff::InputTag, outFile_, AlCaHLTBitMon_QueryRunRegistry::string, tok_hits_, and verbose_.

6  {
7 
8  g4Label = ps.getUntrackedParameter<std::string>("moduleLabel","g4SimHits");
9  hcalHits = ps.getUntrackedParameter<std::string>("HitCollection","HcalHits");
10  outFile_ = ps.getUntrackedParameter<std::string>("outputFile", "hcHit.root");
11  verbose_ = ps.getUntrackedParameter<bool>("Verbose", false);
12  checkHit_= true;
13 
14  tok_hits_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label,hcalHits));
15 
16  edm::LogInfo("HcalSim") << "Module Label: " << g4Label << " Hits: "
17  << hcalHits << " / "<< checkHit_
18  << " Output: " << outFile_;
19 
20 }
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
T getUntrackedParameter(std::string const &, T const &) const
std::string hcalHits
std::string outFile_
std::string g4Label
HcalSimHitStudy::~HcalSimHitStudy ( )

Definition at line 22 of file HcalSimHitStudy.cc.

22 {}

Member Function Documentation

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

Implements edm::stream::EDAnalyzerBase.

Definition at line 87 of file HcalSimHitStudy.cc.

References analyzeHits(), checkHit_, edm::EventID::event(), edm::Event::getByToken(), edm::EventBase::id(), edm::HandleBase::isValid(), LogDebug, edm::EventID::run(), and tok_hits_.

87  {
88 
89  LogDebug("HcalSim") << "Run = " << e.id().run() << " Event = "
90  << e.id().event();
91 
92  std::vector<PCaloHit> caloHits;
94 
95  bool getHits = false;
96  if (checkHit_) {
97  e.getByToken(tok_hits_,hitsHcal);
98  if (hitsHcal.isValid()) getHits = true;
99  }
100 
101  LogDebug("HcalSim") << "HcalValidation: Input flags Hits " << getHits;
102 
103  if (getHits) {
104  caloHits.insert(caloHits.end(),hitsHcal->begin(),hitsHcal->end());
105  LogDebug("HcalSim") << "HcalValidation: Hit buffer "
106  << caloHits.size();
107  analyzeHits (caloHits);
108  }
109 }
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
EventNumber_t event() const
Definition: EventID.h:41
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
void analyzeHits(std::vector< PCaloHit > &)
bool isValid() const
Definition: HandleBase.h:75
edm::EventID id() const
Definition: EventBase.h:60
void HcalSimHitStudy::analyzeHits ( std::vector< PCaloHit > &  hits)
protected

Definition at line 111 of file HcalSimHitStudy.cc.

References TauDecayModes::dec, HLT_25ns14e33_v1_cff::depth, relval_parameters_module::energy, eta(), MonitorElement::Fill(), HcalBarrel, HcalEndcap, HcalForward, HcalOuter, i, LogDebug, meAllNHit_, meBadDetHit_, meBadIdHit_, meBadSubHit_, meDepthHit_, meDetectHit_, meEnergyHit_, meEtaHit_, meHBDepHit_, meHBEneHit2_, meHBEneHit_, meHBEtaHit_, meHBL10Ene_, meHBL10EneP_, meHBNHit_, meHBPhiHit_, meHBTimHit_, meHEDepHit_, meHEEneHit2_, meHEEneHit_, meHEEtaHit_, meHEL10Ene_, meHEL10EneP_, meHENHit_, meHEPhiHit_, meHETimHit_, meHFDepHit_, meHFEneHit2_, meHFEneHit_, meHFEtaHit_, meHFL10Ene_, meHFL10EneP_, meHFNHit_, meHFPhiHit_, meHFTimHit_, meHODepHit_, meHOEneHit2_, meHOEneHit_, meHOEtaHit_, meHOL10Ene_, meHOL10EneP_, meHONHit_, meHOPhiHit_, meHOTimHit_, mePhiHit_, meSubdetHit_, meTimeHit_, meTimeWHit_, phi, and cond::rpcobgas::time.

Referenced by analyze().

111  {
112 
113  int nHit = hits.size();
114  int nHB=0, nHE=0, nHO=0, nHF=0, nBad1=0, nBad2=0, nBad=0;
115  std::vector<double> encontHB(140, 0.);
116  std::vector<double> encontHE(140, 0.);
117  std::vector<double> encontHF(140, 0.);
118  std::vector<double> encontHO(140, 0.);
119  double entotHB = 0, entotHE = 0, entotHF = 0, entotHO = 0;
120 
121  for (int i=0; i<nHit; i++) {
122  double energy = hits[i].energy();
123  double log10en = log10(energy);
124  int log10i = int( (log10en+10.)*10. );
125  double time = hits[i].time();
126  unsigned int id_ = hits[i].id();
127  HcalDetId id = HcalDetId(id_);
128  int det = id.det();
129  int subdet = id.subdet();
130  int depth = id.depth();
131  int eta = id.ieta();
132  int phi = id.iphi();
133  LogDebug("HcalSim") << "Hit[" << i << "] ID " << std::hex << id_
134  << std::dec << " Det " << det << " Sub "
135  << subdet << " depth " << depth << " Eta " << eta
136  << " Phi " << phi << " E " << energy << " time "
137  << time;
138  if (det == 4) { // Check DetId.h
139  if (subdet == static_cast<int>(HcalBarrel)) nHB++;
140  else if (subdet == static_cast<int>(HcalEndcap)) nHE++;
141  else if (subdet == static_cast<int>(HcalOuter)) nHO++;
142  else if (subdet == static_cast<int>(HcalForward)) nHF++;
143  else { nBad++; nBad2++;}
144  } else { nBad++; nBad1++;}
145 
146  meDetectHit_->Fill(double(det));
147  if (det == 4) {
148  meSubdetHit_->Fill(double(subdet));
149  meDepthHit_->Fill(double(depth));
150  meEtaHit_->Fill(double(eta));
151  mePhiHit_->Fill(double(phi));
152  meEnergyHit_->Fill(energy);
153  meTimeHit_->Fill(time);
154  meTimeWHit_->Fill(double(time),energy);
155  if (subdet == static_cast<int>(HcalBarrel)) {
156  meHBDepHit_->Fill(double(depth));
157  meHBEtaHit_->Fill(double(eta));
158  meHBPhiHit_->Fill(double(phi));
159  meHBEneHit_->Fill(energy);
160  meHBEneHit2_->Fill(energy);
161  meHBTimHit_->Fill(time);
162  meHBL10Ene_->Fill(log10en);
163  if( log10i >=0 && log10i < 140 ) encontHB[log10i] += energy;
164  entotHB += energy;
165  } else if (subdet == static_cast<int>(HcalEndcap)) {
166  meHEDepHit_->Fill(double(depth));
167  meHEEtaHit_->Fill(double(eta));
168  meHEPhiHit_->Fill(double(phi));
169  meHEEneHit_->Fill(energy);
170  meHEEneHit2_->Fill(energy);
171  meHETimHit_->Fill(time);
172  meHEL10Ene_->Fill(log10en);
173  if( log10i >=0 && log10i < 140 ) encontHE[log10i] += energy;
174  entotHE += energy;
175  } else if (subdet == static_cast<int>(HcalOuter)) {
176  meHODepHit_->Fill(double(depth));
177  meHOEtaHit_->Fill(double(eta));
178  meHOPhiHit_->Fill(double(phi));
179  meHOEneHit_->Fill(energy);
180  meHOEneHit2_->Fill(energy);
181  meHOTimHit_->Fill(time);
182  meHOL10Ene_->Fill(log10en);
183  if( log10i >=0 && log10i < 140 ) encontHO[log10i] += energy;
184  entotHO += energy;
185  } else if (subdet == static_cast<int>(HcalForward)) {
186  meHFDepHit_->Fill(double(depth));
187  meHFEtaHit_->Fill(double(eta));
188  meHFPhiHit_->Fill(double(phi));
189  meHFEneHit_->Fill(energy);
190  meHFEneHit2_->Fill(energy);
191  meHFTimHit_->Fill(time);
192  meHFL10Ene_->Fill(log10en);
193  if( log10i >=0 && log10i < 140 ) encontHF[log10i] += energy;
194  entotHF += energy;
195  }
196  }
197  }
198  if( entotHB != 0 ) for( int i=0; i<140; i++ ) meHBL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHB[i]/entotHB );
199  if( entotHE != 0 ) for( int i=0; i<140; i++ ) meHEL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHE[i]/entotHE );
200  if( entotHF != 0 ) for( int i=0; i<140; i++ ) meHFL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHF[i]/entotHF );
201  if( entotHO != 0 ) for( int i=0; i<140; i++ ) meHOL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHO[i]/entotHO );
202 
203  meAllNHit_->Fill(double(nHit));
204  meBadDetHit_->Fill(double(nBad1));
205  meBadSubHit_->Fill(double(nBad2));
206  meBadIdHit_->Fill(double(nBad));
207  meHBNHit_->Fill(double(nHB));
208  meHENHit_->Fill(double(nHE));
209  meHONHit_->Fill(double(nHO));
210  meHFNHit_->Fill(double(nHF));
211 
212  LogDebug("HcalSim") << "HcalSimHitStudy::analyzeHits: HB " << nHB
213  << " HE " << nHE << " HO " << nHO << " HF " << nHF
214  << " Bad " << nBad << " All " << nHit;
215 
216 }
#define LogDebug(id)
MonitorElement * mePhiHit_
int i
Definition: DBlmapReader.cc:9
MonitorElement * meHEDepHit_
MonitorElement * meHOL10Ene_
MonitorElement * meHOPhiHit_
MonitorElement * meDetectHit_
MonitorElement * meHFNHit_
MonitorElement * meHFEtaHit_
MonitorElement * meHFL10EneP_
MonitorElement * meHFEneHit2_
MonitorElement * meHOL10EneP_
MonitorElement * meHBNHit_
MonitorElement * meHFTimHit_
T eta() const
MonitorElement * meHBEneHit2_
MonitorElement * meHONHit_
void Fill(long long x)
MonitorElement * meAllNHit_
MonitorElement * meHBL10Ene_
MonitorElement * meBadDetHit_
MonitorElement * meHFL10Ene_
MonitorElement * meHETimHit_
MonitorElement * meTimeHit_
MonitorElement * meHFDepHit_
MonitorElement * meHBTimHit_
MonitorElement * meTimeWHit_
MonitorElement * meHBDepHit_
MonitorElement * meBadSubHit_
MonitorElement * meHENHit_
MonitorElement * meBadIdHit_
MonitorElement * meHEL10EneP_
MonitorElement * meEnergyHit_
MonitorElement * meHOEtaHit_
MonitorElement * meEtaHit_
MonitorElement * meHBEneHit_
MonitorElement * meHOTimHit_
MonitorElement * meHBL10EneP_
MonitorElement * meHEEneHit2_
MonitorElement * meHFPhiHit_
MonitorElement * meDepthHit_
MonitorElement * meHBEtaHit_
MonitorElement * meHODepHit_
MonitorElement * meHEL10Ene_
MonitorElement * meHEEneHit_
MonitorElement * meHEPhiHit_
MonitorElement * meSubdetHit_
MonitorElement * meHBPhiHit_
MonitorElement * meHFEneHit_
MonitorElement * meHEEtaHit_
MonitorElement * meHOEneHit_
MonitorElement * meHOEneHit2_
Definition: DDAxes.h:10
void HcalSimHitStudy::bookHistograms ( DQMStore::IBooker ib,
edm::Run const &  run,
edm::EventSetup const &  es 
)
virtual

Implements DQMEDAnalyzer.

Definition at line 24 of file HcalSimHitStudy.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::bookProfile(), checkHit_, meAllNHit_, meBadDetHit_, meBadIdHit_, meBadSubHit_, meDepthHit_, meDetectHit_, meEnergyHit_, meEtaHit_, meHBDepHit_, meHBEneHit2_, meHBEneHit_, meHBEtaHit_, meHBL10Ene_, meHBL10EneP_, meHBNHit_, meHBPhiHit_, meHBTimHit_, meHEDepHit_, meHEEneHit2_, meHEEneHit_, meHEEtaHit_, meHEL10Ene_, meHEL10EneP_, meHENHit_, meHEPhiHit_, meHETimHit_, meHFDepHit_, meHFEneHit2_, meHFEneHit_, meHFEtaHit_, meHFL10Ene_, meHFL10EneP_, meHFNHit_, meHFPhiHit_, meHFTimHit_, meHODepHit_, meHOEneHit2_, meHOEneHit_, meHOEtaHit_, meHOL10Ene_, meHOL10EneP_, meHONHit_, meHOPhiHit_, meHOTimHit_, mePhiHit_, meSubdetHit_, meTimeHit_, meTimeWHit_, and DQMStore::IBooker::setCurrentFolder().

25 {
26  ib.setCurrentFolder("HcalHitsV/HcalSimHitsTask");
27 
28  //Histograms for Hits
29  if (checkHit_) {
30  meAllNHit_ = ib.book1D("Hit01","Number of Hits in HCal",1000,0.,1000.);
31  meBadDetHit_= ib.book1D("Hit02","Hits with wrong Det", 100,0.,100.);
32  meBadSubHit_= ib.book1D("Hit03","Hits with wrong Subdet",100,0.,100.);
33  meBadIdHit_ = ib.book1D("Hit04","Hits with wrong ID", 100,0.,100.);
34  meHBNHit_ = ib.book1D("Hit05","Number of Hits in HB",1000,0.,1000.);
35  meHENHit_ = ib.book1D("Hit06","Number of Hits in HE",1000,0.,1000.);
36  meHONHit_ = ib.book1D("Hit07","Number of Hits in HO",1000,0.,1000.);
37  meHFNHit_ = ib.book1D("Hit08","Number of Hits in HF",1000,0.,1000.);
38  meDetectHit_= ib.book1D("Hit09","Detector ID", 50,0.,50.);
39  meSubdetHit_= ib.book1D("Hit10","Subdetectors in HCal", 50,0.,50.);
40  meDepthHit_ = ib.book1D("Hit11","Depths in HCal", 20,0.,20.);
41  meEtaHit_ = ib.book1D("Hit12","Eta in HCal", 100,-50.,50.);
42  mePhiHit_ = ib.book1D("Hit13","Phi in HCal", 100,0.,100.);
43  meEnergyHit_= ib.book1D("Hit14","Energy in HCal", 100,0.,1.);
44  meTimeHit_ = ib.book1D("Hit15","Time in HCal", 100,0.,400.);
45  meTimeWHit_ = ib.book1D("Hit16","Time in HCal (E wtd)", 100,0.,400.);
46  meHBDepHit_ = ib.book1D("Hit17","Depths in HB", 20,0.,20.);
47  meHEDepHit_ = ib.book1D("Hit18","Depths in HE", 20,0.,20.);
48  meHODepHit_ = ib.book1D("Hit19","Depths in HO", 20,0.,20.);
49  meHFDepHit_ = ib.book1D("Hit20","Depths in HF", 20,0.,20.);
50  meHBEtaHit_ = ib.book1D("Hit21","Eta in HB", 100,-50.,50.);
51  meHEEtaHit_ = ib.book1D("Hit22","Eta in HE", 100,-50.,50.);
52  meHOEtaHit_ = ib.book1D("Hit23","Eta in HO", 100,-50.,50.);
53  meHFEtaHit_ = ib.book1D("Hit24","Eta in HF", 100,-50.,50.);
54  meHBPhiHit_ = ib.book1D("Hit25","Phi in HB", 100,0.,100.);
55  meHEPhiHit_ = ib.book1D("Hit26","Phi in HE", 100,0.,100.);
56  meHOPhiHit_ = ib.book1D("Hit27","Phi in HO", 100,0.,100.);
57  meHFPhiHit_ = ib.book1D("Hit28","Phi in HF", 100,0.,100.);
58  meHBEneHit_ = ib.book1D("Hit29","Energy in HB", 100,0.,1.);
59  meHEEneHit_ = ib.book1D("Hit30","Energy in HE", 100,0.,1.);
60  meHOEneHit_ = ib.book1D("Hit31","Energy in HO", 100,0.,1.);
61  meHFEneHit_ = ib.book1D("Hit32","Energy in HF", 100,0.,100.);
62  meHBTimHit_ = ib.book1D("Hit33","Time in HB", 100,0.,400.);
63  meHETimHit_ = ib.book1D("Hit34","Time in HE", 100,0.,400.);
64  meHOTimHit_ = ib.book1D("Hit35","Time in HO", 100,0.,400.);
65  meHFTimHit_ = ib.book1D("Hit36","Time in HF", 100,0.,400.);
66  meHBEneHit2_ = ib.book1D("Hit37","Energy in HB 2", 100,0.,0.0001);
67  meHEEneHit2_ = ib.book1D("Hit38","Energy in HE 2", 100,0.,0.0001);
68  meHOEneHit2_ = ib.book1D("Hit39","Energy in HO 2", 100,0.,0.0001);
69  meHFEneHit2_ = ib.book1D("Hit40","Energy in HF 2", 100,0.,0.0001);
70  meHBL10Ene_ = ib.book1D("Hit41","Log10Energy in HB", 140, -10., 4. );
71  meHEL10Ene_ = ib.book1D("Hit42","Log10Energy in HE", 140, -10., 4. );
72  meHFL10Ene_ = ib.book1D("Hit43","Log10Energy in HF", 140, -10., 4. );
73  meHOL10Ene_ = ib.book1D("Hit44","Log10Energy in HO", 140, -10., 4. );
74  meHBL10EneP_ = ib.bookProfile("Hit45","Log10Energy in HB vs Hit contribution", 140, -10., 4., 100, 0., 1. );
75  meHEL10EneP_ = ib.bookProfile("Hit46","Log10Energy in HE vs Hit contribution", 140, -10., 4., 100, 0., 1. );
76  meHFL10EneP_ = ib.bookProfile("Hit47","Log10Energy in HF vs Hit contribution", 140, -10., 4., 100, 0., 1. );
77  meHOL10EneP_ = ib.bookProfile("Hit48","Log10Energy in HO vs Hit contribution", 140, -10., 4., 100, 0., 1. );
78  }
79 
80 }
MonitorElement * mePhiHit_
MonitorElement * meHEDepHit_
MonitorElement * meHOL10Ene_
MonitorElement * meHOPhiHit_
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
MonitorElement * meDetectHit_
MonitorElement * meHFNHit_
MonitorElement * meHFEtaHit_
MonitorElement * meHFL10EneP_
MonitorElement * meHFEneHit2_
MonitorElement * meHOL10EneP_
MonitorElement * meHBNHit_
MonitorElement * meHFTimHit_
MonitorElement * meHBEneHit2_
MonitorElement * meHONHit_
MonitorElement * meAllNHit_
MonitorElement * meHBL10Ene_
MonitorElement * meBadDetHit_
MonitorElement * meHFL10Ene_
MonitorElement * meHETimHit_
MonitorElement * meTimeHit_
MonitorElement * meHFDepHit_
MonitorElement * meHBTimHit_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * meTimeWHit_
MonitorElement * meHBDepHit_
MonitorElement * meBadSubHit_
MonitorElement * meHENHit_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * meBadIdHit_
MonitorElement * meHEL10EneP_
MonitorElement * meEnergyHit_
MonitorElement * meHOEtaHit_
MonitorElement * meEtaHit_
MonitorElement * meHBEneHit_
MonitorElement * meHOTimHit_
MonitorElement * meHBL10EneP_
MonitorElement * meHEEneHit2_
MonitorElement * meHFPhiHit_
MonitorElement * meDepthHit_
MonitorElement * meHBEtaHit_
MonitorElement * meHODepHit_
MonitorElement * meHEL10Ene_
MonitorElement * meHEEneHit_
MonitorElement * meHEPhiHit_
MonitorElement * meSubdetHit_
MonitorElement * meHBPhiHit_
MonitorElement * meHFEneHit_
MonitorElement * meHEEtaHit_
MonitorElement * meHOEneHit_
MonitorElement * meHOEneHit2_

Member Data Documentation

bool HcalSimHitStudy::checkHit_
private

Definition at line 46 of file HcalSimHitStudy.h.

Referenced by analyze(), bookHistograms(), and HcalSimHitStudy().

std::string HcalSimHitStudy::g4Label
private

Definition at line 45 of file HcalSimHitStudy.h.

Referenced by HcalSimHitStudy().

std::string HcalSimHitStudy::hcalHits
private

Definition at line 45 of file HcalSimHitStudy.h.

Referenced by HcalSimHitStudy().

MonitorElement* HcalSimHitStudy::meAllNHit_
private

Definition at line 50 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meBadDetHit_
private

Definition at line 50 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meBadIdHit_
private

Definition at line 50 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meBadSubHit_
private

Definition at line 50 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meDepthHit_
private

Definition at line 52 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meDetectHit_
private

Definition at line 52 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meEnergyHit_
private

Definition at line 53 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meEtaHit_
private

Definition at line 52 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBDepHit_
private

Definition at line 54 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBEneHit2_
private

Definition at line 59 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBEneHit_
private

Definition at line 57 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBEtaHit_
private

Definition at line 55 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBL10Ene_
private

Definition at line 60 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBL10EneP_
private

Definition at line 61 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBNHit_
private

Definition at line 51 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBPhiHit_
private

Definition at line 56 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::meHBTimHit_
private

Definition at line 58 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHEDepHit_
private

Definition at line 54 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHEEneHit2_
private

Definition at line 59 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHEEneHit_
private

Definition at line 57 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHEEtaHit_
private

Definition at line 55 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHEL10Ene_
private

Definition at line 60 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHEL10EneP_
private

Definition at line 61 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHENHit_
private

Definition at line 51 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHEPhiHit_
private

Definition at line 56 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHETimHit_
private

Definition at line 58 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFDepHit_
private

Definition at line 54 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFEneHit2_
private

Definition at line 59 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFEneHit_
private

Definition at line 57 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFEtaHit_
private

Definition at line 55 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFL10Ene_
private

Definition at line 60 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFL10EneP_
private

Definition at line 61 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFNHit_
private

Definition at line 51 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFPhiHit_
private

Definition at line 56 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHFTimHit_
private

Definition at line 58 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHODepHit_
private

Definition at line 54 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHOEneHit2_
private

Definition at line 59 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHOEneHit_
private

Definition at line 57 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHOEtaHit_
private

Definition at line 55 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHOL10Ene_
private

Definition at line 60 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHOL10EneP_
private

Definition at line 61 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHONHit_
private

Definition at line 51 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHOPhiHit_
private

Definition at line 56 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meHOTimHit_
private

Definition at line 58 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement* HcalSimHitStudy::mePhiHit_
private

Definition at line 53 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meSubdetHit_
private

Definition at line 52 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meTimeHit_
private

Definition at line 53 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

MonitorElement * HcalSimHitStudy::meTimeWHit_
private

Definition at line 53 of file HcalSimHitStudy.h.

Referenced by analyzeHits(), and bookHistograms().

std::string HcalSimHitStudy::outFile_
private

Definition at line 45 of file HcalSimHitStudy.h.

Referenced by HcalSimHitStudy().

edm::EDGetTokenT<edm::PCaloHitContainer> HcalSimHitStudy::tok_hits_
private

Definition at line 48 of file HcalSimHitStudy.h.

Referenced by analyze(), and HcalSimHitStudy().

bool HcalSimHitStudy::verbose_
private

Definition at line 46 of file HcalSimHitStudy.h.

Referenced by HcalSimHitStudy().