CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
HGCalHitCalibration Class Reference
Inheritance diagram for HGCalHitCalibration:
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

 HGCalHitCalibration (const edm::ParameterSet &)
 
 ~HGCalHitCalibration () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &ev, edm::EventSetup const &es) final
 
virtual void analyze (edm::Event const &, edm::EventSetup const &)
 
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)
 

Static Public Member Functions

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

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void fillWithRecHits (std::map< DetId, const HGCRecHit * > &, DetId, unsigned int, float, int &, float &)
 

Private Attributes

int algo_
 
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
 
int debug_
 
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electrons_
 
std::array< float, layers_Energy_layer_calib_
 
std::array< float, layers_Energy_layer_calib_fraction_
 
std::map< int, MonitorElement * > h_EoP_CPene_calib_fraction_
 
std::map< int, MonitorElement * > hgcal_ele_EoP_CPene_calib_fraction_
 
std::map< int, MonitorElement * > hgcal_EoP_CPene_calib_fraction_
 
std::map< int, MonitorElement * > hgcal_photon_EoP_CPene_calib_fraction_
 
edm::EDGetTokenT< std::vector< reco::PFCluster > > hgcalMultiClusters_
 
MonitorElementLayerOccupancy_
 
edm::EDGetTokenT< std::vector< reco::Photon > > photons_
 
bool rawRecHits_
 
edm::EDGetTokenT< HGCRecHitCollectionrecHitsBH_
 
edm::EDGetTokenT< HGCRecHitCollectionrecHitsEE_
 
edm::EDGetTokenT< HGCRecHitCollectionrecHitsFH_
 
hgcal::RecHitTools recHitTools_
 

Static Private Attributes

static int layers_ = 60
 

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
 
- 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 HGCalHitCalibration.cc.

Constructor & Destructor Documentation

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

Definition at line 82 of file HGCalHitCalibration.cc.

References algo_, caloTruthProducer_cfi::caloParticles, caloParticles_, gamEcalExtractorBlocks_cff::detector, electrons_cff::electrons, electrons_, Energy_layer_calib_, Energy_layer_calib_fraction_, edm::ParameterSet::getParameter(), hgcalMultiClusters_, nano_cff::photons, photons_, recHitsBH_, EcalDeadCellBoundaryEnergyFilter_cfi::recHitsEE, recHitsEE_, recHitsFH_, and AlCaHLTBitMon_QueryRunRegistry::string.

83  : rawRecHits_(iConfig.getParameter<bool>("rawRecHits")),
84  debug_(iConfig.getParameter<int>("debug")) {
85  auto detector = iConfig.getParameter<std::string>("detector");
86  auto recHitsEE = iConfig.getParameter<edm::InputTag>("recHitsEE");
87  auto recHitsFH = iConfig.getParameter<edm::InputTag>("recHitsFH");
88  auto recHitsBH = iConfig.getParameter<edm::InputTag>("recHitsBH");
89  auto caloParticles = iConfig.getParameter<edm::InputTag>("caloParticles");
90  auto hgcalMultiClusters = iConfig.getParameter<edm::InputTag>("hgcalMultiClusters");
91  auto electrons = iConfig.getParameter<edm::InputTag>("electrons");
92  auto photons = iConfig.getParameter<edm::InputTag>("photons");
93  if (detector == "all") {
94  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
95  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
96  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
97  algo_ = 1;
98  } else if (detector == "EM") {
99  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
100  algo_ = 2;
101  } else if (detector == "HAD") {
102  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
103  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
104  algo_ = 3;
105  }
106  caloParticles_ = consumes<std::vector<CaloParticle> >(caloParticles);
107  hgcalMultiClusters_ = consumes<std::vector<reco::PFCluster> >(hgcalMultiClusters);
108  electrons_ = consumes<std::vector<reco::GsfElectron> >(electrons);
109  photons_ = consumes<std::vector<reco::Photon> >(photons);
110 
111  // size should be HGC layers 52 is enough
112  Energy_layer_calib_.fill(0.);
114 }
T getParameter(std::string const &) const
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
edm::EDGetTokenT< std::vector< reco::Photon > > photons_
std::array< float, layers_ > Energy_layer_calib_fraction_
edm::EDGetTokenT< std::vector< reco::PFCluster > > hgcalMultiClusters_
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electrons_
std::array< float, layers_ > Energy_layer_calib_
HGCalHitCalibration::~HGCalHitCalibration ( )
override

Definition at line 116 of file HGCalHitCalibration.cc.

116  {
117  // do anything here that needs to be done at desctruction time
118  // (e.g. close files, deallocate resources etc.)
119 }

Member Function Documentation

void HGCalHitCalibration::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 182 of file HGCalHitCalibration.cc.

References a, algo_, b, caloTruthProducer_cfi::caloParticles, caloParticles_, fastPrimaryVertexProducer_cfi::clusters, reco::CaloCluster::correctedEnergy(), debug_, reco::deltaR2(), DetId::det(), electrons_, reco::PFCluster::energy(), SimCluster::energy(), Energy_layer_calib_, Energy_layer_calib_fraction_, HcalObjRepresent::Fill(), fillWithRecHits(), DetId::Forward, edm::Event::getByToken(), hgcal::RecHitTools::getEventSetup(), hgcal::RecHitTools::getLayerWithOffset(), hgcal::RecHitTools::getSiThickness(), h_EoP_CPene_calib_fraction_, DetId::Hcal, HcalEndcap, hgcal_ele_EoP_CPene_calib_fraction_, hgcal_EoP_CPene_calib_fraction_, hgcal_photon_EoP_CPene_calib_fraction_, hgcalMultiClusters_, HGCEE, HGCHEB, HGCHEF, SimCluster::hits_and_fractions(), mps_fire::i, IfLogTrace, edm::HandleBase::isValid(), muons2muons_cfi::photon, photons_, rawRecHits_, recHitsBH_, recHitsEE_, recHitsFH_, recHitTools_, SimCluster::simEnergy(), edm::SortedCollection< T, SORT >::size(), and DetId::subdetId().

183  {
184  using namespace edm;
185 
186  recHitTools_.getEventSetup(iSetup);
187 
188  Handle<HGCRecHitCollection> recHitHandleEE;
189  Handle<HGCRecHitCollection> recHitHandleFH;
190  Handle<HGCRecHitCollection> recHitHandleBH;
191 
192  Handle<std::vector<CaloParticle> > caloParticleHandle;
193  iEvent.getByToken(caloParticles_, caloParticleHandle);
194  const std::vector<CaloParticle>& caloParticles = *caloParticleHandle;
195 
196  Handle<std::vector<reco::PFCluster> > hgcalMultiClustersHandle;
197  iEvent.getByToken(hgcalMultiClusters_, hgcalMultiClustersHandle);
198 
199  Handle<std::vector<reco::GsfElectron> > PFElectronHandle;
200  iEvent.getByToken(electrons_, PFElectronHandle);
201 
202  Handle<std::vector<reco::Photon> > PFPhotonHandle;
203  iEvent.getByToken(photons_, PFPhotonHandle);
204 
205  // make a map detid-rechit
206  std::map<DetId, const HGCRecHit*> hitmap;
207  switch (algo_) {
208  case 1: {
209  iEvent.getByToken(recHitsEE_, recHitHandleEE);
210  iEvent.getByToken(recHitsFH_, recHitHandleFH);
211  iEvent.getByToken(recHitsBH_, recHitHandleBH);
212  const auto& rechitsEE = *recHitHandleEE;
213  const auto& rechitsFH = *recHitHandleFH;
214  const auto& rechitsBH = *recHitHandleBH;
215  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
216  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
217  }
218  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
219  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
220  }
221  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
222  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
223  }
224  break;
225  }
226  case 2: {
227  iEvent.getByToken(recHitsEE_, recHitHandleEE);
228  const HGCRecHitCollection& rechitsEE = *recHitHandleEE;
229  for (unsigned int i = 0; i < rechitsEE.size(); i++) {
230  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
231  }
232  break;
233  }
234  case 3: {
235  iEvent.getByToken(recHitsFH_, recHitHandleFH);
236  iEvent.getByToken(recHitsBH_, recHitHandleBH);
237  const auto& rechitsFH = *recHitHandleFH;
238  const auto& rechitsBH = *recHitHandleBH;
239  for (unsigned int i = 0; i < rechitsFH.size(); i++) {
240  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
241  }
242  for (unsigned int i = 0; i < rechitsBH.size(); i++) {
243  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
244  }
245  break;
246  }
247  default:
248  assert(false);
249  break;
250  }
251 
252  // loop over caloParticles
253  int seedDet = 0;
254  float seedEnergy = 0.;
255  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
256  << "Number of caloParticles: " << caloParticles.size() << std::endl;
257  for (const auto& it_caloPart : caloParticles) {
258  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
259  Energy_layer_calib_.fill(0.);
261 
262  seedDet = 0;
263  seedEnergy = 0.;
264  for (const auto& it_sc : simClusterRefVector) {
265  const SimCluster& simCluster = (*(it_sc));
266  IfLogTrace(debug_ > 1, "HGCalHitCalibration")
267  << ">>> SC.energy(): " << simCluster.energy()
268  << " SC.simEnergy(): " << simCluster.simEnergy()
269  << std::endl;
270  const std::vector<std::pair<uint32_t, float> >& hits_and_fractions =
271  simCluster.hits_and_fractions();
272 
273  // loop over hits
274  for (const auto& it_haf : hits_and_fractions) {
275  unsigned int hitlayer = recHitTools_.getLayerWithOffset(it_haf.first);
276  DetId hitid = (it_haf.first);
277  // dump raw RecHits and match
278  if (rawRecHits_) {
279  if ((hitid.det() == DetId::Forward &&
280  (hitid.subdetId() == HGCEE || hitid.subdetId() == HGCHEF ||
281  hitid.subdetId() == HGCHEB)) ||
282  (hitid.det() == DetId::Hcal && hitid.subdetId() == HcalEndcap))
283  fillWithRecHits(hitmap, hitid, hitlayer, it_haf.second, seedDet,
284  seedEnergy);
285  }
286  } // end simHit
287  } // end simCluster
288 
289  auto sumCalibRecHitCalib_fraction = std::accumulate(Energy_layer_calib_fraction_.begin(),
290  Energy_layer_calib_fraction_.end(), 0.);
291  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
292  << ">>> MC Energy: " << it_caloPart.energy()
293  << " reco energy: " << sumCalibRecHitCalib_fraction
294  << std::endl;
295  if (h_EoP_CPene_calib_fraction_.count(seedDet))
296  h_EoP_CPene_calib_fraction_[seedDet]->Fill(sumCalibRecHitCalib_fraction /
297  it_caloPart.energy());
298 
299  // Loop on reconstructed SC.
300  const auto& clusters = *hgcalMultiClustersHandle;
301  float total_energy = 0.;
302  float max_dR2 = 0.0025;
303  auto closest = std::min_element(clusters.begin(),
304  clusters.end(),
305  [&](const reco::PFCluster &a,
306  const reco::PFCluster &b) {
307  auto dR2_a = reco::deltaR2(it_caloPart, a);
308  auto dR2_b = reco::deltaR2(it_caloPart, b);
309  auto ERatio_a = a.correctedEnergy()/it_caloPart.energy();
310  auto ERatio_b = b.correctedEnergy()/it_caloPart.energy();
311  // If both clusters are within 0.0025, mark as first (min) the
312  // element with the highest ratio against the SimCluster
313  if (dR2_a < max_dR2 && dR2_b < max_dR2)
314  return ERatio_a > ERatio_b;
315  return dR2_a < dR2_b;
316  });
317  if (closest != clusters.end()
318  && reco::deltaR2(*closest, it_caloPart) < 0.01) {
319  total_energy = closest->correctedEnergy();
320  seedDet = recHitTools_.getSiThickness(closest->seed());
321  if (hgcal_EoP_CPene_calib_fraction_.count(seedDet)) {
322  hgcal_EoP_CPene_calib_fraction_[seedDet]->Fill(total_energy /
323  it_caloPart.energy());
324  }
325  }
326 
327  auto closest_fcn = [&](auto const &a, auto const&b){
328  auto dR2_a = reco::deltaR2(it_caloPart, a);
329  auto dR2_b = reco::deltaR2(it_caloPart, b);
330  auto ERatio_a = a.energy()/it_caloPart.energy();
331  auto ERatio_b = b.energy()/it_caloPart.energy();
332  // If both clusters are within 0.0025, mark as first (min) the
333  // element with the highest ratio against the SimCluster
334  if (dR2_a < max_dR2 && dR2_b < max_dR2)
335  return ERatio_a > ERatio_b;
336  return dR2_a < dR2_b;
337  };
338  // ELECTRONS in HGCAL
339  if (PFElectronHandle.isValid()) {
340  auto const & ele = (*PFElectronHandle);
341  auto closest = std::min_element(ele.begin(),
342  ele.end(),
343  closest_fcn);
344  if (closest != ele.end()
345  && closest->superCluster()->seed()->seed().det() == DetId::Forward
346  && reco::deltaR2(*closest, it_caloPart) < 0.01) {
347  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
348  if (hgcal_ele_EoP_CPene_calib_fraction_.count(seedDet)) {
349  hgcal_ele_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() /
350  it_caloPart.energy());
351  }
352  }
353  }
354 
355  // PHOTONS in HGCAL
356  if (PFPhotonHandle.isValid()) {
357  auto const & photon = (*PFPhotonHandle);
358  auto closest = std::min_element(photon.begin(),
359  photon.end(),
360  closest_fcn);
361  if (closest != photon.end()
362  && closest->superCluster()->seed()->seed().det() == DetId::Forward
363  && reco::deltaR2(*closest, it_caloPart) < 0.01) {
364  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
365  if (hgcal_photon_EoP_CPene_calib_fraction_.count(seedDet)) {
366  hgcal_photon_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() /
367  it_caloPart.energy());
368  }
369  }
370  }
371  } // end caloparticle
372 }
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
edm::EDGetTokenT< std::vector< reco::Photon > > photons_
std::array< float, layers_ > Energy_layer_calib_fraction_
edm::EDGetTokenT< std::vector< reco::PFCluster > > hgcalMultiClusters_
std::map< int, MonitorElement * > hgcal_EoP_CPene_calib_fraction_
double correctedEnergy() const
Definition: CaloCluster.h:125
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:210
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
std::map< int, MonitorElement * > hgcal_ele_EoP_CPene_calib_fraction_
#define IfLogTrace(cond, cat)
std::map< int, MonitorElement * > hgcal_photon_EoP_CPene_calib_fraction_
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:61
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
Monte Carlo truth information used for tracking validation.
Definition: SimCluster.h:28
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:111
bool isValid() const
Definition: HandleBase.h:74
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:38
double energy() const
cluster energy
Definition: PFCluster.h:82
hgcal::RecHitTools recHitTools_
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:225
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:103
Definition: DetId.h:18
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:179
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electrons_
std::map< int, MonitorElement * > h_EoP_CPene_calib_fraction_
double b
Definition: hdecay.h:120
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
HLT enums.
void fillWithRecHits(std::map< DetId, const HGCRecHit * > &, DetId, unsigned int, float, int &, float &)
size_type size() const
double a
Definition: hdecay.h:121
std::array< float, layers_ > Energy_layer_calib_
Detector det() const
get the detector field from this detid
Definition: DetId.h:36
void HGCalHitCalibration::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &   
)
overrideprivatevirtual

Implements DQMEDAnalyzer.

Definition at line 121 of file HGCalHitCalibration.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::cd(), h_EoP_CPene_calib_fraction_, hgcal_ele_EoP_CPene_calib_fraction_, hgcal_EoP_CPene_calib_fraction_, hgcal_photon_EoP_CPene_calib_fraction_, LayerOccupancy_, layers_, and DQMStore::IBooker::setCurrentFolder().

123  {
124  ibooker.cd();
125  ibooker.setCurrentFolder("HGCalHitCalibration");
127  ibooker.book1D("h_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
129  ibooker.book1D("h_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
131  ibooker.book1D("h_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
133  ibooker.book1D("hgcal_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
135  ibooker.book1D("hgcal_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
137  ibooker.book1D("hgcal_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
139  ibooker.book1D("hgcal_ele_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
141  ibooker.book1D("hgcal_ele_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
143  ibooker.book1D("hgcal_ele_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
145  ibooker.book1D("hgcal_photon_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
147  ibooker.book1D("hgcal_photon_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
149  ibooker.book1D("hgcal_photon_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
150  LayerOccupancy_ = ibooker.book1D("LayerOccupancy", "", layers_, 0., (float)layers_);
151 }
std::map< int, MonitorElement * > hgcal_EoP_CPene_calib_fraction_
std::map< int, MonitorElement * > hgcal_ele_EoP_CPene_calib_fraction_
std::map< int, MonitorElement * > hgcal_photon_EoP_CPene_calib_fraction_
MonitorElement * LayerOccupancy_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
std::map< int, MonitorElement * > h_EoP_CPene_calib_fraction_
void HGCalHitCalibration::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 376 of file HGCalHitCalibration.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), DEFINE_FWK_MODULE, and AlCaHLTBitMon_QueryRunRegistry::string.

377  {
379  desc.add<int>("debug", 0);
380  desc.add<bool>("rawRecHits", true);
381  desc.add<std::string>("detector", "all");
382  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
383  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
384  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
385  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
386  desc.add<edm::InputTag>("hgcalMultiClusters", edm::InputTag("particleFlowClusterHGCalFromMultiCl"));
387  desc.add<edm::InputTag>("electrons", edm::InputTag("ecalDrivenGsfElectronsFromMultiCl"));
388  desc.add<edm::InputTag>("photons", edm::InputTag("photonsFromMultiCl"));
389  descriptions.add("hgcalHitCalibration", desc);
390 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void HGCalHitCalibration::fillWithRecHits ( std::map< DetId, const HGCRecHit * > &  hitmap,
DetId  hitid,
unsigned int  hitlayer,
float  fraction,
int &  seedDet,
float &  seedEnergy 
)
private

Definition at line 153 of file HGCalHitCalibration.cc.

References debug_, TauDecayModes::dec, DetId::det(), Energy_layer_calib_fraction_, MonitorElement::Fill(), DetId::Forward, dedxEstimators_cff::fraction, hgcal::RecHitTools::getLayerWithOffset(), hgcal::RecHitTools::getSiThickness(), IfLogTrace, LayerOccupancy_, DetId::rawId(), recHitTools_, and DetId::subdetId().

Referenced by analyze().

156  {
157  if (hitmap.find(hitid) == hitmap.end()) {
158  // Hit was not reconstructed
159  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
160  << ">>> Failed to find detid " << std::hex << hitid.rawId() << std::dec
161  << " Det " << hitid.det()
162  << " Subdet " << hitid.subdetId() << std::endl;
163  return;
164  }
165  if (hitid.det() != DetId::Forward) {
166  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
167  << ">>> Wrong type of detid " << std::hex << hitid.rawId() << std::dec
168  << " Det " << hitid.det()
169  << " Subdet " << hitid.subdetId() << std::endl;
170  return;
171  }
172  unsigned int layer = recHitTools_.getLayerWithOffset(hitid);
173  assert(hitlayer == layer);
174  Energy_layer_calib_fraction_[layer] += hitmap[hitid]->energy() * fraction;
175  LayerOccupancy_->Fill(layer);
176  if (seedEnergy < hitmap[hitid]->energy()) {
177  seedEnergy = hitmap[hitid]->energy();
178  seedDet = recHitTools_.getSiThickness(hitid);
179  }
180 }
std::array< float, layers_ > Energy_layer_calib_fraction_
#define IfLogTrace(cond, cat)
void Fill(long long x)
uint32_t rawId() const
get the raw id
Definition: DetId.h:44
MonitorElement * LayerOccupancy_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:38
hgcal::RecHitTools recHitTools_
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:103
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:179
Detector det() const
get the detector field from this detid
Definition: DetId.h:36

Member Data Documentation

int HGCalHitCalibration::algo_
private

Definition at line 66 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

edm::EDGetTokenT<std::vector<CaloParticle> > HGCalHitCalibration::caloParticles_
private

Definition at line 61 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

int HGCalHitCalibration::debug_
private

Definition at line 68 of file HGCalHitCalibration.cc.

Referenced by analyze(), and fillWithRecHits().

edm::EDGetTokenT<std::vector<reco::GsfElectron> > HGCalHitCalibration::electrons_
private

Definition at line 63 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

std::array<float, layers_> HGCalHitCalibration::Energy_layer_calib_
private

Definition at line 78 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

std::array<float, layers_> HGCalHitCalibration::Energy_layer_calib_fraction_
private

Definition at line 79 of file HGCalHitCalibration.cc.

Referenced by analyze(), fillWithRecHits(), and HGCalHitCalibration().

std::map<int, MonitorElement*> HGCalHitCalibration::h_EoP_CPene_calib_fraction_
private

Definition at line 71 of file HGCalHitCalibration.cc.

Referenced by analyze(), and bookHistograms().

std::map<int, MonitorElement*> HGCalHitCalibration::hgcal_ele_EoP_CPene_calib_fraction_
private

Definition at line 73 of file HGCalHitCalibration.cc.

Referenced by analyze(), and bookHistograms().

std::map<int, MonitorElement*> HGCalHitCalibration::hgcal_EoP_CPene_calib_fraction_
private

Definition at line 72 of file HGCalHitCalibration.cc.

Referenced by analyze(), and bookHistograms().

std::map<int, MonitorElement*> HGCalHitCalibration::hgcal_photon_EoP_CPene_calib_fraction_
private

Definition at line 74 of file HGCalHitCalibration.cc.

Referenced by analyze(), and bookHistograms().

edm::EDGetTokenT<std::vector<reco::PFCluster> > HGCalHitCalibration::hgcalMultiClusters_
private

Definition at line 62 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

MonitorElement* HGCalHitCalibration::LayerOccupancy_
private

Definition at line 75 of file HGCalHitCalibration.cc.

Referenced by bookHistograms(), and fillWithRecHits().

int HGCalHitCalibration::layers_ = 60
staticprivate

Definition at line 77 of file HGCalHitCalibration.cc.

Referenced by bookHistograms().

edm::EDGetTokenT<std::vector<reco::Photon> > HGCalHitCalibration::photons_
private

Definition at line 64 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

bool HGCalHitCalibration::rawRecHits_
private

Definition at line 67 of file HGCalHitCalibration.cc.

Referenced by analyze().

edm::EDGetTokenT<HGCRecHitCollection> HGCalHitCalibration::recHitsBH_
private

Definition at line 60 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

edm::EDGetTokenT<HGCRecHitCollection> HGCalHitCalibration::recHitsEE_
private

Definition at line 58 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

edm::EDGetTokenT<HGCRecHitCollection> HGCalHitCalibration::recHitsFH_
private

Definition at line 59 of file HGCalHitCalibration.cc.

Referenced by analyze(), and HGCalHitCalibration().

hgcal::RecHitTools HGCalHitCalibration::recHitTools_
private

Definition at line 69 of file HGCalHitCalibration.cc.

Referenced by analyze(), and fillWithRecHits().