CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MonitorTrackResiduals.cc
Go to the documentation of this file.
13 
14 template <TrackerType pixel_or_strip>
16  : conf_(iConfig),
17  tkDetMapToken_{esConsumes<TkDetMap, TrackerTopologyRcd, edm::Transition::BeginRun>()},
18  trackerTopologyRunToken_{esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()},
19  trackerGeometryToken_{esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()},
20  trackerTopologyEventToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
21  m_cacheID_(0),
23  iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"), consumesCollector(), *this)),
24  avalidator_(iConfig, consumesCollector()) {
25  applyVertexCut_ = conf_.getUntrackedParameter<bool>("VertexCut", true);
26  ModOn = conf_.getParameter<bool>("Mod_On");
27  offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(std::string("offlinePrimaryVertices"));
28 }
29 
30 template <TrackerType pixel_or_strip>
33  delete genTriggerEventFlag_;
34 }
35 
36 template <TrackerType pixel_or_strip>
38  const edm::Run &run,
39  const edm::EventSetup &iSetup) {
40  unsigned long long cacheID = iSetup.get<TrackerDigiGeometryRecord>().cacheIdentifier();
41  if (m_cacheID_ != cacheID) {
42  m_cacheID_ = cacheID;
43  this->createMEs(ibooker, iSetup);
44  }
45  std::string topFolderName_ = "SiStrip";
47  folder_organizer.setSiStripFolderName(topFolderName_);
48  const TkDetMap *tkDetMap = &iSetup.getData(tkDetMapToken_);
50  std::make_unique<TkHistoMap>(tkDetMap, ibooker, topFolderName_, "TkHMap_ResidualsMean", 0.0, true);
51 }
52 
53 template <TrackerType pixel_or_strip>
55  // Initialize the GenericTriggerEventFlag
56  if (genTriggerEventFlag_->on())
57  genTriggerEventFlag_->initRun(run, iSetup);
58 }
59 
60 template <TrackerType pixel_or_strip>
62  uint32_t ModuleID, const TrackerTopology *tTopo) {
63  std::string subdet = "";
64  int32_t layer = 0;
65  auto id = DetId(ModuleID);
66  switch (id.subdetId()) {
67  // Pixel Barrel, Endcap
68  case 1:
69  subdet = "BPIX";
70  layer = tTopo->pxbLayer(id);
71  break;
72  case 2:
73  subdet = "FPIX";
74  layer = tTopo->pxfDisk(id) * (tTopo->pxfSide(ModuleID) == 1 ? -1 : +1);
75  break;
76  // Strip TIB, TID, TOB, TEC
77  case 3:
78  subdet = "TIB";
79  layer = tTopo->tibLayer(id);
80  break;
81  case 4:
82  subdet = "TID";
83  layer = tTopo->tidWheel(id) * (tTopo->tidSide(ModuleID) == 1 ? -1 : +1);
84  break;
85  case 5:
86  subdet = "TOB";
87  layer = tTopo->tobLayer(id);
88  break;
89  case 6:
90  subdet = "TEC";
91  layer = tTopo->tecWheel(id) * (tTopo->tecSide(ModuleID) == 1 ? -1 : +1);
92  break;
93  default:
94  // TODO: Fail loudly.
95  subdet = "UNKNOWN";
96  layer = 0;
97  }
98  return std::make_pair(subdet, layer);
99 }
100 
101 template <TrackerType pixel_or_strip>
103  // Retrieve tracker topology and geometry
104  const TrackerTopology *const tTopo = &iSetup.getData(trackerTopologyRunToken_);
105  const TrackerGeometry *TG = &iSetup.getData(trackerGeometryToken_);
106 
107  Parameters = conf_.getParameter<edm::ParameterSet>("TH1ResModules");
108  int32_t i_residuals_Nbins = Parameters.getParameter<int32_t>("Nbinx");
109  double d_residual_xmin = Parameters.getParameter<double>("xmin");
110  double d_residual_xmax = Parameters.getParameter<double>("xmax");
111  Parameters = conf_.getParameter<edm::ParameterSet>("TH1NormResModules");
112  int32_t i_normres_Nbins = Parameters.getParameter<int32_t>("Nbinx");
113  double d_normres_xmin = Parameters.getParameter<double>("xmin");
114  double d_normres_xmax = Parameters.getParameter<double>("xmax");
115 
116  // use SistripHistoId for producing histogram id (and title)
118 
119  SiStripFolderOrganizer strip_organizer;
120  auto pixel_organizer = SiPixelFolderOrganizer(false);
121 
122  // Collect list of modules from Tracker Geometry
123  // book histo per each detector module
124  auto ids = TG->detIds(); // or detUnitIds?
125  for (DetId id : ids) {
126  auto ModuleID = id.rawId();
127  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
128  if (isPixel != (pixel_or_strip == TRACKERTYPE_PIXEL))
129  continue;
130 
131  // Book module histogramms?
132  if (ModOn) {
133  switch (id.subdetId()) {
134  case 1:
135  pixel_organizer.setModuleFolder(ibooker, ModuleID, 0);
136  break;
137  case 2:
138  pixel_organizer.setModuleFolder(ibooker, ModuleID, 0);
139  break;
140  default:
141  strip_organizer.setDetectorFolder(ModuleID, tTopo);
142  }
143  {
144  // this sounds strip specific but also works for pixel
145  std::string hid = hidmanager.createHistoId("HitResidualsX", "det", ModuleID);
146  std::string normhid = hidmanager.createHistoId("NormalizedHitResidualsX", "det", ModuleID);
147  auto &histos = m_ModuleResiduals[std::make_pair("", ModuleID)];
148  histos.x.base = ibooker.book1D(hid, hid, i_residuals_Nbins, d_residual_xmin, d_residual_xmax);
149  histos.x.base->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
150  histos.x.normed = ibooker.book1D(normhid, normhid, i_normres_Nbins, d_normres_xmin, d_normres_xmax);
151  histos.x.normed->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
152  }
153  {
154  std::string hid = hidmanager.createHistoId("HitResidualsY", "det", ModuleID);
155  std::string normhid = hidmanager.createHistoId("NormalizedHitResidualsY", "det", ModuleID);
156  auto &histos = m_ModuleResiduals[std::make_pair("", ModuleID)];
157  histos.y.base = ibooker.book1D(hid, hid, i_residuals_Nbins, d_residual_xmin, d_residual_xmax);
158  histos.y.base->setAxisTitle("(y_{pred} - y_{rec})' [cm]");
159  histos.y.normed = ibooker.book1D(normhid, normhid, i_normres_Nbins, d_normres_xmin, d_normres_xmax);
160  histos.y.normed->setAxisTitle("(y_{pred} - y_{rec})'/#sigma");
161  }
162  }
163 
164  auto subdetandlayer = findSubdetAndLayer(ModuleID, tTopo);
165  if (m_SubdetLayerResiduals.find(subdetandlayer) == m_SubdetLayerResiduals.end()) {
166  // add new histograms
167  auto &histos = m_SubdetLayerResiduals[subdetandlayer];
168  switch (id.subdetId()) {
169  // Pixel Barrel, Endcap
170  // We can't use the folder organizer here (SiPixelActionExecutor.cc#1638
171  // does the same)
172  case 1:
173  ibooker.setCurrentFolder("Pixel/Barrel");
174  break;
175  case 2:
176  ibooker.setCurrentFolder("Pixel/Endcap");
177  break;
178  // All strip
179  default:
180  strip_organizer.setLayerFolder(ModuleID, tTopo, subdetandlayer.second);
181  }
182 
183  auto isBarrel = subdetandlayer.first.find("B") != std::string::npos;
184 
185  auto xy = std::vector<std::pair<HistoPair &, const char *>>{std::make_pair(std::ref(histos.x), "X"),
186  std::make_pair(std::ref(histos.y), "Y")};
187  for (auto &histopair : xy) {
188  // book histogramms on layer level, check for barrel/pixel only for
189  // correct labeling
190 
191  // Skip the Y plots for strips.
192  if (!isPixel && histopair.second[0] == 'Y')
193  continue;
194 
195  std::string histoname = isPixel ? ( // Pixel name
196  Form("HitResiduals%s_%s%d",
197  histopair.second,
198  isBarrel ? "L" : (subdetandlayer.second > 0 ? "Dp" : "Dm"),
199  std::abs(subdetandlayer.second)))
200  : (Form("HitResiduals_%s__%s__%d", // Strip TODO: We use a
201  // legacy name.
202  subdetandlayer.first.c_str(),
203  isBarrel ? "Layer" : "wheel",
204  std::abs(subdetandlayer.second)));
205 
206  std::string histotitle = Form("HitResiduals %s on %s%s full %s %d",
207  histopair.second,
208  subdetandlayer.first.c_str(),
209  isBarrel ? "" : (subdetandlayer.second > 0 ? "+" : "-"),
210  isBarrel ? "Layer" : (isPixel ? "Disk" : "Wheel"),
211  std::abs(subdetandlayer.second));
212 
213  std::string normhistoname = Form("Normalized%s", histoname.c_str());
214  std::string normhistotitle = Form("Normalized%s", histotitle.c_str());
215 
216  // std::cout << "##### Booking: " << ibooker.pwd() << " title " <<
217  // histoname << std::endl;
218 
219  histopair.first.base =
220  ibooker.book1D(histoname.c_str(), histotitle.c_str(), i_residuals_Nbins, d_residual_xmin, d_residual_xmax);
221  histopair.first.base->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
222 
223  histopair.first.normed = ibooker.book1D(
224  normhistoname.c_str(), normhistotitle.c_str(), i_normres_Nbins, d_normres_xmin, d_normres_xmax);
225  histopair.first.normed->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
226  }
227  }
228  } // end loop over activeDets
229 }
230 
231 template <TrackerType pixel_or_strip>
233  auto vtracks = std::vector<TrackerValidationVariables::AVTrackStruct>();
234  // Filter out events if Trigger Filtering is requested
235  if (genTriggerEventFlag_->on() && !genTriggerEventFlag_->accept(iEvent, iSetup))
236  return;
237 
239  if (applyVertexCut_) {
240  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
241  if (!vertices.isValid() || vertices->empty())
242  return;
243  }
244 
245  // Retrieve tracker topology from geometry
246  const TrackerTopology *const tTopo = &iSetup.getData(trackerTopologyEventToken_);
247 
249  iEvent,
250  iSetup,
251  // tell the validator to only look at good tracks
252  [&](const reco::Track &track) -> bool {
253  return (!applyVertexCut_ ||
254  (track.pt() > 0.75 && abs(track.dxy(vertices->at(0).position())) < 5 * track.dxyError()));
255  },
256  vtracks);
257 
258  for (auto &track : vtracks) {
259  for (auto &it : track.hits) {
260  uint RawId = it.rawDetId;
261 
262  auto id = DetId(RawId);
263  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
264  if (isPixel != (pixel_or_strip == TRACKERTYPE_PIXEL))
265  continue;
266 
267  if (ModOn) {
268  auto &mod_histos = m_ModuleResiduals[std::make_pair("", RawId)];
269  mod_histos.x.base->Fill(it.resXprime);
270  mod_histos.x.normed->Fill(it.resXprime / it.resXprimeErr);
271  mod_histos.y.base->Fill(it.resYprime);
272  mod_histos.y.normed->Fill(it.resYprime / it.resYprimeErr);
273  }
274 
275  auto subdetandlayer = findSubdetAndLayer(RawId, tTopo);
276  auto histos = m_SubdetLayerResiduals[subdetandlayer];
277  // fill if its error is not zero
278  if (it.resXprimeErr != 0 && histos.x.base) {
279  histos.x.base->Fill(it.resXprime);
280  histos.x.normed->Fill(it.resXprime / it.resXprimeErr);
281  if (!isPixel)
282  tkhisto_ResidualsMean->fill(RawId, it.resXprime);
283  }
284  if (it.resYprimeErr != 0 && histos.y.base) {
285  histos.y.base->Fill(it.resYprime);
286  histos.y.normed->Fill(it.resYprime / it.resYprimeErr);
287  }
288  }
289  }
290 }
291 
T getUntrackedParameter(std::string const &, T const &) const
edm::ESGetToken< TkDetMap, TrackerTopologyRcd > tkDetMapToken_
void setSiStripFolderName(std::string name)
unsigned int tibLayer(const DetId &id) const
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyEventToken_
const DetIdContainer & detIds() const override
Returm a vector of all GeomDet DetIds (including those of GeomDetUnits)
double dxyError() const
error on dxy
Definition: TrackBase.h:769
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int pxfDisk(const DetId &id) const
TrackerValidationVariables avalidator_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyRunToken_
unsigned int tidWheel(const DetId &id) const
void setLayerFolder(uint32_t rawdetid, const TrackerTopology *tTopo, int32_t layer=0, bool ring_flag=false)
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
edm::EDGetTokenT< reco::VertexCollection > offlinePrimaryVerticesToken_
constexpr std::array< uint8_t, layerIndexSize > layer
void setDetectorFolder(uint32_t rawdetid, const TrackerTopology *tTopo)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
std::unique_ptr< TkHistoMap > tkhisto_ResidualsMean
int iEvent
Definition: GenABIO.cc:224
unsigned int tidSide(const DetId &id) const
vector< ParameterSet > Parameters
bool accept(const edm::Event &event, const edm::EventSetup &setup)
To be called from analyze/filter() methods.
double pt() const
track transverse momentum
Definition: TrackBase.h:637
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorTrackResidualsBase(const edm::ParameterSet &)
bool isValid() const
Definition: HandleBase.h:70
unsigned int pxbLayer(const DetId &id) const
SiStripHistoId hidmanager
Basic2DVector< T > xy() const
Definition: DetId.h:17
void fillTrackQuantities(const edm::Event &, const edm::EventSetup &, std::vector< AVTrackStruct > &v_avtrackout)
GenericTriggerEventFlag * genTriggerEventFlag_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
SiStripFolderOrganizer folder_organizer
void createMEs(DQMStore::IBooker &, const edm::EventSetup &)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeometryToken_
T get() const
Definition: EventSetup.h:82
bool isPixel(HitType hitType)
unsigned int pxfSide(const DetId &id) const
void initRun(const edm::Run &run, const edm::EventSetup &setup)
To be called from beginRun() methods.
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
void analyze(const edm::Event &, const edm::EventSetup &) override
std::pair< std::string, int32_t > findSubdetAndLayer(uint32_t ModuleID, const TrackerTopology *tTopo)
unsigned int tecWheel(const DetId &id) const
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
unsigned int tobLayer(const DetId &id) const
Definition: Run.h:45
unsigned int tecSide(const DetId &id) const
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)