{ "cells": [ { "cell_type": "markdown", "id": "a03a511d", "metadata": {}, "source": [ "# gr012_polar\n", "\n", "Since TGraphPolar is a TGraphErrors, it is painted with\n", "[TGraphPainter](https://root.cern/doc/master/classTGraphPainter.html) options.\n", "\n", "With GetPolargram we retrieve the polar axis to format it; see the\n", "[TGraphPolargram documentation](https://root.cern/doc/master/classTGraphPolargram.html)\n", "\n", "\n", "\n", "**Author:** Olivier Couet \n", "This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Tuesday, May 19, 2026 at 08:37 PM." ] }, { "cell_type": "markdown", "id": "03893c4e", "metadata": {}, "source": [ "Illustrates how to use TGraphPolar" ] }, { "cell_type": "code", "execution_count": null, "id": "2d391d10", "metadata": { "collapsed": false }, "outputs": [], "source": [ "TCanvas * CPol = new TCanvas(\"CPol\",\"TGraphPolar Examples\",1200,600);\n", "CPol->Divide(2,1);" ] }, { "cell_type": "markdown", "id": "2dfd0062", "metadata": {}, "source": [ "Left-side pad. Two graphs without errors" ] }, { "cell_type": "code", "execution_count": null, "id": "f71a019a", "metadata": { "collapsed": false }, "outputs": [], "source": [ "CPol->cd(1);\n", "Double_t xmin=0;\n", "Double_t xmax=TMath::Pi()*2;\n", "\n", "Double_t x[1000];\n", "Double_t y[1000];\n", "Double_t xval1[20];\n", "Double_t yval1[20];" ] }, { "cell_type": "markdown", "id": "5227ee23", "metadata": {}, "source": [ "Graph 1 to be drawn with line and fill" ] }, { "cell_type": "code", "execution_count": null, "id": "9d95f49f", "metadata": { "collapsed": false }, "outputs": [], "source": [ "TF1 * fplot = new TF1(\"fplot\",\"cos(2*x)*cos(20*x)\",xmin,xmax);\n", "for (Int_t ipt = 0; ipt < 1000; ipt++){\n", " x[ipt] = ipt*(xmax-xmin)/1000 + xmin;\n", " y[ipt] = fplot->Eval(x[ipt]);\n", "}\n", "TGraphPolar * grP = new TGraphPolar(1000,x,y);\n", "grP->SetLineColor(2);\n", "grP->SetLineWidth(2);\n", "grP->SetFillStyle(3012);\n", "grP->SetFillColor(2);\n", "grP->Draw(\"AFL\");" ] }, { "cell_type": "markdown", "id": "82a9d507", "metadata": {}, "source": [ "Graph 2 to be drawn superposed over graph 1, with curve and polymarker" ] }, { "cell_type": "code", "execution_count": null, "id": "e87abaa2", "metadata": { "collapsed": false }, "outputs": [], "source": [ "for (Int_t ipt = 0; ipt < 20; ipt++){\n", " xval1[ipt] = x[1000/20*ipt];\n", " yval1[ipt] = y[1000/20*ipt];\n", "}\n", "TGraphPolar * grP1 = new TGraphPolar(20,xval1,yval1);\n", "grP1->SetMarkerStyle(29);\n", "grP1->SetMarkerSize(2);\n", "grP1->SetMarkerColor(4);\n", "grP1->SetLineColor(4);\n", "grP1->Draw(\"CP\");" ] }, { "cell_type": "markdown", "id": "3bc32ed5", "metadata": {}, "source": [ "To format the polar axis, we retrieve the TGraphPolargram.\n", "First update the canvas, otherwise GetPolargram returns 0" ] }, { "cell_type": "code", "execution_count": null, "id": "f0282315", "metadata": { "collapsed": false }, "outputs": [], "source": [ "CPol->Update();\n", "if (grP1->GetPolargram()) {\n", " grP1->GetPolargram()->SetTextColor(8);\n", " grP1->GetPolargram()->SetRangePolar(-TMath::Pi(),TMath::Pi());\n", " grP1->GetPolargram()->SetNdivPolar(703);\n", " grP1->GetPolargram()->SetToRadian(); // tell ROOT that the x and xval1 are in radians\n", "}" ] }, { "cell_type": "markdown", "id": "b9048ade", "metadata": {}, "source": [ "Right-side pad. One graph with errors" ] }, { "cell_type": "code", "execution_count": null, "id": "e4c4fc86", "metadata": { "collapsed": false }, "outputs": [], "source": [ "CPol->cd(2);\n", "Double_t x2[30];\n", "Double_t y2[30];\n", "Double_t ex[30];\n", "Double_t ey[30];\n", "for (Int_t ipt = 0; ipt < 30; ipt++){\n", " x2[ipt] = x[1000/30*ipt];\n", " y2[ipt] = 1.2 + 0.4*sin(TMath::Pi()*2*ipt/30);\n", " ex[ipt] = 0.2+0.1*cos(2*TMath::Pi()/30*ipt);\n", " ey[ipt] = 0.2;\n", "}" ] }, { "cell_type": "markdown", "id": "0e1522b5", "metadata": {}, "source": [ "Grah to be drawn with polymarker and errors" ] }, { "cell_type": "code", "execution_count": null, "id": "dce8a94e", "metadata": { "collapsed": false }, "outputs": [], "source": [ "TGraphPolar * grPE = new TGraphPolar(30,x2,y2,ex,ey);\n", "grPE->SetMarkerStyle(22);\n", "grPE->SetMarkerSize(1.5);\n", "grPE->SetMarkerColor(5);\n", "grPE->SetLineColor(6);\n", "grPE->SetLineWidth(2);\n", "grPE->Draw(\"EP\");" ] }, { "cell_type": "markdown", "id": "11056d07", "metadata": {}, "source": [ "To format the polar axis, we retrieve the TGraphPolargram.\n", "First update the canvas, otherwise GetPolargram returns 0" ] }, { "cell_type": "code", "execution_count": null, "id": "dad844e9", "metadata": { "collapsed": false }, "outputs": [], "source": [ "CPol->Update();\n", "if (grPE->GetPolargram()) {\n", " grPE->GetPolargram()->SetTextSize(0.03);\n", " grPE->GetPolargram()->SetTwoPi();\n", " grPE->GetPolargram()->SetToRadian(); // tell ROOT that the x2 values are in radians\n", "}" ] }, { "cell_type": "markdown", "id": "87fcd096", "metadata": {}, "source": [ "Draw all canvases " ] }, { "cell_type": "code", "execution_count": null, "id": "63847db4", "metadata": { "collapsed": false }, "outputs": [], "source": [ "gROOT->GetListOfCanvases()->Draw()" ] } ], "metadata": { "kernelspec": { "display_name": "ROOT C++", "language": "c++", "name": "root" }, "language_info": { "codemirror_mode": "text/x-c++src", "file_extension": ".C", "mimetype": " text/x-c++src", "name": "c++" } }, "nbformat": 4, "nbformat_minor": 5 }