CMS 3D CMS Logo

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

Classes

struct  hitsinfo
 

Public Member Functions

 HGCalSimHitStudy (const edm::ParameterSet &)
 
 ~HGCalSimHitStudy () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::WatchRuns, edm::one::SharedResources >
 EDAnalyzer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const 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
 
bool wantsStreamRuns () const
 
 ~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
 
ESProxyIndex const * esGetTokenIndices (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::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)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices 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)
 

Protected Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
void beginJob () override
 
void beginRun (edm::Run const &, edm::EventSetup const &) override
 
void endRun (edm::Run const &, edm::EventSetup const &) override
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Private Member Functions

void analyzeHits (int, const std::string &, const std::vector< PCaloHit > &)
 

Private Attributes

const std::vector< std::string > caloHitSources_
 
const double etamax_
 
const double etamin_
 
std::vector< TH1D * > h_C1_
 
std::vector< TH1D * > h_C2_
 
std::vector< TH1D * > h_E_
 
std::vector< TH2D * > h_EtaPhi_
 
std::vector< TH2D * > h_EtFiZm_
 
std::vector< TH2D * > h_EtFiZp_
 
std::vector< TH1D * > h_LayerZm_
 
std::vector< TH1D * > h_LayerZp_
 
std::vector< TH1D * > h_Ly_
 
std::vector< TH2D * > h_RZ_
 
std::vector< TH1D * > h_T_
 
std::vector< TH1D * > h_W1_
 
std::vector< TH1D * > h_W2_
 
std::vector< TH2D * > h_XY_
 
const HcalDDDRecConstantshcons_
 
std::vector< bool > heRebuild_
 
std::vector< const HGCalDDDConstants * > hgcons_
 
const bool ifLayer_
 
const bool ifNose_
 
std::vector< int > layerFront_
 
std::vector< int > layers_
 
const std::vector< std::string > nameDetectors_
 
const int nbinEta_
 
const int nbinR_
 
const int nbinZ_
 
const int nLayers_
 
const double rmax_
 
const double rmin_
 
std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > tok_hits_
 
const int verbosity_
 
const double zmax_
 
const double zmin_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 

Detailed Description

Definition at line 49 of file HGCalSimHitStudy.cc.

Constructor & Destructor Documentation

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

Definition at line 91 of file HGCalSimHitStudy.cc.

References caloHitSources_, heRebuild_, TFileService::kSharedResource, Skims_PA_cff::name, nameDetectors_, source, and tok_hits_.

92  : nameDetectors_(iConfig.getParameter<std::vector<std::string> >("detectorNames")),
93  caloHitSources_(iConfig.getParameter<std::vector<std::string> >("caloHitSources")),
94  rmin_(iConfig.getUntrackedParameter<double>("rMin", 0.0)),
95  rmax_(iConfig.getUntrackedParameter<double>("rMax", 3000.0)),
96  zmin_(iConfig.getUntrackedParameter<double>("zMin", 3000.0)),
97  zmax_(iConfig.getUntrackedParameter<double>("zMax", 6000.0)),
98  etamin_(iConfig.getUntrackedParameter<double>("etaMin", 1.0)),
99  etamax_(iConfig.getUntrackedParameter<double>("etaMax", 3.0)),
100  nbinR_(iConfig.getUntrackedParameter<int>("nBinR", 300)),
101  nbinZ_(iConfig.getUntrackedParameter<int>("nBinZ", 300)),
102  nbinEta_(iConfig.getUntrackedParameter<int>("nBinEta", 200)),
103  nLayers_(iConfig.getUntrackedParameter<int>("layers", 50)),
104  verbosity_(iConfig.getUntrackedParameter<int>("verbosity", 0)),
105  ifNose_(iConfig.getUntrackedParameter<bool>("ifNose", false)),
106  ifLayer_(iConfig.getUntrackedParameter<bool>("ifLayer", false)) {
107  usesResource(TFileService::kSharedResource);
108 
109  for (auto const& name : nameDetectors_) {
110  if (name == "HCal")
111  heRebuild_.emplace_back(true);
112  else
113  heRebuild_.emplace_back(false);
114  }
115  for (auto const& source : caloHitSources_) {
116  tok_hits_.emplace_back(consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", source)));
117  }
118 }
static const std::string kSharedResource
Definition: TFileService.h:76
const std::vector< std::string > nameDetectors_
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > tok_hits_
const double etamin_
const double zmin_
const double rmin_
std::vector< bool > heRebuild_
const double etamax_
const double zmax_
const double rmax_
static std::string const source
Definition: EdmProvDump.cc:47
const std::vector< std::string > caloHitSources_
HGCalSimHitStudy::~HGCalSimHitStudy ( )
inlineoverride

Definition at line 61 of file HGCalSimHitStudy.cc.

References analyze(), beginJob(), beginRun(), and fillDescriptions().

61 {}

Member Function Documentation

void HGCalSimHitStudy::analyze ( edm::Event const &  ,
edm::EventSetup const &   
)
overrideprotectedvirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 142 of file HGCalSimHitStudy.cc.

References analyzeHits(), edm::Event::getByToken(), HcalEndcap, hcons_, heRebuild_, hit::id, edm::HandleBase::isValid(), dqmdumpme::k, nameDetectors_, edm::Handle< T >::product(), DetId::rawId(), HcalHitRelabeller::relabel(), HcalDetId::subdet(), tok_hits_, and verbosity_.

Referenced by ~HGCalSimHitStudy().

142  {
143  //Now the hits
144  for (unsigned int k = 0; k < tok_hits_.size(); ++k) {
145  edm::Handle<edm::PCaloHitContainer> theCaloHitContainers;
146  iEvent.getByToken(tok_hits_[k], theCaloHitContainers);
147  if (theCaloHitContainers.isValid()) {
148  if (verbosity_ > 0)
149  edm::LogVerbatim("HGCalValidation") << " PcalohitItr = " << theCaloHitContainers->size();
150  std::vector<PCaloHit> caloHits;
151  if (heRebuild_[k]) {
152  for (auto const& hit : *(theCaloHitContainers.product())) {
153  unsigned int id = hit.id();
155  if (hid.subdet() != static_cast<int>(HcalEndcap)) {
156  caloHits.emplace_back(hit);
157  caloHits.back().setID(hid.rawId());
158  if (verbosity_ > 0)
159  edm::LogVerbatim("HGCalValidation") << "Hit[" << caloHits.size() << "] " << hid;
160  }
161  }
162  } else {
163  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
164  }
165  analyzeHits(k, nameDetectors_[k], caloHits);
166  } else if (verbosity_ > 0) {
167  edm::LogVerbatim("HGCalValidation") << "PCaloHitContainer does not "
168  << "exist for " << nameDetectors_[k];
169  }
170  }
171 }
const std::vector< std::string > nameDetectors_
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > tok_hits_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
int iEvent
Definition: GenABIO.cc:224
void analyzeHits(int, const std::string &, const std::vector< PCaloHit > &)
std::vector< bool > heRebuild_
bool isValid() const
Definition: HandleBase.h:70
unsigned int id
const HcalDDDRecConstants * hcons_
T const * product() const
Definition: Handle.h:69
DetId relabel(const uint32_t testId) const
void HGCalSimHitStudy::analyzeHits ( int  ,
const std::string &  ,
const std::vector< PCaloHit > &   
)
private

Definition at line 173 of file HGCalSimHitStudy.cc.

References funct::abs(), hgcalTopologyTester_cfi::cell2, HFNoseDetId::cellU(), HGCSiliconDetId::cellU(), HFNoseDetId::cellV(), HGCSiliconDetId::cellV(), funct::cos(), HcalDetId::depth(), DetId::det(), HCALHighEnergyHPDFilter_cfi::energy, HGCalSimHitStudy::hitsinfo::energy, HGCalSimHitStudy::hitsinfo::eta, fastmath::etaphi(), HcalDDDRecConstants::getEtaPhi(), HcalDDDRecConstants::getRZ(), h_C1_, h_C2_, h_E_, h_EtaPhi_, h_EtFiZm_, h_EtFiZp_, h_LayerZm_, h_LayerZp_, h_Ly_, h_RZ_, h_T_, h_W1_, h_W2_, h_XY_, hcons_, heRebuild_, HGCalGeometryMode::Hexagon8, HGCalGeometryMode::Hexagon8Full, hgcons_, triggerObjects_cff::id, hit::id, HGCScintillatorDetId::ieta(), HcalDetId::ietaAbs(), ifLayer_, ifNose_, HGCScintillatorDetId::iphi(), HcalDetId::iphi(), HGCScintillatorDetId::layer(), HFNoseDetId::layer(), HGCalSimHitStudy::hitsinfo::layer, HGCSiliconDetId::layer(), layerFront_, HGCalSimHitStudy::hitsinfo::phi, hgcalTopologyTester_cfi::sector2, funct::sin(), HcalDetId::subdet(), DetId::subdetId(), HGCalSimHitStudy::hitsinfo::time, ntuplemaker::time, HGCalGeometryMode::Trapezoid, HGCScintillatorDetId::type(), HFNoseDetId::type(), HGCSiliconDetId::type(), HGCalTestNumbering::unpackHexagonIndex(), verbosity_, HFNoseDetId::waferU(), HGCSiliconDetId::waferU(), HFNoseDetId::waferV(), HGCSiliconDetId::waferV(), geometryCSVtoXML::xy, HGCScintillatorDetId::zside(), HFNoseDetId::zside(), HGCSiliconDetId::zside(), HcalDetId::zside(), and ecaldqm::zside().

Referenced by analyze(), and endRun().

173  {
174  if (verbosity_ > 0)
175  edm::LogVerbatim("HGCalValidation") << name << " with " << hits.size() << " PcaloHit elements";
176 
177  std::map<uint32_t, hitsinfo> map_hits;
178  map_hits.clear();
179 
180  unsigned int nused(0);
181  for (auto const& hit : hits) {
182  double energy = hit.energy();
183  double time = hit.time();
184  uint32_t id = hit.id();
185  int cell, sector, sector2(0), layer, zside;
186  int subdet(0), cell2(0), type(0);
187  HepGeom::Point3D<float> gcoord;
188  if (heRebuild_[ih]) {
189  HcalDetId detId = HcalDetId(id);
190  subdet = detId.subdet();
191  cell = detId.ietaAbs();
192  sector = detId.iphi();
193  layer = detId.depth();
194  zside = detId.zside();
195  std::pair<double, double> etaphi = hcons_->getEtaPhi(subdet, zside * cell, sector);
196  double rz = hcons_->getRZ(subdet, zside * cell, layer);
197  if (verbosity_ > 2)
198  edm::LogVerbatim("HGCalValidation") << "i/p " << subdet << ":" << zside << ":" << cell << ":" << sector << ":"
199  << layer << " o/p " << etaphi.first << ":" << etaphi.second << ":" << rz;
200  gcoord = HepGeom::Point3D<float>(rz * cos(etaphi.second) / cosh(etaphi.first),
201  rz * sin(etaphi.second) / cosh(etaphi.first),
202  rz * tanh(etaphi.first));
203  } else {
204  std::pair<float, float> xy;
205  if (ifNose_) {
206  HFNoseDetId detId = HFNoseDetId(id);
207  subdet = detId.subdetId();
208  cell = detId.cellU();
209  cell2 = detId.cellV();
210  sector = detId.waferU();
211  sector2 = detId.waferV();
212  type = detId.type();
213  layer = detId.layer();
214  zside = detId.zside();
215  xy = hgcons_[ih]->locateCell(layer, sector, sector2, cell, cell2, false, true);
216  h_W2_[ih]->Fill(sector2);
217  h_C2_[ih]->Fill(cell2);
218 
219  } else if ((hgcons_[ih]->geomMode() == HGCalGeometryMode::Hexagon8) ||
220  (hgcons_[ih]->geomMode() == HGCalGeometryMode::Hexagon8Full)) {
221  HGCSiliconDetId detId = HGCSiliconDetId(id);
222  subdet = static_cast<int>(detId.det());
223  cell = detId.cellU();
224  cell2 = detId.cellV();
225  sector = detId.waferU();
226  sector2 = detId.waferV();
227  type = detId.type();
228  layer = detId.layer();
229  zside = detId.zside();
230  xy = hgcons_[ih]->locateCell(layer, sector, sector2, cell, cell2, false, true);
231  h_W2_[ih]->Fill(sector2);
232  h_C2_[ih]->Fill(cell2);
233  } else if (hgcons_[ih]->geomMode() == HGCalGeometryMode::Trapezoid) {
235  subdet = static_cast<int>(detId.det());
236  sector = detId.ieta();
237  cell = detId.iphi();
238  type = detId.type();
239  layer = detId.layer();
240  zside = detId.zside();
241  xy = hgcons_[ih]->locateCellTrap(layer, sector, cell, false);
242  } else {
243  HGCalTestNumbering::unpackHexagonIndex(id, subdet, zside, layer, sector, type, cell);
244  xy = hgcons_[ih]->locateCell(cell, layer, sector, false);
245  }
246  double zp = hgcons_[ih]->waferZ(layer, false);
247  if (zside < 0)
248  zp = -zp;
249  double xp = (zp < 0) ? -xy.first : xy.first;
250  gcoord = HepGeom::Point3D<float>(xp, xy.second, zp);
251  if (verbosity_ > 2)
252  edm::LogVerbatim("HGCalValidation")
253  << "i/p " << subdet << ":" << zside << ":" << layer << ":" << sector << ":" << sector2 << ":" << cell << ":"
254  << cell2 << " o/p " << xy.first << ":" << xy.second << ":" << zp;
255  }
256  nused++;
257  double tof = (gcoord.mag() * CLHEP::mm) / CLHEP::c_light;
258  if (verbosity_ > 1)
259  edm::LogVerbatim("HGCalValidation")
260  << "Detector " << name << " zside = " << zside << " layer = " << layer << " type = " << type
261  << " wafer = " << sector << ":" << sector2 << " cell = " << cell << ":" << cell2 << " positon = " << gcoord
262  << " energy = " << energy << " time = " << time << ":" << tof;
263  time -= tof;
264  if (time < 0)
265  time = 0;
266  hitsinfo hinfo;
267  if (map_hits.count(id) != 0) {
268  hinfo = map_hits[id];
269  } else {
270  hinfo.layer = layer + layerFront_[ih];
271  hinfo.phi = gcoord.getPhi();
272  hinfo.eta = gcoord.getEta();
273  hinfo.time = time;
274  }
275  hinfo.energy += energy;
276  map_hits[id] = hinfo;
277 
278  //Fill in histograms
279  h_RZ_[0]->Fill(std::abs(gcoord.z()), gcoord.rho());
280  h_RZ_[ih + 1]->Fill(std::abs(gcoord.z()), gcoord.rho());
281  if (ifLayer_) {
282  if (hinfo.layer <= static_cast<int>(h_XY_.size()))
283  h_XY_[hinfo.layer-1]->Fill(gcoord.x(), gcoord.y());
284  } else {
285  h_EtaPhi_[0]->Fill(std::abs(hinfo.eta), hinfo.phi);
286  h_EtaPhi_[ih + 1]->Fill(std::abs(hinfo.eta), hinfo.phi);
287  }
288  h_Ly_[ih]->Fill(layer);
289  h_W1_[ih]->Fill(sector);
290  h_C1_[ih]->Fill(cell);
291  }
292  if (verbosity_ > 0)
293  edm::LogVerbatim("HGCalValidation") << name << " with " << map_hits.size() << ":" << nused << " detector elements"
294  << " being hit";
295 
296  for (auto const& hit : map_hits) {
297  hitsinfo hinfo = hit.second;
298  if (verbosity_ > 1)
299  edm::LogVerbatim("HGCalValidation")
300  << " ---------------------- eta = " << hinfo.eta << " phi = " << hinfo.phi << " layer = " << hinfo.layer
301  << " E = " << hinfo.energy << " T = " << hinfo.time;
302  h_E_[0]->Fill(hinfo.energy);
303  h_E_[ih + 1]->Fill(hinfo.energy);
304  h_T_[0]->Fill(hinfo.time);
305  h_T_[ih + 1]->Fill(hinfo.time);
306  if (hinfo.eta > 0) {
307  if (!ifLayer_) {
308  h_EtFiZp_[0]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
309  h_EtFiZp_[ih + 1]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
310  }
311  h_LayerZp_[0]->Fill(hinfo.layer, hinfo.energy);
312  h_LayerZp_[ih + 1]->Fill(hinfo.layer, hinfo.energy);
313  } else {
314  if (!ifLayer_) {
315  h_EtFiZm_[0]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
316  h_EtFiZm_[ih + 1]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
317  }
318  h_LayerZm_[0]->Fill(hinfo.layer, hinfo.energy);
319  h_LayerZm_[ih + 1]->Fill(hinfo.layer, hinfo.energy);
320  }
321  }
322 }
std::pair< T, T > etaphi(T x, T y, T z)
Definition: FastMath.h:162
type
Definition: HCALResponse.h:21
std::vector< TH2D * > h_XY_
std::vector< TH1D * > h_LayerZm_
std::vector< TH2D * > h_RZ_
std::vector< TH1D * > h_T_
std::vector< TH1D * > h_W1_
int zside() const
get the z-side of the cell (1/-1)
Definition: HFNoseDetId.h:52
std::pair< double, double > getEtaPhi(const int &subdet, const int &ieta, const int &iphi) const
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
std::vector< TH2D * > h_EtFiZm_
std::vector< int > layerFront_
std::vector< TH1D * > h_C1_
int type() const
get the type
Definition: HFNoseDetId.h:49
int cellU() const
get the cell #&#39;s in u,v or in x,y
Definition: HFNoseDetId.h:58
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:141
int waferU() const
int cellV() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
int type() const
get the type
int zside() const
get the z-side of the cell (1/-1)
int zside(DetId const &)
int cellU() const
get the cell #&#39;s in u,v or in x,y
int depth() const
get the tower depth
Definition: HcalDetId.h:164
int cellV() const
Definition: HFNoseDetId.h:59
int waferV() const
Definition: HFNoseDetId.h:77
int layer() const
get the layer #
Definition: HFNoseDetId.h:55
std::vector< TH1D * > h_LayerZp_
int type() const
get the type
int waferU() const
Definition: HFNoseDetId.h:74
int layer() const
get the layer #
double getRZ(const int &subdet, const int &ieta, const int &depth) const
int waferV() const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< TH2D * > h_EtFiZp_
int iphi() const
get the phi index
std::vector< TH1D * > h_C2_
std::vector< bool > heRebuild_
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
unsigned int id
const HcalDDDRecConstants * hcons_
int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
int layer() const
get the layer #
int zside() const
get the z-side of the cell (1/-1)
std::vector< TH1D * > h_Ly_
std::vector< TH1D * > h_E_
std::vector< const HGCalDDDConstants * > hgcons_
std::vector< TH1D * > h_W2_
std::vector< TH2D * > h_EtaPhi_
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
void HGCalSimHitStudy::beginJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 350 of file HGCalSimHitStudy.cc.

References DEFINE_FWK_MODULE, etamax_, etamin_, h_C1_, h_C2_, h_E_, h_EtaPhi_, h_EtFiZm_, h_EtFiZp_, h_LayerZm_, h_LayerZp_, h_Ly_, h_RZ_, h_T_, h_W1_, h_W2_, h_XY_, heRebuild_, ifLayer_, M_PI, TFileService::make(), Skims_PA_cff::name, nameDetectors_, nbinEta_, nbinR_, nbinZ_, nLayers_, rmax_, rmin_, overlapproblemtsosanalyzer_cfi::title, zmax_, and zmin_.

Referenced by ~HGCalSimHitStudy().

350  {
352 
353  std::ostringstream name, title;
354  for (unsigned int ih = 0; ih <= nameDetectors_.size(); ++ih) {
355  name.str("");
356  title.str("");
357  if (ih == 0) {
358  name << "RZ_AllDetectors";
359  title << "R vs Z for All Detectors";
360  } else {
361  name << "RZ_" << nameDetectors_[ih - 1];
362  title << "R vs Z for " << nameDetectors_[ih - 1];
363  }
364  h_RZ_.emplace_back(
365  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinZ_, zmin_, zmax_, nbinR_, rmin_, rmax_));
366  if (ifLayer_) {
367  if (ih == 0) {
368  for (int ly = 0; ly < nLayers_; ++ly) {
369  name.str("");
370  title.str("");
371  name << "XY_L" << (ly + 1);
372  title << "Y vs X at Layer " << (ly + 1);
373  h_XY_.emplace_back(
374  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinR_, -rmax_, rmax_, nbinR_, -rmax_, rmax_));
375  }
376  }
377  } else {
378  name.str("");
379  title.str("");
380  if (ih == 0) {
381  name << "EtaPhi_AllDetectors";
382  title << "#phi vs #eta for All Detectors";
383  } else {
384  name << "EtaPhi_" << nameDetectors_[ih - 1];
385  title << "#phi vs #eta for " << nameDetectors_[ih - 1];
386  }
387  h_EtaPhi_.emplace_back(
388  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinEta_, etamin_, etamax_, 200, -M_PI, M_PI));
389  name.str("");
390  title.str("");
391  if (ih == 0) {
392  name << "EtFiZp_AllDetectors";
393  title << "#phi vs #eta (+z) for All Detectors";
394  } else {
395  name << "EtFiZp_" << nameDetectors_[ih - 1];
396  title << "#phi vs #eta (+z) for " << nameDetectors_[ih - 1];
397  }
398  h_EtFiZp_.emplace_back(
399  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinEta_, etamin_, etamax_, 200, -M_PI, M_PI));
400  name.str("");
401  title.str("");
402  if (ih == 0) {
403  name << "EtFiZm_AllDetectors";
404  title << "#phi vs #eta (-z) for All Detectors";
405  } else {
406  name << "EtFiZm_" << nameDetectors_[ih - 1];
407  title << "#phi vs #eta (-z) for " << nameDetectors_[ih - 1];
408  }
409  h_EtFiZm_.emplace_back(
410  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinEta_, etamin_, etamax_, 200, -M_PI, M_PI));
411  }
412  name.str("");
413  title.str("");
414  if (ih == 0) {
415  name << "LayerZp_AllDetectors";
416  title << "Energy vs Layer (+z) for All Detectors";
417  } else {
418  name << "LayerZp_" << nameDetectors_[ih - 1];
419  title << "Energy vs Layer (+z) for " << nameDetectors_[ih - 1];
420  }
421  h_LayerZp_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 60, 0.0, 60.0));
422  name.str("");
423  title.str("");
424  if (ih == 0) {
425  name << "LayerZm_AllDetectors";
426  title << "Energy vs Layer (-z) for All Detectors";
427  } else {
428  name << "LayerZm_" << nameDetectors_[ih - 1];
429  title << "Energy vs Layer (-z) for " << nameDetectors_[ih - 1];
430  }
431  h_LayerZm_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 60, 0.0, 60.0));
432 
433  name.str("");
434  title.str("");
435  if (ih == 0) {
436  name << "E_AllDetectors";
437  title << "Energy Deposit for All Detectors";
438  } else {
439  name << "E_" << nameDetectors_[ih - 1];
440  title << "Energy Deposit for " << nameDetectors_[ih - 1];
441  }
442  h_E_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 1000, 0.0, 1.0));
443 
444  name.str("");
445  title.str("");
446  if (ih == 0) {
447  name << "T_AllDetectors";
448  title << "Time for All Detectors";
449  } else {
450  name << "T_" << nameDetectors_[ih - 1];
451  title << "Time for " << nameDetectors_[ih - 1];
452  }
453  h_T_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 1000, 0.0, 200.0));
454  }
455 
456  for (unsigned int ih = 0; ih < nameDetectors_.size(); ++ih) {
457  name.str("");
458  title.str("");
459  name << "LY_" << nameDetectors_[ih];
460  title << "Layer number for " << nameDetectors_[ih];
461  h_Ly_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 200, 0, 100));
462  if ((nameDetectors_[ih] == "HGCalHEScintillatorSensitive") || heRebuild_[ih]) {
463  name.str("");
464  title.str("");
465  name << "IR_" << nameDetectors_[ih];
466  title << "Radius index for " << nameDetectors_[ih];
467  h_W1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 200, -50, 50));
468  name.str("");
469  title.str("");
470  name << "FI_" << nameDetectors_[ih];
471  title << "#phi index for " << nameDetectors_[ih];
472  h_C1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 720, 0, 360));
473  } else {
474  name.str("");
475  title.str("");
476  name << "WU_" << nameDetectors_[ih];
477  title << "u index of wafers for " << nameDetectors_[ih];
478  h_W1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 200, -50, 50));
479  name.str("");
480  title.str("");
481  name << "WV_" << nameDetectors_[ih];
482  title << "v index of wafers for " << nameDetectors_[ih];
483  h_W2_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 100, -50, 50));
484  name.str("");
485  title.str("");
486  name << "CU_" << nameDetectors_[ih];
487  title << "u index of cells for " << nameDetectors_[ih];
488  h_C1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 100, 0, 50));
489  name.str("");
490  title.str("");
491  name << "CV_" << nameDetectors_[ih];
492  title << "v index of cells for " << nameDetectors_[ih];
493  h_C2_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 100, 0, 50));
494  }
495  }
496 }
const std::vector< std::string > nameDetectors_
std::vector< TH2D * > h_XY_
std::vector< TH1D * > h_LayerZm_
std::vector< TH2D * > h_RZ_
std::vector< TH1D * > h_T_
std::vector< TH1D * > h_W1_
std::vector< TH2D * > h_EtFiZm_
std::vector< TH1D * > h_C1_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
const double etamin_
const double zmin_
std::vector< TH1D * > h_LayerZp_
std::vector< TH2D * > h_EtFiZp_
std::vector< TH1D * > h_C2_
const double rmin_
std::vector< bool > heRebuild_
#define M_PI
const double etamax_
const double zmax_
const double rmax_
std::vector< TH1D * > h_Ly_
std::vector< TH1D * > h_E_
std::vector< TH1D * > h_W2_
std::vector< TH2D * > h_EtaPhi_
void HGCalSimHitStudy::beginRun ( edm::Run const &  ,
edm::EventSetup const &   
)
overrideprotected

Definition at line 325 of file HGCalSimHitStudy.cc.

References edm::EventSetup::get(), HcalDDDRecConstants::getMaxDepth(), hcons_, heRebuild_, hgcons_, dqmdumpme::k, layerFront_, layers_, nameDetectors_, and verbosity_.

Referenced by ~HGCalSimHitStudy().

325  {
326  for (unsigned int k = 0; k < nameDetectors_.size(); ++k) {
327  if (heRebuild_[k]) {
329  iSetup.get<HcalRecNumberingRecord>().get(pHRNDC);
330  hcons_ = &(*pHRNDC);
331  layers_.emplace_back(hcons_->getMaxDepth(1));
332  hgcons_.emplace_back(nullptr);
333  layerFront_.emplace_back(40);
334  } else {
336  iSetup.get<IdealGeometryRecord>().get(nameDetectors_[k], pHGDC);
337  hgcons_.emplace_back(&(*pHGDC));
338  layers_.emplace_back(hgcons_.back()->layers(false));
339  if (k == 0)
340  layerFront_.emplace_back(0);
341  else
342  layerFront_.emplace_back(28);
343  }
344  if (verbosity_ > 0)
345  edm::LogVerbatim("HGCalValidation") << nameDetectors_[k] << " defined with " << layers_.back() << " Layers with "
346  << layerFront_.back() << " in front";
347  }
348 }
const std::vector< std::string > nameDetectors_
std::vector< int > layerFront_
std::vector< bool > heRebuild_
std::vector< int > layers_
const HcalDDDRecConstants * hcons_
int getMaxDepth(const int &type) const
std::vector< const HGCalDDDConstants * > hgcons_
void HGCalSimHitStudy::endRun ( edm::Run const &  ,
edm::EventSetup const &   
)
inlineoverrideprotected

Definition at line 68 of file HGCalSimHitStudy.cc.

References analyzeHits(), and AlCaHLTBitMon_QueryRunRegistry::string.

68 {}
void HGCalSimHitStudy::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 120 of file HGCalSimHitStudy.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), names, and CalibrationSummaryClient_cfi::sources.

Referenced by ~HGCalSimHitStudy().

120  {
122  std::vector<std::string> names = {"HGCalEESensitive", "HGCalHESiliconSensitive", "Hcal"};
123  std::vector<std::string> sources = {"HGCHitsEE", "HGCHitsHEfront", "HcalHits"};
124  desc.add<std::vector<std::string> >("detectorNames", names);
125  desc.add<std::vector<std::string> >("caloHitSources", sources);
126  desc.addUntracked<double>("rMin", 0.0);
127  desc.addUntracked<double>("rMax", 3000.0);
128  desc.addUntracked<double>("zMin", 3000.0);
129  desc.addUntracked<double>("zMax", 6000.0);
130  desc.addUntracked<double>("etaMin", 1.0);
131  desc.addUntracked<double>("etaMax", 3.0);
132  desc.addUntracked<int>("nBinR", 300);
133  desc.addUntracked<int>("nBinZ", 300);
134  desc.addUntracked<int>("nBinEta", 200);
135  desc.addUntracked<int>("layers", 50);
136  desc.addUntracked<int>("verbosity", 0);
137  desc.addUntracked<bool>("ifNose", false);
138  desc.addUntracked<bool>("ifLayer", false);
139  descriptions.add("hgcalSimHitStudy", desc);
140 }
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const std::string names[nVars_]
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

Member Data Documentation

const std::vector<std::string> HGCalSimHitStudy::caloHitSources_
private

Definition at line 74 of file HGCalSimHitStudy.cc.

Referenced by HGCalSimHitStudy().

const double HGCalSimHitStudy::etamax_
private

Definition at line 76 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

const double HGCalSimHitStudy::etamin_
private

Definition at line 76 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_C1_
private

Definition at line 88 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_C2_
private

Definition at line 88 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_E_
private

Definition at line 87 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH2D*> HGCalSimHitStudy::h_EtaPhi_
private

Definition at line 86 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH2D*> HGCalSimHitStudy::h_EtFiZm_
private

Definition at line 86 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH2D*> HGCalSimHitStudy::h_EtFiZp_
private

Definition at line 86 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_LayerZm_
private

Definition at line 87 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_LayerZp_
private

Definition at line 87 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_Ly_
private

Definition at line 88 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH2D*> HGCalSimHitStudy::h_RZ_
private

Definition at line 86 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_T_
private

Definition at line 87 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_W1_
private

Definition at line 88 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH1D*> HGCalSimHitStudy::h_W2_
private

Definition at line 88 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

std::vector<TH2D*> HGCalSimHitStudy::h_XY_
private

Definition at line 86 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

const HcalDDDRecConstants* HGCalSimHitStudy::hcons_
private

Definition at line 80 of file HGCalSimHitStudy.cc.

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

std::vector<bool> HGCalSimHitStudy::heRebuild_
private

Definition at line 81 of file HGCalSimHitStudy.cc.

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

std::vector<const HGCalDDDConstants*> HGCalSimHitStudy::hgcons_
private

Definition at line 79 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginRun().

const bool HGCalSimHitStudy::ifLayer_
private

Definition at line 78 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginJob().

const bool HGCalSimHitStudy::ifNose_
private

Definition at line 78 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits().

std::vector<int> HGCalSimHitStudy::layerFront_
private

Definition at line 83 of file HGCalSimHitStudy.cc.

Referenced by analyzeHits(), and beginRun().

std::vector<int> HGCalSimHitStudy::layers_
private

Definition at line 83 of file HGCalSimHitStudy.cc.

Referenced by beginRun().

const std::vector<std::string> HGCalSimHitStudy::nameDetectors_
private

Definition at line 74 of file HGCalSimHitStudy.cc.

Referenced by analyze(), beginJob(), beginRun(), and HGCalSimHitStudy().

const int HGCalSimHitStudy::nbinEta_
private

Definition at line 77 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

const int HGCalSimHitStudy::nbinR_
private

Definition at line 77 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

const int HGCalSimHitStudy::nbinZ_
private

Definition at line 77 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

const int HGCalSimHitStudy::nLayers_
private

Definition at line 77 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

const double HGCalSimHitStudy::rmax_
private

Definition at line 75 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

const double HGCalSimHitStudy::rmin_
private

Definition at line 75 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

std::vector<edm::EDGetTokenT<edm::PCaloHitContainer> > HGCalSimHitStudy::tok_hits_
private

Definition at line 82 of file HGCalSimHitStudy.cc.

Referenced by analyze(), and HGCalSimHitStudy().

const int HGCalSimHitStudy::verbosity_
private

Definition at line 77 of file HGCalSimHitStudy.cc.

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

const double HGCalSimHitStudy::zmax_
private

Definition at line 75 of file HGCalSimHitStudy.cc.

Referenced by beginJob().

const double HGCalSimHitStudy::zmin_
private

Definition at line 75 of file HGCalSimHitStudy.cc.

Referenced by beginJob().