1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
|
// -*- C++ -*-
//
// Package: SiPixelPhase1TrackResiduals
// Class: SiPixelPhase1TrackResiduals
//
// Original Author: Marcel Schneider
#include "DQM/SiPixelPhase1Common/interface/SiPixelPhase1Base.h"
#include "Alignment/OfflineValidation/interface/TrackerValidationVariables.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
#include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
#include "Geometry/CommonTopologies/interface/PixelTopology.h"
#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
namespace {
class SiPixelPhase1TrackResiduals final : public SiPixelPhase1Base {
enum {
RESIDUAL_X,
RESIDUAL_Y,
RESONEDGE_X,
RESONEDGE_Y,
RESOTHERBAD_X,
RESOTHERBAD_Y,
NormRes_X,
NormRes_Y,
DRNR_X,
DRNR_Y
};
public:
explicit SiPixelPhase1TrackResiduals(const edm::ParameterSet& conf);
void analyze(const edm::Event&, const edm::EventSetup&) override;
private:
TrackerValidationVariables validator;
const bool applyVertexCut_;
const edm::InputTag vertexInputTag_;
const edm::EDGetTokenT<reco::VertexCollection> offlinePrimaryVerticesToken_;
const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomToken_;
};
SiPixelPhase1TrackResiduals::SiPixelPhase1TrackResiduals(const edm::ParameterSet& iConfig)
: SiPixelPhase1Base(iConfig),
validator(iConfig, consumesCollector()),
applyVertexCut_(iConfig.getUntrackedParameter<bool>("VertexCut", true)),
vertexInputTag_(iConfig.getParameter<edm::InputTag>("vertices")),
offlinePrimaryVerticesToken_(consumes<reco::VertexCollection>(vertexInputTag_)),
trackerGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()) {}
void SiPixelPhase1TrackResiduals::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
if (!checktrigger(iEvent, iSetup, DCS))
return;
edm::ESHandle<TrackerGeometry> tracker = iSetup.getHandle(trackerGeomToken_);
assert(tracker.isValid());
edm::Handle<reco::VertexCollection> vertices;
if (applyVertexCut_) {
iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
if (!vertices.isValid()) {
edm::LogWarning("SiPixelPhase1TrackResiduals")
<< "Vertex collection: " << vertexInputTag_ << " is not valid! Skipping the event.";
return;
}
if (vertices->empty()) {
edm::LogWarning("SiPixelPhase1TrackResiduals")
<< "Size of vertex collection: " << vertexInputTag_ << " is zero! Skipping the event.";
return;
}
}
std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
validator.fillTrackQuantities(
iEvent,
iSetup,
// tell the validator to only look at good tracks
[&](const reco::Track& track) -> bool {
return (!applyVertexCut_ ||
(track.pt() > 0.75 && std::abs(track.dxy(vertices->at(0).position())) < 5 * track.dxyError()));
},
vtracks);
for (auto& track : vtracks) {
for (auto& it : track.hits) {
auto id = DetId(it.rawDetId);
auto isPixel = id.subdetId() == PixelSubdetector::PixelBarrel || id.subdetId() == PixelSubdetector::PixelEndcap;
if (!isPixel)
continue;
histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
if (it.resXprimeErr != 0 && it.resYprimeErr != 0) {
histo[NormRes_X].fill(it.resXprime / it.resXprimeErr, id, &iEvent);
histo[NormRes_Y].fill(it.resYprime / it.resYprimeErr, id, &iEvent);
histo[DRNR_X].fill(it.resXprime / it.resXprimeErr, id, &iEvent);
histo[DRNR_Y].fill(it.resYprime / it.resYprimeErr, id, &iEvent);
}
if (it.isOnEdgePixel) {
histo[RESONEDGE_X].fill(it.resXprime, id, &iEvent);
histo[RESONEDGE_Y].fill(it.resYprime, id, &iEvent);
}
if (it.isOtherBadPixel) {
histo[RESOTHERBAD_X].fill(it.resXprime, id, &iEvent);
histo[RESOTHERBAD_Y].fill(it.resYprime, id, &iEvent);
}
}
}
}
} // namespace
DEFINE_FWK_MODULE(SiPixelPhase1TrackResiduals);
|