/usr/share/doc/mathgl/mathgl_en/mathgl_en_18.html is in mathgl-doc-en 2.1.3.1-4ubuntu3.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 | <!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html401/loose.dtd">
<html>
<!-- This manual is for MathGL (version 2.1.2), a collection of classes and routines for scientific plotting. Please report any errors in this manual to mathgl.abalakin@gmail.org.
Copyright C 2008-2012 Alexey A. Balakin.
Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.2
or any later version published by the Free Software Foundation;
with no Invariant Sections, no Front-Cover Texts, and no Back-Cover
Texts. A copy of the license is included in the section entitled "GNU
Free Documentation License."
-->
<!-- Created on December 22, 2013 by texi2html 1.82
texi2html was written by:
Lionel Cons <Lionel.Cons@cern.ch> (original author)
Karl Berry <karl@freefriends.org>
Olaf Bachmann <obachman@mathematik.uni-kl.de>
and many others.
Maintained by: Many creative people.
Send bugs and suggestions to <texi2html-bug@nongnu.org>
-->
<head>
<title>MathGL 2.1.2: 2.9 Hints</title>
<meta name="description" content="MathGL 2.1.2: 2.9 Hints">
<meta name="keywords" content="MathGL 2.1.2: 2.9 Hints">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="texi2html 1.82">
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.smallquotation {font-size: smaller}
pre.display {font-family: serif}
pre.format {font-family: serif}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: serif; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: serif; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.roman {font-family:serif; font-weight:normal;}
span.sansserif {font-family:sans-serif; font-weight:normal;}
ul.toc {list-style: none}
-->
</style>
</head>
<body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
<a name="Hints"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="mathgl_en_17.html#Dew-sample" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#g_t_0060_0060Compound_0027_0027-graphics" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Hints-1"></a>
<h2 class="section">2.9 Hints</h2>
<p>In this section I’ve included some small hints and advices for the improving of the quality of plots and for the demonstration of some non-trivial features of MathGL library. In contrast to previous examples I showed mostly the idea but not the whole drawing function.
</p>
<table class="menu" border="0" cellspacing="0">
<tr><td align="left" valign="top"><a href="#g_t_0060_0060Compound_0027_0027-graphics">2.9.1 “Compound” graphics</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Transparency-and-lighting">2.9.2 Transparency and lighting</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Types-of-transparency">2.9.3 Types of transparency</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Axis-projection">2.9.4 Axis projection</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Adding-fog">2.9.5 Adding fog</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Several-light-sources">2.9.6 Several light sources</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Using-primitives">2.9.7 Using primitives</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#STFA-sample">2.9.8 STFA sample</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Mapping-visualization">2.9.9 Mapping visualization</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Making-histogram">2.9.10 Making histogram</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Nonlinear-fitting-sample">2.9.11 Nonlinear fitting hints</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#PDE-solving-hints">2.9.12 PDE solving hints</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#MGL-parser-using">2.9.13 MGL parser using</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Using-options">2.9.14 Using options</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#g_t_0060_0060Templates_0027_0027">2.9.15 “Templates”</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Stereo-image">2.9.16 Stereo image</a></td><td> </td><td align="left" valign="top">
</td></tr>
<tr><td align="left" valign="top"><a href="#Reduce-memory-usage">2.9.17 Reduce memory usage</a></td><td> </td><td align="left" valign="top">
</td></tr>
</table>
<hr size="6">
<a name="g_t_0060_0060Compound_0027_0027-graphics"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Hints" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Transparency-and-lighting" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="g_t_0060_0060Compound_0027_0027-graphics-1"></a>
<h3 class="subsection">2.9.1 “Compound” graphics</h3>
<p>As I noted above, MathGL functions (except the special one, like Clf()) do not erase the previous plotting but just add the new one. It allows one to draw “compound” plots easily. For example, popular Matlab command <code>surfc</code> can be emulated in MathGL by 2 calls:
</p><pre class="verbatim"> Surf(a);
Cont(a, "_"); // draw contours at bottom
</pre><p>Here <var>a</var> is 2-dimensional data for the plotting, <code>-1</code> is the value of z-coordinate at which the contour should be plotted (at the bottom in this example). Analogously, one can draw density plot instead of contour lines and so on.
</p>
<p>Another nice plot is contour lines plotted directly on the surface:
</p><pre class="verbatim"> Light(true); // switch on light for the surface
Surf(a, "BbcyrR"); // select 'jet' colormap for the surface
Cont(a, "y"); // and yellow color for contours
</pre><p>The possible difficulties arise in black&white case, when the color of the surface can be close to the color of a contour line. In that case I may suggest the following code:
</p><pre class="verbatim"> Light(true); // switch on light for the surface
Surf(a, "kw"); // select 'gray' colormap for the surface
CAxis(-1,0); // first draw for darker surface colors
Cont(a, "w"); // white contours
CAxis(0,1); // now draw for brighter surface colors
Cont(a, "k"); // black contours
CAxis(-1,1); // return color range to original state
</pre><p>The idea is to divide the color range on 2 parts (dark and bright) and to select the contrasting color for contour lines for each of part.
</p>
<p>Similarly, one can plot flow thread over density plot of vector field amplitude (this is another amusing plot from Matlab) and so on. The list of compound graphics can be prolonged but I hope that the general idea is clear.
</p>
<p>Just for illustration I put here following sample code:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a,b,d; mgls_prepare2v(&a,&b); d = a;
for(int i=0;i<a.nx*a.ny;i++) d.a[i] = hypot(a.a[i],b.a[i]);
mglData c; mgls_prepare3d(&c);
mglData v(10); v.Fill(-0.5,1);
gr->SubPlot(2,2,1,""); gr->Title("Flow + Dens");
gr->Flow(a,b,"br"); gr->Dens(d,"BbcyrR"); gr->Box();
gr->SubPlot(2,2,0); gr->Title("Surf + Cont"); gr->Rotate(50,60);
gr->Light(true); gr->Surf(a); gr->Cont(a,"y"); gr->Box();
gr->SubPlot(2,2,2); gr->Title("Mesh + Cont"); gr->Rotate(50,60);
gr->Box(); gr->Mesh(a); gr->Cont(a,"_");
gr->SubPlot(2,2,3); gr->Title("Surf3 + ContF3");gr->Rotate(50,60);
gr->Box(); gr->ContF3(v,c,"z",0); gr->ContF3(v,c,"x"); gr->ContF3(v,c);
gr->SetCutBox(mglPoint(0,-1,-1), mglPoint(1,0,1.1));
gr->ContF3(v,c,"z",c.nz-1); gr->Surf3(-0.5,c);
return 0;
}
</pre>
<img src="../png/combined.png" alt="Example of “combined” plots">
<hr size="6">
<a name="Transparency-and-lighting"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#g_t_0060_0060Compound_0027_0027-graphics" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Types-of-transparency" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Transparency-and-lighting-1"></a>
<h3 class="subsection">2.9.2 Transparency and lighting</h3>
<p>Here I want to show how transparency and lighting both and separately change the look of a surface. So, there is code and picture for that:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a; mgls_prepare2d(&a);
gr->SubPlot(2,2,0); gr->Title("default"); gr->Rotate(50,60);
gr->Box(); gr->Surf(a);
gr->SubPlot(2,2,1); gr->Title("light on"); gr->Rotate(50,60);
gr->Box(); gr->Light(true); gr->Surf(a);
gr->SubPlot(2,2,3); gr->Title("alpha on; light on"); gr->Rotate(50,60);
gr->Box(); gr->Alpha(true); gr->Surf(a);
gr->SubPlot(2,2,2); gr->Title("alpha on"); gr->Rotate(50,60);
gr->Box(); gr->Light(false); gr->Surf(a);
return 0;
}
</pre>
<img src="../png/alpha.png" alt="Example of transparency and lightings">
<hr size="6">
<a name="Types-of-transparency"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Transparency-and-lighting" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Axis-projection" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Types-of-transparency-1"></a>
<h3 class="subsection">2.9.3 Types of transparency</h3>
<p>MathGL library has advanced features for setting and handling the surface transparency. The simplest way to add transparency is the using of function <a href="mathgl_en_31.html#alpha">alpha</a>. As a result, all further surfaces (and isosurfaces, density plots and so on) become transparent. However, their look can be additionally improved.
</p>
<p>The value of transparency can be different from surface to surface. To do it just use <code>SetAlphaDef</code> before the drawing of the surface, or use option <code>alpha</code> (see <a href="mathgl_en_27.html#Command-options">Command options</a>). If its value is close to 0 then the surface becomes more and more transparent. Contrary, if its value is close to 1 then the surface becomes practically non-transparent.
</p>
<p>Also you can change the way how the light goes through overlapped surfaces. The function <code>SetTranspType</code> defines it. By default the usual transparency is used (‘<samp>0</samp>’) – surfaces below is less visible than the upper ones. A “glass-like” transparency (‘<samp>1</samp>’) has a different look – each surface just decreases the background light (the surfaces are commutable in this case).
</p>
<p>A “neon-like” transparency (‘<samp>2</samp>’) has more interesting look. In this case a surface is the light source (like a lamp on the dark background) and just adds some intensity to the color. At this, the library sets automatically the black color for the background and changes the default line color to white.
</p>
<p>As example I shall show several plots for different types of transparency. The code is the same except the values of <code>SetTranspType</code> function:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
gr->Alpha(true); gr->Light(true);
mglData a; mgls_prepare2d(&a);
gr->SetTranspType(0); gr->Clf();
gr->SubPlot(2,2,0); gr->Rotate(50,60); gr->Surf(a); gr->Box();
gr->SubPlot(2,2,1); gr->Rotate(50,60); gr->Dens(a); gr->Box();
gr->SubPlot(2,2,2); gr->Rotate(50,60); gr->Cont(a); gr->Box();
gr->SubPlot(2,2,3); gr->Rotate(50,60); gr->Axial(a); gr->Box();
return 0;
}
</pre>
<img src="../png/type0.png" alt="Example of SetTranspType(0).">
<img src="../png/type1.png" alt="Example of SetTranspType(1).">
<img src="../png/type2.png" alt="Example of SetTranspType(2).">
<hr size="6">
<a name="Axis-projection"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Types-of-transparency" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Adding-fog" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Axis-projection-1"></a>
<h3 class="subsection">2.9.4 Axis projection</h3>
<p>You can easily make 3D plot and draw its x-,y-,z-projections (like in CAD) by using <a href="mathgl_en_32.html#ternary">ternary</a> function with arguments: 4 for Cartesian, 5 for Ternary and 6 for Quaternary coordinates. The sample code is:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
gr->SetRanges(0,1,0,1,0,1);
mglData x(50),y(50),z(50),rx(10),ry(10), a(20,30);
a.Modify("30*x*y*(1-x-y)^2*(x+y<1)");
x.Modify("0.25*(1+cos(2*pi*x))");
y.Modify("0.25*(1+sin(2*pi*x))");
rx.Modify("rnd"); ry.Modify("(1-v)*rnd",rx);
z.Modify("x");
gr->Title("Projection sample");
gr->Ternary(4);
gr->Rotate(50,60); gr->Light(true);
gr->Plot(x,y,z,"r2"); gr->Surf(a,"#");
gr->Axis(); gr->Grid(); gr->Box();
gr->Label('x',"X",1); gr->Label('y',"Y",1); gr->Label('z',"Z",1);
}
</pre>
<img src="../png/projection.png" alt="Example of axis projections">
<img src="../png/projection5.png" alt="Example of ternary axis projections">
<hr size="6">
<a name="Adding-fog"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Axis-projection" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Several-light-sources" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Adding-fog-1"></a>
<h3 class="subsection">2.9.5 Adding fog</h3>
<p>MathGL can add a fog to the image. Its switching on is rather simple – just use <a href="mathgl_en_31.html#fog">fog</a> function. There is the only feature – fog is applied for whole image. Not to particular subplot. The sample code is:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a; mgls_prepare2d(&a);
gr->Title("Fog sample");
gr->Light(true); gr->Rotate(50,60); gr->Fog(1); gr->Box();
gr->Surf(a);
return 0;
}
</pre>
<img src="../png/fog.png" alt="Example of Fog().">
<hr size="6">
<a name="Several-light-sources"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Adding-fog" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Using-primitives" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Several-light-sources-1"></a>
<h3 class="subsection">2.9.6 Several light sources</h3>
<p>In contrast to the most of other programs, MathGL supports several (up to 10) light sources. Moreover, the color each of them can be different: white (this is usual), yellow, red, cyan, green and so on. The use of several light sources may be interesting for the highlighting of some peculiarities of the plot or just to make an amusing picture. Note, each light source can be switched on/off individually. The sample code is:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a; mgls_prepare2d(&a);
gr->Title("Several light sources");
gr->Rotate(50,60); gr->Light(true);
gr->AddLight(1,mglPoint(0,1,0),'c');
gr->AddLight(2,mglPoint(1,0,0),'y');
gr->AddLight(3,mglPoint(0,-1,0),'m');
gr->Box(); gr->Surf(a,"h");
return 0;
}
</pre>
<img src="../png/several_light.png" alt="Example of several light sources.">
<hr size="6">
<a name="Using-primitives"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Several-light-sources" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#STFA-sample" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Using-primitives-1"></a>
<h3 class="subsection">2.9.7 Using primitives</h3>
<p>MathGL provide a set of functions for drawing primitives (see <a href="mathgl_en_35.html#Primitives">Primitives</a>). Primitives are low level object, which used by most of plotting functions. Picture below demonstrate some of commonly used primitives.
</p>
<img src="../png/primitives.png" alt="Primitives in MathGL.">
<p>Generally, you can create arbitrary new kind of plot using primitives. For example, MathGL don’t provide any special functions for drawing molecules. However, you can do it using only one type of primitives <a href="mathgl_en_35.html#drop">drop</a>. The sample code is:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
gr->Alpha(true); gr->Light(true);
gr->SubPlot(2,2,0,""); gr->Title("Methane, CH_4");
gr->StartGroup("Methane");
gr->Rotate(60,120);
gr->Sphere(mglPoint(0,0,0),0.25,"k");
gr->Drop(mglPoint(0,0,0),mglPoint(0,0,1),0.35,"h",1,2);
gr->Sphere(mglPoint(0,0,0.7),0.25,"g");
gr->Drop(mglPoint(0,0,0),mglPoint(-0.94,0,-0.33),0.35,"h",1,2);
gr->Sphere(mglPoint(-0.66,0,-0.23),0.25,"g");
gr->Drop(mglPoint(0,0,0),mglPoint(0.47,0.82,-0.33),0.35,"h",1,2);
gr->Sphere(mglPoint(0.33,0.57,-0.23),0.25,"g");
gr->Drop(mglPoint(0,0,0),mglPoint(0.47,-0.82,-0.33),0.35,"h",1,2);
gr->Sphere(mglPoint(0.33,-0.57,-0.23),0.25,"g");
gr->EndGroup();
gr->SubPlot(2,2,1,""); gr->Title("Water, H_{2}O");
gr->StartGroup("Water");
gr->Rotate(60,100);
gr->StartGroup("Water_O");
gr->Sphere(mglPoint(0,0,0),0.25,"r");
gr->EndGroup();
gr->StartGroup("Water_Bond_1");
gr->Drop(mglPoint(0,0,0),mglPoint(0.3,0.5,0),0.3,"m",1,2);
gr->EndGroup();
gr->StartGroup("Water_H_1");
gr->Sphere(mglPoint(0.3,0.5,0),0.25,"g");
gr->EndGroup();
gr->StartGroup("Water_Bond_2");
gr->Drop(mglPoint(0,0,0),mglPoint(0.3,-0.5,0),0.3,"m",1,2);
gr->EndGroup();
gr->StartGroup("Water_H_2");
gr->Sphere(mglPoint(0.3,-0.5,0),0.25,"g");
gr->EndGroup();
gr->EndGroup();
gr->SubPlot(2,2,2,""); gr->Title("Oxygen, O_2");
gr->StartGroup("Oxygen");
gr->Rotate(60,120);
gr->Drop(mglPoint(0,0.5,0),mglPoint(0,-0.3,0),0.3,"m",1,2);
gr->Sphere(mglPoint(0,0.5,0),0.25,"r");
gr->Drop(mglPoint(0,-0.5,0),mglPoint(0,0.3,0),0.3,"m",1,2);
gr->Sphere(mglPoint(0,-0.5,0),0.25,"r");
gr->EndGroup();
gr->SubPlot(2,2,3,""); gr->Title("Ammonia, NH_3");
gr->StartGroup("Ammonia");
gr->Rotate(60,120);
gr->Sphere(mglPoint(0,0,0),0.25,"b");
gr->Drop(mglPoint(0,0,0),mglPoint(0.33,0.57,0),0.32,"n",1,2);
gr->Sphere(mglPoint(0.33,0.57,0),0.25,"g");
gr->Drop(mglPoint(0,0,0),mglPoint(0.33,-0.57,0),0.32,"n",1,2);
gr->Sphere(mglPoint(0.33,-0.57,0),0.25,"g");
gr->Drop(mglPoint(0,0,0),mglPoint(-0.65,0,0),0.32,"n",1,2);
gr->Sphere(mglPoint(-0.65,0,0),0.25,"g");
gr->EndGroup();
return 0;
}
</pre>
<img src="../png/molecule.png" alt="Example of molecules drawing.">
<p>Moreover, some of special plots can be more easily produced by primitives rather than by specialized function. For example, Venn diagram can be produced by <code>Error</code> plot:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
double xx[3]={-0.3,0,0.3}, yy[3]={0.3,-0.3,0.3}, ee[3]={0.7,0.7,0.7};
mglData x(3,xx), y(3,yy), e(3,ee);
gr->Title("Venn-like diagram"); gr->Alpha(true);
gr->Error(x,y,e,e,"!rgb@#o");
return 0;
}
</pre><p>You see that you have to specify and fill 3 data arrays. The same picture can be produced by just 3 calls of <a href="mathgl_en_35.html#circle">circle</a> function:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
gr->Title("Venn-like diagram"); gr->Alpha(true);
gr->Circle(mglPoint(-0.3,0.3),0.7,"rr@");
gr->Circle(mglPoint(0,-0.3),0.7,"gg@");
gr->Circle(mglPoint( 0.3,0.3),0.7,"bb@");
return 0;
}
</pre><p>Of course, the first variant is more suitable if you need to plot a lot of circles. But for few ones the usage of primitives looks easy.
</p>
<img src="../png/venn.png" alt="Example of Venn diagram.">
<hr size="6">
<a name="STFA-sample"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Using-primitives" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Mapping-visualization" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="STFA-sample-1"></a>
<h3 class="subsection">2.9.8 STFA sample</h3>
<p>Short-time Fourier Analysis (<a href="mathgl_en_42.html#stfa">stfa</a>) is one of informative method for analyzing long rapidly oscillating 1D data arrays. It is used to determine the sinusoidal frequency and phase content of local sections of a signal as it changes over time.
</p>
<p>MathGL can find and draw STFA result. Just to show this feature I give following sample. Initial data arrays is 1D arrays with step-like frequency. Exactly this you can see at bottom on the STFA plot. The sample code is:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a(2000), b(2000);
gr->Fill(a,"cos(50*pi*x)*(x<-.5)+cos(100*pi*x)*(x<0)*(x>-.5)+\
cos(200*pi*x)*(x<.5)*(x>0)+cos(400*pi*x)*(x>.5)");
gr->SubPlot(1, 2, 0,"<_"); gr->Title("Initial signal");
gr->Plot(a);
gr->Axis();
gr->Label('x', "\\i t");
gr->SubPlot(1, 2, 1,"<_"); gr->Title("STFA plot");
gr->STFA(a, b, 64);
gr->Axis();
gr->Label('x', "\\i t");
gr->Label('y', "\\omega", 0);
return 0;
}
</pre>
<img src="../png/stfa.png" alt="Example of STFA().">
<hr size="6">
<a name="Mapping-visualization"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#STFA-sample" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Making-histogram" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Mapping-visualization-1"></a>
<h3 class="subsection">2.9.9 Mapping visualization</h3>
<p>Sometime ago I worked with mapping and have a question about its visualization. Let me remember you that mapping is some transformation rule for one set of number to another one. The 1d mapping is just an ordinary function – it takes a number and transforms it to another one. The 2d mapping (which I used) is a pair of functions which take 2 numbers and transform them to another 2 ones. Except general plots (like <a href="mathgl_en_42.html#surfc">surfc</a>, <a href="mathgl_en_42.html#surfa">surfa</a>) there is a special plot – Arnold diagram. It shows the area which is the result of mapping of some initial area (usually square).
</p>
<p>I tried to make such plot in <a href="mathgl_en_42.html#map">map</a>. It shows the set of points or set of faces, which final position is the result of mapping. At this, the color gives information about their initial position and the height describes Jacobian value of the transformation. Unfortunately, it looks good only for the simplest mapping but for the real multivalent quasi-chaotic mapping it produces a confusion. So, use it if you like :).
</p>
<p>The sample code for mapping visualization is:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a(50, 40), b(50, 40);
gr->Puts(mglPoint(0, 0), "\\to", ":C", -1.4);
gr->SetRanges(-1,1,-1,1,-2,2);
gr->SubPlot(2, 1, 0);
gr->Fill(a,"x"); gr->Fill(b,"y");
gr->Puts(mglPoint(0, 1.1), "\\{x, y\\}", ":C", -2); gr->Box();
gr->Map(a, b, "brgk");
gr->SubPlot(2, 1, 1);
gr->Fill(a,"(x^3+y^3)/2"); gr->Fill(b,"(x-y)/2");
gr->Puts(mglPoint(0, 1.1), "\\{\\frac{x^3+y^3}{2}, \\frac{x-y}{2}\\}", ":C", -2);
gr->Box();
gr->Map(a, b, "brgk");
return 0;
}
</pre>
<img src="../png/map.png" alt="Example of Map().">
<hr size="6">
<a name="Making-histogram"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Mapping-visualization" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Nonlinear-fitting-sample" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Making-histogram-1"></a>
<h3 class="subsection">2.9.10 Making histogram</h3>
<p>Using the <a href="mathgl_en_57.html#hist">hist</a> function(s) for making regular distributions is one of useful fast methods to process and plot irregular data. <code>Hist</code> can be used to find some momentum of set of points by specifying weight function. It is possible to create not only 1D distributions but also 2D and 3D ones. Below I place the simplest sample code which demonstrate <a href="mathgl_en_57.html#hist">hist</a> usage:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData x(10000), y(10000), z(10000); gr->Fill(x,"2*rnd-1");
gr->Fill(y,"2*rnd-1"); gr->Fill(z,"exp(-6*(v^2+w^2))",x,y);
mglData xx=gr->Hist(x,z), yy=gr->Hist(y,z); xx.Norm(0,1);
yy.Norm(0,1);
gr->MultiPlot(3,3,3,2,2,""); gr->SetRanges(-1,1,-1,1,0,1);
gr->Box(); gr->Dots(x,y,z,"wyrRk");
gr->MultiPlot(3,3,0,2,1,""); gr->SetRanges(-1,1,0,1);
gr->Box(); gr->Bars(xx);
gr->MultiPlot(3,3,5,1,2,""); gr->SetRanges(0,1,-1,1);
gr->Box(); gr->Barh(yy);
gr->SubPlot(3,3,2);
gr->Puts(mglPoint(0.5,0.5),"Hist and\nMultiPlot\nsample","a",-6);
return 0;
}
</pre>
<img src="../png/hist.png" alt="Example of Hist().">
<hr size="6">
<a name="Nonlinear-fitting-sample"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Making-histogram" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#PDE-solving-hints" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Nonlinear-fitting-hints"></a>
<h3 class="subsection">2.9.11 Nonlinear fitting hints</h3>
<p>Nonlinear fitting is rather simple. All that you need is the data to fit, the approximation formula and the list of coefficients to fit (better with its initial guess values). Let me demonstrate it on the following simple example. First, let us use sin function with some random noise:
</p><pre class="verbatim"> mglData rnd(100), in(100); //data to be fitted and ideal data
gr->Fill(rnd,"0.4*rnd+0.1+sin(2*pi*x)");
gr->Fill(in,"0.3+sin(2*pi*x)");
</pre><p>and plot it to see that data we will fit
</p><pre class="verbatim"> gr->Title("Fitting sample");
gr->SetRange('y',-2,2); gr->Box(); gr->Plot(rnd, ". ");
gr->Axis(); gr->Plot(in, "b");
gr->Puts(mglPoint(0, 2.2), "initial: y = 0.3+sin(2\\pi x)", "b");
</pre>
<p>The next step is the fitting itself. For that let me specify an initial values <var>ini</var> for coefficients ‘<samp>abc</samp>’ and do the fitting for approximation formula ‘<samp>a+b*sin(c*x)</samp>’
</p><pre class="verbatim"> mreal ini[3] = {1,1,3};
mglData Ini(3,ini);
mglData res = gr->Fit(rnd, "a+b*sin(c*x)", "abc", Ini);
</pre><p>Now display it
</p><pre class="verbatim"> gr->Plot(res, "r");
gr->Puts(mglPoint(-0.9, -1.3), "fitted:", "r:L");
gr->PutsFit(mglPoint(0, -1.8), "y = ", "r");
</pre>
<p>NOTE! the fitting results may have strong dependence on initial values for coefficients due to algorithm features. The problem is that in general case there are several local "optimums" for coefficients and the program returns only first found one! There are no guaranties that it will be the best. Try for example to set <code>ini[3] = {0, 0, 0}</code> in the code above.
</p>
<p>The full sample code for nonlinear fitting is:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData rnd(100), in(100);
gr->Fill(rnd,"0.4*rnd+0.1+sin(2*pi*x)");
gr->Fill(in,"0.3+sin(2*pi*x)");
mreal ini[3] = {1,1,3};
mglData Ini(3,ini);
mglData res = gr->Fit(rnd, "a+b*sin(c*x)", "abc", Ini);
gr->Title("Fitting sample");
gr->SetRange('y',-2,2); gr->Box(); gr->Plot(rnd, ". ");
gr->Axis(); gr->Plot(res, "r"); gr->Plot(in, "b");
gr->Puts(mglPoint(-0.9, -1.3), "fitted:", "r:L");
gr->PutsFit(mglPoint(0, -1.8), "y = ", "r");
gr->Puts(mglPoint(0, 2.2), "initial: y = 0.3+sin(2\\pi x)", "b");
return 0;
}
</pre>
<img src="../png/fit.png" alt="Example of nonlinear fitting.">
<hr size="6">
<a name="PDE-solving-hints"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Nonlinear-fitting-sample" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#MGL-parser-using" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="PDE-solving-hints-1"></a>
<h3 class="subsection">2.9.12 PDE solving hints</h3>
<p>Solving of Partial Differential Equations (PDE, including beam tracing) and ray tracing (or finding particle trajectory) are more or less common task. So, MathGL have several functions for that. There are <code>mglRay()</code> for ray tracing, <code>mglPDE()</code> for PDE solving, <code>mglQO2d()</code> for beam tracing in 2D case (see <a href="mathgl_en_62.html#Global-functions">Global functions</a>). Note, that these functions take “Hamiltonian” or equations as string values. And I don’t plan now to allow one to use user-defined functions. There are 2 reasons: the complexity of corresponding interface; and the basic nature of used methods which are good for samples but may not good for serious scientific calculations.
</p>
<p>The ray tracing can be done by <code>mglRay()</code> function. Really ray tracing equation is Hamiltonian equation for 3D space. So, the function can be also used for finding a particle trajectory (i.e. solve Hamiltonian ODE) for 1D, 2D or 3D cases. The function have a set of arguments. First of all, it is Hamiltonian which defined the media (or the equation) you are planning to use. The Hamiltonian is defined by string which may depend on coordinates ‘<samp>x</samp>’, ‘<samp>y</samp>’, ‘<samp>z</samp>’, time ‘<samp>t</samp>’ (for particle dynamics) and momentums ‘<samp>p</samp>’=<em>p_x</em>, ‘<samp>q</samp>’=<em>p_y</em>, ‘<samp>v</samp>’=<em>p_z</em>. Next, you have to define the initial conditions for coordinates and momentums at ‘<samp>t</samp>’=0 and set the integrations step (default is 0.1) and its duration (default is 10). The Runge-Kutta method of 4-th order is used for integration.
</p><pre class="verbatim"> const char *ham = "p^2+q^2-x-1+i*0.5*(y+x)*(y>-x)";
mglData r = mglRay(ham, mglPoint(-0.7, -1), mglPoint(0, 0.5), 0.02, 2);
</pre><p>This example calculate the reflection from linear layer (media with Hamiltonian ‘<samp>p^2+q^2-x-1</samp>’=<em>p_x^2+p_y^2-x-1</em>). This is parabolic curve. The resulting array have 7 columns which contain data for {x,y,z,p,q,v,t}.
</p>
<p>The solution of PDE is a bit more complicated. As previous you have to specify the equation as pseudo-differential operator <em>\hat H(x, \nabla)</em> which is called sometime as “Hamiltonian” (for example, in beam tracing). As previously, it is defined by string which may depend on coordinates ‘<samp>x</samp>’, ‘<samp>y</samp>’, ‘<samp>z</samp>’ (but not time!), momentums ‘<samp>p</samp>’=<em>(d/dx)/i k_0</em>, ‘<samp>q</samp>’=<em>(d/dy)/i k_0</em> and field amplitude ‘<samp>u</samp>’=<em>|u|</em>. The evolutionary coordinate is ‘<samp>z</samp>’ in all cases. So that, the equation look like <em>du/dz = ik_0 H(x,y,\hat p, \hat q, |u|)[u]</em>. Dependence on field amplitude ‘<samp>u</samp>’=<em>|u|</em> allows one to solve nonlinear problems too. For example, for nonlinear Shrodinger equation you may set <code>ham="p^2 + q^2 - u^2"</code>. Also you may specify imaginary part for wave absorption, like <code>ham = "p^2 + i*x*(x>0)"</code>, but only if dependence on variable ‘<samp>i</samp>’ is linear (i.e. <em>H = Hre+i*Him</em>).
</p>
<p>Next step is specifying the initial conditions at ‘<samp>z</samp>’=<code>Min.z</code>. The function need 2 arrays for real and for imaginary part. Note, that coordinates x,y,z are supposed to be in specified range [Min, Max]. So, the data arrays should have corresponding scales. Finally, you may set the integration step and parameter k0=<em>k_0</em>. Also keep in mind, that internally the 2 times large box is used (for suppressing numerical reflection from boundaries) and the equation should well defined even in this extended range.
</p>
<p>Final comment is concerning the possible form of pseudo-differential operator <em>H</em>. At this moment, simplified form of operator <em>H</em> is supported – all “mixed” terms (like ‘<samp>x*p</samp>’->x*d/dx) are excluded. For example, in 2D case this operator is effectively <em>H = f(p,z) + g(x,z,u)</em>. However commutable combinations (like ‘<samp>x*q</samp>’->x*d/dy) are allowed for 3D case.
</p>
<p>So, for example let solve the equation for beam deflected from linear layer and absorbed later. The operator will have the form ‘<samp>"p^2+q^2-x-1+i*0.5*(z+x)*(z>-x)"</samp>’ that correspond to equation <em>ik_0 \partial_z u + \Delta u + x \cdot u + i (x+z)/2 \cdot u = 0</em>. This is typical equation for Electron Cyclotron (EC) absorption in magnetized plasmas. For initial conditions let me select the beam with plane phase front <em>exp(-48*(x+0.7)^2)</em>. The corresponding code looks like this:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a,re(128),im(128);
gr->Fill(re,"exp(-48*(x+0.7)^2)");
a = gr->PDE("p^2+q^2-x-1+i*0.5*(z+x)*(z>-x)", re, im, 0.01, 30);
a.Transpose("yxz");
gr->SubPlot(1,1,0,"<_"); gr->Title("PDE solver");
gr->SetRange('c',0,1); gr->Dens(a,"wyrRk");
gr->Axis(); gr->Label('x', "\\i x"); gr->Label('y', "\\i z");
gr->FPlot("-x", "k|");
gr->Puts(mglPoint(0, 0.85), "absorption: (x+z)/2 for x+z>0");
gr->Puts(mglPoint(0,1.1),"Equation: ik_0\\partial_zu + \\Delta u + x\\cdot u + i \\frac{x+z}{2}\\cdot u = 0");
return 0;
}
</pre>
<img src="../png/pde.png" alt="Example of PDE solving.">
<p>The last example is example of beam tracing. Beam tracing equation is special kind of PDE equation written in coordinates accompanied to a ray. Generally this is the same parameters and limitation as for PDE solving but the coordinates are defined by the ray and by parameter of grid width <var>w</var> in direction transverse the ray. So, you don’t need to specify the range of coordinates. <strong>BUT</strong> there is limitation. The accompanied coordinates are well defined only for smooth enough rays, i.e. then the ray curvature <em>K</em> (which is defined as <em>1/K^2 = (|\ddot r|^2 |\dot r|^2 - (\ddot r, \dot r)^2)/|\dot r|^6</em>) is much large then the grid width: <em>K>>w</em>. So, you may receive incorrect results if this condition will be broken.
</p>
<p>You may use following code for obtaining the same solution as in previous example:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData r, xx, yy, a, im(128), re(128);
const char *ham = "p^2+q^2-x-1+i*0.5*(y+x)*(y>-x)";
r = mglRay(ham, mglPoint(-0.7, -1), mglPoint(0, 0.5), 0.02, 2);
gr->SubPlot(1,1,0,"<_"); gr->Title("Beam and ray tracing");
gr->Plot(r.SubData(0), r.SubData(1), "k");
gr->Axis(); gr->Label('x', "\\i x"); gr->Label('y', "\\i z");
// now start beam tracing
gr->Fill(re,"exp(-48*x^2)");
a = mglQO2d(ham, re, im, r, xx, yy, 1, 30);
gr->SetRange('c',0, 1);
gr->Dens(xx, yy, a, "wyrRk");
gr->FPlot("-x", "k|");
gr->Puts(mglPoint(0, 0.85), "absorption: (x+y)/2 for x+y>0");
gr->Puts(mglPoint(0.7, -0.05), "central ray");
return 0;
}
</pre>
<img src="../png/qo2d.png" alt="Example of beam tracing.">
<hr size="6">
<a name="MGL-parser-using"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#PDE-solving-hints" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Using-options" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="MGL-parser-using-1"></a>
<h3 class="subsection">2.9.13 MGL parser using</h3>
<p>Sometimes you may prefer to use MGL scripts in yours code. It is simpler (especially in comparison with C/Fortran interfaces) and provide faster way to plot the data with annotations, labels and so on. Class <code>mglParse</code> (see section <a href="mathgl_en_68.html#mglParse-class">mglParse class</a> parse MGL scripts in C++. It have also the corresponding interface for C/Fortran.
</p>
<p>The key function here is <code>mglParse::Parse()</code> (or <code>mgl_parse()</code> for C/Fortran) which execute one command per string. At this the detailed information about the possible errors or warnings is passed as function value. Or you may execute the whole script as long string with lines separated by ‘<samp>\n</samp>’. Functions <code>mglParse::Execute()</code> and <code>mgl_parse_text()</code> perform it. Also you may set the values of parameters ‘<samp>$0</samp>’...‘<samp>$9</samp>’ for the script by functions <code>mglParse::AddParam()</code> or <code>mgl_add_param()</code>, allow/disable picture resizing, check “once” status and so on. The usage is rather straight-forward.
</p>
<p>The only non-obvious thing is data transition between script and yours program. There are 2 stages: add or find variable; and set data to variable. In C++ you may use functions <code>mglParse::AddVar()</code> and <code>mglParse::FindVar()</code> which return pointer to <code>mglData</code>. In C/Fortran the corresponding functions are <code>mgl_add_var()</code>, <code>mgl_find_var()</code>. This data pointer is valid until next <code>Parse()</code> or <code>Execute()</code> call. Note, you <strong>must not delete or free</strong> the data obtained from these functions!
</p>
<p>So, some simple example at the end. Here I define a data array, create variable, put data into it and plot it. The C++ code looks like this:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
gr->Title("MGL parser sample");
mreal a[100]; // let a_i = sin(4*pi*x), x=0...1
for(int i=0;i<100;i++)a[i]=sin(4*M_PI*i/99);
mglParse *parser = new mglParse;
mglData *d = parser->AddVar("dat");
d->Set(a,100); // set data to variable
parser->Execute(gr, "plot dat; xrange 0 1\nbox\naxis");
// you may break script at any line do something
// and continue after that
parser->Execute(gr, "xlabel 'x'\nylabel 'y'\nbox");
// also you may use cycles or conditions in script
parser->Execute(gr, "for $0 -1 1 0.1\nif $0<0\n"
"line 0 0 -1 $0 'r':else:line 0 0 -1 $0 'g'\n"
"endif\nnext");
delete parser;
return 0;
}
</pre><p>The code in C/Fortran looks practically the same:
</p><pre class="verbatim">int sample(HMGL gr)
{
mgl_title(gr, "MGL parser sample", "", -2);
double a[100]; // let a_i = sin(4*pi*x), x=0...1
int i;
for(i=0;i<100;i++) a[i]=sin(4*M_PI*i/99);
HMPR parser = mgl_create_parser();
HMDT d = mgl_parser_add_var(parser, "dat");
mgl_data_set_double(d,a,100,1,1); // set data to variable
mgl_parse_text(gr, parser, "plot dat; xrange 0 1\nbox\naxis");
// you may break script at any line do something
// and continue after that
mgl_parse_text(gr, parser, "xlabel 'x'\nylabel 'y'");
// also you may use cycles or conditions in script
mgl_parse_text(gr, parser, "for $0 -1 1 0.1\nif $0<0\n"
"line 0 0 -1 $0 'r':else:line 0 0 -1 $0 'g'\n"
"endif\nnext");
mgl_write_png(gr, "test.png", ""); // don't forgot to save picture
return 0;
}
</pre>
<img src="../png/parser.png" alt="Example of MGL script parsing.">
<hr size="6">
<a name="Using-options"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#MGL-parser-using" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#g_t_0060_0060Templates_0027_0027" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Using-options-1"></a>
<h3 class="subsection">2.9.14 Using options</h3>
<p><a href="mathgl_en_27.html#Command-options">Command options</a> allow the easy setup of the selected plot by changing global settings only for this plot. Often, options are used for specifying the range of automatic variables (coordinates). However, options allows easily change plot transparency, numbers of line or faces to be drawn, or add legend entries. The sample function for options usage is:
</p><pre class="verbatim">void template(mglGraph *gr)
{
mglData a(31,41);
gr->Fill(a,"-pi*x*exp(-(y+1)^2-4*x^2)");
gr->SubPlot(2,2,0); gr->Title("Options for coordinates");
gr->Alpha(true); gr->Light(true);
gr->Rotate(40,60); gr->Box();
gr->Surf(a,"r","yrange 0 1"); gr->Surf(a,"b","yrange 0 -1");
if(mini) return;
gr->SubPlot(2,2,1); gr->Title("Option 'meshnum'");
gr->Rotate(40,60); gr->Box();
gr->Mesh(a,"r","yrange 0 1"); gr->Mesh(a,"b","yrange 0 -1; meshnum 5");
gr->SubPlot(2,2,2); gr->Title("Option 'alpha'");
gr->Rotate(40,60); gr->Box();
gr->Surf(a,"r","yrange 0 1; alpha 0.7");
gr->Surf(a,"b","yrange 0 -1; alpha 0.3");
gr->SubPlot(2,2,3,"<_"); gr->Title("Option 'legend'");
gr->FPlot("x^3","r","legend 'y = x^3'");
gr->FPlot("cos(pi*x)","b","legend 'y = cos \\pi x'");
gr->Box(); gr->Axis(); gr->Legend(2,"");
}
</pre>
<img src="../png/mirror.png" alt="Example of options usage.">
<hr size="6">
<a name="g_t_0060_0060Templates_0027_0027"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Using-options" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Stereo-image" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="g_t_0060_0060Templates_0027_0027-1"></a>
<h3 class="subsection">2.9.15 “Templates”</h3>
<p>As I have noted before, the change of settings will influence only for the further plotting commands. This allows one to create “template” function which will contain settings and primitive drawing for often used plots. Correspondingly one may call this template-function for drawing simplification.
</p>
<p>For example, let one has a set of points (experimental or numerical) and wants to compare it with theoretical law (for example, with exponent law <em>\exp(-x/2), x \in [0, 20]</em>). The template-function for this task is:
</p><pre class="verbatim">void template(mglGraph *gr)
{
mglData law(100); // create the law
law.Modify("exp(-10*x)");
gr->SetRanges(0,20, 0.0001,1);
gr->SetFunc(0,"lg(y)",0);
gr->Plot(law,"r2");
gr->Puts(mglPoint(10,0.2),"Theoretical law: e^x","r:L");
gr->Label('x',"x val."); gr->Label('y',"y val.");
gr->Axis(); gr->Grid("xy","g;"); gr->Box();
}
</pre><p>At this, one will only write a few lines for data drawing:
</p><pre class="verbatim"> template(gr); // apply settings and default drawing from template
mglData dat("fname.dat"); // load the data
// and draw it (suppose that data file have 2 columns)
gr->Plot(dat.SubData(0),dat.SubData(1),"bx ");
</pre><p>A template-function can also contain settings for font, transparency, lightning, color scheme and so on.
</p>
<p>I understand that this is obvious thing for any professional programmer, but I several times receive suggestion about “templates” ... So, I decide to point out it here.
</p>
<hr size="6">
<a name="Stereo-image"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#g_t_0060_0060Templates_0027_0027" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="#Reduce-memory-usage" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Stereo-image-1"></a>
<h3 class="subsection">2.9.16 Stereo image</h3>
<p>One can easily create stereo image in MathGL. Stereo image can be produced by making two subplots with slightly different rotation angles. The corresponding code looks like this:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
mglData a; mgls_prepare2d(&a);
gr->Light(true);
gr->SubPlot(2,1,0); gr->Rotate(50,60+1);
gr->Box(); gr->Surf(a);
gr->SubPlot(2,1,1); gr->Rotate(50,60-1);
gr->Box(); gr->Surf(a);
return 0;
}
</pre>
<img src="../png/stereo.png" alt="Example of stereo image.">
<hr size="6">
<a name="Reduce-memory-usage"></a>
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Stereo-image" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_19.html#FAQ" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en.html#Top" title="Cover (top) of document">Top</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_toc.html#SEC_Contents" title="Table of contents">Contents</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_88.html#Index" title="Index">Index</a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_abt.html#SEC_About" title="About (help)"> ? </a>]</td>
</tr></table>
<a name="Reduce-memory-usage-1"></a>
<h3 class="subsection">2.9.17 Reduce memory usage</h3>
<p>By default MathGL save all primitives in memory, rearrange it and only later draw them on bitmaps. Usually, this speed up drawing, but may require a lot of memory for plots which contain a lot of faces (like <a href="mathgl_en_41.html#cloud">cloud</a>, <a href="mathgl_en_43.html#dew">dew</a>). You can use <a href="mathgl_en_34.html#quality">quality</a> function for setting to use direct drawing on bitmap and bypassing keeping any primitives in memory. This function also allow you to decrease the quality of the resulting image but increase the speed of the drawing.
</p>
<p>The code for lowest memory usage looks like this:
</p><pre class="verbatim">int sample(mglGraph *gr)
{
gr->SetQuality(6); // firstly, set to draw directly on bitmap
for(i=0;i<1000;i++)
gr->Sphere(mglPoint(mgl_rnd()*2-1,mgl_rnd()*2-1),0.05);
return 0;
}
</pre>
<hr size="6">
<table cellpadding="1" cellspacing="1" border="0">
<tr><td valign="middle" align="left">[<a href="#Stereo-image" title="Previous section in reading order"> < </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_19.html#FAQ" title="Next section in reading order"> > </a>]</td>
<td valign="middle" align="left"> </td>
<td valign="middle" align="left">[<a href="mathgl_en_9.html#Examples" title="Beginning of this chapter or previous chapter"> << </a>]</td>
<td valign="middle" align="left">[<a href="#Hints" title="Up section"> Up </a>]</td>
<td valign="middle" align="left">[<a href="mathgl_en_20.html#General-concepts" title="Next chapter"> >> </a>]</td>
</tr></table>
<p>
<font size="-1">
This document was generated by <em>Build Daemon user</em> on <em>December 22, 2013</em> using <a href="http://www.nongnu.org/texi2html/"><em>texi2html 1.82</em></a>.
</font>
<br>
</p>
</body>
</html>
|