CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/Fireworks/Tracks/plugins/FWTrackResidualDetailView.cc

Go to the documentation of this file.
00001 
00002 #include "TVector3.h"
00003 #include <TH2.h>
00004 #include <TLine.h>
00005 #include <TLatex.h>
00006 #include <TPaveText.h>
00007 #include <TCanvas.h>
00008 #include <TEveWindow.h>
00009 
00010 #include "DataFormats/TrackReco/interface/Track.h"
00011 #include "DataFormats/TrackReco/interface/HitPattern.h"
00012 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00013 
00014 #include "Fireworks/Core/interface/FWDetailView.h"
00015 #include "Fireworks/Core/interface/FWGeometry.h"
00016 #include "Fireworks/Core/interface/FWEventItem.h"
00017 #include "Fireworks/Core/interface/FWModelId.h"
00018 #include "Fireworks/Core/interface/fwLog.h"
00019 #include "Fireworks/Tracks/plugins/FWTrackResidualDetailView.h"
00020 
00021 using reco::Track;
00022 using reco::TrackBase;
00023 using reco::HitPattern;
00024 using reco::TrackResiduals;
00025 
00026 const char* FWTrackResidualDetailView::m_det_tracker_str[6]={"PXB","PXF","TIB","TID","TOB","TEC"};
00027 
00028 FWTrackResidualDetailView::FWTrackResidualDetailView ():
00029    m_ndet(0),
00030    m_nhits(0),
00031    m_resXFill(3007),
00032    m_resXCol(kGreen-9),
00033    m_resYFill(3006),
00034    m_resYCol(kWhite),
00035    m_stereoFill(3004),
00036    m_stereoCol(kCyan-9),
00037    m_invalidFill(3001),
00038    m_invalidCol(kRed)
00039 {
00040    memset(m_det,        0, sizeof(m_det));
00041    memset(res,          0, sizeof(res));
00042    memset(hittype,      0, sizeof(hittype));
00043    memset(stereo,       0, sizeof(stereo));
00044    memset(substruct,    0, sizeof(substruct));
00045    memset(subsubstruct, 0, sizeof(subsubstruct));
00046    memset(m_detector,   0, sizeof(m_detector));
00047 }
00048 
00049 FWTrackResidualDetailView::~FWTrackResidualDetailView ()
00050 {
00051 }
00052 
00053 void
00054 FWTrackResidualDetailView::prepareData(const FWModelId &id, const reco::Track* track)
00055 {
00056    HitPattern hitpat = track->hitPattern();
00057    TrackResiduals residuals = track->residuals();
00058 
00059    const FWGeometry *geom = id.item()->getGeom();
00060    assert(geom != 0);
00061    m_nhits=hitpat.numberOfHits();
00062    for (int i = 0; i < m_nhits; ++i) {
00063       //        printf("there are %d hits in the pattern, %d in the vector, this is %u\n",
00064       //               m_nhits, track->recHitsEnd() - track->recHitsBegin(), (*(track->recHitsBegin() + i))->geographicalId().rawId());
00065       hittype[i] = 0x3 & hitpat.getHitPattern(i);
00066       stereo[i] = 0x1 & hitpat.getHitPattern(i) >> 2;
00067       subsubstruct[i] = 0xf & hitpat.getHitPattern(i) >> 3;
00068       substruct[i] = 0x7 & hitpat.getHitPattern(i) >> 7;
00069       m_detector[i] = 0x01 & hitpat.getHitPattern(i) >> 10;
00070       if ((*(track->recHitsBegin() + i))->isValid()) {
00071          res[0][i] = getSignedResidual(geom,
00072                                        (*(track->recHitsBegin() + i))->geographicalId().rawId(),
00073                                        residuals.residualX(i, hitpat));
00074       } else {
00075          res[0][i] = 0;
00076       }
00077       res[1][i] = residuals.residualY(i, hitpat);
00078       // printf("%s, %i\n",m_det_tracker_str[substruct[i]-1],subsubstruct[i]);
00079    }
00080 
00081    m_det[0]=0;
00082    for(int j=0; j < m_nhits-1;) {
00083       int k=j+1;
00084       for(; k<m_nhits ; k++) {
00085          if(substruct[j]==substruct[k]  && subsubstruct[j]==subsubstruct[k]) {
00086             if(k==(m_nhits-1)) j=k;
00087          }
00088          else {
00089             m_ndet++;
00090             j=k;
00091             m_det[m_ndet]=j;
00092             break;
00093          }
00094       }
00095    }
00096    m_ndet++;
00097    m_det[m_ndet]=m_nhits;
00098    // printDebug();
00099 }
00100 
00101 void
00102 FWTrackResidualDetailView::build (const FWModelId &id, const reco::Track* track)
00103 {
00104    if (!track->extra().isAvailable()) {
00105       fwLog(fwlog::kError) << " no track extra information is available.\n";
00106      m_viewCanvas->cd();
00107      TLatex* latex = new TLatex();
00108      latex->SetTextSize(0.1);
00109      latex->DrawLatex(0.1, 0.5, "No track extra information");
00110      return;
00111    }
00112    prepareData(id, track);
00113 
00114    // draw histogram
00115    m_viewCanvas->cd();
00116    m_viewCanvas->SetHighLightColor(-1);
00117    TH2F* h_res = new TH2F("h_resx","h_resx",10,-5.5,5.5,m_ndet,0,m_ndet);
00118    TPad* padX = new TPad("pad1","pad1", 0.2, 0., 0.8, 0.99);
00119    padX->SetBorderMode(0);
00120    padX->SetLeftMargin(0.2);
00121    padX->Draw();
00122    padX->cd();
00123    padX->SetFrameLineWidth(0);
00124    padX->Modified();
00125    h_res->SetDirectory(0);
00126    h_res->SetStats(kFALSE);
00127    h_res->SetTitle("");
00128    h_res->SetXTitle("residual");
00129    h_res->GetXaxis()->SetTickLength(0);
00130    h_res->GetYaxis()->SetTickLength(0);
00131    h_res->GetXaxis()->SetNdivisions(20);
00132    h_res->GetYaxis()->SetLabelSize(0.06);
00133    h_res->Draw();
00134    padX->SetGridy();
00135 
00136    float larray[9]={0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.5, 4.5, 5.5};
00137    float larray2[8];
00138    for(int l=0; l<8; l++) {
00139       float diff2=(larray[l+1]-larray[l])/2;
00140       larray2[l]=larray[l]+diff2;
00141       //  printf("(%.1f,%.1f),",larray[i],larray[i+1]);
00142    }
00143 
00144    int resi[2][64];
00145    for(int l=0; l<m_nhits; l++) {
00146       for(int k=0; k<8; k++) {
00147          if(fabs(res[0][l])==larray2[k])
00148             resi[0][l]=k;
00149          if(fabs(res[1][l])==larray2[k])
00150             resi[1][l]=k;
00151       }
00152    }
00153 
00154    TLine* lines[17];
00155    for(int l=0; l<17; l++) {
00156       int ix=l%9;
00157       int sign=1;
00158       sign = (l>8) ? -1 : 1;
00159       lines[l] = new TLine(sign*larray[ix],0,sign*larray[ix],m_ndet);
00160       if(l!=9)
00161          lines[l]->SetLineStyle(3);
00162       padX->cd();
00163       lines[l]->Draw();
00164    }
00165 
00166    float width=0.25;
00167    int filltype;
00168    Color_t color;
00169    Double_t box[4];
00170    padX->cd();
00171 
00172    for(int h=0; h<2; h++) {
00173       float height1=0;
00174       for(int j=0; j<m_ndet; j++) {
00175          // take only X res and Y pixel residals
00176          if (strcmp(m_det_tracker_str[substruct[m_det[j]]-1], "PXB") && h)
00177             continue;
00178 
00179          char det_str2[256];
00180          snprintf(det_str2, 255, "%s/%i",m_det_tracker_str[substruct[m_det[j]]-1],subsubstruct[m_det[j]]);
00181          h_res->GetYaxis()->SetBinLabel(j+1, det_str2);
00182 
00183          int diff=m_det[j+1]-m_det[j];
00184          int k=0;
00185          width=1.0/diff;
00186 
00187          for(int l=m_det[j]; l<(m_det[j]+diff); l++) {
00188             //      g->SetPoint(l,resx[l],j+0.5);
00189             //  printf("%i, %f %f %f\n",l,resx[l],sign*larray[resxi[l]],sign*larray[resxi[l]+1]);
00190             int sign = (res[h][l]<0) ? -1 : 1;
00191             box[0] = (hittype[l]==0) ? sign*larray[resi[h][l]] : -5.5;
00192             box[2] = (hittype[l]==0) ? sign*larray[resi[h][l]+1] : 5.5;
00193             box[1] = height1+width*k;
00194             box[3] = height1+width*(k+1);
00195 
00196             if(stereo[l]==1) {
00197                color    = m_stereoCol;
00198                filltype = m_stereoFill;
00199             }
00200             else if(hittype[l]!=0) {
00201                color    = m_invalidCol;
00202                filltype = m_invalidFill;
00203             }
00204             else {
00205                filltype = h ? m_resYFill : m_resXFill;
00206                color    = h ? m_resYCol  : m_resXCol;
00207             }
00208 
00209             drawCanvasBox(box, color, filltype, h<1);
00210             k++;
00211          }
00212          height1 +=1;
00213       }
00214    }
00215    
00216    //  title
00217    const char* res_str= "residuals in Si detector local x-y coord.";
00218    TPaveText *pt = new TPaveText(0.0,0.91, 1,0.99,"blNDC");
00219    pt->SetBorderSize(0);
00220    pt->SetFillColor(kWhite);
00221    pt->AddText(res_str);
00222    pt->Draw();
00223    
00224    m_viewCanvas->cd();
00225    m_viewCanvas->SetEditable(kFALSE);
00226 
00227    setTextInfo(id, track);
00228 }
00229 
00230 double
00231 FWTrackResidualDetailView::getSignedResidual (const FWGeometry *geom, unsigned int id, double resX)
00232 {
00233    double local1[3] = { 0, 0, 0 };
00234    double local2[3] = { resX, 0, 0 };
00235    double global1[3], global2[3];
00236    const TGeoMatrix *m = geom->getMatrix(id);
00237    assert(m != 0);
00238    m->LocalToMaster(local1, global1);
00239    m->LocalToMaster(local2, global2);
00240    TVector3 g1 = global1;
00241    TVector3 g2 = global2;
00242    if (g2.DeltaPhi(g1) > 0)
00243       return resX;
00244    else return -resX;
00245 }
00246 
00247 void
00248 FWTrackResidualDetailView::printDebug()
00249 {
00250    for(int i=0; i<m_ndet; i++)
00251    {
00252       std::cout <<"<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<" << std::endl;
00253       std::cout << "idx " << i << " det[idx] " <<  m_det[i] << std::endl;
00254       std::cout << "m_det idx " << m_det[i] <<   std::endl;
00255       std::cout << "m_det_tracker_str idx " << substruct[m_det[i]]-1  << std::endl;
00256       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]]);
00257    }
00258 }
00259 
00260 void
00261 FWTrackResidualDetailView::setTextInfo(const FWModelId &/*id*/, const reco::Track* /*track*/)
00262 {
00263    m_infoCanvas->cd();
00264 
00265    char mytext[256];
00266    Double_t fontsize = 0.07;
00267    TLatex* latex = new TLatex();
00268    latex->SetTextSize(fontsize);
00269    latex->Draw();
00270    Double_t x0 = 0.02;
00271    Double_t y = 0.95;
00272 
00273    // summary
00274    int nvalid=0;
00275    int npix=0;
00276    int nstrip=0;
00277    for(int i=0; i<m_nhits; i++)
00278    {
00279       if(hittype[i]==0) nvalid++;
00280       if(substruct[i]<3) npix++;
00281       else nstrip++;
00282    }
00283 
00284    latex->SetTextSize(fontsize);
00285    Double_t boxH = 0.25*fontsize;
00286    
00287    double yStep = 0.04;
00288 
00289    latex->DrawLatex(x0, y, "Residual:");
00290    y-= yStep;
00291    latex->DrawLatex(x0, y, "sgn(#hat{X}#bullet#hat{#phi}) #times #frac{X_{hit} - X_{traj}}{#sqrt{#sigma^{2}_{hit} + #sigma^{2}_{traj}}}" );
00292    y-= 2.5*yStep;
00293    snprintf(mytext, 255, "layers hit: %i", m_ndet);
00294    latex->DrawLatex(x0, y, mytext);
00295    y -= yStep;
00296    snprintf(mytext, 255,"valid Si hits: %i", nvalid);
00297    latex->DrawLatex(x0, y, mytext);
00298    y -= yStep;
00299    snprintf(mytext, 255,"total Si hits: %i", m_nhits);
00300    latex->DrawLatex(x0, y, mytext);
00301    y -= yStep;
00302    snprintf(mytext, 255,"valid Si pixel hits: %i", npix);
00303    latex->DrawLatex(x0, y, mytext);
00304    y -= yStep;
00305    snprintf(mytext, 255, "valid Si strip hits: %i", nstrip);
00306    latex->DrawLatex(x0, y, mytext);
00307    
00308 
00309    Double_t pos[4];
00310    pos[0] = 0.4;
00311    pos[2] = 0.55;
00312 
00313    y -= yStep*2;
00314    latex->DrawLatex(x0, y, "X hit");
00315    pos[1] = y; pos[3] = pos[1] + boxH;
00316    drawCanvasBox(pos, m_resXCol, m_resXFill);
00317 
00318    y -=  yStep;
00319    latex->DrawLatex(x0, y, "Y hit");
00320    pos[1] = y; pos[3] = pos[1] + boxH;
00321    drawCanvasBox(pos, m_resYCol, m_resYFill, 0);
00322 
00323    y -= yStep;
00324    latex->DrawLatex(x0, y, "stereo hit");
00325    pos[1] = y; pos[3] = pos[1] + boxH;
00326    drawCanvasBox(pos, m_stereoCol, m_stereoFill);
00327 
00328    y -= yStep;
00329    latex->DrawLatex(x0, y, "invalid hit");
00330    pos[1] = y; pos[3] = pos[1] + boxH;
00331    drawCanvasBox(pos, m_invalidCol, m_invalidFill);
00332 }
00333 
00334 REGISTER_FWDETAILVIEW(FWTrackResidualDetailView, Residuals);