Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-07-02 22:49:44

0001 
0002 #include "TVector3.h"
0003 #include "TH2.h"
0004 #include "TLine.h"
0005 #include "TLatex.h"
0006 #include "TPaveText.h"
0007 #include "TCanvas.h"
0008 #include "TEveWindow.h"
0009 
0010 #include "DataFormats/TrackReco/interface/Track.h"
0011 #include "DataFormats/TrackReco/interface/HitPattern.h"
0012 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0013 
0014 #include "Fireworks/Core/interface/FWDetailView.h"
0015 #include "Fireworks/Core/interface/FWGeometry.h"
0016 #include "Fireworks/Core/interface/FWEventItem.h"
0017 #include "Fireworks/Core/interface/FWModelId.h"
0018 #include "Fireworks/Core/interface/fwLog.h"
0019 #include "Fireworks/Tracks/plugins/FWTrackResidualDetailView.h"
0020 
0021 using reco::HitPattern;
0022 using reco::Track;
0023 using reco::TrackBase;
0024 using reco::TrackResiduals;
0025 
0026 const char* FWTrackResidualDetailView::m_det_tracker_str[6] = {"PXB", "PXF", "TIB", "TID", "TOB", "TEC"};
0027 
0028 FWTrackResidualDetailView::FWTrackResidualDetailView()
0029     : m_ndet(0),
0030       m_nhits(0),
0031       m_resXFill(3007),
0032       m_resXCol(kGreen - 9),
0033       m_resYFill(3006),
0034       m_resYCol(kWhite),
0035       m_stereoFill(3004),
0036       m_stereoCol(kCyan - 9),
0037       m_invalidFill(3001),
0038       m_invalidCol(kRed) {}
0039 
0040 FWTrackResidualDetailView::~FWTrackResidualDetailView() {}
0041 
0042 void FWTrackResidualDetailView::prepareData(const FWModelId& id, const reco::Track* track) {
0043   auto const& residuals = track->residuals();
0044 
0045   const FWGeometry* geom = id.item()->getGeom();
0046   assert(geom != nullptr);
0047 
0048   const HitPattern& hitpat = track->hitPattern();
0049   m_nhits = hitpat.numberOfAllHits(reco::HitPattern::TRACK_HITS);
0050   hittype.reserve(m_nhits);
0051   stereo.reserve(m_nhits);
0052   subsubstruct.reserve(m_nhits);
0053   substruct.reserve(m_nhits);
0054   m_detector.reserve(m_nhits);
0055   res[0].reserve(m_nhits);
0056   res[1].reserve(m_nhits);
0057 
0058   for (int i = 0; i < m_nhits; ++i) {
0059     //printf("there are %d hits in the pattern, %d in the vector, this is %u\n",
0060     //        m_nhits, track->recHitsEnd() - track->recHitsBegin(), (*(track->recHitsBegin() + i))->geographicalId().rawId());
0061     uint32_t pattern = hitpat.getHitPattern(reco::HitPattern::TRACK_HITS, i);
0062     hittype.push_back(HitPattern::getHitType(pattern));
0063     stereo.push_back(HitPattern::getSide(pattern));
0064     subsubstruct.push_back(HitPattern::getSubSubStructure(pattern));
0065     substruct.push_back(HitPattern::getSubStructure(pattern));
0066     m_detector.push_back(HitPattern::getSubDetector(pattern));
0067     if ((*(track->recHitsBegin() + i))->isValid()) {
0068       res[0].push_back(
0069           getSignedResidual(geom, (*(track->recHitsBegin() + i))->geographicalId().rawId(), residuals.pullX(i)));
0070     } else {
0071       res[0].push_back(0);
0072     }
0073     res[1].push_back(residuals.pullY(i));
0074     // printf("%s, %i\n",m_det_tracker_str[substruct[i]-1],subsubstruct[i]);
0075   }
0076 
0077   m_det[0] = 0;
0078   for (int j = 0; j < m_nhits - 1;) {
0079     int k = j + 1;
0080     for (; k < m_nhits; k++) {
0081       if (substruct[j] == substruct[k] && subsubstruct[j] == subsubstruct[k]) {
0082         if (k == (m_nhits - 1))
0083           j = k;
0084       } else {
0085         m_ndet++;
0086         j = k;
0087         m_det[m_ndet] = j;
0088         break;
0089       }
0090     }
0091   }
0092   m_ndet++;
0093   m_det[m_ndet] = m_nhits;
0094   // printDebug();
0095 }
0096 
0097 void FWTrackResidualDetailView::build(const FWModelId& id, const reco::Track* track) {
0098   if (!track->extra().isAvailable()) {
0099     fwLog(fwlog::kError) << " no track extra information is available.\n";
0100     m_viewCanvas->cd();
0101     TLatex* latex = new TLatex();
0102     latex->SetTextSize(0.1);
0103     latex->DrawLatex(0.1, 0.5, "No track extra information");
0104     return;
0105   }
0106   prepareData(id, track);
0107 
0108   // draw histogram
0109   m_viewCanvas->cd();
0110   m_viewCanvas->SetHighLightColor(-1);
0111   TH2F* h_res = new TH2F("h_resx", "h_resx", 10, -5.5, 5.5, m_ndet, 0, m_ndet);
0112   TPad* padX = new TPad("pad1", "pad1", 0.2, 0., 0.8, 0.99);
0113   padX->SetBorderMode(0);
0114   padX->SetLeftMargin(0.2);
0115   padX->Draw();
0116   padX->cd();
0117   padX->SetFrameLineWidth(0);
0118   padX->Modified();
0119   h_res->SetDirectory(nullptr);
0120   h_res->SetStats(kFALSE);
0121   h_res->SetTitle("");
0122   h_res->SetXTitle("residual");
0123   h_res->GetXaxis()->SetTickLength(0);
0124   h_res->GetYaxis()->SetTickLength(0);
0125   h_res->GetXaxis()->SetNdivisions(20);
0126   h_res->GetYaxis()->SetLabelSize(0.06);
0127   h_res->Draw();
0128   padX->SetGridy();
0129 
0130   float larray[9] = {0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.5, 4.5, 5.5};
0131   float larray2[8];
0132   for (int l = 0; l < 8; l++) {
0133     float diff2 = (larray[l + 1] - larray[l]) / 2;
0134     larray2[l] = larray[l] + diff2;
0135     //  printf("(%.1f,%.1f),",larray[i],larray[i+1]);
0136   }
0137 
0138   int resi[2][64];
0139   for (int l = 0; l < m_nhits; l++) {
0140     for (int k = 0; k < 8; k++) {
0141       if (fabs(res[0][l]) == larray2[k])
0142         resi[0][l] = k;
0143       if (fabs(res[1][l]) == larray2[k])
0144         resi[1][l] = k;
0145     }
0146   }
0147 
0148   TLine* lines[17];
0149   for (int l = 0; l < 17; l++) {
0150     int ix = l % 9;
0151     int sign = 1;
0152     sign = (l > 8) ? -1 : 1;
0153     lines[l] = new TLine(sign * larray[ix], 0, sign * larray[ix], m_ndet);
0154     if (l != 9)
0155       lines[l]->SetLineStyle(3);
0156     padX->cd();
0157     lines[l]->Draw();
0158   }
0159 
0160   float width = 0.25;
0161   int filltype;
0162   Color_t color;
0163   Double_t box[4];
0164   padX->cd();
0165 
0166   for (int h = 0; h < 2; h++) {
0167     float height1 = 0;
0168     for (int j = 0; j < m_ndet; j++) {
0169       // take only X res and Y pixel residals
0170       if (strcmp(m_det_tracker_str[substruct[m_det[j]] - 1], "PXB") && h)
0171         continue;
0172 
0173       char det_str2[256];
0174       snprintf(det_str2, 255, "%s/%i", m_det_tracker_str[substruct[m_det[j]] - 1], subsubstruct[m_det[j]]);
0175       h_res->GetYaxis()->SetBinLabel(j + 1, det_str2);
0176 
0177       int diff = m_det[j + 1] - m_det[j];
0178       int k = 0;
0179       width = 1.0 / diff;
0180 
0181       for (int l = m_det[j]; l < (m_det[j] + diff); l++) {
0182         //      g->SetPoint(l,resx[l],j+0.5);
0183         //  printf("%i, %f %f %f\n",l,resx[l],sign*larray[resxi[l]],sign*larray[resxi[l]+1]);
0184         int sign = (res[h][l] < 0) ? -1 : 1;
0185         box[0] = (hittype[l] == 0) ? sign * larray[resi[h][l]] : -5.5;
0186         box[2] = (hittype[l] == 0) ? sign * larray[resi[h][l] + 1] : 5.5;
0187         box[1] = height1 + width * k;
0188         box[3] = height1 + width * (k + 1);
0189 
0190         if (stereo[l] == 1) {
0191           color = m_stereoCol;
0192           filltype = m_stereoFill;
0193         } else if (hittype[l] != 0) {
0194           color = m_invalidCol;
0195           filltype = m_invalidFill;
0196         } else {
0197           filltype = h ? m_resYFill : m_resXFill;
0198           color = h ? m_resYCol : m_resXCol;
0199         }
0200 
0201         drawCanvasBox(box, color, filltype, h < 1);
0202         k++;
0203       }
0204       height1 += 1;
0205     }
0206   }
0207 
0208   //  title
0209   const char* res_str = "residuals in Si detector local x-y coord.";
0210   TPaveText* pt = new TPaveText(0.0, 0.91, 1, 0.99, "blNDC");
0211   pt->SetBorderSize(0);
0212   pt->SetFillColor(kWhite);
0213   pt->AddText(res_str);
0214   pt->Draw();
0215 
0216   m_viewCanvas->cd();
0217   m_viewCanvas->SetEditable(kFALSE);
0218 
0219   setTextInfo(id, track);
0220 }
0221 
0222 double FWTrackResidualDetailView::getSignedResidual(const FWGeometry* geom, unsigned int id, double resX) {
0223   double local1[3] = {0, 0, 0};
0224   double local2[3] = {resX, 0, 0};
0225   double global1[3], global2[3];
0226   const TGeoMatrix* m = geom->getMatrix(id);
0227   assert(m != nullptr);
0228   m->LocalToMaster(local1, global1);
0229   m->LocalToMaster(local2, global2);
0230   TVector3 g1 = global1;
0231   TVector3 g2 = global2;
0232   if (g2.DeltaPhi(g1) > 0)
0233     return resX;
0234   else
0235     return -resX;
0236 }
0237 
0238 void FWTrackResidualDetailView::printDebug() {
0239   for (int i = 0; i < m_ndet; i++) {
0240     std::cout << "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<" << std::endl;
0241     std::cout << "idx " << i << " det[idx] " << m_det[i] << std::endl;
0242     std::cout << "m_det idx " << m_det[i] << std::endl;
0243     std::cout << "m_det_tracker_str idx " << substruct[m_det[i]] - 1 << std::endl;
0244     printf("m_det[idx] %i m_det_tracker_str %s substruct %i\n",
0245            m_det[i],
0246            m_det_tracker_str[substruct[m_det[i]] - 1],
0247            subsubstruct[m_det[i]]);
0248   }
0249 }
0250 
0251 void FWTrackResidualDetailView::setTextInfo(const FWModelId& /*id*/, const reco::Track* /*track*/) {
0252   m_infoCanvas->cd();
0253 
0254   char mytext[256];
0255   Double_t fontsize = 0.07;
0256   TLatex* latex = new TLatex();
0257   latex->SetTextSize(fontsize);
0258   latex->Draw();
0259   Double_t x0 = 0.02;
0260   Double_t y = 0.95;
0261 
0262   // summary
0263   int nvalid = 0;
0264   int npix = 0;
0265   int nstrip = 0;
0266   for (int i = 0; i < m_nhits; i++) {
0267     if (hittype[i] == 0)
0268       nvalid++;
0269     if (substruct[i] < 3)
0270       npix++;
0271     else
0272       nstrip++;
0273   }
0274 
0275   latex->SetTextSize(fontsize);
0276   Double_t boxH = 0.25 * fontsize;
0277 
0278   double yStep = 0.04;
0279 
0280   latex->DrawLatex(x0, y, "Residual:");
0281   y -= yStep;
0282   latex->DrawLatex(
0283       x0,
0284       y,
0285       "sgn(#hat{X}#bullet#hat{#phi}) #times #frac{X_{hit} - X_{traj}}{#sqrt{#sigma^{2}_{hit} + #sigma^{2}_{traj}}}");
0286   y -= 2.5 * yStep;
0287   snprintf(mytext, 255, "layers hit: %i", m_ndet);
0288   latex->DrawLatex(x0, y, mytext);
0289   y -= yStep;
0290   snprintf(mytext, 255, "valid Si hits: %i", nvalid);
0291   latex->DrawLatex(x0, y, mytext);
0292   y -= yStep;
0293   snprintf(mytext, 255, "total Si hits: %i", m_nhits);
0294   latex->DrawLatex(x0, y, mytext);
0295   y -= yStep;
0296   snprintf(mytext, 255, "valid Si pixel hits: %i", npix);
0297   latex->DrawLatex(x0, y, mytext);
0298   y -= yStep;
0299   snprintf(mytext, 255, "valid Si strip hits: %i", nstrip);
0300   latex->DrawLatex(x0, y, mytext);
0301 
0302   Double_t pos[4];
0303   pos[0] = 0.4;
0304   pos[2] = 0.55;
0305 
0306   y -= yStep * 2;
0307   latex->DrawLatex(x0, y, "X hit");
0308   pos[1] = y;
0309   pos[3] = pos[1] + boxH;
0310   drawCanvasBox(pos, m_resXCol, m_resXFill);
0311 
0312   y -= yStep;
0313   latex->DrawLatex(x0, y, "Y hit");
0314   pos[1] = y;
0315   pos[3] = pos[1] + boxH;
0316   drawCanvasBox(pos, m_resYCol, m_resYFill, false);
0317 
0318   y -= yStep;
0319   latex->DrawLatex(x0, y, "stereo hit");
0320   pos[1] = y;
0321   pos[3] = pos[1] + boxH;
0322   drawCanvasBox(pos, m_stereoCol, m_stereoFill);
0323 
0324   y -= yStep;
0325   latex->DrawLatex(x0, y, "invalid hit");
0326   pos[1] = y;
0327   pos[3] = pos[1] + boxH;
0328   drawCanvasBox(pos, m_invalidCol, m_invalidFill);
0329 }
0330 
0331 REGISTER_FWDETAILVIEW(FWTrackResidualDetailView, Residuals);