CMS 3D CMS Logo

SiStripMonitorApproximateCluster.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: DQM/SiStripMonitorApproximateCluster
4 // Class: SiStripMonitorApproximateCluster
5 //
11 //
12 // Original Author: Marco Musich
13 // Created: Thu, 08 Dec 2022 20:51:10 GMT
14 //
15 //
16 
17 #include <string>
18 
19 // user include files
33 
34 //
35 // class declaration
36 //
37 
38 namespace siStripRawPrime {
40  public:
42  : h_barycenter_{nullptr}, h_width_{nullptr}, h_avgCharge_{nullptr}, h_isSaturated_{nullptr}, isBooked_{false} {}
43 
44  void fill(const SiStripApproximateCluster& cluster) {
45  h_barycenter_->Fill(cluster.barycenter());
46  h_width_->Fill(cluster.width());
47  h_avgCharge_->Fill(cluster.avgCharge());
48  h_isSaturated_->Fill(cluster.isSaturated() ? 1 : -1);
49  }
50 
52  ibook.setCurrentFolder(folder);
54  ibook.book1D("clusterBarycenter", "cluster barycenter;cluster barycenter;#clusters", 7680., 0., 7680.);
55  h_width_ = ibook.book1D("clusterWidth", "cluster width;cluster width;#clusters", 128, -0.5, 127.5);
56  h_avgCharge_ =
57  ibook.book1D("clusterAvgCharge", "average strip charge;average strip charge;#clusters", 256, -0.5, 255.5);
58  h_isSaturated_ = ibook.book1D("clusterSaturation", "cluster saturation;is saturated?;#clusters", 3, -1.5, 1.5);
59  h_isSaturated_->getTH1F()->GetXaxis()->SetBinLabel(1, "Not saturated");
60  h_isSaturated_->getTH1F()->GetXaxis()->SetBinLabel(3, "Saturated");
61  isBooked_ = true;
62  }
63 
64  bool isBooked() { return isBooked_; }
65 
66  private:
71  bool isBooked_;
72  };
73 } // namespace siStripRawPrime
74 
76 public:
78  ~SiStripMonitorApproximateCluster() override = default;
79 
80  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
81 
82 private:
83  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
84 
85  void analyze(const edm::Event&, const edm::EventSetup&) override;
86 
87  // ------------ member data ------------
91 
95 
96  // for comparisons
103 
104  // Event Data
108 
109  // Event Setup Data
111 };
112 
113 //
114 // constructors and destructor
115 //
117  : folder_(iConfig.getParameter<std::string>("folder")),
118  compareClusters_(iConfig.getParameter<bool>("compareClusters")),
119  // Poducer name of input StripClusterCollection
120  approxClustersToken_(consumes<edmNew::DetSetVector<SiStripApproximateCluster>>(
121  iConfig.getParameter<edm::InputTag>("ApproxClustersProducer"))) {
123  if (compareClusters_) {
125  consumes<edmNew::DetSetVector<SiStripCluster>>(iConfig.getParameter<edm::InputTag>("ClustersProducer"));
126  }
127  stripClusterCollection_ = nullptr;
128 }
129 
130 //
131 // member functions
132 //
133 
134 // ------------ method called for each event ------------
136  using namespace edm;
137 
138  const auto& tkGeom = &iSetup.getData(tkGeomToken_);
139  const auto tkDets = tkGeom->dets();
140 
141  // get collection of DetSetVector of clusters from Event
143  iEvent.getByToken(approxClustersToken_, approx_cluster_detsetvector);
144  if (!approx_cluster_detsetvector.isValid()) {
145  edm::LogError("SiStripMonitorApproximateCluster")
146  << "SiStripApproximate cluster collection is not valid!" << std::endl;
147 
148  // if approximate clusters collection not available, then early return
149  return;
150  }
151 
152  // if requested to perform the comparison
153  if (compareClusters_) {
154  // get collection of DetSetVector of clusters from Event
156  iEvent.getByToken(stripClustersToken_, cluster_detsetvector);
157  if (!cluster_detsetvector.isValid()) {
158  edm::LogError("SiStripMonitorApproximateCluster")
159  << "Requested to perform comparison, but regular SiStrip cluster collection is not valid!" << std::endl;
160  return;
161  } else {
162  stripClusterCollection_ = cluster_detsetvector.product();
163  }
164  }
165 
166  int nApproxClusters{0};
167  const edmNew::DetSetVector<SiStripApproximateCluster>* clusterCollection = approx_cluster_detsetvector.product();
168 
169  for (const auto& detClusters : *clusterCollection) {
170  edmNew::DetSet<SiStripCluster> strip_clusters_detset;
171  const auto& detid = detClusters.detId(); // get the detid of the current detset
172 
173  // starts here comaparison with regular clusters
174  if (compareClusters_) {
176  stripClusterCollection_->find(detid); // search clusters of same detid
177  strip_clusters_detset = (*isearch);
178  }
179 
180  for (const auto& cluster : detClusters) {
181  nApproxClusters++;
182 
183  // fill the full cluster collection
184  allClusters.fill(cluster);
185 
186  if (compareClusters_ && !strip_clusters_detset.empty()) {
187  // build the converted cluster for the matching
188  uint16_t nStrips{0};
189  auto det = std::find_if(tkDets.begin(), tkDets.end(), [detid](auto& elem) -> bool {
190  return (elem->geographicalId().rawId() == detid);
191  });
192  const StripTopology& p = dynamic_cast<const StripGeomDetUnit*>(*det)->specificTopology();
193  nStrips = p.nstrips() - 1;
194 
195  const auto convertedCluster = SiStripCluster(cluster, nStrips);
196 
197  float distance{9999.};
198  const SiStripCluster* closestCluster{nullptr};
199  for (const auto& stripCluster : strip_clusters_detset) {
200  // by construction the approximated cluster width has same
201  // size as the original cluster
202 
203  if (cluster.width() != stripCluster.size()) {
204  continue;
205  }
206 
207  float deltaBarycenter = convertedCluster.barycenter() - stripCluster.barycenter();
208  if (std::abs(deltaBarycenter) < distance) {
209  closestCluster = &stripCluster;
210  distance = std::abs(deltaBarycenter);
211  }
212  }
213 
214  // Matching criteria:
215  // - if exists a closest cluster in the DetId
216  // - the size coincides with the original one
217  if (closestCluster) {
218  // comparisong plots
219  h_deltaBarycenter_->Fill(closestCluster->barycenter() - convertedCluster.barycenter());
220  h_deltaSize_->Fill(closestCluster->size() - convertedCluster.size());
221  h_deltaCharge_->Fill(closestCluster->charge() - convertedCluster.charge());
222  h_deltaFirstStrip_->Fill(closestCluster->firstStrip() - convertedCluster.firstStrip());
223  h_deltaEndStrip_->Fill(closestCluster->endStrip() - convertedCluster.endStrip());
224  // monitoring plots
225  matchedClusters.fill(cluster);
226  h_isMatched_->Fill(1);
227  } else {
228  // monitoring plots
229  unMatchedClusters.fill(cluster);
230  h_isMatched_->Fill(-1);
231  }
232  } // if we're doing the comparison cluster by cluster
233 
234  } // loop on clusters in a detset
235  } // loop on the detset vector
236  h_nclusters_->Fill(nApproxClusters);
237 }
238 
240  edm::Run const& run,
241  edm::EventSetup const& iSetup) {
242  ibook.setCurrentFolder(folder_);
243  h_nclusters_ = ibook.book1D("numberOfClusters", "total N. of clusters;N. of clusters;#clusters", 500., 0., 500000.);
244  allClusters.book(ibook, folder_);
245 
246  // for comparisons
247  if (compareClusters_) {
248  // book monitoring for matche and unmatched clusters separately
249  matchedClusters.book(ibook, fmt::format("{}/MatchedClusters", folder_));
250  unMatchedClusters.book(ibook, fmt::format("{}/UnmatchedClusters", folder_));
251 
252  ibook.setCurrentFolder(fmt::format("{}/ClusterComparisons", folder_));
254  ibook.book1D("deltaBarycenter", "#Delta barycenter;#Delta barycenter;cluster pairs", 101, -50.5, 50.5);
255  h_deltaSize_ = ibook.book1D("deltaSize", "#Delta size;#Delta size;cluster pairs", 201, -100.5, 100.5);
256  h_deltaCharge_ = ibook.book1D("deltaCharge", "#Delta charge;#Delta charge;cluster pairs", 401, -200.5, 200.5);
257 
259  ibook.book1D("deltaFirstStrip", "#Delta FirstStrip; #Delta firstStrip;cluster pairs", 101, -50.5, 50.5);
261  ibook.book1D("deltaEndStrip", "#Delta EndStrip; #Delta endStrip; cluster pairs", 101, -50.5, 50.5);
262 
263  h_isMatched_ = ibook.book1D("isClusterMatched", "cluster matching;is matched?;#clusters", 3, -1.5, 1.5);
264  h_isMatched_->getTH1F()->GetXaxis()->SetBinLabel(1, "Not matched");
265  h_isMatched_->getTH1F()->GetXaxis()->SetBinLabel(3, "Matched");
266  }
267 }
268 
269 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
272  desc.setComment("Monitor SiStripApproximateCluster collection and compare with regular SiStrip clusters");
273  desc.add<bool>("compareClusters", false)->setComment("if true, will compare with regualr Strip clusters");
274  desc.add<edm::InputTag>("ApproxClustersProducer", edm::InputTag("hltSiStripClusters2ApproxClusters"))
275  ->setComment("approxmate clusters collection");
276  desc.add<edm::InputTag>("ClustersProducer", edm::InputTag("hltSiStripClusterizerForRawPrime"))
277  ->setComment("regular clusters collection");
278  desc.add<std::string>("folder", "SiStripApproximateClusters")->setComment("Top Level Folder");
279  descriptions.addWithDefaultLabel(desc);
280 }
281 
282 // define this as a plug-in
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
siStripRawPrime::monitorApproxCluster unMatchedClusters
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::string folder_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
T const * product() const
Definition: Handle.h:70
Log< level::Error, false > LogError
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool empty() const
Definition: DetSetNew.h:67
id_type detId() const
Definition: DetSetNew.h:63
nStrips
1.2 is to make the matching window safely the two nearest strips 0.35 is the size of an ME0 chamber i...
void Fill(long long x)
void fill(const SiStripApproximateCluster &cluster)
int iEvent
Definition: GenABIO.cc:224
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
SiStripMonitorApproximateCluster(const edm::ParameterSet &)
~SiStripMonitorApproximateCluster() override=default
siStripRawPrime::monitorApproxCluster allClusters
const edmNew::DetSetVector< SiStripCluster > * stripClusterCollection_
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
const_iterator find(id_type i, bool update=false) const
edm::EDGetTokenT< edmNew::DetSetVector< SiStripApproximateCluster > > approxClustersToken_
float barycenter() const
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
void book(dqm::implementation::DQMStore::IBooker &ibook, const std::string &folder)
Definition: Run.h:45
siStripRawPrime::monitorApproxCluster matchedClusters