CMS 3D CMS Logo

BtlRecHitsValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Validation/MtdValidation
4 // Class: BtlRecHitsValidation
5 //
14 #include <string>
15 
20 
23 
28 
32 
37 
40 
41 struct MTDHit {
42  float energy;
43  float time;
44  float x_local;
45  float y_local;
46  float z_local;
47 };
48 
50 public:
51  explicit BtlRecHitsValidation(const edm::ParameterSet&);
52  ~BtlRecHitsValidation() override;
53 
54  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
55 
56 private:
57  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
58 
59  void analyze(const edm::Event&, const edm::EventSetup&) override;
60 
61  // ------------ member data ------------
62 
64  const float hitMinEnergy_;
65 
68 
69  // --- histograms declaration
70 
72 
75 
77 
83 
91 
96 };
97 
98 // ------------ constructor and destructor --------------
100  : folder_(iConfig.getParameter<std::string>("folder")),
101  hitMinEnergy_(iConfig.getParameter<double>("hitMinimumEnergy")) {
102  btlRecHitsToken_ = consumes<FTLRecHitCollection>(iConfig.getParameter<edm::InputTag>("inputTag"));
103  btlSimHitsToken_ = consumes<CrossingFrame<PSimHit> >(iConfig.getParameter<edm::InputTag>("simHitsTag"));
104 }
105 
107 
108 // ------------ method called for each event ------------
110  using namespace edm;
111  using namespace geant_units::operators;
112 
113  edm::ESHandle<MTDGeometry> geometryHandle;
114  iSetup.get<MTDDigiGeometryRecord>().get(geometryHandle);
115  const MTDGeometry* geom = geometryHandle.product();
116 
117  edm::ESHandle<MTDTopology> topologyHandle;
118  iSetup.get<MTDTopologyRcd>().get(topologyHandle);
119  const MTDTopology* topology = topologyHandle.product();
120 
121  auto btlRecHitsHandle = makeValid(iEvent.getHandle(btlRecHitsToken_));
122  auto btlSimHitsHandle = makeValid(iEvent.getHandle(btlSimHitsToken_));
123  MixCollection<PSimHit> btlSimHits(btlSimHitsHandle.product());
124 
125  // --- Loop over the BLT SIM hits
126  std::unordered_map<uint32_t, MTDHit> m_btlSimHits;
127  for (auto const& simHit : btlSimHits) {
128  // --- Use only hits compatible with the in-time bunch-crossing
129  if (simHit.tof() < 0 || simHit.tof() > 25.)
130  continue;
131 
132  DetId id = simHit.detUnitId();
133 
134  auto simHitIt = m_btlSimHits.emplace(id.rawId(), MTDHit()).first;
135 
136  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
137  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
138 
139  // --- Get the time of the first SIM hit in the cell
140  if ((simHitIt->second).time == 0 || simHit.tof() < (simHitIt->second).time) {
141  (simHitIt->second).time = simHit.tof();
142 
143  auto hit_pos = simHit.entryPoint();
144  (simHitIt->second).x_local = hit_pos.x();
145  (simHitIt->second).y_local = hit_pos.y();
146  (simHitIt->second).z_local = hit_pos.z();
147  }
148 
149  } // simHit loop
150 
151  // --- Loop over the BLT RECO hits
152 
153  unsigned int n_reco_btl = 0;
154 
155  for (const auto& recHit : *btlRecHitsHandle) {
156  BTLDetId detId = recHit.id();
157  DetId geoId = detId.geographicalId(static_cast<BTLDetId::CrysLayout>(topology->getMTDTopologyMode()));
158  const MTDGeomDet* thedet = geom->idToDet(geoId);
159  if (thedet == nullptr)
160  throw cms::Exception("BtlRecHitsValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
161  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
162  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
163  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
164 
165  Local3DPoint local_point(0., 0., 0.);
166  local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
167  const auto& global_point = thedet->toGlobal(local_point);
168 
169  meHitEnergy_->Fill(recHit.energy());
170  meHitTime_->Fill(recHit.time());
171 
172  meOccupancy_->Fill(global_point.z(), global_point.phi());
173 
174  meHitX_->Fill(global_point.x());
175  meHitY_->Fill(global_point.y());
176  meHitZ_->Fill(global_point.z());
177  meHitPhi_->Fill(global_point.phi());
178  meHitEta_->Fill(global_point.eta());
179 
180  meHitTvsE_->Fill(recHit.energy(), recHit.time());
181  meHitEvsPhi_->Fill(global_point.phi(), recHit.energy());
182  meHitEvsEta_->Fill(global_point.eta(), recHit.energy());
183  meHitEvsZ_->Fill(global_point.z(), recHit.energy());
184  meHitTvsPhi_->Fill(global_point.phi(), recHit.time());
185  meHitTvsEta_->Fill(global_point.eta(), recHit.time());
186  meHitTvsZ_->Fill(global_point.z(), recHit.time());
187 
188  // Resolution histograms
189  if (m_btlSimHits.count(detId.rawId()) == 1 && m_btlSimHits[detId.rawId()].energy > hitMinEnergy_) {
190  float time_res = recHit.time() - m_btlSimHits[detId.rawId()].time;
191  float energy_res = recHit.energy() - m_btlSimHits[detId.rawId()].energy;
192 
193  meTimeRes_->Fill(time_res);
194  meEnergyRes_->Fill(energy_res);
195 
196  meTresvsE_->Fill(m_btlSimHits[detId.rawId()].energy, time_res);
197  meEresvsE_->Fill(m_btlSimHits[detId.rawId()].energy, energy_res);
198  }
199 
200  n_reco_btl++;
201 
202  } // recHit loop
203 
204  if (n_reco_btl > 0)
205  meNhits_->Fill(log10(n_reco_btl));
206 }
207 
208 // ------------ method for histogram booking ------------
210  edm::Run const& run,
211  edm::EventSetup const& iSetup) {
212  ibook.setCurrentFolder(folder_);
213 
214  // --- histograms booking
215 
216  meNhits_ = ibook.book1D("BtlNhits", "Number of BTL RECO hits;log_{10}(N_{RECO})", 100, 0., 5.25);
217 
218  meHitEnergy_ = ibook.book1D("BtlHitEnergy", "BTL RECO hits energy;E_{RECO} [MeV]", 100, 0., 20.);
219  meHitTime_ = ibook.book1D("BtlHitTime", "BTL RECO hits ToA;ToA_{RECO} [ns]", 100, 0., 25.);
220 
221  meOccupancy_ = ibook.book2D(
222  "BtlOccupancy", "BTL RECO hits occupancy;Z_{RECO} [cm]; #phi_{RECO} [rad]", 65, -260., 260., 126, -3.15, 3.15);
223 
224  meHitX_ = ibook.book1D("BtlHitX", "BTL RECO hits X;X_{RECO} [cm]", 60, -120., 120.);
225  meHitY_ = ibook.book1D("BtlHitY", "BTL RECO hits Y;Y_{RECO} [cm]", 60, -120., 120.);
226  meHitZ_ = ibook.book1D("BtlHitZ", "BTL RECO hits Z;Z_{RECO} [cm]", 100, -260., 260.);
227  meHitPhi_ = ibook.book1D("BtlHitPhi", "BTL RECO hits #phi;#phi_{RECO} [rad]", 126, -3.15, 3.15);
228  meHitEta_ = ibook.book1D("BtlHitEta", "BTL RECO hits #eta;#eta_{RECO}", 100, -1.55, 1.55);
229 
230  meHitTvsE_ =
231  ibook.bookProfile("BtlHitTvsE", "BTL RECO ToA vs energy;E_{RECO} [MeV];ToA_{RECO} [ns]", 50, 0., 20., 0., 100.);
232  meHitEvsPhi_ = ibook.bookProfile(
233  "BtlHitEvsPhi", "BTL RECO energy vs #phi;#phi_{RECO} [rad];E_{RECO} [MeV]", 50, -3.15, 3.15, 0., 100.);
234  meHitEvsEta_ = ibook.bookProfile(
235  "BtlHitEvsEta", "BTL RECO energy vs #eta;#eta_{RECO};E_{RECO} [MeV]", 50, -1.55, 1.55, 0., 100.);
236  meHitEvsZ_ =
237  ibook.bookProfile("BtlHitEvsZ", "BTL RECO energy vs Z;Z_{RECO} [cm];E_{RECO} [MeV]", 50, -260., 260., 0., 100.);
238  meHitTvsPhi_ = ibook.bookProfile(
239  "BtlHitTvsPhi", "BTL RECO ToA vs #phi;#phi_{RECO} [rad];ToA_{RECO} [ns]", 50, -3.15, 3.15, 0., 100.);
240  meHitTvsEta_ =
241  ibook.bookProfile("BtlHitTvsEta", "BTL RECO ToA vs #eta;#eta_{RECO};ToA_{RECO} [ns]", 50, -1.6, 1.6, 0., 100.);
242  meHitTvsZ_ =
243  ibook.bookProfile("BtlHitTvsZ", "BTL RECO ToA vs Z;Z_{RECO} [cm];ToA_{RECO} [ns]", 50, -260., 260., 0., 100.);
244 
245  meTimeRes_ = ibook.book1D("BtlTimeRes", "BTL time resolution;T_{RECO} - T_{SIM} [ns]", 100, -0.5, 0.5);
246  meEnergyRes_ = ibook.book1D("BtlEnergyRes", "BTL energy resolution;E_{RECO} - E_{SIM} [MeV]", 100, -0.5, 0.5);
247 
248  meTresvsE_ = ibook.bookProfile(
249  "BtlTresvsE", "BTL time resolution vs E;E_{SIM} [MeV];T_{RECO}-T_{SIM} [ns]", 50, 0., 20., 0., 100.);
250  meEresvsE_ = ibook.bookProfile(
251  "BtlEresvsE", "BTL energy resolution vs E;E_{SIM} [MeV];E_{RECO}-E_{SIM} [MeV]", 50, 0., 20., 0., 100.);
252 }
253 
254 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
257 
258  desc.add<std::string>("folder", "MTD/BTL/RecHits");
259  desc.add<edm::InputTag>("inputTag", edm::InputTag("mtdRecHits", "FTLBarrel"));
260  desc.add<edm::InputTag>("simHitsTag", edm::InputTag("mix", "g4SimHitsFastTimerHitsBarrel"));
261  desc.add<double>("hitMinimumEnergy", 1.); // [MeV]
262 
263  descriptions.add("btlRecHits", desc);
264 }
265 
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
T getParameter(std::string const &) const
MonitorElement * meHitEnergy_
const std::string folder_
std::string folder_
constexpr NumType convertUnitsTo(long double desiredUnits, NumType val)
Definition: GeantUnits.h:87
CaloTopology const * topology(0)
MonitorElement * meHitEvsZ_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
virtual const Topology & topology() const
Definition: GeomDet.cc:67
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
MonitorElement * meHitPhi_
int getMTDTopologyMode() const
Definition: MTDTopology.h:73
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
BtlRecHitsValidation(const edm::ParameterSet &)
MonitorElement * meEnergyRes_
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * meHitTime_
int nrows() const override
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * meOccupancy_
MonitorElement * meHitEta_
void Fill(long long x)
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:547
const MTDGeomDet * idToDet(DetId) const override
Definition: MTDGeometry.cc:160
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
BTLDetId geographicalId(CrysLayout lay) const
Definition: BTLDetId.cc:163
MonitorElement * meHitTvsPhi_
MonitorElement * meEresvsE_
MonitorElement * meTimeRes_
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, char const *option="s")
Definition: DQMStore.cc:333
MonitorElement * meHitEvsPhi_
virtual const PixelTopology & specificTopology() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
LocalPoint pixelToModuleLocalPoint(const LocalPoint &plp, int row, int col) const
Definition: DetId.h:17
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
example_stream void bookHistograms(DQMStore::IBooker &,@example_stream edm::Run const &,@example_stream edm::EventSetup const &) override
edm::EDGetTokenT< FTLRecHitCollection > btlRecHitsToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
MonitorElement * meTresvsE_
HLT enums.
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Definition: DQMStore.cc:266
T get() const
Definition: EventSetup.h:73
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:18
edm::EDGetTokenT< CrossingFrame< PSimHit > > btlSimHitsToken_
int column(unsigned nrows=16) const
Definition: BTLDetId.h:110
MonitorElement * meHitTvsE_
auto makeValid(const U &iOtherHandleType) noexcept(false)
Definition: ValidHandle.h:52
MonitorElement * meHitEvsEta_
T const * product() const
Definition: ESHandle.h:86
MonitorElement * meHitTvsEta_
MonitorElement * meHitTvsZ_
Definition: Run.h:45
int row(unsigned nrows=16) const
Definition: BTLDetId.h:105