CMS 3D CMS Logo

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