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
34 
35 //
36 // class declaration
37 //
38 
39 namespace siStripRawPrime {
41  public:
43  : h_barycenter_{nullptr}, h_width_{nullptr}, h_avgCharge_{nullptr}, h_isSaturated_{nullptr}, isBooked_{false} {}
44 
45  void fill(const SiStripApproximateCluster& cluster) {
46  h_barycenter_->Fill(cluster.barycenter());
47  h_width_->Fill(cluster.width());
48  h_avgCharge_->Fill(cluster.avgCharge());
49  h_isSaturated_->Fill(cluster.isSaturated() ? 1 : -1);
50  }
51 
53  ibook.setCurrentFolder(folder);
55  ibook.book1D("clusterBarycenter", "cluster barycenter;cluster barycenter;#clusters", 7680., 0., 7680.);
56  h_width_ = ibook.book1D("clusterWidth", "cluster width;cluster width;#clusters", 128, -0.5, 127.5);
57  h_avgCharge_ =
58  ibook.book1D("clusterAvgCharge", "average strip charge;average strip charge;#clusters", 256, -0.5, 255.5);
59  h_isSaturated_ = ibook.book1D("clusterSaturation", "cluster saturation;is saturated?;#clusters", 3, -1.5, 1.5);
60  h_isSaturated_->getTH1F()->GetXaxis()->SetBinLabel(1, "Not saturated");
61  h_isSaturated_->getTH1F()->GetXaxis()->SetBinLabel(3, "Saturated");
62  isBooked_ = true;
63  }
64 
65  bool isBooked() { return isBooked_; }
66 
67  private:
72  bool isBooked_;
73  };
74 } // namespace siStripRawPrime
75 
77 public:
79  ~SiStripMonitorApproximateCluster() override = default;
80 
81  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
82 
83 private:
84  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
85 
86  void analyze(const edm::Event&, const edm::EventSetup&) override;
87 
88  // ------------ member data ------------
92 
96 
97  // for comparisons
104 
105  // Event Data
109 
110  // Event Setup Data
112 };
113 
114 //
115 // constructors and destructor
116 //
118  : folder_(iConfig.getParameter<std::string>("folder")),
119  compareClusters_(iConfig.getParameter<bool>("compareClusters")),
120  // Poducer name of input StripClusterCollection
121  approxClustersToken_(consumes<SiStripApproximateClusterCollection>(
122  iConfig.getParameter<edm::InputTag>("ApproxClustersProducer"))) {
124  if (compareClusters_) {
126  consumes<edmNew::DetSetVector<SiStripCluster>>(iConfig.getParameter<edm::InputTag>("ClustersProducer"));
127  }
128  stripClusterCollection_ = nullptr;
129 }
130 
131 //
132 // member functions
133 //
134 
135 // ------------ method called for each event ------------
137  using namespace edm;
138 
139  const auto& tkGeom = &iSetup.getData(tkGeomToken_);
140  const auto tkDets = tkGeom->dets();
141 
142  // get collection of DetSetVector of clusters from Event
143  edm::Handle<SiStripApproximateClusterCollection> approx_cluster_detsetvector;
144  iEvent.getByToken(approxClustersToken_, approx_cluster_detsetvector);
145  if (!approx_cluster_detsetvector.isValid()) {
146  edm::LogError("SiStripMonitorApproximateCluster")
147  << "SiStripApproximate cluster collection is not valid!" << std::endl;
148 
149  // if approximate clusters collection not available, then early return
150  return;
151  }
152 
153  // if requested to perform the comparison
154  if (compareClusters_) {
155  // get collection of DetSetVector of clusters from Event
157  iEvent.getByToken(stripClustersToken_, cluster_detsetvector);
158  if (!cluster_detsetvector.isValid()) {
159  edm::LogError("SiStripMonitorApproximateCluster")
160  << "Requested to perform comparison, but regular SiStrip cluster collection is not valid!" << std::endl;
161  return;
162  } else {
163  stripClusterCollection_ = cluster_detsetvector.product();
164  }
165  }
166 
167  int nApproxClusters{0};
168  const SiStripApproximateClusterCollection* clusterCollection = approx_cluster_detsetvector.product();
169 
170  for (const auto& detClusters : *clusterCollection) {
171  edmNew::DetSet<SiStripCluster> strip_clusters_detset;
172  const auto& detid = detClusters.id(); // get the detid of the current detset
173 
174  // starts here comaparison with regular clusters
175  if (compareClusters_) {
177  stripClusterCollection_->find(detid); // search clusters of same detid
178  strip_clusters_detset = (*isearch);
179  }
180 
181  for (const auto& cluster : detClusters) {
182  nApproxClusters++;
183 
184  // fill the full cluster collection
185  allClusters.fill(cluster);
186 
187  if (compareClusters_ && !strip_clusters_detset.empty()) {
188  // build the converted cluster for the matching
189  uint16_t nStrips{0};
190  auto det = std::find_if(tkDets.begin(), tkDets.end(), [detid](auto& elem) -> bool {
191  return (elem->geographicalId().rawId() == detid);
192  });
193  const StripTopology& p = dynamic_cast<const StripGeomDetUnit*>(*det)->specificTopology();
194  nStrips = p.nstrips() - 1;
195 
196  const auto convertedCluster = SiStripCluster(cluster, nStrips);
197 
198  float distance{9999.};
199  const SiStripCluster* closestCluster{nullptr};
200  for (const auto& stripCluster : strip_clusters_detset) {
201  // by construction the approximated cluster width has same
202  // size as the original cluster
203 
204  if (cluster.width() != stripCluster.size()) {
205  continue;
206  }
207 
208  float deltaBarycenter = convertedCluster.barycenter() - stripCluster.barycenter();
209  if (std::abs(deltaBarycenter) < distance) {
210  closestCluster = &stripCluster;
211  distance = std::abs(deltaBarycenter);
212  }
213  }
214 
215  // Matching criteria:
216  // - if exists a closest cluster in the DetId
217  // - the size coincides with the original one
218  if (closestCluster) {
219  // comparisong plots
220  h_deltaBarycenter_->Fill(closestCluster->barycenter() - convertedCluster.barycenter());
221  h_deltaSize_->Fill(closestCluster->size() - convertedCluster.size());
222  h_deltaCharge_->Fill(closestCluster->charge() - convertedCluster.charge());
223  h_deltaFirstStrip_->Fill(closestCluster->firstStrip() - convertedCluster.firstStrip());
224  h_deltaEndStrip_->Fill(closestCluster->endStrip() - convertedCluster.endStrip());
225  // monitoring plots
226  matchedClusters.fill(cluster);
227  h_isMatched_->Fill(1);
228  } else {
229  // monitoring plots
230  unMatchedClusters.fill(cluster);
231  h_isMatched_->Fill(-1);
232  }
233  } // if we're doing the comparison cluster by cluster
234 
235  } // loop on clusters in a detset
236  } // loop on the detset vector
237 
238  h_nclusters_->Fill(nApproxClusters);
239 }
240 
242  edm::Run const& run,
243  edm::EventSetup const& iSetup) {
244  ibook.setCurrentFolder(folder_);
245  h_nclusters_ = ibook.book1D("numberOfClusters", "total N. of clusters;N. of clusters;#clusters", 500., 0., 500000.);
246  allClusters.book(ibook, folder_);
247 
248  // for comparisons
249  if (compareClusters_) {
250  // book monitoring for matche and unmatched clusters separately
251  matchedClusters.book(ibook, fmt::format("{}/MatchedClusters", folder_));
252  unMatchedClusters.book(ibook, fmt::format("{}/UnmatchedClusters", folder_));
253 
254  ibook.setCurrentFolder(fmt::format("{}/ClusterComparisons", folder_));
256  ibook.book1D("deltaBarycenter", "#Delta barycenter;#Delta barycenter;cluster pairs", 101, -50.5, 50.5);
257  h_deltaSize_ = ibook.book1D("deltaSize", "#Delta size;#Delta size;cluster pairs", 201, -100.5, 100.5);
258  h_deltaCharge_ = ibook.book1D("deltaCharge", "#Delta charge;#Delta charge;cluster pairs", 401, -200.5, 200.5);
259 
261  ibook.book1D("deltaFirstStrip", "#Delta FirstStrip; #Delta firstStrip;cluster pairs", 101, -50.5, 50.5);
263  ibook.book1D("deltaEndStrip", "#Delta EndStrip; #Delta endStrip; cluster pairs", 101, -50.5, 50.5);
264 
265  h_isMatched_ = ibook.book1D("isClusterMatched", "cluster matching;is matched?;#clusters", 3, -1.5, 1.5);
266  h_isMatched_->getTH1F()->GetXaxis()->SetBinLabel(1, "Not matched");
267  h_isMatched_->getTH1F()->GetXaxis()->SetBinLabel(3, "Matched");
268  }
269 }
270 
271 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
274  desc.setComment("Monitor SiStripApproximateCluster collection and compare with regular SiStrip clusters");
275  desc.add<bool>("compareClusters", false)->setComment("if true, will compare with regualr Strip clusters");
276  desc.add<edm::InputTag>("ApproxClustersProducer", edm::InputTag("hltSiStripClusters2ApproxClusters"))
277  ->setComment("approxmate clusters collection");
278  desc.add<edm::InputTag>("ClustersProducer", edm::InputTag("hltSiStripClusterizerForRawPrime"))
279  ->setComment("regular clusters collection");
280  desc.add<std::string>("folder", "SiStripApproximateClusters")->setComment("Top Level Folder");
281  descriptions.addWithDefaultLabel(desc);
282 }
283 
284 // 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:303
id_type id() const
Definition: DetSetNew.h:61
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
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
edm::EDGetTokenT< SiStripApproximateClusterCollection > approxClustersToken_
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
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