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 {}
 
const edm::EDGetTokenT< SiStripApproximateClusterCollectionapproxClustersToken_
 
bool compareClusters_
 
std::string folder_
 
MonitorElementh2_CompareBarycenter_ {nullptr}
 
MonitorElementh2_CompareCharge_ {nullptr}
 
MonitorElementh2_CompareEndStrip_ {nullptr}
 
MonitorElementh2_CompareFirstStrip_ {nullptr}
 
MonitorElementh2_CompareSize_ {nullptr}
 
MonitorElementh_deltaBarycenter_ {nullptr}
 
MonitorElementh_deltaCharge_ {nullptr}
 
MonitorElementh_deltaEndStrip_ {nullptr}
 
MonitorElementh_deltaFirstStrip_ {nullptr}
 
MonitorElementh_deltaSize_ {nullptr}
 
MonitorElementh_isMatched_ {nullptr}
 
MonitorElementh_nclusters_
 
dqm::reco::MonitorElementh_numberFEDErrors_
 
siStripRawPrime::monitorApproxCluster matchedClusters {}
 
const edmNew::DetSetVector< SiStripCluster > * stripClusterCollection_
 
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
 
const edm::EDGetTokenT< DetIdVectorstripFEDErrorsToken_
 
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 101 of file SiStripMonitorApproximateCluster.cc.

Constructor & Destructor Documentation

◆ SiStripMonitorApproximateCluster()

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

Definition at line 153 of file SiStripMonitorApproximateCluster.cc.

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

154  : folder_(iConfig.getParameter<std::string>("folder")),
155  compareClusters_(iConfig.getParameter<bool>("compareClusters")),
156  // Poducer name of input StripClusterCollection
158  consumes<SiStripApproximateClusterCollection>(iConfig.getParameter<edm::InputTag>("ApproxClustersProducer"))),
159  stripFEDErrorsToken_(consumes<DetIdVector>(iConfig.getParameter<edm::InputTag>("SiStripFEDErrorVector"))) {
161  if (compareClusters_) {
163  consumes<edmNew::DetSetVector<SiStripCluster>>(iConfig.getParameter<edm::InputTag>("ClustersProducer"));
164  }
165  stripClusterCollection_ = nullptr;
166 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::EDGetTokenT< DetIdVector > stripFEDErrorsToken_
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
const edm::EDGetTokenT< SiStripApproximateClusterCollection > approxClustersToken_
const edmNew::DetSetVector< SiStripCluster > * stripClusterCollection_

◆ ~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 173 of file SiStripMonitorApproximateCluster.cc.

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

173  {
174  using namespace edm;
175 
176  const auto& tkGeom = &iSetup.getData(tkGeomToken_);
177  const auto tkDets = tkGeom->dets();
178 
179  // get SiStripApproximateClusterCollection from Event
180  edm::Handle<SiStripApproximateClusterCollection> approx_cluster_detsetvector;
181  iEvent.getByToken(approxClustersToken_, approx_cluster_detsetvector);
182  if (!approx_cluster_detsetvector.isValid()) {
183  edm::LogError("SiStripMonitorApproximateCluster")
184  << "SiStripApproximate cluster collection is not valid!" << std::endl;
185 
186  // if approximate clusters collection not available, then early return
187  return;
188  }
189 
190  // get the DetIdVector of the SiStrip FED Errors
191  edm::Handle<DetIdVector> sistripFEDErrorsVector;
192  iEvent.getByToken(stripFEDErrorsToken_, sistripFEDErrorsVector);
193  if (!sistripFEDErrorsVector.isValid()) {
194  edm::LogError("SiStripMonitorApproximateCluster")
195  << " DetIdVector collection of SiStrip FED errors is not valid!" << std::endl;
196 
197  // if approximate clusters collection not available, then early return
198  return;
199  }
200 
201  // if requested to perform the comparison
202  if (compareClusters_) {
203  // get collection of DetSetVector of clusters from Event
205  iEvent.getByToken(stripClustersToken_, cluster_detsetvector);
206  if (!cluster_detsetvector.isValid()) {
207  edm::LogError("SiStripMonitorApproximateCluster")
208  << "Requested to perform comparison, but regular SiStrip cluster collection is not valid!" << std::endl;
209  return;
210  } else {
211  stripClusterCollection_ = cluster_detsetvector.product();
212  }
213  }
214 
215  int nApproxClusters{0};
216  const SiStripApproximateClusterCollection* clusterCollection = approx_cluster_detsetvector.product();
217 
218  for (const auto& detClusters : *clusterCollection) {
219  edmNew::DetSet<SiStripCluster> strip_clusters_detset;
220  const auto& detid = detClusters.id(); // get the detid of the current detset
221 
222  // starts here comaparison with regular clusters
223  if (compareClusters_) {
225  edm::LogWarning("SiStripMonitorApproximateCluster")
226  << "Input SiStrip Cluster collecction was empty, skipping event! " << std::endl;
227  return;
228  }
229 
231  stripClusterCollection_->find(detid); // search clusters of same detid
232 
233  // protect against a missing match
234  if (isearch != stripClusterCollection_->end()) {
235  strip_clusters_detset = (*isearch);
236  } else {
237  edm::LogWarning("SiStripMonitorApproximateCluster")
238  << "No edmNew::DetSet<SiStripCluster> was found for detid " << detid << std::endl;
239  }
240  }
241 
242  for (const auto& cluster : detClusters) {
243  nApproxClusters++;
244 
245  // fill the full cluster collection
246  allClusters.fill(cluster);
247 
248  if (compareClusters_ && !strip_clusters_detset.empty()) {
249  // build the converted cluster for the matching
250  uint16_t nStrips{0};
251  auto det = std::find_if(tkDets.begin(), tkDets.end(), [detid](auto& elem) -> bool {
252  return (elem->geographicalId().rawId() == detid);
253  });
254  const StripTopology& p = dynamic_cast<const StripGeomDetUnit*>(*det)->specificTopology();
255  nStrips = p.nstrips() - 1;
256 
257  const auto convertedCluster = SiStripCluster(cluster, nStrips);
258 
259  float distance{9999.};
260  const SiStripCluster* closestCluster{nullptr};
261  for (const auto& stripCluster : strip_clusters_detset) {
262  // by construction the approximated cluster width has same
263  // size as the original cluster
264 
265  if (cluster.width() != stripCluster.size()) {
266  continue;
267  }
268 
269  float deltaBarycenter = convertedCluster.barycenter() - stripCluster.barycenter();
270  if (std::abs(deltaBarycenter) < distance) {
271  closestCluster = &stripCluster;
272  distance = std::abs(deltaBarycenter);
273  }
274  }
275 
276  // Matching criteria:
277  // - if exists a closest cluster in the DetId
278  // - the size coincides with the original one
279  if (closestCluster) {
280  // comparisong plots
281  h_deltaBarycenter_->Fill(closestCluster->barycenter() - convertedCluster.barycenter());
282  h_deltaSize_->Fill(closestCluster->size() - convertedCluster.size());
283  h_deltaCharge_->Fill(closestCluster->charge() - convertedCluster.charge());
284  h_deltaFirstStrip_->Fill(closestCluster->firstStrip() - convertedCluster.firstStrip());
285  h_deltaEndStrip_->Fill(closestCluster->endStrip() - convertedCluster.endStrip());
286 
287  h2_CompareBarycenter_->Fill(closestCluster->barycenter(), convertedCluster.barycenter());
288  h2_CompareSize_->Fill(closestCluster->size(), convertedCluster.size());
289  h2_CompareCharge_->Fill(closestCluster->charge(), convertedCluster.charge());
290  h2_CompareFirstStrip_->Fill(closestCluster->firstStrip(), convertedCluster.firstStrip());
291  h2_CompareEndStrip_->Fill(closestCluster->endStrip(), convertedCluster.endStrip());
292 
293  // monitoring plots
294  matchedClusters.fill(cluster);
295  h_isMatched_->Fill(1);
296  } else {
297  // monitoring plots
298  unMatchedClusters.fill(cluster);
299  h_isMatched_->Fill(-1);
300  }
301  } // if we're doing the comparison cluster by cluster
302 
303  } // loop on clusters in a detset
304  } // loop on the detset vector
305 
306  h_nclusters_->Fill(nApproxClusters);
307  h_numberFEDErrors_->Fill((*sistripFEDErrorsVector).size());
308 }
id_type id() const
Definition: DetSetNew.h:61
siStripRawPrime::monitorApproxCluster unMatchedClusters
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::EDGetTokenT< DetIdVector > stripFEDErrorsToken_
T const * product() const
Definition: Handle.h:70
Log< level::Error, false > LogError
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...
const_iterator end(bool update=false) const
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 edm::EDGetTokenT< SiStripApproximateClusterCollection > approxClustersToken_
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
Log< level::Warning, false > LogWarning
siStripRawPrime::monitorApproxCluster matchedClusters

◆ bookHistograms()

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

Implements DQMEDAnalyzer.

Definition at line 310 of file SiStripMonitorApproximateCluster.cc.

References allClusters, siStripRawPrime::monitorApproxCluster::book(), dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2I(), compareClusters_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), f, folder_, dqm-mbProfile::format, h2_CompareBarycenter_, h2_CompareCharge_, h2_CompareEndStrip_, h2_CompareFirstStrip_, h2_CompareSize_, h_deltaBarycenter_, h_deltaCharge_, h_deltaEndStrip_, h_deltaFirstStrip_, h_deltaSize_, h_isMatched_, h_nclusters_, h_numberFEDErrors_, matchedClusters, dqm::impl::MonitorElement::setBinLabel(), dqm::implementation::NavigatorBase::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, sistrip::STRIPS_PER_APV, and unMatchedClusters.

312  {
313  ibook.setCurrentFolder(folder_);
314  h_nclusters_ = ibook.book1D("numberOfClusters", "total N. of clusters;N. of clusters;# events", 500., 0., 500000.);
315  h_numberFEDErrors_ = ibook.book1D(
316  "numberOfFEDErrors", "number of SiStrip modules in FED Error;N. of Modules in Error; #events", 100, 0, 1000);
317 
318  allClusters.book(ibook, folder_);
319 
320  // for comparisons
321  if (compareClusters_) {
322  // book monitoring for matche and unmatched clusters separately
323  matchedClusters.book(ibook, fmt::format("{}/MatchedClusters", folder_));
324  unMatchedClusters.book(ibook, fmt::format("{}/UnmatchedClusters", folder_));
325 
326  // 1D histograms
327  ibook.setCurrentFolder(fmt::format("{}/ClusterComparisons", folder_));
329  ibook.book1D("deltaBarycenter", "#Delta barycenter;#Delta barycenter;cluster pairs", 101, -50.5, 50.5);
330  h_deltaSize_ = ibook.book1D("deltaSize", "#Delta size;#Delta size;cluster pairs", 201, -100.5, 100.5);
331  h_deltaCharge_ = ibook.book1D("deltaCharge", "#Delta charge;#Delta charge;cluster pairs", 401, -200.5, 200.5);
332 
334  ibook.book1D("deltaFirstStrip", "#Delta FirstStrip; #Delta firstStrip;cluster pairs", 101, -50.5, 50.5);
336  ibook.book1D("deltaEndStrip", "#Delta EndStrip; #Delta endStrip; cluster pairs", 101, -50.5, 50.5);
337 
338  // geometric constants
339  constexpr int maxNStrips = 6 * sistrip::STRIPS_PER_APV;
340  constexpr float minStrip = -0.5f;
341  constexpr float maxStrip = maxNStrips - 0.5f;
342 
343  // cluster constants
344  constexpr float maxCluSize = 100;
345  constexpr float maxCluCharge = 3000;
346 
347  // 2D histograms (use TH2I for counts to limit memory allocation)
348  std::string toRep = "SiStrip Cluster Barycenter";
349  h2_CompareBarycenter_ = ibook.book2I("compareBarycenter",
350  fmt::sprintf("; %s;Approx %s", toRep, toRep),
351  maxNStrips,
352  minStrip,
353  maxStrip,
354  maxNStrips,
355  minStrip,
356  maxStrip);
357 
358  toRep = "SiStrip Cluster Size";
359  h2_CompareSize_ = ibook.book2I("compareSize",
360  fmt::sprintf("; %s;Approx %s", toRep, toRep),
361  maxCluSize,
362  -0.5f,
363  maxCluSize - 0.5f,
364  maxCluSize,
365  -0.5f,
366  maxCluSize - 0.5f);
367 
368  toRep = "SiStrip Cluster Charge";
369  h2_CompareCharge_ = ibook.book2I("compareCharge",
370  fmt::sprintf("; %s;Approx %s", toRep, toRep),
371  (maxCluCharge / 5),
372  -0.5f,
373  maxCluCharge - 0.5f,
374  (maxCluCharge / 5),
375  -0.5f,
376  maxCluCharge - 0.5f);
377 
378  toRep = "SiStrip Cluster First Strip number";
379  h2_CompareFirstStrip_ = ibook.book2I("compareFirstStrip",
380  fmt::sprintf("; %s;Approx %s", toRep, toRep),
381  maxNStrips,
382  minStrip,
383  maxStrip,
384  maxNStrips,
385  minStrip,
386  maxStrip);
387 
388  toRep = "SiStrip Cluster Last Strip number";
389  h2_CompareEndStrip_ = ibook.book2I("compareLastStrip",
390  fmt::sprintf("; %s;Approx %s", toRep, toRep),
391  maxNStrips,
392  minStrip,
393  maxStrip,
394  maxNStrips,
395  minStrip,
396  maxStrip);
397 
398  h_isMatched_ = ibook.book1D("isClusterMatched", "cluster matching;is matched?;#clusters", 3, -1.5, 1.5);
399  h_isMatched_->setBinLabel(1, "Not matched");
400  h_isMatched_->setBinLabel(3, "Matched");
401  }
402 }
siStripRawPrime::monitorApproxCluster unMatchedClusters
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
double f[11][100]
virtual void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
siStripRawPrime::monitorApproxCluster allClusters
static const uint16_t STRIPS_PER_APV
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * book2I(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:296
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 405 of file SiStripMonitorApproximateCluster.cc.

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

405  {
407  desc.setComment("Monitor SiStripApproximateCluster collection and compare with regular SiStrip clusters");
408  desc.add<bool>("compareClusters", false)->setComment("if true, will compare with regualr Strip clusters");
409  desc.add<edm::InputTag>("ApproxClustersProducer", edm::InputTag("hltSiStripClusters2ApproxClusters"))
410  ->setComment("approxmate clusters collection");
411  desc.add<edm::InputTag>("SiStripFEDErrorVector", edm::InputTag("hltSiStripRawToDigi"))
412  ->setComment("SiStrip FED Errors collection");
413  desc.add<edm::InputTag>("ClustersProducer", edm::InputTag("hltSiStripClusterizerForRawPrime"))
414  ->setComment("regular clusters collection");
415  desc.add<std::string>("folder", "SiStripApproximateClusters")->setComment("Top Level Folder");
416  descriptions.addWithDefaultLabel(desc);
417 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ allClusters

siStripRawPrime::monitorApproxCluster SiStripMonitorApproximateCluster::allClusters {}
private

Definition at line 118 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ approxClustersToken_

const edm::EDGetTokenT<SiStripApproximateClusterCollection> SiStripMonitorApproximateCluster::approxClustersToken_
private

Definition at line 140 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze().

◆ compareClusters_

bool SiStripMonitorApproximateCluster::compareClusters_
private

◆ folder_

std::string SiStripMonitorApproximateCluster::folder_
private

Definition at line 114 of file SiStripMonitorApproximateCluster.cc.

Referenced by bookHistograms().

◆ h2_CompareBarycenter_

MonitorElement* SiStripMonitorApproximateCluster::h2_CompareBarycenter_ {nullptr}
private

Definition at line 133 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h2_CompareCharge_

MonitorElement* SiStripMonitorApproximateCluster::h2_CompareCharge_ {nullptr}
private

Definition at line 135 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h2_CompareEndStrip_

MonitorElement* SiStripMonitorApproximateCluster::h2_CompareEndStrip_ {nullptr}
private

Definition at line 137 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h2_CompareFirstStrip_

MonitorElement* SiStripMonitorApproximateCluster::h2_CompareFirstStrip_ {nullptr}
private

Definition at line 136 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h2_CompareSize_

MonitorElement* SiStripMonitorApproximateCluster::h2_CompareSize_ {nullptr}
private

Definition at line 134 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaBarycenter_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaBarycenter_ {nullptr}
private

Definition at line 127 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaCharge_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaCharge_ {nullptr}
private

Definition at line 129 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaEndStrip_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaEndStrip_ {nullptr}
private

Definition at line 131 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaFirstStrip_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaFirstStrip_ {nullptr}
private

Definition at line 130 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_deltaSize_

MonitorElement* SiStripMonitorApproximateCluster::h_deltaSize_ {nullptr}
private

Definition at line 128 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_isMatched_

MonitorElement* SiStripMonitorApproximateCluster::h_isMatched_ {nullptr}
private

Definition at line 126 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_nclusters_

MonitorElement* SiStripMonitorApproximateCluster::h_nclusters_
private

Definition at line 116 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ h_numberFEDErrors_

dqm::reco::MonitorElement* SiStripMonitorApproximateCluster::h_numberFEDErrors_
private

Definition at line 123 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().

◆ matchedClusters

siStripRawPrime::monitorApproxCluster SiStripMonitorApproximateCluster::matchedClusters {}
private

Definition at line 119 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

◆ stripFEDErrorsToken_

const edm::EDGetTokenT<DetIdVector> SiStripMonitorApproximateCluster::stripFEDErrorsToken_
private

Definition at line 141 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze().

◆ tkGeomToken_

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

◆ unMatchedClusters

siStripRawPrime::monitorApproxCluster SiStripMonitorApproximateCluster::unMatchedClusters {}
private

Definition at line 120 of file SiStripMonitorApproximateCluster.cc.

Referenced by analyze(), and bookHistograms().