TMarker3DBox.cxx

Go to the documentation of this file.
00001 // @(#)root/g3d:$Id: TMarker3DBox.cxx 20882 2007-11-19 11:31:26Z rdm $
00002 // Author: Rene Brun , Olivier Couet 31/10/97
00003 
00004 
00005 /*************************************************************************
00006  * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers.               *
00007  * All rights reserved.                                                  *
00008  *                                                                       *
00009  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00010  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00011  *************************************************************************/
00012 
00013 #include "Riostream.h"
00014 #include "TROOT.h"
00015 #include "TView.h"
00016 #include "TMarker3DBox.h"
00017 #include "TVirtualPad.h"
00018 #include "TH1.h"
00019 #include "TH3.h"
00020 #include "TBuffer3D.h"
00021 #include "TBuffer3DTypes.h"
00022 #include "TVirtualViewer3D.h"
00023 #include "TGeometry.h"
00024 #include "TClass.h"
00025 #include "TMath.h"
00026 
00027 #include <assert.h>
00028 
00029 ClassImp(TMarker3DBox)
00030 
00031 
00032 //______________________________________________________________________________
00033 // Marker3DBox is a special 3-D marker designed for event display.
00034 // It has the following parameters:
00035 //    fX;               X coordinate of the center of the box
00036 //    fY;               Y coordinate of the center of the box
00037 //    fZ;               Z coordinate of the center of the box
00038 //    fDx;              half length in X
00039 //    fDy;              half length in Y
00040 //    fDz;              half length in Z
00041 //    fTheta;           Angle of box z axis with respect to main Z axis
00042 //    fPhi;             Angle of box x axis with respect to main Xaxis
00043 //    fRefObject;       A reference to an object
00044 // Begin_Html <P ALIGN=CENTER> <IMG SRC="gif/Marker3DBox.gif"> </P> End_Html
00045 
00046 
00047 //______________________________________________________________________________
00048 TMarker3DBox::TMarker3DBox()
00049 {
00050    // Marker3DBox  default constructor
00051 
00052    fRefObject = 0;
00053    fDx        = 1;
00054    fDy        = 1;
00055    fDz        = 1;
00056    fX         = 0;
00057    fY         = 0;
00058    fZ         = 0;
00059    fTheta     = 0;
00060    fPhi       = 0;
00061    SetBit(kTemporary,kFALSE);
00062 }
00063 
00064 
00065 //______________________________________________________________________________
00066 TMarker3DBox::TMarker3DBox( Float_t x, Float_t y, Float_t z,
00067                             Float_t dx, Float_t dy, Float_t dz,
00068                             Float_t theta, Float_t phi)
00069               :TAttLine(1,1,1), TAttFill(1,0)
00070 {
00071    // Marker3DBox normal constructor
00072 
00073    fDx        = dx;
00074    fDy        = dy;
00075    fDz        = dz;
00076    fX         = x;
00077    fY         = y;
00078    fZ         = z;
00079    fTheta     = theta;
00080    fPhi       = phi;
00081    fRefObject = 0;
00082    SetBit(kTemporary,kFALSE);
00083 }
00084 
00085 //______________________________________________________________________________
00086 TMarker3DBox::TMarker3DBox(const TMarker3DBox& m3d) :
00087  TObject(m3d),
00088  TAttLine(m3d),
00089  TAttFill(m3d),
00090  TAtt3D(m3d),
00091  fX(m3d.fX),
00092  fY(m3d.fY),
00093  fZ(m3d.fZ),
00094  fDx(m3d.fDx),
00095  fDy(m3d.fDy),
00096  fDz(m3d.fDz),
00097  fTheta(m3d.fTheta),
00098  fPhi(m3d.fPhi),
00099  fRefObject(m3d.fRefObject)
00100 {
00101    //copy constructor
00102 }
00103 
00104 //______________________________________________________________________________
00105 TMarker3DBox& TMarker3DBox::operator=(const TMarker3DBox& m3d)
00106 {
00107    //assignement operator
00108    if(this!=&m3d) {
00109       TObject::operator=(m3d);
00110       TAttLine::operator=(m3d);
00111       TAttFill::operator=(m3d);
00112       TAtt3D::operator=(m3d);
00113       fX=m3d.fX;
00114       fY=m3d.fY;
00115       fZ=m3d.fZ;
00116       fDx=m3d.fDx;
00117       fDy=m3d.fDy;
00118       fDz=m3d.fDz;
00119       fTheta=m3d.fTheta;
00120       fPhi=m3d.fPhi;
00121       fRefObject=m3d.fRefObject;
00122    }
00123    return *this;
00124 }
00125 
00126 //______________________________________________________________________________
00127 TMarker3DBox::~TMarker3DBox()
00128 {
00129    // Marker3DBox shape default destructor
00130 }
00131 
00132 
00133 //______________________________________________________________________________
00134 Int_t TMarker3DBox::DistancetoPrimitive(Int_t px, Int_t py)
00135 {
00136    // Compute distance from point px,py to a Marker3DBox
00137    //
00138    // Compute the closest distance of approach from point px,py to each corner
00139    // point of the Marker3DBox.
00140 
00141    const Int_t numPoints = 8;
00142    Int_t dist = 9999;
00143    Double_t points[3*numPoints];
00144 
00145    TView *view = gPad->GetView();
00146    if (!view) return dist;
00147    const Int_t seg1[12] = {0,1,2,3,4,5,6,7,0,1,2,3};
00148    const Int_t seg2[12] = {1,2,3,0,5,6,7,4,4,5,6,7};
00149 
00150    SetPoints(points);
00151 
00152    Int_t i, i1, i2, dsegment;
00153    Double_t x1,y1,x2,y2;
00154    Double_t xndc[3];
00155    for (i = 0; i < 12; i++) {
00156       i1 = 3*seg1[i];
00157       view->WCtoNDC(&points[i1], xndc);
00158       x1 = xndc[0];
00159       y1 = xndc[1];
00160 
00161       i2 = 3*seg2[i];
00162       view->WCtoNDC(&points[i2], xndc);
00163       x2 = xndc[0];
00164       y2 = xndc[1];
00165       dsegment = DistancetoLine(px,py,x1,y1,x2,y2);
00166       if (dsegment < dist) dist = dsegment;
00167    }
00168    if (dist < 5) {
00169       gPad->SetCursor(kCross);
00170       if (fRefObject) {gPad->SetSelected(fRefObject); return 0;}
00171    }
00172    return dist;
00173 }
00174 
00175 
00176 //______________________________________________________________________________
00177 void TMarker3DBox::ExecuteEvent(Int_t event, Int_t px, Int_t py)
00178 {
00179    // Execute action corresponding to one event
00180    //
00181    // This member function must be implemented to realize the action
00182    // corresponding to the mouse click on the object in the window
00183 
00184    if (gPad->GetView()) gPad->GetView()->ExecuteRotateView(event, px, py);
00185 }
00186 
00187 
00188 //______________________________________________________________________________
00189 void TMarker3DBox::Paint(Option_t * /* option */ )
00190 {
00191    // Paint marker 3D box.
00192 
00193    static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
00194 
00195    buffer.ClearSectionsValid();
00196 
00197    // Section kCore
00198 
00199    // If we are just a temporary object then no 'real object' to
00200    // pass to viewer
00201    if (TestBit(kTemporary)) {
00202       buffer.fID = 0;
00203    } else {
00204       buffer.fID           = this;
00205    }
00206    buffer.fColor        = GetLineColor();
00207    buffer.fTransparency = 0;
00208    buffer.fLocalFrame   = kFALSE;
00209    buffer.SetSectionsValid(TBuffer3D::kCore);
00210 
00211    // We fill kCore and kRawSizes on first pass and try with viewer
00212    Int_t reqSections = gPad->GetViewer3D()->AddObject(buffer);
00213    if (reqSections == TBuffer3D::kNone) {
00214       return;
00215    }
00216 
00217    if (reqSections & TBuffer3D::kRawSizes) {
00218       Int_t nbPnts = 8;
00219       Int_t nbSegs = 12;
00220       Int_t nbPols = 6;
00221       if (!buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
00222          return;
00223       }
00224       buffer.SetSectionsValid(TBuffer3D::kRawSizes);
00225    }
00226 
00227    if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
00228       // Points
00229       SetPoints(buffer.fPnts);
00230 
00231       // Transform points
00232       if (gGeometry && !buffer.fLocalFrame) {
00233          Double_t dlocal[3];
00234          Double_t dmaster[3];
00235          for (UInt_t j=0; j<buffer.NbPnts(); j++) {
00236             dlocal[0] = buffer.fPnts[3*j];
00237             dlocal[1] = buffer.fPnts[3*j+1];
00238             dlocal[2] = buffer.fPnts[3*j+2];
00239             gGeometry->Local2Master(&dlocal[0],&dmaster[0]);
00240             buffer.fPnts[3*j]   = dmaster[0];
00241             buffer.fPnts[3*j+1] = dmaster[1];
00242             buffer.fPnts[3*j+2] = dmaster[2];
00243          }
00244       }
00245 
00246       // Basic colors: 0, 1, ... 8
00247       Int_t c = (((GetLineColor()) %8) -1) * 4;
00248       if (c < 0) c = 0;
00249 
00250       // Segments
00251       buffer.fSegs[ 0] = c   ; buffer.fSegs[ 1] = 0 ; buffer.fSegs[ 2] = 1;
00252       buffer.fSegs[ 3] = c+1 ; buffer.fSegs[ 4] = 1 ; buffer.fSegs[ 5] = 2;
00253       buffer.fSegs[ 6] = c+1 ; buffer.fSegs[ 7] = 2 ; buffer.fSegs[ 8] = 3;
00254       buffer.fSegs[ 9] = c   ; buffer.fSegs[10] = 3 ; buffer.fSegs[11] = 0;
00255       buffer.fSegs[12] = c+2 ; buffer.fSegs[13] = 4 ; buffer.fSegs[14] = 5;
00256       buffer.fSegs[15] = c+2 ; buffer.fSegs[16] = 5 ; buffer.fSegs[17] = 6;
00257       buffer.fSegs[18] = c+3 ; buffer.fSegs[19] = 6 ; buffer.fSegs[20] = 7;
00258       buffer.fSegs[21] = c+3 ; buffer.fSegs[22] = 7 ; buffer.fSegs[23] = 4;
00259       buffer.fSegs[24] = c   ; buffer.fSegs[25] = 0 ; buffer.fSegs[26] = 4;
00260       buffer.fSegs[27] = c+2 ; buffer.fSegs[28] = 1 ; buffer.fSegs[29] = 5;
00261       buffer.fSegs[30] = c+1 ; buffer.fSegs[31] = 2 ; buffer.fSegs[32] = 6;
00262       buffer.fSegs[33] = c+3 ; buffer.fSegs[34] = 3 ; buffer.fSegs[35] = 7;
00263 
00264       // Polygons
00265       buffer.fPols[ 0] = c   ; buffer.fPols[ 1] = 4 ; buffer.fPols[ 2] = 0;
00266       buffer.fPols[ 3] = 9   ; buffer.fPols[ 4] = 4 ; buffer.fPols[ 5] = 8;
00267       buffer.fPols[ 6] = c+1 ; buffer.fPols[ 7] = 4 ; buffer.fPols[ 8] = 1;
00268       buffer.fPols[ 9] = 10  ; buffer.fPols[10] = 5 ; buffer.fPols[11] = 9;
00269       buffer.fPols[12] = c   ; buffer.fPols[13] = 4 ; buffer.fPols[14] = 2;
00270       buffer.fPols[15] = 11  ; buffer.fPols[16] = 6 ; buffer.fPols[17] = 10;
00271       buffer.fPols[18] = c+1 ; buffer.fPols[19] = 4 ; buffer.fPols[20] = 3;
00272       buffer.fPols[21] = 8   ; buffer.fPols[22] = 7 ; buffer.fPols[23] = 11;
00273       buffer.fPols[24] = c+2 ; buffer.fPols[25] = 4 ; buffer.fPols[26] = 0;
00274       buffer.fPols[27] = 3   ; buffer.fPols[28] = 2 ; buffer.fPols[29] = 1;
00275       buffer.fPols[30] = c+3 ; buffer.fPols[31] = 4 ; buffer.fPols[32] = 4;
00276       buffer.fPols[33] = 5   ; buffer.fPols[34] = 6 ; buffer.fPols[35] = 7;
00277 
00278       buffer.SetSectionsValid(TBuffer3D::kRaw);
00279 
00280       TAttLine::Modify();
00281       TAttFill::Modify();
00282    }
00283 
00284    gPad->GetViewer3D()->AddObject(buffer);
00285 }
00286 
00287 
00288 //______________________________________________________________________________
00289 void TMarker3DBox::PaintH3(TH1 *h, Option_t *option)
00290 {
00291    // Paint 3-d histogram h with marker3dboxes
00292 
00293    Int_t bin,ix,iy,iz;
00294    Double_t xmin,xmax,ymin,ymax,zmin,zmax,wmin,wmax,w;
00295    TAxis *xaxis = h->GetXaxis();
00296    TAxis *yaxis = h->GetYaxis();
00297    TAxis *zaxis = h->GetZaxis();
00298 
00299    //compute min and max of all cells
00300    wmin = wmax = 0;
00301    for (iz=zaxis->GetFirst();iz<=zaxis->GetLast();iz++) {
00302       for (iy=yaxis->GetFirst();iy<=yaxis->GetLast();iy++) {
00303          for (ix=xaxis->GetFirst();ix<=xaxis->GetLast();ix++) {
00304             bin = h->GetBin(ix,iy,iz);
00305             w = h->GetBinContent(bin);
00306             if (w > wmax) wmax = w;
00307             if (w < wmin) wmin = w;
00308          }
00309       }
00310    }
00311 
00312    //Create or modify 3-d view object
00313    TView *view = gPad->GetView();
00314    if (!view) {
00315       gPad->Range(-1,-1,1,1);
00316       view = TView::CreateView(1,0,0);
00317    }
00318    view->SetRange(xaxis->GetBinLowEdge(xaxis->GetFirst()),
00319                   yaxis->GetBinLowEdge(yaxis->GetFirst()),
00320                   zaxis->GetBinLowEdge(zaxis->GetFirst()),
00321                   xaxis->GetBinUpEdge(xaxis->GetLast()),
00322                   yaxis->GetBinUpEdge(yaxis->GetLast()),
00323                   zaxis->GetBinUpEdge(zaxis->GetLast()));
00324 
00325    view->PadRange(gPad->GetFrameFillColor());
00326 
00327    //Draw TMarker3DBox with size proportional to cell content
00328    TMarker3DBox m3;
00329    m3.SetBit(kTemporary,kTRUE);
00330    m3.SetRefObject(h);
00331    m3.SetDirection(0,0);
00332    m3.SetLineColor(h->GetMarkerColor());
00333    Double_t scale;
00334    for (ix=xaxis->GetFirst();ix<=xaxis->GetLast();ix++) {
00335       xmin = h->GetXaxis()->GetBinLowEdge(ix);
00336       xmax = xmin + h->GetXaxis()->GetBinWidth(ix);
00337       for (iy=yaxis->GetFirst();iy<=yaxis->GetLast();iy++) {
00338          ymin = h->GetYaxis()->GetBinLowEdge(iy);
00339          ymax = ymin + h->GetYaxis()->GetBinWidth(iy);
00340          for (iz=zaxis->GetFirst();iz<=zaxis->GetLast();iz++) {
00341             zmin = h->GetZaxis()->GetBinLowEdge(iz);
00342             zmax = zmin + h->GetZaxis()->GetBinWidth(iz);
00343             bin = h->GetBin(ix,iy,iz);
00344             w = h->GetBinContent(bin);
00345             if (w == 0) continue;
00346             scale = (w-wmin)/(wmax-wmin);
00347             m3.SetPosition(0.5*(xmin+xmax),0.5*(ymin+ymax),0.5*(zmin+zmax));
00348             m3.SetSize(scale*(xmax-xmin),scale*(ymax-ymin),scale*(zmax-zmin));
00349             m3.Paint(option);
00350          }
00351       }
00352    }
00353 }
00354 
00355 
00356 //______________________________________________________________________________
00357 void TMarker3DBox::SavePrimitive(ostream &out, Option_t * /*= ""*/)
00358 {
00359    // Save primitive as a C++ statement(s) on output stream out
00360 
00361    out<<"   "<<endl;
00362    if (gROOT->ClassSaved(TMarker3DBox::Class())) {
00363       out<<"   ";
00364    } else {
00365       out<<"   TMarker3DBox *";
00366    }
00367    out<<"marker3DBox = new TMarker3DBox("<<fX<<","
00368                                          <<fY<<","
00369                                          <<fZ<<","
00370                                          <<fDx<<","
00371                                          <<fDy<<","
00372                                          <<fDz<<","
00373                                          <<fTheta<<","
00374                                          <<fPhi<<");"<<endl;
00375 
00376    SaveLineAttributes(out,"marker3DBox",1,1,1);
00377    SaveFillAttributes(out,"marker3DBox",1,0);
00378 
00379    out<<"   marker3DBox->Draw();"<<endl;
00380 }
00381 
00382 
00383 //______________________________________________________________________________
00384 void TMarker3DBox::SetDirection(Float_t theta, Float_t phi)
00385 {
00386    // Set direction.
00387 
00388    fTheta = theta;
00389    fPhi   = phi;
00390 }
00391 
00392 
00393 //______________________________________________________________________________
00394 void TMarker3DBox::SetSize(Float_t dx, Float_t dy, Float_t dz)
00395 {
00396    // Set size.
00397 
00398    fDx = dx;
00399    fDy = dy;
00400    fDz = dz;
00401 }
00402 
00403 
00404 //______________________________________________________________________________
00405 void TMarker3DBox::SetPosition(Float_t x, Float_t y, Float_t z)
00406 {
00407    // Set position.
00408 
00409    fX  = x;
00410    fY  = y;
00411    fZ  = z;
00412 }
00413 
00414 
00415 //______________________________________________________________________________
00416 void TMarker3DBox::SetPoints(Double_t *points) const
00417 {
00418    // Set points.
00419 
00420    if (points) {
00421       points[ 0] = -fDx ; points[ 1] = -fDy ; points[ 2] = -fDz;
00422       points[ 3] = -fDx ; points[ 4] =  fDy ; points[ 5] = -fDz;
00423       points[ 6] =  fDx ; points[ 7] =  fDy ; points[ 8] = -fDz;
00424       points[ 9] =  fDx ; points[10] = -fDy ; points[11] = -fDz;
00425       points[12] = -fDx ; points[13] = -fDy ; points[14] =  fDz;
00426       points[15] = -fDx ; points[16] =  fDy ; points[17] =  fDz;
00427       points[18] =  fDx ; points[19] =  fDy ; points[20] =  fDz;
00428       points[21] =  fDx ; points[22] = -fDy ; points[23] =  fDz;
00429 
00430       Double_t x, y, z;
00431       const Double_t kPI = TMath::Pi();
00432       Double_t theta  = fTheta*kPI/180;
00433       Double_t phi    = fPhi*kPI/180;
00434       Double_t sinth = TMath::Sin(theta);
00435       Double_t costh = TMath::Cos(theta);
00436       Double_t sinfi = TMath::Sin(phi);
00437       Double_t cosfi = TMath::Cos(phi);
00438 
00439       // Matrix to convert from fruit frame to master frame
00440       Double_t m[9];
00441       m[0] =  costh * cosfi;       m[1] = -sinfi;          m[2] = sinth*cosfi;
00442       m[3] =  costh * sinfi;       m[4] =  cosfi;          m[5] = sinth*sinfi;
00443       m[6] = -sinth;               m[7] =  0;              m[8] = costh;
00444       for (Int_t i = 0; i < 8; i++) {
00445          x = points[3*i];
00446          y = points[3*i+1];
00447          z = points[3*i+2];
00448 
00449          points[3*i]   = fX + m[0] * x + m[1] * y + m[2] * z;
00450          points[3*i+1] = fY + m[3] * x + m[4] * y + m[5] * z;
00451          points[3*i+2] = fZ + m[6] * x + m[7] * y + m[8] * z;
00452       }
00453    }
00454 }
00455 
00456 
00457 //______________________________________________________________________________
00458 void TMarker3DBox::Streamer(TBuffer &R__b)
00459 {
00460    // Stream an object of class TMarker3DBox.
00461 
00462    if (R__b.IsReading()) {
00463       UInt_t R__s, R__c;
00464       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
00465       if (R__v > 1) {
00466          R__b.ReadClassBuffer(TMarker3DBox::Class(), this, R__v, R__s, R__c);
00467          return;
00468       }
00469       //====process old versions before automatic schema evolution
00470       TObject::Streamer(R__b);
00471       TAttLine::Streamer(R__b);
00472       TAttFill::Streamer(R__b);
00473       if (R__b.GetVersionOwner() > 22300) {
00474          TAtt3D::Streamer(R__b);
00475       } else {
00476          TAtt3D::Streamer(R__b);
00477       }
00478       R__b >> fX;
00479       R__b >> fY;
00480       R__b >> fZ;
00481       R__b >> fDx;
00482       R__b >> fDy;
00483       R__b >> fDz;
00484       R__b >> fTheta;
00485       R__b >> fPhi;
00486       R__b >> fRefObject;
00487       R__b.CheckByteCount(R__s, R__c, TMarker3DBox::IsA());
00488       //====end of old versions
00489 
00490    } else {
00491       R__b.WriteClassBuffer(TMarker3DBox::Class(),this);
00492    }
00493 }

Generated on Tue Jul 5 14:16:40 2011 for ROOT_528-00b_version by  doxygen 1.5.1