TCrown.cxx

Go to the documentation of this file.
00001 // @(#)root/graf:$Id: TCrown.cxx 27688 2009-03-04 09:48:43Z couet $
00002 // Author: Rene Brun   108/08/2002
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include "Riostream.h"
00013 #include "TROOT.h"
00014 #include "TMath.h"
00015 #include "TCrown.h"
00016 #include "TVirtualPad.h"
00017 
00018 ClassImp(TCrown)
00019 
00020 
00021 //______________________________________________________________________________
00022 /* Begin_Html
00023 <center><h2>TCrown : to draw crown</h2></center>
00024 A crown is specified with the position of its centre, its inner/outer radius
00025 a minimum and maximum angle. The attributes of the outline line are given via
00026 TAttLine. The attributes of the fill area are given via TAttFill.
00027 <p> Example:
00028 End_Html
00029 Begin_Macro(source)
00030 {
00031    TCanvas *c1 = new TCanvas("c1","c1",400,400);
00032    TCrown cr1(.5,.5,.3,.4);
00033    cr1->SetLineStyle(2);
00034    cr1->SetLineWidth(4);
00035    cr1.Draw();
00036    TCrown cr2(.5,.5,.2,.3,45,315);
00037    cr2.SetFillColor(38);
00038    cr2.SetFillStyle(3010);
00039    cr2.Draw();
00040    TCrown cr3(.5,.5,.2,.3,-45,45);
00041    cr3.SetFillColor(50);
00042    cr3.SetFillStyle(3025);
00043    cr3.Draw();
00044    TCrown cr4(.5,.5,.0,.2);
00045    cr4.SetFillColor(4);
00046    cr4.SetFillStyle(3008);
00047    cr4.Draw();
00048    return c1;
00049 }
00050 End_Macro */
00051 
00052 
00053 //______________________________________________________________________________
00054 TCrown::TCrown(): TEllipse()
00055 {
00056    /* Begin_Html
00057    Crown default constructor.
00058    End_Html */
00059 
00060 }
00061 
00062 
00063 //______________________________________________________________________________
00064 TCrown::TCrown(Double_t x1, Double_t y1,Double_t radin, Double_t radout,Double_t phimin,Double_t phimax)
00065       :TEllipse(x1,y1,radin,radout,phimin,phimax,0)
00066 {
00067    /* Begin_Html
00068    Crown normal constructor.
00069    <ul>
00070    <li> x1,y1  : coordinates of centre of crown
00071    <li> radin  : inner crown radius
00072    <li> radout : outer crown radius
00073    <li> phimin : min angle in degrees (default is 0)
00074    <li> phimax : max angle in degrees (default is 360)
00075    </ul>
00076    When a crown sector only is drawn, the lines connecting the center
00077    of the crown to the edges are drawn by default. One can specify
00078    the drawing option "only" to not draw these lines.
00079    End_Html */
00080 
00081 }
00082 
00083 
00084 //______________________________________________________________________________
00085 TCrown::TCrown(const TCrown &crown) : TEllipse(crown)
00086 {
00087    /* Begin_Html
00088    Crown copy constructor.
00089    End_Html */
00090 
00091    ((TCrown&)crown).Copy(*this);
00092 }
00093 
00094 
00095 //______________________________________________________________________________
00096 TCrown::~TCrown()
00097 {
00098    /* Begin_Html
00099    Crown default destructor.
00100    End_Html */
00101 }
00102 
00103 
00104 //______________________________________________________________________________
00105 void TCrown::Copy(TObject &crown) const
00106 {
00107    /* Begin_Html
00108    Copy this crown to crown.
00109    End_Html */
00110 
00111    TEllipse::Copy(crown);
00112 }
00113 
00114 
00115 //______________________________________________________________________________
00116 Int_t TCrown::DistancetoPrimitive(Int_t px, Int_t py)
00117 {
00118    /* Begin_Html
00119    Compute distance from point px,py to a crown.
00120    <p>
00121    If crown is filled, return OK if we are inside
00122    otherwise, crown is found if near the crown edges.
00123    End_Html */
00124 
00125    const Double_t kPI = TMath::Pi();
00126    Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px)) - fX1;
00127    Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py)) - fY1;
00128 
00129    Double_t r1 = fR1;
00130    Double_t r2 = fR2;
00131    Double_t r  = TMath::Sqrt(x*x+y*y);
00132 
00133    if (r1>r2) {
00134       r1 = fR2;
00135       r2 = fR1;
00136    }
00137 
00138    Int_t dist = 9999;
00139    if (r > r2) return dist;
00140    if (r < r1) return dist;
00141    if (fPhimax-fPhimin < 360) {
00142       Double_t phi = 180*TMath::ACos(x/r)/kPI;
00143       if (y<0) phi = 360-phi;
00144       Double_t phi1 = fPhimin;
00145       Double_t phi2 = fPhimax;
00146       if (phi1<0) phi1=phi1+360;
00147       if (phi2<0) phi2=phi2+360;
00148       if (phi2<phi1) {
00149          if (phi < phi1 && phi > phi2) return dist;
00150       } else {
00151          if (phi < phi1) return dist;
00152          if (phi > phi2) return dist;
00153       }
00154    }
00155 
00156    if (GetFillColor() && GetFillStyle()) {
00157       return 0;
00158    } else {
00159       if (TMath::Abs(r2-r)/r2 < 0.02) return 0;
00160       if (TMath::Abs(r1-r)/r1 < 0.02) return 0;
00161    }
00162    return dist;
00163 }
00164 
00165 
00166 //______________________________________________________________________________
00167 void TCrown::DrawCrown(Double_t x1, Double_t y1,Double_t radin,Double_t radout,Double_t phimin,Double_t phimax,Option_t *option)
00168 {
00169    /* Begin_Html
00170    Draw this crown with new coordinates.
00171    End_Html */
00172 
00173    TCrown *newcrown = new TCrown(x1, y1, radin, radout, phimin, phimax);
00174    TAttLine::Copy(*newcrown);
00175    TAttFill::Copy(*newcrown);
00176    newcrown->SetBit(kCanDelete);
00177    newcrown->AppendPad(option);
00178 }
00179 
00180 
00181 //______________________________________________________________________________
00182 void TCrown::ExecuteEvent(Int_t event, Int_t px, Int_t py)
00183 {
00184    /* Begin_Html
00185    Execute action corresponding to one event
00186    <p>
00187    For the time being TEllipse::ExecuteEvent is used.
00188    End_Html */
00189 
00190    TEllipse::ExecuteEvent(event,px,py);
00191 }
00192 
00193 
00194 //______________________________________________________________________________
00195 void TCrown::Paint(Option_t *)
00196 {
00197    /* Begin_Html
00198    Paint this crown with its current attributes.
00199    End_Html */
00200 
00201    const Double_t kPI = TMath::Pi();
00202    const Int_t np = 40;
00203    static Double_t x[2*np+3], y[2*np+3];
00204    TAttLine::Modify();
00205    TAttFill::Modify();
00206 
00207    Double_t angle,dx,dy;
00208    Double_t dphi = (fPhimax-fPhimin)*kPI/(180*np);
00209    Double_t ct   = TMath::Cos(kPI*fTheta/180);
00210    Double_t st   = TMath::Sin(kPI*fTheta/180);
00211    Int_t i;
00212    //compute outer points
00213    for (i=0;i<=np;i++) {
00214       angle = fPhimin*kPI/180 + Double_t(i)*dphi;
00215       dx    = fR2*TMath::Cos(angle);
00216       dy    = fR2*TMath::Sin(angle);
00217       x[i]  = fX1 + dx*ct - dy*st;
00218       y[i]  = fY1 + dx*st + dy*ct;
00219    }
00220    //compute inner points
00221    for (i=0;i<=np;i++) {
00222       angle = fPhimin*kPI/180 + Double_t(i)*dphi;
00223       dx    = fR1*TMath::Cos(angle);
00224       dy    = fR1*TMath::Sin(angle);
00225       x[2*np-i+1]  = fX1 + dx*ct - dy*st;
00226       y[2*np-i+1]  = fY1 + dx*st + dy*ct;
00227    }
00228    x[2*np+2]  = x[0];
00229    y[2*np+2]  = y[0];
00230    if (fPhimax-fPhimin >= 360 ) {
00231       // a complete filled crown
00232       if (GetFillColor()  && GetFillStyle()) {
00233          gPad->PaintFillArea(2*np+2,x,y);
00234       }
00235       // a complete empty crown
00236       if (GetLineStyle()) {
00237          gPad->PaintPolyLine(np+1,x,y);
00238          gPad->PaintPolyLine(np+1,&x[np+1],&y[np+1]);
00239       }
00240    } else {
00241       //crown segment
00242       if (GetFillColor()  && GetFillStyle()) gPad->PaintFillArea(2*np+2,x,y);
00243       if (GetLineStyle()) gPad->PaintPolyLine(2*np+3,x,y);
00244    }
00245 }
00246 
00247 
00248 //______________________________________________________________________________
00249 void TCrown::SavePrimitive(ostream &out, Option_t * /*= ""*/)
00250 {
00251    /* Begin_Html
00252    Save primitive as a C++ statement(s) on output stream out.
00253    End_Html */
00254 
00255    out<<"   "<<endl;
00256    if (gROOT->ClassSaved(TCrown::Class())) {
00257       out<<"   ";
00258    } else {
00259       out<<"   TCrown *";
00260    }
00261    out<<"crown = new TCrown("<<fX1<<","<<fY1<<","<<fR1<<","<<fR2
00262       <<","<<fPhimin<<","<<fPhimax<<");"<<endl;
00263 
00264    SaveFillAttributes(out,"crown",0,1001);
00265    SaveLineAttributes(out,"crown",1,1,1);
00266 
00267    if (GetNoEdges()) out<<"   crown->SetNoEdges();"<<endl;
00268 
00269    out<<"   crown->Draw();"<<endl;
00270 }

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