00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042 #include "TGraphPolar.h"
00043 #include "TGraphPolargram.h"
00044 #include "TGaxis.h"
00045 #include "THLimitsFinder.h"
00046 #include "TVirtualPad.h"
00047 #include "TROOT.h"
00048 #include "TLatex.h"
00049 #include "TEllipse.h"
00050 #include "TMath.h"
00051
00052
00053 ClassImp(TGraphPolargram);
00054
00055
00056 TGraphPolargram::TGraphPolargram(const char* name, Double_t rmin, Double_t rmax,
00057 Double_t tmin, Double_t tmax):
00058 TNamed(name,"Polargram")
00059 {
00060
00061 Init();
00062 fNdivRad = 508;
00063 fNdivPol = 508;
00064 fPolarLabels = NULL;
00065 fRwrmax = rmax;
00066 fRwrmin = rmin;
00067 fRwtmin = tmin;
00068 fRwtmax = tmax;
00069 }
00070
00071
00072
00073 TGraphPolargram::TGraphPolargram(const char* name):
00074 TNamed(name,"Polargram")
00075 {
00076
00077
00078 Init();
00079 fNdivRad = 0;
00080 fNdivPol = 0;
00081 fPolarLabels = NULL;
00082 fRwrmax = 1;
00083 fRwrmin = 0;
00084 fRwtmax = 0;
00085 fRwtmin = 0;
00086 }
00087
00088
00089
00090 TGraphPolargram::~TGraphPolargram()
00091 {
00092
00093
00094 if (fPolarLabels != NULL) delete [] fPolarLabels;
00095 }
00096
00097
00098
00099 void TGraphPolargram::ChangeRangePolar(Double_t tmin, Double_t tmax)
00100 {
00101
00102
00103
00104
00105 if (tmin < tmax) {
00106 fRwtmin = tmin;
00107 fRwtmax = tmax;
00108 }
00109 if (gPad) gPad->Modified();
00110 }
00111
00112
00113
00114 Int_t TGraphPolargram::DistancetoPrimitive(Int_t px, Int_t py)
00115 {
00116
00117
00118 Int_t i;
00119 Double_t x = gPad->AbsPixeltoX(px);
00120 Double_t y = gPad->AbsPixeltoY(py);
00121
00122
00123 Double_t rad = TMath::Sqrt(x*x+y*y);
00124 Int_t div = (Int_t)rad*(fNdivRad%100);
00125 Double_t dr = TMath::Min(TMath::Abs(rad-div*1./(fNdivRad%100)),
00126 TMath::Abs(rad-(div+1)*1./(fNdivRad%100)));
00127 Int_t drad = gPad->XtoPixel(dr)-gPad->XtoPixel(0);
00128
00129
00130
00131 Int_t dt = kMaxPixel;
00132 for (i=0; i<(fNdivPol%100); i++) {
00133 Double_t theta = i*2*TMath::Pi()/(fNdivPol%100);
00134
00135
00136 Int_t dthis = DistancetoLine(px,py,0.,0.,TMath::Cos(theta),
00137 TMath::Sin(theta));
00138
00139
00140
00141 if (dthis==9999) {
00142
00143
00144 if (rad>1) {
00145 dthis = (Int_t)TMath::Sqrt(
00146 TMath::Power(px-gPad->XtoPixel(TMath::Cos(theta)),2)+
00147 TMath::Power(py-gPad->YtoPixel(TMath::Sin(theta)),2));
00148 } else {
00149
00150
00151 if (((TMath::Abs(theta-TMath::Pi())<0.1) &&
00152 ((px-gPad->XtoPixel(0))<0)) ||
00153 ((TMath::Abs(theta)<0.1) &&
00154 ((px-gPad->XtoPixel(0))>0))) {
00155 dthis = TMath::Abs(py-gPad->YtoPixel(0.));
00156 }
00157
00158
00159 if (((TMath::Abs(theta-TMath::PiOver2())<0.1) &&
00160 ((py-gPad->YtoPixel(0))>0)) ||
00161 ((TMath::Abs(theta-3*TMath::PiOver2())<0.1) &&
00162 (py-gPad->YtoPixel(0))<0)) {
00163 dthis = TMath::Abs(px-gPad->XtoPixel(0.));
00164 }
00165 if (dthis==9999) {
00166
00167
00168
00169 dthis = (Int_t)TMath::Sqrt(
00170 TMath::Power(px-gPad->XtoPixel(0.),2)+
00171 TMath::Power(py-gPad->YtoPixel(0.),2));
00172 }
00173 }
00174 }
00175
00176
00177 dt = TMath::Min(dthis,dt);
00178 }
00179 return TMath::Min(drad, dt);
00180 }
00181
00182
00183
00184 void TGraphPolargram::Draw(Option_t* options)
00185 {
00186
00187
00188 Paint(options);
00189 AppendPad(options);
00190 }
00191
00192
00193
00194 void TGraphPolargram::ExecuteEvent(Int_t event, Int_t px, Int_t py)
00195 {
00196
00197
00198 Int_t kMaxDiff = 20;
00199 static Int_t d1, d2, d3, px1, py1, px3, py3, px4, py4;
00200 static Bool_t p1, p2, p3, p4, p5, p6, p7, p8;
00201 Double_t px2, py2;
00202 p2 = p3 = p4 = p5 = p6 = p7 = p8 = kFALSE;
00203 if (!gPad->IsEditable()) return;
00204 switch (event) {
00205 case kMouseMotion:
00206 px1 = gPad->XtoAbsPixel(TMath::Cos(GetAngle()));
00207 py1 = gPad->YtoAbsPixel(TMath::Sin(GetAngle()));
00208 d1 = TMath::Abs(px1 - px) + TMath::Abs(py1-py);
00209 p1 = kFALSE;
00210 px2 = gPad->XtoAbsPixel(-1);
00211 py2 = gPad->YtoAbsPixel(1);
00212 d2 = (Int_t)(TMath::Abs(px2 - px) + TMath::Abs(py2 - py)) ;
00213 px3 = gPad->XtoAbsPixel(-1);
00214 py3 = gPad->YtoAbsPixel(-1);
00215 d3 = TMath::Abs(px3 - px) + TMath::Abs(py3 - py) ;
00216
00217 if (d1 < kMaxDiff) {
00218 gPad->SetCursor(kMove);
00219 p1 = kTRUE;
00220 }
00221
00222 if ( d2 < kMaxDiff) {
00223 gPad->SetCursor(kHand);
00224 p7 = kTRUE;
00225 }
00226
00227 if ( d3 < kMaxDiff) {
00228 gPad->SetCursor(kHand);
00229 p8 = kTRUE;
00230 }
00231
00232 if (!p1 && !p7 ) {
00233 p6 = kTRUE;
00234 gPad->SetCursor(kHand);
00235 }
00236 break;
00237
00238 case kButton1Down:
00239
00240 px4 = px;
00241 py4 = py;
00242
00243 case kButton1Motion:
00244 if (p1) {
00245 px2 = gPad->AbsPixeltoX(px);
00246 py2 = gPad->AbsPixeltoY(py);
00247 if ( px2 < 0 && py2 < 0) {p2 = kTRUE;};
00248 if ( px2 < 0 && py2 > 0 ) {p3 = kTRUE;};
00249 if ( px2 > 0 && py2 > 0 ) {p4 = kTRUE;};
00250 if ( px2 > 0 && py2 < 0 ) {p5 = kTRUE;};
00251 px2 = TMath::ACos(TMath::Abs(px2));
00252 py2 = TMath::ASin(TMath::Abs(py2));
00253 if (p2) {
00254 fAxisAngle = TMath::Pi()+(px2+py2)/2;
00255 p2 = kFALSE;
00256 };
00257 if (p3) {
00258 fAxisAngle = TMath::Pi()-(px2+py2)/2;
00259 p3 = kFALSE;
00260 };
00261 if (p4) {
00262 fAxisAngle = (px2+py2)/2;
00263 p4 = kFALSE;
00264 };
00265 if (p5) {
00266 fAxisAngle = -(px2+py2)/2;
00267 p5 = kFALSE;
00268 };
00269 }
00270 break;
00271
00272 case kButton1Up:
00273 Paint();
00274 }
00275 }
00276
00277
00278
00279 Int_t TGraphPolargram::FindAlign(Double_t angle)
00280 {
00281
00282
00283 Double_t pi = TMath::Pi();
00284
00285 while(angle < 0 || angle > 2*pi){
00286 if(angle < 0) angle+=2*pi;
00287 if(angle > 2*pi) angle-=2*pi;
00288 }
00289 if(!TestBit(TGraphPolargram::kLabelOrtho)){
00290 if(angle > 0 && angle < pi/2) return 11;
00291 else if(angle > pi/2 && angle < pi) return 31;
00292 else if(angle > pi && angle < 3*pi/2) return 33;
00293 else if(angle > 3*pi/2 && angle < 2*pi) return 13;
00294 else if(angle == 0 || angle == 2*pi) return 12;
00295 else if(angle == pi/2) return 21;
00296 else if(angle == pi) return 32;
00297 else if(angle == 3*pi/2) return 23;
00298 else return 0;
00299 }
00300 else{
00301 if(angle >= 0 && angle <= pi/2) return 12;
00302 else if((angle > pi/2 && angle <= pi) || (angle > pi && angle <= 3*pi/2)) return 32;
00303 else if(angle > 3*pi/2 && angle <= 2*pi) return 12;
00304 else return 0;
00305 }
00306 }
00307
00308
00309
00310 Double_t TGraphPolargram::FindTextAngle(Double_t angle)
00311 {
00312
00313
00314 Double_t pi = TMath::Pi();
00315 Double_t convraddeg = 180.0/pi;
00316
00317 while(angle < 0 || angle > 2*pi){
00318 if(angle < 0) angle+=2*pi;
00319 if(angle > 2*pi) angle-=2*pi;
00320 }
00321
00322 if(angle >= 0 && angle <= pi/2) return angle*convraddeg;
00323 else if(angle > pi/2 && angle <= pi) return (angle + pi)*convraddeg;
00324 else if(angle > pi && angle <= 3*pi/2) return (angle - pi)*convraddeg;
00325 else if(angle > 3*pi/2 && angle <= 2*pi) return angle*convraddeg;
00326 else return 0;
00327 }
00328
00329
00330
00331 void TGraphPolargram::Init()
00332 {
00333
00334
00335 fAxisAngle = 0;
00336 fCutRadial = 0;
00337 fDegree = kFALSE;
00338 fGrad = kFALSE;
00339 fLineStyle = 3;
00340 fPolarLabelColor = 1;
00341 fPolarLabelFont = 62;
00342 fPolarOffset = 0.04;
00343 fPolarTextSize = 0.04;
00344 fRadialOffset = 0.025;
00345 fRadian = kTRUE;
00346 fRadialLabelColor = 1;
00347 fRadialLabelFont = 62;
00348 fRadialTextSize = 0.035;
00349 fTickpolarSize = 0.02;
00350 }
00351
00352
00353
00354 void TGraphPolargram::Paint(Option_t * chopt)
00355 {
00356
00357
00358 Int_t optionpoldiv, optionraddiv;
00359 Bool_t optionLabels = kTRUE;
00360
00361 TString opt = chopt;
00362 opt.ToUpper();
00363
00364 if(opt.Contains('P')) optionpoldiv=1; else optionpoldiv=0;
00365 if(opt.Contains('R')) optionraddiv=1; else optionraddiv=0;
00366 if(opt.Contains('O')) SetBit(TGraphPolargram::kLabelOrtho);
00367 else ResetBit(TGraphPolargram::kLabelOrtho);
00368 if(!opt.Contains('P') && !opt.Contains('R')) optionpoldiv=optionraddiv=1;
00369 if(opt.Contains('N')) optionLabels = kFALSE;
00370
00371 if(optionraddiv) PaintRadialDivisions(kTRUE);
00372 else PaintRadialDivisions(kFALSE);
00373 if(optionpoldiv) PaintPolarDivisions(optionLabels);
00374 }
00375
00376
00377
00378 void TGraphPolargram::PaintCircle(Double_t x1, Double_t y1, Double_t r,
00379 Double_t phimin, Double_t phimax, Double_t theta)
00380 {
00381
00382
00383
00384 Int_t i;
00385 const Int_t np = 200;
00386 static Double_t x[np+3], y[np+3];
00387
00388
00389
00390
00391 Double_t circ = TMath::Pi()*2*r*(phimax-phimin)/36;
00392 Int_t n = (Int_t)(np*circ/((gPad->GetX2()-gPad->GetX1())+
00393 (gPad->GetY2()-gPad->GetY1())));
00394 if (n < 8) n = 8;
00395 if (n > np) n = np;
00396 Double_t angle,dx,dy;
00397 Double_t dphi = (phimax-phimin)*TMath::Pi()/(180*n);
00398 Double_t ct = TMath::Cos(TMath::Pi()*theta/180);
00399 Double_t st = TMath::Sin(TMath::Pi()*theta/180);
00400 for (i=0; i<=n; i++) {
00401 angle = phimin*TMath::Pi()/180 + Double_t(i)*dphi;
00402 dx = r*TMath::Cos(angle);
00403 dy = r*TMath::Sin(angle);
00404 x[i] = x1 + dx*ct - dy*st;
00405 y[i] = y1 + dx*st + dy*ct;
00406 }
00407 gPad->PaintPolyLine(n+1,x,y);
00408 }
00409
00410
00411
00412 void TGraphPolargram::PaintPolarDivisions(Bool_t optionLabels)
00413 {
00414
00415
00416
00417 Int_t i, j, rnum, rden, first, last;
00418 if (!gPad) return ;
00419
00420 gPad->RangeAxis(-1,-1,1,1);
00421 gPad->Range(-1.25,-1.25,1.25,1.25);
00422 Int_t ndivMajor = fNdivPol%100;
00423 Int_t ndivMinor = fNdivPol/100;
00424
00425 if (!gPad->GetLogy()) {
00426 for (i=0; i<ndivMajor; i++) {
00427 Double_t txtval = fRwtmin + i*(fRwtmax-fRwtmin)/ndivMajor;
00428 Double_t theta = i*2*TMath::Pi()/ndivMajor;
00429 Double_t costheta = TMath::Cos(theta);
00430 Double_t sintheta = TMath::Sin(theta);
00431 Double_t tantheta = TMath::Tan(theta);
00432 Double_t costhetas = (1+fPolarOffset)*costheta;
00433 Double_t sinthetas = (1+fPolarOffset)*sintheta;
00434 Double_t corr = 0.01;
00435
00436 TLatex *textangular = new TLatex();
00437 textangular->SetTextColor(GetPolarColorLabel());
00438 textangular->SetTextFont(GetPolarLabelFont());
00439
00440 const char* form = (char *)" ";
00441 TGaxis axis;
00442 if (TestBit(TGraphPolargram::kLabelOrtho)) {
00443
00444 if(fPolarLabels == NULL && optionLabels){;
00445 if (fRadian) {
00446
00447 ReduceFraction(2*i, ndivMajor, rnum, rden);
00448 if (rnum == 0) form = Form("%d",rnum);
00449 if (rnum == 1 && rden == 1) form = Form("#pi");
00450 if (rnum == 1 && rden != 1) form = Form("#frac{#pi}{%d}",rden);
00451 if (rnum != 1 && rden == 1 && i !=0) form= Form("%d#pi",rnum);
00452 if (rnum != 1 && rden != 1) form = Form("#frac{%d#pi}{%d}",rnum,rden);
00453 textangular->SetTextAlign(FindAlign(theta));
00454 textangular->PaintLatex(costhetas,
00455 sinthetas, FindTextAngle(theta),
00456 GetPolarLabelSize(), form);
00457 } else {
00458
00459 form = Form("%5.3g",txtval);
00460 axis.LabelsLimits(form,first,last);
00461 TString s = Form("%s",form);
00462 if (first != 0) s.Remove(0, first);
00463 textangular->SetTextAlign(FindAlign(theta));
00464 textangular->PaintLatex(costhetas,
00465 sinthetas, FindTextAngle(theta),
00466 GetPolarLabelSize(), s);
00467 }
00468 } else if (fPolarLabels){
00469
00470 textangular->SetTextAlign(FindAlign(theta));
00471 textangular->PaintLatex(costhetas,sinthetas,FindTextAngle(theta),
00472 GetPolarLabelSize(), fPolarLabels[i]);
00473 }
00474 } else {
00475
00476 if(fPolarLabels == NULL && optionLabels){
00477 if (fRadian) {
00478
00479 ReduceFraction(2*i, ndivMajor, rnum, rden);
00480 if (rnum == 0) form = Form("%d",rnum);
00481 if (rnum == 1 && rden == 1) form = Form("#pi");
00482 if (rnum == 1 && rden != 1) form = Form("#frac{#pi}{%d}",rden);
00483 if (rnum != 1 && rden == 1 && i !=0) form = Form("%d#pi",rnum);
00484 if (rnum != 1 && rden != 1) form = Form("#frac{%d#pi}{%d}",rnum,rden);
00485 if(theta >= 3*TMath::Pi()/12.0 && theta < 2*TMath::Pi()/3.0) corr=0.04;
00486 textangular->SetTextAlign(FindAlign(theta));
00487 textangular->PaintLatex(costhetas,corr+sinthetas,0,
00488 GetPolarLabelSize(),form);
00489 } else {
00490
00491 form = Form("%5.3g",txtval);
00492 axis.LabelsLimits(form,first,last);
00493 TString s = Form("%s",form);
00494 if (first != 0) s.Remove(0, first);
00495 if(theta >= 3*TMath::Pi()/12.0 && theta < 2*TMath::Pi()/3.0) corr=0.04;
00496 textangular->SetTextAlign(FindAlign(theta));
00497 textangular->PaintLatex(costhetas,
00498 corr+sinthetas,0,GetPolarLabelSize(),s);
00499 }
00500 } else if (fPolarLabels) {
00501
00502 textangular->SetTextAlign(FindAlign(theta));
00503 textangular->PaintText(costhetas,sinthetas,fPolarLabels[i]);
00504 }
00505 }
00506 TAttLine::Modify();
00507
00508 Bool_t issettickpolar = gPad->GetTicky();
00509
00510 if (issettickpolar) {
00511 if (theta != 0 && theta !=TMath::Pi()) {
00512 gPad->PaintLine((sintheta-GetTickpolarSize())/tantheta,sintheta-GetTickpolarSize(),
00513 (sintheta+GetTickpolarSize())/tantheta,sintheta+GetTickpolarSize());
00514 }
00515 if (theta == 0 || theta ==TMath::Pi()) {
00516 gPad->PaintLine(1-GetTickpolarSize(),0,1+GetTickpolarSize(),0);
00517 gPad->PaintLine(-1+GetTickpolarSize(),0,-1-GetTickpolarSize(),0);
00518 }
00519 }
00520 TAttLine::SetLineStyle(1);
00521 TAttLine::Modify();
00522 gPad->PaintLine(0.,0.,costheta,sintheta);
00523 delete textangular;
00524
00525 Int_t oldLineStyle = GetLineStyle();
00526 TAttLine::SetLineStyle(2);
00527 TAttLine::Modify();
00528 for (j=1; j<ndivMinor; j++) {
00529 Double_t thetamin = theta+j*2*TMath::Pi()/(ndivMajor*ndivMinor);
00530 gPad->PaintLine(0.,0.,TMath::Cos(thetamin),TMath::Sin(thetamin));
00531 }
00532 TAttLine::SetLineStyle(oldLineStyle);
00533 TAttLine::Modify();
00534 }
00535 } else {
00536 Int_t big = (Int_t)fRwtmax;
00537 Int_t test= 1;
00538 while (big >= 10) {
00539 big = big/10;
00540 test++;
00541 }
00542 for (i=1; i<=test; i++) {
00543 Double_t txtval = pow((double)10,(double)(i-1));
00544 Double_t theta = (i-1)*2*TMath::Pi()/(double)(test);
00545 Double_t costheta = TMath::Cos(theta);
00546 Double_t sintheta = TMath::Sin(theta);
00547 Double_t tantheta = TMath::Tan(theta);
00548 Double_t costhetas = (1+fPolarOffset)*costheta;
00549 Double_t sinthetas = (1+fPolarOffset)*sintheta;
00550 Double_t corr = 0.01;
00551
00552 TLatex *textangular = new TLatex();
00553 textangular->SetTextColor(GetPolarColorLabel());
00554 textangular->SetTextFont(GetPolarLabelFont());
00555
00556 const char* form = (char *)" ";
00557 TGaxis axis;
00558
00559 if (TestBit(TGraphPolargram::kLabelOrtho)) {
00560 if(fPolarLabels==NULL && optionLabels){
00561
00562 form = Form("%5.3g",txtval);
00563 axis.LabelsLimits(form,first,last);
00564 TString s = Form("%s",form);
00565 if (first != 0) s.Remove(0, first);
00566 textangular->SetTextAlign(FindAlign(theta));
00567 textangular->PaintLatex(costhetas,
00568 sinthetas, FindTextAngle(theta), GetPolarLabelSize(), s);
00569 }
00570 else if (fPolarLabels){
00571
00572 textangular->SetTextAlign(FindAlign(theta));
00573 textangular->PaintText(costhetas,sinthetas,fPolarLabels[i]);
00574 }
00575
00576 } else {
00577 if(fPolarLabels==NULL && optionLabels){
00578
00579 form = Form("%5.3g",txtval);
00580 axis.LabelsLimits(form,first,last);
00581 TString s = Form("%s",form);
00582 if (first != 0) s.Remove(0, first);
00583 if(theta >= 3*TMath::Pi()/12.0 && theta < 2*TMath::Pi()/3.0) corr=0.04;
00584 textangular->SetTextAlign(FindAlign(theta));
00585 textangular->PaintLatex(costhetas,
00586 corr+sinthetas,0,GetPolarLabelSize(),s);
00587 } else if (fPolarLabels){
00588
00589 textangular->SetTextAlign(FindAlign(theta));
00590 textangular->PaintText(costhetas,sinthetas,fPolarLabels[i]);
00591 }
00592 }
00593
00594 TAttLine::Modify();
00595
00596 Bool_t issettickpolar = gPad->GetTicky();
00597 if (issettickpolar) {
00598 if (theta != 0 && theta !=TMath::Pi()) {
00599 gPad->PaintLine((sintheta-GetTickpolarSize())/tantheta,sintheta-GetTickpolarSize(),
00600 (sintheta+GetTickpolarSize())/tantheta,sintheta+GetTickpolarSize());
00601 }
00602 if (theta == 0 || theta ==TMath::Pi()) {
00603 gPad->PaintLine(1-GetTickpolarSize(),0,1+GetTickpolarSize(),0);
00604 gPad->PaintLine(-1+GetTickpolarSize(),0,-1-GetTickpolarSize(),0);
00605 }
00606 }
00607 TAttLine::SetLineStyle(1);
00608 TAttLine::Modify();
00609 gPad->PaintLine(0.,0.,costheta,sintheta);
00610 delete textangular;
00611
00612 Int_t oldLineStyle = GetLineStyle();
00613 TAttLine::SetLineStyle(2);
00614 TAttLine::Modify();
00615 Double_t a=0;
00616 Double_t b,c,d;
00617 b = TMath::Log(10)*test;
00618 d= 2*TMath::Pi()/(double)test;
00619 for (j=1; j<9; j++) {
00620 a=TMath::Log(j+1)-TMath::Log(j)+a;
00621 c=a/b*6.28+d*(i-1);
00622 gPad->PaintLine(0.,0.,TMath::Cos(c),TMath::Sin(c));
00623 }
00624 TAttLine::SetLineStyle(oldLineStyle);
00625 TAttLine::Modify();
00626 }
00627 }
00628 }
00629
00630
00631
00632 void TGraphPolargram::PaintRadialDivisions(Bool_t drawaxis)
00633 {
00634
00635
00636
00637 static char chopt[8] = "";
00638 Int_t i,j;
00639 Int_t ndiv = TMath::Abs(fNdivRad);
00640 Int_t ndivMajor = ndiv%100;
00641 Int_t ndivMinor = ndiv/100;
00642 Int_t ndivmajor = 0;
00643 Double_t frwrmin = 0., frwrmax = 0., binWidth = 0;
00644
00645 THLimitsFinder::Optimize(fRwrmin,fRwrmax,ndivMajor,frwrmin,
00646 frwrmax, ndivmajor,binWidth,"");
00647
00648 if (!gPad) return ;
00649 if (!gPad->GetLogx()) {
00650 gPad->RangeAxis(-1,-1,1,1);
00651 gPad->Range(-1.25,-1.25,1.25,1.25);
00652 Double_t umin = fRwrmin;
00653 Double_t umax = fRwrmax;
00654 Double_t rmajmin = (frwrmin-fRwrmin)/(fRwrmax-fRwrmin);
00655 Double_t rmajmax = (frwrmax-fRwrmin)/(fRwrmax-fRwrmin);
00656 Double_t dist = (rmajmax-rmajmin)/ndivmajor;
00657 Int_t ndivminor = 0;
00658
00659 chopt[0] = 0;
00660 strncat(chopt, "SDH", 3);
00661 if (fNdivRad < 0) strncat(chopt, "N",1);
00662 if(drawaxis){
00663
00664 TGaxis axis;
00665 axis.SetLabelSize(GetRadialLabelSize());
00666 axis.SetLabelColor(GetRadialColorLabel());
00667 axis.SetLabelFont(GetRadialLabelFont());
00668 axis.SetLabelOffset(GetRadialOffset());
00669 axis.PaintAxis(0, 0, TMath::Cos(GetAngle()), TMath::Sin(GetAngle()),
00670 umin, umax, ndiv, chopt, 0., kFALSE);
00671 }
00672
00673
00674
00675 PaintCircle(0.,0.,1,0.,360,0);
00676
00677 if (fNdivRad>0 ) {
00678 Double_t frwrmini = 0., frwrmaxi = 0., binWidth2 =0;
00679 THLimitsFinder::Optimize(frwrmin,frwrmin+binWidth,ndivMinor,frwrmini,
00680 frwrmaxi, ndivminor,binWidth2,"");
00681 Double_t dist2 = dist/(ndivminor);
00682
00683 for (i=1; i<=ndivmajor+2; i++) {
00684 TAttLine::SetLineStyle(1);
00685 TAttLine::Modify();
00686 PaintCircle(0.,0.,rmajmin,0.,360,0);
00687
00688
00689 TAttLine::SetLineStyle(2);
00690 TAttLine::Modify();
00691 for (j=1; j<ndivminor+1; j++) {
00692 if (rmajmin+j*dist2<=1) PaintCircle(0.,0.,rmajmin+j*dist2,0.,360,0);
00693 }
00694 rmajmin = (frwrmin-fRwrmin)/(fRwrmax-fRwrmin)+(i-1)*dist;
00695 }
00696
00697 } else {
00698
00699
00700 for (i=1; i<=ndivMajor; i++) {
00701 TAttLine::SetLineStyle(1);
00702 TAttLine::Modify();
00703 Double_t rmaj = i*1./ndivMajor;
00704 PaintCircle(0.,0.,rmaj,0.,360,0);
00705
00706
00707 for (j=1; j<ndivMinor; j++) {
00708 TAttLine::SetLineStyle(2);
00709 TAttLine::Modify();
00710 PaintCircle(0.,0.,rmaj- j*1./(ndivMajor*ndivMinor),0.,360,0);
00711 }
00712 }
00713 }
00714 } else {
00715
00716 Int_t big = (Int_t)fRwrmax;
00717 Int_t test= 0;
00718 while (big >= 10) {
00719 big = big/10;
00720 test++;
00721 }
00722 for (i=1; i<=test; i++) {
00723 TAttLine::SetLineStyle(1);
00724 TAttLine::Modify();
00725 Double_t ecart;
00726 ecart = ((double) i)/ ((double) test);
00727 PaintCircle(0.,0.,ecart,0,360,0);
00728 TAttLine::SetLineStyle(GetLineStyle());
00729 TAttLine::Modify();
00730 Double_t a=0;
00731 Double_t b,c,d;
00732 b = TMath::Log(10)*test;
00733 d = 1/(double)test;
00734 for (j=1; j<9; j++) {
00735 a = TMath::Log(j+1)-TMath::Log(j)+a;
00736 c = a/b+d*(i-1);
00737 PaintCircle(0,0.,c,0.,360,0);
00738 }
00739 }
00740 }
00741 TAttLine::SetLineStyle(1);
00742 TAttLine::Modify();
00743 }
00744
00745
00746
00747 void TGraphPolargram::ReduceFraction(Int_t num, Int_t den, Int_t &rnum, Int_t &rden)
00748 {
00749
00750
00751 Int_t a = 0;
00752 Int_t b = 0;
00753 Int_t i = 0;
00754 Int_t j = 0;
00755 a = den;
00756 b = num;
00757 if (b > a) {
00758 j = b;
00759 } else {
00760 j = a;
00761 }
00762 for (i=j; i > 1; i--) {
00763 if ((a % i == 0) && (b % i == 0)) {
00764 a = a/i;
00765 b = b/i;
00766 }
00767 }
00768 rden = a;
00769 rnum = b;
00770 }
00771
00772
00773
00774 void TGraphPolargram::SetAxisAngle(Double_t angle)
00775 {
00776
00777
00778 fAxisAngle = angle/180*TMath::Pi();
00779 }
00780
00781
00782
00783 void TGraphPolargram::SetNdivPolar(Int_t ndiv)
00784 {
00785
00786
00787
00788
00789
00790 if (ndiv > 0)
00791 fNdivPol = ndiv;
00792 if (gPad) gPad->Modified();
00793 }
00794
00795
00796
00797 void TGraphPolargram::SetNdivRadial(Int_t ndiv)
00798 {
00799
00800
00801
00802
00803
00804 fNdivRad = ndiv;
00805 if (gPad) gPad->Modified();
00806 }
00807
00808
00809
00810 void TGraphPolargram::SetPolarLabel(Int_t div, const TString & label)
00811 {
00812
00813 if(fPolarLabels == NULL)
00814 fPolarLabels = new TString[fNdivPol];
00815 fPolarLabels[div]=label;
00816 if (gPad) gPad->Modified();
00817 }
00818
00819
00820
00821 void TGraphPolargram::SetPolarLabelColor(Color_t tcolorangular )
00822 {
00823
00824
00825 fPolarLabelColor = tcolorangular;
00826 }
00827
00828
00829
00830 void TGraphPolargram::SetPolarLabelFont(Font_t tfontangular)
00831 {
00832
00833
00834
00835 fPolarLabelFont = tfontangular;
00836 }
00837
00838
00839
00840 void TGraphPolargram::SetPolarLabelSize(Double_t angularsize )
00841 {
00842
00843
00844 fPolarTextSize = angularsize;
00845 }
00846
00847
00848
00849 void TGraphPolargram::SetPolarOffset(Double_t angularOffset)
00850 {
00851
00852
00853 fPolarOffset = angularOffset;
00854 if (gPad) gPad->Modified();
00855 }
00856
00857
00858
00859 void TGraphPolargram::SetRadialLabelColor(Color_t tcolorradial )
00860 {
00861
00862
00863 fRadialLabelColor = tcolorradial;
00864 }
00865
00866
00867
00868 void TGraphPolargram::SetRadialLabelFont(Font_t tfontradial)
00869 {
00870
00871
00872 fRadialLabelFont = tfontradial;
00873 }
00874
00875
00876
00877 void TGraphPolargram::SetRadialLabelSize(Double_t radialsize )
00878 {
00879
00880
00881 fRadialTextSize = radialsize;
00882 }
00883
00884
00885
00886 void TGraphPolargram::SetRadialOffset(Double_t radialOffset)
00887 {
00888
00889
00890 fRadialOffset = radialOffset;
00891 if (gPad) gPad->Modified();
00892 }
00893
00894
00895
00896 void TGraphPolargram::SetRangePolar(Double_t tmin, Double_t tmax)
00897 {
00898
00899
00900
00901
00902 fDegree = kFALSE;
00903 fGrad = kFALSE;
00904 fRadian = kFALSE;
00905
00906 if (tmin < tmax) {
00907 fRwtmin = tmin;
00908 fRwtmax = tmax;
00909 }
00910 if (gPad) gPad->Modified();
00911 }
00912
00913
00914
00915 void TGraphPolargram::SetRangeRadial(Double_t rmin, Double_t rmax)
00916 {
00917
00918
00919
00920
00921 if (rmin < rmax) {
00922 fRwrmin = rmin;
00923 fRwrmax = rmax;
00924 }
00925 if (gPad) gPad->Modified();
00926 }
00927
00928
00929
00930 void TGraphPolargram::SetTickpolarSize(Double_t tickpolarsize)
00931 {
00932
00933
00934 fTickpolarSize = tickpolarsize;
00935 }
00936
00937
00938
00939 void TGraphPolargram::SetToDegree()
00940 {
00941
00942
00943 fDegree = kTRUE;
00944 fGrad = kFALSE;
00945 fRadian = kFALSE;
00946 ChangeRangePolar(0,360);
00947 }
00948
00949
00950
00951 void TGraphPolargram::SetToGrad()
00952 {
00953
00954
00955 fGrad = kTRUE;
00956 fRadian = kFALSE;
00957 fDegree = kFALSE;
00958 ChangeRangePolar(0,200);
00959 }
00960
00961
00962
00963 void TGraphPolargram::SetToRadian()
00964 {
00965
00966
00967 fRadian = kTRUE;
00968 fGrad = kFALSE;
00969 fDegree = kFALSE;
00970 ChangeRangePolar(0,2*TMath::Pi());
00971 }
00972
00973
00974
00975 void TGraphPolargram::SetTwoPi()
00976 {
00977
00978 SetRangePolar(0,2*TMath::Pi());
00979 }