CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelPhase1CompareVertexSoA.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 // Package: SiPixelPhase1CompareVertexSoA
3 // Class: SiPixelPhase1CompareVertexSoA
4 //
7 //
8 // Author: Suvankar Roy Chowdhury
9 //
17 // DQM Histograming
23 
24 class SiPixelPhase1CompareVertexSoA : public DQMEDAnalyzer {
25 public:
26  using IndToEdm = std::vector<uint16_t>;
28  ~SiPixelPhase1CompareVertexSoA() override = default;
29  void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
30  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
31  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
32 
33 private:
38  const float dzCut_;
47 };
48 
49 //
50 // constructors
51 //
52 
54  : tokenSoAVertexCPU_(consumes<ZVertexHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelVertexSrcCPU"))),
55  tokenSoAVertexGPU_(consumes<ZVertexHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelVertexSrcGPU"))),
56  tokenBeamSpot_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
57  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
58  dzCut_(iConfig.getParameter<double>("dzCut")) {}
59 
60 //
61 // -- Analyze
62 //
64  const auto& vsoaHandleCPU = iEvent.getHandle(tokenSoAVertexCPU_);
65  const auto& vsoaHandleGPU = iEvent.getHandle(tokenSoAVertexGPU_);
66  if (not vsoaHandleCPU or not vsoaHandleGPU) {
67  edm::LogWarning out("SiPixelPhase1CompareTrackSoA");
68  if (not vsoaHandleCPU) {
69  out << "reference (cpu) tracks not found; ";
70  }
71  if (not vsoaHandleGPU) {
72  out << "target (gpu) tracks not found; ";
73  }
74  out << "the comparison will not run.";
75  return;
76  }
77 
78  auto const& vsoaCPU = *vsoaHandleCPU->get();
79  int nVerticesCPU = vsoaCPU.nvFinal;
80  auto const& vsoaGPU = *vsoaHandleGPU->get();
81  int nVerticesGPU = vsoaGPU.nvFinal;
82 
83  auto bsHandle = iEvent.getHandle(tokenBeamSpot_);
84  float x0 = 0., y0 = 0., z0 = 0., dxdz = 0., dydz = 0.;
85  if (!bsHandle.isValid()) {
86  edm::LogWarning("PixelVertexProducer") << "No beamspot found. returning vertexes with (0,0,Z) ";
87  } else {
88  const reco::BeamSpot& bs = *bsHandle;
89  x0 = bs.x0();
90  y0 = bs.y0();
91  z0 = bs.z0();
92  dxdz = bs.dxdz();
93  dydz = bs.dydz();
94  }
95 
96  for (int ivc = 0; ivc < nVerticesCPU; ivc++) {
97  auto sic = vsoaCPU.sortInd[ivc];
98  auto zc = vsoaCPU.zv[sic];
99  auto xc = x0 + dxdz * zc;
100  auto yc = y0 + dydz * zc;
101  zc += z0;
102 
103  auto ndofCPU = vsoaCPU.ndof[sic];
104  auto chi2CPU = vsoaCPU.chi2[sic];
105 
106  const int32_t notFound = -1;
107  int32_t closestVtxidx = notFound;
108  float mindz = dzCut_;
109 
110  for (int ivg = 0; ivg < nVerticesGPU; ivg++) {
111  auto sig = vsoaGPU.sortInd[ivg];
112  auto zgc = vsoaGPU.zv[sig] + z0;
113  auto zDist = std::abs(zc - zgc);
114  //insert some matching condition
115  if (zDist > dzCut_)
116  continue;
117  if (mindz > zDist) {
118  mindz = zDist;
119  closestVtxidx = sig;
120  }
121  }
122  if (closestVtxidx == notFound)
123  continue;
124 
125  auto zg = vsoaGPU.zv[closestVtxidx];
126  auto xg = x0 + dxdz * zg;
127  auto yg = y0 + dydz * zg;
128  zg += z0;
129  auto ndofGPU = vsoaGPU.ndof[closestVtxidx];
130  auto chi2GPU = vsoaGPU.chi2[closestVtxidx];
131 
132  hx_->Fill(xc, xg);
133  hy_->Fill(yc, yg);
134  hz_->Fill(zc, zg);
135  hchi2_->Fill(chi2CPU, chi2GPU);
136  hchi2oNdof_->Fill(chi2CPU / ndofCPU, chi2GPU / ndofGPU);
137  hptv2_->Fill(vsoaCPU.ptv2[sic], vsoaGPU.ptv2[closestVtxidx]);
138  hntrks_->Fill(ndofCPU + 1, ndofGPU + 1);
139  }
140  hnVertex_->Fill(nVerticesCPU, nVerticesGPU);
141 }
142 
143 //
144 // -- Book Histograms
145 //
147  edm::Run const& iRun,
148  edm::EventSetup const& iSetup) {
149  //std::string top_folder = ""//
150  ibooker.cd();
152 
153  // FIXME: all the 2D correlation plots are quite heavy in terms of memory consumption, so a as soon as DQM supports either TH2I or THnSparse
154  // these should be moved to a less resource consuming format
155  hnVertex_ = ibooker.book2D("nVertex", "# of Vertex;CPU;GPU", 101, -0.5, 100.5, 101, -0.5, 100.5);
156  hx_ = ibooker.book2D("vx", "Vertez x;CPU;GPU", 20, -0.1, 0.1, 20, -0.1, 0.1);
157  hy_ = ibooker.book2D("vy", "Vertez y;CPU;GPU", 20, -0.1, 0.1, 20, -0.1, 0.1);
158  hz_ = ibooker.book2D("vz", "Vertez z;CPU;GPU", 30, -30., 30., 30, -30., 30.);
159  hchi2_ = ibooker.book2D("chi2", "Vertex chi-squared;CPU;GPU", 40, 0., 20., 40, 0., 20.);
160  hchi2oNdof_ = ibooker.book2D("chi2oNdof", "Vertex chi-squared/Ndof;CPU;GPU", 40, 0., 20., 40, 0., 20.);
161  hptv2_ = ibooker.book2D("ptsq", "Vertex p_T squared;CPU;GPU", 200, 0., 200., 200, 0., 200.);
162  hntrks_ = ibooker.book2D("ntrk", "#tracks associated;CPU;GPU", 100, -0.5, 99.5, 100, -0.5, 99.5);
163  hntrks_ = ibooker.book2D("ntrk", "#tracks associated;CPU;GPU", 100, -0.5, 99.5, 100, -0.5, 99.5);
164 }
165 
167  // monitorpixelVertexSoA
169  desc.add<edm::InputTag>("pixelVertexSrcCPU", edm::InputTag("pixelVerticesSoA@cpu"));
170  desc.add<edm::InputTag>("pixelVertexSrcGPU", edm::InputTag("pixelVerticesSoA@cuda"));
171  desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
172  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelVertexCompareSoAGPUvsCPU");
173  desc.add<double>("dzCut", 1.);
174  descriptions.addWithDefaultLabel(desc);
175 }
double z0() const
z coordinate
Definition: BeamSpot.h:65
float dydz
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< ZVertexHeterogeneous > tokenSoAVertexGPU_
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
float dxdz
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~SiPixelPhase1CompareVertexSoA() override=default
void Fill(long long x)
double dydz() const
dydz slope
Definition: BeamSpot.h:80
void bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
int iEvent
Definition: GenABIO.cc:224
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double dxdz() const
dxdz slope
Definition: BeamSpot.h:78
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
static const GlobalPoint notFound(0, 0, 0)
const edm::EDGetTokenT< ZVertexHeterogeneous > tokenSoAVertexCPU_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double y0() const
y coordinate
Definition: BeamSpot.h:63
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
SiPixelPhase1CompareVertexSoA(const edm::ParameterSet &)
Log< level::Warning, false > LogWarning
const edm::EDGetTokenT< reco::BeamSpot > tokenBeamSpot_
Definition: Run.h:45
double x0() const
x coordinate
Definition: BeamSpot.h:61