1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191
|
<!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/vec.html" />
<meta charset="utf-8" />
<title>Vectors and Parallel Data — 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="Matrices" href="mat.html" />
<link rel="prev" title="Programming with PETSc" href="programming.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/vec.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="mat.html" title="Matrices"
accesskey="N">next</a> |</li>
<li class="right" >
<a href="programming.html" title="Programming with PETSc"
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="programming.html" accesskey="U">Programming with PETSc</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="#">Vectors and Parallel Data</a><ul>
<li><a class="reference internal" href="#creating-and-assembling-vectors">Creating and Assembling Vectors</a></li>
<li><a class="reference internal" href="#basic-vector-operations">Basic Vector Operations</a></li>
<li><a class="reference internal" href="#indexing-and-ordering">Indexing and Ordering</a><ul>
<li><a class="reference internal" href="#application-orderings">Application Orderings</a></li>
<li><a class="reference internal" href="#local-to-global-mappings">Local to Global Mappings</a></li>
</ul>
</li>
<li><a class="reference internal" href="#structured-grids-using-distributed-arrays">Structured Grids Using Distributed Arrays</a><ul>
<li><a class="reference internal" href="#creating-distributed-arrays">Creating Distributed Arrays</a></li>
<li><a class="reference internal" href="#local-global-vectors-and-scatters">Local/Global Vectors and Scatters</a></li>
<li><a class="reference internal" href="#local-ghosted-work-vectors">Local (Ghosted) Work Vectors</a></li>
<li><a class="reference internal" href="#accessing-the-vector-entries-for-dmda-vectors">Accessing the Vector Entries for DMDA Vectors</a></li>
<li><a class="reference internal" href="#grid-information">Grid Information</a></li>
<li><a class="reference internal" href="#staggered-grids">Staggered Grids</a></li>
</ul>
</li>
<li><a class="reference internal" href="#vectors-related-to-unstructured-grids">Vectors Related to Unstructured Grids</a><ul>
<li><a class="reference internal" href="#index-sets">Index Sets</a></li>
<li><a class="reference internal" href="#scatters-and-gathers">Scatters and Gathers</a></li>
<li><a class="reference internal" href="#scattering-ghost-values">Scattering Ghost Values</a></li>
<li><a class="reference internal" href="#vectors-with-locations-for-ghost-values">Vectors with Locations for Ghost Values</a></li>
</ul>
</li>
</ul>
</li>
</ul>
<h4>Previous topic</h4>
<p class="topless"><a href="programming.html"
title="previous chapter">Programming with PETSc</a></p>
<h4>Next topic</h4>
<p class="topless"><a href="mat.html"
title="next chapter">Matrices</a></p>
<div role="note" aria-label="source link">
<h3>This Page</h3>
<ul class="this-page-menu">
<li><a href="../_sources/manual/vec.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="vectors-and-parallel-data">
<span id="chapter-vectors"></span><h1>Vectors and Parallel Data<a class="headerlink" href="#vectors-and-parallel-data" title="Permalink to this headline">¶</a></h1>
<p>The vector (denoted by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span></code>) is one of the simplest PETSc objects.
Vectors are used to store discrete PDE solutions, right-hand sides for
linear systems, etc. This chapter is organized as follows:</p>
<ul class="simple">
<li><p>(<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span></code>) <a class="reference internal" href="#sec-veccreate"><span class="std std-ref">Creating and Assembling Vectors</span></a> and
<a class="reference internal" href="#sec-vecbasic"><span class="std std-ref">Basic Vector Operations</span></a> - basic usage of vectors</p></li>
<li><p>Section <a class="reference internal" href="#sec-indexingandordering"><span class="std std-ref">Indexing and Ordering</span></a> - management of the
various numberings of degrees of freedom, vertices, cells, etc.</p>
<ul>
<li><p>(<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code>) Mapping between different global numberings</p></li>
<li><p>(<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span></code>) Mapping between local and global
numberings</p></li>
</ul>
</li>
<li><p>(<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DM.html#DM">DM</a></span></code>) <a class="reference internal" href="#sec-struct"><span class="std std-ref">Structured Grids Using Distributed Arrays</span></a> - management of grids</p></li>
<li><p>(<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>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span></code>) <a class="reference internal" href="#sec-unstruct"><span class="std std-ref">Vectors Related to Unstructured Grids</span></a> - management
of vectors related to unstructured grids</p></li>
</ul>
<div class="section" id="creating-and-assembling-vectors">
<span id="sec-veccreate"></span><h2>Creating and Assembling Vectors<a class="headerlink" href="#creating-and-assembling-vectors" title="Permalink to this headline">¶</a></h2>
<p>PETSc currently provides two basic vector types: sequential and parallel
(MPI-based). To create a sequential vector with <code class="docutils literal notranslate"><span class="pre">m</span></code> components, one
can use the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">m</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">x</span><span class="p">);</span>
</pre></div>
</div>
<p>To create a parallel vector one can either specify the number of
components that will be stored on each process or let PETSc decide. The
command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateMPI.html#VecCreateMPI">VecCreateMPI</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">m</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">M</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">x</span><span class="p">);</span>
</pre></div>
</div>
<p>creates a vector distributed over all processes in the communicator,
<code class="docutils literal notranslate"><span class="pre">comm</span></code>, where <code class="docutils literal notranslate"><span class="pre">m</span></code> indicates the number of components to store on the
local process, and <code class="docutils literal notranslate"><span class="pre">M</span></code> is the total number of vector components.
Either the local or global dimension, but not both, can be set to
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span></code> or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DETERMINE.html#PETSC_DETERMINE">PETSC_DETERMINE</a></span></code>, respectively, to indicate that
PETSc should decide or determine it. More generally, one can use the
routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/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/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">v</span><span class="p">);</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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</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">m</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">M</span><span class="p">);</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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</span><span class="p">);</span>
</pre></div>
</div>
<p>which automatically generates the appropriate vector type (sequential or
parallel) over all processes in <code class="docutils literal notranslate"><span class="pre">comm</span></code>. The option <code class="docutils literal notranslate"><span class="pre">-vec_type</span> <span class="pre">mpi</span></code>
can be used in conjunction with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreate.html#VecCreate">VecCreate</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/VecSetFromOptions.html#VecSetFromOptions">VecSetFromOptions</a>()</span></code> to specify the use of MPI vectors even for the
uniprocessor case.</p>
<p>We emphasize that all processes in <code class="docutils literal notranslate"><span class="pre">comm</span></code> <em>must</em> call the vector
creation routines, since these routines are collective over all
processes in the communicator. If you are not familiar with MPI
communicators, see the discussion in <a class="reference internal" href="getting_started.html#sec-writing"><span class="std std-ref">Writing PETSc Programs</span></a> on
page . In addition, if a sequence of <code class="docutils literal notranslate"><span class="pre">VecCreateXXX()</span></code> routines is
used, they must be called in the same order on each process in the
communicator.</p>
<p>One can assign a single value to all components of a vector with the
command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</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">value</span><span class="p">);</span>
</pre></div>
</div>
<p>Assigning values to individual components of the vector is more
complicated, in order to make it possible to write efficient parallel
code. Assigning a set of components is a two-step process: one first
calls</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">indices</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="o">*</span><span class="n">values</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">);</span>
</pre></div>
</div>
<p>any number of times on any or all of the processes. The argument <code class="docutils literal notranslate"><span class="pre">n</span></code>
gives the number of components being set in this insertion. The integer
array <code class="docutils literal notranslate"><span class="pre">indices</span></code> contains the <em>global component indices</em>, and
<code class="docutils literal notranslate"><span class="pre">values</span></code> is the array of values to be inserted. Any process can set
any components of the vector; PETSc ensures that they are automatically
stored in the correct location. Once all of the values have been
inserted with <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>, one must call</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAssemblyBegin.html#VecAssemblyBegin">VecAssemblyBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">);</span>
</pre></div>
</div>
<p>followed by</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAssemblyEnd.html#VecAssemblyEnd">VecAssemblyEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">);</span>
</pre></div>
</div>
<p>to perform any needed message passing of nonlocal components. In order
to allow the overlap of communication and calculation, the user’s code
can perform any series of other actions between these two calls while
the messages are in transition.</p>
<p>Example usage 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> may be found in
<code class="docutils literal notranslate"><span class="pre">$PETSC_DIR/src/vec/vec/tutorials/ex2.c</span></code> or <code class="docutils literal notranslate"><span class="pre">ex2f.F</span></code>.</p>
<p>Often, rather than inserting elements in a vector, one may wish to add
values. This process is also done with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">indices</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="o">*</span><span class="n">values</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span><span class="p">);</span>
</pre></div>
</div>
<p>Again one must call the assembly routines <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAssemblyBegin.html#VecAssemblyBegin">VecAssemblyBegin</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/VecAssemblyEnd.html#VecAssemblyEnd">VecAssemblyEnd</a>()</span></code> after all of the values have been added. Note that
addition and insertion calls to <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> <em>cannot</em> be mixed.
Instead, one must add and insert vector elements in phases, with
intervening calls to the assembly routines. This phased assembly
procedure overcomes the nondeterministic behavior that would occur if
two different processes generated values for the same location, with one
process adding while the other is inserting its value. (In this case the
addition and insertion actions could be performed in either order, thus
resulting in different values at the particular location. Since PETSc
does not allow the simultaneous use of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</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/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span></code> this nondeterministic behavior will not occur in PETSc.)</p>
<p>You can call <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetValues.html#VecGetValues">VecGetValues</a>()</span></code> to pull local values from a vector (but
not off-process values), an alternative method for extracting some
components of a vector are the vector scatter routines. See
<a class="reference internal" href="#sec-scatter"><span class="std std-ref">Scatters and Gathers</span></a> for details; see also below for
<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>.</p>
<p>One can examine a vector with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">v</span><span class="p">);</span>
</pre></div>
</div>
<p>To print the vector to the screen, one can use the viewer
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_WORLD.html#PETSC_VIEWER_STDOUT_WORLD">PETSC_VIEWER_STDOUT_WORLD</a></span></code>, which ensures that parallel vectors are
printed correctly to <code class="docutils literal notranslate"><span class="pre">stdout</span></code>. To display the vector in an X-window,
one can use the default X-windows viewer <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_DRAW_WORLD.html#PETSC_VIEWER_DRAW_WORLD">PETSC_VIEWER_DRAW_WORLD</a></span></code>, or
one can create a viewer with the routine <code class="docutils literal notranslate"><span class="pre">PetscViewerDrawOpenX()</span></code>. A
variety of viewers are discussed further in
<a class="reference internal" href="other.html#sec-viewers"><span class="std std-ref">Viewers: Looking at PETSc Objects</span></a>.</p>
<p>To create a new vector of the same format as an existing vector, one
uses the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">old</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">new</span><span class="p">);</span>
</pre></div>
</div>
<p>To create several new vectors of the same format as an existing vector,
one uses the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/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/Vec.html#Vec">Vec</a></span> <span class="n">old</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">**</span><span class="n">new</span><span class="p">);</span>
</pre></div>
</div>
<p>This routine creates an array of pointers to vectors. The two routines
are very useful because they allow one to write library code that does
not depend on the particular format of the vectors being used. Instead,
the subroutines can automatically correctly create work vectors based on
the specified existing vector. As discussed in
<a class="reference internal" href="fortran.html#sec-fortvecd"><span class="std std-ref">Duplicating Multiple Vectors</span></a>, the Fortran interface for
<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.</p>
<p>When a vector is no longer needed, it should be destroyed with the
command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">x</span><span class="p">);</span>
</pre></div>
</div>
<p>To destroy an array of vectors, use the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/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/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">**</span><span class="n">vecs</span><span class="p">);</span>
</pre></div>
</div>
<p>Note that the Fortran interface for <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> differs
slightly, as described in <a class="reference internal" href="fortran.html#sec-fortvecd"><span class="std std-ref">Duplicating Multiple Vectors</span></a>.</p>
<p>It is also possible to create vectors that use an array provided by the
user, rather than having PETSc internally allocate the array space. Such
vectors can be created with the routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateSeqWithArray.html#VecCreateSeqWithArray">VecCreateSeqWithArray</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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">bs</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="o">*</span><span class="n">array</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">V</span><span class="p">);</span>
</pre></div>
</div>
<p>and</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateMPIWithArray.html#VecCreateMPIWithArray">VecCreateMPIWithArray</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">bs</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">N</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="o">*</span><span class="n">array</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">vv</span><span class="p">);</span>
</pre></div>
</div>
<p>Note that here one must provide the value <code class="docutils literal notranslate"><span class="pre">n</span></code>; it cannot be
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span></code> and the user is responsible for providing enough space
in the array; <code class="docutils literal notranslate"><span class="pre">n*sizeof(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a>)</span></code>.</p>
</div>
<div class="section" id="basic-vector-operations">
<span id="sec-vecbasic"></span><h2>Basic Vector Operations<a class="headerlink" href="#basic-vector-operations" title="Permalink to this headline">¶</a></h2>
<div class="docutils container" id="fig-vectorops">
<table class="docutils align-default" id="id3">
<caption><span class="caption-number">Table 1 </span><span class="caption-text">PETSc Vector Operations</span><a class="headerlink" href="#id3" title="Permalink to this table">¶</a></caption>
<colgroup>
<col style="width: 63%" />
<col style="width: 37%" />
</colgroup>
<thead>
<tr class="row-odd"><th class="head"><p><strong>Function Name</strong></p></th>
<th class="head"><p><strong>Operation</strong></p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAXPY.html#VecAXPY">VecAXPY</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">a,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x);</span></code></p></td>
<td><p><span class="math">\(y = y + a*x\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAYPX.html#VecAYPX">VecAYPX</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">a,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x);</span></code></p></td>
<td><p><span class="math">\(y = x + a*y\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecWAXPY.html#VecWAXPY">VecWAXPY</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span>  <span class="pre">w,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">a,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y);</span></code></p></td>
<td><p><span class="math">\(w = a*x + y\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAXPBY.html#VecAXPBY">VecAXPBY</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">a,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">b,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x);</span></code></p></td>
<td><p><span class="math">\(y = a*x + b*y\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScale.html#VecScale">VecScale</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">a);</span></code></p></td>
<td><p><span class="math">\(x = a*x\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDot.html#VecDot">VecDot</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r = \bar{x}^T*y\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecTDot.html#VecTDot">VecTDot</a>(</span>
<span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r = x'*y\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></span> <span class="pre">type,</span>  <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r = ||x||_{type}\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSum.html#VecSum">VecSum</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r = \sum x_{i}\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCopy.html#VecCopy">VecCopy</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y);</span></code></p></td>
<td><p><span class="math">\(y = x\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSwap.html#VecSwap">VecSwap</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y);</span></code></p></td>
<td><p><span class="math">\(y = x\)</span> while
<span class="math">\(x = y\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecPointwiseMult.html#VecPointwiseMult">VecPointwiseMult</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">w,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y);</span></code></p></td>
<td><p><span class="math">\(w_{i} = x_{i}*y_{i}\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecPointwiseDivide.html#VecPointwiseDivide">VecPointwiseDivide</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">w,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y);</span></code></p></td>
<td><p><span class="math">\(w_{i} = x_{i}/y_{i}\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMDot.html#VecMDot">VecMDot</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="pre">n,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y[],<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r[i] = \bar{x}^T*y[i]\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMTDot.html#VecMTDot">VecMTDot</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="pre">n,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y[],<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r[i] = x^T*y[i]\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMAXPY.html#VecMAXPY">VecMAXPY</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">y,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="pre">n,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">*a,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x[]);</span></code></p></td>
<td><p><span class="math">\(y = y + \sum_i a_{i}*x[i]\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMax.html#VecMax">VecMax</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="pre">*idx,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r = \max x_{i}\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMin.html#VecMin">VecMin</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="pre">*idx,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="pre">*r);</span></code></p></td>
<td><p><span class="math">\(r = \min x_{i}\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAbs.html#VecAbs">VecAbs</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x);</span></code></p></td>
<td><p><span class="math">\(x_i = |x_{i}|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecReciprocal.html#VecReciprocal">VecReciprocal</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x);</span></code></p></td>
<td><p><span class="math">\(x_i = 1/x_{i}\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecShift.html#VecShift">VecShift</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">s);</span></code></p></td>
<td><p><span class="math">\(x_i = s + x_{i}\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="pre">x,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="pre">alpha);</span></code></p></td>
<td><p><span class="math">\(x_i = \alpha\)</span></p></td>
</tr>
</tbody>
</table>
</div>
<p>As listed in the table, we have chosen certain
basic vector operations to support within the PETSc vector library.
These operations were selected because they often arise in application
codes. The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></span></code> argument to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a>()</span></code> is one of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_1.html#NORM_1">NORM_1</a></span></code>,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_2.html#NORM_2">NORM_2</a></span></code>, or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_INFINITY.html#NORM_INFINITY">NORM_INFINITY</a></span></code>. The 1-norm is <span class="math">\(\sum_i |x_{i}|\)</span>,
the 2-norm is <span class="math">\(( \sum_{i} x_{i}^{2})^{1/2}\)</span> and the infinity norm
is <span class="math">\(\max_{i} |x_{i}|\)</span>.</p>
<p>For parallel vectors that are distributed across the processes by
ranges, it is possible to determine a process’s local range with the
routine</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetOwnershipRange.html#VecGetOwnershipRange">VecGetOwnershipRange</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">vec</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">low</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">high</span><span class="p">);</span>
</pre></div>
</div>
<p>The argument <code class="docutils literal notranslate"><span class="pre">low</span></code> indicates the first component owned by the local
process, while <code class="docutils literal notranslate"><span class="pre">high</span></code> specifies <em>one more than</em> the last owned by the
local process. This command is useful, for instance, in assembling
parallel vectors.</p>
<p>On occasion, the user needs to access the actual elements of the vector.
The routine <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> returns a pointer to the elements local to
the process:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/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/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</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="o">**</span><span class="n">array</span><span class="p">);</span>
</pre></div>
</div>
<p>When access to the array is no longer needed, the user should call</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</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="o">**</span><span class="n">array</span><span class="p">);</span>
</pre></div>
</div>
<p>If the values do not need to be modified, the routines
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</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/VecRestoreArrayRead.html#VecRestoreArrayRead">VecRestoreArrayRead</a>()</span></code> provide read-only
access and should be used instead.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</span><span class="p">,</span> <span class="k">const</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="o">**</span><span class="n">array</span><span class="p">);</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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</span><span class="p">,</span> <span class="k">const</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="o">**</span><span class="n">array</span><span class="p">);</span>
</pre></div>
</div>
<p>Minor differences exist in the Fortran interface for <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>
and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a>()</span></code>, as discussed in
<a class="reference internal" href="fortran.html#sec-fortranarrays"><span class="std std-ref">Array Arguments</span></a>. It is important to note that
<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> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a>()</span></code> do <em>not</em> copy the vector
elements; they merely give users direct access to the vector elements.
Thus, these routines require essentially no time to call and can be used
efficiently.</p>
<p>The number of elements stored locally can be accessed with</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetLocalSize.html#VecGetLocalSize">VecGetLocalSize</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">size</span><span class="p">);</span>
</pre></div>
</div>
<p>The global vector length can be determined by</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetSize.html#VecGetSize">VecGetSize</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">size</span><span class="p">);</span>
</pre></div>
</div>
<p>In addition to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDot.html#VecDot">VecDot</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/VecMDot.html#VecMDot">VecMDot</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/VecNorm.html#VecNorm">VecNorm</a>()</span></code>, PETSc
provides split phase versions of these that allow several independent
inner products and/or norms to share the same communication (thus
improving parallel efficiency). For example, one may have code such as</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDot.html#VecDot">VecDot</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">y</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="o">*</span><span class="n">dot</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMDot.html#VecMDot">VecMDot</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</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">nv</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">y</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="o">*</span><span class="n">dot</span><span class="p">);</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"><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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">norm2</span><span class="p">);</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"><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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_1.html#NORM_1">NORM_1</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">norm1</span><span class="p">);</span>
</pre></div>
</div>
<p>This code works fine, but it performs three separate parallel
communication operations. Instead, one can write</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDotBegin.html#VecDotBegin">VecDotBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">y</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="o">*</span><span class="n">dot</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMDotBegin.html#VecMDotBegin">VecMDotBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</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">nv</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">y</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="o">*</span><span class="n">dot</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNormBegin.html#VecNormBegin">VecNormBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">norm2</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNormBegin.html#VecNormBegin">VecNormBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_1.html#NORM_1">NORM_1</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">norm1</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDotEnd.html#VecDotEnd">VecDotEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">y</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="o">*</span><span class="n">dot</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMDotEnd.html#VecMDotEnd">VecMDotEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</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">nv</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">y</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="o">*</span><span class="n">dot</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNormEnd.html#VecNormEnd">VecNormEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">norm2</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNormEnd.html#VecNormEnd">VecNormEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NormType.html#NormType">NormType</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_1.html#NORM_1">NORM_1</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="o">*</span><span class="n">norm1</span><span class="p">);</span>
</pre></div>
</div>
<p>With this code, the communication is delayed until the first call to
<code class="docutils literal notranslate"><span class="pre">VecxxxEnd()</span></code> at which a single MPI reduction is used to communicate
all the required values. It is required that the calls to the
<code class="docutils literal notranslate"><span class="pre">VecxxxEnd()</span></code> are performed in the same order as the calls to the
<code class="docutils literal notranslate"><span class="pre">VecxxxBegin()</span></code>; however, if you mistakenly make the calls in the
wrong order, PETSc will generate an error informing you of this. There
are additional routines <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecTDotBegin.html#VecTDotBegin">VecTDotBegin</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/VecTDotEnd.html#VecTDotEnd">VecTDotEnd</a>()</span></code>,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMTDotBegin.html#VecMTDotBegin">VecMTDotBegin</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecMTDotEnd.html#VecMTDotEnd">VecMTDotEnd</a>()</span></code>.</p>
<p>Note: these routines use only MPI-1 functionality; they do not allow you
to overlap computation and communication (assuming no threads are
spawned within a MPI process). Once MPI-2 implementations are more
common we’ll improve these routines to allow overlap of inner product
and norm calculations with other calculations. Also currently these
routines only work for the PETSc built in vector types.</p>
</div>
<div class="section" id="indexing-and-ordering">
<span id="sec-indexingandordering"></span><h2>Indexing and Ordering<a class="headerlink" href="#indexing-and-ordering" title="Permalink to this headline">¶</a></h2>
<p>When writing parallel PDE codes, there is extra complexity caused by
having multiple ways of indexing (numbering) and ordering objects such
as vertices and degrees of freedom. For example, a grid generator or
partitioner may renumber the nodes, requiring adjustment of the other
data structures that refer to these objects; see Figure
<a class="reference internal" href="#fig-daao"><span class="std std-ref">Natural Ordering and PETSc Ordering for a 2D Distributed Array (Four
Processes)</span></a>. In addition, local numbering (on a single process)
of objects may be different than the global (cross-process) numbering.
PETSc provides a variety of tools to help to manage the mapping amongst
the various numbering systems. The two most basic are the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code>
(application ordering), which enables mapping between different global
(cross-process) numbering schemes and the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span></code>,
which allows mapping between local (on-process) and global
(cross-process) numbering.</p>
<div class="section" id="application-orderings">
<span id="sec-ao"></span><h3>Application Orderings<a class="headerlink" href="#application-orderings" title="Permalink to this headline">¶</a></h3>
<p>In many applications it is desirable to work with one or more
“orderings” (or numberings) of degrees of freedom, cells, nodes, etc.
Doing so in a parallel environment is complicated by the fact that each
process cannot keep complete lists of the mappings between different
orderings. In addition, the orderings used in the PETSc linear algebra
routines (often contiguous ranges) may not correspond to the “natural”
orderings for the application.</p>
<p>PETSc provides certain utility routines that allow one to deal cleanly
and efficiently with the various orderings. To define a new application
ordering (called an <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> in PETSc), one can call the routine</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasic.html#AOCreateBasic">AOCreateBasic</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="k">const</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">apordering</span><span class="p">[],</span><span class="k">const</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">petscordering</span><span class="p">[],</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="o">*</span><span class="n">ao</span><span class="p">);</span>
</pre></div>
</div>
<p>The arrays <code class="docutils literal notranslate"><span class="pre">apordering</span></code> and <code class="docutils literal notranslate"><span class="pre">petscordering</span></code>, respectively, contain a
list of integers in the application ordering and their corresponding
mapped values in the PETSc ordering. Each process can provide whatever
subset of the ordering it chooses, but multiple processes should never
contribute duplicate values. The argument <code class="docutils literal notranslate"><span class="pre">n</span></code> indicates the number of
local contributed values.</p>
<p>For example, consider a vector of length 5, where node 0 in the
application ordering corresponds to node 3 in the PETSc ordering. In
addition, nodes 1, 2, 3, and 4 of the application ordering correspond,
respectively, to nodes 2, 1, 4, and 0 of the PETSc ordering. We can
write this correspondence as</p>
<div class="math">
\[\{ 0, 1, 2, 3, 4 \} \to \{ 3, 2, 1, 4, 0 \}.
\]</div>
<p>The user can create the PETSc <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> mappings in a number of ways. For
example, if using two processes, one could call</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasic.html#AOCreateBasic">AOCreateBasic</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="mi">2</span><span class="p">,{</span><span class="mi">0</span><span class="p">,</span><span class="mi">3</span><span class="p">},{</span><span class="mi">3</span><span class="p">,</span><span class="mi">4</span><span class="p">},</span><span class="o">&</span><span class="n">ao</span><span class="p">);</span>
</pre></div>
</div>
<p>on the first process and</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasic.html#AOCreateBasic">AOCreateBasic</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="mi">3</span><span class="p">,{</span><span class="mi">1</span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="mi">4</span><span class="p">},{</span><span class="mi">2</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="mi">0</span><span class="p">},</span><span class="o">&</span><span class="n">ao</span><span class="p">);</span>
</pre></div>
</div>
<p>on the other process.</p>
<p>Once the application ordering has been created, it can be used with
either of the commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOPetscToApplication.html#AOPetscToApplication">AOPetscToApplication</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="n">ao</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">indices</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOApplicationToPetsc.html#AOApplicationToPetsc">AOApplicationToPetsc</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="n">ao</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">indices</span><span class="p">);</span>
</pre></div>
</div>
<p>Upon input, the <code class="docutils literal notranslate"><span class="pre">n</span></code>-dimensional array <code class="docutils literal notranslate"><span class="pre">indices</span></code> specifies the
indices to be mapped, while upon output, <code class="docutils literal notranslate"><span class="pre">indices</span></code> contains the mapped
values. Since we, in general, employ a parallel database for the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code>
mappings, it is crucial that all processes that called
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasic.html#AOCreateBasic">AOCreateBasic</a>()</span></code> also call these routines; these routines <em>cannot</em> be
called by just a subset of processes in the MPI communicator that was
used in the call to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasic.html#AOCreateBasic">AOCreateBasic</a>()</span></code>.</p>
<p>An alternative routine to create the application ordering, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code>, is</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasicIS.html#AOCreateBasicIS">AOCreateBasicIS</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">apordering</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">petscordering</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="o">*</span><span class="n">ao</span><span class="p">);</span>
</pre></div>
</div>
<p>where index sets (see <a class="reference internal" href="#sec-indexset"><span class="std std-ref">Index Sets</span></a>) are used
instead of integer arrays.</p>
<p>The mapping routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOPetscToApplicationIS.html#AOPetscToApplicationIS">AOPetscToApplicationIS</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="n">ao</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">indices</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOApplicationToPetscIS.html#AOApplicationToPetscIS">AOApplicationToPetscIS</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="n">ao</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">indices</span><span class="p">);</span>
</pre></div>
</div>
<p>will map index sets (<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> objects) between orderings. Both the
<code class="docutils literal notranslate"><span class="pre">AOXxxToYyy()</span></code> and <code class="docutils literal notranslate"><span class="pre">AOXxxToYyyIS()</span></code> routines can be used regardless
of whether the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> was created with a <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasic.html#AOCreateBasic">AOCreateBasic</a>()</span></code> or
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOCreateBasicIS.html#AOCreateBasicIS">AOCreateBasicIS</a>()</span></code>.</p>
<p>The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> context should be destroyed with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AODestroy.html#AODestroy">AODestroy</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="pre">*ao)</span></code> and
viewed with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AOView.html#AOView">AOView</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="pre">ao,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="pre">viewer)</span></code>.</p>
<p>Although we refer to the two orderings as “PETSc” and “application”
orderings, the user is free to use them both for application orderings
and to maintain relationships among a variety of orderings by employing
several <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> contexts.</p>
<p>The <code class="docutils literal notranslate"><span class="pre">AOxxToxx()</span></code> routines allow negative entries in the input integer
array. These entries are not mapped; they simply remain unchanged. This
functionality enables, for example, mapping neighbor lists that use
negative numbers to indicate nonexistent neighbors due to boundary
conditions, etc.</p>
</div>
<div class="section" id="local-to-global-mappings">
<span id="sec-islocaltoglobalmapping"></span><h3>Local to Global Mappings<a class="headerlink" href="#local-to-global-mappings" title="Permalink to this headline">¶</a></h3>
<p>In many applications one works with a global representation of a vector
(usually on a vector obtained with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateMPI.html#VecCreateMPI">VecCreateMPI</a>()</span></code>) and a local
representation of the same vector that includes ghost points required
for local computation. PETSc provides routines to help map indices from
a local numbering scheme to the PETSc global numbering scheme. This is
done via the following routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMappingCreate.html#ISLocalToGlobalMappingCreate">ISLocalToGlobalMappingCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">bs</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">N</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span><span class="o">*</span> <span class="n">globalnum</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscCopyMode.html#PetscCopyMode">PetscCopyMode</a></span> <span class="n">mode</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span><span class="o">*</span> <span class="n">ctx</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMappingApply.html#ISLocalToGlobalMappingApply">ISLocalToGlobalMappingApply</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">in</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">out</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMappingApplyIS.html#ISLocalToGlobalMappingApplyIS">ISLocalToGlobalMappingApplyIS</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">isin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span><span class="o">*</span> <span class="n">isout</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMappingDestroy.html#ISLocalToGlobalMappingDestroy">ISLocalToGlobalMappingDestroy</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>Here <code class="docutils literal notranslate"><span class="pre">N</span></code> denotes the number of local indices, <code class="docutils literal notranslate"><span class="pre">globalnum</span></code> contains
the global number of each local number, and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span></code>
is the resulting PETSc object that contains the information needed to
apply the mapping with either <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMappingApply.html#ISLocalToGlobalMappingApply">ISLocalToGlobalMappingApply</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/ISLocalToGlobalMappingApplyIS.html#ISLocalToGlobalMappingApplyIS">ISLocalToGlobalMappingApplyIS</a>()</span></code>.</p>
<p>Note that the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span></code> routines serve a different
purpose than the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> routines. In the former case they provide a
mapping from a local numbering scheme (including ghost points) to a
global numbering scheme, while in the latter they provide a mapping
between two global numbering schemes. In fact, many applications may use
both <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</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/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span></code> routines. The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> routines
are first used to map from an application global ordering (that has no
relationship to parallel processing etc.) to the PETSc ordering scheme
(where each process has a contiguous set of indices in the numbering).
Then in order to perform function or Jacobian evaluations locally on
each process, one works with a local numbering scheme that includes
ghost points. The mapping from this local numbering scheme back to the
global PETSc numbering can be handled with the
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span></code> routines.</p>
<p>If one is given a list of block indices in a global numbering, the
routine</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGlobalToLocalMappingApplyBlock.html#ISGlobalToLocalMappingApplyBlock">ISGlobalToLocalMappingApplyBlock</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGlobalToLocalMappingMode.html#ISGlobalToLocalMappingMode">ISGlobalToLocalMappingMode</a></span> <span class="n">type</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">nin</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">idxin</span><span class="p">[],</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">nout</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">idxout</span><span class="p">[]);</span>
</pre></div>
</div>
<p>will provide a new list of indices in the local numbering. Again,
negative values in <code class="docutils literal notranslate"><span class="pre">idxin</span></code> are left unmapped. But, in addition, if
<code class="docutils literal notranslate"><span class="pre">type</span></code> is set to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGlobalToLocalMappingMode.html#ISGlobalToLocalMappingMode">IS_GTOLM_MASK</a></span></code> , then <code class="docutils literal notranslate"><span class="pre">nout</span></code> is set to <code class="docutils literal notranslate"><span class="pre">nin</span></code>
and all global values in <code class="docutils literal notranslate"><span class="pre">idxin</span></code> that are not represented in the local
to global mapping are replaced by -1. When <code class="docutils literal notranslate"><span class="pre">type</span></code> is set to
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGlobalToLocalMappingMode.html#ISGlobalToLocalMappingMode">IS_GTOLM_DROP</a></span></code>, the values in <code class="docutils literal notranslate"><span class="pre">idxin</span></code> that are not represented
locally in the mapping are not included in <code class="docutils literal notranslate"><span class="pre">idxout</span></code>, so that
potentially <code class="docutils literal notranslate"><span class="pre">nout</span></code> is smaller than <code class="docutils literal notranslate"><span class="pre">nin</span></code>. One must pass in an array
long enough to hold all the indices. One can call
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGlobalToLocalMappingApplyBlock.html#ISGlobalToLocalMappingApplyBlock">ISGlobalToLocalMappingApplyBlock</a>()</span></code> with <code class="docutils literal notranslate"><span class="pre">idxout</span></code> equal to <code class="docutils literal notranslate"><span class="pre">NULL</span></code>
to determine the required length (returned in <code class="docutils literal notranslate"><span class="pre">nout</span></code>) and then
allocate the required space and call
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGlobalToLocalMappingApplyBlock.html#ISGlobalToLocalMappingApplyBlock">ISGlobalToLocalMappingApplyBlock</a>()</span></code> a second time to set the values.</p>
<p>Often it is convenient to set elements into a vector using the local
node numbering rather than the global node numbering (e.g., each process
may maintain its own sublist of vertices and elements and number them
locally). To set values into a vector with the local numbering, one must
first call</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetLocalToGlobalMapping.html#VecSetLocalToGlobalMapping">VecSetLocalToGlobalMapping</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>and then call</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValuesLocal.html#VecSetValuesLocal">VecSetValuesLocal</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="k">const</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">indices</span><span class="p">[],</span><span class="k">const</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">values</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>
</pre></div>
</div>
<p>Now the <code class="docutils literal notranslate"><span class="pre">indices</span></code> use the local numbering, rather than the global,
meaning the entries lie in <span class="math">\([0,n)\)</span> where <span class="math">\(n\)</span> is the local
size of the vector.</p>
</div>
</div>
<div class="section" id="structured-grids-using-distributed-arrays">
<span id="sec-struct"></span><h2>Structured Grids Using Distributed Arrays<a class="headerlink" href="#structured-grids-using-distributed-arrays" title="Permalink to this headline">¶</a></h2>
<p>Distributed arrays (DMDAs), which are used in conjunction with PETSc
vectors, are intended for use with <em>logically regular rectangular grids</em>
when communication of nonlocal data is needed before certain local
computations can occur. PETSc distributed arrays are designed only for
the case in which data can be thought of as being stored in a standard
multidimensional array; thus, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code>s are <em>not</em> intended for
parallelizing unstructured grid problems, etc. DAs are intended for
communicating vector (field) information; they are not intended for
storing matrices.</p>
<p>For example, a typical situation one encounters in solving PDEs in
parallel is that, to evaluate a local function, <code class="docutils literal notranslate"><span class="pre">f(x)</span></code>, each process
requires its local portion of the vector <code class="docutils literal notranslate"><span class="pre">x</span></code> as well as its ghost
points (the bordering portions of the vector that are owned by
neighboring processes). Figure <a class="reference internal" href="#fig-ghosts"><span class="std std-ref">Ghost Points for Two Stencil Types on the Seventh Process</span></a> illustrates the
ghost points for the seventh process of a two-dimensional, regular
parallel grid. Each box represents a process; the ghost points for the
seventh process’s local part of a parallel array are shown in gray.</p>
<div class="figure align-default" id="fig-ghosts">
<img alt="Ghost Points for Two Stencil Types on the Seventh Process" src="../_images/ghost.png" />
<p class="caption"><span class="caption-number">Fig. 2 </span><span class="caption-text">Ghost Points for Two Stencil Types on the Seventh Process</span><a class="headerlink" href="#fig-ghosts" title="Permalink to this image">¶</a></p>
</div>
<div class="section" id="creating-distributed-arrays">
<h3>Creating Distributed Arrays<a class="headerlink" href="#creating-distributed-arrays" title="Permalink to this headline">¶</a></h3>
<p>The PETSc <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> object manages the parallel communication required
while working with data stored in regular arrays. The actual data is
stored in appropriately sized vector objects; the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> object only
contains the parallel data layout information and communication
information, however it may be used to create vectors and matrices with
the proper layout.</p>
<p>One creates a distributed array communication data structure in two
dimensions with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDACreate2d.html#DMDACreate2d">DMDACreate2d</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DMBoundaryType</a></span> <span class="n">xperiod</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DMBoundaryType</a></span> <span class="n">yperiod</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDAStencilType</a></span> <span class="n">st</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">M</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">N</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">m</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">dof</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">s</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">lx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">ly</span><span class="p">,</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="o">*</span><span class="n">da</span><span class="p">);</span>
</pre></div>
</div>
<p>The arguments <code class="docutils literal notranslate"><span class="pre">M</span></code> and <code class="docutils literal notranslate"><span class="pre">N</span></code> indicate the global numbers of grid points
in each direction, while <code class="docutils literal notranslate"><span class="pre">m</span></code> and <code class="docutils literal notranslate"><span class="pre">n</span></code> denote the process partition in
each direction; <code class="docutils literal notranslate"><span class="pre">m*n</span></code> must equal the number of processes in the MPI
communicator, <code class="docutils literal notranslate"><span class="pre">comm</span></code>. Instead of specifying the process layout, one
may use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span></code> for <code class="docutils literal notranslate"><span class="pre">m</span></code> and <code class="docutils literal notranslate"><span class="pre">n</span></code> so that PETSc will
determine the partition using MPI. The type of periodicity of the array
is specified by <code class="docutils literal notranslate"><span class="pre">xperiod</span></code> and <code class="docutils literal notranslate"><span class="pre">yperiod</span></code>, which can be
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DM_BOUNDARY_NONE</a></span></code> (no periodicity), <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DM_BOUNDARY_PERIODIC</a></span></code>
(periodic in that direction), <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DM_BOUNDARY_TWIST</a></span></code> (periodic in that
direction, but identified in reverse order), <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DM_BOUNDARY_GHOSTED</a></span></code> ,
or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DM_BOUNDARY_MIRROR</a></span></code>. The argument <code class="docutils literal notranslate"><span class="pre">dof</span></code> indicates the number of
degrees of freedom at each array point, and <code class="docutils literal notranslate"><span class="pre">s</span></code> is the stencil width
(i.e., the width of the ghost point region). The optional arrays <code class="docutils literal notranslate"><span class="pre">lx</span></code>
and <code class="docutils literal notranslate"><span class="pre">ly</span></code> may contain the number of nodes along the x and y axis for
each cell, i.e. the dimension of <code class="docutils literal notranslate"><span class="pre">lx</span></code> is <code class="docutils literal notranslate"><span class="pre">m</span></code> and the dimension of
<code class="docutils literal notranslate"><span class="pre">ly</span></code> is <code class="docutils literal notranslate"><span class="pre">n</span></code>; alternately, <code class="docutils literal notranslate"><span class="pre">NULL</span></code> may be passed in.</p>
<p>Two types of distributed array communication data structures can be
created, as specified by <code class="docutils literal notranslate"><span class="pre">st</span></code>. Star-type stencils that radiate outward
only in the coordinate directions are indicated by
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_STAR</a></span></code>, while box-type stencils are specified by
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_BOX</a></span></code>. For example, for the two-dimensional case,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_STAR</a></span></code> with width 1 corresponds to the standard 5-point
stencil, while <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_BOX</a></span></code> with width 1 denotes the standard
9-point stencil. In both instances the ghost points are identical, the
only difference being that with star-type stencils certain ghost points
are ignored, decreasing substantially the number of messages sent. Note
that the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_STAR</a></span></code> stencils can save interprocess
communication in two and three dimensions.</p>
<p>These DMDA stencils have nothing directly to do with any finite
difference stencils one might chose to use for a discretization; they
only ensure that the correct values are in place for application of a
user-defined finite difference stencil (or any other discretization
technique).</p>
<p>The commands for creating distributed array communication data
structures in one and three dimensions are analogous:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDACreate1d.html#DMDACreate1d">DMDACreate1d</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DMBoundaryType</a></span> <span class="n">xperiod</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">M</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">w</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">s</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">lc</span><span class="p">,</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="o">*</span><span class="n">inra</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDACreate3d.html#DMDACreate3d">DMDACreate3d</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DMBoundaryType</a></span> <span class="n">xperiod</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DMBoundaryType</a></span> <span class="n">yperiod</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DMBoundaryType</a></span> <span class="n">zperiod</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDAStencilType</a></span> <span class="n">stencil_type</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">M</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">N</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">P</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">m</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">p</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">w</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">s</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">lx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">ly</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">lz</span><span class="p">,</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="o">*</span><span class="n">inra</span><span class="p">);</span>
</pre></div>
</div>
<p>The routines to create distributed arrays are collective, so that all
processes in the communicator <code class="docutils literal notranslate"><span class="pre">comm</span></code> must call <code class="docutils literal notranslate"><span class="pre">DACreateXXX()</span></code>.</p>
</div>
<div class="section" id="local-global-vectors-and-scatters">
<h3>Local/Global Vectors and Scatters<a class="headerlink" href="#local-global-vectors-and-scatters" title="Permalink to this headline">¶</a></h3>
<p>Each <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> object defines the layout of two vectors: a distributed
global vector and a local vector that includes room for the appropriate
ghost points. The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> object provides information about the size
and layout of these vectors, but does not internally allocate any
associated storage space for field values. Instead, the user can create
vector objects that use the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> layout information with the
routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateGlobalVector.html#DMCreateGlobalVector">DMCreateGlobalVector</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">g</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateLocalVector.html#DMCreateLocalVector">DMCreateLocalVector</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">l</span><span class="p">);</span>
</pre></div>
</div>
<p>These vectors will generally serve as the building blocks for local and
global PDE solutions, etc. If additional vectors with such layout
information are needed in a code, they can be obtained by duplicating
<code class="docutils literal notranslate"><span class="pre">l</span></code> or <code class="docutils literal notranslate"><span class="pre">g</span></code> via <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a>()</span></code> or <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>.</p>
<p>We emphasize that a distributed array provides the information needed to
communicate the ghost value information between processes. In most
cases, several different vectors can share the same communication
information (or, in other words, can share a given <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code>). The design
of the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> object makes this easy, as each <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> operation may
operate on vectors of the appropriate size, as obtained via
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateLocalVector.html#DMCreateLocalVector">DMCreateLocalVector</a>()</span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateGlobalVector.html#DMCreateGlobalVector">DMCreateGlobalVector</a>()</span></code> or as produced
by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a>()</span></code>. As such, the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> scatter/gather operations
(e.g., <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGlobalToLocalBegin.html#DMGlobalToLocalBegin">DMGlobalToLocalBegin</a>()</span></code>) require vector input/output
arguments, as discussed below.</p>
<p>PETSc currently provides no container for multiple arrays sharing the
same distributed array communication; note, however, that the <code class="docutils literal notranslate"><span class="pre">dof</span></code>
parameter handles many cases of interest.</p>
<p>At certain stages of many applications, there is a need to work on a
local portion of the vector, including the ghost points. This may be
done by scattering a global vector into its local parts by using the
two-stage commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGlobalToLocalBegin.html#DMGlobalToLocalBegin">DMGlobalToLocalBegin</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">g</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n">iora</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGlobalToLocalEnd.html#DMGlobalToLocalEnd">DMGlobalToLocalEnd</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">g</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n">iora</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">);</span>
</pre></div>
</div>
<p>which allow the overlap of communication and computation. Since the
global and local vectors, given by <code class="docutils literal notranslate"><span class="pre">g</span></code> and <code class="docutils literal notranslate"><span class="pre">l</span></code>, respectively, must
be compatible with the distributed array, <code class="docutils literal notranslate"><span class="pre">da</span></code>, they should be
generated by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateGlobalVector.html#DMCreateGlobalVector">DMCreateGlobalVector</a>()</span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateLocalVector.html#DMCreateLocalVector">DMCreateLocalVector</a>()</span></code>
(or be duplicates of such a vector obtained via <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a>()</span></code>). The
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span></code> can be either <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span></code> or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span></code>.</p>
<p>One can scatter the local patches into the distributed vector with the
command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMLocalToGlobal.html#DMLocalToGlobal">DMLocalToGlobal</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n">mode</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">g</span><span class="p">);</span>
</pre></div>
</div>
<p>or the commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMLocalToGlobalBegin.html#DMLocalToGlobalBegin">DMLocalToGlobalBegin</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n">mode</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">g</span><span class="p">);</span>
<span class="cm">/* (Computation to overlap with communication) */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMLocalToGlobalEnd.html#DMLocalToGlobalEnd">DMLocalToGlobalEnd</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n">mode</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">g</span><span class="p">);</span>
</pre></div>
</div>
<p>In general this is used with an <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span></code> of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span></code>,
because if one wishes to insert values into the global vector they
should just access the global vector directly and put in the values.</p>
<p>A third type of distributed array scatter is from a local vector
(including ghost points that contain irrelevant values) to a local
vector with correct ghost point values. This scatter may be done with
the commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMLocalToLocalBegin.html#DMLocalToLocalBegin">DMLocalToLocalBegin</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l1</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n">iora</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l2</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMLocalToLocalEnd.html#DMLocalToLocalEnd">DMLocalToLocalEnd</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l1</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n">iora</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l2</span><span class="p">);</span>
</pre></div>
</div>
<p>Since both local vectors, <code class="docutils literal notranslate"><span class="pre">l1</span></code> and <code class="docutils literal notranslate"><span class="pre">l2</span></code>, must be compatible with the
distributed array, <code class="docutils literal notranslate"><span class="pre">da</span></code>, they should be generated by
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateLocalVector.html#DMCreateLocalVector">DMCreateLocalVector</a>()</span></code> (or be duplicates of such vectors obtained via
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a>()</span></code>). The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span></code> can be either <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span></code> or
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span></code>.</p>
<p>It is possible to directly access the vector scatter contexts (see
below) used in the local-to-global (<code class="docutils literal notranslate"><span class="pre">ltog</span></code>), global-to-local
(<code class="docutils literal notranslate"><span class="pre">gtol</span></code>), and local-to-local (<code class="docutils literal notranslate"><span class="pre">ltol</span></code>) scatters with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetScatter.html#DMDAGetScatter">DMDAGetScatter</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="o">*</span><span class="n">ltog</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="o">*</span><span class="n">gtol</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="o">*</span><span class="n">ltol</span><span class="p">);</span>
</pre></div>
</div>
<p>Most users should not need to use these contexts.</p>
</div>
<div class="section" id="local-ghosted-work-vectors">
<h3>Local (Ghosted) Work Vectors<a class="headerlink" href="#local-ghosted-work-vectors" title="Permalink to this headline">¶</a></h3>
<p>In most applications the local ghosted vectors are only needed during
user “function evaluations”. PETSc provides an easy, light-weight
(requiring essentially no CPU time) way to obtain these work vectors and
return them when they are no longer needed. This is done with the
routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGetLocalVector.html#DMGetLocalVector">DMGetLocalVector</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">l</span><span class="p">);</span>
<span class="p">...</span> <span class="n">use</span> <span class="n">the</span> <span class="n">local</span> <span class="n">vector</span> <span class="n">l</span> <span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMRestoreLocalVector.html#DMRestoreLocalVector">DMRestoreLocalVector</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">l</span><span class="p">);</span>
</pre></div>
</div>
</div>
<div class="section" id="accessing-the-vector-entries-for-dmda-vectors">
<h3>Accessing the Vector Entries for DMDA Vectors<a class="headerlink" href="#accessing-the-vector-entries-for-dmda-vectors" title="Permalink to this headline">¶</a></h3>
<p>PETSc provides an easy way to set values into the DMDA Vectors and
access them using the natural grid indexing. This is done with the
routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
<span class="p">...</span> <span class="n">use</span> <span class="n">the</span> <span class="n">array</span> <span class="n">indexing</span> <span class="n">it</span> <span class="n">with</span> <span class="mi">1</span> <span class="n">or</span> <span class="mi">2</span> <span class="n">or</span> <span class="mi">3</span> <span class="n">dimensions</span> <span class="p">...</span>
<span class="p">...</span> <span class="n">depending</span> <span class="n">on</span> <span class="n">the</span> <span class="n">dimension</span> <span class="n">of</span> <span class="n">the</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span> <span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArrayRead.html#DMDAVecGetArrayRead">DMDAVecGetArrayRead</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
<span class="p">...</span> <span class="n">use</span> <span class="n">the</span> <span class="n">array</span> <span class="n">indexing</span> <span class="n">it</span> <span class="n">with</span> <span class="mi">1</span> <span class="n">or</span> <span class="mi">2</span> <span class="n">or</span> <span class="mi">3</span> <span class="n">dimensions</span> <span class="p">...</span>
<span class="p">...</span> <span class="n">depending</span> <span class="n">on</span> <span class="n">the</span> <span class="n">dimension</span> <span class="n">of</span> <span class="n">the</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span> <span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArrayRead.html#DMDAVecRestoreArrayRead">DMDAVecRestoreArrayRead</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
</pre></div>
</div>
<p>and</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArrayDOF.html#DMDAVecGetArrayDOF">DMDAVecGetArrayDOF</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
<span class="p">...</span> <span class="n">use</span> <span class="n">the</span> <span class="n">array</span> <span class="n">indexing</span> <span class="n">it</span> <span class="n">with</span> <span class="mi">1</span> <span class="n">or</span> <span class="mi">2</span> <span class="n">or</span> <span class="mi">3</span> <span class="n">dimensions</span> <span class="p">...</span>
<span class="p">...</span> <span class="n">depending</span> <span class="n">on</span> <span class="n">the</span> <span class="n">dimension</span> <span class="n">of</span> <span class="n">the</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span> <span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArrayDOF.html#DMDAVecRestoreArrayDOF">DMDAVecRestoreArrayDOF</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArrayDOFRead.html#DMDAVecGetArrayDOFRead">DMDAVecGetArrayDOFRead</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
<span class="p">...</span> <span class="n">use</span> <span class="n">the</span> <span class="n">array</span> <span class="n">indexing</span> <span class="n">it</span> <span class="n">with</span> <span class="mi">1</span> <span class="n">or</span> <span class="mi">2</span> <span class="n">or</span> <span class="mi">3</span> <span class="n">dimensions</span> <span class="p">...</span>
<span class="p">...</span> <span class="n">depending</span> <span class="n">on</span> <span class="n">the</span> <span class="n">dimension</span> <span class="n">of</span> <span class="n">the</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span> <span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArrayDOFRead.html#DMDAVecRestoreArrayDOFRead">DMDAVecRestoreArrayDOFRead</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">l</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">array</span><span class="p">);</span>
</pre></div>
</div>
<p>where <code class="docutils literal notranslate"><span class="pre">array</span></code> is a multidimensional C array with the same dimension as
<code class="docutils literal notranslate"><span class="pre">da</span></code>. The vector <code class="docutils literal notranslate"><span class="pre">l</span></code> can be either a global vector or a local
vector. The <code class="docutils literal notranslate"><span class="pre">array</span></code> is accessed using the usual <em>global</em> indexing on
the entire grid, but the user may <em>only</em> refer to the local and ghost
entries of this array as all other entries are undefined. For example,
for a scalar problem in two dimensions one could use</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="o">**</span><span class="n">f</span><span class="p">,</span><span class="o">**</span><span class="n">u</span><span class="p">;</span>
<span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">local</span><span class="p">,</span><span class="o">&</span><span class="n">u</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">global</span><span class="p">,</span><span class="o">&</span><span class="n">f</span><span class="p">);</span>
<span class="p">...</span>
<span class="n">f</span><span class="p">[</span><span class="n">i</span><span class="p">][</span><span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">u</span><span class="p">[</span><span class="n">i</span><span class="p">][</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="p">...</span>
<span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">local</span><span class="p">,</span><span class="o">&</span><span class="n">u</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">global</span><span class="p">,</span><span class="o">&</span><span class="n">f</span><span class="p">);</span>
</pre></div>
</div>
<p>The recommended approach for multi-component PDEs is to declare a
<code class="docutils literal notranslate"><span class="pre">struct</span></code> representing the fields defined at each node of the grid,
e.g.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="k">typedef</span> <span class="k">struct</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">u</span><span class="p">,</span><span class="n">v</span><span class="p">,</span><span class="n">omega</span><span class="p">,</span><span class="n">temperature</span><span class="p">;</span>
<span class="p">}</span> <span class="n">Node</span><span class="p">;</span>
</pre></div>
</div>
<p>and write residual evaluation using</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">Node</span> <span class="o">**</span><span class="n">f</span><span class="p">,</span><span class="o">**</span><span class="n">u</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">local</span><span class="p">,</span><span class="o">&</span><span class="n">u</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">global</span><span class="p">,</span><span class="o">&</span><span class="n">f</span><span class="p">);</span>
<span class="p">...</span>
<span class="n">f</span><span class="p">[</span><span class="n">i</span><span class="p">][</span><span class="n">j</span><span class="p">].</span><span class="n">omega</span> <span class="o">=</span> <span class="p">...</span>
<span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">local</span><span class="p">,</span><span class="o">&</span><span class="n">u</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">global</span><span class="p">,</span><span class="o">&</span><span class="n">f</span><span class="p">);</span>
</pre></div>
</div>
<p>See
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex5.c.html">SNES Tutorial ex5</a>
for a complete example and see
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex19.c.html">SNES Tutorial ex19</a>
for an example for a multi-component PDE.</p>
</div>
<div class="section" id="grid-information">
<h3>Grid Information<a class="headerlink" href="#grid-information" title="Permalink to this headline">¶</a></h3>
<p>The global indices of the lower left corner of the local portion of the
array as well as the local array size can be obtained with the commands</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetCorners.html#DMDAGetCorners">DMDAGetCorners</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">z</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">m</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">p</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetGhostCorners.html#DMDAGetGhostCorners">DMDAGetGhostCorners</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">z</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">m</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">p</span><span class="p">);</span>
</pre></div>
</div>
<p>The first version excludes any ghost points, while the second version
includes them. The routine <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetGhostCorners.html#DMDAGetGhostCorners">DMDAGetGhostCorners</a>()</span></code> deals with the fact
that subarrays along boundaries of the problem domain have ghost points
only on their interior edges, but not on their boundary edges.</p>
<p>When either type of stencil is used, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_STAR</a></span></code> or
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_BOX</a></span></code>, the local vectors (with the ghost points)
represent rectangular arrays, including the extra corner elements in the
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_STAR</a></span></code> case. This configuration provides simple access to
the elements by employing two- (or three-) dimensional indexing. The
only difference between the two cases is that when <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAStencilType.html#DMDAStencilType">DMDA_STENCIL_STAR</a></span></code>
is used, the extra corner components are <em>not</em> scattered between the
processes and thus contain undefined values that should <em>not</em> be used.</p>
<p>To assemble global stiffness matrices, one can use these global indices
with <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> or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValuesStencil.html#MatSetValuesStencil">MatSetValuesStencil</a>()</span></code>. Alternately, the
global node number of each local node, including the ghost nodes, can be
obtained by calling</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGetLocalToGlobalMapping.html#DMGetLocalToGlobalMapping">DMGetLocalToGlobalMapping</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="o">*</span><span class="n">map</span><span class="p">);</span>
</pre></div>
</div>
<p>followed by</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetLocalToGlobalMapping.html#VecSetLocalToGlobalMapping">VecSetLocalToGlobalMapping</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">v</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="n">map</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetLocalToGlobalMapping.html#MatSetLocalToGlobalMapping">MatSetLocalToGlobalMapping</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="n">A</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="n">rmapping</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</a></span> <span class="n">cmapping</span><span class="p">);</span>
</pre></div>
</div>
<p>Now entries may be added to the vector and matrix using the local
numbering and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValuesLocal.html#VecSetValuesLocal">VecSetValuesLocal</a>()</span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValuesLocal.html#MatSetValuesLocal">MatSetValuesLocal</a>()</span></code>.</p>
<p>Since the global ordering that PETSc uses to manage its parallel vectors
(and matrices) does not usually correspond to the “natural” ordering of
a two- or three-dimensional array, the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code> structure provides an
application ordering <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span></code> (see <a class="reference internal" href="#sec-ao"><span class="std std-ref">Application Orderings</span></a>) that maps
between the natural ordering on a rectangular grid and the ordering
PETSc uses to parallelize. This ordering context can be obtained with
the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetAO.html#DMDAGetAO">DMDAGetAO</a></span><span class="p">(</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">da</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/AO/AO.html#AO">AO</a></span> <span class="o">*</span><span class="n">ao</span><span class="p">);</span>
</pre></div>
</div>
<p>In Figure <a class="reference internal" href="#fig-daao"><span class="std std-ref">Natural Ordering and PETSc Ordering for a 2D Distributed Array (Four
Processes)</span></a> we indicate the orderings for a
two-dimensional distributed array, divided among four processes.</p>
<div class="figure align-default" id="fig-daao">
<img alt="Natural Ordering and PETSc Ordering for a 2D Distributed Array (Four Processes)" src="../_images/danumbering.png" />
<p class="caption"><span class="caption-number">Fig. 3 </span><span class="caption-text">Natural Ordering and PETSc Ordering for a 2D Distributed Array (Four
Processes)</span><a class="headerlink" href="#fig-daao" title="Permalink to this image">¶</a></p>
</div>
<p>The example
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex5.c.html">SNES Tutorial ex5</a>
illustrates the use of a distributed array in the solution of a
nonlinear problem. The analogous Fortran program is
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex5f.F90.html">SNES Tutorial ex5f</a>;
see <a class="reference internal" href="snes.html#chapter-snes"><span class="std std-ref">SNES: Nonlinear Solvers</span></a> for a discussion of the
nonlinear solvers.</p>
</div>
<div class="section" id="staggered-grids">
<h3>Staggered Grids<a class="headerlink" href="#staggered-grids" title="Permalink to this headline">¶</a></h3>
<p>For regular grids with staggered data (living on elements, faces, edges,
and/or vertices), the <code class="docutils literal notranslate"><span class="pre">DMStag</span></code> object is available. It behaves much
like <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a></span></code>; see the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMSTAG/DMSTAG.html#DMSTAG">DMSTAG</a></span></code> manual page for more information.</p>
</div>
</div>
<div class="section" id="vectors-related-to-unstructured-grids">
<span id="sec-unstruct"></span><h2>Vectors Related to Unstructured Grids<a class="headerlink" href="#vectors-related-to-unstructured-grids" title="Permalink to this headline">¶</a></h2>
<div class="section" id="index-sets">
<span id="sec-indexset"></span><h3>Index Sets<a class="headerlink" href="#index-sets" title="Permalink to this headline">¶</a></h3>
<p>To facilitate general vector scatters and gathers used, for example, in
updating ghost points for problems defined on unstructured grids <a class="footnote-reference brackets" href="#id2" id="id1">1</a>,
PETSc employs the concept of an <em>index set</em>, via the <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> class. An
index set, which is a generalization of a set of integer indices, is
used to define scatters, gathers, and similar operations on vectors and
matrices.</p>
<p>The following command creates an index set based on a list of integers:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISCreateGeneral.html#ISCreateGeneral">ISCreateGeneral</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">indices</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscCopyMode.html#PetscCopyMode">PetscCopyMode</a></span> <span class="n">mode</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="o">*</span><span class="n">is</span><span class="p">);</span>
</pre></div>
</div>
<p>When <code class="docutils literal notranslate"><span class="pre">mode</span></code> is <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscCopyMode.html#PetscCopyMode">PETSC_COPY_VALUES</a></span></code>, this routine copies the <code class="docutils literal notranslate"><span class="pre">n</span></code>
indices passed to it by the integer array <code class="docutils literal notranslate"><span class="pre">indices</span></code>. Thus, the user
should be sure to free the integer array <code class="docutils literal notranslate"><span class="pre">indices</span></code> when it is no
longer needed, perhaps directly after the call to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISCreateGeneral.html#ISCreateGeneral">ISCreateGeneral</a>()</span></code>.
The communicator, <code class="docutils literal notranslate"><span class="pre">comm</span></code>, should consist of all processes that will be
using the <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>.</p>
<p>Another standard index set is defined by a starting point (<code class="docutils literal notranslate"><span class="pre">first</span></code>)
and a stride (<code class="docutils literal notranslate"><span class="pre">step</span></code>), and can be created with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISCreateStride.html#ISCreateStride">ISCreateStride</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">first</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">step</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="o">*</span><span class="n">is</span><span class="p">);</span>
</pre></div>
</div>
<p>Index sets can be destroyed with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISDestroy.html#ISDestroy">ISDestroy</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="o">&</span><span class="n">is</span><span class="p">);</span>
</pre></div>
</div>
<p>On rare occasions the user may need to access information directly from
an index set. Several commands assist in this process:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetSize.html#ISGetSize">ISGetSize</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">is</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">size</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISStrideGetInfo.html#ISStrideGetInfo">ISStrideGetInfo</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">is</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">first</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">stride</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/IS/IS.html#IS">IS</a></span> <span class="n">is</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">**</span><span class="n">indices</span><span class="p">);</span>
</pre></div>
</div>
<p>The function <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetIndices.html#ISGetIndices">ISGetIndices</a>()</span></code> returns a pointer to a list of the
indices in the index set. For certain index sets, this may be a
temporary array of indices created specifically for a given routine.
Thus, once the user finishes using the array of indices, the routine</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISRestoreIndices.html#ISRestoreIndices">ISRestoreIndices</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">is</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">**</span><span class="n">indices</span><span class="p">);</span>
</pre></div>
</div>
<p>should be called to ensure that the system can free the space it may
have used to generate the list of indices.</p>
<p>A blocked version of the index sets can be created with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISCreateBlock.html#ISCreateBlock">ISCreateBlock</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">bs</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">indices</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscCopyMode.html#PetscCopyMode">PetscCopyMode</a></span> <span class="n">mode</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="o">*</span><span class="n">is</span><span class="p">);</span>
</pre></div>
</div>
<p>This version is used for defining operations in which each element of
the index set refers to a block of <code class="docutils literal notranslate"><span class="pre">bs</span></code> vector entries. Related
routines analogous to those described above exist as well, including
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISBlockGetIndices.html#ISBlockGetIndices">ISBlockGetIndices</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISBlockGetSize.html#ISBlockGetSize">ISBlockGetSize</a>()</span></code>,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISBlockGetLocalSize.html#ISBlockGetLocalSize">ISBlockGetLocalSize</a>()</span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISGetBlockSize.html#ISGetBlockSize">ISGetBlockSize</a>()</span></code>. See the man pages for
details.</p>
</div>
<div class="section" id="scatters-and-gathers">
<span id="sec-scatter"></span><h3>Scatters and Gathers<a class="headerlink" href="#scatters-and-gathers" title="Permalink to this headline">¶</a></h3>
<p>PETSc vectors have full support for general scatters and gathers. One
can select any subset of the components of a vector to insert or add to
any subset of the components of another vector. We refer to these
operations as <em>generalized scatters</em>, though they are actually a
combination of scatters and gathers.</p>
<p>To copy selected components from one vector to another, one uses the
following set of commands:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterCreate.html#VecScatterCreate">VecScatterCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">ix</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">iy</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterBegin.html#VecScatterBegin">VecScatterBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></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/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterEnd.html#VecScatterEnd">VecScatterEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></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/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterDestroy.html#VecScatterDestroy">VecScatterDestroy</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>Here <code class="docutils literal notranslate"><span class="pre">ix</span></code> denotes the index set of the first vector, while <code class="docutils literal notranslate"><span class="pre">iy</span></code>
indicates the index set of the destination vector. The vectors can be
parallel or sequential. The only requirements are that the number of
entries in the index set of the first vector, <code class="docutils literal notranslate"><span class="pre">ix</span></code>, equals the number
in the destination index set, <code class="docutils literal notranslate"><span class="pre">iy</span></code>, and that the vectors be long
enough to contain all the indices referred to in the index sets. If both
<code class="docutils literal notranslate"><span class="pre">x</span></code> and <code class="docutils literal notranslate"><span class="pre">y</span></code> are parallel, their communicator must have the same set
of processes, but their process order can be different. The argument
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span></code> specifies that the vector elements will be inserted
into the specified locations of the destination vector, overwriting any
existing values. To add the components, rather than insert them, the
user should select the option <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span></code> instead of
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span></code>. One can also use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MAX_VALUES.html#MAX_VALUES">MAX_VALUES</a></span></code> or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MIN_VALUES.html#MIN_VALUES">MIN_VALUES</a></span></code> to
replace destination with the maximal or minimal of its current value and
the scattered values.</p>
<p>To perform a conventional gather operation, the user simply makes the
destination index set, <code class="docutils literal notranslate"><span class="pre">iy</span></code>, be a stride index set with a stride of
one. Similarly, a conventional scatter can be done with an initial
(sending) index set consisting of a stride. The scatter routines are
collective operations (i.e. all processes that own a parallel vector
<em>must</em> call the scatter routines). When scattering from a parallel
vector to sequential vectors, each process has its own sequential vector
that receives values from locations as indicated in its own index set.
Similarly, in scattering from sequential vectors to a parallel vector,
each process has its own sequential vector that makes contributions to
the parallel vector.</p>
<p><em>Caution</em>: When <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span></code> is used, if two different processes
contribute different values to the same component in a parallel vector,
either value may end up being inserted. When <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span></code> is used, the
correct sum is added to the correct location.</p>
<p>In some cases one may wish to “undo” a scatter, that is perform the
scatter backwards, switching the roles of the sender and receiver. This
is done by using</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterBegin.html#VecScatterBegin">VecScatterBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></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/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterEnd.html#VecScatterEnd">VecScatterEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">ctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></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/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span><span class="p">);</span>
</pre></div>
</div>
<p>Note that the roles of the first two arguments to these routines must be
swapped whenever the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span></code> option is used.</p>
<p>Once a <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span></code> object has been created it may be used with any
vectors that have the appropriate parallel data layout. That is, one can
call <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterBegin.html#VecScatterBegin">VecScatterBegin</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/VecScatterEnd.html#VecScatterEnd">VecScatterEnd</a>()</span></code> with different
vectors than used in the call to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterCreate.html#VecScatterCreate">VecScatterCreate</a>()</span></code> as long as they
have the same parallel layout (number of elements on each process are
the same). Usually, these “different” vectors would have been obtained
via calls to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a>()</span></code> from the original vectors used in the
call to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterCreate.html#VecScatterCreate">VecScatterCreate</a>()</span></code>.</p>
<p>There is a PETSc routine that is nearly the opposite 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>, that is, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetValues.html#VecGetValues">VecGetValues</a>()</span></code>, but it can only get
local values from the vector. To get off-process values, the user should
create a new vector where the components are to be stored, and then
perform the appropriate vector scatter. For example, if one desires to
obtain the values of the 100th and 200th entries of a parallel vector,
<code class="docutils literal notranslate"><span class="pre">p</span></code>, one could use a code such as that below. In this example, the
values of the 100th and 200th components are placed in the array values.
In this example each process now has the 100th and 200th component, but
obviously each process could gather any elements it needed, or none by
creating an index set with no entries.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">p</span><span class="p">,</span> <span class="n">x</span><span class="p">;</span> <span class="cm">/* initial vector, destination vector */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">scatter</span><span class="p">;</span> <span class="cm">/* scatter context */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/IS.html#IS">IS</a></span> <span class="n">from</span><span class="p">,</span> <span class="n">to</span><span class="p">;</span> <span class="cm">/* index sets that define the scatter */</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="o">*</span><span class="n">values</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">idx_from</span><span class="p">[]</span> <span class="o">=</span> <span class="p">{</span><span class="mi">100</span><span class="p">,</span><span class="mi">200</span><span class="p">},</span> <span class="n">idx_to</span><span class="p">[]</span> <span class="o">=</span> <span class="p">{</span><span class="mi">0</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/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="mi">2</span><span class="p">,</span><span class="o">&</span><span class="n">x</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISCreateGeneral.html#ISCreateGeneral">ISCreateGeneral</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="mi">2</span><span class="p">,</span><span class="n">idx_from</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscCopyMode.html#PetscCopyMode">PETSC_COPY_VALUES</a></span><span class="p">,</span><span class="o">&</span><span class="n">from</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISCreateGeneral.html#ISCreateGeneral">ISCreateGeneral</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="mi">2</span><span class="p">,</span><span class="n">idx_to</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscCopyMode.html#PetscCopyMode">PETSC_COPY_VALUES</a></span><span class="p">,</span><span class="o">&</span><span class="n">to</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterCreate.html#VecScatterCreate">VecScatterCreate</a></span><span class="p">(</span><span class="n">p</span><span class="p">,</span><span class="n">from</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="n">to</span><span class="p">,</span><span class="o">&</span><span class="n">scatter</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterBegin.html#VecScatterBegin">VecScatterBegin</a></span><span class="p">(</span><span class="n">scatter</span><span class="p">,</span><span class="n">p</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/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterEnd.html#VecScatterEnd">VecScatterEnd</a></span><span class="p">(</span><span class="n">scatter</span><span class="p">,</span><span class="n">p</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/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></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">x</span><span class="p">,</span><span class="o">&</span><span class="n">values</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISDestroy.html#ISDestroy">ISDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">from</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISDestroy.html#ISDestroy">ISDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">to</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterDestroy.html#VecScatterDestroy">VecScatterDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">scatter</span><span class="p">);</span>
</pre></div>
</div>
<p>The scatter comprises two stages, in order to allow overlap of
communication and computation. The introduction of the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span></code>
context allows the communication patterns for the scatter to be computed
once and then reused repeatedly. Generally, even setting up the
communication for a scatter requires communication; hence, it is best to
reuse such information when possible.</p>
</div>
<div class="section" id="scattering-ghost-values">
<h3>Scattering Ghost Values<a class="headerlink" href="#scattering-ghost-values" title="Permalink to this headline">¶</a></h3>
<p>Generalized scatters provide a very general method for managing the
communication of required ghost values for unstructured grid
computations. One scatters the global vector into a local “ghosted” work
vector, performs the computation on the local work vectors, and then
scatters back into the global solution vector. In the simplest case this
may be written as</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterBegin.html#VecScatterBegin">VecScatterBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">scatter</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">localin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterEnd.html#VecScatterEnd">VecScatterEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">scatter</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">localin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span><span class="p">);</span>
<span class="cm">/* For example, do local calculations from localin to localout */</span>
<span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterBegin.html#VecScatterBegin">VecScatterBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">scatter</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">localout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterEnd.html#VecScatterEnd">VecScatterEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatter.html#VecScatter">VecScatter</a></span> <span class="n">scatter</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">localout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span><span class="p">);</span>
</pre></div>
</div>
</div>
<div class="section" id="vectors-with-locations-for-ghost-values">
<h3>Vectors with Locations for Ghost Values<a class="headerlink" href="#vectors-with-locations-for-ghost-values" title="Permalink to this headline">¶</a></h3>
<p>There are two minor drawbacks to the basic approach described above:</p>
<ul class="simple">
<li><p>the extra memory requirement for the local work vector, <code class="docutils literal notranslate"><span class="pre">localin</span></code>,
which duplicates the memory in <code class="docutils literal notranslate"><span class="pre">globalin</span></code>, and</p></li>
<li><p>the extra time required to copy the local values from <code class="docutils literal notranslate"><span class="pre">localin</span></code> to
<code class="docutils literal notranslate"><span class="pre">globalin</span></code>.</p></li>
</ul>
<p>An alternative approach is to allocate global vectors with space
preallocated for the ghost values; this may be done with either</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateGhost.html#VecCreateGhost">VecCreateGhost</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">N</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">nghost</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">ghosts</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">vv</span><span class="p">)</span>
</pre></div>
</div>
<p>or</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateGhostWithArray.html#VecCreateGhostWithArray">VecCreateGhostWithArray</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</a></span> <span class="n">comm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">n</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">N</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">nghost</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="o">*</span><span class="n">ghosts</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="o">*</span><span class="n">array</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">vv</span><span class="p">)</span>
</pre></div>
</div>
<p>Here <code class="docutils literal notranslate"><span class="pre">n</span></code> is the number of local vector entries, <code class="docutils literal notranslate"><span class="pre">N</span></code> is the number of
global entries (or <code class="docutils literal notranslate"><span class="pre">NULL</span></code>) and <code class="docutils literal notranslate"><span class="pre">nghost</span></code> is the number of ghost
entries. The array <code class="docutils literal notranslate"><span class="pre">ghosts</span></code> is of size <code class="docutils literal notranslate"><span class="pre">nghost</span></code> and contains the
global vector location for each local ghost location. Using
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a>()</span></code> or <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> on a ghosted vector will
generate additional ghosted vectors.</p>
<p>In many ways, a ghosted vector behaves just like any other MPI vector
created by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreateMPI.html#VecCreateMPI">VecCreateMPI</a>()</span></code>. The difference is that the ghosted vector
has an additional “local” representation that allows one to access the
ghost locations. This is done through the call to</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostGetLocalForm.html#VecGhostGetLocalForm">VecGhostGetLocalForm</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">g</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">l</span><span class="p">);</span>
</pre></div>
</div>
<p>The vector <code class="docutils literal notranslate"><span class="pre">l</span></code> is a sequential representation of the parallel vector
<code class="docutils literal notranslate"><span class="pre">g</span></code> that shares the same array space (and hence numerical values); but
allows one to access the “ghost” values past “the end of the” array.
Note that one access the entries in <code class="docutils literal notranslate"><span class="pre">l</span></code> using the local numbering of
elements and ghosts, while they are accessed in <code class="docutils literal notranslate"><span class="pre">g</span></code> using the global
numbering.</p>
<p>A common usage of a ghosted vector is given by</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostUpdateBegin.html#VecGhostUpdateBegin">VecGhostUpdateBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostUpdateEnd.html#VecGhostUpdateEnd">VecGhostUpdateEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></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"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostGetLocalForm.html#VecGhostGetLocalForm">VecGhostGetLocalForm</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">localin</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostGetLocalForm.html#VecGhostGetLocalForm">VecGhostGetLocalForm</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">localout</span><span class="p">);</span>
<span class="p">...</span> <span class="n">Do</span> <span class="n">local</span> <span class="n">calculations</span> <span class="n">from</span> <span class="n">localin</span> <span class="n">to</span> <span class="n">localout</span> <span class="p">...</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostRestoreLocalForm.html#VecGhostRestoreLocalForm">VecGhostRestoreLocalForm</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalin</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">localin</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostRestoreLocalForm.html#VecGhostRestoreLocalForm">VecGhostRestoreLocalForm</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="o">*</span><span class="n">localout</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostUpdateBegin.html#VecGhostUpdateBegin">VecGhostUpdateBegin</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostUpdateEnd.html#VecGhostUpdateEnd">VecGhostUpdateEnd</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">globalout</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/InsertMode.html#InsertMode">InsertMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</a></span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/ScatterMode.html#ScatterMode">ScatterMode</a></span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span><span class="p">);</span>
</pre></div>
</div>
<p>The routines <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostUpdateBegin.html#VecGhostUpdateBegin">VecGhostUpdateBegin</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/VecGhostUpdateEnd.html#VecGhostUpdateEnd">VecGhostUpdateEnd</a>()</span></code> are
equivalent to the routines <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecScatterBegin.html#VecScatterBegin">VecScatterBegin</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/VecScatterEnd.html#VecScatterEnd">VecScatterEnd</a>()</span></code>
above except that since they are scattering into the ghost locations,
they do not need to copy the local vector values, which are already in
place. In addition, the user does not have to allocate the local work
vector, since the ghosted vector already has allocated slots to contain
the ghost values.</p>
<p>The input arguments <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</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/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</a></span></code> cause the
ghost values to be correctly updated from the appropriate process. The
arguments <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</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/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span></code> update the “local”
portions of the vector from all the other processes’ ghost values. This
would be appropriate, for example, when performing a finite element
assembly of a load vector. One can also use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MAX_VALUES.html#MAX_VALUES">MAX_VALUES</a></span></code> or
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/MIN_VALUES.html#MIN_VALUES">MIN_VALUES</a></span></code> with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/SCATTER_REVERSE.html#SCATTER_REVERSE">SCATTER_REVERSE</a></span></code>.</p>
<p><a class="reference internal" href="mat.html#sec-partitioning"><span class="std std-ref">Partitioning</span></a> discusses the important topic of
partitioning an unstructured grid.</p>
<dl class="footnote brackets">
<dt class="label" id="id2"><span class="brackets"><a class="fn-backref" href="#id1">1</a></span></dt>
<dd><p>Also see <code class="docutils literal notranslate"><span class="pre">DMPlex</span></code> (<a class="reference internal" href="dmplex.html#chapter-unstructured"><span class="std std-ref">DMPlex: Unstructured Grids in PETSc</span></a>), an abstraction for working with unstructured grids.</p>
</dd>
</dl>
</div>
</div>
</div>
</div>
</div>
</div>
<div class="clearer"></div>
</div>
<div class="related" role="navigation" aria-label="related navigation">
<h3>Navigation</h3>
<ul>
<li class="right" style="margin-right: 10px">
<a href="../genindex.html" title="General Index"
>index</a></li>
<li class="right" >
<a href="mat.html" title="Matrices"
>next</a> |</li>
<li class="right" >
<a href="programming.html" title="Programming with PETSc"
>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="programming.html" >Programming with PETSc</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>
|