CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FWTrackResidualDetailView.cc
Go to the documentation of this file.
1 
2 #include "TVector3.h"
3 #include <TH2.h>
4 #include <TLine.h>
5 #include <TLatex.h>
6 #include <TPaveText.h>
7 #include <TCanvas.h>
8 #include <TEveWindow.h>
9 
13 
20 
21 using reco::Track;
22 using reco::TrackBase;
23 using reco::HitPattern;
25 
26 const char* FWTrackResidualDetailView::m_det_tracker_str[6]={"PXB","PXF","TIB","TID","TOB","TEC"};
27 
29  m_ndet(0),
30  m_nhits(0),
31  m_resXFill(3007),
32  m_resXCol(kGreen-9),
33  m_resYFill(3006),
34  m_resYCol(kWhite),
35  m_stereoFill(3004),
36  m_stereoCol(kCyan-9),
37  m_invalidFill(3001),
38  m_invalidCol(kRed)
39 {
40  memset(m_det, 0, sizeof(m_det));
41  memset(res, 0, sizeof(res));
42  memset(hittype, 0, sizeof(hittype));
43  memset(stereo, 0, sizeof(stereo));
44  memset(substruct, 0, sizeof(substruct));
45  memset(subsubstruct, 0, sizeof(subsubstruct));
46  memset(m_detector, 0, sizeof(m_detector));
47 }
48 
50 {
51 }
52 
53 void
55 {
56  HitPattern hitpat = track->hitPattern();
57  TrackResiduals residuals = track->residuals();
58 
59  const FWGeometry *geom = id.item()->getGeom();
60  assert(geom != 0);
61  m_nhits=hitpat.numberOfHits();
62  for (int i = 0; i < m_nhits; ++i) {
63  // printf("there are %d hits in the pattern, %d in the vector, this is %u\n",
64  // m_nhits, track->recHitsEnd() - track->recHitsBegin(), (*(track->recHitsBegin() + i))->geographicalId().rawId());
65  hittype[i] = 0x3 & hitpat.getHitPattern(i);
66  stereo[i] = 0x1 & hitpat.getHitPattern(i) >> 2;
67  subsubstruct[i] = 0xf & hitpat.getHitPattern(i) >> 3;
68  substruct[i] = 0x7 & hitpat.getHitPattern(i) >> 7;
69  m_detector[i] = 0x01 & hitpat.getHitPattern(i) >> 10;
70  if ((*(track->recHitsBegin() + i))->isValid()) {
71  res[0][i] = getSignedResidual(geom,
72  (*(track->recHitsBegin() + i))->geographicalId().rawId(),
73  residuals.residualX(i, hitpat));
74  } else {
75  res[0][i] = 0;
76  }
77  res[1][i] = residuals.residualY(i, hitpat);
78  // printf("%s, %i\n",m_det_tracker_str[substruct[i]-1],subsubstruct[i]);
79  }
80 
81  m_det[0]=0;
82  for(int j=0; j < m_nhits-1;) {
83  int k=j+1;
84  for(; k<m_nhits ; k++) {
85  if(substruct[j]==substruct[k] && subsubstruct[j]==subsubstruct[k]) {
86  if(k==(m_nhits-1)) j=k;
87  }
88  else {
89  m_ndet++;
90  j=k;
91  m_det[m_ndet]=j;
92  break;
93  }
94  }
95  }
96  m_ndet++;
98  // printDebug();
99 }
100 
101 void
103 {
104  if (!track->extra().isAvailable()) {
105  fwLog(fwlog::kError) << " no track extra information is available.\n";
106  m_viewCanvas->cd();
107  TLatex* latex = new TLatex();
108  latex->SetTextSize(0.1);
109  latex->DrawLatex(0.1, 0.5, "No track extra information");
110  return;
111  }
112  prepareData(id, track);
113 
114  // draw histogram
115  m_viewCanvas->cd();
116  m_viewCanvas->SetHighLightColor(-1);
117  TH2F* h_res = new TH2F("h_resx","h_resx",10,-5.5,5.5,m_ndet,0,m_ndet);
118  TPad* padX = new TPad("pad1","pad1", 0.2, 0., 0.8, 0.99);
119  padX->SetBorderMode(0);
120  padX->SetLeftMargin(0.2);
121  padX->Draw();
122  padX->cd();
123  padX->SetFrameLineWidth(0);
124  padX->Modified();
125  h_res->SetDirectory(0);
126  h_res->SetStats(kFALSE);
127  h_res->SetTitle("");
128  h_res->SetXTitle("residual");
129  h_res->GetXaxis()->SetTickLength(0);
130  h_res->GetYaxis()->SetTickLength(0);
131  h_res->GetXaxis()->SetNdivisions(20);
132  h_res->GetYaxis()->SetLabelSize(0.06);
133  h_res->Draw();
134  padX->SetGridy();
135 
136  float larray[9]={0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.5, 4.5, 5.5};
137  float larray2[8];
138  for(int l=0; l<8; l++) {
139  float diff2=(larray[l+1]-larray[l])/2;
140  larray2[l]=larray[l]+diff2;
141  // printf("(%.1f,%.1f),",larray[i],larray[i+1]);
142  }
143 
144  int resi[2][64];
145  for(int l=0; l<m_nhits; l++) {
146  for(int k=0; k<8; k++) {
147  if(fabs(res[0][l])==larray2[k])
148  resi[0][l]=k;
149  if(fabs(res[1][l])==larray2[k])
150  resi[1][l]=k;
151  }
152  }
153 
154  TLine* lines[17];
155  for(int l=0; l<17; l++) {
156  int ix=l%9;
157  int sign=1;
158  sign = (l>8) ? -1 : 1;
159  lines[l] = new TLine(sign*larray[ix],0,sign*larray[ix],m_ndet);
160  if(l!=9)
161  lines[l]->SetLineStyle(3);
162  padX->cd();
163  lines[l]->Draw();
164  }
165 
166  float width=0.25;
167  int filltype;
168  Color_t color;
169  Double_t box[4];
170  padX->cd();
171 
172  for(int h=0; h<2; h++) {
173  float height1=0;
174  for(int j=0; j<m_ndet; j++) {
175  // take only X res and Y pixel residals
176  if (strcmp(m_det_tracker_str[substruct[m_det[j]]-1], "PXB") && h)
177  continue;
178 
179  char det_str2[256];
180  snprintf(det_str2, 255, "%s/%i",m_det_tracker_str[substruct[m_det[j]]-1],subsubstruct[m_det[j]]);
181  h_res->GetYaxis()->SetBinLabel(j+1, det_str2);
182 
183  int diff=m_det[j+1]-m_det[j];
184  int k=0;
185  width=1.0/diff;
186 
187  for(int l=m_det[j]; l<(m_det[j]+diff); l++) {
188  // g->SetPoint(l,resx[l],j+0.5);
189  // printf("%i, %f %f %f\n",l,resx[l],sign*larray[resxi[l]],sign*larray[resxi[l]+1]);
190  int sign = (res[h][l]<0) ? -1 : 1;
191  box[0] = (hittype[l]==0) ? sign*larray[resi[h][l]] : -5.5;
192  box[2] = (hittype[l]==0) ? sign*larray[resi[h][l]+1] : 5.5;
193  box[1] = height1+width*k;
194  box[3] = height1+width*(k+1);
195 
196  if(stereo[l]==1) {
197  color = m_stereoCol;
198  filltype = m_stereoFill;
199  }
200  else if(hittype[l]!=0) {
201  color = m_invalidCol;
202  filltype = m_invalidFill;
203  }
204  else {
205  filltype = h ? m_resYFill : m_resXFill;
206  color = h ? m_resYCol : m_resXCol;
207  }
208 
209  drawCanvasBox(box, color, filltype, h<1);
210  k++;
211  }
212  height1 +=1;
213  }
214  }
215 
216  // title
217  const char* res_str= "residuals in Si detector local x-y coord.";
218  TPaveText *pt = new TPaveText(0.0,0.91, 1,0.99,"blNDC");
219  pt->SetBorderSize(0);
220  pt->SetFillColor(kWhite);
221  pt->AddText(res_str);
222  pt->Draw();
223 
224  m_viewCanvas->cd();
225  m_viewCanvas->SetEditable(kFALSE);
226 
227  setTextInfo(id, track);
228 }
229 
230 double
231 FWTrackResidualDetailView::getSignedResidual (const FWGeometry *geom, unsigned int id, double resX)
232 {
233  double local1[3] = { 0, 0, 0 };
234  double local2[3] = { resX, 0, 0 };
235  double global1[3], global2[3];
236  const TGeoMatrix *m = geom->getMatrix(id);
237  assert(m != 0);
238  m->LocalToMaster(local1, global1);
239  m->LocalToMaster(local2, global2);
240  TVector3 g1 = global1;
241  TVector3 g2 = global2;
242  if (g2.DeltaPhi(g1) > 0)
243  return resX;
244  else return -resX;
245 }
246 
247 void
249 {
250  for(int i=0; i<m_ndet; i++)
251  {
252  std::cout <<"<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<" << std::endl;
253  std::cout << "idx " << i << " det[idx] " << m_det[i] << std::endl;
254  std::cout << "m_det idx " << m_det[i] << std::endl;
255  std::cout << "m_det_tracker_str idx " << substruct[m_det[i]]-1 << std::endl;
256  printf("m_det[idx] %i m_det_tracker_str %s substruct %i\n",m_det[i], m_det_tracker_str[substruct[m_det[i]]-1], subsubstruct[m_det[i]]);
257  }
258 }
259 
260 void
262 {
263  m_infoCanvas->cd();
264 
265  char mytext[256];
266  Double_t fontsize = 0.07;
267  TLatex* latex = new TLatex();
268  latex->SetTextSize(fontsize);
269  latex->Draw();
270  Double_t x0 = 0.02;
271  Double_t y = 0.95;
272 
273  // summary
274  int nvalid=0;
275  int npix=0;
276  int nstrip=0;
277  for(int i=0; i<m_nhits; i++)
278  {
279  if(hittype[i]==0) nvalid++;
280  if(substruct[i]<3) npix++;
281  else nstrip++;
282  }
283 
284  latex->SetTextSize(fontsize);
285  Double_t boxH = 0.25*fontsize;
286 
287  double yStep = 0.04;
288 
289  latex->DrawLatex(x0, y, "Residual:");
290  y-= yStep;
291  latex->DrawLatex(x0, y, "sgn(#hat{X}#bullet#hat{#phi}) #times #frac{X_{hit} - X_{traj}}{#sqrt{#sigma^{2}_{hit} + #sigma^{2}_{traj}}}" );
292  y-= 2.5*yStep;
293  snprintf(mytext, 255, "layers hit: %i", m_ndet);
294  latex->DrawLatex(x0, y, mytext);
295  y -= yStep;
296  snprintf(mytext, 255,"valid Si hits: %i", nvalid);
297  latex->DrawLatex(x0, y, mytext);
298  y -= yStep;
299  snprintf(mytext, 255,"total Si hits: %i", m_nhits);
300  latex->DrawLatex(x0, y, mytext);
301  y -= yStep;
302  snprintf(mytext, 255,"valid Si pixel hits: %i", npix);
303  latex->DrawLatex(x0, y, mytext);
304  y -= yStep;
305  snprintf(mytext, 255, "valid Si strip hits: %i", nstrip);
306  latex->DrawLatex(x0, y, mytext);
307 
308 
309  Double_t pos[4];
310  pos[0] = 0.4;
311  pos[2] = 0.55;
312 
313  y -= yStep*2;
314  latex->DrawLatex(x0, y, "X hit");
315  pos[1] = y; pos[3] = pos[1] + boxH;
317 
318  y -= yStep;
319  latex->DrawLatex(x0, y, "Y hit");
320  pos[1] = y; pos[3] = pos[1] + boxH;
322 
323  y -= yStep;
324  latex->DrawLatex(x0, y, "stereo hit");
325  pos[1] = y; pos[3] = pos[1] + boxH;
327 
328  y -= yStep;
329  latex->DrawLatex(x0, y, "invalid hit");
330  pos[1] = y; pos[3] = pos[1] + boxH;
332 }
333 
int i
Definition: DBlmapReader.cc:9
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:97
double getSignedResidual(const FWGeometry *geom, unsigned int id, double resX)
double residualY(int i, const HitPattern &) const
bool isAvailable() const
Definition: Ref.h:276
const TGeoMatrix * getMatrix(unsigned int id) const
Definition: FWGeometry.cc:166
int numberOfHits() const
Definition: HitPattern.cc:213
int j
Definition: DBlmapReader.cc:9
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
int k[5][pyjets_maxn]
double residualX(int i, const HitPattern &) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
#define fwLog(_level_)
Definition: fwLog.h:51
void prepareData(const FWModelId &id, const reco::Track *)
virtual void build(const FWModelId &id, const reco::Track *)
#define REGISTER_FWDETAILVIEW(_classname_, _name_)
const TrackResiduals & residuals() const
Definition: Track.h:117
tuple cout
Definition: gather_cfg.py:121
static const char * m_det_tracker_str[]
uint32_t getHitPattern(int position) const
Definition: HitPattern.cc:144
static void drawCanvasBox(Double_t *pos, Color_t fillCol, Int_t fillType=0, bool bg=kTRUE)
virtual void setTextInfo(const FWModelId &id, const reco::Track *)