CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HFPMTHitAnalyzer Class Reference
Inheritance diagram for HFPMTHitAnalyzer:
edm::one::EDAnalyzer< edm::one::WatchRuns, edm::one::SharedResources > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 HFPMTHitAnalyzer (const edm::ParameterSet &)
 
 ~HFPMTHitAnalyzer () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::WatchRuns, edm::one::SharedResources >
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const noexcept final
 
bool wantsGlobalRuns () const noexcept final
 
bool wantsInputProcessBlocks () const noexcept final
 
bool wantsProcessBlocks () const noexcept final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const noexcept
 
bool wantsStreamRuns () const noexcept
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void analyzeHits (std::vector< PCaloHit > &, const std::vector< SimTrack > &)
 
void beginJob () override
 
void beginRun (edm::Run const &, edm::EventSetup const &) override
 
void endJob () override
 
void endRun (edm::Run const &, edm::EventSetup const &) override
 

Private Attributes

int event_no
 
const std::string g4Label_
 
TH1F * h_HFDepHit
 
TH1F * h_HFEta [3]
 
TH1F * h_HFPhi [3]
 
const std::string hcalHits_
 
TH1F * hHF12_time [3]
 
TH1F * hHF12_time_Ewt [3]
 
TH1F * hHF1_time [3]
 
TH1F * hHF1_time_Ewt [3]
 
TH1F * hHF2_time [3]
 
TH1F * hHF2_time_Ewt [3]
 
TH1F * hHF_e_1 [3]
 
TH1F * hHF_e_12 [3]
 
TH1F * hHF_e_2 [3]
 
TH1F * hHF_em_1 [3]
 
TH1F * hHF_em_12 [3]
 
TH1F * hHF_em_2 [3]
 
TH1F * hHF_had_1 [3]
 
TH1F * hHF_had_12 [3]
 
TH1F * hHF_had_2 [3]
 
TH1F * hHF_MC_e
 
const edm::EDGetTokenT< edm::PCaloHitContainertok_calo_
 
const edm::EDGetTokenT< edm::HepMCProducttok_evt_
 
const edm::EDGetTokenT< edm::SimTrackContainertok_track_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 34 of file HFPMTHitAnalyzer.cc.

Constructor & Destructor Documentation

◆ HFPMTHitAnalyzer()

HFPMTHitAnalyzer::HFPMTHitAnalyzer ( const edm::ParameterSet iConfig)
explicit

Definition at line 68 of file HFPMTHitAnalyzer.cc.

References TFileService::kSharedResource.

69  : g4Label_(iConfig.getUntrackedParameter<std::string>("ModuleLabel", "g4SimHits")),
70  hcalHits_(iConfig.getUntrackedParameter<std::string>("HitCollection", "HcalHits")),
71  tok_evt_(consumes<edm::HepMCProduct>(
72  edm::InputTag(iConfig.getUntrackedParameter<std::string>("SourceLabel", "VtxSmeared")))),
73  tok_calo_(consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label_, hcalHits_))),
74  tok_track_(consumes<edm::SimTrackContainer>(edm::InputTag(g4Label_))) {
75  usesResource(TFileService::kSharedResource);
76 }
static const std::string kSharedResource
Definition: TFileService.h:76
T getUntrackedParameter(std::string const &, T const &) const
const edm::EDGetTokenT< edm::PCaloHitContainer > tok_calo_
const std::string hcalHits_
const edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
const edm::EDGetTokenT< edm::SimTrackContainer > tok_track_
const std::string g4Label_

◆ ~HFPMTHitAnalyzer()

HFPMTHitAnalyzer::~HFPMTHitAnalyzer ( )
inlineoverride

Definition at line 37 of file HFPMTHitAnalyzer.cc.

37 {}

Member Function Documentation

◆ analyze()

void HFPMTHitAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 201 of file HFPMTHitAnalyzer.cc.

References analyzeHits(), event_no, edm::HepMCProduct::GetEvent(), hHF_MC_e, iEvent, AlCaHLTBitMon_ParallelJobs::p, tok_calo_, tok_evt_, tok_track_, and gather_cfg::Tracks.

201  {
202  ++event_no;
203  if (event_no % 500 == 0)
204  edm::LogVerbatim("HcalSim") << "Event # " << event_no << " processed.";
205 
206  const edm::Handle<edm::PCaloHitContainer> &hitsHcal = iEvent.getHandle(tok_calo_);
208  const edm::Handle<edm::HepMCProduct> EvtHandle = iEvent.getHandle(tok_evt_);
209  const HepMC::GenEvent *myGenEvent = EvtHandle->GetEvent();
210 
211  float orig_energy = 0;
212  for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
213  ++p) {
214  orig_energy = (*p)->momentum().e();
215  hHF_MC_e->Fill(orig_energy);
216  }
217 
218  std::vector<PCaloHit> caloHits;
219  caloHits.insert(caloHits.end(), hitsHcal->begin(), hitsHcal->end());
220  analyzeHits(caloHits, *Tracks);
221 }
Log< level::Info, true > LogVerbatim
const edm::EDGetTokenT< edm::PCaloHitContainer > tok_calo_
int iEvent
Definition: GenABIO.cc:224
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
const edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
const edm::EDGetTokenT< edm::SimTrackContainer > tok_track_
void analyzeHits(std::vector< PCaloHit > &, const std::vector< SimTrack > &)

◆ analyzeHits()

void HFPMTHitAnalyzer::analyzeHits ( std::vector< PCaloHit > &  hits,
const std::vector< SimTrack > &  tracks1 
)
private

Definition at line 223 of file HFPMTHitAnalyzer.cc.

References TauDecayModes::dec, hcalRecHitTable_cff::depth, hcalRecHitTable_cff::energy, event_no, h_HFDepHit, h_HFEta, h_HFPhi, HcalForward, hHF1_time, hHF1_time_Ewt, hHF2_time, hHF2_time_Ewt, hHF_e_1, hHF_e_12, hHF_e_2, hHF_em_1, hHF_em_12, hHF_em_2, hHF_had_1, hHF_had_12, hHF_had_2, hfClusterShapes_cfi::hits, mps_fire::i, hcalRecHitTable_cff::ieta, hcalRecHitTable_cff::iphi, and hcalRecHitTable_cff::time.

Referenced by analyze().

223  {
224  int nHit = hits.size();
225  float energy1[3], energy2[3], energy12[3];
226  float em1[3], had1[3], em2[3], had2[3], em12[3], had12[3];
227  for (int i = 0; i < 3; i++) {
228  energy1[i] = energy2[i] = energy12[i] = 0;
229  em1[i] = em2[i] = em12[i] = 0;
230  had1[i] = had2[i] = had12[i] = 0;
231  }
232  edm::LogVerbatim("HFShower") << "HFPMTHitAnalyser::Entry " << event_no << " with " << nHit << " hits";
233  for (int i = 0; i < nHit; i++) {
234  double energy = hits[i].energy();
235  double em = hits[i].energyEM();
236  double had = hits[i].energyHad();
237  double time = hits[i].time();
238  uint32_t id_ = hits[i].id();
239  uint16_t pmtHit = hits[i].depth();
240  uint16_t depthX = pmtHit;
241  edm::LogVerbatim("HFShower") << "HFPMTHitAnalyser::Hit " << i << " ID " << std::hex << id_ << std::dec << " PMT "
242  << pmtHit << " E (e|h|t) " << em << " " << had << " " << energy << " Time " << time;
243 
244  HcalDetId id = HcalDetId(id_);
245  int det = id.det();
246  int subdet = id.subdet();
247  int depth = id.depth();
248  if (pmtHit != 0)
249  pmtHit = 1;
250 
251  if (det == 4) {
252  if (subdet == static_cast<int>(HcalForward)) {
253  h_HFDepHit->Fill(double(depth + 10 * depthX));
254  if (depthX > 0) {
255  int ieta = id.ietaAbs();
256  int iphi = id.iphi();
257  if (depth != 1) {
258  ieta += 50;
259  iphi += 50;
260  }
261  if (depthX <= 3) {
262  h_HFEta[depthX - 1]->Fill(double(ieta));
263  h_HFPhi[depthX - 1]->Fill(double(iphi));
264  }
265  }
266  if (depth == 1) {
267  energy1[pmtHit] += energy;
268  energy12[pmtHit] += energy;
269  em1[pmtHit] += em;
270  em12[pmtHit] += em;
271  had1[pmtHit] += had;
272  had12[pmtHit] += had;
273  energy1[2] += energy;
274  energy12[2] += energy;
275  em1[2] += em;
276  em12[2] += em;
277  had1[2] += had;
278  had12[2] += had;
279  hHF1_time[pmtHit]->Fill(time);
280  hHF1_time[2]->Fill(time);
281  hHF1_time_Ewt[pmtHit]->Fill(time, energy);
282  hHF1_time_Ewt[2]->Fill(time, energy);
283  }
284  if (depth == 2) {
285  energy2[pmtHit] += energy;
286  energy12[pmtHit] += energy;
287  em2[pmtHit] += em;
288  em12[pmtHit] += em;
289  had2[pmtHit] += had;
290  had12[pmtHit] += had;
291  energy2[2] += energy;
292  energy12[2] += energy;
293  em2[2] += em;
294  em12[2] += em;
295  had2[2] += had;
296  had12[2] += had;
297  hHF2_time[pmtHit]->Fill(time);
298  hHF2_time[2]->Fill(time);
299  hHF2_time_Ewt[pmtHit]->Fill(time, energy);
300  hHF2_time_Ewt[2]->Fill(time, energy);
301  }
302  }
303  }
304  }
305  for (int i = 0; i < 3; i++) {
306  hHF_e_1[i]->Fill(energy1[i]);
307  hHF_e_2[i]->Fill(energy2[i]);
308  hHF_e_12[i]->Fill(energy12[i]);
309  hHF_em_1[i]->Fill(em1[i]);
310  hHF_em_2[i]->Fill(em2[i]);
311  hHF_em_12[i]->Fill(em12[i]);
312  hHF_had_1[i]->Fill(had1[i]);
313  hHF_had_2[i]->Fill(had2[i]);
314  hHF_had_12[i]->Fill(had12[i]);
315  edm::LogVerbatim("HFShower") << "HFPMTHitAnalyser:: Type " << i << " Energy 1|2| " << energy1[i] << " "
316  << energy2[i] << " " << energy12[i] << " EM Energy 1|2| " << em1[i] << " " << em2[i]
317  << " " << em12[i] << " Had Energy 1|2| " << had1[i] << " " << had2[i] << " "
318  << had12[i];
319  }
320 }
Log< level::Info, true > LogVerbatim

◆ beginJob()

void HFPMTHitAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 86 of file HFPMTHitAnalyzer.cc.

References event_no, compareTotals::fs, h_HFDepHit, h_HFEta, h_HFPhi, hHF12_time, hHF12_time_Ewt, hHF1_time, hHF1_time_Ewt, hHF2_time, hHF2_time_Ewt, hHF_e_1, hHF_e_12, hHF_e_2, hHF_em_1, hHF_em_12, hHF_em_2, hHF_had_1, hHF_had_12, hHF_had_2, hHF_MC_e, mps_fire::i, TFileDirectory::make(), Skims_PA_cff::name, and runGCPTkAlMap::title.

86  {
87  event_no = 0;
88  char name[20], title[120], sub[11];
89 
91  if (!fs.isAvailable())
92  throw cms::Exception("BadConfig") << "TFileService unavailable: "
93  << "please add it to config file";
94 
95  TFileDirectory HFHitsDir = fs->mkdir("HFPMTHits");
96  h_HFDepHit = HFHitsDir.make<TH1F>("Hit20", "Depths in HF", 50, 0., 50.);
97  h_HFDepHit->GetXaxis()->SetTitle("Depths in HF");
98  for (int i = 0; i < 3; i++) {
99  if (i == 0)
100  sprintf(sub, "(PMT)");
101  else if (i == 1)
102  sprintf(sub, "(Bundle)");
103  else
104  sprintf(sub, "(Jungle)");
105  sprintf(name, "Eta%d", i);
106  sprintf(title, "Eta Index of hits in %s", sub);
107  h_HFEta[i] = HFHitsDir.make<TH1F>(name, title, 100, 0, 100.);
108  h_HFEta[i]->GetXaxis()->SetTitle(title);
109  sprintf(name, "Phi%d", i);
110  sprintf(title, "Phi Index of hits in %s", sub);
111  h_HFPhi[i] = HFHitsDir.make<TH1F>(name, title, 100, 0, 100.);
112  h_HFPhi[i]->GetXaxis()->SetTitle(title);
113  }
114  TFileDirectory HFSourcePart = fs->mkdir("HFMCinfo");
115  hHF_MC_e = HFSourcePart.make<TH1F>("MCEnergy", "Energy of Generated Particle", 1000, 0., 500.);
116  hHF_MC_e->GetXaxis()->SetTitle("Energy of Generated Particle");
117 
118  //energy Histograms
119  TFileDirectory HFPCaloHitEnergyDir = fs->mkdir("HFPCaloHitEnergy");
120  for (int i = 0; i < 3; i++) {
121  if (i == 0)
122  sprintf(sub, "(Absorber)");
123  else if (i == 1)
124  sprintf(sub, "(PMT)");
125  else
126  sprintf(sub, "(All)");
127  sprintf(name, "Energy1%d", i);
128  sprintf(title, "Energy in depth 1 %s", sub);
129  hHF_e_1[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 500.);
130  hHF_e_1[i]->GetXaxis()->SetTitle(title);
131  sprintf(name, "Energy2%d", i);
132  sprintf(title, "Energy in depth 2 %s", sub);
133  hHF_e_2[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 500.);
134  hHF_e_2[i]->GetXaxis()->SetTitle(title);
135  sprintf(name, "Energy12%d", i);
136  sprintf(title, "Energy in depths 1,2 %s", sub);
137  hHF_e_12[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 500.);
138  hHF_e_12[i]->GetXaxis()->SetTitle(title);
139  sprintf(name, "Em1%d", i);
140  sprintf(title, "EM energy in depth 1 %s", sub);
141  hHF_em_1[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 500.);
142  hHF_em_1[i]->GetXaxis()->SetTitle(title);
143  sprintf(name, "Em2%d", i);
144  sprintf(title, "EM energy in depth 2 %s", sub);
145  hHF_em_2[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 500.);
146  hHF_em_2[i]->GetXaxis()->SetTitle(title);
147  sprintf(name, "Em12%d", i);
148  sprintf(title, "EM energy in depths 1,2 %s", sub);
149  hHF_em_12[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 500.);
150  hHF_em_12[i]->GetXaxis()->SetTitle(title);
151  sprintf(name, "Had1%d", i);
152  sprintf(title, "Had energy in depth 1 %s", sub);
153  hHF_had_1[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 0.1);
154  hHF_had_1[i]->GetXaxis()->SetTitle(title);
155  sprintf(name, "Had2%d", i);
156  sprintf(title, "Had energy in depth 2 %s", sub);
157  hHF_had_2[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 0.1);
158  hHF_had_2[i]->GetXaxis()->SetTitle(title);
159  sprintf(name, "Had12%d", i);
160  sprintf(title, "Had energy in depths 1,2 %s", sub);
161  hHF_had_12[i] = HFPCaloHitEnergyDir.make<TH1F>(name, title, 1000, 0., 0.1);
162  hHF_had_12[i]->GetXaxis()->SetTitle(title);
163  }
164 
165  //Timing Histograms
166  TFileDirectory HFPCaloHitTimeDir = fs->mkdir("HFPCaloHitTime");
167  for (int i = 0; i < 3; i++) {
168  if (i == 0)
169  sprintf(sub, "(Absorber)");
170  else if (i == 1)
171  sprintf(sub, "(PMT)");
172  else
173  sprintf(sub, "(All)");
174  sprintf(name, "Time1Ewt%d", i);
175  sprintf(title, "Time (energy weighted) in depth 1 %s", sub);
176  hHF1_time_Ewt[i] = HFPCaloHitTimeDir.make<TH1F>(name, title, 400, 0., 400.);
177  hHF1_time_Ewt[i]->GetXaxis()->SetTitle(title);
178  sprintf(name, "Time2Ewt%d", i);
179  sprintf(title, "Time (energy weighted) in depth 2 %s", sub);
180  hHF2_time_Ewt[i] = HFPCaloHitTimeDir.make<TH1F>(name, title, 400, 0., 400.);
181  hHF2_time_Ewt[i]->GetXaxis()->SetTitle(title);
182  sprintf(name, "Time12Ewt%d", i);
183  sprintf(title, "Time (energy weighted) in depths 1,2 %s", sub);
184  hHF12_time_Ewt[i] = HFPCaloHitTimeDir.make<TH1F>(name, title, 400, 0., 400.);
185  hHF12_time_Ewt[i]->GetXaxis()->SetTitle(title);
186  sprintf(name, "Time1%d", i);
187  sprintf(title, "Time in depth 1 %s", sub);
188  hHF1_time[i] = HFPCaloHitTimeDir.make<TH1F>(name, title, 400, 0., 400.);
189  hHF1_time[i]->GetXaxis()->SetTitle(title);
190  sprintf(name, "Time2%d", i);
191  sprintf(title, "Time in depth 2 %s", sub);
192  hHF2_time[i] = HFPCaloHitTimeDir.make<TH1F>(name, title, 400, 0., 400.);
193  hHF2_time[i]->GetXaxis()->SetTitle(title);
194  sprintf(name, "Time12%d", i);
195  sprintf(title, "Time in depths 1,2 %s", sub);
196  hHF12_time[i] = HFPCaloHitTimeDir.make<TH1F>(name, title, 400, 0., 400.);
197  hHF12_time[i]->GetXaxis()->SetTitle(title);
198  }
199 }
TH1F * hHF12_time_Ewt[3]
T * make(const Args &...args) const
make new ROOT object

◆ beginRun()

void HFPMTHitAnalyzer::beginRun ( edm::Run const &  ,
edm::EventSetup const &   
)
inlineoverrideprivate

Definition at line 42 of file HFPMTHitAnalyzer.cc.

42 {}

◆ endJob()

void HFPMTHitAnalyzer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 322 of file HFPMTHitAnalyzer.cc.

322 {}

◆ endRun()

void HFPMTHitAnalyzer::endRun ( edm::Run const &  ,
edm::EventSetup const &   
)
inlineoverrideprivate

Definition at line 45 of file HFPMTHitAnalyzer.cc.

45 {}

◆ fillDescriptions()

void HFPMTHitAnalyzer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 78 of file HFPMTHitAnalyzer.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, and AlCaHLTBitMon_QueryRunRegistry::string.

78  {
80  desc.addUntracked<std::string>("SourceLabel", "generatorSmeared");
81  desc.addUntracked<std::string>("ModuleLabel", "g4SimHits");
82  desc.addUntracked<std::string>("HitCollection", "HcalHits");
83  descriptions.add("HFPMTHitAnalyzer", desc);
84 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ event_no

int HFPMTHitAnalyzer::event_no
private

Definition at line 54 of file HFPMTHitAnalyzer.cc.

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

◆ g4Label_

const std::string HFPMTHitAnalyzer::g4Label_
private

Definition at line 49 of file HFPMTHitAnalyzer.cc.

◆ h_HFDepHit

TH1F* HFPMTHitAnalyzer::h_HFDepHit
private

Definition at line 57 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ h_HFEta

TH1F * HFPMTHitAnalyzer::h_HFEta[3]
private

Definition at line 57 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ h_HFPhi

TH1F * HFPMTHitAnalyzer::h_HFPhi[3]
private

Definition at line 57 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hcalHits_

const std::string HFPMTHitAnalyzer::hcalHits_
private

Definition at line 49 of file HFPMTHitAnalyzer.cc.

◆ hHF12_time

TH1F* HFPMTHitAnalyzer::hHF12_time[3]
private

Definition at line 65 of file HFPMTHitAnalyzer.cc.

Referenced by beginJob().

◆ hHF12_time_Ewt

TH1F * HFPMTHitAnalyzer::hHF12_time_Ewt[3]
private

Definition at line 65 of file HFPMTHitAnalyzer.cc.

Referenced by beginJob().

◆ hHF1_time

TH1F* HFPMTHitAnalyzer::hHF1_time[3]
private

Definition at line 63 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF1_time_Ewt

TH1F * HFPMTHitAnalyzer::hHF1_time_Ewt[3]
private

Definition at line 63 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF2_time

TH1F* HFPMTHitAnalyzer::hHF2_time[3]
private

Definition at line 64 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF2_time_Ewt

TH1F * HFPMTHitAnalyzer::hHF2_time_Ewt[3]
private

Definition at line 64 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_e_1

TH1F* HFPMTHitAnalyzer::hHF_e_1[3]
private

Definition at line 59 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_e_12

TH1F* HFPMTHitAnalyzer::hHF_e_12[3]
private

Definition at line 61 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_e_2

TH1F* HFPMTHitAnalyzer::hHF_e_2[3]
private

Definition at line 60 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_em_1

TH1F * HFPMTHitAnalyzer::hHF_em_1[3]
private

Definition at line 59 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_em_12

TH1F * HFPMTHitAnalyzer::hHF_em_12[3]
private

Definition at line 61 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_em_2

TH1F * HFPMTHitAnalyzer::hHF_em_2[3]
private

Definition at line 60 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_had_1

TH1F * HFPMTHitAnalyzer::hHF_had_1[3]
private

Definition at line 59 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_had_12

TH1F * HFPMTHitAnalyzer::hHF_had_12[3]
private

Definition at line 61 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_had_2

TH1F * HFPMTHitAnalyzer::hHF_had_2[3]
private

Definition at line 60 of file HFPMTHitAnalyzer.cc.

Referenced by analyzeHits(), and beginJob().

◆ hHF_MC_e

TH1F * HFPMTHitAnalyzer::hHF_MC_e
private

Definition at line 57 of file HFPMTHitAnalyzer.cc.

Referenced by analyze(), and beginJob().

◆ tok_calo_

const edm::EDGetTokenT<edm::PCaloHitContainer> HFPMTHitAnalyzer::tok_calo_
private

Definition at line 51 of file HFPMTHitAnalyzer.cc.

Referenced by analyze().

◆ tok_evt_

const edm::EDGetTokenT<edm::HepMCProduct> HFPMTHitAnalyzer::tok_evt_
private

Definition at line 50 of file HFPMTHitAnalyzer.cc.

Referenced by analyze().

◆ tok_track_

const edm::EDGetTokenT<edm::SimTrackContainer> HFPMTHitAnalyzer::tok_track_
private

Definition at line 52 of file HFPMTHitAnalyzer.cc.

Referenced by analyze().