CMS 3D CMS Logo

HGCalShowerSeparation.cc
Go to the documentation of this file.
1 // user include files
2 
6 
10 
12 
14 
17 
24 
29 
35 
37 
38 #include <map>
39 #include <array>
40 #include <string>
41 #include <numeric>
42 
44  public:
46  ~HGCalShowerSeparation() override;
47 
48  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
49 
50  private:
52  edm::EventSetup const&) override;
53  void analyze(const edm::Event&, const edm::EventSetup&) override;
54 
55  void fillWithRecHits(std::map<DetId, const HGCRecHit*>&, DetId, unsigned int,
56  float, int&, float&);
57 
62 
63  int debug_;
66 
78  std::vector<MonitorElement*> profileOnLayer_;
79  std::vector<MonitorElement*> globalProfileOnLayer_;
80  std::vector<MonitorElement*> distanceOnLayer_;
81  std::vector<MonitorElement*> idealDistanceOnLayer_;
82  std::vector<MonitorElement*> idealDeltaXY_;
83  std::vector<MonitorElement*> centers_;
84 
85  static constexpr int layers_ = 52;
86 };
87 
89  : debug_(iConfig.getParameter<int>("debug")),
90  filterOnEnergyAndCaloP_(iConfig.getParameter<bool>("filterOnEnergyAndCaloP")){
91  auto recHitsEE = iConfig.getParameter<edm::InputTag>("recHitsEE");
92  auto recHitsFH = iConfig.getParameter<edm::InputTag>("recHitsFH");
93  auto recHitsBH = iConfig.getParameter<edm::InputTag>("recHitsBH");
94  auto caloParticles = iConfig.getParameter<edm::InputTag>("caloParticles");
95  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
96  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
97  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
98  caloParticles_ = consumes<std::vector<CaloParticle> >(caloParticles);
99 }
100 
102  // do anything here that needs to be done at desctruction time
103  // (e.g. close files, deallocate resources etc.)
104 }
105 
107  edm::Run const& iRun,
108  edm::EventSetup const& /* iSetup */) {
109  ibooker.cd();
110  ibooker.setCurrentFolder("HGCalShowerSeparation");
111  scEnergy_ = ibooker.book1D("SCEnergy", "SCEnergy", 240, 0., 120.);
112  eta1_ = ibooker.book1D("eta1", "eta1", 80, 0., 4.);
113  eta2_ = ibooker.book1D("eta2", "eta2", 80, 0., 4.);
114  energy1_ = ibooker.book1D("energy1", "energy1", 240, 0., 120.);
115  energy2_ = ibooker.book1D("energy2", "energy2", 240, 0., 120.);
116  energytot_ = ibooker.book1D("energytot", "energytot", 200, 100., 200.);
117  showerProfile_ = ibooker.book2D("ShowerProfile", "ShowerProfile",
118  800, -400., 400.,
119  layers_, 0., (float)layers_);
120  layerEnergy_ = ibooker.book2D("LayerEnergy", "LayerEnergy",
121  60, 0., 60.,
122  50, 0., 0.1);
123  layerDistance_ = ibooker.book2D("LayerDistance", "LayerDistance",
124  60, 0., 60.,
125  400, -400., 400.);
126  etaPhi_ = ibooker.book2D("EtaPhi", "EtaPhi",
127  800, -4., 4.,
128  800, -4., 4.);
129  deltaEtaPhi_ = ibooker.book2D("DeltaEtaPhi", "DeltaEtaPhi",
130  100, -0.5, 0.5,
131  100, -0.5, 0.5);
132  for (int i = 0; i < layers_; ++i) {
133  profileOnLayer_.push_back(ibooker.book2D(std::string("ProfileOnLayer_") + std::to_string(i),
134  std::string("ProfileOnLayer_") + std::to_string(i),
135  120, -600., 600.,
136  120, -600., 600.)
137  );
138  globalProfileOnLayer_.push_back(ibooker.book2D(std::string("GlobalProfileOnLayer_") + std::to_string(i),
139  std::string("GlobalProfileOnLayer_") + std::to_string(i),
140  320, -160., 160.,
141  320, -160., 160.)
142  );
143  distanceOnLayer_.push_back(ibooker.book1D(std::string("DistanceOnLayer_") + std::to_string(i),
144  std::string("DistanceOnLayer_") + std::to_string(i),
145  120, -600., 600.)
146  );
147  idealDistanceOnLayer_.push_back(ibooker.book1D(std::string("IdealDistanceOnLayer_") + std::to_string(i),
148  std::string("IdealDistanceOnLayer_") + std::to_string(i),
149  120, -600., 600.)
150  );
151  idealDeltaXY_.push_back(ibooker.book2D(std::string("IdealDeltaXY_") + std::to_string(i),
152  std::string("IdealDeltaXY_") + std::to_string(i),
153  800, -400., 400.,
154  800, -400., 400.)
155  );
156  centers_.push_back(ibooker.book2D(std::string("Centers_") + std::to_string(i),
157  std::string("Centers_") + std::to_string(i),
158  320, -1600., 1600.,
159  320, -1600., 1600.)
160  );
161  }
162 }
163 
165  const edm::EventSetup& iSetup) {
166  using namespace edm;
167 
168  recHitTools_.getEventSetup(iSetup);
169 
170  Handle<HGCRecHitCollection> recHitHandleEE;
171  Handle<HGCRecHitCollection> recHitHandleFH;
172  Handle<HGCRecHitCollection> recHitHandleBH;
173 
174  Handle<std::vector<CaloParticle> > caloParticleHandle;
175  iEvent.getByToken(caloParticles_, caloParticleHandle);
176  const std::vector<CaloParticle>& caloParticles = *caloParticleHandle;
177 
178 
179  iEvent.getByToken(recHitsEE_, recHitHandleEE);
180  iEvent.getByToken(recHitsFH_, recHitHandleFH);
181  iEvent.getByToken(recHitsBH_, recHitHandleBH);
182  const auto& rechitsEE = *recHitHandleEE;
183  const auto& rechitsFH = *recHitHandleFH;
184  const auto& rechitsBH = *recHitHandleBH;
185 
186  std::map<DetId, const HGCRecHit*> hitmap;
187  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
188  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
189  }
190  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
191  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
192  }
193  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
194  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
195  }
196 
197  // loop over caloParticles
198  IfLogTrace(debug_ > 0, "HGCalShowerSeparation")
199  << "Number of caloParticles: " << caloParticles.size() << std::endl;
200  if (caloParticles.size() == 2) {
201  auto eta1 = caloParticles[0].eta();
202  auto phi1 = caloParticles[0].phi();
203  auto theta1 = 2.*atan(exp(-eta1));
204  auto eta2 = caloParticles[1].eta();
205  auto phi2 = caloParticles[1].phi();
206  auto theta2 = 2.*atan(exp(-eta2));
207  eta1_->Fill(eta1);
208  eta2_->Fill(eta2);
209 
210 
211  // Select event only if the sum of the energy of its recHits
212  // is close enough to the gen energy
213  int count = 0;
214  int size = 0;
215  float energy = 0.;
216  float energy_tmp = 0.;
217  for (const auto& it_caloPart : caloParticles) {
218  count++;
219  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
220  size += simClusterRefVector.size();
221  for (const auto& it_sc : simClusterRefVector) {
222  const SimCluster& simCluster = (*(it_sc));
223  const std::vector<std::pair<uint32_t, float> >& hits_and_fractions =
224  simCluster.hits_and_fractions();
225  for (const auto& it_haf : hits_and_fractions) {
226  if (hitmap.count(it_haf.first))
227  energy += hitmap[it_haf.first]->energy()*it_haf.second;
228  } //hits and fractions
229  } // simcluster
230  if (count == 1) {
231  energy1_->Fill(energy);
232  energy_tmp = energy;
233  } else {
234  energy2_->Fill(energy-energy_tmp);
235  }
236  } // caloParticle
237  energytot_->Fill(energy);
238  if (filterOnEnergyAndCaloP_ && (energy < 2.*0.8*80 or size !=2))
239  return;
240 
241  deltaEtaPhi_->Fill(eta1-eta2, phi1-phi2);
242 
243  for (const auto& it_caloPart : caloParticles) {
244  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
245  IfLogTrace(debug_ > 0, "HGCalShowerSeparation")
246  << ">>> " << simClusterRefVector.size() << std::endl;
247  for (const auto& it_sc : simClusterRefVector) {
248  const SimCluster& simCluster = (*(it_sc));
249  if (simCluster.energy() < 80*0.8)
250  continue;
251  scEnergy_->Fill(simCluster.energy());
252  IfLogTrace(debug_ > 1, "HGCalShowerSeparation")
253  << ">>> SC.energy(): " << simCluster.energy()
254  << " SC.simEnergy(): " << simCluster.simEnergy()
255  << std::endl;
256  const std::vector<std::pair<uint32_t, float> >& hits_and_fractions =
257  simCluster.hits_and_fractions();
258 
259  for (const auto& it_haf : hits_and_fractions) {
260  if (!hitmap.count(it_haf.first))
261  continue;
262  unsigned int hitlayer = recHitTools_.getLayerWithOffset(it_haf.first);
263  auto global = recHitTools_.getPosition(it_haf.first);
264  float globalx = global.x();
265  float globaly = global.y();
266  float globalz = global.z();
267  if (globalz == 0)
268  continue;
269  auto rho1 = globalz*tan(theta1);
270  auto rho2 = globalz*tan(theta2);
271  auto x1 = rho1*cos(phi1);
272  auto y1 = rho1*sin(phi1);
273  auto x2 = rho2*cos(phi2);
274  auto y2 = rho2*sin(phi2);
275  auto half_point_x = (x1+x2)/2.;
276  auto half_point_y = (y1+y2)/2.;
277  auto half_point = sqrt((x1-half_point_x)*(x1-half_point_x)+(y1-half_point_y)*(y1-half_point_y));
278  auto d_len = sqrt((x2-x1)*(x2-x1)+(y2-y1)*(y2-y1));
279  auto dn_x = (x2-x1)/d_len;
280  auto dn_y = (y2-y1)/d_len;
281  auto distance = (globalx-x1)*dn_x + (globaly - y1)*dn_y;
282  distance -= half_point;
283  auto idealDistance = sqrt((x1-x2)*(x1-x2) + (y1-y2)*(y1-y2));
284  if (hitmap.count(it_haf.first)) {
285  profileOnLayer_[hitlayer]->Fill(10.*(globalx-half_point_x),
286  10.*(globaly-half_point_y),
287  hitmap[it_haf.first]->energy()*it_haf.second);
288  profileOnLayer_[55]->Fill(10.*(globalx-half_point_x),
289  10.*(globaly-half_point_y),
290  hitmap[it_haf.first]->energy()*it_haf.second);
291  globalProfileOnLayer_[hitlayer]->Fill(globalx,
292  globaly,
293  hitmap[it_haf.first]->energy()*it_haf.second);
294  globalProfileOnLayer_[55]->Fill(globalx,
295  globaly,
296  hitmap[it_haf.first]->energy()*it_haf.second);
297  layerEnergy_->Fill(hitlayer, hitmap[it_haf.first]->energy());
298  layerDistance_->Fill(hitlayer, std::abs(10.*distance), hitmap[it_haf.first]->energy()*it_haf.second);
299  etaPhi_->Fill(global.eta(), global.phi());
300  distanceOnLayer_[hitlayer]->Fill(10.*distance);//,
301  idealDistanceOnLayer_[hitlayer]->Fill(10.*idealDistance);//,
302  idealDeltaXY_[hitlayer]->Fill(10.*(x1-x2), 10.*(y1-y2));//,
303  centers_[hitlayer]->Fill(10.*half_point_x, 10.*half_point_y);//,
304  IfLogTrace(debug_ > 0, "HGCalShowerSeparation")
305  << ">>> " << distance
306  << " " << hitlayer
307  << " " << hitmap[it_haf.first]->energy()*it_haf.second
308  << std::endl;
310  hitlayer,
311  hitmap[it_haf.first]->energy()*it_haf.second);
312  }
313  } // end simHit
314  } // end simCluster
315  } // end caloparticle
316  }
317 }
318 
319 // ------------ method fills 'descriptions' with the allowed parameters for the
320 // module ------------
322  edm::ConfigurationDescriptions& descriptions) {
324  desc.add<int>("debug", 1);
325  desc.add<bool>("filterOnEnergyAndCaloP", false);
326  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
327  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
328  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
329  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
330  descriptions.add("hgcalShowerSeparationDefault", desc);
331 }
332 
333 // define this as a plug-in
size
Write out results.
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
T getParameter(std::string const &) const
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< MonitorElement * > globalProfileOnLayer_
#define IfLogTrace(cond, cat)
std::vector< MonitorElement * > idealDeltaXY_
std::vector< MonitorElement * > profileOnLayer_
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:73
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
Monte Carlo truth information used for tracking validation.
Definition: SimCluster.h:29
T sqrt(T t)
Definition: SSEVec.h:18
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
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
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
std::vector< MonitorElement * > distanceOnLayer_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
ParameterDescriptionBase * add(U const &iLabel, T const &value)
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:196
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: DetId.h:18
hgcal::RecHitTools recHitTools_
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:283
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
void fillWithRecHits(std::map< DetId, const HGCRecHit * > &, DetId, unsigned int, float, int &, float &)
MonitorElement * layerDistance_
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:112
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void analyze(const edm::Event &, const edm::EventSetup &) override
HLT enums.
std::vector< MonitorElement * > centers_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
std::vector< MonitorElement * > idealDistanceOnLayer_
T x() const
Definition: PV3DBase.h:62
MonitorElement * showerProfile_
#define constexpr
Definition: Run.h:45
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
HGCalShowerSeparation(const edm::ParameterSet &)