CMS 3D CMS Logo

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

#include <EcalBarrelSimHitsValidation.h>

Inheritance diagram for EcalBarrelSimHitsValidation:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Types

typedef dqm::legacy::DQMStore DQMStore
 
typedef dqm::legacy::MonitorElement MonitorElement
 
- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Public Member Functions

 EcalBarrelSimHitsValidation (const edm::ParameterSet &ps)
 Constructor. More...
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c) override
 Analyze. More...
 
void bookHistograms (DQMStore::IBooker &ib, edm::Run const &, edm::EventSetup const &c) override
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 

Private Types

typedef std::map< uint32_t, float, std::less< uint32_t > > MapType
 

Private Member Functions

float eCluster2x2 (MapType &themap)
 
float eCluster4x4 (float e33, MapType &themap)
 
virtual float energyInMatrixEB (int nCellInEta, int nCellInPhi, int centralEta, int centralPhi, int centralZ, MapType &themap)
 
bool fillEBMatrix (int nCellInEta, int nCellInPhi, int CentralEta, int CentralPhi, int CentralZ, MapType &fillmap, MapType &themap)
 
std::vector< uint32_t > getIdsAroundMax (int nCellInEta, int nCellInPhi, int centralEta, int centralPhi, int centralZ, MapType &themap)
 
uint32_t getUnitWithMaxEnergy (MapType &themap)
 

Private Attributes

std::string EBHitsCollection
 
edm::EDGetTokenT< edm::PCaloHitContainerEBHitsToken
 
float eRLength [26]
 
std::string g4InfoLabel
 
MonitorElementmeEBcrystalEnergy2_
 
MonitorElementmeEBcrystalEnergy_
 
MonitorElementmeEBe16_
 
MonitorElementmeEBe16oe25_
 
MonitorElementmeEBe1_
 
MonitorElementmeEBe1oe25_
 
MonitorElementmeEBe1oe4_
 
MonitorElementmeEBe1oe9_
 
MonitorElementmeEBe25_
 
MonitorElementmeEBe4_
 
MonitorElementmeEBe4oe9_
 
MonitorElementmeEBe9_
 
MonitorElementmeEBe9oe16_
 
MonitorElementmeEBe9oe25_
 
MonitorElementmeEBhitEnergy2_
 
MonitorElementmeEBhitEnergy_
 
MonitorElementmeEBhitLog10Energy25Norm_
 
MonitorElementmeEBhitLog10Energy_
 
MonitorElementmeEBhitLog10EnergyNorm_
 
MonitorElementmeEBLongitudinalShower_
 
MonitorElementmeEBOccupancy_
 
MonitorElementmenEBCrystals_
 
MonitorElementmenEBHits_
 
int myEntries
 
std::string ValidationCollection
 
edm::EDGetTokenT< PEcalValidInfoValidationCollectionToken
 
bool verbose_
 

Additional Inherited Members

- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

Definition at line 29 of file EcalBarrelSimHitsValidation.h.

Member Typedef Documentation

◆ DQMStore

Definition at line 33 of file EcalBarrelSimHitsValidation.h.

◆ MapType

typedef std::map<uint32_t, float, std::less<uint32_t> > EcalBarrelSimHitsValidation::MapType
private

Definition at line 30 of file EcalBarrelSimHitsValidation.h.

◆ MonitorElement

Definition at line 34 of file EcalBarrelSimHitsValidation.h.

Constructor & Destructor Documentation

◆ EcalBarrelSimHitsValidation()

EcalBarrelSimHitsValidation::EcalBarrelSimHitsValidation ( const edm::ParameterSet ps)

Constructor.

Definition at line 17 of file EcalBarrelSimHitsValidation.cc.

References EBHitsCollection, EBHitsToken, eRLength, g4InfoLabel, edm::ParameterSet::getUntrackedParameter(), ProducerED_cfi::InputTag, myEntries, AlCaHLTBitMon_QueryRunRegistry::string, ValidationCollection, ValidationCollectionToken, and verbose_.

18  : g4InfoLabel(ps.getParameter<std::string>("moduleLabelG4")),
19  EBHitsCollection(ps.getParameter<std::string>("EBHitsCollection")),
20  ValidationCollection(ps.getParameter<std::string>("ValidationCollection")) {
21  EBHitsToken =
22  consumes<edm::PCaloHitContainer>(edm::InputTag(std::string(g4InfoLabel), std::string(EBHitsCollection)));
25  // verbosity switch
26  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
27 
28  myEntries = 0;
29  for (int myStep = 0; myStep < 26; myStep++) {
30  eRLength[myStep] = 0.0;
31  }
32 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::PCaloHitContainer > EBHitsToken
edm::EDGetTokenT< PEcalValidInfo > ValidationCollectionToken

Member Function Documentation

◆ analyze()

void EcalBarrelSimHitsValidation::analyze ( const edm::Event e,
const edm::EventSetup c 
)
overrideprotectedvirtual

Analyze.

Reimplemented from DQMEDAnalyzer.

Definition at line 108 of file EcalBarrelSimHitsValidation.cc.

References nano_mu_digi_cff::bx, PEcalValidInfo::bX0(), MillePedeFileConverter_cfg::e, PEcalValidInfo::eb1x1(), EBHitsToken, eCluster2x2(), eCluster4x4(), hcalRecHitTable_cff::energy, energyInMatrixEB(), eRLength, dqm::impl::MonitorElement::Fill(), fillEBMatrix(), getIdsAroundMax(), getUnitWithMaxEnergy(), mps_fire::i, EBDetId::ietaAbs(), createfilelist::int, EBDetId::iphi(), edm::HandleBase::isValid(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, dqmiolumiharvest::j, LogDebug, genParticles_cff::map, meEBcrystalEnergy2_, meEBcrystalEnergy_, meEBe16_, meEBe16oe25_, meEBe1_, meEBe1oe25_, meEBe1oe4_, meEBe1oe9_, meEBe25_, meEBe4_, meEBe4oe9_, meEBe9_, meEBe9oe16_, meEBe9oe25_, meEBhitEnergy2_, meEBhitEnergy_, meEBhitLog10Energy25Norm_, meEBhitLog10Energy_, meEBhitLog10EnergyNorm_, meEBLongitudinalShower_, meEBOccupancy_, menEBCrystals_, menEBHits_, myEntries, dqm::impl::MonitorElement::Reset(), ValidationCollectionToken, and EBDetId::zside().

108  {
109  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
110 
112  e.getByToken(EBHitsToken, EcalHitsEB);
113 
114  // Do nothing if no Barrel data available
115  if (!EcalHitsEB.isValid())
116  return;
117 
118  edm::Handle<PEcalValidInfo> MyPEcalValidInfo;
119  e.getByToken(ValidationCollectionToken, MyPEcalValidInfo);
120 
121  std::vector<PCaloHit> theEBCaloHits;
122  theEBCaloHits.insert(theEBCaloHits.end(), EcalHitsEB->begin(), EcalHitsEB->end());
123 
124  myEntries++;
125 
126  double EBEnergy_ = 0.;
127  std::map<unsigned int, std::vector<PCaloHit *>, std::less<unsigned int>> CaloHitMap;
128 
129  double eb1 = 0.0;
130  double eb4 = 0.0;
131  double eb9 = 0.0;
132  double eb16 = 0.0;
133  double eb25 = 0.0;
134  std::vector<double> econtr(140, 0.);
135  std::vector<double> econtr25(140, 0.);
136 
137  MapType ebmap;
138  uint32_t nEBHits = 0;
139 
140  for (std::vector<PCaloHit>::iterator isim = theEBCaloHits.begin(); isim != theEBCaloHits.end(); ++isim) {
141  if (isim->time() > 500.) {
142  continue;
143  }
144 
145  CaloHitMap[isim->id()].push_back(&(*isim));
146 
147  EBDetId ebid(isim->id());
148 
149  LogDebug("HitInfo") << " CaloHit " << isim->getName() << "\n"
150  << " DetID = " << isim->id() << " EBDetId = " << ebid.ieta() << " " << ebid.iphi() << "\n"
151  << " Time = " << isim->time() << "\n"
152  << " Track Id = " << isim->geantTrackId() << "\n"
153  << " Energy = " << isim->energy();
154 
155  meEBOccupancy_->Fill(ebid.iphi(), ebid.ieta());
156 
157  uint32_t crystid = ebid.rawId();
158  ebmap[crystid] += isim->energy();
159 
160  EBEnergy_ += isim->energy();
161  nEBHits++;
162  meEBhitEnergy_->Fill(isim->energy());
163  if (isim->energy() > 0) {
164  meEBhitLog10Energy_->Fill(log10(isim->energy()));
165  int log10i = int((log10(isim->energy()) + 10.) * 10.);
166  if (log10i >= 0 && log10i < 140)
167  econtr[log10i] += isim->energy();
168  }
169  meEBhitEnergy2_->Fill(isim->energy());
170  }
171 
172  menEBCrystals_->Fill(ebmap.size());
173  if (meEBcrystalEnergy_) {
174  for (std::map<uint32_t, float, std::less<uint32_t>>::iterator it = ebmap.begin(); it != ebmap.end(); ++it)
175  meEBcrystalEnergy_->Fill((*it).second);
176  }
177  if (meEBcrystalEnergy2_) {
178  for (std::map<uint32_t, float, std::less<uint32_t>>::iterator it = ebmap.begin(); it != ebmap.end(); ++it)
179  meEBcrystalEnergy2_->Fill((*it).second);
180  }
181 
182  menEBHits_->Fill(nEBHits);
183 
184  if (nEBHits > 0) {
185  uint32_t ebcenterid = getUnitWithMaxEnergy(ebmap);
186  EBDetId myEBid(ebcenterid);
187  int bx = myEBid.ietaAbs();
188  int by = myEBid.iphi();
189  int bz = myEBid.zside();
190  eb1 = energyInMatrixEB(1, 1, bx, by, bz, ebmap);
191  meEBe1_->Fill(eb1);
192  eb9 = energyInMatrixEB(3, 3, bx, by, bz, ebmap);
193  meEBe9_->Fill(eb9);
194  eb25 = energyInMatrixEB(5, 5, bx, by, bz, ebmap);
195  meEBe25_->Fill(eb25);
196 
197  std::vector<uint32_t> ids25;
198  ids25 = getIdsAroundMax(5, 5, bx, by, bz, ebmap);
199 
200  for (unsigned i = 0; i < 25; i++) {
201  for (unsigned int j = 0; j < CaloHitMap[ids25[i]].size(); j++) {
202  if (CaloHitMap[ids25[i]][j]->energy() > 0) {
203  int log10i = int((log10(CaloHitMap[ids25[i]][j]->energy()) + 10.) * 10.);
204  if (log10i >= 0 && log10i < 140)
205  econtr25[log10i] += CaloHitMap[ids25[i]][j]->energy();
206  }
207  }
208  }
209 
210  MapType newebmap;
211  if (fillEBMatrix(3, 3, bx, by, bz, newebmap, ebmap)) {
212  eb4 = eCluster2x2(newebmap);
213  meEBe4_->Fill(eb4);
214  }
215  if (fillEBMatrix(5, 5, bx, by, bz, newebmap, ebmap)) {
216  eb16 = eCluster4x4(eb9, newebmap);
217  meEBe16_->Fill(eb16);
218  }
219 
220  if (eb4 > 0.1)
221  meEBe1oe4_->Fill(eb1 / eb4);
222  if (eb9 > 0.1)
223  meEBe1oe9_->Fill(eb1 / eb9);
224  if (eb9 > 0.1)
225  meEBe4oe9_->Fill(eb4 / eb9);
226  if (eb16 > 0.1)
227  meEBe9oe16_->Fill(eb9 / eb16);
228  if (eb25 > 0.1)
229  meEBe1oe25_->Fill(eb1 / eb25);
230  if (eb25 > 0.1)
231  meEBe9oe25_->Fill(eb9 / eb25);
232  if (eb25 > 0.1)
233  meEBe16oe25_->Fill(eb16 / eb25);
234 
235  if (EBEnergy_ != 0) {
236  for (int i = 0; i < 140; i++) {
237  meEBhitLog10EnergyNorm_->Fill(-10. + (float(i) + 0.5) / 10., econtr[i] / EBEnergy_);
238  }
239  }
240 
241  if (eb25 != 0) {
242  for (int i = 0; i < 140; i++) {
243  meEBhitLog10Energy25Norm_->Fill(-10. + (float(i) + 0.5) / 10., econtr25[i] / eb25);
244  }
245  }
246  }
247 
248  if (MyPEcalValidInfo.isValid()) {
249  if (MyPEcalValidInfo->eb1x1() > 0.) {
250  std::vector<float> BX0 = MyPEcalValidInfo->bX0();
252  for (int myStep = 0; myStep < 26; myStep++) {
253  eRLength[myStep] += BX0[myStep];
254  meEBLongitudinalShower_->Fill(float(myStep), eRLength[myStep] / myEntries);
255  }
256  }
257  }
258 }
float eb1x1() const
std::vector< uint32_t > getIdsAroundMax(int nCellInEta, int nCellInPhi, int centralEta, int centralPhi, int centralZ, MapType &themap)
void Fill(long long x)
virtual void Reset()
Remove all data from the ME, keept the empty histogram with all its settings.
float eCluster4x4(float e33, MapType &themap)
FloatVector bX0() const
std::map< uint32_t, float, std::less< uint32_t > > MapType
edm::EDGetTokenT< edm::PCaloHitContainer > EBHitsToken
edm::EDGetTokenT< PEcalValidInfo > ValidationCollectionToken
uint32_t getUnitWithMaxEnergy(MapType &themap)
virtual float energyInMatrixEB(int nCellInEta, int nCellInPhi, int centralEta, int centralPhi, int centralZ, MapType &themap)
Log< level::Info, false > LogInfo
bool isValid() const
Definition: HandleBase.h:70
bool fillEBMatrix(int nCellInEta, int nCellInPhi, int CentralEta, int CentralPhi, int CentralZ, MapType &fillmap, MapType &themap)
#define LogDebug(id)

◆ bookHistograms()

void EcalBarrelSimHitsValidation::bookHistograms ( DQMStore::IBooker ib,
edm::Run const &  ,
edm::EventSetup const &  c 
)
overrideprotectedvirtual

Implements DQMEDAnalyzer.

Definition at line 34 of file EcalBarrelSimHitsValidation.cc.

References timingPdfMaker::histo, cuy::ib, meEBcrystalEnergy2_, meEBcrystalEnergy_, meEBe16_, meEBe16oe25_, meEBe1_, meEBe1oe25_, meEBe1oe4_, meEBe1oe9_, meEBe25_, meEBe4_, meEBe4oe9_, meEBe9_, meEBe9oe16_, meEBe9oe25_, meEBhitEnergy2_, meEBhitEnergy_, meEBhitLog10Energy25Norm_, meEBhitLog10Energy_, meEBhitLog10EnergyNorm_, meEBLongitudinalShower_, meEBOccupancy_, menEBCrystals_, menEBHits_, and AlCaHLTBitMon_QueryRunRegistry::string.

34  {
35  ib.setCurrentFolder("EcalHitsV/EcalSimHitsValidation");
36  ib.setScope(MonitorElementData::Scope::RUN);
37 
38  std::string histo = "EB hits multiplicity";
39  menEBHits_ = ib.book1D(histo, histo, 50, 0., 5000.);
40 
41  histo = "EB crystals multiplicity";
42  menEBCrystals_ = ib.book1D(histo, histo, 200, 0., 2000.);
43 
44  histo = "EB occupancy";
45  meEBOccupancy_ = ib.book2D(histo, histo, 360, 0., 360., 170, -85., 85.);
46 
47  histo = "EB longitudinal shower profile";
48  meEBLongitudinalShower_ = ib.bookProfile(histo, histo, 26, 0., 26., 100, 0., 20000.);
49 
50  histo = "EB hits energy spectrum";
51  meEBhitEnergy_ = ib.book1D(histo, histo, 4000, 0., 400.);
52 
53  histo = "EB hits log10energy spectrum";
54  meEBhitLog10Energy_ = ib.book1D(histo, histo, 140, -10., 4.);
55 
56  histo = "EB hits log10energy spectrum vs normalized energy";
57  meEBhitLog10EnergyNorm_ = ib.bookProfile(histo, histo, 140, -10., 4., 100, 0., 1.);
58 
59  histo = "EB hits log10energy spectrum vs normalized energy25";
60  meEBhitLog10Energy25Norm_ = ib.bookProfile(histo, histo, 140, -10., 4., 100, 0., 1.);
61 
62  histo = "EB hits energy spectrum 2";
63  meEBhitEnergy2_ = ib.book1D(histo, histo, 1000, 0., 0.001);
64 
65  histo = "EB crystal energy spectrum";
66  meEBcrystalEnergy_ = ib.book1D(histo, histo, 5000, 0., 50.);
67 
68  histo = "EB crystal energy spectrum 2";
69  meEBcrystalEnergy2_ = ib.book1D(histo, histo, 1000, 0., 0.001);
70 
71  histo = "EB E1";
72  meEBe1_ = ib.book1D(histo, histo, 400, 0., 400.);
73 
74  histo = "EB E4";
75  meEBe4_ = ib.book1D(histo, histo, 400, 0., 400.);
76 
77  histo = "EB E9";
78  meEBe9_ = ib.book1D(histo, histo, 400, 0., 400.);
79 
80  histo = "EB E16";
81  meEBe16_ = ib.book1D(histo, histo, 400, 0., 400.);
82 
83  histo = "EB E25";
84  meEBe25_ = ib.book1D(histo, histo, 400, 0., 400.);
85 
86  histo = "EB E1oE4";
87  meEBe1oe4_ = ib.book1D(histo, histo, 100, 0.4, 1.1);
88 
89  histo = "EB E1oE9";
90  meEBe1oe9_ = ib.book1D(histo, histo, 100, 0.4, 1.1);
91 
92  histo = "EB E4oE9";
93  meEBe4oe9_ = ib.book1D(histo, histo, 100, 0.4, 1.1);
94 
95  histo = "EB E9oE16";
96  meEBe9oe16_ = ib.book1D(histo, histo, 100, 0.4, 1.1);
97 
98  histo = "EB E1oE25";
99  meEBe1oe25_ = ib.book1D(histo, histo, 100, 0.4, 1.1);
100 
101  histo = "EB E9oE25";
102  meEBe9oe25_ = ib.book1D(histo, histo, 100, 0.4, 1.1);
103 
104  histo = "EB E16oE25";
105  meEBe16oe25_ = ib.book1D(histo, histo, 100, 0.4, 1.1);
106 }
ib
Definition: cuy.py:661

◆ eCluster2x2()

float EcalBarrelSimHitsValidation::eCluster2x2 ( MapType themap)
private

Definition at line 359 of file EcalBarrelSimHitsValidation.cc.

Referenced by analyze().

359  {
360  float E22 = 0.;
361  float e012 = themap[0] + themap[1] + themap[2];
362  float e036 = themap[0] + themap[3] + themap[6];
363  float e678 = themap[6] + themap[7] + themap[8];
364  float e258 = themap[2] + themap[5] + themap[8];
365 
366  if ((e012 > e678 || e012 == e678) && (e036 > e258 || e036 == e258))
367  E22 = themap[0] + themap[1] + themap[3] + themap[4];
368  else if ((e012 > e678 || e012 == e678) && (e036 < e258 || e036 == e258))
369  E22 = themap[1] + themap[2] + themap[4] + themap[5];
370  else if ((e012 < e678 || e012 == e678) && (e036 > e258 || e036 == e258))
371  E22 = themap[3] + themap[4] + themap[6] + themap[7];
372  else if ((e012 < e678 || e012 == e678) && (e036 < e258 || e036 == e258))
373  E22 = themap[4] + themap[5] + themap[7] + themap[8];
374 
375  return E22;
376 }

◆ eCluster4x4()

float EcalBarrelSimHitsValidation::eCluster4x4 ( float  e33,
MapType themap 
)
private

Definition at line 378 of file EcalBarrelSimHitsValidation.cc.

Referenced by analyze().

378  {
379  float E44 = 0.;
380  float e0_4 = themap[0] + themap[1] + themap[2] + themap[3] + themap[4];
381  float e0_20 = themap[0] + themap[5] + themap[10] + themap[15] + themap[20];
382  float e4_24 = themap[4] + themap[9] + themap[14] + themap[19] + themap[24];
383  float e0_24 = themap[20] + themap[21] + themap[22] + themap[23] + themap[24];
384 
385  if ((e0_4 > e0_24 || e0_4 == e0_24) && (e0_20 > e4_24 || e0_20 == e4_24))
386  E44 = e33 + themap[0] + themap[1] + themap[2] + themap[3] + themap[5] + themap[10] + themap[15];
387  else if ((e0_4 > e0_24 || e0_4 == e0_24) && (e0_20 < e4_24 || e0_20 == e4_24))
388  E44 = e33 + themap[1] + themap[2] + themap[3] + themap[4] + themap[9] + themap[14] + themap[19];
389  else if ((e0_4 < e0_24 || e0_4 == e0_24) && (e0_20 > e4_24 || e0_20 == e4_24))
390  E44 = e33 + themap[5] + themap[10] + themap[15] + themap[20] + themap[21] + themap[22] + themap[23];
391  else if ((e0_4 < e0_24 || e0_4 == e0_24) && (e0_20 < e4_24 || e0_20 == e4_24))
392  E44 = e33 + themap[21] + themap[22] + themap[23] + themap[24] + themap[9] + themap[14] + themap[19];
393 
394  return E44;
395 }

◆ energyInMatrixEB()

float EcalBarrelSimHitsValidation::energyInMatrixEB ( int  nCellInEta,
int  nCellInPhi,
int  centralEta,
int  centralPhi,
int  centralZ,
MapType themap 
)
privatevirtual

Definition at line 260 of file EcalBarrelSimHitsValidation.cc.

References funct::abs(), hcalRecHitTable_cff::ieta, hcalRecHitTable_cff::iphi, LogDebug, and DetId::rawId().

Referenced by analyze().

261  {
262  int ncristals = 0;
263  float totalEnergy = 0.;
264 
265  int goBackInEta = nCellInEta / 2;
266  int goBackInPhi = nCellInPhi / 2;
267  int startEta = centralZ * centralEta - goBackInEta;
268  int startPhi = centralPhi - goBackInPhi;
269 
270  for (int ieta = startEta; ieta < startEta + nCellInEta; ieta++) {
271  for (int iphi = startPhi; iphi < startPhi + nCellInPhi; iphi++) {
272  uint32_t index;
273  if (abs(ieta) > 85 || abs(ieta) < 1) {
274  continue;
275  }
276  if (iphi < 1) {
277  index = EBDetId(ieta, iphi + 360).rawId();
278  } else if (iphi > 360) {
279  index = EBDetId(ieta, iphi - 360).rawId();
280  } else {
281  index = EBDetId(ieta, iphi).rawId();
282  }
283 
284  totalEnergy += themap[index];
285  ncristals += 1;
286  }
287  }
288 
289  LogDebug("GeomInfo") << nCellInEta << " x " << nCellInPhi << " EB matrix energy = " << totalEnergy << " for "
290  << ncristals << " crystals";
291  return totalEnergy;
292 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
#define LogDebug(id)

◆ fillEBMatrix()

bool EcalBarrelSimHitsValidation::fillEBMatrix ( int  nCellInEta,
int  nCellInPhi,
int  CentralEta,
int  CentralPhi,
int  CentralZ,
MapType fillmap,
MapType themap 
)
private

Definition at line 325 of file EcalBarrelSimHitsValidation.cc.

References funct::abs(), getUnitWithMaxEnergy(), mps_fire::i, hcalRecHitTable_cff::ieta, hcalRecHitTable_cff::iphi, and DetId::rawId().

Referenced by analyze().

326  {
327  int goBackInEta = nCellInEta / 2;
328  int goBackInPhi = nCellInPhi / 2;
329 
330  int startEta = CentralZ * CentralEta - goBackInEta;
331  int startPhi = CentralPhi - goBackInPhi;
332 
333  int i = 0;
334  for (int ieta = startEta; ieta < startEta + nCellInEta; ieta++) {
335  for (int iphi = startPhi; iphi < startPhi + nCellInPhi; iphi++) {
336  uint32_t index;
337  if (abs(ieta) > 85 || abs(ieta) < 1) {
338  continue;
339  }
340  if (iphi < 1) {
341  index = EBDetId(ieta, iphi + 360).rawId();
342  } else if (iphi > 360) {
343  index = EBDetId(ieta, iphi - 360).rawId();
344  } else {
345  index = EBDetId(ieta, iphi).rawId();
346  }
347  fillmap[i++] = themap[index];
348  }
349  }
350 
351  uint32_t ebcenterid = getUnitWithMaxEnergy(themap);
352 
353  if (fillmap[i / 2] == themap[ebcenterid])
354  return true;
355  else
356  return false;
357 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
uint32_t getUnitWithMaxEnergy(MapType &themap)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57

◆ getIdsAroundMax()

std::vector< uint32_t > EcalBarrelSimHitsValidation::getIdsAroundMax ( int  nCellInEta,
int  nCellInPhi,
int  centralEta,
int  centralPhi,
int  centralZ,
MapType themap 
)
private

Definition at line 294 of file EcalBarrelSimHitsValidation.cc.

References funct::abs(), hcalRecHitTable_cff::ieta, hcalRecHitTable_cff::iphi, and DetId::rawId().

Referenced by analyze().

295  {
296  int ncristals = 0;
297  std::vector<uint32_t> ids(nCellInEta * nCellInPhi);
298 
299  int goBackInEta = nCellInEta / 2;
300  int goBackInPhi = nCellInPhi / 2;
301  int startEta = centralZ * centralEta - goBackInEta;
302  int startPhi = centralPhi - goBackInPhi;
303 
304  for (int ieta = startEta; ieta < startEta + nCellInEta; ieta++) {
305  for (int iphi = startPhi; iphi < startPhi + nCellInPhi; iphi++) {
306  uint32_t index;
307  if (abs(ieta) > 85 || abs(ieta) < 1) {
308  continue;
309  }
310  if (iphi < 1) {
311  index = EBDetId(ieta, iphi + 360).rawId();
312  } else if (iphi > 360) {
313  index = EBDetId(ieta, iphi - 360).rawId();
314  } else {
315  index = EBDetId(ieta, iphi).rawId();
316  }
317  ids[ncristals] = index;
318  ncristals += 1;
319  }
320  }
321 
322  return ids;
323 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57

◆ getUnitWithMaxEnergy()

uint32_t EcalBarrelSimHitsValidation::getUnitWithMaxEnergy ( MapType themap)
private

Definition at line 397 of file EcalBarrelSimHitsValidation.cc.

References LogDebug.

Referenced by analyze(), and fillEBMatrix().

397  {
398  // look for max
399  uint32_t unitWithMaxEnergy = 0;
400  float maxEnergy = 0.;
401 
402  MapType::iterator iter;
403  for (iter = themap.begin(); iter != themap.end(); iter++) {
404  if (maxEnergy < (*iter).second) {
405  maxEnergy = (*iter).second;
406  unitWithMaxEnergy = (*iter).first;
407  }
408  }
409 
410  LogDebug("GeomInfo") << " max energy of " << maxEnergy << " GeV in Unit id " << unitWithMaxEnergy;
411  return unitWithMaxEnergy;
412 }
#define LogDebug(id)

Member Data Documentation

◆ EBHitsCollection

std::string EcalBarrelSimHitsValidation::EBHitsCollection
private

Definition at line 61 of file EcalBarrelSimHitsValidation.h.

Referenced by EcalBarrelSimHitsValidation().

◆ EBHitsToken

edm::EDGetTokenT<edm::PCaloHitContainer> EcalBarrelSimHitsValidation::EBHitsToken
private

Definition at line 64 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and EcalBarrelSimHitsValidation().

◆ eRLength

float EcalBarrelSimHitsValidation::eRLength[26]
private

Definition at line 70 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and EcalBarrelSimHitsValidation().

◆ g4InfoLabel

std::string EcalBarrelSimHitsValidation::g4InfoLabel
private

Definition at line 60 of file EcalBarrelSimHitsValidation.h.

Referenced by EcalBarrelSimHitsValidation().

◆ meEBcrystalEnergy2_

MonitorElement* EcalBarrelSimHitsValidation::meEBcrystalEnergy2_
private

Definition at line 92 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBcrystalEnergy_

MonitorElement* EcalBarrelSimHitsValidation::meEBcrystalEnergy_
private

Definition at line 90 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe16_

MonitorElement* EcalBarrelSimHitsValidation::meEBe16_
private

Definition at line 97 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe16oe25_

MonitorElement* EcalBarrelSimHitsValidation::meEBe16oe25_
private

Definition at line 106 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe1_

MonitorElement* EcalBarrelSimHitsValidation::meEBe1_
private

Definition at line 94 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe1oe25_

MonitorElement* EcalBarrelSimHitsValidation::meEBe1oe25_
private

Definition at line 104 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe1oe4_

MonitorElement* EcalBarrelSimHitsValidation::meEBe1oe4_
private

Definition at line 100 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe1oe9_

MonitorElement* EcalBarrelSimHitsValidation::meEBe1oe9_
private

Definition at line 101 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe25_

MonitorElement* EcalBarrelSimHitsValidation::meEBe25_
private

Definition at line 98 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe4_

MonitorElement* EcalBarrelSimHitsValidation::meEBe4_
private

Definition at line 95 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe4oe9_

MonitorElement* EcalBarrelSimHitsValidation::meEBe4oe9_
private

Definition at line 102 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe9_

MonitorElement* EcalBarrelSimHitsValidation::meEBe9_
private

Definition at line 96 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe9oe16_

MonitorElement* EcalBarrelSimHitsValidation::meEBe9oe16_
private

Definition at line 103 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBe9oe25_

MonitorElement* EcalBarrelSimHitsValidation::meEBe9oe25_
private

Definition at line 105 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBhitEnergy2_

MonitorElement* EcalBarrelSimHitsValidation::meEBhitEnergy2_
private

Definition at line 88 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBhitEnergy_

MonitorElement* EcalBarrelSimHitsValidation::meEBhitEnergy_
private

Definition at line 80 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBhitLog10Energy25Norm_

MonitorElement* EcalBarrelSimHitsValidation::meEBhitLog10Energy25Norm_
private

Definition at line 86 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBhitLog10Energy_

MonitorElement* EcalBarrelSimHitsValidation::meEBhitLog10Energy_
private

Definition at line 82 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBhitLog10EnergyNorm_

MonitorElement* EcalBarrelSimHitsValidation::meEBhitLog10EnergyNorm_
private

Definition at line 84 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBLongitudinalShower_

MonitorElement* EcalBarrelSimHitsValidation::meEBLongitudinalShower_
private

Definition at line 78 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ meEBOccupancy_

MonitorElement* EcalBarrelSimHitsValidation::meEBOccupancy_
private

Definition at line 76 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ menEBCrystals_

MonitorElement* EcalBarrelSimHitsValidation::menEBCrystals_
private

Definition at line 74 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ menEBHits_

MonitorElement* EcalBarrelSimHitsValidation::menEBHits_
private

Definition at line 72 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and bookHistograms().

◆ myEntries

int EcalBarrelSimHitsValidation::myEntries
private

Definition at line 69 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and EcalBarrelSimHitsValidation().

◆ ValidationCollection

std::string EcalBarrelSimHitsValidation::ValidationCollection
private

Definition at line 62 of file EcalBarrelSimHitsValidation.h.

Referenced by EcalBarrelSimHitsValidation().

◆ ValidationCollectionToken

edm::EDGetTokenT<PEcalValidInfo> EcalBarrelSimHitsValidation::ValidationCollectionToken
private

Definition at line 65 of file EcalBarrelSimHitsValidation.h.

Referenced by analyze(), and EcalBarrelSimHitsValidation().

◆ verbose_

bool EcalBarrelSimHitsValidation::verbose_
private

Definition at line 67 of file EcalBarrelSimHitsValidation.h.

Referenced by EcalBarrelSimHitsValidation().