CMS 3D CMS Logo

HGCalShowerSeparation.cc
Go to the documentation of this file.
1 // user include files
2 
6 
9 
11 
13 
20 
25 
28 
30 
31 #include <map>
32 #include <array>
33 #include <string>
34 #include <numeric>
35 
37 public:
39  ~HGCalShowerSeparation() override = default;
40 
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
42 
43 private:
44  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
45  void analyze(const edm::Event&, const edm::EventSetup&) override;
46 
47  void fillWithRecHits(std::unordered_map<DetId, const unsigned int>&, DetId, unsigned int, float, int&, float&);
48 
52  std::vector<edm::InputTag> hits_label;
53  std::vector<edm::EDGetTokenT<HGCRecHitCollection>> hits_token_;
54 
55  int debug_;
58  std::vector<HGCRecHit> hits_;
59 
71  std::vector<MonitorElement*> profileOnLayer_;
72  std::vector<MonitorElement*> globalProfileOnLayer_;
73  std::vector<MonitorElement*> distanceOnLayer_;
74  std::vector<MonitorElement*> idealDistanceOnLayer_;
75  std::vector<MonitorElement*> idealDeltaXY_;
76  std::vector<MonitorElement*> centers_;
77 
78  static constexpr int layers_ = 52;
79 };
80 
83  debug_(iConfig.getParameter<int>("debug")),
84  filterOnEnergyAndCaloP_(iConfig.getParameter<bool>("filterOnEnergyAndCaloP")) {
85  auto hitMapInputTag = iConfig.getParameter<edm::InputTag>("hitMapTag");
86  auto caloParticles = iConfig.getParameter<edm::InputTag>("caloParticles");
87  hitMap_ = consumes<std::unordered_map<DetId, const unsigned int>>(hitMapInputTag);
88  caloParticles_ = consumes<std::vector<CaloParticle>>(caloParticles);
89  auto hits_label_ = iConfig.getParameter<std::vector<edm::InputTag>>("hits");
90 
91  for (auto& label : hits_label_) {
92  hits_token_.push_back(consumes<HGCRecHitCollection>(label));
93  }
94 }
95 
97  ibooker.cd();
98  ibooker.setCurrentFolder("HGCalShowerSeparation");
99  scEnergy_ = ibooker.book1D("SCEnergy", "SCEnergy", 240, 0., 120.);
100  eta1_ = ibooker.book1D("eta1", "eta1", 80, 0., 4.);
101  eta2_ = ibooker.book1D("eta2", "eta2", 80, 0., 4.);
102  energy1_ = ibooker.book1D("energy1", "energy1", 240, 0., 120.);
103  energy2_ = ibooker.book1D("energy2", "energy2", 240, 0., 120.);
104  energytot_ = ibooker.book1D("energytot", "energytot", 200, 100., 200.);
105  showerProfile_ = ibooker.book2D("ShowerProfile", "ShowerProfile", 800, -400., 400., layers_, 0., (float)layers_);
106  layerEnergy_ = ibooker.book2D("LayerEnergy", "LayerEnergy", 60, 0., 60., 50, 0., 0.1);
107  layerDistance_ = ibooker.book2D("LayerDistance", "LayerDistance", 60, 0., 60., 400, -400., 400.);
108  etaPhi_ = ibooker.book2D("EtaPhi", "EtaPhi", 800, -4., 4., 800, -4., 4.);
109  deltaEtaPhi_ = ibooker.book2D("DeltaEtaPhi", "DeltaEtaPhi", 100, -0.5, 0.5, 100, -0.5, 0.5);
110  for (int i = 0; i < layers_; ++i) {
111  profileOnLayer_.push_back(ibooker.book2D(std::string("ProfileOnLayer_") + std::to_string(i),
112  std::string("ProfileOnLayer_") + std::to_string(i),
113  120,
114  -600.,
115  600.,
116  120,
117  -600.,
118  600.));
119  globalProfileOnLayer_.push_back(ibooker.book2D(std::string("GlobalProfileOnLayer_") + std::to_string(i),
120  std::string("GlobalProfileOnLayer_") + std::to_string(i),
121  320,
122  -160.,
123  160.,
124  320,
125  -160.,
126  160.));
127  distanceOnLayer_.push_back(ibooker.book1D(std::string("DistanceOnLayer_") + std::to_string(i),
128  std::string("DistanceOnLayer_") + std::to_string(i),
129  120,
130  -600.,
131  600.));
132  idealDistanceOnLayer_.push_back(ibooker.book1D(std::string("IdealDistanceOnLayer_") + std::to_string(i),
133  std::string("IdealDistanceOnLayer_") + std::to_string(i),
134  120,
135  -600.,
136  600.));
137  idealDeltaXY_.push_back(ibooker.book2D(std::string("IdealDeltaXY_") + std::to_string(i),
138  std::string("IdealDeltaXY_") + std::to_string(i),
139  800,
140  -400.,
141  400.,
142  800,
143  -400.,
144  400.));
145  centers_.push_back(ibooker.book2D(std::string("Centers_") + std::to_string(i),
146  std::string("Centers_") + std::to_string(i),
147  320,
148  -1600.,
149  1600.,
150  320,
151  -1600.,
152  1600.));
153  }
154 }
155 
158 
159  for (auto& token : hits_token_) {
161  iEvent.getByToken(token, hits_handle);
162  hits_.insert(hits_.end(), (*hits_handle).begin(), (*hits_handle).end());
163  }
164 
165  const edm::Handle<std::vector<CaloParticle>>& caloParticleHandle = iEvent.getHandle(caloParticles_);
166  const std::vector<CaloParticle>& caloParticles = *(caloParticleHandle.product());
167 
169  iEvent.getByToken(hitMap_, hitMapHandle);
170  const auto hitmap = *hitMapHandle;
171 
172  // loop over caloParticles
173  IfLogTrace(debug_ > 0, "HGCalShowerSeparation") << "Number of caloParticles: " << caloParticles.size() << std::endl;
174  if (caloParticles.size() == 2) {
175  auto eta1 = caloParticles[0].eta();
176  auto phi1 = caloParticles[0].phi();
177  auto theta1 = 2. * std::atan(exp(-eta1));
178  auto eta2 = caloParticles[1].eta();
179  auto phi2 = caloParticles[1].phi();
180  auto theta2 = 2. * std::atan(exp(-eta2));
181  eta1_->Fill(eta1);
182  eta2_->Fill(eta2);
183 
184  // Select event only if the sum of the energy of its recHits
185  // is close enough to the gen energy
186  int count = 0;
187  int size = 0;
188  float energy = 0.;
189  float energy_tmp = 0.;
190  for (const auto& it_caloPart : caloParticles) {
191  count++;
192  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
193  size += simClusterRefVector.size();
194  for (const auto& it_sc : simClusterRefVector) {
195  const SimCluster& simCluster = (*(it_sc));
196  const std::vector<std::pair<uint32_t, float>>& hits_and_fractions = simCluster.hits_and_fractions();
197  for (const auto& it_haf : hits_and_fractions) {
198  if (hitmap.count(it_haf.first)) {
199  const HGCRecHit* hit = &(hits_[hitmap.at(it_haf.first)]);
200  energy += hit->energy() * it_haf.second;
201  }
202  } //hits and fractions
203  } // simcluster
204  if (count == 1) {
205  energy1_->Fill(energy);
206  energy_tmp = energy;
207  } else {
208  energy2_->Fill(energy - energy_tmp);
209  }
210  } // caloParticle
212  if (filterOnEnergyAndCaloP_ && (energy < 2. * 0.8 * 80 or size != 2))
213  return;
214 
215  deltaEtaPhi_->Fill(eta1 - eta2, phi1 - phi2);
216 
217  for (const auto& it_caloPart : caloParticles) {
218  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
219  IfLogTrace(debug_ > 0, "HGCalShowerSeparation") << ">>> " << simClusterRefVector.size() << std::endl;
220  for (const auto& it_sc : simClusterRefVector) {
221  const SimCluster& simCluster = (*(it_sc));
222  if (simCluster.energy() < 80 * 0.8)
223  continue;
224  scEnergy_->Fill(simCluster.energy());
225  IfLogTrace(debug_ > 1, "HGCalShowerSeparation")
226  << ">>> SC.energy(): " << simCluster.energy() << " SC.simEnergy(): " << simCluster.simEnergy() << std::endl;
227  const std::vector<std::pair<uint32_t, float>>& hits_and_fractions = simCluster.hits_and_fractions();
228 
229  for (const auto& it_haf : hits_and_fractions) {
230  if (!hitmap.count(it_haf.first))
231  continue;
232  unsigned int hitlayer = recHitTools_.getLayerWithOffset(it_haf.first);
233  auto global = recHitTools_.getPosition(it_haf.first);
234  float globalx = global.x();
235  float globaly = global.y();
236  float globalz = global.z();
237  if (globalz == 0)
238  continue;
239  auto rho1 = globalz * tan(theta1);
240  auto rho2 = globalz * tan(theta2);
241  auto x1 = rho1 * cos(phi1);
242  auto y1 = rho1 * sin(phi1);
243  auto x2 = rho2 * cos(phi2);
244  auto y2 = rho2 * sin(phi2);
245  auto half_point_x = (x1 + x2) / 2.;
246  auto half_point_y = (y1 + y2) / 2.;
247  auto half_point = sqrt((x1 - half_point_x) * (x1 - half_point_x) + (y1 - half_point_y) * (y1 - half_point_y));
248  auto d_len = sqrt((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1));
249  auto dn_x = (x2 - x1) / d_len;
250  auto dn_y = (y2 - y1) / d_len;
251  auto distance = (globalx - x1) * dn_x + (globaly - y1) * dn_y;
252  distance -= half_point;
253  auto idealDistance = sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2));
254  if (hitmap.count(it_haf.first)) {
255  const HGCRecHit* hit = &(hits_[hitmap.at(it_haf.first)]);
256  profileOnLayer_[hitlayer]->Fill(
257  10. * (globalx - half_point_x), 10. * (globaly - half_point_y), hit->energy() * it_haf.second);
258  profileOnLayer_[55]->Fill(
259  10. * (globalx - half_point_x), 10. * (globaly - half_point_y), hit->energy() * it_haf.second);
260  globalProfileOnLayer_[hitlayer]->Fill(globalx, globaly, hit->energy() * it_haf.second);
261  globalProfileOnLayer_[55]->Fill(globalx, globaly, hit->energy() * it_haf.second);
262  layerEnergy_->Fill(hitlayer, hit->energy());
263  layerDistance_->Fill(hitlayer, std::abs(10. * distance), hit->energy() * it_haf.second);
264  etaPhi_->Fill(global.eta(), global.phi());
265  distanceOnLayer_[hitlayer]->Fill(10. * distance); //,
266  idealDistanceOnLayer_[hitlayer]->Fill(10. * idealDistance); //,
267  idealDeltaXY_[hitlayer]->Fill(10. * (x1 - x2), 10. * (y1 - y2)); //,
268  centers_[hitlayer]->Fill(10. * half_point_x, 10. * half_point_y); //,
269  IfLogTrace(debug_ > 0, "HGCalShowerSeparation")
270  << ">>> " << distance << " " << hitlayer << " " << hit->energy() * it_haf.second << std::endl;
271  showerProfile_->Fill(10. * distance, hitlayer, hit->energy() * it_haf.second);
272  }
273  } // end simHit
274  } // end simCluster
275  } // end caloparticle
276  }
277 }
278 
279 // ------------ method fills 'descriptions' with the allowed parameters for the
280 // module ------------
283  desc.add<int>("debug", 1);
284  desc.add<bool>("filterOnEnergyAndCaloP", false);
285  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
286  desc.add<edm::InputTag>("hitMapTag", edm::InputTag("recHitMapProducer", "hgcalRecHitMap"));
287  desc.add<std::vector<edm::InputTag>>("hits",
288  {edm::InputTag("HGCalRecHit", "HGCEERecHits"),
289  edm::InputTag("HGCalRecHit", "HGCHEFRecHits"),
290  edm::InputTag("HGCalRecHit", "HGCHEBRecHits")});
291  descriptions.add("hgcalShowerSeparationDefault", desc);
292 }
293 
294 // define this as a plug-in
size
Write out results.
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
#define IfLogTrace(cond, cat)
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tok_geom_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T const * product() const
Definition: Handle.h:70
std::vector< edm::InputTag > hits_label
static constexpr int layers_
std::vector< MonitorElement * > globalProfileOnLayer_
~HGCalShowerSeparation() override=default
std::vector< MonitorElement * > idealDeltaXY_
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:188
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:241
static std::string to_string(const XMLCh *ch)
std::vector< MonitorElement * > profileOnLayer_
void Fill(long long x)
char const * label
T x() const
Definition: PV3DBase.h:59
int iEvent
Definition: GenABIO.cc:224
Monte Carlo truth information used for tracking validation.
Definition: SimCluster.h:33
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:140
T sqrt(T t)
Definition: SSEVec.h:23
void fillWithRecHits(std::unordered_map< DetId, const unsigned int > &, DetId, unsigned int, float, int &, float &)
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
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: DetId.h:17
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
hgcal::RecHitTools recHitTools_
MonitorElement * layerDistance_
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void analyze(const edm::Event &, const edm::EventSetup &) override
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79
std::vector< edm::EDGetTokenT< HGCRecHitCollection > > hits_token_
std::vector< HGCRecHit > hits_
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:108
edm::EDGetTokenT< std::unordered_map< DetId, const unsigned int > > hitMap_
std::vector< MonitorElement * > centers_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
std::vector< MonitorElement * > idealDistanceOnLayer_
MonitorElement * showerProfile_
Definition: Run.h:45
HGCalShowerSeparation(const edm::ParameterSet &)
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381