CMS 3D CMS Logo

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

#include <HcalSimHitsValidation.h>

Inheritance diagram for HcalSimHitsValidation:
DQMEDAnalyzer edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void analyze (edm::Event const &ev, edm::EventSetup const &c) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void endJob () override
 
 HcalSimHitsValidation (edm::ParameterSet const &conf)
 
 ~HcalSimHitsValidation () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &ev, edm::EventSetup const &es) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) override
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
 DQMEDAnalyzer (DQMEDAnalyzer const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer &&)=delete
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) override
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) override
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) override
 
 ~DQMEDAnalyzer () override=default
 
- Public Member Functions inherited from edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns >
 EDProducer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () 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
 
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::vector< ModuleDescription const * > &modules, 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
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

double dPhiWsign (double phi1, double phi2)
 
double dR (double eta1, double phi1, double eta2, double phi2)
 
double phi12 (double phi1, double en1, double phi2, double en2)
 

Private Attributes

bool auxPlots_
 
std::string ebHits_
 
std::string eeHits_
 
std::vector< MonitorElement * > emean_vs_ieta_HB
 
std::vector< MonitorElement * > emean_vs_ieta_HE
 
std::vector< MonitorElement * > emean_vs_ieta_HF
 
MonitorElementemean_vs_ieta_HO
 
std::string g4Label_
 
edm::ESHandle< CaloGeometrygeometry
 
std::string hcalHits_
 
const HcalDDDRecConstantshcons
 
int maxDepthHB_
 
int maxDepthHE_
 
int maxDepthHF_
 
int maxDepthHO_
 
MonitorElementmeEnConeEtaProfile
 
MonitorElementmeEnConeEtaProfile_E
 
MonitorElementmeEnConeEtaProfile_EH
 
std::vector< MonitorElement * > meSimHitsEnergyHB
 
std::vector< MonitorElement * > meSimHitsEnergyHE
 
std::vector< MonitorElement * > meSimHitsEnergyHF
 
MonitorElementmeSimHitsEnergyHO
 
int nevtot
 
std::vector< MonitorElement * > Nhb
 
std::vector< MonitorElement * > Nhe
 
std::vector< MonitorElement * > Nhf
 
MonitorElementNho
 
std::vector< MonitorElement * > occupancy_vs_ieta_HB
 
std::vector< MonitorElement * > occupancy_vs_ieta_HE
 
std::vector< MonitorElement * > occupancy_vs_ieta_HF
 
MonitorElementoccupancy_vs_ieta_HO
 
std::string outputFile_
 
bool testNumber_
 
edm::EDGetTokenT< edm::PCaloHitContainertok_ecalEB_
 
edm::EDGetTokenT< edm::PCaloHitContainertok_ecalEE_
 
edm::EDGetTokenT< edm::HepMCProducttok_evt_
 
edm::EDGetTokenT< edm::PCaloHitContainertok_hcal_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- 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)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

Definition at line 43 of file HcalSimHitsValidation.h.

Constructor & Destructor Documentation

HcalSimHitsValidation::HcalSimHitsValidation ( edm::ParameterSet const &  conf)

Definition at line 5 of file HcalSimHitsValidation.cc.

References auxPlots_, ebHits_, eeHits_, g4Label_, edm::ParameterSet::getUntrackedParameter(), hcalHits_, nevtot, outputFile_, AlCaHLTBitMon_QueryRunRegistry::string, testNumber_, tok_ecalEB_, tok_ecalEE_, tok_evt_, and tok_hcal_.

5  {
6  // DQM ROOT output
7  outputFile_ = conf.getUntrackedParameter<std::string>("outputFile", "myfile.root");
8  testNumber_ = conf.getUntrackedParameter<bool>("TestNumber",false);
9  auxPlots_ = conf.getUntrackedParameter<bool>("auxiliaryPlots",false);
10 
11  // register for data access
12  g4Label_ = conf.getUntrackedParameter<std::string>("ModuleLabel","g4SimHits");
13  hcalHits_ = conf.getUntrackedParameter<std::string>("HcalHitCollection","HcalHits");
14  ebHits_ = conf.getUntrackedParameter<std::string>("EBHitCollection","EcalHitsEB");
15  eeHits_ = conf.getUntrackedParameter<std::string>("EEHitCollection","EcalHitsEE");
16 
17  tok_evt_ = consumes<edm::HepMCProduct>(edm::InputTag("generatorSmeared"));
18  tok_hcal_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label_,hcalHits_));
19  tok_ecalEB_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label_,ebHits_));
20  tok_ecalEE_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label_,eeHits_));
21 
22  if ( !outputFile_.empty() ) { edm::LogInfo("OutputInfo") << " Hcal SimHit Task histograms will be saved to '" << outputFile_.c_str() << "'";
23  } else {
24  edm::LogInfo("OutputInfo") << " Hcal SimHit Task histograms will NOT be saved";
25  }
26 
27  nevtot = 0;
28 
29 }
edm::EDGetTokenT< edm::PCaloHitContainer > tok_ecalEB_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_ecalEE_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hcal_
edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
HcalSimHitsValidation::~HcalSimHitsValidation ( )
override

Definition at line 32 of file HcalSimHitsValidation.cc.

32 { }

Member Function Documentation

void HcalSimHitsValidation::analyze ( edm::Event const &  ev,
edm::EventSetup const &  c 
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 253 of file HcalSimHitsValidation.cc.

References funct::abs(), auxPlots_, gather_cfg::cout, egammaForCoreTracking_cff::depth, dR(), ebHits_, eeHits_, emean_vs_ieta_HB, emean_vs_ieta_HE, emean_vs_ieta_HF, emean_vs_ieta_HO, ALCARECOTkAlBeamHalo_cff::etaMax, MonitorElement::Fill(), edm::EventSetup::get(), edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), hcons, edm::HandleBase::isValid(), MuonErrorMatrixAnalyzer_cfi::maxPt, meEnConeEtaProfile, meEnConeEtaProfile_E, meEnConeEtaProfile_EH, meSimHitsEnergyHB, meSimHitsEnergyHE, meSimHitsEnergyHF, meSimHitsEnergyHO, nevtot, npart, occupancy_vs_ieta_HB, occupancy_vs_ieta_HE, occupancy_vs_ieta_HF, occupancy_vs_ieta_HO, AlCaHLTBitMon_ParallelJobs::p, edm::Handle< T >::product(), EnergyCorrector::pt, alignCSCRings::r, HcalHitRelabeller::relabel(), HcalDetId::subdet(), testNumber_, tok_ecalEB_, tok_ecalEE_, tok_evt_, and tok_hcal_.

253  {
254 
255  using namespace edm;
256  using namespace std;
257 
258  //===========================================================================
259  // Getting SimHits
260  //===========================================================================
261 
262  double phi_MC = -999.; // phi of initial particle from HepMC
263  double eta_MC = -999.; // eta of initial particle from HepMC
264 
266  ev.getByToken(tok_evt_,evtMC); // generator in late 310_preX
267  if (!evtMC.isValid()) {
268  std::cout << "no HepMCProduct found" << std::endl;
269  }
270 
271  // MC particle with highest pt is taken as a direction reference
272  double maxPt = -99999.;
273  int npart = 0;
274 
275  const HepMC::GenEvent * myGenEvent = evtMC->GetEvent();
276  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
277  p != myGenEvent->particles_end(); ++p ) {
278  double phip = (*p)->momentum().phi();
279  double etap = (*p)->momentum().eta();
280  double pt = (*p)->momentum().perp();
281  if(pt > maxPt) {npart++; maxPt = pt; phi_MC = phip; eta_MC = etap; }
282  }
283 
284  double partR = 0.3;
285 
286 
287  //Hcal SimHits
288 
289  //Approximate calibration constants
290  const float calib_HB = 120.;
291  const float calib_HE = 190.;
292  const float calib_HF1 = 1.0/0.383;
293  const float calib_HF2 = 1.0/0.368;
294 
296  ev.getByToken(tok_hcal_,hcalHits);
297  const PCaloHitContainer * SimHitResult = hcalHits.product () ;
298 
299  float eta_diff;
300  float etaMax = 9999;
301  int ietaMax = 0;
302 
303  double HcalCone = 0;
304 
305  c.get<CaloGeometryRecord>().get (geometry);
306 
307  for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResult->begin () ; SimHits != SimHitResult->end(); ++SimHits) {
308  HcalDetId cell;
309  if (testNumber_) cell = HcalHitRelabeller::relabel(SimHits->id(),hcons);
310  else cell = HcalDetId(SimHits->id());
311 
312  auto cellGeometry = geometry->getSubdetectorGeometry(cell)->getGeometry(cell);
313  double etaS = cellGeometry->getPosition().eta () ;
314  double phiS = cellGeometry->getPosition().phi () ;
315  double en = SimHits->energy();
316 
317  int sub = cell.subdet();
318  int depth = cell.depth();
319  double ieta = cell.ieta();
320 
321  //Energy in Cone
322  double r = dR(eta_MC, phi_MC, etaS, phiS);
323 
324  if (r < partR){
325  eta_diff = std::abs(eta_MC - etaS);
326  if(eta_diff < etaMax) {
327  etaMax = eta_diff;
328  ietaMax = cell.ieta();
329  }
330  //Approximation of calibration
331  if (sub == 1) HcalCone += en*calib_HB;
332  else if (sub == 2) HcalCone += en*calib_HE;
333  else if (sub == 4 && (depth == 1 || depth == 3)) HcalCone += en*calib_HF1;
334  else if (sub == 4 && (depth == 2 || depth == 4)) HcalCone += en*calib_HF2;
335  }
336 
337  if (auxPlots_) {
338 
339  //HB
340  if (sub == 1){
341  meSimHitsEnergyHB[0]->Fill(en);
342  meSimHitsEnergyHB[depth]->Fill(en);
343 
344  emean_vs_ieta_HB[0]->Fill(double(ieta), en);
345  emean_vs_ieta_HB[depth]->Fill(double(ieta), en);
346 
347  occupancy_vs_ieta_HB[0]->Fill(double(ieta));
348  occupancy_vs_ieta_HB[depth]->Fill(double(ieta));
349  }
350  //HE
351  if (sub == 2){
352  meSimHitsEnergyHE[0]->Fill(en);
353  meSimHitsEnergyHE[depth]->Fill(en);
354 
355  emean_vs_ieta_HE[0]->Fill(double(ieta), en);
356  emean_vs_ieta_HE[depth]->Fill(double(ieta), en);
357 
358  occupancy_vs_ieta_HE[0]->Fill(double(ieta));
359  occupancy_vs_ieta_HE[depth]->Fill(double(ieta));
360  }
361  //HO
362  if (sub == 3){
363  meSimHitsEnergyHO->Fill(en);
364 
365  emean_vs_ieta_HO->Fill(double(ieta), en);
366 
367  occupancy_vs_ieta_HO->Fill(double(ieta));
368  }
369  //HF
370  if (sub == 4){
371  meSimHitsEnergyHF[0]->Fill(en);
372  meSimHitsEnergyHF[depth]->Fill(en);
373 
374  emean_vs_ieta_HF[0]->Fill(double(ieta), en);
375  emean_vs_ieta_HF[depth]->Fill(double(ieta), en);
376 
377  occupancy_vs_ieta_HF[0]->Fill(double(ieta));
378  occupancy_vs_ieta_HF[depth]->Fill(double(ieta));
379  }
380 
381  } // auxPlots_
382 
383  } //Loop over SimHits
384 
385  //Ecal EB SimHits
386  double EcalCone = 0;
387 
388  if (!ebHits_.empty()){
390  ev.getByToken(tok_ecalEB_,ecalEBHits);
391  const PCaloHitContainer * SimHitResultEB = ecalEBHits.product () ;
392 
393  for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResultEB->begin () ; SimHits != SimHitResultEB->end(); ++SimHits) {
394 
395  EBDetId EBid = EBDetId(SimHits->id());
396 
397  auto cellGeometry = geometry->getSubdetectorGeometry(EBid)->getGeometry(EBid);
398  double etaS = cellGeometry->getPosition().eta () ;
399  double phiS = cellGeometry->getPosition().phi () ;
400  double en = SimHits->energy();
401 
402  double r = dR(eta_MC, phi_MC, etaS, phiS);
403 
404  if (r < partR) EcalCone += en;
405  }
406  } // ebHits_
407 
408  //Ecal EE SimHits
409  if (!eeHits_.empty()){
411  ev.getByToken(tok_ecalEE_,ecalEEHits);
412  const PCaloHitContainer * SimHitResultEE = ecalEEHits.product () ;
413 
414  for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResultEE->begin () ; SimHits != SimHitResultEE->end(); ++SimHits) {
415 
416  EEDetId EEid = EEDetId(SimHits->id());
417 
418  auto cellGeometry = geometry->getSubdetectorGeometry(EEid)->getGeometry(EEid) ;
419  double etaS = cellGeometry->getPosition().eta () ;
420  double phiS = cellGeometry->getPosition().phi () ;
421  double en = SimHits->energy();
422 
423  double r = dR(eta_MC, phi_MC, etaS, phiS);
424 
425  if (r < partR) EcalCone += en;
426  }
427  } // eeHits_
428 
429  if (ietaMax != 0){ //If ietaMax == 0, there were no good HCAL SimHits
430  meEnConeEtaProfile ->Fill(double(ietaMax), HcalCone);
431  meEnConeEtaProfile_E ->Fill(double(ietaMax), EcalCone);
432  meEnConeEtaProfile_EH ->Fill(double(ietaMax), HcalCone+EcalCone);
433  }
434 
435  nevtot++;
436 }
MonitorElement * meEnConeEtaProfile
std::vector< MonitorElement * > meSimHitsEnergyHB
std::vector< PCaloHit > PCaloHitContainer
double dR(double eta1, double phi1, double eta2, double phi2)
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
std::vector< MonitorElement * > emean_vs_ieta_HE
MonitorElement * meSimHitsEnergyHO
double npart
Definition: HydjetWrapper.h:49
bool ev
std::vector< MonitorElement * > occupancy_vs_ieta_HF
edm::EDGetTokenT< edm::PCaloHitContainer > tok_ecalEB_
void Fill(long long x)
MonitorElement * meEnConeEtaProfile_E
std::vector< MonitorElement * > emean_vs_ieta_HF
edm::EDGetTokenT< edm::PCaloHitContainer > tok_ecalEE_
MonitorElement * emean_vs_ieta_HO
std::vector< MonitorElement * > emean_vs_ieta_HB
std::vector< MonitorElement * > meSimHitsEnergyHF
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hcal_
bool isValid() const
Definition: HandleBase.h:74
const HcalDDDRecConstants * hcons
std::vector< MonitorElement * > occupancy_vs_ieta_HE
std::vector< MonitorElement * > meSimHitsEnergyHE
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
T const * product() const
Definition: Handle.h:81
MonitorElement * meEnConeEtaProfile_EH
HLT enums.
MonitorElement * occupancy_vs_ieta_HO
DetId relabel(const uint32_t testId) const
edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
std::vector< MonitorElement * > occupancy_vs_ieta_HB
void HcalSimHitsValidation::bookHistograms ( DQMStore::IBooker ib,
edm::Run const &  run,
edm::EventSetup const &  es 
)
overridevirtual

Implements DQMEDAnalyzer.

Definition at line 35 of file HcalSimHitsValidation.cc.

References auxPlots_, DQMStore::IBooker::book1D(), DQMStore::IBooker::bookProfile(), egammaForCoreTracking_cff::depth, emean_vs_ieta_HB, emean_vs_ieta_HE, emean_vs_ieta_HF, emean_vs_ieta_HO, edm::EventSetup::get(), HcalDDDRecConstants::getEtaRange(), HcalDDDRecConstants::getMaxDepth(), HcalDDDRecConstants::getNPhi(), hcons, trackerHits::histo, hcalTTPDigis_cfi::iEtaMax, createfilelist::int, maxDepthHB_, maxDepthHE_, maxDepthHF_, maxDepthHO_, meEnConeEtaProfile, meEnConeEtaProfile_E, meEnConeEtaProfile_EH, meSimHitsEnergyHB, meSimHitsEnergyHE, meSimHitsEnergyHF, meSimHitsEnergyHO, Nhb, Nhe, Nhf, Nho, occupancy_vs_ieta_HB, occupancy_vs_ieta_HE, occupancy_vs_ieta_HF, occupancy_vs_ieta_HO, and DQMStore::IBooker::setCurrentFolder().

36 {
38  es.get<HcalRecNumberingRecord>().get( pHRNDC );
39  hcons = &(*pHRNDC);
44 
45  //Get Phi segmentation from geometry, use the max phi number so that all iphi values are included.
46 
47  int NphiMax = hcons->getNPhi(0);
48 
49  NphiMax = (hcons->getNPhi(1) > NphiMax ? hcons->getNPhi(1) : NphiMax);
50  NphiMax = (hcons->getNPhi(2) > NphiMax ? hcons->getNPhi(2) : NphiMax);
51  NphiMax = (hcons->getNPhi(3) > NphiMax ? hcons->getNPhi(3) : NphiMax);
52 
53  //Center the iphi bins on the integers
54  //float iphi_min = 0.5;
55  //float iphi_max = NphiMax + 0.5;
56  //int iphi_bins = (int) (iphi_max - iphi_min);
57 
58  int iEtaHBMax = hcons->getEtaRange(0).second;
59  int iEtaHEMax = hcons->getEtaRange(1).second;
60  int iEtaHFMax = hcons->getEtaRange(2).second;
61  int iEtaHOMax = hcons->getEtaRange(3).second;
62 
63  //Retain classic behavior, all plots have same ieta range.
64  //Comment out code to allow each subdetector to have its on range
65 
66  int iEtaMax = (iEtaHBMax > iEtaHEMax ? iEtaHBMax : iEtaHEMax);
67  iEtaMax = (iEtaMax > iEtaHFMax ? iEtaMax : iEtaHFMax);
68  iEtaMax = (iEtaMax > iEtaHOMax ? iEtaMax : iEtaHOMax);
69 
70  iEtaHBMax = iEtaMax;
71  iEtaHEMax = iEtaMax;
72  iEtaHFMax = iEtaMax;
73  iEtaHOMax = iEtaMax;
74 
75  //Give an empty bin around the subdet ieta range to make it clear that all ieta rings have been included
76  float ieta_min_HB = -iEtaHBMax - 1.5;
77  float ieta_max_HB = iEtaHBMax + 1.5;
78  int ieta_bins_HB = (int) (ieta_max_HB - ieta_min_HB);
79 
80  float ieta_min_HE = -iEtaHEMax - 1.5;
81  float ieta_max_HE = iEtaHEMax + 1.5;
82  int ieta_bins_HE = (int) (ieta_max_HE - ieta_min_HE);
83 
84  float ieta_min_HF = -iEtaHFMax - 1.5;
85  float ieta_max_HF = iEtaHFMax + 1.5;
86  int ieta_bins_HF = (int) (ieta_max_HF - ieta_min_HF);
87 
88  float ieta_min_HO = -iEtaHOMax - 1.5;
89  float ieta_max_HO = iEtaHOMax + 1.5;
90  int ieta_bins_HO = (int) (ieta_max_HO - ieta_min_HO);
91 
92  Char_t histo[200];
93 
94  ib.setCurrentFolder("HcalHitsV/HcalSimHitTask");
95 
96  if (auxPlots_) {
97 
98  // General counters
99  for(int depth = 0; depth <= maxDepthHB_; depth++){
100  if(depth == 0){ sprintf (histo, "N_HB" ); }
101  else { sprintf (histo, "N_HB%d", depth ); }
102 
103  Nhb.push_back( ib.book1D(histo, histo, 2600,0.,2600.) );
104  }
105  for(int depth = 0; depth <= maxDepthHE_; depth++){
106  if(depth == 0){ sprintf (histo, "N_HE" ); }
107  else { sprintf (histo, "N_HE%d", depth ); }
108 
109  Nhe.push_back( ib.book1D(histo, histo, 2600,0.,2600.) );
110  }
111 
112  sprintf (histo, "N_HO" );
113  Nho = ib.book1D(histo, histo, 2200,0.,2200.);
114 
115  for(int depth = 0; depth <= maxDepthHF_; depth++){
116  if(depth == 0){ sprintf (histo, "N_HF" ); }
117  else { sprintf (histo, "N_HF%d", depth ); }
118 
119  Nhf.push_back( ib.book1D(histo, histo, 1800,0.,1800.) );
120  }
121 
122  //Mean energy vs iEta TProfiles
123  for(int depth = 0; depth <= maxDepthHB_; depth++){
124  if(depth == 0){ sprintf (histo, "emean_vs_ieta_HB" ); }
125  else { sprintf (histo, "emean_vs_ieta_HB%d", depth ); }
126 
127  emean_vs_ieta_HB.push_back( ib.bookProfile(histo, histo, ieta_bins_HB, ieta_min_HB, ieta_max_HB, -10., 2000., " ") );
128  }
129  for(int depth = 0; depth <= maxDepthHE_; depth++){
130  if(depth == 0){ sprintf (histo, "emean_vs_ieta_HE" ); }
131  else { sprintf (histo, "emean_vs_ieta_HE%d", depth ); }
132 
133  emean_vs_ieta_HE.push_back( ib.bookProfile(histo, histo, ieta_bins_HE, ieta_min_HE, ieta_max_HE, -10., 2000., " ") );
134  }
135 
136  sprintf (histo, "emean_vs_ieta_HO" );
137  emean_vs_ieta_HO = ib.bookProfile(histo, histo, ieta_bins_HO, ieta_min_HO, ieta_max_HO, -10., 2000., " ");
138 
139  for(int depth = 0; depth <= maxDepthHF_; depth++){
140  if(depth == 0){ sprintf (histo, "emean_vs_ieta_HF" ); }
141  else { sprintf (histo, "emean_vs_ieta_HF%d", depth ); }
142 
143  emean_vs_ieta_HF.push_back( ib.bookProfile(histo, histo, ieta_bins_HF, ieta_min_HF, ieta_max_HF, -10., 2000., " ") );
144  }
145 
146  //Occupancy vs. iEta TH1Fs
147  for(int depth = 0; depth <= maxDepthHB_; depth++){
148  if(depth == 0){ sprintf (histo, "occupancy_vs_ieta_HB" ); }
149  else { sprintf (histo, "occupancy_vs_ieta_HB%d", depth ); }
150 
151  occupancy_vs_ieta_HB.push_back( ib.book1D(histo, histo, ieta_bins_HB, ieta_min_HB, ieta_max_HB) );
152  }
153  for(int depth = 0; depth <= maxDepthHE_; depth++){
154  if(depth == 0){ sprintf (histo, "occupancy_vs_ieta_HE" ); }
155  else { sprintf (histo, "occupancy_vs_ieta_HE%d", depth ); }
156 
157  occupancy_vs_ieta_HE.push_back( ib.book1D(histo, histo, ieta_bins_HE, ieta_min_HE, ieta_max_HE) );
158  }
159 
160  sprintf (histo, "occupancy_vs_ieta_HO" );
161  occupancy_vs_ieta_HO = ib.book1D(histo, histo, ieta_bins_HO, ieta_min_HO, ieta_max_HO);
162 
163  for(int depth = 0; depth <= maxDepthHF_; depth++){
164  if(depth == 0){ sprintf (histo, "occupancy_vs_ieta_HF" ); }
165  else { sprintf (histo, "occupancy_vs_ieta_HF%d", depth ); }
166 
167  occupancy_vs_ieta_HF.push_back( ib.book1D(histo, histo, ieta_bins_HF, ieta_min_HF, ieta_max_HF) );
168  }
169 
170  //Energy spectra
171  for(int depth = 0; depth <= maxDepthHB_; depth++){
172  if(depth == 0){ sprintf (histo, "HcalSimHitTask_energy_of_simhits_HB" ); }
173  else { sprintf (histo, "HcalSimHitTask_energy_of_simhits_HB%d", depth ); }
174 
175  meSimHitsEnergyHB.push_back( ib.book1D(histo, histo, 510 , -0.1 , 5.) );
176  }
177  for(int depth = 0; depth <= maxDepthHE_; depth++){
178  if(depth == 0){ sprintf (histo, "HcalSimHitTask_energy_of_simhits_HE" ); }
179  else { sprintf (histo, "HcalSimHitTask_energy_of_simhits_HE%d", depth ); }
180 
181  meSimHitsEnergyHE.push_back( ib.book1D(histo, histo, 510 , -0.1 , 5.) );
182  }
183 
184  sprintf (histo, "HcalSimHitTask_energy_of_simhits_HO" );
185  meSimHitsEnergyHO = ib.book1D(histo, histo, 510 , -0.1 , 5.);
186 
187  for(int depth = 0; depth <= maxDepthHF_; depth++){
188  if(depth == 0){ sprintf (histo, "HcalSimHitTask_energy_of_simhits_HF" ); }
189  else { sprintf (histo, "HcalSimHitTask_energy_of_simhits_HF%d", depth ); }
190 
191  meSimHitsEnergyHF.push_back( ib.book1D(histo, histo, 1010 , -5. , 500.) );
192  }
193 
194  } // auxPlots_
195 
196  //Energy in Cone
197  sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths");
198  meEnConeEtaProfile = ib.bookProfile(histo, histo, ieta_bins_HF, ieta_min_HF, ieta_max_HF, -10., 200., " ");
199 
200  sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths_E");
201  meEnConeEtaProfile_E = ib.bookProfile(histo, histo, ieta_bins_HF, ieta_min_HF, ieta_max_HF, -10., 200., " ");
202 
203  sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths_EH");
204  meEnConeEtaProfile_EH = ib.bookProfile(histo, histo, ieta_bins_HF, ieta_min_HF, ieta_max_HF, -10., 200., " ");
205 
206 
207 }
MonitorElement * meEnConeEtaProfile
std::vector< MonitorElement * > meSimHitsEnergyHB
std::vector< MonitorElement * > emean_vs_ieta_HE
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:160
MonitorElement * meSimHitsEnergyHO
std::vector< MonitorElement * > occupancy_vs_ieta_HF
std::vector< MonitorElement * > Nhe
std::pair< int, int > getEtaRange(const int &i) const
MonitorElement * meEnConeEtaProfile_E
std::vector< MonitorElement * > emean_vs_ieta_HF
std::vector< MonitorElement * > Nhb
MonitorElement * emean_vs_ieta_HO
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
std::vector< MonitorElement * > emean_vs_ieta_HB
std::vector< MonitorElement * > meSimHitsEnergyHF
const HcalDDDRecConstants * hcons
std::vector< MonitorElement * > occupancy_vs_ieta_HE
std::vector< MonitorElement * > meSimHitsEnergyHE
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
int getMaxDepth(const int &type) const
MonitorElement * meEnConeEtaProfile_EH
std::vector< MonitorElement * > Nhf
MonitorElement * occupancy_vs_ieta_HO
std::vector< MonitorElement * > occupancy_vs_ieta_HB
int getNPhi(const int &type) const
double HcalSimHitsValidation::dPhiWsign ( double  phi1,
double  phi2 
)
private

Definition at line 465 of file HcalSimHitsValidation.cc.

References DEFINE_FWK_MODULE, PI, and tmp.

465  {
466  // clockwise phi2 w.r.t phi1 means "+" phi distance
467  // anti-clockwise phi2 w.r.t phi1 means "-" phi distance
468 
469  double PI = 3.1415926535898;
470  double a1 = phi1; double a2 = phi2;
471  double tmp = a2 - a1;
472  if( a1*a2 < 0.) {
473  if(a1 > 0.5 * PI) tmp += 2.*PI;
474  if(a2 > 0.5 * PI) tmp -= 2.*PI;
475  }
476  return tmp;
477 
478 }
#define PI
Definition: QcdUeDQM.h:36
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double HcalSimHitsValidation::dR ( double  eta1,
double  phi1,
double  eta2,
double  phi2 
)
private

Definition at line 439 of file HcalSimHitsValidation.cc.

References PI, mathSSE::sqrt(), and tmp.

Referenced by analyze().

439  {
440  double PI = 3.1415926535898;
441  double deltaphi= phi1 - phi2;
442  if( phi2 > phi1 ) { deltaphi= phi2 - phi1;}
443  if(deltaphi > PI) { deltaphi = 2.*PI - deltaphi;}
444  double deltaeta = eta2 - eta1;
445  double tmp = sqrt(deltaeta* deltaeta + deltaphi*deltaphi);
446  return tmp;
447 }
T sqrt(T t)
Definition: SSEVec.h:18
#define PI
Definition: QcdUeDQM.h:36
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
void HcalSimHitsValidation::endJob ( void  )
overridevirtual

Reimplemented from edm::one::EDProducerBase.

Definition at line 210 of file HcalSimHitsValidation.cc.

References funct::abs(), auxPlots_, egammaForCoreTracking_cff::depth, MonitorElement::getBinContent(), mps_fire::i, maxDepthHB_, maxDepthHE_, maxDepthHF_, nevtot, occupancy_vs_ieta_HB, occupancy_vs_ieta_HE, occupancy_vs_ieta_HF, occupancy_vs_ieta_HO, and MonitorElement::setBinContent().

210  {
211 
212  if (auxPlots_){
213 
214  for (int i = 1; i <= occupancy_vs_ieta_HB[0]->getNbinsX(); i++){
215 
216  int ieta = i - 43; // -41 -1, 1 41
217 
218  float phi_factor;
219 
220  if (std::abs(ieta) <= 20) phi_factor = 72.;
221  else if (std::abs(ieta) < 40) phi_factor = 36.;
222  else phi_factor = 18.;
223 
224  float cnorm;
225 
226  //Occupancy vs. iEta TH1Fs
227  for(int depth = 0; depth <= maxDepthHB_; depth++){
228  cnorm = occupancy_vs_ieta_HB[depth]->getBinContent(i) / (phi_factor * nevtot);
229  occupancy_vs_ieta_HB[depth]->setBinContent(i, cnorm);
230  }
231  for(int depth = 0; depth <= maxDepthHE_; depth++){
232  cnorm = occupancy_vs_ieta_HE[depth]->getBinContent(i) / (phi_factor * nevtot);
233  occupancy_vs_ieta_HE[depth]->setBinContent(i, cnorm);
234  }
235 
236  cnorm = occupancy_vs_ieta_HO->getBinContent(i) / (phi_factor * nevtot);
238 
239  for(int depth = 0; depth <= maxDepthHF_; depth++){
240  cnorm = occupancy_vs_ieta_HF[depth]->getBinContent(i) / (phi_factor * nevtot);
241  occupancy_vs_ieta_HF[depth]->setBinContent(i, cnorm);
242  }
243 
244  }
245 
246  }
247 
248  // let's see if this breaks anything
249  //if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
250 }
void setBinContent(int binx, double content)
set content of bin (1-D)
std::vector< MonitorElement * > occupancy_vs_ieta_HF
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< MonitorElement * > occupancy_vs_ieta_HE
double getBinContent(int binx) const
get content of bin (1-D)
MonitorElement * occupancy_vs_ieta_HO
std::vector< MonitorElement * > occupancy_vs_ieta_HB
double HcalSimHitsValidation::phi12 ( double  phi1,
double  en1,
double  phi2,
double  en2 
)
private

Definition at line 449 of file HcalSimHitsValidation.cc.

References PI, and tmp.

449  {
450  // weighted mean value of phi1 and phi2
451 
452  double tmp;
453  double PI = 3.1415926535898;
454  double a1 = phi1; double a2 = phi2;
455 
456  if( a1 > 0.5*PI && a2 < 0.) a2 += 2*PI;
457  if( a2 > 0.5*PI && a1 < 0.) a1 += 2*PI;
458  tmp = (a1 * en1 + a2 * en2)/(en1 + en2);
459  if(tmp > PI) tmp -= 2.*PI;
460 
461  return tmp;
462 
463 }
#define PI
Definition: QcdUeDQM.h:36
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100

Member Data Documentation

bool HcalSimHitsValidation::auxPlots_
private

Definition at line 73 of file HcalSimHitsValidation.h.

Referenced by analyze(), bookHistograms(), endJob(), and HcalSimHitsValidation().

std::string HcalSimHitsValidation::ebHits_
private

Definition at line 61 of file HcalSimHitsValidation.h.

Referenced by analyze(), and HcalSimHitsValidation().

std::string HcalSimHitsValidation::eeHits_
private

Definition at line 61 of file HcalSimHitsValidation.h.

Referenced by analyze(), and HcalSimHitsValidation().

std::vector<MonitorElement*> HcalSimHitsValidation::emean_vs_ieta_HB
private

Definition at line 86 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::emean_vs_ieta_HE
private

Definition at line 87 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::emean_vs_ieta_HF
private

Definition at line 89 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* HcalSimHitsValidation::emean_vs_ieta_HO
private

Definition at line 88 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

std::string HcalSimHitsValidation::g4Label_
private

Definition at line 61 of file HcalSimHitsValidation.h.

Referenced by HcalSimHitsValidation().

edm::ESHandle<CaloGeometry> HcalSimHitsValidation::geometry
private
std::string HcalSimHitsValidation::hcalHits_
private

Definition at line 61 of file HcalSimHitsValidation.h.

Referenced by HcalSimHitsValidation().

const HcalDDDRecConstants* HcalSimHitsValidation::hcons
private

Definition at line 68 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

int HcalSimHitsValidation::maxDepthHB_
private

Definition at line 69 of file HcalSimHitsValidation.h.

Referenced by bookHistograms(), and endJob().

int HcalSimHitsValidation::maxDepthHE_
private

Definition at line 69 of file HcalSimHitsValidation.h.

Referenced by bookHistograms(), and endJob().

int HcalSimHitsValidation::maxDepthHF_
private

Definition at line 70 of file HcalSimHitsValidation.h.

Referenced by bookHistograms(), and endJob().

int HcalSimHitsValidation::maxDepthHO_
private

Definition at line 70 of file HcalSimHitsValidation.h.

Referenced by bookHistograms().

MonitorElement* HcalSimHitsValidation::meEnConeEtaProfile
private

Definition at line 97 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* HcalSimHitsValidation::meEnConeEtaProfile_E
private

Definition at line 98 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* HcalSimHitsValidation::meEnConeEtaProfile_EH
private

Definition at line 99 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::meSimHitsEnergyHB
private

Definition at line 102 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::meSimHitsEnergyHE
private

Definition at line 103 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::meSimHitsEnergyHF
private

Definition at line 105 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* HcalSimHitsValidation::meSimHitsEnergyHO
private

Definition at line 104 of file HcalSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

int HcalSimHitsValidation::nevtot
private

Definition at line 110 of file HcalSimHitsValidation.h.

Referenced by analyze(), endJob(), and HcalSimHitsValidation().

std::vector<MonitorElement*> HcalSimHitsValidation::Nhb
private

Definition at line 76 of file HcalSimHitsValidation.h.

Referenced by bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::Nhe
private

Definition at line 77 of file HcalSimHitsValidation.h.

Referenced by bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::Nhf
private

Definition at line 79 of file HcalSimHitsValidation.h.

Referenced by bookHistograms().

MonitorElement* HcalSimHitsValidation::Nho
private

Definition at line 78 of file HcalSimHitsValidation.h.

Referenced by bookHistograms().

std::vector<MonitorElement*> HcalSimHitsValidation::occupancy_vs_ieta_HB
private

Definition at line 91 of file HcalSimHitsValidation.h.

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

std::vector<MonitorElement*> HcalSimHitsValidation::occupancy_vs_ieta_HE
private

Definition at line 92 of file HcalSimHitsValidation.h.

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

std::vector<MonitorElement*> HcalSimHitsValidation::occupancy_vs_ieta_HF
private

Definition at line 94 of file HcalSimHitsValidation.h.

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

MonitorElement* HcalSimHitsValidation::occupancy_vs_ieta_HO
private

Definition at line 93 of file HcalSimHitsValidation.h.

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

std::string HcalSimHitsValidation::outputFile_
private

Definition at line 60 of file HcalSimHitsValidation.h.

Referenced by HcalSimHitsValidation().

bool HcalSimHitsValidation::testNumber_
private

Definition at line 72 of file HcalSimHitsValidation.h.

Referenced by analyze(), and HcalSimHitsValidation().

edm::EDGetTokenT<edm::PCaloHitContainer> HcalSimHitsValidation::tok_ecalEB_
private

Definition at line 65 of file HcalSimHitsValidation.h.

Referenced by analyze(), and HcalSimHitsValidation().

edm::EDGetTokenT<edm::PCaloHitContainer> HcalSimHitsValidation::tok_ecalEE_
private

Definition at line 66 of file HcalSimHitsValidation.h.

Referenced by analyze(), and HcalSimHitsValidation().

edm::EDGetTokenT<edm::HepMCProduct> HcalSimHitsValidation::tok_evt_
private

Definition at line 63 of file HcalSimHitsValidation.h.

Referenced by analyze(), and HcalSimHitsValidation().

edm::EDGetTokenT<edm::PCaloHitContainer> HcalSimHitsValidation::tok_hcal_
private

Definition at line 64 of file HcalSimHitsValidation.h.

Referenced by analyze(), and HcalSimHitsValidation().