CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
CaloParticleValidation Class Reference
Inheritance diagram for CaloParticleValidation:
DQMGlobalEDAnalyzer< Histograms_CaloParticleValidation > edm::global::EDAnalyzer< edm::RunCache< Histograms_CaloParticleValidation >, Args... > edm::global::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 CaloParticleValidation (const edm::ParameterSet &)
 
 ~CaloParticleValidation () override
 
- Public Member Functions inherited from edm::global::EDAnalyzer< edm::RunCache< Histograms_CaloParticleValidation >, Args... >
 EDAnalyzer ()=default
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () 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
 
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::global::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

void bookHistograms (DQMStore::ConcurrentBooker &, edm::Run const &, edm::EventSetup const &, Histograms_CaloParticleValidation &) const override
 
void dqmAnalyze (edm::Event const &, edm::EventSetup const &, Histograms_CaloParticleValidation const &) const override
 

Private Attributes

edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
 
std::string folder_
 
std::vector< int > particles_to_monitor_
 
edm::EDGetTokenT< HGCRecHitCollectionrecHitsBH_
 
edm::EDGetTokenT< HGCRecHitCollectionrecHitsEE_
 
edm::EDGetTokenT< HGCRecHitCollectionrecHitsFH_
 
edm::EDGetTokenT< reco::PFCandidateCollectionsimPFCandidates_
 
edm::EDGetTokenT< std::vector< reco::SuperCluster > > simPFClusters_
 
edm::EDGetTokenT< std::vector< SimVertex > > simVertices_
 

Additional Inherited Members

- Public Types inherited from edm::global::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
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)
 

Detailed Description

Definition at line 56 of file CaloParticleValidation.cc.

Constructor & Destructor Documentation

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

Definition at line 98 of file CaloParticleValidation.cc.

99  : folder_(iConfig.getParameter<std::string>("folder")),
100  particles_to_monitor_(iConfig.getParameter<std::vector<int> >("particles_to_monitor")),
101  simVertices_(consumes<std::vector<SimVertex>>(iConfig.getParameter<edm::InputTag>("simVertices"))),
102  caloParticles_(consumes<std::vector<CaloParticle> >(iConfig.getParameter<edm::InputTag>("caloParticles"))),
103  simPFClusters_(consumes<std::vector<reco::SuperCluster>>(iConfig.getParameter<edm::InputTag>("simPFClusters"))),
104  simPFCandidates_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("simPFCandidates"))),
105  recHitsEE_(consumes<HGCRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsEE"))),
106  recHitsFH_(consumes<HGCRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsFH"))),
107  recHitsBH_(consumes<HGCRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsBH")))
108 {
109  //now do what ever initialization is needed
110 }
T getParameter(std::string const &) const
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
edm::EDGetTokenT< std::vector< reco::SuperCluster > > simPFClusters_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< std::vector< SimVertex > > simVertices_
edm::EDGetTokenT< reco::PFCandidateCollection > simPFCandidates_
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
std::vector< int > particles_to_monitor_
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
CaloParticleValidation::~CaloParticleValidation ( )
override

Definition at line 113 of file CaloParticleValidation.cc.

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

Member Function Documentation

void CaloParticleValidation::bookHistograms ( DQMStore::ConcurrentBooker ibook,
edm::Run const &  run,
edm::EventSetup const &  iSetup,
Histograms_CaloParticleValidation histos 
) const
overrideprivatevirtual

Implements DQMGlobalEDAnalyzer< Histograms_CaloParticleValidation >.

Definition at line 230 of file CaloParticleValidation.cc.

References DQMStore::ConcurrentBooker::book1D(), DQMStore::ConcurrentBooker::book2D(), reco::PFCandidate::egamma_HF, folder_, reco::PFCandidate::h, trackerHits::histo, PFRecoTauDiscriminationByIsolation_cfi::offset, particles_to_monitor_, and DQMStore::IBooker::setCurrentFolder().

234 {
235  for (auto const particle : particles_to_monitor_) {
236  ibook.setCurrentFolder(folder_ + "CaloParticles/" + std::to_string(particle));
237  auto & histo = histos[particle];
238  histo.eta_ = ibook.book1D("Eta", "Eta", 80, -4., 4.);
239  histo.energy_ = ibook.book1D("Energy", "Energy", 250, 0., 500.);
240  histo.pt_ = ibook.book1D("Pt", "Pt", 100, 0., 100.);
241  histo.nSimClusters_ = ibook.book1D("NSimClusters", "NSimClusters", 100, 0., 100.);
242  histo.nHitInSimClusters_ = ibook.book1D("NHitInSimClusters", "NHitInSimClusters", 100, 0., 100.);
243  histo.selfEnergy_ = ibook.book1D("SelfEnergy", "SelfEnergy", 250, 0., 500.);
244  histo.energyDifference_ = ibook.book1D("EnergyDifference", "(Energy-SelfEnergy)/Energy", 300, -5., 1.);
245  histo.eta_Zorigin_map_ = ibook.book2D("Eta vs Zorigin", "Eta vs Zorigin", 80, -4., 4., 1100, -550., 550.);
246  }
247  int offset = 100000;
248  ibook.setCurrentFolder(folder_ + "PFCandidates");
249  histos[offset].pfcandidateType_ = ibook.book1D("PFCandidateType", "PFCandidateType", 10, 0, 10);
250  histos[offset].pfcandidate_vect_sum_pt_ = ibook.book1D("PFCandidatePtVectSum", "PFCandidatePtVectSum", 200, 0., 200.);
252  ibook.setCurrentFolder(folder_ + "PFCandidates/" + std::to_string(type));
253  auto & histo = histos[offset + type];
254  histo.pfcandidateEnergy_ = ibook.book1D("PFCandidateEnergy", "PFCandidateEnergy", 250, 0., 250.);
255  histo.pfcandidatePt_ = ibook.book1D("PFCandidatePt", "PFCandidatePt", 250, 0., 250.);
256  histo.pfcandidateEta_ = ibook.book1D("PFCandidateEta", "PFCandidateEta", 100, -5., 5.);
257  histo.pfcandidatePhi_ = ibook.book1D("PFCandidatePhi", "PFCandidatePhi", 100, -4., 4.);
258  histo.pfcandidateElementsInBlocks_ = ibook.book1D("PFCandidateElements", "PFCandidateElements", 20, 0., 20.);
259  }
260  // Folder '0' is meant to be cumulative, with no connection to pdgId
261  ibook.setCurrentFolder(folder_ + std::to_string(0));
262  histos[0].simPFSuperClusterSize_ = ibook.book1D("SimPFSuperClusterSize", "SimPFSuperClusterSize", 40, 0., 40.);
263  histos[0].simPFSuperClusterEnergy_ = ibook.book1D("SimPFSuperClusterEnergy", "SimPFSuperClusterEnergy", 250, 0., 500.);
264 }
type
Definition: HCALResponse.h:21
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
ConcurrentMonitorElement book2D(Args &&...args)
Definition: DQMStore.h:163
ConcurrentMonitorElement book1D(Args &&...args)
Definition: DQMStore.h:160
std::vector< int > particles_to_monitor_
void CaloParticleValidation::dqmAnalyze ( edm::Event const &  iEvent,
edm::EventSetup const &  iSetup,
Histograms_CaloParticleValidation const &  histos 
) const
overrideprivate

Definition at line 127 of file CaloParticleValidation.cc.

References caloTruthProducer_cfi::caloParticles, caloParticles_, edm::Event::getByToken(), trackerHits::histo, SimCluster::hits_and_fractions(), mps_fire::i, LogDebug, PFRecoTauDiscriminationByIsolation_cfi::offset, or, position, recHitsBH_, recHitsEE_, recHitsFH_, SimDataFormats::CaloAnalysis::sc, trackerHits::simHits, simPFCandidates_, simPFClusters_, simVertices_, and mathSSE::sqrt().

129 {
130  using namespace edm;
131 
132  Handle<HGCRecHitCollection> recHitHandleEE;
133  Handle<HGCRecHitCollection> recHitHandleFH;
134  Handle<HGCRecHitCollection> recHitHandleBH;
135  // make a map detid-rechit
136 
137  iEvent.getByToken(recHitsEE_, recHitHandleEE);
138  iEvent.getByToken(recHitsFH_, recHitHandleFH);
139  iEvent.getByToken(recHitsBH_, recHitHandleBH);
140  const auto& rechitsEE = *recHitHandleEE;
141  const auto& rechitsFH = *recHitHandleFH;
142  const auto& rechitsBH = *recHitHandleBH;
143  std::map<DetId, const HGCRecHit*> hitmap;
144  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
145  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
146  }
147  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
148  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
149  }
150  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
151  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
152  }
153 
154  Handle<std::vector<SimVertex>> simVerticesHandle;
155  iEvent.getByToken(simVertices_, simVerticesHandle);
156  std::vector<SimVertex> const & simVertices = *simVerticesHandle;
157 
158  Handle<std::vector<CaloParticle> > caloParticleHandle;
159  iEvent.getByToken(caloParticles_, caloParticleHandle);
160  std::vector<CaloParticle> const & caloParticles = *caloParticleHandle;
161 
162  Handle<std::vector<reco::SuperCluster>> simPFClustersHandle;
163  iEvent.getByToken(simPFClusters_, simPFClustersHandle);
164  std::vector<reco::SuperCluster> const & simPFClusters = *simPFClustersHandle;
165 
166  Handle<reco::PFCandidateCollection> simPFCandidatesHandle;
167  iEvent.getByToken(simPFCandidates_, simPFCandidatesHandle);
168  reco::PFCandidateCollection const & simPFCandidates = *simPFCandidatesHandle;
169 
170  for (auto const caloParticle : caloParticles) {
171  if (caloParticle.g4Tracks()[0].eventId().event() != 0 or caloParticle.g4Tracks()[0].eventId().bunchCrossing() != 0) {
172  LogDebug("CaloParticleValidation") << "Excluding CaloParticles from event: "
173  << caloParticle.g4Tracks()[0].eventId().event()
174  << " with BX: " << caloParticle.g4Tracks()[0].eventId().bunchCrossing() << std::endl;
175  continue;
176  }
177  int id = caloParticle.pdgId();
178  if (histos.count(id)) {
179  auto & histo = histos.at(id);
180  histo.eta_.fill(caloParticle.eta());
181  histo.pt_.fill(caloParticle.pt());
182  histo.energy_.fill(caloParticle.energy());
183  histo.nSimClusters_.fill(caloParticle.simClusters().size());
184  // Find the corresponding vertex.
185  histo.eta_Zorigin_map_.fill(
186  simVertices.at(caloParticle.g4Tracks()[0].vertIndex()).position().z(), caloParticle.eta());
187  int simHits = 0;
188  float energy = 0.;
189  for (auto const sc : caloParticle.simClusters()) {
190  simHits += sc->hits_and_fractions().size();
191  for (auto const h_and_f : sc->hits_and_fractions()) {
192  if (hitmap.count(h_and_f.first))
193  energy += hitmap[h_and_f.first]->energy() * h_and_f.second;
194  }
195  }
196  histo.nHitInSimClusters_.fill((float)simHits);
197  histo.selfEnergy_.fill(energy);
198  histo.energyDifference_.fill(1.- energy/caloParticle.energy());
199  }
200  }
201 
202  // simPFSuperClusters
203  for (auto const sc : simPFClusters) {
204  histos.at(0).simPFSuperClusterSize_.fill((float)sc.clustersSize());
205  histos.at(0).simPFSuperClusterEnergy_.fill(sc.rawEnergy());
206  }
207 
208  // simPFCandidates
209  int offset = 100000;
210  double ptx_tot = 0.;
211  double pty_tot = 0.;
212  for (auto const pfc : simPFCandidates) {
213  size_t type = offset + pfc.particleId();
214  ptx_tot += pfc.px();
215  pty_tot += pfc.py();
216  histos.at(offset).pfcandidateType_.fill(type - offset);
217  auto & histo = histos.at(type);
218  histo.pfcandidateEnergy_.fill(pfc.energy());
219  histo.pfcandidatePt_.fill(pfc.pt());
220  histo.pfcandidateEta_.fill(pfc.eta());
221  histo.pfcandidatePhi_.fill(pfc.phi());
222  histo.pfcandidateElementsInBlocks_.fill(pfc.elementsInBlocks().size());
223  }
224  auto & histo = histos.at(offset);
225  histo.pfcandidate_vect_sum_pt_.fill(std::sqrt(ptx_tot*ptx_tot+pty_tot*pty_tot));
226 }
#define LogDebug(id)
type
Definition: HCALResponse.h:21
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
edm::EDGetTokenT< std::vector< reco::SuperCluster > > simPFClusters_
edm::EDGetTokenT< std::vector< SimVertex > > simVertices_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::PFCandidateCollection > simPFCandidates_
T sqrt(T t)
Definition: SSEVec.h:18
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
HLT enums.
static int position[264][3]
Definition: ReadPGInfo.cc:509
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
void CaloParticleValidation::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 268 of file CaloParticleValidation.cc.

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

268  {
269  //The following says we do not know what parameters are allowed so do no validation
270  // Please change this to state exactly what you do use, even if it is no parameters
272  desc.add<std::string>("folder", "HGCAL/"); // Please keep the trailing '/'
273  desc.add<std::vector<int> > ("particles_to_monitor", {11, -11, 13, -13, 22, 111, 211, -211, 321, -321});
274  desc.add<edm::InputTag>("simVertices", edm::InputTag("g4SimHits"));
275  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
276  desc.add<edm::InputTag>("simPFClusters", edm::InputTag("simPFProducer", "perfect"));
277  desc.add<edm::InputTag>("simPFCandidates", edm::InputTag("simPFProducer"));
278  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit","HGCEERecHits"));
279  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit","HGCHEFRecHits"));
280  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit","HGCHEBRecHits"));
281  descriptions.add("caloparticlevalidationDefault", desc);
282 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

Member Data Documentation

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

Definition at line 79 of file CaloParticleValidation.cc.

Referenced by dqmAnalyze().

std::string CaloParticleValidation::folder_
private

Definition at line 75 of file CaloParticleValidation.cc.

Referenced by bookHistograms().

std::vector<int> CaloParticleValidation::particles_to_monitor_
private

Definition at line 76 of file CaloParticleValidation.cc.

Referenced by bookHistograms().

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

Definition at line 84 of file CaloParticleValidation.cc.

Referenced by dqmAnalyze().

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

Definition at line 82 of file CaloParticleValidation.cc.

Referenced by dqmAnalyze().

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

Definition at line 83 of file CaloParticleValidation.cc.

Referenced by dqmAnalyze().

edm::EDGetTokenT<reco::PFCandidateCollection> CaloParticleValidation::simPFCandidates_
private

Definition at line 81 of file CaloParticleValidation.cc.

Referenced by dqmAnalyze().

edm::EDGetTokenT<std::vector<reco::SuperCluster> > CaloParticleValidation::simPFClusters_
private

Definition at line 80 of file CaloParticleValidation.cc.

Referenced by dqmAnalyze().

edm::EDGetTokenT<std::vector<SimVertex> > CaloParticleValidation::simVertices_
private

Definition at line 78 of file CaloParticleValidation.cc.

Referenced by dqmAnalyze().