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 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190
|
<!DOCTYPE html>
<html xmlns="http://www.w3.org/1999/xhtml">
<head> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/sphinx_docs/html/manual/other.html" />
<meta charset="utf-8" />
<title>Other PETSc Features — PETSc 3.14.5 documentation</title>
<link rel="stylesheet" href="../_static/sphinxdoc.css" type="text/css" />
<link rel="stylesheet" href="../_static/pygments.css" type="text/css" />
<link rel="stylesheet" type="text/css" href="../_static/graphviz.css" />
<link rel="stylesheet" type="text/css" href="https://cdn.jsdelivr.net/npm/katex@0.12.0/dist/katex.min.css" />
<link rel="stylesheet" type="text/css" href="../_static/katex-math.css" />
<script id="documentation_options" data-url_root="../" src="../_static/documentation_options.js"></script>
<script src="../_static/jquery.js"></script>
<script src="../_static/underscore.js"></script>
<script src="../_static/doctools.js"></script>
<script src="../_static/language_data.js"></script>
<script src="https://cdn.jsdelivr.net/npm/katex@0.12.0/dist/katex.min.js"></script>
<script src="https://cdn.jsdelivr.net/npm/katex@0.12.0/dist/contrib/auto-render.min.js"></script>
<script src="../_static/katex_autorenderer.js"></script>
<link rel="shortcut icon" href="../_static/PETSc_RGB-logo.png"/>
<link rel="index" title="Index" href="../genindex.html" />
<link rel="search" title="Search" href="../search.html" />
<link rel="next" title="Unimportant and Advanced Features of Matrices and Solvers" href="advanced.html" />
<link rel="prev" title="Hints for Performance Tuning" href="performance.html" />
</head><body>
<div id="version" align=right><b>petsc-3.14.5 2021-03-03</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.14.5 v3.14.5 docs/sphinx_docs/html/manual/other.html "><small>Report Typos and Errors</small></a></div>
<div class="related" role="navigation" aria-label="related navigation">
<h3>Navigation</h3>
<ul>
<li class="right" style="margin-right: 10px">
<a href="../genindex.html" title="General Index"
accesskey="I">index</a></li>
<li class="right" >
<a href="advanced.html" title="Unimportant and Advanced Features of Matrices and Solvers"
accesskey="N">next</a> |</li>
<li class="right" >
<a href="performance.html" title="Hints for Performance Tuning"
accesskey="P">previous</a> |</li>
<li class="nav-item nav-item-0"><a href="../index.html">PETSc 3.14.5 documentation</a> »</li>
<li class="nav-item nav-item-1"><a href="index.html" >PETSc Users Manual</a> »</li>
<li class="nav-item nav-item-2"><a href="additional.html" accesskey="U">Additional Information</a> »</li>
</ul>
</div>
<div class="sphinxsidebar" role="navigation" aria-label="main navigation">
<div class="sphinxsidebarwrapper">
<p class="logo"><a href="../index.html">
<img class="logo" src="../_static/PETSc-TAO_RGB.svg" alt="Logo"/>
</a></p>
<h3><a href="../index.html">Table of Contents</a></h3>
<ul>
<li><a class="reference internal" href="#">Other PETSc Features</a><ul>
<li><a class="reference internal" href="#petsc-on-a-process-subset">PETSc on a process subset</a></li>
<li><a class="reference internal" href="#runtime-options">Runtime Options</a><ul>
<li><a class="reference internal" href="#the-options-database">The Options Database</a></li>
<li><a class="reference internal" href="#options-prefixes">Options Prefixes</a></li>
<li><a class="reference internal" href="#user-defined-petscoptions">User-Defined PetscOptions</a></li>
<li><a class="reference internal" href="#keeping-track-of-options">Keeping Track of Options</a></li>
</ul>
</li>
<li><a class="reference internal" href="#viewers-looking-at-petsc-objects">Viewers: Looking at PETSc Objects</a><ul>
<li><a class="reference internal" href="#viewing-from-options">Viewing From Options</a></li>
<li><a class="reference internal" href="#using-viewers-to-check-load-imbalance">Using Viewers to Check Load Imbalance</a></li>
</ul>
</li>
<li><a class="reference internal" href="#using-saws-with-petsc">Using SAWs with PETSc</a></li>
<li><a class="reference internal" href="#debugging">Debugging</a></li>
<li><a class="reference internal" href="#error-handling">Error Handling</a></li>
<li><a class="reference internal" href="#numbers">Numbers</a></li>
<li><a class="reference internal" href="#parallel-communication">Parallel Communication</a></li>
<li><a class="reference internal" href="#graphics">Graphics</a><ul>
<li><a class="reference internal" href="#windows-as-petscviewers">Windows as PetscViewers</a></li>
<li><a class="reference internal" href="#simple-petscdrawing">Simple PetscDrawing</a></li>
<li><a class="reference internal" href="#line-graphs">Line Graphs</a></li>
<li><a class="reference internal" href="#graphical-convergence-monitor">Graphical Convergence Monitor</a></li>
<li><a class="reference internal" href="#disabling-graphics-at-compile-time">Disabling Graphics at Compile Time</a></li>
</ul>
</li>
<li><a class="reference internal" href="#emacs-users">Emacs Users</a><ul>
<li><a class="reference internal" href="#tags">Tags</a></li>
</ul>
</li>
<li><a class="reference internal" href="#vs-code-users">VS Code Users</a></li>
<li><a class="reference internal" href="#vi-and-vim-users">Vi and Vim Users</a></li>
<li><a class="reference internal" href="#eclipse-users">Eclipse Users</a></li>
<li><a class="reference internal" href="#qt-creator-users">Qt Creator Users</a><ul>
<li><a class="reference internal" href="#creating-a-project">Creating a Project</a></li>
</ul>
</li>
<li><a class="reference internal" href="#visual-studio-users">Visual Studio Users</a></li>
<li><a class="reference internal" href="#xcode-users-the-apple-gui-development-system">XCode Users (The Apple GUI Development System)</a><ul>
<li><a class="reference internal" href="#mac-os-x">Mac OS X</a></li>
<li><a class="reference internal" href="#iphone-ipad-ios">iPhone/iPad iOS</a></li>
</ul>
</li>
</ul>
</li>
</ul>
<h4>Previous topic</h4>
<p class="topless"><a href="performance.html"
title="previous chapter">Hints for Performance Tuning</a></p>
<h4>Next topic</h4>
<p class="topless"><a href="advanced.html"
title="next chapter">Unimportant and Advanced Features of Matrices and Solvers</a></p>
<div role="note" aria-label="source link">
<h3>This Page</h3>
<ul class="this-page-menu">
<li><a href="../_sources/manual/other.rst.txt"
rel="nofollow">Show Source</a></li>
</ul>
</div>
<div id="searchbox" style="display: none" role="search">
<h3 id="searchlabel">Quick search</h3>
<div class="searchformwrapper">
<form class="search" action="../search.html" method="get">
<input type="text" name="q" aria-labelledby="searchlabel" />
<input type="submit" value="Go" />
</form>
</div>
</div>
<script>$('#searchbox').show(0);</script>
</div>
</div>
<div class="document">
<div class="documentwrapper">
<div class="bodywrapper">
<div class="body" role="main">
<div class="section" id="other-petsc-features">
<h1>Other PETSc Features<a class="headerlink" href="#other-petsc-features" title="Permalink to this headline">¶</a></h1>
<div class="section" id="petsc-on-a-process-subset">
<h2>PETSc on a process subset<a class="headerlink" href="#petsc-on-a-process-subset" title="Permalink to this headline">¶</a></h2>
<p>Users who wish to employ PETSc routines on only a subset of processes
within a larger parallel job, or who wish to use a “master” process to
coordinate the work of “slave” PETSc processes, should specify an
alternative communicator for <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span></code> by directly setting
its value, for example to an existing <code class="docutils literal notranslate"><span class="pre">MPI_COMM_WORLD</span></code>,</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="o">=</span><span class="n">MPI_COMM_WORLD</span><span class="p">;</span> <span class="cm">/* To use an existing MPI_COMM_WORLD */</span>
</pre></div>
</div>
<p><em>before</em> calling <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a>()</span></code>, but, obviously, after calling
<code class="docutils literal notranslate"><span class="pre"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Init.html#MPI_Init">MPI_Init</a>()</span></code>.</p>
</div>
<div class="section" id="runtime-options">
<span id="sec-options"></span><h2>Runtime Options<a class="headerlink" href="#runtime-options" title="Permalink to this headline">¶</a></h2>
<p>Allowing the user to modify parameters and options easily at runtime is
very desirable for many applications. PETSc provides a simple mechanism
to enable such customization. To print a list of available options for a
given program, simply specify the option <code class="docutils literal notranslate"><span class="pre">-help</span></code> at
runtime, e.g.,</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">mpiexec</span> <span class="o">-</span><span class="n">n</span> <span class="mi">1</span> <span class="p">.</span><span class="o">/</span><span class="n">ex1</span> <span class="o">-</span><span class="n">help</span>
</pre></div>
</div>
<p>Note that all runtime options correspond to particular PETSc routines
that can be explicitly called from within a program to set compile-time
defaults. For many applications it is natural to use a combination of
compile-time and runtime choices. For example, when solving a linear
system, one could explicitly specify use of the Krylov subspace
technique BiCGStab by calling</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetType.html#KSPSetType">KSPSetType</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPBCGS.html#KSPBCGS">KSPBCGS</a></span><span class="p">);</span>
</pre></div>
</div>
<p>One could then override this choice at runtime with the option</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="o">-</span><span class="n">ksp_type</span> <span class="n">tfqmr</span>
</pre></div>
</div>
<p>to select the Transpose-Free QMR algorithm. (See
<a class="reference internal" href="ksp.html#chapter-ksp"><span class="std std-ref">KSP: Linear System Solvers</span></a> for details.)</p>
<p>The remainder of this section discusses details of runtime options.</p>
<div class="section" id="the-options-database">
<span id="the-options-database-1"></span><h3>The Options Database<a class="headerlink" href="#the-options-database" title="Permalink to this headline">¶</a></h3>
<p>Each PETSc process maintains a database of option names and values
(stored as text strings). This database is generated with the command
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a>()</span></code>, which is listed below in its C/C++ and Fortran
variants, respectively:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a></span><span class="p">(</span><span class="kt">int</span> <span class="o">*</span><span class="n">argc</span><span class="p">,</span><span class="kt">char</span> <span class="o">***</span><span class="n">args</span><span class="p">,</span><span class="k">const</span> <span class="kt">char</span> <span class="o">*</span><span class="n">file</span><span class="p">,</span><span class="k">const</span> <span class="kt">char</span> <span class="o">*</span><span class="n">help</span><span class="p">);</span> <span class="cm">/* C */</span>
</pre></div>
</div>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a></span><span class="p">(</span><span class="kt">character </span><span class="k">file</span><span class="p">,</span><span class="kt">integer </span><span class="n">ierr</span><span class="p">)</span> <span class="c">! Fortran</span>
</pre></div>
</div>
<p>The arguments <code class="docutils literal notranslate"><span class="pre">argc</span></code> and <code class="docutils literal notranslate"><span class="pre">args</span></code> (in the C/C++ version only) are the
addresses of usual command line arguments, while the <code class="docutils literal notranslate"><span class="pre">file</span></code> is a name
of a file that can contain additional options. By default this file is
called <code class="docutils literal notranslate"><span class="pre">.petscrc</span></code> in the user’s home directory. The user can also
specify options via the environmental variable <code class="docutils literal notranslate"><span class="pre">PETSC_OPTIONS</span></code>. The
options are processed in the following order:</p>
<ol class="arabic simple">
<li><p>file</p></li>
<li><p>environmental variable</p></li>
<li><p>command line</p></li>
</ol>
<p>Thus, the command line options supersede the environmental variable
options, which in turn supersede the options file.</p>
<p>The file format for specifying options is</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>-optionname possible_value
-anotheroptionname possible_value
...
</pre></div>
</div>
<p>All of the option names must begin with a dash (-) and have no
intervening spaces. Note that the option values cannot have intervening
spaces either, and tab characters cannot be used between the option
names and values. The user can employ any naming convention. For
uniformity throughout PETSc, we employ the format
<code class="docutils literal notranslate"><span class="pre">-[prefix_]package_option</span></code> (for instance, <code class="docutils literal notranslate"><span class="pre">-ksp_type</span></code>,
<code class="docutils literal notranslate"><span class="pre">-mat_view</span> <span class="pre">::info</span></code>, or <code class="docutils literal notranslate"><span class="pre">-mg_levels_ksp_type</span></code>).</p>
<p>Users can specify an alias for any option name (to avoid typing the
sometimes lengthy default name) by adding an alias to the <code class="docutils literal notranslate"><span class="pre">.petscrc</span></code>
file in the format</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>alias -newname -oldname
</pre></div>
</div>
<p>For example,</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>alias -kspt -ksp_type
alias -sd -start_in_debugger
</pre></div>
</div>
<p>Comments can be placed in the <code class="docutils literal notranslate"><span class="pre">.petscrc</span></code> file by using <code class="docutils literal notranslate"><span class="pre">#</span></code> in the
first column of a line.</p>
</div>
<div class="section" id="options-prefixes">
<h3>Options Prefixes<a class="headerlink" href="#options-prefixes" title="Permalink to this headline">¶</a></h3>
<p>Options prefixes allow specific objects to be controlled from the
options database. For instance, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html#PCMG">PCMG</a></span></code> gives prefixes to its nested
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span></code> objects; one may control the coarse grid solver by adding the
<code class="docutils literal notranslate"><span class="pre">mg_coarse</span></code> prefix, for example <code class="docutils literal notranslate"><span class="pre">-mg_coarse_ksp_type</span> <span class="pre">preonly</span></code>. One
may also use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetOptionsPrefix.html#KSPSetOptionsPrefix">KSPSetOptionsPrefix</a>()</span></code>,<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMSetOptionsPrefix.html#DMSetOptionsPrefix">DMSetOptionsPrefix</a>()</span></code> ,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetOptionsPrefix.html#SNESSetOptionsPrefix">SNESSetOptionsPrefix</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetOptionsPrefix.html#TSSetOptionsPrefix">TSSetOptionsPrefix</a>()</span></code>, and similar
functions to assign custom prefixes, useful for applications with
multiple or nested solvers.</p>
</div>
<div class="section" id="user-defined-petscoptions">
<h3>User-Defined PetscOptions<a class="headerlink" href="#user-defined-petscoptions" title="Permalink to this headline">¶</a></h3>
<p>Any subroutine in a PETSc program can add entries to the database with
the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsSetValue.html#PetscOptionsSetValue">PetscOptionsSetValue</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">value</span><span class="p">);</span>
</pre></div>
</div>
<p>though this is rarely done. To locate options in the database, one
should use the commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">pre</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">flg</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">pre</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">value</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">flg</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetReal.html#PetscOptionsGetReal">PetscOptionsGetReal</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">pre</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">value</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">flg</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetString.html#PetscOptionsGetString">PetscOptionsGetString</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">pre</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">value</span><span class="p">,</span><span class="kt">int</span> <span class="n">maxlen</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">flg</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetStringArray.html#PetscOptionsGetStringArray">PetscOptionsGetStringArray</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">pre</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="kt">char</span> <span class="o">**</span><span class="n">values</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">nmax</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">flg</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetIntArray.html#PetscOptionsGetIntArray">PetscOptionsGetIntArray</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">pre</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="kt">int</span> <span class="o">*</span><span class="n">value</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">nmax</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">flg</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetRealArray.html#PetscOptionsGetRealArray">PetscOptionsGetRealArray</a></span><span class="p">(</span><span class="n">PetscOptions</span> <span class="n">options</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">pre</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">name</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">value</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">nmax</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">flg</span><span class="p">);</span>
</pre></div>
</div>
<p>All of these routines set <code class="docutils literal notranslate"><span class="pre">flg=<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_TRUE.html#PETSC_TRUE">PETSC_TRUE</a></span></code> if the corresponding option
was found, <code class="docutils literal notranslate"><span class="pre">flg=<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</a></span></code> if it was not found. The optional
argument <code class="docutils literal notranslate"><span class="pre">pre</span></code> indicates that the true name of the option is the given
name (with the dash “-” removed) prepended by the prefix <code class="docutils literal notranslate"><span class="pre">pre</span></code>.
Usually <code class="docutils literal notranslate"><span class="pre">pre</span></code> should be set to <code class="docutils literal notranslate"><span class="pre">NULL</span></code> (or <code class="docutils literal notranslate"><span class="pre">PETSC_NULL_CHARACTER</span></code>
for Fortran); its purpose is to allow someone to rename all the options
in a package without knowing the names of the individual options. For
example, when using block Jacobi preconditioning, the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span></code>
methods used on the individual blocks can be controlled via the options
<code class="docutils literal notranslate"><span class="pre">-sub_ksp_type</span></code> and <code class="docutils literal notranslate"><span class="pre">-sub_pc_type</span></code>.</p>
</div>
<div class="section" id="keeping-track-of-options">
<h3>Keeping Track of Options<a class="headerlink" href="#keeping-track-of-options" title="Permalink to this headline">¶</a></h3>
<p>One useful means of keeping track of user-specified runtime options is
use of <code class="docutils literal notranslate"><span class="pre">-options_view</span></code>, which prints to <code class="docutils literal notranslate"><span class="pre">stdout</span></code> during
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFinalize.html#PetscFinalize">PetscFinalize</a>()</span></code> a table of all runtime options that the user has
specified. A related option is <code class="docutils literal notranslate"><span class="pre">-options_left</span></code>, which prints the
options table and indicates any options that have <em>not</em> been requested
upon a call to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFinalize.html#PetscFinalize">PetscFinalize</a>()</span></code>. This feature is useful to check
whether an option has been activated for a particular PETSc object (such
as a solver or matrix format), or whether an option name may have been
accidentally misspelled.</p>
</div>
</div>
<div class="section" id="viewers-looking-at-petsc-objects">
<span id="sec-viewers"></span><h2>Viewers: Looking at PETSc Objects<a class="headerlink" href="#viewers-looking-at-petsc-objects" title="Permalink to this headline">¶</a></h2>
<p>PETSc employs a consistent scheme for examining, printing, and saving
objects through commands of the form</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">XXXView</span><span class="p">(</span><span class="n">XXX</span> <span class="n">obj</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">);</span>
</pre></div>
</div>
<p>Here <code class="docutils literal notranslate"><span class="pre">obj</span></code> is any PETSc object of type <code class="docutils literal notranslate"><span class="pre">XXX</span></code>, where <code class="docutils literal notranslate"><span class="pre">XXX</span></code> is
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>, etc. There are several predefined viewers.</p>
<ul class="simple">
<li><p>Passing in a zero (<code class="docutils literal notranslate"><span class="pre">0</span></code>) for the viewer causes the object to be
printed to the screen; this is useful when viewing an object in a
debugger but should be avoided in source code.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_SELF.html#PETSC_VIEWER_STDOUT_SELF">PETSC_VIEWER_STDOUT_SELF</a></span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_WORLD.html#PETSC_VIEWER_STDOUT_WORLD">PETSC_VIEWER_STDOUT_WORLD</a></span></code> causes
the object to be printed to the screen.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_DRAW_SELF.html#PETSC_VIEWER_DRAW_SELF">PETSC_VIEWER_DRAW_SELF</a></span></code> <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_DRAW_WORLD.html#PETSC_VIEWER_DRAW_WORLD">PETSC_VIEWER_DRAW_WORLD</a></span></code> causes the
object to be drawn in a default X window.</p></li>
<li><p>Passing in a viewer obtained by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerDrawOpen.html#PetscViewerDrawOpen">PetscViewerDrawOpen</a>()</span></code> causes the
object to be displayed graphically. See
<a class="reference internal" href="#sec-graphics"><span class="std std-ref">Graphics</span></a> for more on PETSc’s graphics support.</p></li>
<li><p>To save an object to a file in ASCII format, the user creates the
viewer object with the command
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerASCIIOpen.html#PetscViewerASCIIOpen">PetscViewerASCIIOpen</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="pre">comm,</span> <span class="pre">char*</span> <span class="pre">file,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="pre">*viewer)</span></code>.
This object is analogous to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_SELF.html#PETSC_VIEWER_STDOUT_SELF">PETSC_VIEWER_STDOUT_SELF</a></span></code> (for a
communicator of <code class="docutils literal notranslate"><span class="pre">MPI_COMM_SELF</span></code>) and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_WORLD.html#PETSC_VIEWER_STDOUT_WORLD">PETSC_VIEWER_STDOUT_WORLD</a></span></code>
(for a parallel communicator).</p></li>
<li><p>To save an object to a file in binary format, the user creates the
viewer object with the command
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerBinaryOpen.html#PetscViewerBinaryOpen">PetscViewerBinaryOpen</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="pre">comm,char*</span> <span class="pre">file,PetscViewerBinaryType</span> <span class="pre">type,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="pre">*viewer)</span></code>.
Details of binary I/O are discussed below.</p></li>
<li><p>Vector and matrix objects can be passed to a running MATLAB process
with a viewer created by
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerSocketOpen.html#PetscViewerSocketOpen">PetscViewerSocketOpen</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="pre">comm,char</span> <span class="pre">*machine,int</span> <span class="pre">port,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="pre">*viewer)</span></code>.
For more, see <a class="reference internal" href="matlab.html#sec-matlabsocket"><span class="std std-ref">Sending Data to an Interactive MATLAB Session</span></a>.</p></li>
</ul>
<p>The user can control the format of ASCII printed objects with viewers
created by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerASCIIOpen.html#PetscViewerASCIIOpen">PetscViewerASCIIOpen</a>()</span></code> by calling</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerPushFormat.html#PetscViewerPushFormat">PetscViewerPushFormat</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PetscViewerFormat</a></span> <span class="n">format</span><span class="p">);</span>
</pre></div>
</div>
<p>Formats include <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PETSC_VIEWER_DEFAULT</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PETSC_VIEWER_ASCII_MATLAB</a></span></code>,
and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PETSC_VIEWER_ASCII_IMPL</a></span></code>. The implementation-specific format,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PETSC_VIEWER_ASCII_IMPL</a></span></code>, displays the object in the most natural way
for a particular implementation.</p>
<p>The routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerPushFormat.html#PetscViewerPushFormat">PetscViewerPushFormat</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PetscViewerFormat</a></span> <span class="n">format</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerPopFormat.html#PetscViewerPopFormat">PetscViewerPopFormat</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">);</span>
</pre></div>
</div>
<p>allow one to temporarily change the format of a viewer.</p>
<p>As discussed above, one can output PETSc objects in binary format by
first opening a binary viewer with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerBinaryOpen.html#PetscViewerBinaryOpen">PetscViewerBinaryOpen</a>()</span></code> and then
using <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatView.html#MatView">MatView</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a>()</span></code>, etc. The corresponding routines for
input of a binary object have the form <code class="docutils literal notranslate"><span class="pre">XXXLoad()</span></code>. In particular,
matrix and vector binary input is handled by the following routines:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatLoad.html#MatLoad">MatLoad</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatType.html#MatType">MatType</a></span> <span class="n">outtype</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="o">*</span><span class="n">newmat</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecLoad.html#VecLoad">VecLoad</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecType.html#VecType">VecType</a></span> <span class="n">outtype</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">newvec</span><span class="p">);</span>
</pre></div>
</div>
<p>These routines generate parallel matrices and vectors if the viewer’s
communicator has more than one process. The particular matrix and vector
formats are determined from the options database; see the manual pages
for details.</p>
<p>One can provide additional information about matrix data for matrices
stored on disk by providing an optional file <code class="docutils literal notranslate"><span class="pre">matrixfilename.info</span></code>,
where <code class="docutils literal notranslate"><span class="pre">matrixfilename</span></code> is the name of the file containing the matrix.
The format of the optional file is the same as the <code class="docutils literal notranslate"><span class="pre">.petscrc</span></code> file and
can (currently) contain the following:</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>-matload_block_size <bs>
</pre></div>
</div>
<p>The block size indicates the size of blocks to use if the matrix is read
into a block oriented data structure (for example, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATMPIBAIJ.html#MATMPIBAIJ">MATMPIBAIJ</a></span></code>). The
diagonal information <code class="docutils literal notranslate"><span class="pre">s1,s2,s3,...</span></code> indicates which (block) diagonals
in the matrix have nonzero values.</p>
<div class="section" id="viewing-from-options">
<span id="sec-viewfromoptions"></span><h3>Viewing From Options<a class="headerlink" href="#viewing-from-options" title="Permalink to this headline">¶</a></h3>
<p>Command-line options provide a particularly convenient way to view PETSc
objects. All options of the form <code class="docutils literal notranslate"><span class="pre">-xxx_view</span></code> accept
colon(<code class="docutils literal notranslate"><span class="pre">:</span></code>)-separated compound arguments which specify a viewer type,
format, and/or destination (e.g. file name or socket) if appropriate.
For example, to quickly export a binary file containing a matrix, one
may use <code class="docutils literal notranslate"><span class="pre">-mat_view</span> <span class="pre">binary:matrix.out</span></code>, or to output to a
MATLAB-compatible ASCII file, one may use
<code class="docutils literal notranslate"><span class="pre">-mat_view</span> <span class="pre">ascii:matrix.m:ascii_matlab</span></code>. See the
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscOptionsGetViewer.html#PetscOptionsGetViewer">PetscOptionsGetViewer</a>()</span></code> man page for full details, as well as the
<code class="docutils literal notranslate"><span class="pre">XXXViewFromOptions()</span></code> man pages (for instance,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetFromOptions.html#PetscDrawSetFromOptions">PetscDrawSetFromOptions</a>()</span></code>) for many other convenient command-line
options.</p>
</div>
<div class="section" id="using-viewers-to-check-load-imbalance">
<h3>Using Viewers to Check Load Imbalance<a class="headerlink" href="#using-viewers-to-check-load-imbalance" title="Permalink to this headline">¶</a></h3>
<p>The PetscViewer format <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PETSC_VIEWER_LOAD_BALANCE</a></span></code> will cause certain
objects to display simple measures of their imbalance. For example</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>-n 4 ./ex32 -ksp_view_mat ::load_balance
</pre></div>
</div>
<p>will display</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>Nonzeros: Min 162 avg 168 max 174
</pre></div>
</div>
<p>indicating that one process has 162 nonzero entries in the matrix, the
average number of nonzeros per process is 168 and the maximum number of
nonzeros is 174. Similar for vectors one can see the load balancing
with, for example,</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>-n 4 ./ex32 -ksp_view_rhs ::load_balance
</pre></div>
</div>
<p>The measurements of load balancing can also be done within the program
with calls to the appropriate object viewer with the viewer format
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerFormat.html#PetscViewerFormat">PETSC_VIEWER_LOAD_BALANCE</a></span></code>.</p>
</div>
</div>
<div class="section" id="using-saws-with-petsc">
<span id="sec-saws"></span><h2>Using SAWs with PETSc<a class="headerlink" href="#using-saws-with-petsc" title="Permalink to this headline">¶</a></h2>
<p>The Scientific Application Web server, SAWs <a class="footnote-reference brackets" href="#id5" id="id1">8</a>, allows one to monitor
running PETSc applications from a browser. <code class="docutils literal notranslate"><span class="pre">./configure</span></code> PETSc with
the additional option <code class="docutils literal notranslate"><span class="pre">--download-saws</span></code>. Options to use SAWs include</p>
<ul class="simple">
<li><p><code class="docutils literal notranslate"><span class="pre">-saws_options</span></code> - allows setting values in the PETSc options
database via the browser (works only on one process).</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-stack_view</span> <span class="pre">saws</span></code> - allows monitoring the current stack frame that
PETSc is in; refresh to see the new location.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_monitor_saws,</span> <span class="pre">-ksp_monitor_saws</span></code> - monitor the solvers’
iterations from the web browser.</p></li>
</ul>
<p>For each of these you need to point your browser to
<code class="docutils literal notranslate"><span class="pre">http://hostname:8080</span></code>, for example <code class="docutils literal notranslate"><span class="pre">http://localhost:8080</span></code>. Options
that control behavior of SAWs include</p>
<ul class="simple">
<li><p><code class="docutils literal notranslate"><span class="pre">-saws_log</span> <span class="pre">filename</span></code> - log all SAWs actions in a file.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-saws_https</span> <span class="pre">certfile</span></code> - use HTTPS instead of HTTP with a
certificate.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-saws_port_auto_select</span></code> - have SAWs pick a port number instead of
using 8080.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-saws_port</span> <span class="pre">port</span></code> - use <code class="docutils literal notranslate"><span class="pre">port</span></code> instead of 8080.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-saws_root</span> <span class="pre">rootdirectory</span></code> - local directory to which the SAWs
browser will have read access.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-saws_local</span></code> - use the local file system to obtain the SAWS
javascript files (they much be in <code class="docutils literal notranslate"><span class="pre">rootdirectory/js</span></code>).</p></li>
</ul>
<p>Also see the manual pages for <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscSAWsBlock.html#PetscSAWsBlock">PetscSAWsBlock</a></span></code>,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSAWsTakeAccess.html#PetscObjectSAWsTakeAccess">PetscObjectSAWsTakeAccess</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSAWsGrantAccess.html#PetscObjectSAWsGrantAccess">PetscObjectSAWsGrantAccess</a></span></code>,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSAWsSetBlock.html#PetscObjectSAWsSetBlock">PetscObjectSAWsSetBlock</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscStackSAWsGrantAccess.html#PetscStackSAWsGrantAccess">PetscStackSAWsGrantAccess</a></span></code>
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscStackSAWsTakeAccess.html#PetscStackSAWsTakeAccess">PetscStackSAWsTakeAccess</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPMonitorSAWs.html#KSPMonitorSAWs">KSPMonitorSAWs</a></span></code>, and
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESMonitorSAWs.html#SNESMonitorSAWs">SNESMonitorSAWs</a></span></code>.</p>
</div>
<div class="section" id="debugging">
<span id="sec-debugging"></span><h2>Debugging<a class="headerlink" href="#debugging" title="Permalink to this headline">¶</a></h2>
<p>PETSc programs may be debugged using one of the two options below.</p>
<ul class="simple">
<li><p><code class="docutils literal notranslate"><span class="pre">-start_in_debugger</span></code> <code class="docutils literal notranslate"><span class="pre">[noxterm,dbx,xxgdb,xdb,xldb,lldb]</span></code>
<code class="docutils literal notranslate"><span class="pre">[-display</span> <span class="pre">name]</span></code> - start all processes in debugger</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-on_error_attach_debugger</span></code> <code class="docutils literal notranslate"><span class="pre">[noxterm,dbx,xxgdb,xdb,xldb,lldb]</span></code>
<code class="docutils literal notranslate"><span class="pre">[-display</span> <span class="pre">name]</span></code> - start debugger only on encountering an error</p></li>
</ul>
<p>Note that, in general, debugging MPI programs cannot be done in the
usual manner of starting the programming in the debugger (because then
it cannot set up the MPI communication and remote processes).</p>
<p>By default the GNU debugger <code class="docutils literal notranslate"><span class="pre">gdb</span></code> is used when <code class="docutils literal notranslate"><span class="pre">-start_in_debugger</span></code>
or <code class="docutils literal notranslate"><span class="pre">-on_error_attach_debugger</span></code> is specified. To employ either
<code class="docutils literal notranslate"><span class="pre">xxgdb</span></code> or the common UNIX debugger <code class="docutils literal notranslate"><span class="pre">dbx</span></code>, one uses command line
options as indicated above. On HP-UX machines the debugger <code class="docutils literal notranslate"><span class="pre">xdb</span></code>
should be used instead of <code class="docutils literal notranslate"><span class="pre">dbx</span></code>; on RS/6000 machines the <code class="docutils literal notranslate"><span class="pre">xldb</span></code>
debugger is supported as well. On OS X systems with XCode tools,
<code class="docutils literal notranslate"><span class="pre">lldb</span></code> is available. By default, the debugger will be started in a new
xterm (to enable running separate debuggers on each process), unless the
option <code class="docutils literal notranslate"><span class="pre">noxterm</span></code> is used. In order to handle the MPI startup phase,
the debugger command <code class="docutils literal notranslate"><span class="pre">cont</span></code> should be used to continue execution of
the program within the debugger. Rerunning the program through the
debugger requires terminating the first job and restarting the
processor(s); the usual <code class="docutils literal notranslate"><span class="pre">run</span></code> option in the debugger will not
correctly handle the MPI startup and should not be used. Not all
debuggers work on all machines, the user may have to experiment to find
one that works correctly.</p>
<p>You can select a subset of the processes to be debugged (the rest just
run without the debugger) with the option</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>-debugger_ranks rank1,rank2,...
</pre></div>
</div>
<p>where you simply list the ranks you want the debugger to run with.</p>
</div>
<div class="section" id="error-handling">
<h2>Error Handling<a class="headerlink" href="#error-handling" title="Permalink to this headline">¶</a></h2>
<p>Errors are handled through the routine <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscError.html#PetscError">PetscError</a>()</span></code>. This routine
checks a stack of error handlers and calls the one on the top. If the
stack is empty, it selects <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscTraceBackErrorHandler.html#PetscTraceBackErrorHandler">PetscTraceBackErrorHandler</a>()</span></code>, which tries
to print a traceback. A new error handler can be put on the stack with</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPushErrorHandler.html#PetscPushErrorHandler">PetscPushErrorHandler</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">HandlerFunction</span><span class="p">)(</span><span class="kt">int</span> <span class="n">line</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">dir</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">file</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">message</span><span class="p">,</span><span class="kt">int</span> <span class="n">number</span><span class="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">),</span><span class="kt">void</span> <span class="o">*</span><span class="n">HandlerContext</span><span class="p">)</span>
</pre></div>
</div>
<p>The arguments to <code class="docutils literal notranslate"><span class="pre">HandlerFunction()</span></code> are the line number where the
error occurred, the file in which the error was detected, the
corresponding directory, the error message, the error integer, and the
<code class="docutils literal notranslate"><span class="pre">HandlerContext.</span></code> The routine</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPopErrorHandler.html#PetscPopErrorHandler">PetscPopErrorHandler</a></span><span class="p">()</span>
</pre></div>
</div>
<p>removes the last error handler and discards it.</p>
<p>PETSc provides two additional error handlers besides
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscTraceBackErrorHandler.html#PetscTraceBackErrorHandler">PetscTraceBackErrorHandler</a>()</span></code>:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscAbortErrorHandler.html#PetscAbortErrorHandler">PetscAbortErrorHandler</a></span><span class="p">()</span>
<span class="n">PetscAttachErrorHandler</span><span class="p">()</span>
</pre></div>
</div>
<p>The function <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscAbortErrorHandler.html#PetscAbortErrorHandler">PetscAbortErrorHandler</a>()</span></code> calls abort on encountering an
error, while <code class="docutils literal notranslate"><span class="pre">PetscAttachErrorHandler()</span></code> attaches a debugger to the
running process if an error is detected. At runtime, these error
handlers can be set with the options <code class="docutils literal notranslate"><span class="pre">-on_error_abort</span></code> or
<code class="docutils literal notranslate"><span class="pre">-on_error_attach_debugger</span></code> <code class="docutils literal notranslate"><span class="pre">[noxterm,</span> <span class="pre">dbx,</span> <span class="pre">xxgdb,</span> <span class="pre">xldb]</span></code>
<code class="docutils literal notranslate"><span class="pre">[-display</span> <span class="pre">DISPLAY]</span></code>.</p>
<p>All PETSc calls can be traced (useful for determining where a program is
hanging without running in the debugger) with the option</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>-log_trace [filename]
</pre></div>
</div>
<p>where <code class="docutils literal notranslate"><span class="pre">filename</span></code> is optional. By default the traces are printed to the
screen. This can also be set with the command
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Profiling/PetscLogTraceBegin.html#PetscLogTraceBegin">PetscLogTraceBegin</a>(FILE*)</span></code>.</p>
<p>It is also possible to trap signals by using the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPushSignalHandler.html#PetscPushSignalHandler">PetscPushSignalHandler</a></span><span class="p">(</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">Handler</span><span class="p">)(</span><span class="kt">int</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="p">),</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>The default handler <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscSignalHandlerDefault.html#PetscSignalHandlerDefault">PetscSignalHandlerDefault</a>()</span></code> calls
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscError.html#PetscError">PetscError</a>()</span></code> and then terminates. In general, a signal in PETSc
indicates a catastrophic failure. Any error handler that the user
provides should try to clean up only before exiting. By default all
PETSc programs use the default signal handler, although the user can
turn this off at runtime with the option <code class="docutils literal notranslate"><span class="pre">-no_signal_handler</span></code> .</p>
<p>There is a separate signal handler for floating-point exceptions. The
option <code class="docutils literal notranslate"><span class="pre">-fp_trap</span></code> turns on the floating-point trap at runtime, and the
routine</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscSetFPTrap.html#PetscSetFPTrap">PetscSetFPTrap</a></span><span class="p">(</span><span class="n">PetscFPTrap</span> <span class="n">flag</span><span class="p">);</span>
</pre></div>
</div>
<p>can be used in-line. A <code class="docutils literal notranslate"><span class="pre">flag</span></code> of <code class="docutils literal notranslate"><span class="pre">PETSC_FP_TRAP_ON</span></code> indicates that
floating-point exceptions should be trapped, while a value of
<code class="docutils literal notranslate"><span class="pre">PETSC_FP_TRAP_OFF</span></code> (the default) indicates that they should be
ignored. Note that on certain machines, in particular the IBM RS/6000,
trapping is very expensive.</p>
<p>A small set of macros is used to make the error handling lightweight.
These macros are used throughout the PETSc libraries and can be employed
by the application programmer as well. When an error is first detected,
one should set it by calling</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/SETERRQ.html#SETERRQ">SETERRQ</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">flag</span><span class="p">,,</span><span class="kt">char</span> <span class="o">*</span><span class="n">message</span><span class="p">);</span>
</pre></div>
</div>
<p>The user should check the return codes for all PETSc routines (and
possibly user-defined routines as well) with</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">ierr</span> <span class="o">=</span> <span class="n">PetscRoutine</span><span class="p">(...);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">);</span>
</pre></div>
</div>
<p>Likewise, all memory allocations should be checked with</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscMalloc1.html#PetscMalloc1">PetscMalloc1</a></span><span class="p">(</span><span class="n">n</span><span class="p">,</span> <span class="o">&</span><span class="n">ptr</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
</pre></div>
</div>
<p>If this procedure is followed throughout all of the user’s libraries and
codes, any error will by default generate a clean traceback of the
location of the error.</p>
<p>Note that the macro <code class="docutils literal notranslate"><span class="pre">PETSC_FUNCTION_NAME</span></code> is used to keep track of
routine names during error tracebacks. Users need not worry about this
macro in their application codes; however, users can take advantage of
this feature if desired by setting this macro before each user-defined
routine that may call <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/SETERRQ.html#SETERRQ">SETERRQ</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a>()</span></code>. A simple example of
usage is given below.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">MyRoutine1</span><span class="p">()</span> <span class="p">{</span>
<span class="cm">/* Declarations Here */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="cm">/* code here */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
</pre></div>
</div>
</div>
<div class="section" id="numbers">
<span id="sec-complex"></span><h2>Numbers<a class="headerlink" href="#numbers" title="Permalink to this headline">¶</a></h2>
<p>PETSc supports the use of complex numbers in application programs
written in C, C++, and Fortran. To do so, we employ either the C99
<code class="docutils literal notranslate"><span class="pre">complex</span></code> type or the C++ versions of the PETSc libraries in which the
basic “scalar” datatype, given in PETSc codes by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span></code>, is
defined as <code class="docutils literal notranslate"><span class="pre">complex</span></code> (or <code class="docutils literal notranslate"><span class="pre">complex<double></span></code> for machines using
templated complex class libraries). To work with complex numbers, the
user should run <code class="docutils literal notranslate"><span class="pre">./configure</span></code> with the additional option
<code class="docutils literal notranslate"><span class="pre">--with-scalar-type=complex</span></code>. The
<a class="reference external" href="https://www.mcs.anl.gov/petsc/documentation/installation.html">installation instructions</a>
provide detailed instructions for installing PETSc. You can use
<code class="docutils literal notranslate"><span class="pre">--with-clanguage=c</span></code> (the default) to use the C99 complex numbers or
<code class="docutils literal notranslate"><span class="pre">--with-clanguage=c++</span></code> to use the C++ complex type <a class="footnote-reference brackets" href="#id6" id="id2">9</a>.</p>
<p>Recall that each variant of the PETSc libraries is stored in a different
directory, given by <code class="docutils literal notranslate"><span class="pre">${PETSC_DIR}/lib/${PETSC_ARCH}</span></code></p>
<p>according to the architecture. Thus, the libraries for complex numbers
are maintained separately from those for real numbers. When using any of
the complex numbers versions of PETSc, <em>all</em> vector and matrix elements
are treated as complex, even if their imaginary components are zero. Of
course, one can elect to use only the real parts of the complex numbers
when using the complex versions of the PETSc libraries; however, when
working <em>only</em> with real numbers in a code, one should use a version of
PETSc for real numbers for best efficiency.</p>
<p>The program
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/ksp/ksp/tutorials/ex11.c.html">KSP Tutorial ex11</a>
solves a linear system with a complex coefficient matrix. Its Fortran
counterpart is
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/ksp/ksp/tutorials/ex11f.F90.html">KSP Tutorial ex11f</a>.</p>
</div>
<div class="section" id="parallel-communication">
<h2>Parallel Communication<a class="headerlink" href="#parallel-communication" title="Permalink to this headline">¶</a></h2>
<p>When used in a message-passing environment, all communication within
PETSc is done through MPI, the message-passing interface standard
<span id="id3">[<a class="reference internal" href="getting_started.html#id201"><span>For94</span></a>]</span>. Any file that includes <code class="docutils literal notranslate"><span class="pre">petscsys.h</span></code> (or
any other PETSc include file) can freely use any MPI routine.</p>
</div>
<div class="section" id="graphics">
<span id="sec-graphics"></span><h2>Graphics<a class="headerlink" href="#graphics" title="Permalink to this headline">¶</a></h2>
<p>The PETSc graphics library is not intended to compete with high-quality
graphics packages. Instead, it is intended to be easy to use
interactively with PETSc programs. We urge users to generate their
publication-quality graphics using a professional graphics package. If a
user wants to hook certain packages into PETSc, he or she should send a
message to
<a class="reference external" href="mailto:petsc-maint%40mcs.anl.gov">petsc-maint<span>@</span>mcs<span>.</span>anl<span>.</span>gov</a>; we
will see whether it is reasonable to try to provide direct interfaces.</p>
<div class="section" id="windows-as-petscviewers">
<h3>Windows as PetscViewers<a class="headerlink" href="#windows-as-petscviewers" title="Permalink to this headline">¶</a></h3>
<p>For drawing predefined PETSc objects such as matrices and vectors, one
may first create a viewer using the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerDrawOpen.html#PetscViewerDrawOpen">PetscViewerDrawOpen</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">display</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">title</span><span class="p">,</span><span class="kt">int</span> <span class="n">x</span><span class="p">,</span><span class="kt">int</span> <span class="n">y</span><span class="p">,</span><span class="kt">int</span> <span class="n">w</span><span class="p">,</span><span class="kt">int</span> <span class="n">h</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="o">*</span><span class="n">viewer</span><span class="p">);</span>
</pre></div>
</div>
<p>This viewer may be passed to any of the <code class="docutils literal notranslate"><span class="pre">XXXView()</span></code> routines.
Alternately, one may use command-line options to quickly specify viewer
formats, including <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span></code>-based ones; see
<a class="reference internal" href="#sec-viewfromoptions"><span class="std std-ref">Viewing From Options</span></a>.</p>
<p>To draw directly into the viewer, one must obtain the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span></code>
object with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerDrawGetDraw.html#PetscViewerDrawGetDraw">PetscViewerDrawGetDraw</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="o">*</span><span class="n">draw</span><span class="p">);</span>
</pre></div>
</div>
<p>Then one can call any of the <code class="docutils literal notranslate"><span class="pre">PetscDrawXXX</span></code> commands on the <code class="docutils literal notranslate"><span class="pre">draw</span></code>
object. If one obtains the <code class="docutils literal notranslate"><span class="pre">draw</span></code> object in this manner, one does not
call the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawOpenX.html#PetscDrawOpenX">PetscDrawOpenX</a>()</span></code> command discussed below.</p>
<p>Predefined viewers, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_DRAW_WORLD.html#PETSC_VIEWER_DRAW_WORLD">PETSC_VIEWER_DRAW_WORLD</a></span></code> and
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_DRAW_SELF.html#PETSC_VIEWER_DRAW_SELF">PETSC_VIEWER_DRAW_SELF</a></span></code>, may be used at any time. Their initial use
will cause the appropriate window to be created.</p>
<p>Implementations using OpenGL, TikZ, and other formats may be selected
with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetType.html#PetscDrawSetType">PetscDrawSetType</a>()</span></code>. PETSc can also produce movies; see
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetSaveMovie.html#PetscDrawSetSaveMovie">PetscDrawSetSaveMovie</a>()</span></code>, and note that command-line options can also
be convenient; see the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetFromOptions.html#PetscDrawSetFromOptions">PetscDrawSetFromOptions</a>()</span></code> man page.</p>
<p>By default, PETSc drawing tools employ a private colormap, which
remedies the problem of poor color choices for contour plots due to an
external program’s mangling of the colormap. Unfortunately, this may
cause flashing of colors as the mouse is moved between the PETSc windows
and other windows. Alternatively, a shared colormap can be used via the
option <code class="docutils literal notranslate"><span class="pre">-draw_x_shared_colormap</span></code>.</p>
</div>
<div class="section" id="simple-petscdrawing">
<h3>Simple PetscDrawing<a class="headerlink" href="#simple-petscdrawing" title="Permalink to this headline">¶</a></h3>
<p>With the default format, one can open a window that is not associated
with a viewer directly under the X11 Window System or OpenGL with the
command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawCreate.html#PetscDrawCreate">PetscDrawCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">display</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">title</span><span class="p">,</span><span class="kt">int</span> <span class="n">x</span><span class="p">,</span><span class="kt">int</span> <span class="n">y</span><span class="p">,</span><span class="kt">int</span> <span class="n">w</span><span class="p">,</span><span class="kt">int</span> <span class="n">h</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="o">*</span><span class="n">win</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetFromOptions.html#PetscDrawSetFromOptions">PetscDrawSetFromOptions</a></span><span class="p">(</span><span class="n">win</span><span class="p">);</span>
</pre></div>
</div>
<p>All drawing routines are performed relative to the window’s coordinate
system and viewport. By default, the drawing coordinates are from
<code class="docutils literal notranslate"><span class="pre">(0,0)</span></code> to <code class="docutils literal notranslate"><span class="pre">(1,1)</span></code>, where <code class="docutils literal notranslate"><span class="pre">(0,0)</span></code> indicates the lower left corner
of the window. The application program can change the window coordinates
with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetCoordinates.html#PetscDrawSetCoordinates">PetscDrawSetCoordinates</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xl</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">yl</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xr</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">yr</span><span class="p">);</span>
</pre></div>
</div>
<p>By default, graphics will be drawn in the entire window. To restrict the
drawing to a portion of the window, one may use the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetViewPort.html#PetscDrawSetViewPort">PetscDrawSetViewPort</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xl</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">yl</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xr</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">yr</span><span class="p">);</span>
</pre></div>
</div>
<p>These arguments, which indicate the fraction of the window in which the
drawing should be done, must satisfy
<span class="math">\(0 \leq {\tt xl} \leq {\tt xr} \leq 1\)</span> and
<span class="math">\(0 \leq {\tt yl} \leq {\tt yr} \leq 1.\)</span></p>
<p>To draw a line, one uses the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLine.html#PetscDrawLine">PetscDrawLine</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xl</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">yl</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xr</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">yr</span><span class="p">,</span><span class="kt">int</span> <span class="n">cl</span><span class="p">);</span>
</pre></div>
</div>
<p>The argument <code class="docutils literal notranslate"><span class="pre">cl</span></code> indicates the color (which is an integer between 0
and 255) of the line. A list of predefined colors may be found in
<code class="docutils literal notranslate"><span class="pre">include/petscdraw.h</span></code> and includes <code class="docutils literal notranslate"><span class="pre">PETSC_DRAW_BLACK</span></code>,
<code class="docutils literal notranslate"><span class="pre">PETSC_DRAW_RED</span></code>, <code class="docutils literal notranslate"><span class="pre">PETSC_DRAW_BLUE</span></code> etc.</p>
<p>To ensure that all graphics actually have been displayed, one should use
the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawFlush.html#PetscDrawFlush">PetscDrawFlush</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">);</span>
</pre></div>
</div>
<p>When displaying by using double buffering, which is set with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetDoubleBuffer.html#PetscDrawSetDoubleBuffer">PetscDrawSetDoubleBuffer</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">);</span>
</pre></div>
</div>
<p><em>all</em> processes must call</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawFlush.html#PetscDrawFlush">PetscDrawFlush</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">);</span>
</pre></div>
</div>
<p>in order to swap the buffers. From the options database one may use
<code class="docutils literal notranslate"><span class="pre">-draw_pause</span></code> <code class="docutils literal notranslate"><span class="pre">n</span></code>, which causes the PETSc application to pause <code class="docutils literal notranslate"><span class="pre">n</span></code>
seconds at each <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawPause.html#PetscDrawPause">PetscDrawPause</a>()</span></code>. A time of <code class="docutils literal notranslate"><span class="pre">-1</span></code> indicates that
the application should pause until receiving mouse input from the user.</p>
<p>Text can be drawn with commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawString.html#PetscDrawString">PetscDrawString</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">y</span><span class="p">,</span><span class="kt">int</span> <span class="n">color</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">text</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawStringVertical.html#PetscDrawStringVertical">PetscDrawStringVertical</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">y</span><span class="p">,</span><span class="kt">int</span> <span class="n">color</span><span class="p">,</span><span class="k">const</span> <span class="kt">char</span> <span class="o">*</span><span class="n">text</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawStringCentered.html#PetscDrawStringCentered">PetscDrawStringCentered</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">y</span><span class="p">,</span><span class="kt">int</span> <span class="n">color</span><span class="p">,</span><span class="k">const</span> <span class="kt">char</span> <span class="o">*</span><span class="n">text</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawStringBoxed.html#PetscDrawStringBoxed">PetscDrawStringBoxed</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">draw</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">sxl</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">syl</span><span class="p">,</span><span class="kt">int</span> <span class="n">sc</span><span class="p">,</span><span class="kt">int</span> <span class="n">bc</span><span class="p">,</span><span class="k">const</span> <span class="kt">char</span> <span class="n">text</span><span class="p">[],</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">w</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">h</span><span class="p">);</span>
</pre></div>
</div>
<p>The user can set the text font size or determine it with the commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawStringSetSize.html#PetscDrawStringSetSize">PetscDrawStringSetSize</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">width</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">height</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawStringGetSize.html#PetscDrawStringGetSize">PetscDrawStringGetSize</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">width</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">height</span><span class="p">);</span>
</pre></div>
</div>
</div>
<div class="section" id="line-graphs">
<h3>Line Graphs<a class="headerlink" href="#line-graphs" title="Permalink to this headline">¶</a></h3>
<p>PETSc includes a set of routines for manipulating simple two-dimensional
graphs. These routines, which begin with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxisDraw.html#PetscDrawAxisDraw">PetscDrawAxisDraw</a>()</span></code>, are
usually not used directly by the application programmer. Instead, the
programmer employs the line graph routines to draw simple line graphs.
As shown in the <a class="reference internal" href="#listing-draw-test-ex3"><span class="std std-ref">listing below</span></a>, line
graphs are created with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGCreate.html#PetscDrawLGCreate">PetscDrawLGCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">win</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">ncurves</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>The argument <code class="docutils literal notranslate"><span class="pre">ncurves</span></code> indicates how many curves are to be drawn.
Points can be added to each of the curves with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGAddPoint.html#PetscDrawLGAddPoint">PetscDrawLGAddPoint</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">y</span><span class="p">);</span>
</pre></div>
</div>
<p>The arguments <code class="docutils literal notranslate"><span class="pre">x</span></code> and <code class="docutils literal notranslate"><span class="pre">y</span></code> are arrays containing the next point value
for each curve. Several points for each curve may be added with</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGAddPoints.html#PetscDrawLGAddPoints">PetscDrawLGAddPoints</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">**</span><span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">**</span><span class="n">y</span><span class="p">);</span>
</pre></div>
</div>
<p>The line graph is drawn (or redrawn) with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGDraw.html#PetscDrawLGDraw">PetscDrawLGDraw</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>A line graph that is no longer needed can be destroyed with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGDestroy.html#PetscDrawLGDestroy">PetscDrawLGDestroy</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>To plot new curves, one can reset a linegraph with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGReset.html#PetscDrawLGReset">PetscDrawLGReset</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>The line graph automatically determines the range of values to display
on the two axes. The user can change these defaults with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGSetLimits.html#PetscDrawLGSetLimits">PetscDrawLGSetLimits</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xmin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xmax</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">ymin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">ymax</span><span class="p">);</span>
</pre></div>
</div>
<p>It is also possible to change the display of the axes and to label them.
This procedure is done by first obtaining the axes context with the
command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGGetAxis.html#PetscDrawLGGetAxis">PetscDrawLGGetAxis</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxis.html#PetscDrawAxis">PetscDrawAxis</a></span> <span class="o">*</span><span class="n">axis</span><span class="p">);</span>
</pre></div>
</div>
<p>One can set the axes’ colors and labels, respectively, by using the
commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxisSetColors.html#PetscDrawAxisSetColors">PetscDrawAxisSetColors</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxis.html#PetscDrawAxis">PetscDrawAxis</a></span> <span class="n">axis</span><span class="p">,</span><span class="kt">int</span> <span class="n">axis_lines</span><span class="p">,</span><span class="kt">int</span> <span class="n">ticks</span><span class="p">,</span><span class="kt">int</span> <span class="n">text</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxisSetLabels.html#PetscDrawAxisSetLabels">PetscDrawAxisSetLabels</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxis.html#PetscDrawAxis">PetscDrawAxis</a></span> <span class="n">axis</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">top</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">x</span><span class="p">,</span><span class="kt">char</span> <span class="o">*</span><span class="n">y</span><span class="p">);</span>
</pre></div>
</div>
<p>It is possible to turn off all graphics with the option <code class="docutils literal notranslate"><span class="pre">-nox</span></code>. This
will prevent any windows from being opened or any drawing actions to be
done. This is useful for running large jobs when the graphics overhead
is too large, or for timing.</p>
<p>The full example, <a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/sys/classes/draw/tests/ex3.c.html">Draw Test ex3</a>,
follows.</p>
<div class="admonition-listing-src-classes-draw-tests-ex3-c admonition" id="snes-ex1">
<span id="listing-draw-test-ex3"></span><p class="admonition-title">Listing: <code class="docutils literal notranslate"><span class="pre">src/classes/draw/tests/ex3.c</span></code></p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span>
<span class="k">static</span> <span class="kt">char</span> <span class="n">help</span><span class="p">[]</span> <span class="o">=</span> <span class="s">"Plots a simple line graph.</span><span class="se">\n</span><span class="s">"</span><span class="p">;</span>
<span class="cp">#if defined(PETSC_APPLE_FRAMEWORK)</span>
<span class="cp">#import <PETSc/petscsys.h></span>
<span class="cp">#import <PETSc/petscdraw.h></span>
<span class="cp">#else</span>
<span class="cp">#include</span> <span class="cpf"><petscsys.h></span><span class="cp"></span>
<span class="cp">#include</span> <span class="cpf"><petscdraw.h></span><span class="cp"></span>
<span class="cp">#endif</span>
<span class="kt">int</span> <span class="nf">main</span><span class="p">(</span><span class="kt">int</span> <span class="n">argc</span><span class="p">,</span><span class="kt">char</span> <span class="o">**</span><span class="n">argv</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a></span> <span class="n">draw</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLG.html#PetscDrawLG">PetscDrawLG</a></span> <span class="n">lg</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxis.html#PetscDrawAxis">PetscDrawAxis</a></span> <span class="n">axis</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span> <span class="o">=</span> <span class="mi">15</span><span class="p">,</span><span class="n">i</span><span class="p">,</span><span class="n">x</span> <span class="o">=</span> <span class="mi">0</span><span class="p">,</span><span class="n">y</span> <span class="o">=</span> <span class="mi">0</span><span class="p">,</span><span class="n">width</span> <span class="o">=</span> <span class="mi">400</span><span class="p">,</span><span class="n">height</span> <span class="o">=</span> <span class="mi">300</span><span class="p">,</span><span class="n">nports</span> <span class="o">=</span> <span class="mi">1</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">useports</span><span class="p">,</span><span class="n">flg</span><span class="p">;</span>
<span class="k">const</span> <span class="kt">char</span> <span class="o">*</span><span class="n">xlabel</span><span class="p">,</span><span class="o">*</span><span class="n">ylabel</span><span class="p">,</span><span class="o">*</span><span class="n">toplabel</span><span class="p">,</span><span class="o">*</span><span class="n">legend</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xd</span><span class="p">,</span><span class="n">yd</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawViewPorts.html#PetscDrawViewPorts">PetscDrawViewPorts</a></span> <span class="o">*</span><span class="n">ports</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n">toplabel</span> <span class="o">=</span> <span class="s">"Top Label"</span><span class="p">;</span> <span class="n">xlabel</span> <span class="o">=</span> <span class="s">"X-axis Label"</span><span class="p">;</span> <span class="n">ylabel</span> <span class="o">=</span> <span class="s">"Y-axis Label"</span><span class="p">;</span> <span class="n">legend</span> <span class="o">=</span> <span class="s">"Legend"</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a></span><span class="p">(</span><span class="o">&</span><span class="n">argc</span><span class="p">,</span><span class="o">&</span><span class="n">argv</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="n">help</span><span class="p">);</span><span class="k">if</span> <span class="p">(</span><span class="n">ierr</span><span class="p">)</span> <span class="k">return</span> <span class="n">ierr</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-x"</span><span class="p">,</span><span class="o">&</span><span class="n">x</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-y"</span><span class="p">,</span><span class="o">&</span><span class="n">y</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-width"</span><span class="p">,</span><span class="o">&</span><span class="n">width</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-height"</span><span class="p">,</span><span class="o">&</span><span class="n">height</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-n"</span><span class="p">,</span><span class="o">&</span><span class="n">n</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-nports"</span><span class="p">,</span><span class="o">&</span><span class="n">nports</span><span class="p">,</span><span class="o">&</span><span class="n">useports</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-nolegend"</span><span class="p">,</span><span class="o">&</span><span class="n">flg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">flg</span><span class="p">)</span> <span class="n">legend</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-notoplabel"</span><span class="p">,</span><span class="o">&</span><span class="n">flg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">flg</span><span class="p">)</span> <span class="n">toplabel</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-noxlabel"</span><span class="p">,</span><span class="o">&</span><span class="n">flg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">flg</span><span class="p">)</span> <span class="n">xlabel</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-noylabel"</span><span class="p">,</span><span class="o">&</span><span class="n">flg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">flg</span><span class="p">)</span> <span class="n">ylabel</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-nolabels"</span><span class="p">,</span><span class="o">&</span><span class="n">flg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">flg</span><span class="p">)</span> <span class="p">{</span><span class="n">toplabel</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;</span> <span class="n">xlabel</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;</span> <span class="n">ylabel</span> <span class="o">=</span> <span class="nb">NULL</span><span class="p">;}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawCreate.html#PetscDrawCreate">PetscDrawCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="s">"Title"</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">y</span><span class="p">,</span><span class="n">width</span><span class="p">,</span><span class="n">height</span><span class="p">,</span><span class="o">&</span><span class="n">draw</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawSetFromOptions.html#PetscDrawSetFromOptions">PetscDrawSetFromOptions</a></span><span class="p">(</span><span class="n">draw</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">useports</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawViewPortsCreate.html#PetscDrawViewPortsCreate">PetscDrawViewPortsCreate</a></span><span class="p">(</span><span class="n">draw</span><span class="p">,</span><span class="n">nports</span><span class="p">,</span><span class="o">&</span><span class="n">ports</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawViewPortsSet.html#PetscDrawViewPortsSet">PetscDrawViewPortsSet</a></span><span class="p">(</span><span class="n">ports</span><span class="p">,</span><span class="mi">0</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGCreate.html#PetscDrawLGCreate">PetscDrawLGCreate</a></span><span class="p">(</span><span class="n">draw</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="o">&</span><span class="n">lg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGSetUseMarkers.html#PetscDrawLGSetUseMarkers">PetscDrawLGSetUseMarkers</a></span><span class="p">(</span><span class="n">lg</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_TRUE.html#PETSC_TRUE">PETSC_TRUE</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGGetAxis.html#PetscDrawLGGetAxis">PetscDrawLGGetAxis</a></span><span class="p">(</span><span class="n">lg</span><span class="p">,</span><span class="o">&</span><span class="n">axis</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxisSetColors.html#PetscDrawAxisSetColors">PetscDrawAxisSetColors</a></span><span class="p">(</span><span class="n">axis</span><span class="p">,</span><span class="n">PETSC_DRAW_BLACK</span><span class="p">,</span><span class="n">PETSC_DRAW_RED</span><span class="p">,</span><span class="n">PETSC_DRAW_BLUE</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawAxisSetLabels.html#PetscDrawAxisSetLabels">PetscDrawAxisSetLabels</a></span><span class="p">(</span><span class="n">axis</span><span class="p">,</span><span class="n">toplabel</span><span class="p">,</span><span class="n">xlabel</span><span class="p">,</span><span class="n">ylabel</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGSetLegend.html#PetscDrawLGSetLegend">PetscDrawLGSetLegend</a></span><span class="p">(</span><span class="n">lg</span><span class="p">,</span><span class="o">&</span><span class="n">legend</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGSetFromOptions.html#PetscDrawLGSetFromOptions">PetscDrawLGSetFromOptions</a></span><span class="p">(</span><span class="n">lg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">for</span> <span class="p">(</span><span class="n">i</span><span class="o">=</span><span class="mi">0</span><span class="p">;</span> <span class="n">i</span><span class="o"><=</span><span class="n">n</span><span class="p">;</span> <span class="n">i</span><span class="o">++</span><span class="p">)</span> <span class="p">{</span>
<span class="n">xd</span> <span class="o">=</span> <span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span><span class="p">)(</span><span class="n">i</span> <span class="o">-</span> <span class="mi">5</span><span class="p">);</span> <span class="n">yd</span> <span class="o">=</span> <span class="n">xd</span><span class="o">*</span><span class="n">xd</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGAddPoint.html#PetscDrawLGAddPoint">PetscDrawLGAddPoint</a></span><span class="p">(</span><span class="n">lg</span><span class="p">,</span><span class="o">&</span><span class="n">xd</span><span class="p">,</span><span class="o">&</span><span class="n">yd</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGDraw.html#PetscDrawLGDraw">PetscDrawLGDraw</a></span><span class="p">(</span><span class="n">lg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGSave.html#PetscDrawLGSave">PetscDrawLGSave</a></span><span class="p">(</span><span class="n">lg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawViewPortsDestroy.html#PetscDrawViewPortsDestroy">PetscDrawViewPortsDestroy</a></span><span class="p">(</span><span class="n">ports</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawLGDestroy.html#PetscDrawLGDestroy">PetscDrawLGDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">lg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDrawDestroy.html#PetscDrawDestroy">PetscDrawDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">draw</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFinalize.html#PetscFinalize">PetscFinalize</a></span><span class="p">();</span>
<span class="k">return</span> <span class="n">ierr</span><span class="p">;</span>
<span class="p">}</span>
</pre></div>
</div>
</div>
</div>
<div class="section" id="graphical-convergence-monitor">
<h3>Graphical Convergence Monitor<a class="headerlink" href="#graphical-convergence-monitor" title="Permalink to this headline">¶</a></h3>
<p>For both the linear and nonlinear solvers default routines allow one to
graphically monitor convergence of the iterative method. These are
accessed via the command line with <code class="docutils literal notranslate"><span class="pre">-ksp_monitor_lg_residualnorm</span></code> and
<code class="docutils literal notranslate"><span class="pre">-snes_monitor_lg_residualnorm</span></code>. See also
<a class="reference internal" href="ksp.html#sec-kspmonitor"><span class="std std-ref">Convergence Monitoring</span></a> and <a class="reference internal" href="snes.html#sec-snesmonitor"><span class="std std-ref">Convergence Monitoring</span></a>.</p>
<p>The two functions used are <code class="docutils literal notranslate"><span class="pre">KSPMonitorLGResidualNorm()</span></code> and
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPMonitorLGResidualNormCreate.html#KSPMonitorLGResidualNormCreate">KSPMonitorLGResidualNormCreate</a>()</span></code>. These can easily be modified to
serve specialized needs.</p>
</div>
<div class="section" id="disabling-graphics-at-compile-time">
<h3>Disabling Graphics at Compile Time<a class="headerlink" href="#disabling-graphics-at-compile-time" title="Permalink to this headline">¶</a></h3>
<p>To disable all X-window-based graphics, run <code class="docutils literal notranslate"><span class="pre">./configure</span></code> with the
additional option <code class="docutils literal notranslate"><span class="pre">--with-x=0</span></code></p>
</div>
</div>
<div class="section" id="emacs-users">
<span id="sec-emacs"></span><h2>Emacs Users<a class="headerlink" href="#emacs-users" title="Permalink to this headline">¶</a></h2>
<p>Many PETSc developers use Emacs, which can be used as a “simple” text editor or a comprehensive development environment.
For a more integrated development environment, we recommend using <a class="reference external" href="https://emacs-lsp.github.io/lsp-mode/">lsp-mode</a> (or <a class="reference external" href="https://github.com/joaotavora/eglot">eglot</a>) with <a class="reference external" href="https://clangd.llvm.org/">clangd</a>.
The most convenient way to teach clangd what compilation flags to use is to install <a class="reference external" href="https://github.com/rizsotto/Bear">Bear</a> (“build ear”) and run:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">bear</span> <span class="n">make</span> <span class="o">-</span><span class="n">B</span>
</pre></div>
</div>
<p>which will do a complete rebuild (<code class="docutils literal notranslate"><span class="pre">-B</span></code>) of PETSc and capture the compilation commands in a file named <code class="docutils literal notranslate"><span class="pre">compile_commands.json</span></code>, which will be automatically picked up by clangd.
You can use the same procedure when building examples or your own project.
It can also be used with any other editor that supports clangd, including VS Code and Vim.
When lsp-mode is accompanied by <a class="reference external" href="https://www.flycheck.org/en/latest/">flycheck</a>, Emacs will provide real-time feedback and syntax checking, along with refactoring tools provided by clangd.</p>
<p>The easiest way to install packages in recent Emacs is to use the “Options” menu to select “Manage Emacs Packages”.</p>
<div class="section" id="tags">
<h3>Tags<a class="headerlink" href="#tags" title="Permalink to this headline">¶</a></h3>
<p>It is sometimes useful to cross-reference tags across projects.
Regardless of whether you use lsp-mode, it can be useful to use <a class="reference external" href="https://www.gnu.org/software/global/">GNU Global</a> (install <code class="docutils literal notranslate"><span class="pre">gtags</span></code>) to provide reverse lookups (e.g. find all call sites
for a given function) across all projects you might work on/browse.
Tags for PETSc can be generated by running <code class="docutils literal notranslate"><span class="pre">make</span> <span class="pre">allgtags</span></code> from <code class="docutils literal notranslate"><span class="pre">${PETSC_DIR}</span></code>, or one can generate tags for all projects by running a command such as</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>find $PETSC_DIR/{include,src,tutorials,$PETSC_ARCH/include} any/other/paths \
-regex '.*\.\(cc\|hh\|cpp\|cxx\|C\|hpp\|c\|h\|cu\)$' \
| grep -v ftn-auto | gtags -f -
</pre></div>
</div>
<p>from your home directory or wherever you keep source code. If you are
making large changes, it is useful to either set this up to run as a
cron job or to make a convenient alias so that refreshing is easy. Then
add the following to <code class="docutils literal notranslate"><span class="pre">~/.emacs</span></code> to enable gtags and specify key bindings.</p>
<div class="highlight-emacs notranslate"><div class="highlight"><pre><span></span><span class="p">(</span><span class="nb">when</span> <span class="p">(</span><span class="nb">require</span> <span class="ss">'gtags</span><span class="p">)</span>
<span class="p">(</span><span class="nv">global-set-key</span> <span class="p">(</span><span class="nv">kbd</span> <span class="s">"C-c f"</span><span class="p">)</span> <span class="ss">'gtags-find-file</span><span class="p">)</span>
<span class="p">(</span><span class="nv">global-set-key</span> <span class="p">(</span><span class="nv">kbd</span> <span class="s">"C-c ."</span><span class="p">)</span> <span class="ss">'gtags-find-tag</span><span class="p">)</span>
<span class="p">(</span><span class="nv">global-set-key</span> <span class="p">(</span><span class="nv">kbd</span> <span class="s">"C-c r"</span><span class="p">)</span> <span class="ss">'gtags-find-rtag</span><span class="p">)</span>
<span class="p">(</span><span class="nv">global-set-key</span> <span class="p">(</span><span class="nv">kbd</span> <span class="s">"C-c ,"</span><span class="p">)</span> <span class="ss">'gtags-pop-stack</span><span class="p">))</span>
<span class="p">(</span><span class="nv">add-hook</span> <span class="ss">'c-mode-common-hook</span>
<span class="o">'</span><span class="p">(</span><span class="nb">lambda</span> <span class="p">()</span> <span class="p">(</span><span class="nv">gtags-mode</span> <span class="no">t</span><span class="p">)))</span> <span class="c1">; Or add to existing hook</span>
</pre></div>
</div>
<p>A more basic alternative to the GNU Global (<code class="docutils literal notranslate"><span class="pre">gtags</span></code>) approach that does not require adding packages is to use
the builtin <code class="docutils literal notranslate"><span class="pre">etags</span></code> feature. First, run <code class="docutils literal notranslate"><span class="pre">make</span> <span class="pre">alletags</span></code> from the
PETSc home directory to generate the file <code class="docutils literal notranslate"><span class="pre">${PETSC_DIR}/TAGS</span></code>, and
then from within Emacs, run</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>M-x visit-tags-table
</pre></div>
</div>
<p>where <code class="docutils literal notranslate"><span class="pre">M</span></code> denotes the Emacs Meta key, and enter the name of the
<code class="docutils literal notranslate"><span class="pre">TAGS</span></code> file. Then the command <code class="docutils literal notranslate"><span class="pre">M-.</span></code> will cause Emacs to find the
file and line number where a desired PETSc function is defined. Any
string in any of the PETSc files can be found with the command <code class="docutils literal notranslate"><span class="pre">M-x</span> <span class="pre">tags-search</span></code>.
To find repeated occurrences, one can simply use <code class="docutils literal notranslate"><span class="pre">M-,</span></code> to find the next occurrence.</p>
</div>
</div>
<div class="section" id="vs-code-users">
<h2>VS Code Users<a class="headerlink" href="#vs-code-users" title="Permalink to this headline">¶</a></h2>
<p><a class="reference external" href="https://code.visualstudio.com/">VS Code</a> (unlike <a class="reference internal" href="#sec-visual-studio"><span class="std std-ref">Visual Studio Users</span></a>, described below) is an open source editor with a rich extension ecosystem.
It has <a class="reference external" href="https://marketplace.visualstudio.com/items?itemName=llvm-vs-code-extensions.vscode-clangd">excellent integration</a> with clangd and will automatically pick up <code class="docutils literal notranslate"><span class="pre">compile_commands.json</span></code> as produced by a command such as <code class="docutils literal notranslate"><span class="pre">bear</span> <span class="pre">make</span> <span class="pre">-B</span></code> (see <a class="reference internal" href="#sec-emacs"><span class="std std-ref">Emacs Users</span></a>).
If you have no prior attachment to a specific code editor, we recommend trying VS Code.</p>
</div>
<div class="section" id="vi-and-vim-users">
<h2>Vi and Vim Users<a class="headerlink" href="#vi-and-vim-users" title="Permalink to this headline">¶</a></h2>
<p>See the <a class="reference internal" href="#sec-emacs"><span class="std std-ref">Emacs Users</span></a> discussion above for configuration of clangd, which provides integrated development environment.</p>
<p>If users develop application codes using Vi or Vim the <code class="docutils literal notranslate"><span class="pre">tags</span></code> feature
can be used to search PETSc files quickly and efficiently. To use this
feature, one should first check if the file, <code class="docutils literal notranslate"><span class="pre">${PETSC_DIR}/CTAGS</span></code>
exists. If this file is not present, it should be generated by running
<code class="docutils literal notranslate"><span class="pre">make</span> <span class="pre">alletags</span></code> from the PETSc home directory. Once the file
exists, from Vi/Vim the user should issue the command</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>:set tags=CTAGS
</pre></div>
</div>
<p>from the <code class="docutils literal notranslate"><span class="pre">PETSC_DIR</span></code> directory and enter the name of the <code class="docutils literal notranslate"><span class="pre">CTAGS</span></code>
file. Then the command “tag functionname” will cause Vi/Vim to find the
file and line number where a desired PETSc function is defined. See, <a class="reference external" href="http://www.yolinux.com/TUTORIALS/LinuxTutorialAdvanced_vi.html">online tutorials</a>
for additional Vi/Vim options that allow searches, etc. It is also
possible to use GNU Global with Vim; see the description for Emacs
above.</p>
</div>
<div class="section" id="eclipse-users">
<h2>Eclipse Users<a class="headerlink" href="#eclipse-users" title="Permalink to this headline">¶</a></h2>
<p>If you are interested in developing code that uses PETSc from Eclipse or
developing PETSc in Eclipse and have knowledge of how to do indexing and
build libraries in Eclipse, please contact us at
<a class="reference external" href="mailto:petsc-dev%40mcs.anl.gov">petsc-dev<span>@</span>mcs<span>.</span>anl<span>.</span>gov</a>.</p>
<p>One way to index and build PETSc in Eclipse is as follows.</p>
<ol class="arabic simple">
<li><p>Open
“File<span class="math">\(\rightarrow\)</span>Import<span class="math">\(\rightarrow\)</span>Git<span class="math">\(\rightarrow\)</span>Projects
from Git”. In the next two panels, you can either add your existing
local repository or download PETSc from Bitbucket by providing the
URL. Most Eclipse distributions come with Git support. If not,
install the EGit plugin. When importing the project, select the
wizard “Import as general project”.</p></li>
<li><p>Right-click on the project (or the “File” menu on top) and select
“New <span class="math">\(\rightarrow\)</span> Convert to a C/C++ Project (Adds C/C++
Nature)”. In the setting window, choose “C Project” and specify the
project type as “Shared Library”.</p></li>
<li><p>Right-click on the C project and open the “Properties” panel. Under
“C/C++ Build <span class="math">\(\rightarrow\)</span> Builder Settings”, set the Build
directory to <code class="docutils literal notranslate"><span class="pre">PETSC_DIR</span></code> and make sure “Generate Makefiles
automatically” is unselected. Under the section “C/C++
General<span class="math">\(\rightarrow\)</span>Paths and Symbols”, add the PETSc paths
to “Includes”.</p></li>
</ol>
<blockquote>
<div><div class="highlight-none notranslate"><div class="highlight"><pre><span></span> ${PETSC_DIR}/include
${PETSC_DIR}/${PETSC_ARCH}/include
Under the section “C/C++ General\ :math:`\rightarrow`\ index”, choose
“Use active build configuration”.
</pre></div>
</div>
</div></blockquote>
<ol class="arabic simple">
<li><p>Configure PETSc normally outside Eclipse to generate a makefile and
then build the project in Eclipse. The source code will be parsed by
Eclipse.</p></li>
</ol>
<p>If you launch Eclipse from the Dock on Mac OS X, <code class="docutils literal notranslate"><span class="pre">.bashrc</span></code> will not be
loaded (a known OS X behavior, for security reasons). This will be a
problem if you set the environment variables <code class="docutils literal notranslate"><span class="pre">PETSC_DIR</span></code> and
<code class="docutils literal notranslate"><span class="pre">PETSC_ARCH</span></code> in <code class="docutils literal notranslate"><span class="pre">.bashrc</span></code>. A solution which involves replacing the
executable can be found at
<code class="docutils literal notranslate"><span class="pre">`/questions/829749/launch-mac-eclipse-with-environment-variables-set</span></code> </questions/829749/launch-mac-eclipse-with-environment-variables-set>`__.
Alternatively, you can add <code class="docutils literal notranslate"><span class="pre">PETSC_DIR</span></code> and <code class="docutils literal notranslate"><span class="pre">PETSC_ARCH</span></code> manually
under “Properties <span class="math">\(\rightarrow\)</span> C/C++ Build <span class="math">\(\rightarrow\)</span>
Environment”.</p>
<p>To allow an Eclipse code to compile with the PETSc include files and
link with the PETSc libraries, a PETSc user has suggested the following.</p>
<ol class="arabic simple">
<li><p>Right-click on your C project and select “Properties
<span class="math">\(\rightarrow\)</span> C/C++ Build <span class="math">\(\rightarrow\)</span> Settings”</p></li>
<li><p>A new window on the righthand side appears with various settings
options. Select “Includes” and add the required PETSc paths,</p></li>
</ol>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>${PETSC_DIR}/include
${PETSC_DIR}/${PETSC_ARCH}/include
</pre></div>
</div>
<ol class="arabic simple">
<li><p>Select “Libraries” under the header Linker and set the library search
path:</p></li>
</ol>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span> ${PETSC_DIR}/${PETSC_ARCH}/lib
and the libraries, for example
</pre></div>
</div>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>m, petsc, stdc++, mpichxx, mpich, lapack, blas, gfortran, dl, rt,gcc_s, pthread, X11
</pre></div>
</div>
<p>Another PETSc user has provided the following steps to build an Eclipse
index for PETSc that can be used with their own code, without compiling
PETSc source into their project.</p>
<ol class="arabic simple">
<li><p>In the user project source directory, create a symlink to the PETSC
<code class="docutils literal notranslate"><span class="pre">src/</span></code> directory.</p></li>
<li><p>Refresh the project explorer in Eclipse, so the new symlink is
followed.</p></li>
<li><p>Right-click on the project in the project explorer, and choose “Index
<span class="math">\(\rightarrow\)</span> Rebuild”. The index should now be built.</p></li>
<li><p>Right-click on the PETSc symlink in the project explorer, and choose
“Exclude from build…” to make sure Eclipse does not try to compile
PETSc with the project.</p></li>
</ol>
<p>For further examples of using Eclipse with a PETSc-based application,
see the documentation for LaMEM <a class="footnote-reference brackets" href="#id7" id="id4">10</a>.</p>
</div>
<div class="section" id="qt-creator-users">
<h2>Qt Creator Users<a class="headerlink" href="#qt-creator-users" title="Permalink to this headline">¶</a></h2>
<p>This information was provided by Mohammad Mirzadeh. The Qt Creator IDE
is part of the Qt SDK, developed for cross-platform GUI programming
using C++. It is available under GPL v3, LGPL v2 and a commercial
license and may be obtained, either as part of the Qt SDK or as
stand-alone software. It supports
automatic makefile generation using cross-platform <code class="docutils literal notranslate"><span class="pre">qmake</span></code> and
<code class="docutils literal notranslate"><span class="pre">cmake</span></code> build systems as well as allowing one to import projects based
on existing, possibly hand-written, makefiles. Qt Creator has a visual
debugger using GDB and LLDB (on Linux and OS X) or Microsoft’s CDB (on
Windows) as backends. It also has an interface to Valgrind’s “memcheck”
and “callgrind” tools to detect memory leaks and profile code. It has
built-in support for a variety of version control systems including git,
mercurial, and subversion. Finally, Qt Creator comes fully equipped with
auto-completion, function look-up, and code refactoring tools. This
enables one to easily browse source files, find relevant functions, and
refactor them across an entire project.</p>
<div class="section" id="creating-a-project">
<h3>Creating a Project<a class="headerlink" href="#creating-a-project" title="Permalink to this headline">¶</a></h3>
<p>When using Qt Creator with <code class="docutils literal notranslate"><span class="pre">qmake</span></code>, one needs a <code class="docutils literal notranslate"><span class="pre">.pro</span></code> file. This
configuration file tells Qt Creator about all build/compile options and
locations of source files. One may start with a blank <code class="docutils literal notranslate"><span class="pre">.pro</span></code> file and
fill in configuration options as needed. For example:</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span># The name of the application executable
TARGET = ex1
# There are two ways to add PETSc functionality
# 1-Manual: Set all include path and libs required by PETSc
PETSC_INCLUDE = path/to/petsc_includes # e.g. obtained via running `make getincludedirs'
PETSC_LIBS = path/to/petsc_libs # e.g. obtained via running `make getlinklibs'
INCLUDEPATH += $$PETSC_INCLUDES
LIBS += $$PETSC_LIBS
# 2-Automatic: Use the PKGCONFIG funtionality
# NOTE: PETSc.pc must be in the pkgconfig path. You might need to adjust PKG_CONFIG_PATH
CONFIG += link_pkgconfig
PKGCONFIG += PETSc
# Set appropriate compiler and its flags
QMAKE_CC = path/to/mpicc
QMAKE_CXX = path/to/mpicxx # if this is a cpp project
QMAKE_LINK = path/to/mpicxx # if this is a cpp project
QMAKE_CFLAGS += -O3 # add extra flags here
QMAKE_CXXFLAGS += -O3
QMAKE_LFLAGS += -O3
# Add all files that must be compiled
SOURCES += ex1.c source1.c source2.cpp
HEADERS += source1.h source2.h
# OTHER_FILES are ignored during compilation but will be shown in file panel in Qt Creator
OTHER_FILES += \
path/to/resource_file \
path/to/another_file
</pre></div>
</div>
<p>In this example, keywords include:</p>
<ul class="simple">
<li><p><code class="docutils literal notranslate"><span class="pre">TARGET</span></code>: The name of the application executable.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">INCLUDEPATH</span></code>: Used at compile time to point to required include
files. Essentially, it is used as an <code class="docutils literal notranslate"><span class="pre">-I</span> <span class="pre">\$\$INCLUDEPATH</span></code> flag for
the compiler. This should include all application-specific header
files and those related to PETSc (which may be found via
<code class="docutils literal notranslate"><span class="pre">make</span> <span class="pre">getincludedirs</span></code>).</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">LIBS</span></code>: Defines all required external libraries to link with the
application. To get PETSc’s linking libraries, use
<code class="docutils literal notranslate"><span class="pre">make</span> <span class="pre">getlinklibs</span></code>.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">CONFIG</span></code>: Configuration options to be used by <code class="docutils literal notranslate"><span class="pre">qmake</span></code>. Here, the
option <code class="docutils literal notranslate"><span class="pre">link_pkgconfig</span></code> instructs <code class="docutils literal notranslate"><span class="pre">qmake</span></code> to internally use
<code class="docutils literal notranslate"><span class="pre">pkgconfig</span></code> to resolve <code class="docutils literal notranslate"><span class="pre">INCLUDEPATH</span></code> and <code class="docutils literal notranslate"><span class="pre">LIBS</span></code> variables.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">PKGCONFIG</span></code>: Name of the configuration file (the <code class="docutils literal notranslate"><span class="pre">.pc</span></code> file –
here <code class="docutils literal notranslate"><span class="pre">PETSc.pc</span></code>) to be passed to <code class="docutils literal notranslate"><span class="pre">pkgconfig</span></code>. Note that for this
functionality to work, <code class="docutils literal notranslate"><span class="pre">PETSc.pc</span></code> must be in path which might
require adjusting the <code class="docutils literal notranslate"><span class="pre">PKG_CONFIG_PATH</span></code> enviroment variable. For
more information see
<a class="reference external" href="https://doc.qt.io/qtcreator/creator-build-settings.html">the Qt Creator documentation</a>.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">QMAKE_CC</span></code> and <code class="docutils literal notranslate"><span class="pre">QMAKE_CXX</span></code>: Define which C/C++ compilers use.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">QMAKE_LINK</span></code>: Defines the proper linker to be used. Relevant if
compiling C++ projects.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">QMAKE_CFLAGS</span></code>, <code class="docutils literal notranslate"><span class="pre">QMAKE_CXXFLAGS</span></code> and <code class="docutils literal notranslate"><span class="pre">QMAKE_LFLAGS</span></code>: Set the
corresponding compile and linking flags.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">SOURCES</span></code>: Source files to be compiled.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">HEADERS</span></code>: Header files required by the application.</p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">OTHER_FILES</span></code>: Other files to include (source, header, or any other
extension). Note that none of the source files placed here are
compiled.</p></li>
</ul>
<p>More options can be included in a <code class="docutils literal notranslate"><span class="pre">.pro</span></code> file; see
<a class="reference external" href="https://doc.qt.io/qt-5/qmake-project-files.html">https://doc.qt.io/qt-5/qmake-project-files.html</a>. Once the <code class="docutils literal notranslate"><span class="pre">.pro</span></code> file
is generated, the user can simply open it via Qt Creator. Upon opening,
one has the option to create two different build options, debug and
release, and switch between the two. For more information on using the
Qt Creator interface and other more advanced aspects of the IDE, refer
to <a class="reference external" href="https://www.qt.io/qt-features-libraries-apis-tools-and-ide/">https://www.qt.io/qt-features-libraries-apis-tools-and-ide/</a></p>
</div>
</div>
<div class="section" id="visual-studio-users">
<span id="sec-visual-studio"></span><h2>Visual Studio Users<a class="headerlink" href="#visual-studio-users" title="Permalink to this headline">¶</a></h2>
<p>To use PETSc from MS Visual Studio, one would have to compile a PETSc
example with its corresponding makefile and then transcribe all compiler
and linker options used in this build into a Visual Studio project file,
in the appropriate format in Visual Studio project settings.</p>
</div>
<div class="section" id="xcode-users-the-apple-gui-development-system">
<h2>XCode Users (The Apple GUI Development System)<a class="headerlink" href="#xcode-users-the-apple-gui-development-system" title="Permalink to this headline">¶</a></h2>
<div class="section" id="mac-os-x">
<h3>Mac OS X<a class="headerlink" href="#mac-os-x" title="Permalink to this headline">¶</a></h3>
<p>Follow the instructions in <code class="docutils literal notranslate"><span class="pre">$PETSC_DIR/systems/Apple/OSX/bin/makeall</span></code>
to build the PETSc framework and documentation suitable for use in
XCode.</p>
<p>You can then use the PETSc framework in
<code class="docutils literal notranslate"><span class="pre">$PETSC_DIR/arch-osx/PETSc.framework</span></code> in the usual manner for Apple
frameworks. See the examples in
<code class="docutils literal notranslate"><span class="pre">$PETSC_DIR/systems/Apple/OSX/examples</span></code>. When working in XCode, things
like function name completion should work for all PETSc functions as
well as MPI functions. You must also link against the Apple
<code class="docutils literal notranslate"><span class="pre">Accelerate.framework</span></code>.</p>
</div>
<div class="section" id="iphone-ipad-ios">
<h3>iPhone/iPad iOS<a class="headerlink" href="#iphone-ipad-ios" title="Permalink to this headline">¶</a></h3>
<p>Follow the instructions in
<code class="docutils literal notranslate"><span class="pre">$PETSC_DIR/systems/Apple/iOS/bin/iosbuilder.py</span></code> to build the PETSc
library for use on the iPhone/iPad.</p>
<p>You can then use the PETSc static library in
<code class="docutils literal notranslate"><span class="pre">$PETSC_DIR/arch-osx/libPETSc.a</span></code> in the usual manner for Apple
libraries inside your iOS XCode projects; see the examples in
<code class="docutils literal notranslate"><span class="pre">$PETSC_DIR/systems/Apple/iOS/examples</span></code>. You must also link against
the Apple <code class="docutils literal notranslate"><span class="pre">Accelerate.framework</span></code>.</p>
<dl class="footnote brackets">
<dt class="label" id="id5"><span class="brackets"><a class="fn-backref" href="#id1">8</a></span></dt>
<dd><p><a class="reference external" href="https://bitbucket.org/saws/saws/wiki/Home">Saws wiki on Bitbucket</a></p>
</dd>
<dt class="label" id="id6"><span class="brackets"><a class="fn-backref" href="#id2">9</a></span></dt>
<dd><p>Note that this option is not required to use PETSc with C++</p>
</dd>
<dt class="label" id="id7"><span class="brackets"><a class="fn-backref" href="#id4">10</a></span></dt>
<dd><p><code class="docutils literal notranslate"><span class="pre">doc/</span></code> at <a class="reference external" href="https://bitbucket.org/bkaus/lamem">https://bitbucket.org/bkaus/lamem</a></p>
</dd>
</dl>
</div>
</div>
</div>
</div>
</div>
</div>
<div class="clearer"></div>
</div>
<div class="related" role="navigation" aria-label="related navigation">
<h3>Navigation</h3>
<ul>
<li class="right" style="margin-right: 10px">
<a href="../genindex.html" title="General Index"
>index</a></li>
<li class="right" >
<a href="advanced.html" title="Unimportant and Advanced Features of Matrices and Solvers"
>next</a> |</li>
<li class="right" >
<a href="performance.html" title="Hints for Performance Tuning"
>previous</a> |</li>
<li class="nav-item nav-item-0"><a href="../index.html">PETSc 3.14.5 documentation</a> »</li>
<li class="nav-item nav-item-1"><a href="index.html" >PETSc Users Manual</a> »</li>
<li class="nav-item nav-item-2"><a href="additional.html" >Additional Information</a> »</li>
</ul>
</div>
<div class="footer" role="contentinfo">
© Copyright 1991-2021, UChicago Argonne, LLC and the PETSc Development Team.
Created using <a href="http://sphinx-doc.org/">Sphinx</a> 2.4.4.
</div>
</body>
</html>
|