CMS 3D CMS Logo

CaloParticleValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Class: CaloParticleValidation
4 // Original Author: Marco Rovere
5 // Created: Thu, 18 Jan 2018 15:54:55 GMT
6 //
7 //
8 
9 #include <string>
10 #include <unordered_map>
11 
12 // user include files
15 
18 
20 
28 
29 //
30 // class declaration
31 //
32 
39  ConcurrentMonitorElement selfEnergy_; // this is the sum of the energy associated to all recHits linked to all SimClusters
40  ConcurrentMonitorElement energyDifference_; // This contains (energy-selfEnergy)/energy
50 };
51 
52 
53 using Histograms_CaloParticleValidation = std::unordered_map<int, Histogram_CaloParticleSingle>;
54 
55 class CaloParticleValidation : public DQMGlobalEDAnalyzer<Histograms_CaloParticleValidation> {
56  public:
58  ~CaloParticleValidation() override;
59 
60  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
61 
62 
63  private:
65  edm::Run const&,
66  edm::EventSetup const&,
67  Histograms_CaloParticleValidation&) const override;
68 
69  void dqmAnalyze(edm::Event const&,
70  edm::EventSetup const&,
71  Histograms_CaloParticleValidation const&) const override;
72 
73  // ----------member data ---------------------------
75  std::vector<int> particles_to_monitor_;
76 
84 };
85 
86 //
87 // constants, enums and typedefs
88 //
89 
90 //
91 // static data member definitions
92 //
93 
94 //
95 // constructors and destructor
96 //
98  : folder_(iConfig.getParameter<std::string>("folder")),
99  particles_to_monitor_(iConfig.getParameter<std::vector<int> >("particles_to_monitor")),
100  simVertices_(consumes<std::vector<SimVertex>>(iConfig.getParameter<edm::InputTag>("simVertices"))),
101  caloParticles_(consumes<std::vector<CaloParticle> >(iConfig.getParameter<edm::InputTag>("caloParticles"))),
102  simPFClusters_(consumes<std::vector<reco::SuperCluster>>(iConfig.getParameter<edm::InputTag>("simPFClusters"))),
103  simPFCandidates_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("simPFCandidates"))),
104  recHitsEE_(consumes<HGCRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsEE"))),
105  recHitsFH_(consumes<HGCRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsFH"))),
106  recHitsBH_(consumes<HGCRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsBH")))
107 {
108  //now do what ever initialization is needed
109 }
110 
111 
113 {
114  // do anything here that needs to be done at desctruction time
115  // (e.g. close files, deallocate resources etc.)
116 }
117 
118 
119 //
120 // member functions
121 //
122 
123 // ------------ method called for each event ------------
124 
125 void
128 {
129  using namespace edm;
130 
131  Handle<HGCRecHitCollection> recHitHandleEE;
132  Handle<HGCRecHitCollection> recHitHandleFH;
133  Handle<HGCRecHitCollection> recHitHandleBH;
134  // make a map detid-rechit
135 
136  iEvent.getByToken(recHitsEE_, recHitHandleEE);
137  iEvent.getByToken(recHitsFH_, recHitHandleFH);
138  iEvent.getByToken(recHitsBH_, recHitHandleBH);
139  const auto& rechitsEE = *recHitHandleEE;
140  const auto& rechitsFH = *recHitHandleFH;
141  const auto& rechitsBH = *recHitHandleBH;
142  std::map<DetId, const HGCRecHit*> hitmap;
143  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
144  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
145  }
146  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
147  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
148  }
149  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
150  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
151  }
152 
153  Handle<std::vector<SimVertex>> simVerticesHandle;
154  iEvent.getByToken(simVertices_, simVerticesHandle);
155  std::vector<SimVertex> const & simVertices = *simVerticesHandle;
156 
157  Handle<std::vector<CaloParticle> > caloParticleHandle;
158  iEvent.getByToken(caloParticles_, caloParticleHandle);
159  std::vector<CaloParticle> const & caloParticles = *caloParticleHandle;
160 
161  Handle<std::vector<reco::SuperCluster>> simPFClustersHandle;
162  iEvent.getByToken(simPFClusters_, simPFClustersHandle);
163  std::vector<reco::SuperCluster> const & simPFClusters = *simPFClustersHandle;
164 
165  Handle<reco::PFCandidateCollection> simPFCandidatesHandle;
166  iEvent.getByToken(simPFCandidates_, simPFCandidatesHandle);
167  reco::PFCandidateCollection const & simPFCandidates = *simPFCandidatesHandle;
168 
169  for (auto const caloParticle : caloParticles) {
170  int id = caloParticle.pdgId();
171  if (histos.count(id)) {
172  auto & histo = histos.at(id);
173  histo.eta_.fill(caloParticle.eta());
174  histo.pt_.fill(caloParticle.pt());
175  histo.energy_.fill(caloParticle.energy());
176  histo.nSimClusters_.fill(caloParticle.simClusters().size());
177  // Find the corresponding vertex.
178  histo.eta_Zorigin_map_.fill(
179  simVertices.at(caloParticle.g4Tracks()[0].vertIndex()).position().z(), caloParticle.eta());
180  int simHits = 0;
181  float energy = 0.;
182  for (auto const sc : caloParticle.simClusters()) {
183  simHits += sc->hits_and_fractions().size();
184  for (auto const h_and_f : sc->hits_and_fractions()) {
185  if (hitmap.count(h_and_f.first))
186  energy += hitmap[h_and_f.first]->energy() * h_and_f.second;
187  }
188  }
189  histo.nHitInSimClusters_.fill((float)simHits);
190  histo.selfEnergy_.fill(energy);
191  histo.energyDifference_.fill(1.- energy/caloParticle.energy());
192  }
193  }
194 
195  // simPFSuperClusters
196  for (auto const sc : simPFClusters) {
197  histos.at(0).simPFSuperClusterSize_.fill((float)sc.clustersSize());
198  histos.at(0).simPFSuperClusterEnergy_.fill(sc.rawEnergy());
199  }
200 
201  // simPFCandidates
202  int offset = 100000;
203  for (auto const pfc : simPFCandidates) {
204  size_t type = offset + pfc.particleId();
205  histos.at(offset).pfcandidateType_.fill(type - offset);
206  auto & histo = histos.at(type);
207  histo.pfcandidateEnergy_.fill(pfc.energy());
208  histo.pfcandidatePt_.fill(pfc.pt());
209  histo.pfcandidateEta_.fill(pfc.eta());
210  histo.pfcandidatePhi_.fill(pfc.phi());
211  histo.pfcandidateElementsInBlocks_.fill(pfc.elementsInBlocks().size());
212  }
213 }
214 
215 
216 void
218  edm::Run const & run,
219  edm::EventSetup const & iSetup,
221 {
222  for (auto const particle : particles_to_monitor_) {
223  ibook.setCurrentFolder(folder_ + "CaloParticles/" + std::to_string(particle));
224  auto & histo = histos[particle];
225  histo.eta_ = ibook.book1D("Eta", "Eta", 80, -4., 4.);
226  histo.energy_ = ibook.book1D("Energy", "Energy", 250, 0., 500.);
227  histo.pt_ = ibook.book1D("Pt", "Pt", 100, 0., 100.);
228  histo.nSimClusters_ = ibook.book1D("NSimClusters", "NSimClusters", 100, 0., 100.);
229  histo.nHitInSimClusters_ = ibook.book1D("NHitInSimClusters", "NHitInSimClusters", 100, 0., 100.);
230  histo.selfEnergy_ = ibook.book1D("SelfEnergy", "SelfEnergy", 250, 0., 500.);
231  histo.energyDifference_ = ibook.book1D("EnergyDifference", "(Energy-SelfEnergy)/Energy", 300, -5., 1.);
232  histo.eta_Zorigin_map_ = ibook.book2D("Eta vs Zorigin", "Eta vs Zorigin", 80, -4., 4., 1100, -550., 550.);
233  }
234  int offset = 100000;
235  ibook.setCurrentFolder(folder_ + "PFCandidates");
236  histos[offset].pfcandidateType_ = ibook.book1D("PFCandidateType", "PFCandidateType", 10, 0, 10);
238  ibook.setCurrentFolder(folder_ + "PFCandidates/" + std::to_string(type));
239  auto & histo = histos[offset + type];
240  histo.pfcandidateEnergy_ = ibook.book1D("PFCandidateEnergy", "PFCandidateEnergy", 250, 0., 250.);
241  histo.pfcandidatePt_ = ibook.book1D("PFCandidatePt", "PFCandidatePt", 250, 0., 250.);
242  histo.pfcandidateEta_ = ibook.book1D("PFCandidateEta", "PFCandidateEta", 100, -5., 5.);
243  histo.pfcandidatePhi_ = ibook.book1D("PFCandidatePhi", "PFCandidatePhi", 100, -4., 4.);
244  histo.pfcandidateElementsInBlocks_ = ibook.book1D("PFCandidateElements", "PFCandidateElements", 20, 0., 20.);
245  }
246  // Folder '0' is meant to be cumulative, with no connection to pdgId
247  ibook.setCurrentFolder(folder_ + std::to_string(0));
248  histos[0].simPFSuperClusterSize_ = ibook.book1D("SimPFSuperClusterSize", "SimPFSuperClusterSize", 40, 0., 40.);
249  histos[0].simPFSuperClusterEnergy_ = ibook.book1D("SimPFSuperClusterEnergy", "SimPFSuperClusterEnergy", 250, 0., 500.);
250 }
251 
252 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
253 void
255  //The following says we do not know what parameters are allowed so do no validation
256  // Please change this to state exactly what you do use, even if it is no parameters
258  desc.add<std::string>("folder", "HGCAL/"); // Please keep the trailing '/'
259  desc.add<std::vector<int> > ("particles_to_monitor", {11, -11, 13, 22, 111, 211, -211, 321, -321});
260  desc.add<edm::InputTag>("simVertices", edm::InputTag("g4SimHits"));
261  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
262  desc.add<edm::InputTag>("simPFClusters", edm::InputTag("simPFProducer", "perfect"));
263  desc.add<edm::InputTag>("simPFCandidates", edm::InputTag("simPFProducer"));
264  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit","HGCEERecHits"));
265  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit","HGCHEFRecHits"));
266  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit","HGCHEBRecHits"));
267  descriptions.add("caloparticlevalidationDefault", desc);
268 }
269 
270 //define this as a plug-in
type
Definition: HCALResponse.h:21
std::unordered_map< int, Histogram_CaloParticleSingle > Histograms_CaloParticleValidation
ConcurrentMonitorElement simPFSuperClusterEnergy_
ConcurrentMonitorElement eta_Zorigin_map_
ConcurrentMonitorElement energyDifference_
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
ConcurrentMonitorElement pfcandidatePhi_
ConcurrentMonitorElement pfcandidateElementsInBlocks_
CaloParticleValidation(const edm::ParameterSet &)
ConcurrentMonitorElement nHitInSimClusters_
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:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual example_global void dqmAnalyze(edm::Event const &,@example_global edm::EventSetup const &,@example_global Histograms___class__ const &) const override
edm::EDGetTokenT< std::vector< reco::SuperCluster > > simPFClusters_
edm::EDGetTokenT< std::vector< SimVertex > > simVertices_
int iEvent
Definition: GenABIO.cc:230
ConcurrentMonitorElement pfcandidateEta_
edm::EDGetTokenT< reco::PFCandidateCollection > simPFCandidates_
ConcurrentMonitorElement book2D(Args &&...args)
Definition: DQMStore.h:244
ConcurrentMonitorElement pfcandidateType_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual example_stream void bookHistograms(DQMStore::IBooker &,@example_stream edm::Run const &,@example_stream edm::EventSetup const &) override
ConcurrentMonitorElement book1D(Args &&...args)
Definition: DQMStore.h:223
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, Histograms_CaloParticleValidation const &) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
std::string folder_
ConcurrentMonitorElement pfcandidateEnergy_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
ConcurrentMonitorElement simPFSuperClusterSize_
static int position[264][3]
Definition: ReadPGInfo.cc:509
ConcurrentMonitorElement nSimClusters_
ConcurrentMonitorElement pfcandidatePt_
void bookHistograms(DQMStore::ConcurrentBooker &, edm::Run const &, edm::EventSetup const &, Histograms_CaloParticleValidation &) const override
ConcurrentMonitorElement energy_
std::vector< int > particles_to_monitor_
ConcurrentMonitorElement selfEnergy_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
Definition: Run.h:44