CMS 3D CMS Logo

CaloParticleDebugger.cc
Go to the documentation of this file.
1 //
2 // Original Author: Marco Rovere
3 // Created: Fri, 10 Nov 2017 14:39:18 GMT
4 //
5 //
6 //
7 // system include files
8 #include <memory>
9 #include <iostream>
10 #include <numeric>
11 #include <algorithm>
12 
13 // user include files
17 
22 
24 
33 
42 
43 //
44 // class declaration
45 //
46 
48 public:
49  explicit CaloParticleDebugger(const edm::ParameterSet&);
50  ~CaloParticleDebugger() override;
51 
52  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
53 
54 private:
55  void beginJob() override;
56  void analyze(const edm::Event&, const edm::EventSetup&) override;
57  void endJob() override;
58  void fillSimHits(std::map<int, float>&, const edm::Event&, const edm::EventSetup&);
65  std::vector<edm::InputTag> collectionTags_;
72  std::vector<edm::EDGetTokenT<std::vector<PCaloHit>>> collectionTagsToken_;
73  int geometryType_ = 0;
74  // ----------member data ---------------------------
75 };
76 
77 //
78 // constants, enums and typedefs
79 //
80 
81 //
82 // static data member definitions
83 //
84 
85 //
86 // constructors and destructor
87 //
89  : simTracks_(iConfig.getParameter<edm::InputTag>("simTracks")),
90  genParticles_(iConfig.getParameter<edm::InputTag>("genParticles")),
91  simVertices_(iConfig.getParameter<edm::InputTag>("simVertices")),
92  trackingParticles_(iConfig.getParameter<edm::InputTag>("trackingParticles")),
93  caloParticles_(iConfig.getParameter<edm::InputTag>("caloParticles")),
94  simClusters_(iConfig.getParameter<edm::InputTag>("simClusters")),
95  collectionTags_(iConfig.getParameter<std::vector<edm::InputTag>>("collectionTags")) {
97  simTracksToken_ = iC.consumes<std::vector<SimTrack>>(simTracks_);
98  genParticlesToken_ = iC.consumes<std::vector<reco::GenParticle>>(genParticles_);
99  simVerticesToken_ = iC.consumes<std::vector<SimVertex>>(simVertices_);
100  trackingParticlesToken_ = iC.consumes<std::vector<TrackingParticle>>(trackingParticles_);
101  caloParticlesToken_ = iC.consumes<std::vector<CaloParticle>>(caloParticles_);
102  simClustersToken_ = iC.consumes<std::vector<SimCluster>>(simClusters_);
103  for (auto const& collectionTag : collectionTags_) {
104  collectionTagsToken_.push_back(iC.consumes<std::vector<PCaloHit>>(collectionTag));
105  }
106 }
107 
109 
110 //
111 // member functions
112 //
113 
114 // ------------ method called for each event ------------
116  using namespace edm;
117  using std::begin;
118  using std::end;
119  using std::iota;
120  using std::sort;
121 
125  edm::Handle<std::vector<TrackingParticle>> trackingParticlesH;
128 
129  iEvent.getByToken(simTracksToken_, simTracksH);
130  auto const& tracks = *simTracksH.product();
131  std::vector<int> sorted_tracks_idx(tracks.size());
132  iota(begin(sorted_tracks_idx), end(sorted_tracks_idx), 0);
133  sort(begin(sorted_tracks_idx), end(sorted_tracks_idx), [&tracks](int i, int j) {
134  return tracks[i].momentum().eta() < tracks[j].momentum().eta();
135  });
136 
137  iEvent.getByToken(genParticlesToken_, genParticlesH);
138  auto const& genParticles = *genParticlesH.product();
139  std::vector<int> sorted_genParticles_idx(genParticles.size());
140  iota(begin(sorted_genParticles_idx), end(sorted_genParticles_idx), 0);
141  sort(begin(sorted_genParticles_idx), end(sorted_genParticles_idx), [&genParticles](int i, int j) {
142  return genParticles[i].momentum().eta() < genParticles[j].momentum().eta();
143  });
144 
145  iEvent.getByToken(simVerticesToken_, simVerticesH);
146  auto const& vertices = *simVerticesH.product();
147  std::vector<int> sorted_vertices_idx(vertices.size());
148  iota(begin(sorted_vertices_idx), end(sorted_vertices_idx), 0);
149  sort(begin(sorted_vertices_idx), end(sorted_vertices_idx), [&vertices](int i, int j) {
150  return vertices[i].vertexId() < vertices[j].vertexId();
151  });
152 
153  iEvent.getByToken(trackingParticlesToken_, trackingParticlesH);
154  auto const& trackingpart = *trackingParticlesH.product();
155  std::vector<int> sorted_tp_idx(trackingpart.size());
156  iota(begin(sorted_tp_idx), end(sorted_tp_idx), 0);
157  sort(begin(sorted_tp_idx), end(sorted_tp_idx), [&trackingpart](int i, int j) {
158  return trackingpart[i].eta() < trackingpart[j].eta();
159  });
160 
161  iEvent.getByToken(caloParticlesToken_, caloParticlesH);
162  auto const& calopart = *caloParticlesH.product();
163  std::vector<int> sorted_cp_idx(calopart.size());
164  iota(begin(sorted_cp_idx), end(sorted_cp_idx), 0);
165  sort(begin(sorted_cp_idx), end(sorted_cp_idx), [&calopart](int i, int j) {
166  return calopart[i].eta() < calopart[j].eta();
167  });
168 
169  iEvent.getByToken(simClustersToken_, simClustersH);
170  auto const& simclusters = *simClustersH.product();
171  std::vector<int> sorted_simcl_idx(simclusters.size());
172  iota(begin(sorted_simcl_idx), end(sorted_simcl_idx), 0);
173  sort(begin(sorted_simcl_idx), end(sorted_simcl_idx), [&simclusters](int i, int j) {
174  return simclusters[i].eta() < simclusters[j].eta();
175  });
176 
177  // Let's first fill in hits information
178  std::map<int, float> detIdToTotalSimEnergy;
179  fillSimHits(detIdToTotalSimEnergy, iEvent, iSetup);
180 
181  int idx = 0;
182 
183  std::map<int, int> trackid_to_track_index;
184  std::cout << "Printing SimTracks information" << std::endl;
185  std::cout << "IDX\tTrackId\tPDGID\tMOMENTUM(x,y,z,E)\tVertexIdx\tGenPartIdx" << std::endl;
186  for (auto i : sorted_tracks_idx) {
187  auto const& t = tracks[i];
188  std::cout << idx << "\t" << t.trackId() << "\t" << t << std::endl;
189  trackid_to_track_index[t.trackId()] = idx;
190  idx++;
191  }
192 
193  std::cout << "Printing GenParticles information" << std::endl;
194  std::cout << "IDX\tPDGID\tMOMENTUM(x,y,z)\tVertex(x,y,z)" << std::endl;
195  for (auto i : sorted_genParticles_idx) {
196  auto const& gp = genParticles[i];
197  std::cout << i << "\t" << gp.pdgId() << "\t" << gp.momentum() << "\t" << gp.vertex() << std::endl;
198  }
199 
200  std::cout << "Printing SimVertex information" << std::endl;
201  std::cout << "IDX\tPOSITION(x,y,z)\tPARENT_INDEX\tVERTEX_ID" << std::endl;
202  for (auto i : sorted_vertices_idx) {
203  auto const& v = vertices[i];
204  std::cout << i << "\t" << v << std::endl;
205  }
206  std::cout << "Printing TrackingParticles information" << std::endl;
207  for (auto i : sorted_tp_idx) {
208  auto const& tp = trackingpart[i];
209  std::cout << i << "\t" << tp << std::endl;
210  }
211 
212  std::cout << "Printing CaloParticles information" << std::endl;
213  idx = 0;
214  for (auto i : sorted_cp_idx) {
215  auto const& cp = calopart[i];
216  std::cout << "\n\n"
217  << idx++ << " |Eta|: " << std::abs(cp.momentum().eta()) << "\tType: " << cp.pdgId()
218  << "\tEnergy: " << cp.energy() << "\tIdx: " << cp.g4Tracks()[0].trackId()
219  << std::endl; // << cp << std::endl;
220  double total_sim_energy = 0.;
221  double total_cp_energy = 0.;
222  std::cout << "--> Overall simclusters's size: " << cp.simClusters().size() << std::endl;
223  // All the next mess just to print the simClusters ordered
224  auto const& simcs = cp.simClusters();
225  std::vector<int> sorted_sc_idx(simcs.size());
226  iota(begin(sorted_sc_idx), end(sorted_sc_idx), 0);
227  sort(begin(sorted_sc_idx), end(sorted_sc_idx), [&simcs](int i, int j) {
228  return simcs[i]->momentum().eta() < simcs[j]->momentum().eta();
229  });
230  for (auto i : sorted_sc_idx) {
231  std::cout << *(simcs[i]);
232  }
233 
234  for (auto const& sc : cp.simClusters()) {
235  for (auto const& cl : sc->hits_and_fractions()) {
236  total_sim_energy += detIdToTotalSimEnergy[cl.first] * cl.second;
237  total_cp_energy += cp.energy() * cl.second;
238  }
239  }
240  std::cout << "--> Overall SC energy (sum using sim energies): " << total_sim_energy << std::endl;
241  std::cout << "--> Overall SC energy (sum using CaloP energies): " << total_cp_energy << std::endl;
242  }
243 
244  idx = 0;
245  std::cout << "Printing SimClusters information" << std::endl;
246  for (auto i : sorted_simcl_idx) {
247  auto const& simcl = simclusters[i];
248  std::cout << "\n\n"
249  << idx++ << " |Eta|: " << std::abs(simcl.momentum().eta()) << "\tType: " << simcl.pdgId()
250  << "\tEnergy: " << simcl.energy() << "\tKey: " << i << std::endl; // << simcl << std::endl;
251  double total_sim_energy = 0.;
252  std::cout << "--> Overall simclusters's size: " << simcl.numberOfRecHits() << std::endl;
253  for (auto const& cl : simcl.hits_and_fractions()) {
254  total_sim_energy += detIdToTotalSimEnergy[cl.first] * cl.second;
255  }
256  std::cout << simcl << std::endl;
257  std::cout << "--> Overall SimCluster energy (sum using sim energies): " << total_sim_energy << std::endl;
258  }
259 }
260 
261 // ------------ method called once each job just before starting event loop ------------
263 
264 // ------------ method called once each job just after ending the event loop ------------
266 
267 void CaloParticleDebugger::fillSimHits(std::map<int, float>& detIdToTotalSimEnergy,
268  const edm::Event& iEvent,
269  const edm::EventSetup& iSetup) {
270  // Taken needed quantities from the EventSetup
272  iSetup.get<CaloGeometryRecord>().get(geom);
273  const HGCalGeometry *eegeom = nullptr, *fhgeom = nullptr, *bhgeomnew = nullptr;
274  const HcalGeometry* bhgeom = nullptr;
275  const HGCalDDDConstants* hgddd[3];
276  const HGCalTopology* hgtopo[3];
277  const HcalDDDRecConstants* hcddd = nullptr;
278  eegeom =
280  //check if it's the new geometry
281  if (eegeom) {
282  geometryType_ = 1;
283  fhgeom = static_cast<const HGCalGeometry*>(
285  bhgeomnew = static_cast<const HGCalGeometry*>(
287  } else {
288  geometryType_ = 0;
289  eegeom = static_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward, HGCEE));
290  fhgeom = static_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward, HGCHEF));
291  bhgeom = static_cast<const HcalGeometry*>(geom->getSubdetectorGeometry(DetId::Hcal, HcalEndcap));
292  }
293  hgtopo[0] = &(eegeom->topology());
294  hgtopo[1] = &(fhgeom->topology());
295  if (bhgeomnew)
296  hgtopo[2] = &(bhgeomnew->topology());
297 
298  for (unsigned i = 0; i < 3; ++i) {
299  if (hgtopo[i])
300  hgddd[i] = &(hgtopo[i]->dddConstants());
301  }
302 
303  if (bhgeom)
304  hcddd = bhgeom->topology().dddConstants();
305 
306  // loop over the collections
307  int token = 0;
308  for (auto const& collectionTag : collectionTags_) {
310  const bool isHcal = (collectionTag.instance().find("HcalHits") != std::string::npos);
311  iEvent.getByToken(collectionTagsToken_[token++], hSimHits);
312  for (auto const& simHit : *hSimHits) {
313  DetId id(0);
314  const uint32_t simId = simHit.id();
315  if (geometryType_ == 1) {
316  //no test numbering in new geometry
317  id = simId;
318  } else if (isHcal) {
319  assert(hcddd);
320  HcalDetId hid = HcalHitRelabeller::relabel(simId, hcddd);
321  if (hid.subdet() == HcalEndcap)
322  id = hid;
323  } else {
324  int subdet, layer, cell, sec, subsec, zp;
325  HGCalTestNumbering::unpackHexagonIndex(simId, subdet, zp, layer, sec, subsec, cell);
326  const HGCalDDDConstants* ddd = hgddd[subdet - 3];
327  std::pair<int, int> recoLayerCell = ddd->simToReco(cell, layer, sec, hgtopo[subdet - 3]->detectorType());
328  cell = recoLayerCell.first;
329  layer = recoLayerCell.second;
330  // skip simhits with bad barcodes or non-existant layers
331  if (layer == -1 || simHit.geantTrackId() == 0)
332  continue;
333  id = HGCalDetId((ForwardSubdetector)subdet, zp, layer, subsec, sec, cell);
334  }
335 
336  if (DetId(0) == id)
337  continue;
338 
339  detIdToTotalSimEnergy[id.rawId()] += simHit.energy();
340  }
341  } // end of loop over InputTags
342 }
343 
344 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
347  desc.add<edm::InputTag>("simTracks", edm::InputTag("g4SimHits"));
348  desc.add<edm::InputTag>("genParticles", edm::InputTag("genParticles"));
349  desc.add<edm::InputTag>("simVertices", edm::InputTag("g4SimHits"));
350  desc.add<edm::InputTag>("trackingParticles", edm::InputTag("mix", "MergedTrackTruth"));
351  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
352  desc.add<edm::InputTag>("simClusters", edm::InputTag("mix", "MergedCaloTruth"));
353  desc.add<std::vector<edm::InputTag>>("collectionTags",
354  {edm::InputTag("g4SimHits", "HGCHitsEE"),
355  edm::InputTag("g4SimHits", "HGCHitsHEfront"),
356  edm::InputTag("g4SimHits", "HcalHits")});
357  descriptions.add("caloParticleDebugger", desc);
358 }
359 
360 // define this as a plug-in
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:168
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:146
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< std::vector< TrackingParticle > > trackingParticlesToken_
const std::vector< SimTrack > & g4Tracks() const
Definition: CaloParticle.h:74
ForwardSubdetector
edm::EDGetTokenT< std::vector< SimCluster > > simClustersToken_
const HcalTopology & topology() const
Definition: HcalGeometry.h:117
CaloParticleDebugger(const edm::ParameterSet &)
std::pair< int, int > simToReco(int cell, int layer, int mod, bool half) const
int iEvent
Definition: GenABIO.cc:224
std::vector< edm::EDGetTokenT< std::vector< PCaloHit > > > collectionTagsToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const SimClusterRefVector & simClusters() const
Definition: CaloParticle.h:72
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
const HGCalTopology & topology() const
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:98
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< std::vector< SimTrack > > simTracksToken_
#define end
Definition: vmac.h:39
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticlesToken_
Definition: DetId.h:18
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
math::XYZVectorF momentum() const
spatial momentum vector
Definition: CaloParticle.h:88
T const * product() const
Definition: Handle.h:74
const HGCalDDDConstants & dddConstants() const
void fillSimHits(std::map< int, float > &, const edm::Event &, const edm::EventSetup &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int pdgId() const
PDG ID.
Definition: CaloParticle.h:43
edm::EDGetTokenT< std::vector< reco::GenParticle > > genParticlesToken_
#define begin
Definition: vmac.h:32
HLT enums.
T get() const
Definition: EventSetup.h:71
void analyze(const edm::Event &, const edm::EventSetup &) override
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
DetId relabel(const uint32_t testId) const
edm::InputTag trackingParticles_
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:181
std::vector< edm::InputTag > collectionTags_
edm::EDGetTokenT< std::vector< SimVertex > > simVerticesToken_