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
|
<!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/fortran.html" />
<meta charset="utf-8" />
<title>PETSc for Fortran Users — 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="Using MATLAB with PETSc" href="matlab.html" />
<link rel="prev" title="Additional Information" href="additional.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/fortran.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="matlab.html" title="Using MATLAB with PETSc"
accesskey="N">next</a> |</li>
<li class="right" >
<a href="additional.html" title="Additional Information"
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="#">PETSc for Fortran Users</a><ul>
<li><a class="reference internal" href="#c-vs-fortran-interfaces">C vs. Fortran Interfaces</a><ul>
<li><a class="reference internal" href="#fortran-include-files">Fortran Include Files</a></li>
<li><a class="reference internal" href="#error-checking">Error Checking</a></li>
<li><a class="reference internal" href="#calling-fortran-routines-from-c-and-c-routines-from-fortran">Calling Fortran Routines from C (and C Routines from Fortran)</a></li>
<li><a class="reference internal" href="#passing-null-pointers">Passing Null Pointers</a></li>
<li><a class="reference internal" href="#duplicating-multiple-vectors">Duplicating Multiple Vectors</a></li>
<li><a class="reference internal" href="#matrix-vector-and-is-indices">Matrix, Vector and IS Indices</a></li>
<li><a class="reference internal" href="#setting-routines">Setting Routines</a></li>
<li><a class="reference internal" href="#compiling-and-linking-fortran-programs">Compiling and Linking Fortran Programs</a></li>
<li><a class="reference internal" href="#routines-with-different-fortran-interfaces">Routines with Different Fortran Interfaces</a></li>
</ul>
</li>
<li><a class="reference internal" href="#sample-fortran-programs">Sample Fortran Programs</a><ul>
<li><a class="reference internal" href="#array-arguments">Array Arguments</a></li>
</ul>
</li>
</ul>
</li>
</ul>
<h4>Previous topic</h4>
<p class="topless"><a href="additional.html"
title="previous chapter">Additional Information</a></p>
<h4>Next topic</h4>
<p class="topless"><a href="matlab.html"
title="next chapter">Using MATLAB with PETSc</a></p>
<div role="note" aria-label="source link">
<h3>This Page</h3>
<ul class="this-page-menu">
<li><a href="../_sources/manual/fortran.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="petsc-for-fortran-users">
<span id="chapter-fortran"></span><h1>PETSc for Fortran Users<a class="headerlink" href="#petsc-for-fortran-users" title="Permalink to this headline">¶</a></h1>
<p>Most of the functionality of PETSc can be obtained by people who program
purely in Fortran.</p>
<div class="section" id="c-vs-fortran-interfaces">
<h2>C vs. Fortran Interfaces<a class="headerlink" href="#c-vs-fortran-interfaces" title="Permalink to this headline">¶</a></h2>
<p>Only a few differences exist between the C and Fortran PETSc interfaces,
are due to Fortran syntax differences. All Fortran routines have the
same names as the corresponding C versions, and PETSc command line
options are fully supported. The routine arguments follow the usual
Fortran conventions; the user need not worry about passing pointers or
values. The calling sequences for the Fortran version are in most cases
identical to the C version, except for the error checking variable
discussed in <a class="reference internal" href="#sec-fortran-errors"><span class="std std-ref">Error Checking</span></a> and a few routines
listed in <a class="reference internal" href="#sec-fortran-exceptions"><span class="std std-ref">Routines with Different Fortran Interfaces</span></a>.</p>
<div class="section" id="fortran-include-files">
<span id="sec-fortran-includes"></span><h3>Fortran Include Files<a class="headerlink" href="#fortran-include-files" title="Permalink to this headline">¶</a></h3>
<p>The Fortran include files for PETSc are located in the directory
<code class="docutils literal notranslate"><span class="pre">${PETSC_DIR}/include/petsc/finclude</span></code> and should be used via
statements such as the following:</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="cp">#include <petsc/finclude/petscXXX.h></span>
</pre></div>
</div>
<p>for example,</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="cp">#include <petsc/finclude/petscksp.h></span>
</pre></div>
</div>
<p>You must also use the appropriate Fortran module which is done with</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="k">use </span><span class="n">petscXXX</span>
</pre></div>
</div>
<p>for example,</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="k">use </span><span class="n">petscksp</span>
</pre></div>
</div>
</div>
<div class="section" id="error-checking">
<span id="sec-fortran-errors"></span><h3>Error Checking<a class="headerlink" href="#error-checking" title="Permalink to this headline">¶</a></h3>
<p>In the Fortran version, each PETSc routine has as its final argument an
integer error variable, in contrast to the C convention of providing the
error variable as the routine’s return value. The error code is set to
be nonzero if an error has been detected; otherwise, it is zero. For
example, the Fortran and C variants of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSolve.html#KSPSolve">KSPSolve</a>()</span></code> are given,
respectively, below, where <code class="docutils literal notranslate"><span class="pre">ierr</span></code> denotes the error variable:</p>
<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/KSP/KSPSolve.html#KSPSolve">KSPSolve</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="n">b</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span> <span class="c">! Fortran</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/KSP/KSPSolve.html#KSPSolve">KSPSolve</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="n">b</span><span class="p">,</span><span class="n">x</span><span class="p">);</span> <span class="o">/*</span> <span class="n">C</span> <span class="o">*/</span>
</pre></div>
</div>
<p>Fortran programmers can check these error codes with <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>(ierr)</span></code>,
which terminates all processes when an error is encountered. Likewise,
one can set error codes within Fortran programs by using
<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>(comm,p,'</span> <span class="pre">',ierr)</span></code>, which again terminates all processes upon
detection of an error. Note that complete error tracebacks with
<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> and <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>, as described in
<a class="reference internal" href="getting_started.html#sec-simple"><span class="std std-ref">Simple PETSc Examples</span></a> for C routines, are <em>not</em> directly supported for
Fortran routines; however, Fortran programmers can easily use the error
codes in writing their own tracebacks. For example, one could use code
such as the following:</p>
<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/KSP/KSPSolve.html#KSPSolve">KSPSolve</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="n">b</span><span class="p">,</span><span class="n">x</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">ierr</span> <span class="p">.</span><span class="n">ne</span><span class="p">.</span> <span class="mi">0</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> print</span><span class="o">*</span><span class="p">,</span> <span class="s1">'Error in routine ...'</span>
<span class="k">return</span>
<span class="k">end if</span>
</pre></div>
</div>
</div>
<div class="section" id="calling-fortran-routines-from-c-and-c-routines-from-fortran">
<h3>Calling Fortran Routines from C (and C Routines from Fortran)<a class="headerlink" href="#calling-fortran-routines-from-c-and-c-routines-from-fortran" title="Permalink to this headline">¶</a></h3>
<p>Different machines have different methods of naming Fortran routines
called from C (or C routines called from Fortran). Most Fortran
compilers change all the capital letters in Fortran routines to
lowercase. On some machines, the Fortran compiler appends an underscore
to the end of each Fortran routine name; for example, the Fortran
routine <code class="docutils literal notranslate"><span class="pre">Dabsc()</span></code> would be called from C with <code class="docutils literal notranslate"><span class="pre">dabsc_()</span></code>. Other
machines change all the letters in Fortran routine names to capitals.</p>
<p>PETSc provides two macros (defined in C/C++) to help write portable code
that mixes C/C++ and Fortran. They are <code class="docutils literal notranslate"><span class="pre">PETSC_HAVE_FORTRAN_UNDERSCORE</span></code>
and <code class="docutils literal notranslate"><span class="pre">PETSC_HAVE_FORTRAN_CAPS</span></code> , which are defined in the file
<code class="docutils literal notranslate"><span class="pre">${PETSC_DIR}/${PETSC_ARCH}/include/petscconf.h</span></code>. The macros are used,
for example, as follows:</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="cp">#if defined(PETSC_HAVE_FORTRAN_CAPS)</span>
<span class="cp">#define dabsc_ DMDABSC</span>
<span class="cp">#elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)</span>
<span class="cp">#define dabsc_ dabsc</span>
<span class="cp">#endif</span>
<span class="p">.....</span>
<span class="n">dabsc_</span><span class="p">(</span> <span class="p">&</span><span class="n">n</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="o">/*</span> <span class="k">call </span><span class="n">the</span> <span class="n">Fortran</span> <span class="k">function</span> <span class="o">*/</span>
</pre></div>
</div>
</div>
<div class="section" id="passing-null-pointers">
<h3>Passing Null Pointers<a class="headerlink" href="#passing-null-pointers" title="Permalink to this headline">¶</a></h3>
<p>In several PETSc C functions, one has the option of passing a NULL (0)
argument (for example, the fifth argument of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateSeqAIJ.html#MatCreateSeqAIJ">MatCreateSeqAIJ</a>()</span></code>).
From Fortran, users <em>must</em> pass <code class="docutils literal notranslate"><span class="pre">PETSC_NULL_XXX</span></code> to indicate a null
argument (where <code class="docutils literal notranslate"><span class="pre">XXX</span></code> is <code class="docutils literal notranslate"><span class="pre">INTEGER</span></code>, <code class="docutils literal notranslate"><span class="pre">DOUBLE</span></code>, <code class="docutils literal notranslate"><span class="pre">CHARACTER</span></code>, or
<code class="docutils literal notranslate"><span class="pre">SCALAR</span></code> depending on the type of argument required); passing 0 from
Fortran will crash the code. Note that the C convention of passing NULL
(or 0) <em>cannot</em> be used. For example, when no options prefix is desired
in the routine <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a>()</span></code>, one must use the following
command in Fortran:</p>
<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/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="n">PETSC_NULL_OPTIONS</span><span class="p">,</span><span class="n">PETSC_NULL_CHARACTER</span><span class="p">,</span><span class="n">PETSC_NULL_CHARACTER</span><span class="p">,</span><span class="s1">'-name'</span><span class="p">,</span><span class="n">N</span><span class="p">,</span><span class="n">flg</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
</pre></div>
</div>
<p>This Fortran requirement is inconsistent with C, where the user can
employ <code class="docutils literal notranslate"><span class="pre">NULL</span></code> for all null arguments.</p>
</div>
<div class="section" id="duplicating-multiple-vectors">
<span id="sec-fortvecd"></span><h3>Duplicating Multiple Vectors<a class="headerlink" href="#duplicating-multiple-vectors" title="Permalink to this headline">¶</a></h3>
<p>The Fortran interface to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicateVecs.html#VecDuplicateVecs">VecDuplicateVecs</a>()</span></code> differs slightly from
the C/C++ variant because Fortran does not allow conventional arrays to
be returned in routine arguments. To create <code class="docutils literal notranslate"><span class="pre">n</span></code> vectors of the same
format as an existing vector, the user must declare a vector array,
<code class="docutils literal notranslate"><span class="pre">v_new</span></code> of size <code class="docutils literal notranslate"><span class="pre">n</span></code>. Then, after <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicateVecs.html#VecDuplicateVecs">VecDuplicateVecs</a>()</span></code> has been
called, <code class="docutils literal notranslate"><span class="pre">v_new</span></code> will contain (pointers to) the new PETSc vector
objects. When finished with the vectors, the user should destroy them by
calling <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroyVecs.html#VecDestroyVecs">VecDestroyVecs</a>()</span></code>. For example, the following code fragment
duplicates <code class="docutils literal notranslate"><span class="pre">v_old</span></code> to form two new vectors, <code class="docutils literal notranslate"><span class="pre">v_new(1)</span></code> and
<code class="docutils literal notranslate"><span class="pre">v_new(2)</span></code>.</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></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="n">v_old</span><span class="p">,</span> <span class="n">v_new</span><span class="p">(</span><span class="mi">2</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">ierr</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">alpha</span>
<span class="p">....</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicateVecs.html#VecDuplicateVecs">VecDuplicateVecs</a></span><span class="p">(</span><span class="n">v_old</span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="n">v_new</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="n">alpha</span> <span class="o">=</span> <span class="mf">4.3</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n">v_new</span><span class="p">(</span><span class="mi">1</span><span class="p">),</span><span class="n">alpha</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="n">alpha</span> <span class="o">=</span> <span class="mf">6.0</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n">v_new</span><span class="p">(</span><span class="mi">2</span><span class="p">),</span><span class="n">alpha</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="p">....</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroyVecs.html#VecDestroyVecs">VecDestroyVecs</a></span><span class="p">(</span><span class="mi">2</span><span class="p">,</span> <span class="p">&</span><span class="n">v_new</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
</pre></div>
</div>
</div>
<div class="section" id="matrix-vector-and-is-indices">
<h3>Matrix, Vector and IS Indices<a class="headerlink" href="#matrix-vector-and-is-indices" title="Permalink to this headline">¶</a></h3>
<p>All matrices, vectors and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span></code> in PETSc use zero-based indexing,
regardless of whether C or Fortran is being used. The interface
routines, such as <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a>()</span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</a>()</span></code>, always use
zero indexing. See <a class="reference internal" href="mat.html#sec-matoptions"><span class="std std-ref">Basic Matrix Operations</span></a> for further
details.</p>
</div>
<div class="section" id="setting-routines">
<h3>Setting Routines<a class="headerlink" href="#setting-routines" title="Permalink to this headline">¶</a></h3>
<p>When a function pointer is passed as an argument to a PETSc function,
such as the test in <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetConvergenceTest.html#KSPSetConvergenceTest">KSPSetConvergenceTest</a>()</span></code>, it is assumed that this
pointer references a routine written in the same language as the PETSc
interface function that was called. For instance, if
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetConvergenceTest.html#KSPSetConvergenceTest">KSPSetConvergenceTest</a>()</span></code> is called from C, the test argument is
assumed to be a C function. Likewise, if it is called from Fortran, the
test is assumed to be written in Fortran.</p>
</div>
<div class="section" id="compiling-and-linking-fortran-programs">
<span id="sec-fortcompile"></span><h3>Compiling and Linking Fortran Programs<a class="headerlink" href="#compiling-and-linking-fortran-programs" title="Permalink to this headline">¶</a></h3>
<p>See <a class="reference internal" href="getting_started.html#sec-writing-application-codes"><span class="std std-ref">Writing Application Codes with PETSc</span></a>.</p>
</div>
<div class="section" id="routines-with-different-fortran-interfaces">
<span id="sec-fortran-exceptions"></span><h3>Routines with Different Fortran Interfaces<a class="headerlink" href="#routines-with-different-fortran-interfaces" title="Permalink to this headline">¶</a></h3>
<p>The following Fortran routines differ slightly from their C
counterparts; see the manual pages and previous discussion in this
chapter for details:</p>
<div class="highlight-fortran 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="nb">char</span> <span class="o">*</span><span class="n">filename</span><span class="p">,</span><span class="nb">int </span><span class="n">ierr</span><span class="p">)</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscError.html#PetscError">PetscError</a></span><span class="p">(</span><span class="n">MPI_COMM</span><span class="p">,</span><span class="nb">int </span><span class="n">err</span><span class="p">,</span><span class="nb">char</span> <span class="o">*</span><span class="n">message</span><span class="p">,</span><span class="nb">int </span><span class="n">ierr</span><span class="p">)</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDenseGetArray.html#MatDenseGetArray">MatDenseGetArray</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetIndices.html#ISGetIndices">ISGetIndices</a></span><span class="p">(),</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicateVecs.html#VecDuplicateVecs">VecDuplicateVecs</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroyVecs.html#VecDestroyVecs">VecDestroyVecs</a></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>
</pre></div>
</div>
<p>The following functions are not supported in Fortran:</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFClose.html#PetscFClose">PetscFClose</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFOpen.html#PetscFOpen">PetscFOpen</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFPrintf.html#PetscFPrintf">PetscFPrintf</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">()</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> <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/Profiling/PetscInfo.html#PetscInfo">PetscInfo</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscSetDebugger.html#PetscSetDebugger">PetscSetDebugger</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrays.html#VecGetArrays">VecGetArrays</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrays.html#VecRestoreArrays">VecRestoreArrays</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerASCIIGetPointer.html#PetscViewerASCIIGetPointer">PetscViewerASCIIGetPointer</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerBinaryGetDescriptor.html#PetscViewerBinaryGetDescriptor">PetscViewerBinaryGetDescriptor</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerStringOpen.html#PetscViewerStringOpen">PetscViewerStringOpen</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerStringSPrintf.html#PetscViewerStringSPrintf">PetscViewerStringSPrintf</a></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>
</pre></div>
</div>
<p>PETSc includes some support for direct use of Fortran90 pointers.
Current routines include:</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayF90.html#VecGetArrayF90">VecGetArrayF90</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayF90.html#VecRestoreArrayF90">VecRestoreArrayF90</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayReadF90.html#VecGetArrayReadF90">VecGetArrayReadF90</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayReadF90.html#VecRestoreArrayReadF90">VecRestoreArrayReadF90</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicateVecsF90.html#VecDuplicateVecsF90">VecDuplicateVecsF90</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroyVecsF90.html#VecDestroyVecsF90">VecDestroyVecsF90</a></span><span class="p">()</span>
<span class="n">DMDAVecGetArrayF90</span><span class="p">(),</span> <span class="n">DMDAVecGetArrayReadF90</span><span class="p">(),</span> <span class="n">ISLocalToGlobalMappingGetIndicesF90</span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDenseGetArrayF90.html#MatDenseGetArrayF90">MatDenseGetArrayF90</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDenseRestoreArrayF90.html#MatDenseRestoreArrayF90">MatDenseRestoreArrayF90</a></span><span class="p">()</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetIndicesF90.html#ISGetIndicesF90">ISGetIndicesF90</a></span><span class="p">(),</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISRestoreIndicesF90.html#ISRestoreIndicesF90">ISRestoreIndicesF90</a></span><span class="p">()</span>
</pre></div>
</div>
<p>See the manual pages for details and pointers to example programs.</p>
</div>
</div>
<div class="section" id="sample-fortran-programs">
<span id="sec-fortran-examples"></span><h2>Sample Fortran Programs<a class="headerlink" href="#sample-fortran-programs" title="Permalink to this headline">¶</a></h2>
<p>Sample programs that illustrate the PETSc interface for Fortran are
given below, corresponding to
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/vec/vec/tests/ex19f.F.html">Vec Test ex19f</a>,
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/vec/vec/tutorials/ex4f.F.html">Vec Tutorial ex4f</a>,
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/sys/classes/draw/tests/ex5f.F.html">Draw Test ex5f</a>,
and
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex1f.F90.html">SNES Tutorial ex1f</a>,
respectively. We also refer Fortran programmers to the C examples listed
throughout the manual, since PETSc usage within the two languages
differs only slightly.</p>
<div class="admonition-listing-src-vec-vec-tests-ex19f-f admonition" id="vec-test-ex19f">
<p class="admonition-title">Listing: <code class="docutils literal notranslate"><span class="pre">src/vec/vec/tests/ex19f.F</span></code></p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="c">!</span>
<span class="c">!</span>
<span class="k">program </span><span class="n">main</span>
<span class="cp">#include <petsc/finclude/petscvec.h></span>
<span class="k">use </span><span class="n">petscvec</span>
<span class="k">implicit none</span>
<span class="c">!</span>
<span class="c">! This example demonstrates basic use of the PETSc Fortran interface</span>
<span class="c">! to vectors.</span>
<span class="c">!</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="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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">flg</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">one</span><span class="p">,</span><span class="n">two</span><span class="p">,</span><span class="n">three</span><span class="p">,</span><span class="n">dot</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">norm</span><span class="p">,</span><span class="n">rdot</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="n">x</span><span class="p">,</span><span class="n">y</span><span class="p">,</span><span class="n">w</span>
<span class="n">PetscOptions</span> <span class="k">options</span>
<span class="k"> </span><span class="n">n</span> <span class="o">=</span> <span class="mi">20</span>
<span class="n">one</span> <span class="o">=</span> <span class="mf">1.0</span>
<span class="n">two</span> <span class="o">=</span> <span class="mf">2.0</span>
<span class="n">three</span> <span class="o">=</span> <span class="mf">3.0</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="n">PETSC_NULL_CHARACTER</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">ierr</span> <span class="p">.</span><span class="n">ne</span><span class="p">.</span> <span class="mi">0</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> print</span><span class="o">*</span><span class="p">,</span><span class="s1">'Unable to initialize PETSc'</span>
<span class="k">stop</span>
<span class="k"> endif</span>
<span class="k"> call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsCreate.html#PetscOptionsCreate">PetscOptionsCreate</a></span><span class="p">(</span><span class="k">options</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="k">options</span><span class="p">,</span><span class="n">PETSC_NULL_CHARACTER</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="s1">'-n'</span><span class="p">,</span><span class="n">n</span><span class="p">,</span><span class="n">flg</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsDestroy.html#PetscOptionsDestroy">PetscOptionsDestroy</a></span><span class="p">(</span><span class="k">options</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Create a vector, then duplicate it</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreate.html#VecCreate">VecCreate</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="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetSizes.html#VecSetSizes">VecSetSizes</a></span><span class="p">(</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/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="n">n</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetFromOptions.html#VecSetFromOptions">VecSetFromOptions</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></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">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">w</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">one</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n">two</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDot.html#VecDot">VecDot</a></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">dot</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="n">rdot</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscRealPart.html#PetscRealPart">PetscRealPart</a></span><span class="p">(</span><span class="n">dot</span><span class="p">)</span>
<span class="k">write</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span><span class="mi">100</span><span class="p">)</span> <span class="n">rdot</span>
<span class="mi">100</span> <span class="k">format</span><span class="p">(</span><span class="s1">'Result of inner product '</span><span class="p">,</span><span class="n">f10</span><span class="p">.</span><span class="mi">4</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScale.html#VecScale">VecScale</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">two</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</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/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="n">norm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">write</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span><span class="mi">110</span><span class="p">)</span> <span class="n">norm</span>
<span class="mi">110</span> <span class="k">format</span><span class="p">(</span><span class="s1">'Result of scaling '</span><span class="p">,</span><span class="n">f10</span><span class="p">.</span><span class="mi">4</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCopy.html#VecCopy">VecCopy</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">w</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</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/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="n">norm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">write</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span><span class="mi">120</span><span class="p">)</span> <span class="n">norm</span>
<span class="mi">120</span> <span class="k">format</span><span class="p">(</span><span class="s1">'Result of copy '</span><span class="p">,</span><span class="n">f10</span><span class="p">.</span><span class="mi">4</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAXPY.html#VecAXPY">VecAXPY</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n">three</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="n">norm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">write</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span><span class="mi">130</span><span class="p">)</span> <span class="n">norm</span>
<span class="mi">130</span> <span class="k">format</span><span class="p">(</span><span class="s1">'Result of axpy '</span><span class="p">,</span><span class="n">f10</span><span class="p">.</span><span class="mi">4</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n">w</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="n">ierr</span><span class="p">)</span>
<span class="k">end</span>
<span class="c">!/*TEST</span>
<span class="c">!</span>
<span class="c">! test:</span>
<span class="c">!</span>
<span class="c">!TEST*/</span>
</pre></div>
</div>
</div>
<div class="admonition-listing-src-vec-vec-tutorials-ex4f-f admonition" id="vec-ex4f">
<span id="listing-vec-ex4f"></span><p class="admonition-title">Listing: <code class="docutils literal notranslate"><span class="pre">src/vec/vec/tutorials/ex4f.F</span></code></p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="c">!</span>
<span class="c">!</span>
<span class="c">! Description: Illustrates the use of <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</a>() to set</span>
<span class="c">! multiple values at once; demonstrates <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>().</span>
<span class="c">!</span>
<span class="c">!/*T</span>
<span class="c">! Concepts: vectors^assembling;</span>
<span class="c">! Concepts: vectors^arrays of vectors;</span>
<span class="c">! Processors: 1</span>
<span class="c">!T*/</span>
<span class="c">! -----------------------------------------------------------------------</span>
<span class="k">program </span><span class="n">main</span>
<span class="cp">#include <petsc/finclude/petscvec.h></span>
<span class="k">use </span><span class="n">petscvec</span>
<span class="k">implicit none</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Macro definitions</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">!</span>
<span class="c">! Macros to make clearer the process of setting values in vectors and</span>
<span class="c">! getting values from vectors.</span>
<span class="c">!</span>
<span class="c">! - The element xx_a(ib) is element ib+1 in the vector x</span>
<span class="c">! - Here we add 1 to the base array index to facilitate the use of</span>
<span class="c">! conventional Fortran 1-based array indexing.</span>
<span class="c">!</span>
<span class="cp">#define xx_a(ib) xx_v(xx_i + (ib))</span>
<span class="cp">#define yy_a(ib) yy_v(yy_i + (ib))</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Beginning of program</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">xwork</span><span class="p">(</span><span class="mi">6</span><span class="p">)</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">xx_v</span><span class="p">(</span><span class="mi">1</span><span class="p">),</span><span class="n">yy_v</span><span class="p">(</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/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">i</span><span class="p">,</span><span class="n">n</span><span class="p">,</span><span class="nb">loc</span><span class="p">(</span><span class="mi">6</span><span class="p">),</span><span class="n">isix</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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOffset.html#PetscOffset">PetscOffset</a></span> <span class="n">xx_i</span><span class="p">,</span><span class="n">yy_i</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="n">x</span><span class="p">,</span><span class="n">y</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="n">PETSC_NULL_CHARACTER</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">ierr</span> <span class="p">.</span><span class="n">ne</span><span class="p">.</span> <span class="mi">0</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> print</span><span class="o">*</span><span class="p">,</span><span class="s1">'<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a> failed'</span>
<span class="k">stop</span>
<span class="k"> endif</span>
<span class="k"> </span><span class="n">n</span> <span class="o">=</span> <span class="mi">6</span>
<span class="n">isix</span> <span class="o">=</span> <span class="mi">6</span>
<span class="c">! Create initial vector and duplicate it</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateSeq.html#VecCreateSeq">VecCreateSeq</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_SELF.html#PETSC_COMM_SELF">PETSC_COMM_SELF</a></span><span class="p">,</span><span class="n">n</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></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">ierr</span><span class="p">)</span>
<span class="c">! Fill work arrays with vector entries and locations. Note that</span>
<span class="c">! the vector indices are 0-based in PETSc (for both Fortran and</span>
<span class="c">! C vectors)</span>
<span class="k">do </span><span class="mi">10</span> <span class="n">i</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">n</span>
<span class="nb">loc</span><span class="p">(</span><span class="n">i</span><span class="p">)</span> <span class="o">=</span> <span class="n">i</span><span class="o">-</span><span class="mi">1</span>
<span class="n">xwork</span><span class="p">(</span><span class="n">i</span><span class="p">)</span> <span class="o">=</span> <span class="mi">1</span><span class="mf">0.0</span><span class="o">*</span><span class="kt">real</span><span class="p">(</span><span class="n">i</span><span class="p">)</span>
<span class="mi">10</span> <span class="k">continue</span>
<span class="c">! Set vector values. Note that we set multiple entries at once.</span>
<span class="c">! Of course, usually one would create a work array that is the</span>
<span class="c">! natural size for a particular problem (not one that is as long</span>
<span class="c">! as the full vector).</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">isix</span><span class="p">,</span><span class="nb">loc</span><span class="p">,</span><span class="n">xwork</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Assemble vector</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAssemblyBegin.html#VecAssemblyBegin">VecAssemblyBegin</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAssemblyEnd.html#VecAssemblyEnd">VecAssemblyEnd</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! View vector</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSetName.html#PetscObjectSetName">PetscObjectSetName</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span> <span class="s1">'initial vector:'</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a></span><span class="p">(</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/Viewer/PETSC_VIEWER_STDOUT_SELF.html#PETSC_VIEWER_STDOUT_SELF">PETSC_VIEWER_STDOUT_SELF</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCopy.html#VecCopy">VecCopy</a></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">ierr</span><span class="p">)</span>
<span class="c">! Get a pointer to vector data.</span>
<span class="c">! - For default PETSc vectors, <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>() returns a pointer to</span>
<span class="c">! the data array. Otherwise, the routine is implementation dependent.</span>
<span class="c">! - You MUST call <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a>() when you no longer need access to</span>
<span class="c">! the array.</span>
<span class="c">! - Note that the Fortran interface to <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>() differs from the</span>
<span class="c">! C version. See the users manual for details.</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">xx_v</span><span class="p">,</span><span class="n">xx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n">yy_v</span><span class="p">,</span><span class="n">yy_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Modify vector data</span>
<span class="k">do </span><span class="mi">30</span> <span class="n">i</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">n</span>
<span class="n">xx_a</span><span class="p">(</span><span class="n">i</span><span class="p">)</span> <span class="o">=</span> <span class="mi">10</span><span class="mf">0.0</span><span class="o">*</span><span class="kt">real</span><span class="p">(</span><span class="n">i</span><span class="p">)</span>
<span class="n">yy_a</span><span class="p">(</span><span class="n">i</span><span class="p">)</span> <span class="o">=</span> <span class="mi">100</span><span class="mf">0.0</span><span class="o">*</span><span class="kt">real</span><span class="p">(</span><span class="n">i</span><span class="p">)</span>
<span class="mi">30</span> <span class="k">continue</span>
<span class="c">! Restore vectors</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">xx_v</span><span class="p">,</span><span class="n">xx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n">yy_v</span><span class="p">,</span><span class="n">yy_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! View vectors</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSetName.html#PetscObjectSetName">PetscObjectSetName</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span> <span class="s1">'new vector 1:'</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a></span><span class="p">(</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/Viewer/PETSC_VIEWER_STDOUT_SELF.html#PETSC_VIEWER_STDOUT_SELF">PETSC_VIEWER_STDOUT_SELF</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSetName.html#PetscObjectSetName">PetscObjectSetName</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span> <span class="s1">'new vector 2:'</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n"><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><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Free work space. All PETSc objects should be destroyed when they</span>
<span class="c">! are no longer needed.</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="n">ierr</span><span class="p">)</span>
<span class="k">end</span>
<span class="c">!/*TEST</span>
<span class="c">!</span>
<span class="c">! test:</span>
<span class="c">!</span>
<span class="c">!TEST*/</span>
</pre></div>
</div>
</div>
<div class="admonition-listing-src-sys-classes-draw-tests-ex5f-f admonition" id="draw-test-ex5f">
<p class="admonition-title">Listing: <code class="docutils literal notranslate"><span class="pre">src/sys/classes/draw/tests/ex5f.F</span></code></p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="c">!</span>
<span class="c">!</span>
<span class="k">program </span><span class="n">main</span>
<span class="cp">#include <petsc/finclude/petscsys.h></span>
<span class="cp">#include <petsc/finclude/petscdraw.h></span>
<span class="k">use </span><span class="n">petscsys</span>
<span class="k">implicit none</span>
<span class="c">!</span>
<span class="c">! This example demonstrates basic use of the Fortran interface for</span>
<span class="c">! <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Draw/PetscDraw.html#PetscDraw">PetscDraw</a> routines.</span>
<span class="c">!</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="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="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="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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">flg</span>
<span class="kt">integer </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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">xd</span><span class="p">,</span><span class="n">yd</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">ten</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">i</span><span class="p">,</span><span class="n">n</span><span class="p">,</span><span class="n">w</span><span class="p">,</span><span class="n">h</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">one</span>
<span class="n">n</span> <span class="o">=</span> <span class="mi">15</span>
<span class="n">x</span> <span class="o">=</span> <span class="mi">0</span>
<span class="n">y</span> <span class="o">=</span> <span class="mi">0</span>
<span class="n">w</span> <span class="o">=</span> <span class="mi">400</span>
<span class="n">h</span> <span class="o">=</span> <span class="mi">300</span>
<span class="n">ten</span> <span class="o">=</span> <span class="mi">1</span><span class="mf">0.0</span>
<span class="n">one</span> <span class="o">=</span> <span class="mi">1</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="n">PETSC_NULL_CHARACTER</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">ierr</span> <span class="p">.</span><span class="n">ne</span><span class="p">.</span> <span class="mi">0</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> print</span><span class="o">*</span><span class="p">,</span><span class="s1">'Unable to initialize PETSc'</span>
<span class="k">stop</span>
<span class="k"> endif</span>
<span class="c">! GetInt requires a <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a> so have to do this ugly setting</span>
<span class="k">call </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">PETSC_NULL_OPTIONS</span><span class="p">,</span><span class="n">PETSC_NULL_CHARACTER</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="s1">'-width'</span><span class="p">,</span><span class="n">w</span><span class="p">,</span> <span class="n">flg</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="n">width</span> <span class="o">=</span> <span class="nb">int</span><span class="p">(</span><span class="n">w</span><span class="p">)</span>
<span class="k">call </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">PETSC_NULL_OPTIONS</span><span class="p">,</span><span class="n">PETSC_NULL_CHARACTER</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="s1">'-height'</span><span class="p">,</span><span class="n">h</span><span class="p">,</span><span class="n">flg</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="n">height</span> <span class="o">=</span> <span class="nb">int</span><span class="p">(</span><span class="n">h</span><span class="p">)</span>
<span class="k">call </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">PETSC_NULL_OPTIONS</span><span class="p">,</span><span class="n">PETSC_NULL_CHARACTER</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="s1">'-n'</span><span class="p">,</span><span class="n">n</span><span class="p">,</span><span class="n">flg</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="n">PETSC_NULL_CHARACTER</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="n">PETSC_NULL_CHARACTER</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="n">draw</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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">ierr</span><span class="p">)</span>
<span class="k">call </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="n">one</span><span class="p">,</span><span class="n">lg</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="n">axis</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="p">&</span>
<span class="p">&</span> <span class="n">PETSC_DRAW_BLUE</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="s1">'toplabel'</span><span class="p">,</span><span class="s1">'xlabel'</span><span class="p">,</span><span class="s1">'ylabel'</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="n">ierr</span><span class="p">)</span>
<span class="k">do </span><span class="mi">10</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">n</span><span class="o">-</span><span class="mi">1</span>
<span class="n">xd</span> <span class="o">=</span> <span class="kt">real</span><span class="p">(</span><span class="n">i</span><span class="p">)</span> <span class="o">-</span> <span class="mf">5.0</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="k">call </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="n">xd</span><span class="p">,</span><span class="n">yd</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="mi">10</span> <span class="k">continue</span>
<span class="k"> call </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">ierr</span><span class="p">)</span>
<span class="k">call </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">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscSleep.html#PetscSleep">PetscSleep</a></span><span class="p">(</span><span class="n">ten</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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">lg</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="n">draw</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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="n">ierr</span><span class="p">)</span>
<span class="k">end</span>
<span class="c">!/*TEST</span>
<span class="c">!</span>
<span class="c">! build:</span>
<span class="c">! requires: x</span>
<span class="c">!</span>
<span class="c">! test:</span>
<span class="c">! output_file: output/ex1_1.out</span>
<span class="c">!</span>
<span class="c">!TEST*/</span>
</pre></div>
</div>
</div>
<div class="admonition-listing-src-snes-tutorials-ex1f-f90 admonition" id="snes-ex1f">
<p class="admonition-title">Listing: <code class="docutils literal notranslate"><span class="pre">src/snes/tutorials/ex1f.F90</span></code></p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="c">!</span>
<span class="c">!</span>
<span class="c">! Description: Uses the Newton method to solve a two-variable system.</span>
<span class="c">!</span>
<span class="c">!!/*T</span>
<span class="c">! Concepts: <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a>^basic uniprocessor example</span>
<span class="c">! Processors: 1</span>
<span class="c">!T*/</span>
<span class="k">program </span><span class="n">main</span>
<span class="cp">#include <petsc/finclude/petsc.h></span>
<span class="k">use </span><span class="n">petsc</span>
<span class="k">implicit none</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Variable declarations</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">!</span>
<span class="c">! Variables:</span>
<span class="c">! snes - nonlinear solver</span>
<span class="c">! ksp - linear solver</span>
<span class="c">! pc - preconditioner context</span>
<span class="c">! ksp - Krylov subspace method context</span>
<span class="c">! x, r - solution, residual vectors</span>
<span class="c">! J - Jacobian matrix</span>
<span class="c">! its - iterations for convergence</span>
<span class="c">!</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span> <span class="n">pc</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span> <span class="n">ksp</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="n">x</span><span class="p">,</span><span class="n">r</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="n">J</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">linesearch</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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">its</span><span class="p">,</span><span class="n">i2</span><span class="p">,</span><span class="n">i20</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscMPIInt.html#PetscMPIInt">PetscMPIInt</a></span> <span class="n">size</span><span class="p">,</span><span class="n">rank</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">pfive</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">tol</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">setls</span>
<span class="cp">#if defined(PETSC_USE_LOG)</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="cp">#endif</span>
<span class="kt">double precision </span><span class="n">threshold</span><span class="p">,</span><span class="n">oldthreshold</span>
<span class="c">! Note: Any user-defined Fortran routines (such as FormJacobian)</span>
<span class="c">! MUST be declared as external.</span>
<span class="k">external </span><span class="n">FormFunction</span><span class="p">,</span> <span class="n">FormJacobian</span><span class="p">,</span> <span class="n">MyLineSearch</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Macro definitions</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">!</span>
<span class="c">! Macros to make clearer the process of setting values in vectors and</span>
<span class="c">! getting values from vectors. These vectors are used in the routines</span>
<span class="c">! FormFunction() and FormJacobian().</span>
<span class="c">! - The element lx_a(ib) is element ib in the vector x</span>
<span class="c">!</span>
<span class="cp">#define lx_a(ib) lx_v(lx_i + (ib))</span>
<span class="cp">#define lf_a(ib) lf_v(lf_i + (ib))</span>
<span class="c">!</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Beginning of program</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</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="n">PETSC_NULL_CHARACTER</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">ierr</span> <span class="p">.</span><span class="n">ne</span><span class="p">.</span> <span class="mi">0</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> print</span><span class="o">*</span><span class="p">,</span><span class="s1">'Unable to initialize PETSc'</span>
<span class="k">stop</span>
<span class="k"> endif</span>
<span class="k"> call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Profiling/PetscLogNestedBegin.html#PetscLogNestedBegin">PetscLogNestedBegin</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span><span class="n">CHKERRA</span><span class="p">(</span><span class="n">ierr</span><span class="p">)</span>
<span class="n">threshold</span> <span class="o">=</span> <span class="mf">1.0</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Profiling/PetscLogSetThreshold.html#PetscLogSetThreshold">PetscLogSetThreshold</a></span><span class="p">(</span><span class="n">threshold</span><span class="p">,</span><span class="n">oldthreshold</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! dummy test of logging a reduction</span>
<span class="cp">#if defined(PETSC_USE_LOG)</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n">PetscAReduce</span><span class="p">()</span>
<span class="cp">#endif</span>
<span class="k">call </span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Comm_size.html#MPI_Comm_size">MPI_Comm_size</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="n">size</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Comm_rank.html#MPI_Comm_rank">MPI_Comm_rank</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="n">rank</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">size</span> <span class="p">.</span><span class="n">ne</span><span class="p">.</span> <span class="mi">1</span><span class="p">)</span> <span class="k">then</span><span class="p">;</span> <span class="n">SETERRA</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_SELF.html#PETSC_COMM_SELF">PETSC_COMM_SELF</a></span><span class="p">,</span><span class="n">PETSC_ERR_WRONG_MPI_SIZE</span><span class="p">,</span><span class="s1">'Uniprocessor example'</span><span class="p">);</span> <span class="k">endif</span>
<span class="k"> </span><span class="n">i2</span> <span class="o">=</span> <span class="mi">2</span>
<span class="n">i20</span> <span class="o">=</span> <span class="mi">20</span>
<span class="c">! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Create nonlinear solver context</span>
<span class="c">! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESCreate.html#SNESCreate">SNESCreate</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="n">snes</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Create matrix and vector data structures; set corresponding routines</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Create vectors for solution and nonlinear function</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateSeq.html#VecCreateSeq">VecCreateSeq</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_SELF.html#PETSC_COMM_SELF">PETSC_COMM_SELF</a></span><span class="p">,</span><span class="n">i2</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">r</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Create Jacobian matrix data structure</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_SELF.html#PETSC_COMM_SELF">PETSC_COMM_SELF</a></span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetSizes.html#MatSetSizes">MatSetSizes</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="n">i2</span><span class="p">,</span><span class="n">i2</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetFromOptions.html#MatSetFromOptions">MatSetFromOptions</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetUp.html#MatSetUp">MatSetUp</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Set function evaluation routine and vector</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html#SNESSetFunction">SNESSetFunction</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">r</span><span class="p">,</span><span class="n">FormFunction</span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Set Jacobian matrix data structure and Jacobian evaluation routine</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">FormJacobian</span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Customize nonlinear solver; set runtime options</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Set linear solver defaults for this problem. By extracting the</span>
<span class="c">! <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>, <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>, and <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> contexts from the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context, we can then</span>
<span class="c">! directly call any <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>, <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>, and <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> routines to set various options.</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetKSP.html#SNESGetKSP">SNESGetKSP</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">ksp</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGetPC.html#KSPGetPC">KSPGetPC</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="n">pc</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSetType.html#PCSetType">PCSetType</a></span><span class="p">(</span><span class="n">pc</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCNONE.html#PCNONE">PCNONE</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="n">tol</span> <span class="o">=</span> <span class="mf">1.e-4</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetTolerances.html#KSPSetTolerances">KSPSetTolerances</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="n">tol</span><span class="p">,</span><span class="n">PETSC_DEFAULT_REAL</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="n">PETSC_DEFAULT_REAL</span><span class="p">,</span><span class="n">i20</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Set <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> runtime options, e.g.,</span>
<span class="c">! -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc></span>
<span class="c">! These options will override those specified above as long as</span>
<span class="c">! <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFromOptions.html#SNESSetFromOptions">SNESSetFromOptions</a>() is called _after_ any other customization</span>
<span class="c">! routines.</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFromOptions.html#SNESSetFromOptions">SNESSetFromOptions</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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">PETSC_NULL_OPTIONS</span><span class="p">,</span><span class="n">PETSC_NULL_CHARACTER</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="s1">'-setls'</span><span class="p">,</span><span class="n">setls</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">setls</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetLineSearch.html#SNESGetLineSearch">SNESGetLineSearch</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span> <span class="n">linesearch</span><span class="p">,</span> <span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetType.html#SNESLineSearchSetType">SNESLineSearchSetType</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="s1">'shell'</span><span class="p">,</span> <span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchShellSetUserFunc.html#SNESLineSearchShellSetUserFunc">SNESLineSearchShellSetUserFunc</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="n">MyLineSearch</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="mi">0</span><span class="p">,</span> <span class="n">ierr</span><span class="p">)</span>
<span class="k">endif</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Evaluate initial guess; then solve nonlinear system</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Note: The user should initialize the vector, x, with the initial guess</span>
<span class="c">! for the nonlinear solver prior to calling <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSolve.html#SNESSolve">SNESSolve</a>(). In particular,</span>
<span class="c">! to employ an initial guess of zero, the user should explicitly set</span>
<span class="c">! this vector to zero by calling <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a>().</span>
<span class="n">pfive</span> <span class="o">=</span> <span class="mf">0.5</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">pfive</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSolve.html#SNESSolve">SNESSolve</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">PETSC_NULL_VEC</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! View solver converged reason; we could instead use the option -snes_converged_reason</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESConvergedReasonView.html#SNESConvergedReasonView">SNESConvergedReasonView</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n"><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><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetIterationNumber.html#SNESGetIterationNumber">SNESGetIterationNumber</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">its</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">rank</span> <span class="p">.</span><span class="n">eq</span><span class="p">.</span> <span class="mi">0</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> write</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span><span class="mi">100</span><span class="p">)</span> <span class="n">its</span>
<span class="k">endif</span>
<span class="k"> </span><span class="mi">100</span> <span class="k">format</span><span class="p">(</span><span class="s1">'Number of <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> iterations = '</span><span class="p">,</span><span class="n">i5</span><span class="p">)</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="c">! Free work space. All PETSc objects should be destroyed when they</span>
<span class="c">! are no longer needed.</span>
<span class="c">! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n">r</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDestroy.html#MatDestroy">MatDestroy</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESDestroy.html#SNESDestroy">SNESDestroy</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="cp">#if defined(PETSC_USE_LOG)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerASCIIOpen.html#PetscViewerASCIIOpen">PetscViewerASCIIOpen</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="s1">'filename.xml'</span><span class="p">,</span><span class="n">viewer</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </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">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">PETSC_VIEWER_ASCII_XML</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Profiling/PetscLogView.html#PetscLogView">PetscLogView</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerDestroy.html#PetscViewerDestroy">PetscViewerDestroy</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="cp">#endif</span>
<span class="k">call </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="n">ierr</span><span class="p">)</span>
<span class="k">end</span>
<span class="c">!</span>
<span class="c">! ------------------------------------------------------------------------</span>
<span class="c">!</span>
<span class="c">! FormFunction - Evaluates nonlinear function, F(x).</span>
<span class="c">!</span>
<span class="c">! Input Parameters:</span>
<span class="c">! snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="c">! x - input vector</span>
<span class="c">! dummy - optional user-defined context (not used here)</span>
<span class="c">!</span>
<span class="c">! Output Parameter:</span>
<span class="c">! f - function vector</span>
<span class="c">!</span>
<span class="k">subroutine </span><span class="n">FormFunction</span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">f</span><span class="p">,</span><span class="n">dummy</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">use </span><span class="n">petscsnes</span>
<span class="k">implicit none</span>
<span class="k"> </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</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="n">x</span><span class="p">,</span><span class="n">f</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="kt">integer </span><span class="n">dummy</span><span class="p">(</span><span class="o">*</span><span class="p">)</span>
<span class="c">! Declarations for use with local arrays</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">lx_v</span><span class="p">(</span><span class="mi">2</span><span class="p">),</span><span class="n">lf_v</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOffset.html#PetscOffset">PetscOffset</a></span> <span class="n">lx_i</span><span class="p">,</span><span class="n">lf_i</span>
<span class="c">! Get pointers to vector data.</span>
<span class="c">! - For default PETSc vectors, <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>() returns a pointer to</span>
<span class="c">! the data array. Otherwise, the routine is implementation dependent.</span>
<span class="c">! - You MUST call <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a>() when you no longer need access to</span>
<span class="c">! the array.</span>
<span class="c">! - Note that the Fortran interface to <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>() differs from the</span>
<span class="c">! C version. See the Fortran chapter of the users manual for details.</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">lx_v</span><span class="p">,</span><span class="n">lx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(</span><span class="n">f</span><span class="p">,</span><span class="n">lf_v</span><span class="p">,</span><span class="n">lf_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Compute function</span>
<span class="n">lf_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span> <span class="o">=</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="n">lx_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span> <span class="p">&</span>
<span class="p">&</span> <span class="o">+</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="n">lx_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span> <span class="o">-</span> <span class="mf">3.0</span>
<span class="n">lf_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span> <span class="o">=</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span><span class="o">*</span><span class="n">lx_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span> <span class="p">&</span>
<span class="p">&</span> <span class="o">+</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span><span class="o">*</span><span class="n">lx_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span> <span class="o">-</span> <span class="mf">6.0</span>
<span class="c">! Restore vectors</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayRead.html#VecRestoreArrayRead">VecRestoreArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">lx_v</span><span class="p">,</span><span class="n">lx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n">f</span><span class="p">,</span><span class="n">lf_v</span><span class="p">,</span><span class="n">lf_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">return</span>
<span class="k"> end</span>
<span class="c">! ---------------------------------------------------------------------</span>
<span class="c">!</span>
<span class="c">! FormJacobian - Evaluates Jacobian matrix.</span>
<span class="c">!</span>
<span class="c">! Input Parameters:</span>
<span class="c">! snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="c">! x - input vector</span>
<span class="c">! dummy - optional user-defined context (not used here)</span>
<span class="c">!</span>
<span class="c">! Output Parameters:</span>
<span class="c">! A - Jacobian matrix</span>
<span class="c">! B - optionally different preconditioning matrix</span>
<span class="c">! flag - flag indicating matrix structure</span>
<span class="c">!</span>
<span class="k">subroutine </span><span class="n">FormJacobian</span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">X</span><span class="p">,</span><span class="n">jac</span><span class="p">,</span><span class="n">B</span><span class="p">,</span><span class="n">dummy</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">use </span><span class="n">petscsnes</span>
<span class="k">implicit none</span>
<span class="k"> </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</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="n">X</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="n">jac</span><span class="p">,</span><span class="n">B</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">A</span><span class="p">(</span><span class="mi">4</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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">idx</span><span class="p">(</span><span class="mi">2</span><span class="p">),</span><span class="n">i2</span>
<span class="kt">integer </span><span class="n">dummy</span><span class="p">(</span><span class="o">*</span><span class="p">)</span>
<span class="c">! Declarations for use with local arrays</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">lx_v</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOffset.html#PetscOffset">PetscOffset</a></span> <span class="n">lx_i</span>
<span class="c">! Get pointer to vector data</span>
<span class="n">i2</span> <span class="o">=</span> <span class="mi">2</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">lx_v</span><span class="p">,</span><span class="n">lx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Compute Jacobian entries and insert into matrix.</span>
<span class="c">! - Since this is such a small problem, we set all entries for</span>
<span class="c">! the matrix at once.</span>
<span class="c">! - Note that <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a>() uses 0-based row and column numbers</span>
<span class="c">! in Fortran as well as in C (as set here in the array idx).</span>
<span class="n">idx</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span> <span class="o">=</span> <span class="mi">0</span>
<span class="n">idx</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span> <span class="o">=</span> <span class="mi">1</span>
<span class="n">A</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span> <span class="o">=</span> <span class="mf">2.0</span><span class="o">*</span><span class="n">lx_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span> <span class="o">+</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span>
<span class="n">A</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span> <span class="o">=</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span>
<span class="n">A</span><span class="p">(</span><span class="mi">3</span><span class="p">)</span> <span class="o">=</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span>
<span class="n">A</span><span class="p">(</span><span class="mi">4</span><span class="p">)</span> <span class="o">=</span> <span class="n">lx_a</span><span class="p">(</span><span class="mi">1</span><span class="p">)</span> <span class="o">+</span> <span class="mf">2.0</span><span class="o">*</span><span class="n">lx_a</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="n">i2</span><span class="p">,</span><span class="n">idx</span><span class="p">,</span><span class="n">i2</span><span class="p">,</span><span class="n">idx</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Restore vector</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayRead.html#VecRestoreArrayRead">VecRestoreArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">lx_v</span><span class="p">,</span><span class="n">lx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="c">! Assemble matrix</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</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">B</span> <span class="p">.</span><span class="n">ne</span><span class="p">.</span> <span class="n">jac</span><span class="p">)</span> <span class="k">then</span>
<span class="k"> call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">endif</span>
<span class="k"> return</span>
<span class="k"> end</span>
<span class="k"> subroutine </span><span class="n">MyLineSearch</span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="n">lctx</span><span class="p">,</span> <span class="n">ierr</span><span class="p">)</span>
<span class="k">use </span><span class="n">petscsnes</span>
<span class="k">implicit none</span>
<span class="k"> </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">linesearch</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span>
<span class="kt">integer </span><span class="n">lctx</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="n">x</span><span class="p">,</span> <span class="n">f</span><span class="p">,</span><span class="n">g</span><span class="p">,</span> <span class="n">y</span><span class="p">,</span> <span class="n">w</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">ynorm</span><span class="p">,</span><span class="n">gnorm</span><span class="p">,</span><span class="n">xnorm</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">flag</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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">mone</span>
<span class="n">mone</span> <span class="o">=</span> <span class="o">-</span><span class="mf">1.0</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchGetSNES.html#SNESLineSearchGetSNES">SNESLineSearchGetSNES</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="n">snes</span><span class="p">,</span> <span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchGetVecs.html#SNESLineSearchGetVecs">SNESLineSearchGetVecs</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="n">x</span><span class="p">,</span> <span class="n">f</span><span class="p">,</span> <span class="n">y</span><span class="p">,</span> <span class="n">w</span><span class="p">,</span> <span class="n">g</span><span class="p">,</span> <span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="n">ynorm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAXPY.html#VecAXPY">VecAXPY</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">mone</span><span class="p">,</span><span class="n">y</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESComputeFunction.html#SNESComputeFunction">SNESComputeFunction</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">f</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</span><span class="n">f</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="n">gnorm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</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/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="n">xnorm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="n">ynorm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetNorms.html#SNESLineSearchSetNorms">SNESLineSearchSetNorms</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="n">xnorm</span><span class="p">,</span> <span class="n">gnorm</span><span class="p">,</span> <span class="n">ynorm</span><span class="p">,</span> <span class="p">&</span>
<span class="p">&</span> <span class="n">ierr</span><span class="p">)</span>
<span class="n">flag</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</a></span>
<span class="k">return</span>
<span class="k"> end</span>
<span class="c">!/*TEST</span>
<span class="c">!</span>
<span class="c">! test:</span>
<span class="c">! args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short</span>
<span class="c">! requires: !single</span>
<span class="c">!</span>
<span class="c">!TEST*/</span>
</pre></div>
</div>
</div>
<div class="section" id="array-arguments">
<span id="sec-fortranarrays"></span><h3>Array Arguments<a class="headerlink" href="#array-arguments" title="Permalink to this headline">¶</a></h3>
<p>This material is no longer relevent since one should use
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayF90.html#VecGetArrayF90">VecGetArrayF90</a>()</span></code> and the other routines that utilize Fortran
pointers, instead of the code below, but it is included for historical
reasons and because many of the Fortran examples still utilize the old
approach.</p>
<p>Since Fortran 77 does not allow arrays to be returned in routine
arguments, all PETSc routines that return arrays, such as
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDenseGetArray.html#MatDenseGetArray">MatDenseGetArray</a>()</span></code>, and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetIndices.html#ISGetIndices">ISGetIndices</a>()</span></code>, are
defined slightly differently in Fortran than in C. Instead of returning
the array itself, these routines accept as input a user-specified array
of dimension one and return an integer index to the actual array used
for data storage within PETSc. The Fortran interface for several
routines is as follows:</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">xx_v</span><span class="p">(</span><span class="mi">1</span><span class="p">),</span> <span class="n">aa_v</span><span class="p">(</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/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</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">ss_v</span><span class="p">(</span><span class="mi">1</span><span class="p">),</span> <span class="n">dd_v</span><span class="p">(</span><span class="mi">1</span><span class="p">),</span> <span class="n">nloc</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOffset.html#PetscOffset">PetscOffset</a></span> <span class="n">ss_i</span><span class="p">,</span> <span class="n">xx_i</span><span class="p">,</span> <span class="n">aa_i</span><span class="p">,</span> <span class="n">dd_i</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="n">x</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="n">A</span>
<span class="k"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a> </span><span class="n">s</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DM.html#DM">DM</a></span> <span class="n">d</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">xx_v</span><span class="p">,</span><span class="n">xx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDenseGetArray.html#MatDenseGetArray">MatDenseGetArray</a></span><span class="p">(</span><span class="n">A</span><span class="p">,</span><span class="n">aa_v</span><span class="p">,</span><span class="n">aa_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetIndices.html#ISGetIndices">ISGetIndices</a></span><span class="p">(</span><span class="n">s</span><span class="p">,</span><span class="n">ss_v</span><span class="p">,</span><span class="n">ss_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
</pre></div>
</div>
<p>To access array elements directly, both the user-specified array and the
integer index <em>must</em> then be used together. For example, the following
Fortran program fragment illustrates directly setting the values of a
vector array instead of using <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</a>()</span></code>. Note the (optional)
use of the preprocessor <code class="docutils literal notranslate"><span class="pre">#define</span></code> statement to enable array
manipulations in the conventional Fortran manner.</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="cp">#define xx_a(ib) xx_v(xx_i + (ib))</span>
<span class="kt">double precision </span><span class="n">xx_v</span><span class="p">(</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/PetscOffset.html#PetscOffset">PetscOffset</a></span> <span class="n">xx_i</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="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">i</span><span class="p">,</span> <span class="n">n</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="n">x</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">xx_v</span><span class="p">,</span><span class="n">xx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetLocalSize.html#VecGetLocalSize">VecGetLocalSize</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">n</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
<span class="k">do </span><span class="mi">10</span><span class="p">,</span> <span class="n">i</span><span class="o">=</span><span class="mi">1</span><span class="p">,</span><span class="n">n</span>
<span class="n">xx_a</span><span class="p">(</span><span class="n">i</span><span class="p">)</span> <span class="o">=</span> <span class="mi">3</span><span class="o">*</span><span class="n">i</span> <span class="o">+</span> <span class="mi">1</span>
<span class="mi">10</span> <span class="k">continue</span>
<span class="k"> call </span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">xx_v</span><span class="p">,</span><span class="n">xx_i</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
</pre></div>
</div>
<p><a class="reference internal" href="#listing-vec-ex4f"><span class="std std-ref">The Vec ex4f Tutorial listed above</span></a> contains an example of
using <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>()</span></code> within a Fortran routine.</p>
<p>Since in this case the array is accessed directly from Fortran, indexing
begins with 1, not 0 (unless the array is declared as <code class="docutils literal notranslate"><span class="pre">xx_v(0:1)</span></code>).
This is different from the use of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</a>()</span></code> where, indexing
always starts with 0.</p>
<p><em>Note</em>: If using <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDenseGetArray.html#MatDenseGetArray">MatDenseGetArray</a>()</span></code>, or
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetIndices.html#ISGetIndices">ISGetIndices</a>()</span></code>, from Fortran, the user <em>must not</em> compile the
Fortran code with options to check for “array entries out of bounds”
(e.g., on the IBM RS/6000 this is done with the <code class="docutils literal notranslate"><span class="pre">-C</span></code> compiler option,
so never use the <code class="docutils literal notranslate"><span class="pre">-C</span></code> option with this).</p>
</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="matlab.html" title="Using MATLAB with PETSc"
>next</a> |</li>
<li class="right" >
<a href="additional.html" title="Additional Information"
>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>
|