CMS 3D CMS Logo

MonitorTrackResiduals.cc
Go to the documentation of this file.
18 
19 template<TrackerType pixel_or_strip>
21  : conf_(iConfig), m_cacheID_(0)
22  , genTriggerEventFlag_(new GenericTriggerEventFlag(iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"), consumesCollector(), *this))
23  , avalidator_(iConfig, consumesCollector()) {
24  ModOn = conf_.getParameter<bool>("Mod_On");
25  offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(std::string("offlinePrimaryVertices"));
26 }
27 
28 template<TrackerType pixel_or_strip>
31 }
32 
33 template<TrackerType pixel_or_strip>
35 {
36  unsigned long long cacheID = iSetup.get<TrackerDigiGeometryRecord>().cacheIdentifier();
37  if (m_cacheID_ != cacheID) {
38  m_cacheID_ = cacheID;
39  this->createMEs( ibooker , iSetup);
40  }
41  std::string topFolderName_ = "SiStrip";
42  SiStripFolderOrganizer folder_organizer;
43  folder_organizer.setSiStripFolderName(topFolderName_);
44  tkhisto_ResidualsMean = std::unique_ptr<TkHistoMap>(new TkHistoMap(ibooker , topFolderName_ ,"TkHMap_ResidualsMean", 0.0,true));
45 }
46 
47 template<TrackerType pixel_or_strip>
49 
50  // Initialize the GenericTriggerEventFlag
51  if ( genTriggerEventFlag_->on() ) genTriggerEventFlag_->initRun( run, iSetup );
52 }
53 
54 template<TrackerType pixel_or_strip>
55 std::pair<std::string, int32_t> MonitorTrackResidualsBase<pixel_or_strip>::findSubdetAndLayer(uint32_t ModuleID, const TrackerTopology* tTopo) {
56  std::string subdet = "";
57  int32_t layer = 0;
58  auto id = DetId(ModuleID);
59  switch (id.subdetId()) {
60  // Pixel Barrel, Endcap
61  case 1:
62  subdet = "BPIX";
63  layer = tTopo->pxbLayer(id);
64  break;
65  case 2:
66  subdet = "FPIX";
67  layer = tTopo->pxfDisk(id) * ( tTopo->pxfSide(ModuleID)==1 ? -1 : +1);
68  break;
69  // Strip TIB, TID, TOB, TEC
70  case 3:
71  subdet = "TIB";
72  layer = tTopo->tibLayer(id);
73  break;
74  case 4:
75  subdet = "TID";
76  layer = tTopo->tidWheel(id) * ( tTopo->tidSide(ModuleID)==1 ? -1 : +1);
77  break;
78  case 5:
79  subdet = "TOB";
80  layer = tTopo->tobLayer(id);
81  break;
82  case 6:
83  subdet = "TEC";
84  layer = tTopo->tecWheel(id) * ( tTopo->tecSide(ModuleID)==1 ? -1 : +1);
85  break;
86  default:
87  // TODO: Fail loudly.
88  subdet = "UNKNOWN";
89  layer = 0;
90  }
91  return std::make_pair(subdet, layer);
92 }
93 
94 
95 template<TrackerType pixel_or_strip>
97 
98  //Retrieve tracker topology and geometry
100  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
101  const TrackerTopology* const tTopo = tTopoHandle.product();
102 
104  iSetup.get<TrackerDigiGeometryRecord>().get(TG);
105 
106  Parameters = conf_.getParameter<edm::ParameterSet>("TH1ResModules");
107  int32_t i_residuals_Nbins = Parameters.getParameter<int32_t>("Nbinx");
108  double d_residual_xmin = Parameters.getParameter<double>("xmin");
109  double d_residual_xmax = Parameters.getParameter<double>("xmax");
110  Parameters = conf_.getParameter<edm::ParameterSet>("TH1NormResModules");
111  int32_t i_normres_Nbins = Parameters.getParameter<int32_t>("Nbinx");
112  double d_normres_xmin = Parameters.getParameter<double>("xmin");
113  double d_normres_xmax = Parameters.getParameter<double>("xmax");
114 
115 
116  // use SistripHistoId for producing histogram id (and title)
117  SiStripHistoId hidmanager;
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  {
127  auto ModuleID = id.rawId();
128  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
129  if (isPixel != (pixel_or_strip == TRACKERTYPE_PIXEL)) continue;
130 
131  // Book module histogramms?
132  if (ModOn) {
133  switch (id.subdetId()) {
134  case 1: pixel_organizer.setModuleFolder(ibooker, ModuleID, 0); break;
135  case 2: pixel_organizer.setModuleFolder(ibooker, ModuleID, 0); break;
136  default: strip_organizer.setDetectorFolder(ModuleID,tTopo);
137  }
138  {
139  // this sounds strip specific but also works for pixel
140  std::string hid = hidmanager.createHistoId("HitResidualsX","det",ModuleID);
141  std::string normhid = hidmanager.createHistoId("NormalizedHitResidualsX","det",ModuleID);
142  auto& histos = m_ModuleResiduals[std::make_pair("", ModuleID)];
143  histos.x.base = ibooker.book1D(hid, hid, i_residuals_Nbins,d_residual_xmin,d_residual_xmax);
144  histos.x.base->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
145  histos.x.normed = ibooker.book1D(normhid, normhid, i_normres_Nbins,d_normres_xmin,d_normres_xmax);
146  histos.x.normed->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
147  }{
148  std::string hid = hidmanager.createHistoId("HitResidualsY","det",ModuleID);
149  std::string normhid = hidmanager.createHistoId("NormalizedHitResidualsY","det",ModuleID);
150  auto& histos = m_ModuleResiduals[std::make_pair("", ModuleID)];
151  histos.y.base = ibooker.book1D(hid, hid, i_residuals_Nbins,d_residual_xmin,d_residual_xmax);
152  histos.y.base->setAxisTitle("(y_{pred} - y_{rec})' [cm]");
153  histos.y.normed = ibooker.book1D(normhid, normhid, i_normres_Nbins,d_normres_xmin,d_normres_xmax);
154  histos.y.normed->setAxisTitle("(y_{pred} - y_{rec})'/#sigma");
155  }
156  }
157 
158  auto subdetandlayer = findSubdetAndLayer(ModuleID, tTopo);
159  if(m_SubdetLayerResiduals.find(subdetandlayer) == m_SubdetLayerResiduals.end()) {
160  // add new histograms
161  auto& histos = m_SubdetLayerResiduals[subdetandlayer];
162  switch (id.subdetId()) {
163  // Pixel Barrel, Endcap
164  // We can't use the folder organizer here (SiPixelActionExecutor.cc#1638 does the same)
165  case 1: ibooker.setCurrentFolder("Pixel/Barrel"); break;
166  case 2: ibooker.setCurrentFolder("Pixel/Endcap"); break;
167  // All strip
168  default: strip_organizer.setLayerFolder(ModuleID,tTopo,subdetandlayer.second);
169  }
170 
171  auto isBarrel = subdetandlayer.first.find("B") != std::string::npos;
172 
173  auto xy = std::vector<std::pair<HistoPair&, const char*> >
174  { std::make_pair(std::ref(histos.x), "X"),
175  std::make_pair(std::ref(histos.y), "Y") };
176  for (auto& histopair : xy) {
177  // book histogramms on layer level, check for barrel/pixel only for correct labeling
178 
179  // Skip the Y plots for strips.
180  if (!isPixel && histopair.second[0] == 'Y') continue;
181 
182  std::string histoname = isPixel ? ( // Pixel name
183  Form("HitResiduals%s_%s%d",
184  histopair.second,
185  isBarrel ? "L" : (subdetandlayer.second > 0 ? "Dp" : "Dm"),
186  std::abs(subdetandlayer.second)))
187  : (Form("HitResiduals_%s__%s__%d", // Strip TODO: We use a legacy name.
188  subdetandlayer.first.c_str(),
189  isBarrel ? "Layer" : "wheel",
190  std::abs(subdetandlayer.second)));
191 
192  std::string histotitle = Form("HitResiduals %s on %s%s full %s %d",
193  histopair.second,
194  subdetandlayer.first.c_str(),
195  isBarrel ? "" : (subdetandlayer.second > 0 ? "+" : "-"),
196  isBarrel ? "Layer" : (isPixel ? "Disk" : "Wheel"),
197  std::abs(subdetandlayer.second));
198 
199 
200  std::string normhistoname = Form("Normalized%s", histoname.c_str());
201  std::string normhistotitle = Form("Normalized%s", histotitle.c_str());
202 
203  //std::cout << "##### Booking: " << ibooker.pwd() << " title " << histoname << std::endl;
204 
205  histopair.first.base =
206  ibooker.book1D(histoname.c_str(),histotitle.c_str(),
207  i_residuals_Nbins,d_residual_xmin,d_residual_xmax);
208  histopair.first.base->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
209 
210  histopair.first.normed =
211  ibooker.book1D(normhistoname.c_str(),normhistotitle.c_str(),
212  i_normres_Nbins,d_normres_xmin,d_normres_xmax);
213  histopair.first.normed->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
214  }
215  }
216  } // end loop over activeDets
217 }
218 
219 template<TrackerType pixel_or_strip>
221 
222  auto vtracks = std::vector<TrackerValidationVariables::AVTrackStruct>();
223  // Filter out events if Trigger Filtering is requested
224  if (genTriggerEventFlag_->on()&& ! genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
225 
227  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
228  if (!vertices.isValid() || vertices->size() == 0) return;
229  const auto primaryVertex = vertices->at(0);
230 
231  //Retrieve tracker topology from geometry
232  edm::ESHandle<TrackerTopology> tTopoHandle;
233  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
234  const TrackerTopology* const tTopo = tTopoHandle.product();
235 
236  avalidator_.fillTrackQuantities(iEvent, iSetup,
237  // tell the validator to only look at good tracks
238  [&](const reco::Track& track) -> bool {
239  return track.pt() > 0.75
240  && abs( track.dxy(primaryVertex.position()) ) < 5*track.dxyError();
241  }, vtracks);
242 
243  for (auto& track : vtracks) {
244  for (auto& it : track.hits) {
245  uint RawId = it.rawDetId;
246 
247  auto id = DetId(RawId);
248  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
249  if (isPixel != (pixel_or_strip == TRACKERTYPE_PIXEL)) continue;
250 
251 
252  if (ModOn) {
253  auto& mod_histos = m_ModuleResiduals[std::make_pair("",RawId)];
254  mod_histos.x.base->Fill(it.resXprime);
255  mod_histos.x.normed->Fill(it.resXprime/it.resXprimeErr);
256  mod_histos.y.base->Fill(it.resYprime);
257  mod_histos.y.normed->Fill(it.resYprime/it.resYprimeErr);
258  }
259 
260  auto subdetandlayer = findSubdetAndLayer(RawId, tTopo);
261  auto histos = m_SubdetLayerResiduals[subdetandlayer];
262  // fill if its error is not zero
263  if(it.resXprimeErr != 0 && histos.x.base) {
264  histos.x.base->Fill(it.resXprime);
265  histos.x.normed->Fill(it.resXprime/it.resXprimeErr);
266  if(!isPixel) tkhisto_ResidualsMean->fill(RawId,it.resXprime);
267  }
268  if(it.resYprimeErr != 0 && histos.y.base) {
269  histos.y.base->Fill(it.resYprime);
270  histos.y.normed->Fill(it.resYprime/it.resYprimeErr);
271  }
272  }
273  }
274 }
275 
276 
277 
280 
T getParameter(std::string const &) const
void setSiStripFolderName(std::string name)
unsigned int tibLayer(const DetId &id) const
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
double dxyError() const
error on dxy
Definition: TrackBase.h:796
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
unsigned int pxfDisk(const DetId &id) const
TrackerValidationVariables avalidator_
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_
void setDetectorFolder(uint32_t rawdetid, const TrackerTopology *tTopo)
std::unique_ptr< TkHistoMap > tkhisto_ResidualsMean
int iEvent
Definition: GenABIO.cc:230
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:621
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
const DetIdContainer & detIds() const override
Returm a vector of all GeomDet DetIds (including those of GeomDetUnits)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorTrackResidualsBase(const edm::ParameterSet &)
bool isValid() const
Definition: HandleBase.h:74
unsigned int pxbLayer(const DetId &id) const
Definition: DetId.h:18
void fillTrackQuantities(const edm::Event &, const edm::EventSetup &, std::vector< AVTrackStruct > &v_avtrackout)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
GenericTriggerEventFlag * genTriggerEventFlag_
const T & get() const
Definition: EventSetup.h:55
void createMEs(DQMStore::IBooker &, const edm::EventSetup &)
HLT enums.
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.
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:591
virtual 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 setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
T const * product() const
Definition: ESHandle.h:86
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
unsigned int tobLayer(const DetId &id) const
Definition: Run.h:42
unsigned int tecSide(const DetId &id) const