CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
SiStripMonitorApproximateCluster Class Reference

#include <DQM/SiStripMonitorApproximateCluster/plugins/SiStripMonitorApproximateCluster.cc>

Inheritance diagram for SiStripMonitorApproximateCluster:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Member Functions

 SiStripMonitorApproximateCluster (const edm::ParameterSet &)
 
 ~SiStripMonitorApproximateCluster () override=default
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 

Private Attributes

siStripRawPrime::monitorApproxCluster allClusters {}
 
edm::EDGetTokenT< edmNew::DetSetVector< SiStripApproximateCluster > > approxClustersToken_
 
bool compareClusters_
 
std::string folder_
 
MonitorElementh_deltaBarycenter_ {nullptr}
 
MonitorElementh_deltaCharge_ {nullptr}
 
MonitorElementh_deltaEndStrip_ {nullptr}
 
MonitorElementh_deltaFirstStrip_ {nullptr}
 
MonitorElementh_deltaSize_ {nullptr}
 
MonitorElementh_isMatched_ {nullptr}
 
MonitorElementh_nclusters_
 
siStripRawPrime::monitorApproxCluster matchedClusters {}
 
const edmNew::DetSetVector< SiStripCluster > * stripClusterCollection_
 
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
 
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtkGeomToken_
 
siStripRawPrime::monitorApproxCluster unMatchedClusters {}
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

Description: Monitor SiStripApproximateClusters and on-demand compare properties with original SiStripClusters

Definition at line 75 of file SiStripMonitorApproximateCluster.cc.

Constructor & Destructor Documentation

◆ SiStripMonitorApproximateCluster()

SiStripMonitorApproximateCluster::SiStripMonitorApproximateCluster ( const edm::ParameterSet iConfig)
explicit

Definition at line 116 of file SiStripMonitorApproximateCluster.cc.

References compareClusters_, deDxTools::esConsumes(), edm::ParameterSet::getParameter(), stripClusterCollection_, stripClustersToken_, and tkGeomToken_.

117  : folder_(iConfig.getParameter<std::string>("folder")),
118  compareClusters_(iConfig.getParameter<bool>("compareClusters")),
119  // Poducer name of input StripClusterCollection
121  iConfig.getParameter<edm::InputTag>("ApproxClustersProducer"))) {
123  if (compareClusters_) {
125  consumes<edmNew::DetSetVector<SiStripCluster>>(iConfig.getParameter<edm::InputTag>("ClustersProducer"));
126  }
127  stripClusterCollection_ = nullptr;
128 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
const edmNew::DetSetVector< SiStripCluster > * stripClusterCollection_
edm::EDGetTokenT< edmNew::DetSetVector< SiStripApproximateCluster > > approxClustersToken_

◆ ~SiStripMonitorApproximateCluster()

SiStripMonitorApproximateCluster::~SiStripMonitorApproximateCluster ( )
overridedefault

Member Function Documentation

◆ analyze()

void SiStripMonitorApproximateCluster::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 135 of file SiStripMonitorApproximateCluster.cc.

References funct::abs(), allClusters, approxClustersToken_, SiStripCluster::barycenter(), ALCARECOSiPixelCalSingleMuonDQM_cff::clusterCollection, compareClusters_, edmNew::DetSet< T >::detId(), HLT_2023v12_cff::distance, edmNew::DetSet< T >::empty(), siStripRawPrime::monitorApproxCluster::fill(), dqm::impl::MonitorElement::Fill(), edmNew::DetSetVector< T >::find(), edm::EventSetup::getData(), h_deltaBarycenter_, h_deltaCharge_, h_deltaEndStrip_, h_deltaFirstStrip_, h_deltaSize_, h_isMatched_, h_nclusters_, iEvent, edm::HandleBase::isValid(), matchedClusters, me0TriggerPseudoDigis_cff::nStrips, AlCaHLTBitMon_ParallelJobs::p, edm::Handle< T >::product(), StripGeomDetUnit::specificTopology(), TrackingMonitor_cfi::stripCluster, stripClusterCollection_, stripClustersToken_, tkGeomToken_, and unMatchedClusters.

135  {
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 }
siStripRawPrime::monitorApproxCluster unMatchedClusters
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
T const * product() const
Definition: Handle.h:70
Log< level::Error, false > LogError
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
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
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
siStripRawPrime::monitorApproxCluster matchedClusters

◆ bookHistograms()

void SiStripMonitorApproximateCluster::bookHistograms ( DQMStore::IBooker ibook,
edm::Run const &  run,
edm::EventSetup const &  iSetup 
)
overrideprivatevirtual

Implements DQMEDAnalyzer.

Definition at line 239 of file SiStripMonitorApproximateCluster.cc.

References allClusters, siStripRawPrime::monitorApproxCluster::book(), dqm::implementation::IBooker::book1D(), compareClusters_, folder_, dqm-mbProfile::format, dqm::impl::MonitorElement::getTH1F(), h_deltaBarycenter_, h_deltaCharge_, h_deltaEndStrip_, h_deltaFirstStrip_, h_deltaSize_, h_isMatched_, h_nclusters_, matchedClusters, dqm::implementation::NavigatorBase::setCurrentFolder(), and unMatchedClusters.

241  {
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 }
siStripRawPrime::monitorApproxCluster unMatchedClusters
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
siStripRawPrime::monitorApproxCluster allClusters
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)
siStripRawPrime::monitorApproxCluster matchedClusters

◆ fillDescriptions()

void SiStripMonitorApproximateCluster::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 270 of file SiStripMonitorApproximateCluster.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

270  {
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 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ allClusters

siStripRawPrime::monitorApproxCluster SiStripMonitorApproximateCluster::allClusters {}
private

Definition at line 92 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ approxClustersToken_

edm::EDGetTokenT<edmNew::DetSetVector<SiStripApproximateCluster> > SiStripMonitorApproximateCluster::approxClustersToken_
private

Definition at line 105 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze().

◆ compareClusters_

bool SiStripMonitorApproximateCluster::compareClusters_
private

◆ folder_

std::string SiStripMonitorApproximateCluster::folder_
private

Definition at line 88 of file SiStripMonitorApproximateCluster.cc.

Referenced by bookHistograms().

◆ h_deltaBarycenter_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaBarycenter_ {nullptr}
private

Definition at line 98 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaCharge_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaCharge_ {nullptr}
private

Definition at line 100 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaEndStrip_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaEndStrip_ {nullptr}
private

Definition at line 102 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaFirstStrip_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaFirstStrip_ {nullptr}
private

Definition at line 101 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaSize_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaSize_ {nullptr}
private

Definition at line 99 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_isMatched_

MonitorElement* SiStripMonitorApproximateCluster::h_isMatched_ {nullptr}
private

Definition at line 97 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_nclusters_

MonitorElement* SiStripMonitorApproximateCluster::h_nclusters_
private

Definition at line 90 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ matchedClusters

siStripRawPrime::monitorApproxCluster SiStripMonitorApproximateCluster::matchedClusters {}
private

Definition at line 93 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ stripClusterCollection_

const edmNew::DetSetVector<SiStripCluster>* SiStripMonitorApproximateCluster::stripClusterCollection_
private

◆ stripClustersToken_

edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster> > SiStripMonitorApproximateCluster::stripClustersToken_
private

◆ tkGeomToken_

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> SiStripMonitorApproximateCluster::tkGeomToken_
private

◆ unMatchedClusters

siStripRawPrime::monitorApproxCluster SiStripMonitorApproximateCluster::unMatchedClusters {}
private

Definition at line 94 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().