CMS 3D CMS Logo

SiPixelCompareVertices.cc
Go to the documentation of this file.
1 // TODO: change file name to SiPixelCompareVerticesSoA.cc when CUDA code is removed
2 
3 // -*- C++ -*-
4 // Package: SiPixelCompareVertices
5 // Class: SiPixelCompareVertices
6 //
9 //
10 // Author: Suvankar Roy Chowdhury
11 //
12 
25 
26 // TODO: change class name to SiPixelCompareVerticesSoA when CUDA code is removed
28 public:
29  using IndToEdm = std::vector<uint16_t>;
31  ~SiPixelCompareVertices() override = default;
32  void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
33  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
34  // analyzeSeparate is templated to accept distinct types of SoAs
35  // The default use case is to use vertices from Alpaka reconstructed on CPU and GPU;
36  template <typename U, typename V>
37  void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event& iEvent);
38  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
39 
40 private:
41  // these two are both on Host but originally they have been produced on Host or on Device
46  const float dzCut_;
58 };
59 
60 //
61 // constructors
62 //
63 
65  : tokenSoAVertexReferenceSoA_(
66  consumes<ZVertexHost>(iConfig.getParameter<edm::InputTag>("pixelVertexReferenceSoA"))),
67  tokenSoAVertexTargetSoA_(consumes<ZVertexHost>(iConfig.getParameter<edm::InputTag>("pixelVertexTargetSoA"))),
68  tokenBeamSpot_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
69  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
70  dzCut_(iConfig.getParameter<double>("dzCut")) {}
71 
72 template <typename U, typename V>
73 void SiPixelCompareVertices::analyzeSeparate(U tokenRef, V tokenTar, const edm::Event& iEvent) {
74  const auto& vsoaHandleRef = iEvent.getHandle(tokenRef);
75  const auto& vsoaHandleTar = iEvent.getHandle(tokenTar);
76 
77  if (not vsoaHandleRef or not vsoaHandleTar) {
78  edm::LogWarning out("SiPixelCompareVertices");
79  if (not vsoaHandleRef) {
80  out << "reference vertices not found; ";
81  }
82  if (not vsoaHandleTar) {
83  out << "target vertices not found; ";
84  }
85  out << "the comparison will not run.";
86  return;
87  }
88 
89  auto const& vsoaRef = *vsoaHandleRef;
90  int nVerticesRef = vsoaRef.view().nvFinal();
91  auto const& vsoaTar = *vsoaHandleTar;
92  int nVerticesTar = vsoaTar.view().nvFinal();
93 
94  auto bsHandle = iEvent.getHandle(tokenBeamSpot_);
95  float x0 = 0., y0 = 0., z0 = 0., dxdz = 0., dydz = 0.;
96  if (not bsHandle.isValid()) {
97  edm::LogWarning("SiPixelCompareVertices") << "No beamspot found, returning vertexes with (0,0,Z).";
98  } else {
99  const reco::BeamSpot& bs = *bsHandle;
100  x0 = bs.x0();
101  y0 = bs.y0();
102  z0 = bs.z0();
103  dxdz = bs.dxdz();
104  dydz = bs.dydz();
105  }
106 
107  for (int ivc = 0; ivc < nVerticesRef; ivc++) {
108  auto sic = vsoaRef.view()[ivc].sortInd();
109  auto zc = vsoaRef.view()[sic].zv();
110  auto xc = x0 + dxdz * zc;
111  auto yc = y0 + dydz * zc;
112  zc += z0;
113 
114  auto ndofRef = vsoaRef.template view<reco::ZVertexTracksSoA>()[sic].ndof();
115  auto chi2Ref = vsoaRef.view()[sic].chi2();
116 
117  const int32_t notFound = -1;
118  int32_t closestVtxidx = notFound;
119  float mindz = dzCut_;
120 
121  for (int ivg = 0; ivg < nVerticesTar; ivg++) {
122  auto sig = vsoaTar.view()[ivg].sortInd();
123  auto zgc = vsoaTar.view()[sig].zv() + z0;
124  auto zDist = std::abs(zc - zgc);
125  //insert some matching condition
126  if (zDist > dzCut_)
127  continue;
128  if (mindz > zDist) {
129  mindz = zDist;
130  closestVtxidx = sig;
131  }
132  }
133  if (closestVtxidx == notFound)
134  continue;
135 
136  auto zg = vsoaTar.view()[closestVtxidx].zv();
137  auto xg = x0 + dxdz * zg;
138  auto yg = y0 + dydz * zg;
139  zg += z0;
140  auto ndofTar = vsoaTar.template view<reco::ZVertexTracksSoA>()[closestVtxidx].ndof();
141  auto chi2Tar = vsoaTar.view()[closestVtxidx].chi2();
142 
143  hx_->Fill(xc - x0, xg - x0);
144  hy_->Fill(yc - y0, yg - y0);
145  hz_->Fill(zc, zg);
146  hxdiff_->Fill(xc - xg);
147  hydiff_->Fill(yc - yg);
148  hzdiff_->Fill(zc - zg);
149  hchi2_->Fill(chi2Ref, chi2Tar);
150  hchi2oNdof_->Fill(chi2Ref / ndofRef, chi2Tar / ndofTar);
151  hptv2_->Fill(vsoaRef.view()[sic].ptv2(), vsoaTar.view()[closestVtxidx].ptv2());
152  hntrks_->Fill(ndofRef + 1, ndofTar + 1);
153  }
154  hnVertex_->Fill(nVerticesRef, nVerticesTar);
155 }
156 
157 //
158 // -- Analyze
159 //
161  // The default use case is to use vertices from Alpaka reconstructed on CPU and GPU;
162  // The function is left templated if any other cases need to be added
164 }
165 
166 //
167 // -- Book Histograms
168 //
170  edm::Run const& iRun,
171  edm::EventSetup const& iSetup) {
172  ibooker.cd();
174 
175  // 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
176  // these should be moved to a less resource consuming format
177  hnVertex_ = ibooker.book2I("nVertex", "# of Vertices;Reference;Target", 101, -0.5, 100.5, 101, -0.5, 100.5);
178  hx_ = ibooker.book2I("vx", "Vertez x - Beamspot x;Reference;Target", 50, -0.1, 0.1, 50, -0.1, 0.1);
179  hy_ = ibooker.book2I("vy", "Vertez y - Beamspot y;Reference;Target", 50, -0.1, 0.1, 50, -0.1, 0.1);
180  hz_ = ibooker.book2I("vz", "Vertez z;Reference;Target", 30, -30., 30., 30, -30., 30.);
181  hchi2_ = ibooker.book2I("chi2", "Vertex chi-squared;Reference;Target", 40, 0., 20., 40, 0., 20.);
182  hchi2oNdof_ = ibooker.book2I("chi2oNdof", "Vertex chi-squared/Ndof;Reference;Target", 40, 0., 20., 40, 0., 20.);
183  hptv2_ = ibooker.book2I("ptsq", "Vertex #sum (p_{T})^{2};Reference;Target", 200, 0., 200., 200, 0., 200.);
184  hntrks_ = ibooker.book2I("ntrk", "#tracks associated;Reference;Target", 100, -0.5, 99.5, 100, -0.5, 99.5);
185  hntrks_ = ibooker.book2I("ntrk", "#tracks associated;Reference;Target", 100, -0.5, 99.5, 100, -0.5, 99.5);
186  hxdiff_ = ibooker.book1D("vxdiff", ";Vertex x difference (Reference - Target);#entries", 100, -0.001, 0.001);
187  hydiff_ = ibooker.book1D("vydiff", ";Vertex y difference (Reference - Target);#entries", 100, -0.001, 0.001);
188  hzdiff_ = ibooker.book1D("vzdiff", ";Vertex z difference (Reference - Target);#entries", 100, -2.5, 2.5);
189 }
190 
192  // monitorpixelVertexSoA
194  desc.add<edm::InputTag>("pixelVertexReferenceSoA", edm::InputTag("pixelVerticesAlpakaSerial"));
195  desc.add<edm::InputTag>("pixelVertexTargetSoA", edm::InputTag("pixelVerticesAlpaka"));
196  desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
197  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelVertexCompareSoADeviceVSHost");
198  desc.add<double>("dzCut", 1.);
199  descriptions.addWithDefaultLabel(desc);
200 }
201 
202 // TODO: change module name to SiPixelCompareVerticesSoA when CUDA code is removed
float dydz
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const std::string topFolderName_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
float dxdz
SiPixelCompareVertices(const edm::ParameterSet &)
const edm::EDGetTokenT< reco::BeamSpot > tokenBeamSpot_
std::vector< uint16_t > IndToEdm
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t V
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::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
const edm::EDGetTokenT< ZVertexHost > tokenSoAVertexTargetSoA_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event &iEvent)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< ZVertexHost > tokenSoAVertexReferenceSoA_
static const GlobalPoint notFound(0, 0, 0)
fixed size matrix
HLT enums.
void bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
Log< level::Warning, false > LogWarning
~SiPixelCompareVertices() override=default
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:305
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: Run.h:45