Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Added pixel cluster ratio histograms #18071

Merged
merged 1 commit into from Apr 5, 2017
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
7 changes: 5 additions & 2 deletions DQM/SiPixelPhase1Clusters/interface/SiPixelPhase1Clusters.h
Expand Up @@ -10,6 +10,7 @@

#include "DQM/SiPixelPhase1Common/interface/SiPixelPhase1Base.h"
#include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h"

class SiPixelPhase1Clusters : public SiPixelPhase1Base {
enum {
Expand All @@ -26,7 +27,8 @@ class SiPixelPhase1Clusters : public SiPixelPhase1Base {
POSITION_YZ,
SIZE_VS_ETA,
READOUT_CHARGE,
READOUT_NCLUSTERS
READOUT_NCLUSTERS,
PIXEL_TO_STRIP_RATIO
};
// Uncomment to add trigger event flag enumerators
// Make sure enum corresponds correctly with flags defined in _cfi.py file
Expand All @@ -40,7 +42,8 @@ class SiPixelPhase1Clusters : public SiPixelPhase1Base {
void analyze(const edm::Event&, const edm::EventSetup&);

private:
edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > srcToken_;
edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > pixelSrcToken_;
edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster> > stripSrcToken_;
};

#endif
26 changes: 24 additions & 2 deletions DQM/SiPixelPhase1Clusters/python/SiPixelPhase1Clusters_cfi.py
Expand Up @@ -211,6 +211,26 @@
)
)

SiPixelPhase1ClustersPixelToStripRatio = DefaultHistoDigiCluster.clone(
name = "cluster_ratio",
title = "Pixel to Strip clusters ratio",

xlabel = "ratio",
dimensions = 1,

specs = VPSet(
Specification().groupBy("PXAll").save(100, 0, 1),
Specification().groupBy("PXAll/Lumisection")
.reduce("MEAN")
.groupBy("PXAll", "EXTEND_X")
.save(),
Specification().groupBy("PXAll/BX")
.reduce("MEAN")
.groupBy("PXAll", "EXTEND_X")
.save(),
)
)

SiPixelPhase1ClustersConf = cms.VPSet(
SiPixelPhase1ClustersCharge,
SiPixelPhase1ClustersSize,
Expand All @@ -225,7 +245,8 @@
SiPixelPhase1ClustersPositionYZ,
SiPixelPhase1ClustersSizeVsEta,
SiPixelPhase1ClustersReadoutCharge,
SiPixelPhase1ClustersReadoutNClusters
SiPixelPhase1ClustersReadoutNClusters,
SiPixelPhase1ClustersPixelToStripRatio
)

## Uncomment to add trigger event flag settings
Expand All @@ -236,7 +257,8 @@
)

SiPixelPhase1ClustersAnalyzer = cms.EDAnalyzer("SiPixelPhase1Clusters",
src = cms.InputTag("siPixelClusters"),
pixelSrc = cms.InputTag("siPixelClusters"),
stripSrc = cms.InputTag("siStripClusters"),
histograms = SiPixelPhase1ClustersConf,
geometry = SiPixelPhase1Geometry,
triggerflag = SiPixelPhase1ClustersTriggers,
Expand Down
22 changes: 17 additions & 5 deletions DQM/SiPixelPhase1Clusters/src/SiPixelPhase1Clusters.cc
Expand Up @@ -21,13 +21,25 @@
SiPixelPhase1Clusters::SiPixelPhase1Clusters(const edm::ParameterSet& iConfig) :
SiPixelPhase1Base(iConfig)
{
srcToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("src"));
pixelSrcToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelSrc"));

stripSrcToken_ = consumes<edmNew::DetSetVector<SiStripCluster>>(iConfig.getParameter<edm::InputTag>("stripSrc"));
}

void SiPixelPhase1Clusters::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
edm::Handle<edmNew::DetSetVector<SiPixelCluster>> input;
iEvent.getByToken(srcToken_, input);
if (!input.isValid()) return;
edm::Handle<edmNew::DetSetVector<SiPixelCluster>> inputPixel;
iEvent.getByToken(pixelSrcToken_, inputPixel);
if (!inputPixel.isValid()) return;

edm::Handle<edmNew::DetSetVector<SiStripCluster>> inputStrip;
iEvent.getByToken(stripSrcToken_, inputStrip);
if (inputStrip.isValid())
{
if (inputStrip.product()->data().size())
{
histo[PIXEL_TO_STRIP_RATIO].fill((double)inputPixel.product()->data().size() / (double) inputStrip.product()->data().size(), DetId(0), &iEvent);
}
}

bool hasClusters=false;

Expand All @@ -36,7 +48,7 @@ void SiPixelPhase1Clusters::analyze(const edm::Event& iEvent, const edm::EventSe
assert(tracker.isValid());

edmNew::DetSetVector<SiPixelCluster>::const_iterator it;
for (it = input->begin(); it != input->end(); ++it) {
for (it = inputPixel->begin(); it != inputPixel->end(); ++it) {
auto id = DetId(it->detId());

const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*> ( tracker->idToDet(id) );
Expand Down