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 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733 1734 1735 1736 1737 1738 1739 1740 1741 1742 1743 1744 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 1776 1777 1778 1779 1780 1781 1782 1783 1784 1785 1786 1787 1788 1789 1790 1791 1792 1793 1794 1795 1796 1797 1798 1799 1800 1801 1802 1803 1804 1805 1806 1807 1808 1809 1810 1811 1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822 1823 1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835 1836 1837 1838 1839 1840 1841 1842 1843 1844 1845 1846 1847 1848 1849 1850 1851 1852 1853 1854 1855 1856 1857 1858 1859 1860 1861 1862 1863 1864 1865 1866 1867 1868 1869 1870 1871 1872 1873 1874 1875 1876 1877 1878 1879 1880 1881 1882 1883 1884 1885 1886 1887 1888 1889 1890 1891 1892 1893 1894 1895 1896 1897 1898 1899 1900 1901 1902 1903 1904 1905 1906 1907 1908 1909 1910 1911 1912 1913 1914 1915 1916 1917 1918 1919 1920 1921 1922 1923 1924 1925 1926 1927 1928 1929 1930 1931 1932 1933 1934 1935 1936 1937 1938 1939 1940 1941 1942 1943 1944 1945 1946 1947 1948 1949 1950 1951 1952 1953 1954 1955 1956 1957 1958 1959 1960 1961 1962 1963 1964 1965 1966 1967 1968 1969 1970 1971 1972 1973 1974 1975 1976 1977 1978 1979 1980 1981 1982 1983 1984 1985 1986 1987 1988 1989 1990 1991 1992 1993 1994 1995 1996 1997 1998 1999 2000 2001 2002 2003 2004 2005 2006 2007 2008 2009 2010 2011 2012 2013 2014 2015 2016 2017 2018 2019 2020 2021 2022 2023 2024 2025 2026 2027 2028 2029 2030 2031 2032 2033 2034 2035 2036 2037 2038 2039 2040 2041 2042 2043 2044 2045 2046 2047 2048 2049 2050 2051 2052 2053 2054 2055 2056 2057 2058 2059 2060 2061 2062 2063 2064 2065 2066 2067 2068 2069 2070 2071 2072 2073 2074 2075 2076 2077 2078 2079 2080 2081 2082 2083 2084 2085 2086 2087 2088 2089 2090 2091 2092 2093 2094 2095 2096 2097 2098 2099 2100 2101 2102 2103 2104 2105 2106 2107 2108 2109 2110 2111 2112 2113 2114 2115 2116 2117 2118 2119 2120 2121 2122 2123 2124 2125 2126 2127 2128 2129 2130 2131 2132 2133 2134 2135 2136 2137 2138 2139 2140 2141 2142 2143 2144 2145 2146 2147 2148 2149 2150 2151 2152 2153 2154 2155 2156 2157 2158 2159 2160 2161 2162 2163 2164 2165 2166 2167 2168 2169 2170 2171 2172 2173 2174 2175 2176 2177 2178 2179 2180 2181 2182 2183 2184 2185 2186 2187 2188 2189 2190 2191 2192 2193 2194 2195 2196 2197 2198 2199 2200 2201 2202 2203 2204 2205 2206 2207 2208 2209 2210 2211 2212 2213 2214 2215 2216 2217 2218 2219 2220 2221 2222 2223 2224 2225 2226 2227 2228 2229 2230 2231 2232 2233 2234 2235 2236 2237 2238 2239 2240 2241 2242 2243 2244 2245 2246 2247 2248 2249 2250 2251 2252 2253 2254 2255 2256 2257 2258 2259 2260 2261 2262 2263 2264 2265 2266 2267 2268 2269 2270 2271 2272 2273 2274 2275 2276 2277 2278 2279 2280 2281 2282 2283 2284 2285 2286 2287 2288 2289 2290 2291 2292 2293 2294 2295 2296 2297 2298
|
<!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/snes.html" />
<meta charset="utf-8" />
<title>SNES: Nonlinear Solvers — 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="TS: Scalable ODE and DAE Solvers" href="ts.html" />
<link rel="prev" title="KSP: Linear System Solvers" href="ksp.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/snes.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="ts.html" title="TS: Scalable ODE and DAE Solvers"
accesskey="N">next</a> |</li>
<li class="right" >
<a href="ksp.html" title="KSP: Linear System Solvers"
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="#">SNES: Nonlinear Solvers</a><ul>
<li><a class="reference internal" href="#basic-snes-usage">Basic SNES Usage</a><ul>
<li><a class="reference internal" href="#nonlinear-function-evaluation">Nonlinear Function Evaluation</a></li>
<li><a class="reference internal" href="#jacobian-evaluation">Jacobian Evaluation</a></li>
</ul>
</li>
<li><a class="reference internal" href="#the-nonlinear-solvers">The Nonlinear Solvers</a><ul>
<li><a class="reference internal" href="#line-search-newton">Line Search Newton</a></li>
<li><a class="reference internal" href="#trust-region-methods">Trust Region Methods</a></li>
<li><a class="reference internal" href="#nonlinear-krylov-methods">Nonlinear Krylov Methods</a></li>
<li><a class="reference internal" href="#quasi-newton-methods">Quasi-Newton Methods</a></li>
<li><a class="reference internal" href="#the-full-approximation-scheme">The Full Approximation Scheme</a></li>
<li><a class="reference internal" href="#nonlinear-additive-schwarz">Nonlinear Additive Schwarz</a></li>
</ul>
</li>
<li><a class="reference internal" href="#general-options">General Options</a><ul>
<li><a class="reference internal" href="#convergence-tests">Convergence Tests</a></li>
<li><a class="reference internal" href="#convergence-monitoring">Convergence Monitoring</a></li>
<li><a class="reference internal" href="#checking-accuracy-of-derivatives">Checking Accuracy of Derivatives</a></li>
</ul>
</li>
<li><a class="reference internal" href="#inexact-newton-like-methods">Inexact Newton-like Methods</a></li>
<li><a class="reference internal" href="#matrix-free-methods">Matrix-Free Methods</a></li>
<li><a class="reference internal" href="#finite-difference-jacobian-approximations">Finite Difference Jacobian Approximations</a></li>
<li><a class="reference internal" href="#variational-inequalities">Variational Inequalities</a></li>
<li><a class="reference internal" href="#nonlinear-preconditioning">Nonlinear Preconditioning</a></li>
</ul>
</li>
</ul>
<h4>Previous topic</h4>
<p class="topless"><a href="ksp.html"
title="previous chapter">KSP: Linear System Solvers</a></p>
<h4>Next topic</h4>
<p class="topless"><a href="ts.html"
title="next chapter">TS: Scalable ODE and DAE Solvers</a></p>
<div role="note" aria-label="source link">
<h3>This Page</h3>
<ul class="this-page-menu">
<li><a href="../_sources/manual/snes.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="snes-nonlinear-solvers">
<span id="chapter-snes"></span><h1>SNES: Nonlinear Solvers<a class="headerlink" href="#snes-nonlinear-solvers" title="Permalink to this headline">¶</a></h1>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>This chapter is being cleaned up by Jed Brown. Contributions are welcome.</p>
</div>
<p>The solution of large-scale nonlinear problems pervades many facets of
computational science and demands robust and flexible solution
strategies. The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> library of PETSc provides a powerful suite of
data-structure-neutral numerical routines for such problems. Built on
top of the linear solvers and data structures discussed in preceding
chapters, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> enables the user to easily customize the nonlinear
solvers according to the application at hand. Also, the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>
interface is <em>identical</em> for the uniprocess and parallel cases; the only
difference in the parallel version is that each process typically forms
only its local contribution to various matrices and vectors.</p>
<p>The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> class includes methods for solving systems of nonlinear
equations of the form</p>
<div class="math" id="equation-fx0">
<span class="eqno">(3)<a class="headerlink" href="#equation-fx0" title="Permalink to this equation">¶</a></span>\[\mathbf{F}(\mathbf{x}) = 0,\]</div>
<p>where <span class="math">\(\mathbf{F}: \, \Re^n \to \Re^n\)</span>. Newton-like methods provide the
core of the package, including both line search and trust region
techniques. A suite of nonlinear Krylov methods and methods based upon
problem decomposition are also included. The solvers are discussed
further in <a class="reference internal" href="#sec-nlsolvers"><span class="std std-ref">The Nonlinear Solvers</span></a>. Following the PETSc design
philosophy, the interfaces to the various solvers are all virtually
identical. In addition, the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> software is completely flexible, so
that the user can at runtime change any facet of the solution process.</p>
<p>PETSc’s default method for solving the nonlinear equation is Newton’s
method. The general form of the <span class="math">\(n\)</span>-dimensional Newton’s method
for solving <a class="reference internal" href="#equation-fx0">(3)</a> is</p>
<div class="math" id="equation-newton">
<span class="eqno">(4)<a class="headerlink" href="#equation-newton" title="Permalink to this equation">¶</a></span>\[\mathbf{x}_{k+1} = \mathbf{x}_k - \mathbf{J}(\mathbf{x}_k)^{-1} \mathbf{F}(\mathbf{x}_k), \;\; k=0,1, \ldots,\]</div>
<p>where <span class="math">\(\mathbf{x}_0\)</span> is an initial approximation to the solution and
<span class="math">\(\mathbf{J}(\mathbf{x}_k) = \mathbf{F}'(\mathbf{x}_k)\)</span>, the Jacobian, is nonsingular at each
iteration. In practice, the Newton iteration <a class="reference internal" href="#equation-newton">(4)</a> is
implemented by the following two steps:</p>
<div class="math">
\[\begin{aligned}
1. & \text{(Approximately) solve} & \mathbf{J}(\mathbf{x}_k) \Delta \mathbf{x}_k &= -\mathbf{F}(\mathbf{x}_k). \\
2. & \text{Update} & \mathbf{x}_{k+1} &\gets \mathbf{x}_k + \Delta \mathbf{x}_k.
\end{aligned}\]</div>
<p>Other defect-correction algorithms can be implemented by using different
choices for <span class="math">\(J(\mathbf{x}_k)\)</span>.</p>
<div class="section" id="basic-snes-usage">
<span id="sec-snesusage"></span><h2>Basic SNES Usage<a class="headerlink" href="#basic-snes-usage" title="Permalink to this headline">¶</a></h2>
<p>In the simplest usage of the nonlinear solvers, the user must merely
provide a C, C++, or Fortran routine to evaluate the nonlinear function
<a class="reference internal" href="#equation-fx0">(3)</a>. The corresponding Jacobian matrix
can be approximated with finite differences. For codes that are
typically more efficient and accurate, the user can provide a routine to
compute the Jacobian; details regarding these application-provided
routines are discussed below. To provide an overview of the use of the
nonlinear solvers, browse the concrete example in <a class="reference external" href="#snes-ex1">ex1.c</a> or skip ahead to the discussion.</p>
<div class="admonition-listing-src-snes-tutorials-ex1-c admonition" id="snes-ex1">
<p class="admonition-title">Listing: <code class="docutils literal notranslate"><span class="pre">src/snes/tutorials/ex1.c</span></code></p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span>
<span class="k">static</span> <span class="kt">char</span> <span class="n">help</span><span class="p">[]</span> <span class="o">=</span> <span class="s">"Newton's method for a two-variable system, sequential.</span><span class="se">\n\n</span><span class="s">"</span><span class="p">;</span>
<span class="cm">/*T</span>
<span class="cm"> Concepts: <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a>^basic example</span>
<span class="cm">T*/</span>
<span class="cm">/*</span>
<span class="cm"> Include "petscsnes.h" so that we can use <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> solvers. Note that this</span>
<span class="cm"> file automatically includes:</span>
<span class="cm"> petscsys.h - base PETSc routines petscvec.h - vectors</span>
<span class="cm"> petscmat.h - matrices</span>
<span class="cm"> petscis.h - index sets petscksp.h - Krylov subspace methods</span>
<span class="cm"> petscviewer.h - viewers petscpc.h - preconditioners</span>
<span class="cm"> petscksp.h - linear solvers</span>
<span class="cm">*/</span>
<span class="cm">/*F</span>
<span class="cm">This examples solves either</span>
<span class="cm">\begin{equation}</span>
<span class="cm"> F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{x^2_0 + x_0 x_1 - 3}{x_0 x_1 + x^2_1 - 6}</span>
<span class="cm">\end{equation}</span>
<span class="cm">or if the {\tt -hard} options is given</span>
<span class="cm">\begin{equation}</span>
<span class="cm"> F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1}</span>
<span class="cm">\end{equation}</span>
<span class="cm">F*/</span>
<span class="cp">#include</span> <span class="cpf"><petscsnes.h></span><span class="cp"></span>
<span class="cm">/*</span>
<span class="cm"> User-defined routines</span>
<span class="cm">*/</span>
<span class="k">extern</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">FormJacobian1</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="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="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="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="k">extern</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">FormFunction1</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="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="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="k">extern</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">FormJacobian2</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="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="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="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="k">extern</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">FormFunction2</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="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="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="kt">int</span> <span class="nf">main</span><span class="p">(</span><span class="kt">int</span> <span class="n">argc</span><span class="p">,</span><span class="kt">char</span> <span class="o">**</span><span class="n">argv</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">;</span> <span class="cm">/* nonlinear solver context */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span> <span class="n">ksp</span><span class="p">;</span> <span class="cm">/* linear solver context */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span> <span class="n">pc</span><span class="p">;</span> <span class="cm">/* preconditioner context */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n">r</span><span class="p">;</span> <span class="cm">/* solution, residual vectors */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="n">J</span><span class="p">;</span> <span class="cm">/* Jacobian matrix */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscMPIInt.html#PetscMPIInt">PetscMPIInt</a></span> <span class="n">size</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">pfive</span> <span class="o">=</span> <span class="mf">.5</span><span class="p">,</span><span class="o">*</span><span class="n">xx</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">flg</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a></span><span class="p">(</span><span class="o">&</span><span class="n">argc</span><span class="p">,</span><span class="o">&</span><span class="n">argv</span><span class="p">,(</span><span class="kt">char</span><span class="o">*</span><span class="p">)</span><span class="mi">0</span><span class="p">,</span><span class="n">help</span><span class="p">);</span><span class="k">if</span> <span class="p">(</span><span class="n">ierr</span><span class="p">)</span> <span class="k">return</span> <span class="n">ierr</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Comm_size.html#MPI_Comm_size">MPI_Comm_size</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">size</span> <span class="o">></span> <span class="mi">1</span><span class="p">)</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/SETERRQ.html#SETERRQ">SETERRQ</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="n">PETSC_ERR_SUP</span><span class="p">,</span><span class="s">"Example is only for sequential runs"</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Create nonlinear solver context</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESCreate.html#SNESCreate">SNESCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="o">&</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Create matrix and vector data structures; set corresponding routines</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="cm">/*</span>
<span class="cm"> Create vectors for solution and nonlinear function</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCreate.html#VecCreate">VecCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetSizes.html#VecSetSizes">VecSetSizes</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="mi">2</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetFromOptions.html#VecSetFromOptions">VecSetFromOptions</a></span><span class="p">(</span><span class="n">x</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">r</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Create Jacobian matrix data structure</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="o">&</span><span class="n">J</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetSizes.html#MatSetSizes">MatSetSizes</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="mi">2</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetFromOptions.html#MatSetFromOptions">MatSetFromOptions</a></span><span class="p">(</span><span class="n">J</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetUp.html#MatSetUp">MatSetUp</a></span><span class="p">(</span><span class="n">J</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-hard"</span><span class="p">,</span><span class="o">&</span><span class="n">flg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="o">!</span><span class="n">flg</span><span class="p">)</span> <span class="p">{</span>
<span class="cm">/*</span>
<span class="cm"> Set function evaluation routine and vector.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html#SNESSetFunction">SNESSetFunction</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">r</span><span class="p">,</span><span class="n">FormFunction1</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set Jacobian matrix data structure and Jacobian evaluation routine</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">FormJacobian1</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span> <span class="k">else</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html#SNESSetFunction">SNESSetFunction</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">r</span><span class="p">,</span><span class="n">FormFunction2</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">FormJacobian2</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Customize nonlinear solver; set runtime options</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="cm">/*</span>
<span class="cm"> Set linear solver defaults for this problem. By extracting the</span>
<span class="cm"> <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a> and <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> contexts from the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context, we can then</span>
<span class="cm"> directly call any <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a> and <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> routines to set various options.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetKSP.html#SNESGetKSP">SNESGetKSP</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">ksp</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGetPC.html#KSPGetPC">KSPGetPC</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="o">&</span><span class="n">pc</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSetType.html#PCSetType">PCSetType</a></span><span class="p">(</span><span class="n">pc</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCNONE.html#PCNONE">PCNONE</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetTolerances.html#KSPSetTolerances">KSPSetTolerances</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="mf">1.</span><span class="n">e</span><span class="mi">-4</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DEFAULT.html#PETSC_DEFAULT">PETSC_DEFAULT</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DEFAULT.html#PETSC_DEFAULT">PETSC_DEFAULT</a></span><span class="p">,</span><span class="mi">20</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> runtime options, e.g.,</span>
<span class="cm"> -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc></span>
<span class="cm"> These options will override those specified above as long as</span>
<span class="cm"> <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFromOptions.html#SNESSetFromOptions">SNESSetFromOptions</a>() is called _after_ any other customization</span>
<span class="cm"> routines.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFromOptions.html#SNESSetFromOptions">SNESSetFromOptions</a></span><span class="p">(</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Evaluate initial guess; then solve nonlinear system</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="k">if</span> <span class="p">(</span><span class="o">!</span><span class="n">flg</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">pfive</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span> <span class="k">else</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/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">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="mf">2.0</span><span class="p">;</span> <span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="mf">3.0</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Note: The user should initialize the vector, x, with the initial guess</span>
<span class="cm"> for the nonlinear solver prior to calling <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSolve.html#SNESSolve">SNESSolve</a>(). In particular,</span>
<span class="cm"> to employ an initial guess of zero, the user should explicitly set</span>
<span class="cm"> this vector to zero by calling <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a>().</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSolve.html#SNESSolve">SNESSolve</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="nb">NULL</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">flg</span><span class="p">)</span> <span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">f</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_WORLD.html#PETSC_VIEWER_STDOUT_WORLD">PETSC_VIEWER_STDOUT_WORLD</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetFunction.html#SNESGetFunction">SNESGetFunction</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">f</span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="mi">0</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</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">r</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_WORLD.html#PETSC_VIEWER_STDOUT_WORLD">PETSC_VIEWER_STDOUT_WORLD</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Free work space. All PETSc objects should be destroyed when they</span>
<span class="cm"> are no longer needed.</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span> <span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">r</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDestroy.html#MatDestroy">MatDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">J</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span> <span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESDestroy.html#SNESDestroy">SNESDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFinalize.html#PetscFinalize">PetscFinalize</a></span><span class="p">();</span>
<span class="k">return</span> <span class="n">ierr</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> FormFunction1 - Evaluates nonlinear function, F(x).</span>
<span class="cm"> Input Parameters:</span>
<span class="cm">. snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="cm">. x - input vector</span>
<span class="cm">. ctx - optional user-defined context</span>
<span class="cm"> Output Parameter:</span>
<span class="cm">. f - function vector</span>
<span class="cm"> */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormFunction1</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">f</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="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">xx</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">ff</span><span class="p">;</span>
<span class="cm">/*</span>
<span class="cm"> Get pointers to vector data.</span>
<span class="cm"> - For default PETSc vectors, <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>() returns a pointer to</span>
<span class="cm"> the data array. Otherwise, the routine is implementation dependent.</span>
<span class="cm"> - You MUST call <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a>() when you no longer need access to</span>
<span class="cm"> the array.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(</span><span class="n">f</span><span class="p">,</span><span class="o">&</span><span class="n">ff</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* Compute function */</span>
<span class="n">ff</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">+</span> <span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">-</span> <span class="mf">3.0</span><span class="p">;</span>
<span class="n">ff</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">+</span> <span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">-</span> <span class="mf">6.0</span><span class="p">;</span>
<span class="cm">/* Restore vectors */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayRead.html#VecRestoreArrayRead">VecRestoreArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n">f</span><span class="p">,</span><span class="o">&</span><span class="n">ff</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">return</span> <span class="mi">0</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> FormJacobian1 - Evaluates Jacobian matrix.</span>
<span class="cm"> Input Parameters:</span>
<span class="cm">. snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="cm">. x - input vector</span>
<span class="cm">. dummy - optional user-defined context (not used here)</span>
<span class="cm"> Output Parameters:</span>
<span class="cm">. jac - Jacobian matrix</span>
<span class="cm">. B - optionally different preconditioning matrix</span>
<span class="cm">. flag - flag indicating matrix structure</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormJacobian1</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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/Mat/Mat.html#Mat">Mat</a></span> <span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="n">B</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">dummy</span><span class="p">)</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">xx</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">A</span><span class="p">[</span><span class="mi">4</span><span class="p">];</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="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</span><span class="p">[</span><span class="mi">2</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="cm">/*</span>
<span class="cm"> Get pointer to vector data</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Compute Jacobian entries and insert into matrix.</span>
<span class="cm"> - Since this is such a small problem, we set all entries for</span>
<span class="cm"> the matrix at once.</span>
<span class="cm"> */</span>
<span class="n">A</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="mf">2.0</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">+</span> <span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">];</span> <span class="n">A</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">];</span>
<span class="n">A</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">];</span> <span class="n">A</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">+</span> <span class="mf">2.0</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">];</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="n">idx</span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="n">idx</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Restore vector</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayRead.html#VecRestoreArrayRead">VecRestoreArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Assemble matrix</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">jac</span> <span class="o">!=</span> <span class="n">B</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="k">return</span> <span class="mi">0</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormFunction2</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">f</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">dummy</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="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">xx</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">ff</span><span class="p">;</span>
<span class="cm">/*</span>
<span class="cm"> Get pointers to vector data.</span>
<span class="cm"> - For default PETSc vectors, <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a>() returns a pointer to</span>
<span class="cm"> the data array. Otherwise, the routine is implementation dependent.</span>
<span class="cm"> - You MUST call <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a>() when you no longer need access to</span>
<span class="cm"> the array.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArray.html#VecGetArray">VecGetArray</a></span><span class="p">(</span><span class="n">f</span><span class="p">,</span><span class="o">&</span><span class="n">ff</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Compute function</span>
<span class="cm"> */</span>
<span class="n">ff</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="n">PetscSinScalar</span><span class="p">(</span><span class="mf">3.0</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">])</span> <span class="o">+</span> <span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">];</span>
<span class="n">ff</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="mi">1</span><span class="p">];</span>
<span class="cm">/*</span>
<span class="cm"> Restore vectors</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayRead.html#VecRestoreArrayRead">VecRestoreArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</a></span><span class="p">(</span><span class="n">f</span><span class="p">,</span><span class="o">&</span><span class="n">ff</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">return</span> <span class="mi">0</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormJacobian2</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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/Mat/Mat.html#Mat">Mat</a></span> <span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="n">B</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">dummy</span><span class="p">)</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">xx</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">A</span><span class="p">[</span><span class="mi">4</span><span class="p">];</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="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</span><span class="p">[</span><span class="mi">2</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="cm">/*</span>
<span class="cm"> Get pointer to vector data</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Compute Jacobian entries and insert into matrix.</span>
<span class="cm"> - Since this is such a small problem, we set all entries for</span>
<span class="cm"> the matrix at once.</span>
<span class="cm"> */</span>
<span class="n">A</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="mf">3.0</span><span class="o">*</span><span class="n">PetscCosScalar</span><span class="p">(</span><span class="mf">3.0</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">])</span> <span class="o">+</span> <span class="mf">1.0</span><span class="p">;</span> <span class="n">A</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="mf">0.0</span><span class="p">;</span>
<span class="n">A</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span> <span class="o">=</span> <span class="mf">0.0</span><span class="p">;</span> <span class="n">A</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span> <span class="o">=</span> <span class="mf">1.0</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="n">idx</span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="n">idx</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Restore vector</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecRestoreArrayRead.html#VecRestoreArrayRead">VecRestoreArrayRead</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Assemble matrix</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a></span><span class="p">(</span><span class="n">B</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">jac</span> <span class="o">!=</span> <span class="n">B</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="k">return</span> <span class="mi">0</span><span class="p">;</span>
<span class="p">}</span>
</pre></div>
</div>
</div>
<p>To create a <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> solver, one must first call <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESCreate.html#SNESCreate">SNESCreate</a>()</span></code> as
follows:</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/SNES/SNESCreate.html#SNESCreate">SNESCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/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/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">snes</span><span class="p">);</span>
</pre></div>
</div>
<p>The user must then set routines for evaluating the residual function <a class="reference internal" href="#equation-fx0">(3)</a> and its associated Jacobian matrix, as
discussed in the following sections.</p>
<p>To choose a nonlinear solution method, the user can either 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/SNES/SNESSetType.html#SNESSetType">SNESSetType</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESType.html#SNESType">SNESType</a></span> <span class="n">method</span><span class="p">);</span>
</pre></div>
</div>
<p>or use the option <code class="docutils literal notranslate"><span class="pre">-snes_type</span> <span class="pre"><method></span></code>, where details regarding the
available methods are presented in <a class="reference internal" href="#sec-nlsolvers"><span class="std std-ref">The Nonlinear Solvers</span></a>. The
application code can take complete control of the linear and nonlinear
techniques used in the Newton-like method 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/SNES/SNESSetFromOptions.html#SNESSetFromOptions">SNESSetFromOptions</a></span><span class="p">(</span><span class="n">snes</span><span class="p">);</span>
</pre></div>
</div>
<p>This routine provides an interface to the PETSc options database, so
that at runtime the user can select a particular nonlinear solver, set
various parameters and customized routines (e.g., specialized line
search variants), prescribe the convergence tolerance, and set
monitoring routines. With this routine the user can also control all
linear solver options in the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span></code>, and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span></code> modules, as discussed
in <a class="reference internal" href="ksp.html#chapter-ksp"><span class="std std-ref">KSP: Linear System Solvers</span></a>.</p>
<p>After having set these routines and options, the user solves the problem
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/SNES/SNESSolve.html#SNESSolve">SNESSolve</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">b</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>where <code class="docutils literal notranslate"><span class="pre">x</span></code> should be initialized to the initial guess before calling and contains the solution on return.
In particular, to employ an initial guess of
zero, the user should explicitly set this vector to zero by calling
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecZeroEntries.html#VecZeroEntries">VecZeroEntries</a>(x)</span></code>. Finally, after solving the nonlinear system (or several
systems), the user should destroy the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> context 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/SNES/SNESDestroy.html#SNESDestroy">SNESDestroy</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">snes</span><span class="p">);</span>
</pre></div>
</div>
<div class="section" id="nonlinear-function-evaluation">
<span id="sec-snesfunction"></span><h3>Nonlinear Function Evaluation<a class="headerlink" href="#nonlinear-function-evaluation" title="Permalink to this headline">¶</a></h3>
<p>When solving a system of nonlinear equations, the user must provide a
a residual function <a class="reference internal" href="#equation-fx0">(3)</a>, which is set 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/SNES/SNESSetFunction.html#SNESSetFunction">SNESSetFunction</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">f</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">FormFunction</span><span class="p">)(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">f</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">),</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>The argument <code class="docutils literal notranslate"><span class="pre">f</span></code> is an optional vector for storing the solution; pass <code class="docutils literal notranslate"><span class="pre">NULL</span></code> to have the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> allocate it for you.
The argument <code class="docutils literal notranslate"><span class="pre">ctx</span></code> is an optional user-defined context, which can
store any private, application-specific data required by the function
evaluation routine; <code class="docutils literal notranslate"><span class="pre">NULL</span></code> should be used if such information is not
needed. In C and C++, a user-defined context is merely a structure in
which various objects can be stashed; in Fortran a user context can be
an integer array that contains both parameters and pointers to PETSc
objects.
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex5.c.html">SNES Tutorial ex5</a>
and
<a class="reference external" href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/tutorials/ex5f.F90.html">SNES Tutorial ex5f</a>
give examples of user-defined application contexts in C and Fortran,
respectively.</p>
</div>
<div class="section" id="jacobian-evaluation">
<span id="sec-snesjacobian"></span><h3>Jacobian Evaluation<a class="headerlink" href="#jacobian-evaluation" title="Permalink to this headline">¶</a></h3>
<p>The user must also specify a routine to form some approximation of the
Jacobian matrix, <code class="docutils literal notranslate"><span class="pre">A</span></code>, at the current iterate, <code class="docutils literal notranslate"><span class="pre">x</span></code>, as is typically
done 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/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">Amat</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">Pmat</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">FormJacobian</span><span class="p">)(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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/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/Mat/Mat.html#Mat">Mat</a></span> <span class="n">B</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">),</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>The arguments of the routine <code class="docutils literal notranslate"><span class="pre">FormJacobian()</span></code> are the current iterate,
<code class="docutils literal notranslate"><span class="pre">x</span></code>; the (approximate) Jacobian matrix, <code class="docutils literal notranslate"><span class="pre">Amat</span></code>; the matrix from
which the preconditioner is constructed, <code class="docutils literal notranslate"><span class="pre">Pmat</span></code> (which is usually the
same as <code class="docutils literal notranslate"><span class="pre">Amat</span></code>); and an optional user-defined Jacobian context,
<code class="docutils literal notranslate"><span class="pre">ctx</span></code>, for application-specific data. Note that the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> solvers
are all data-structure neutral, so the full range of PETSc matrix
formats (including “matrix-free” methods) can be used.
<a class="reference internal" href="mat.html#chapter-matrices"><span class="std std-ref">Matrices</span></a> discusses information regarding
available matrix formats and options, while <a class="reference internal" href="#sec-nlmatrixfree"><span class="std std-ref">Matrix-Free Methods</span></a> focuses on matrix-free methods in
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>. We briefly touch on a few details of matrix usage that are
particularly important for efficient use of the nonlinear solvers.</p>
<p>A common usage paradigm is to assemble the problem Jacobian in the
preconditioner storage <code class="docutils literal notranslate"><span class="pre">B</span></code>, rather than <code class="docutils literal notranslate"><span class="pre">A</span></code>. In the case where they
are identical, as in many simulations, this makes no difference.
However, it allows us to check the analytic Jacobian we construct in
<code class="docutils literal notranslate"><span class="pre">FormJacobian()</span></code> by passing the <code class="docutils literal notranslate"><span class="pre">-snes_mf_operator</span></code> flag. This
causes PETSc to approximate the Jacobian using finite differencing of
the function evaluation (discussed in <a class="reference internal" href="#sec-fdmatrix"><span class="std std-ref">Finite Difference Jacobian Approximations</span></a>),
and the analytic Jacobian becomes merely the preconditioner. Even if the
analytic Jacobian is incorrect, it is likely that the finite difference
approximation will converge, and thus this is an excellent method to
verify the analytic Jacobian. Moreover, if the analytic Jacobian is
incomplete (some terms are missing or approximate),
<code class="docutils literal notranslate"><span class="pre">-snes_mf_operator</span></code> may be used to obtain the exact solution, where
the Jacobian approximation has been transferred to the preconditioner.</p>
<p>One such approximate Jacobian comes from “Picard linearization” which
writes the nonlinear system as</p>
<div class="math">
\[\mathbf{F}(\mathbf{x}) := \mathbf{A}(\mathbf{x}) \mathbf{x} - \mathbf{b} = 0
\]</div>
<p>where <span class="math">\(\mathbf{A}(\mathbf{x})\)</span> usually contains the lower-derivative parts of the
equation. For example, the nonlinear diffusion problem</p>
<div class="math">
\[- \nabla\cdot(\kappa(u) \nabla u) = 0
\]</div>
<p>would be linearized as</p>
<div class="math">
\[A(u) v \simeq -\nabla\cdot(\kappa(u) \nabla v).
\]</div>
<p>Usually this linearization is simpler to implement than Newton and the
linear problems are somewhat easier to solve. In addition to using
<code class="docutils literal notranslate"><span class="pre">-snes_mf_operator</span></code> with this approximation to the Jacobian, the
Picard iterative procedure can be performed by defining <span class="math">\(\mathbf{J}(\mathbf{x})\)</span>
to be <span class="math">\(\mathbf{A}(\mathbf{x})\)</span>. Sometimes this iteration exhibits better global
convergence than Newton linearization.</p>
<p>During successive calls to <code class="docutils literal notranslate"><span class="pre">FormJacobian()</span></code>, the user can either
insert new matrix contexts or reuse old ones, depending on the
application requirements. For many sparse matrix formats, reusing the
old space (and merely changing the matrix elements) is more efficient;
however, if the matrix structure completely changes, creating an
entirely new matrix context may be preferable. Upon subsequent calls to
the <code class="docutils literal notranslate"><span class="pre">FormJacobian()</span></code> routine, the user may wish to reinitialize the
matrix entries to zero by calling <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatZeroEntries.html#MatZeroEntries">MatZeroEntries</a>()</span></code>. See
<a class="reference internal" href="mat.html#sec-othermat"><span class="std std-ref">Other Matrix Operations</span></a> for details on the reuse of the matrix
context.</p>
<p>The directory <code class="docutils literal notranslate"><span class="pre">${PETSC_DIR}/src/snes/tutorials</span></code> provides a variety of
examples.</p>
</div>
</div>
<div class="section" id="the-nonlinear-solvers">
<span id="sec-nlsolvers"></span><h2>The Nonlinear Solvers<a class="headerlink" href="#the-nonlinear-solvers" title="Permalink to this headline">¶</a></h2>
<p>As summarized in Table <a class="reference internal" href="#tab-snesdefaults"><span class="std std-ref">PETSc Nonlinear Solvers</span></a>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> includes
several Newton-like nonlinear solvers based on line search techniques
and trust region methods. Also provided are several nonlinear Krylov
methods, as well as nonlinear methods involving decompositions of the
problem.</p>
<p>Each solver may have associated with it a set of options, which can be
set with routines and options database commands provided for this
purpose. A complete list can be found by consulting the manual pages or
by running a program with the <code class="docutils literal notranslate"><span class="pre">-help</span></code> option; we discuss just a few in
the sections below.</p>
<table class="docutils align-default" id="tab-snesdefaults">
<caption><span class="caption-number">Table 7 </span><span class="caption-text">PETSc Nonlinear Solvers</span><a class="headerlink" href="#tab-snesdefaults" title="Permalink to this table">¶</a></caption>
<colgroup>
<col style="width: 25%" />
<col style="width: 25%" />
<col style="width: 25%" />
<col style="width: 25%" />
</colgroup>
<thead>
<tr class="row-odd"><th class="head"><p>Method</p></th>
<th class="head"><p>SNESType</p></th>
<th class="head"><p>Options Name</p></th>
<th class="head"><p>Default Line Search</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Line Search Newton</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNEWTONLS.html#SNESNEWTONLS">SNESNEWTONLS</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">newtonls</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBT.html#SNESLINESEARCHBT">SNESLINESEARCHBT</a></span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Trust region Newton</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNEWTONTR.html#SNESNEWTONTR">SNESNEWTONTR</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">newtontr</span></code></p></td>
<td><p>—</p></td>
</tr>
<tr class="row-even"><td><p>Nonlinear Richardson</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNRICHARDSON.html#SNESNRICHARDSON">SNESNRICHARDSON</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">nrichardson</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHL2.html#SNESLINESEARCHL2">SNESLINESEARCHL2</a></span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Nonlinear CG</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNCG.html#SNESNCG">SNESNCG</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">ncg</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHCP.html#SNESLINESEARCHCP">SNESLINESEARCHCP</a></span></code></p></td>
</tr>
<tr class="row-even"><td><p>Nonlinear GMRES</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNGMRES.html#SNESNGMRES">SNESNGMRES</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">ngmres</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHL2.html#SNESLINESEARCHL2">SNESLINESEARCHL2</a></span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Quasi-Newton</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESQN.html#SNESQN">SNESQN</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">qn</span></code></p></td>
<td><p>see <a class="reference internal" href="#tab-qndefaults"><span class="std std-ref">PETSc quasi-Newton solvers</span></a></p></td>
</tr>
<tr class="row-even"><td><p>Full Approximation Scheme</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNESFAS/SNESFAS.html#SNESFAS">SNESFAS</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">fas</span></code></p></td>
<td><p>—</p></td>
</tr>
<tr class="row-odd"><td><p>Nonlinear ASM</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNASM.html#SNESNASM">SNESNASM</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">nasm</span></code></p></td>
<td><p>–</p></td>
</tr>
<tr class="row-even"><td><p>ASPIN</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESASPIN.html#SNESASPIN">SNESASPIN</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">aspin</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBT.html#SNESLINESEARCHBT">SNESLINESEARCHBT</a></span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Nonlinear Gauss-Seidel</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNGS.html#SNESNGS">SNESNGS</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">ngs</span></code></p></td>
<td><p>–</p></td>
</tr>
<tr class="row-even"><td><p>Anderson Mixing</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESANDERSON.html#SNESANDERSON">SNESANDERSON</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">anderson</span></code></p></td>
<td><p>–</p></td>
</tr>
<tr class="row-odd"><td><p>Newton with constraints (1)</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESVINEWTONRSLS.html#SNESVINEWTONRSLS">SNESVINEWTONRSLS</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">vinewtonrsls</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBT.html#SNESLINESEARCHBT">SNESLINESEARCHBT</a></span></code></p></td>
</tr>
<tr class="row-even"><td><p>Newton with constraints (2)</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESVINEWTONSSLS.html#SNESVINEWTONSSLS">SNESVINEWTONSSLS</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">vinewtonssls</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBT.html#SNESLINESEARCHBT">SNESLINESEARCHBT</a></span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Multi-stage Smoothers</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESMS.html#SNESMS">SNESMS</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">ms</span></code></p></td>
<td><p>–</p></td>
</tr>
<tr class="row-even"><td><p>Composite</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESCOMPOSITE.html#SNESCOMPOSITE">SNESCOMPOSITE</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">composite</span></code></p></td>
<td><p>–</p></td>
</tr>
<tr class="row-odd"><td><p>Linear solve only</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESKSPONLY.html#SNESKSPONLY">SNESKSPONLY</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">ksponly</span></code></p></td>
<td><p>–</p></td>
</tr>
<tr class="row-even"><td><p>Python Shell</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">SNESPYTHON</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">python</span></code></p></td>
<td><p>–</p></td>
</tr>
<tr class="row-odd"><td><p>Shell (user-defined)</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSHELL.html#SNESSHELL">SNESSHELL</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">shell</span></code></p></td>
<td><p>–</p></td>
</tr>
</tbody>
</table>
<div class="section" id="line-search-newton">
<h3>Line Search Newton<a class="headerlink" href="#line-search-newton" title="Permalink to this headline">¶</a></h3>
<p>The method <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNEWTONLS.html#SNESNEWTONLS">SNESNEWTONLS</a></span></code> (<code class="docutils literal notranslate"><span class="pre">-snes_type</span> <span class="pre">newtonls</span></code>) provides a
line search Newton method for solving systems of nonlinear equations. By
default, this technique employs cubic backtracking
<span id="id1">[<a class="reference internal" href="#id299"><span>DS83</span></a>]</span>. Alternative line search techniques are
listed in Table <a class="reference internal" href="#tab-linesearches"><span class="std std-ref">PETSc Line Search Methods</span></a>.</p>
<table class="docutils align-default" id="tab-linesearches">
<caption><span class="caption-number">Table 8 </span><span class="caption-text">PETSc Line Search Methods</span><a class="headerlink" href="#tab-linesearches" title="Permalink to this table">¶</a></caption>
<colgroup>
<col style="width: 34%" />
<col style="width: 39%" />
<col style="width: 27%" />
</colgroup>
<thead>
<tr class="row-odd"><th class="head"><p><strong>Line Search</strong></p></th>
<th class="head"><p><strong>SNESLineSearchType</strong></p></th>
<th class="head"><p><strong>Options Name</strong></p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Backtracking</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBT.html#SNESLINESEARCHBT">SNESLINESEARCHBT</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">bt</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>(damped) step</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBASIC.html#SNESLINESEARCHBASIC">SNESLINESEARCHBASIC</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">basic</span></code></p></td>
</tr>
<tr class="row-even"><td><p>L2-norm Minimization</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHL2.html#SNESLINESEARCHL2">SNESLINESEARCHL2</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">l2</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Critical point</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHCP.html#SNESLINESEARCHCP">SNESLINESEARCHCP</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">cp</span></code></p></td>
</tr>
<tr class="row-even"><td><p>Shell</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHSHELL.html#SNESLINESEARCHSHELL">SNESLINESEARCHSHELL</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">shell</span></code></p></td>
</tr>
</tbody>
</table>
<p>Every <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> has a line search context of type <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span></code> that
may be retrieved 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/SNES/SNESGetLineSearch.html#SNESGetLineSearch">SNESGetLineSearch</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="o">*</span><span class="n">ls</span><span class="p">);.</span>
</pre></div>
</div>
<p>There are several default options for the line searches. The order of
polynomial approximation may be set with <code class="docutils literal notranslate"><span class="pre">-snes_linesearch_order</span></code> 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/SNES/SNESLineSearchSetOrder.html#SNESLineSearchSetOrder">SNESLineSearchSetOrder</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">ls</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">order</span><span class="p">);</span>
</pre></div>
</div>
<p>for instance, 2 for quadratic or 3 for cubic. Sometimes, it may not be
necessary to monitor the progress of the nonlinear iteration. In this
case, <code class="docutils literal notranslate"><span class="pre">-snes_linesearch_norms</span></code> 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/SNES/SNESLineSearchSetComputeNorms.html#SNESLineSearchSetComputeNorms">SNESLineSearchSetComputeNorms</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">ls</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">norms</span><span class="p">);</span>
</pre></div>
</div>
<p>may be used to turn off function, step, and solution norm computation at
the end of the linesearch.</p>
<p>The default line search for the line search Newton method,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBT.html#SNESLINESEARCHBT">SNESLINESEARCHBT</a></span></code> involves several parameters, which are set to
defaults that are reasonable for many applications. The user can
override the defaults by using the following options:</p>
<ul class="simple">
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_linesearch_alpha</span> <span class="pre"><alpha></span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_linesearch_maxstep</span> <span class="pre"><max></span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_linesearch_minlambda</span> <span class="pre"><tol></span></code></p></li>
</ul>
<p>Besides the backtracking linesearch, there are <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHL2.html#SNESLINESEARCHL2">SNESLINESEARCHL2</a></span></code>,
which uses a polynomial secant minimization of <span class="math">\(||F(x)||_2\)</span>, and
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHCP.html#SNESLINESEARCHCP">SNESLINESEARCHCP</a></span></code>, which minimizes <span class="math">\(F(x) \cdot Y\)</span> where
<span class="math">\(Y\)</span> is the search direction. These are both potentially iterative
line searches, which may be used to find a better-fitted steplength in
the case where a single secant search is not sufficient. The number of
iterations may be set with <code class="docutils literal notranslate"><span class="pre">-snes_linesearch_max_it</span></code>. In addition, the
convergence criteria of the iterative line searches may be set using
function tolerances <code class="docutils literal notranslate"><span class="pre">-snes_linesearch_rtol</span></code> and
<code class="docutils literal notranslate"><span class="pre">-snes_linesearch_atol</span></code>, and steplength tolerance
<code class="docutils literal notranslate"><span class="pre">snes_linesearch_ltol</span></code>.</p>
<p>Custom line search types may either be defined using
<code class="docutils literal notranslate"><span class="pre">SNESLineSearchShell</span></code>, or by creating a custom user line search type
in the model of the preexisting ones and register it 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/SNES/SNESLineSearchRegister.html#SNESLineSearchRegister">SNESLineSearchRegister</a></span><span class="p">(</span><span class="k">const</span> <span class="kt">char</span> <span class="n">sname</span><span class="p">[],</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">function</span><span class="p">)(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span><span class="p">));.</span>
</pre></div>
</div>
</div>
<div class="section" id="trust-region-methods">
<h3>Trust Region Methods<a class="headerlink" href="#trust-region-methods" title="Permalink to this headline">¶</a></h3>
<p>The trust region method in <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> for solving systems of nonlinear
equations, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNEWTONTR.html#SNESNEWTONTR">SNESNEWTONTR</a></span></code> (<code class="docutils literal notranslate"><span class="pre">-snes_type</span> <span class="pre">newtontr</span></code>), is taken from the
MINPACK project <span id="id2">[<a class="reference internal" href="#id291"><span>MoreSGH84</span></a>]</span>. Several parameters can be
set to control the variation of the trust region size during the
solution process. In particular, the user can control the initial trust
region radius, computed by</p>
<div class="math">
\[\Delta = \Delta_0 \| F_0 \|_2,
\]</div>
<p>by setting <span class="math">\(\Delta_0\)</span> via the option <code class="docutils literal notranslate"><span class="pre">-snes_tr_delta0</span> <span class="pre"><delta0></span></code>.</p>
</div>
<div class="section" id="nonlinear-krylov-methods">
<h3>Nonlinear Krylov Methods<a class="headerlink" href="#nonlinear-krylov-methods" title="Permalink to this headline">¶</a></h3>
<p>A number of nonlinear Krylov methods are provided, including Nonlinear
Richardson, conjugate gradient, GMRES, and Anderson Mixing. These
methods are described individually below. They are all instrumental to
PETSc’s nonlinear preconditioning.</p>
<p><strong>Nonlinear Richardson.</strong> The nonlinear Richardson iteration merely
takes the form of a line search-damped fixed-point iteration of the form</p>
<div class="math">
\[\mathbf{x}_{k+1} = \mathbf{x}_k - \lambda \mathbf{F}(\mathbf{x}_k), \;\; k=0,1, \ldots,\]</div>
<p>where the default linesearch is <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHL2.html#SNESLINESEARCHL2">SNESLINESEARCHL2</a></span></code>. This simple solver
is mostly useful as a nonlinear smoother, or to provide line search
stabilization to an inner method.</p>
<p><strong>Nonlinear Conjugate Gradients.</strong> Nonlinear CG is equivalent to linear
CG, but with the steplength determined by line search
(<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHCP.html#SNESLINESEARCHCP">SNESLINESEARCHCP</a></span></code> by default). Five variants (Fletcher-Reed,
Hestenes-Steifel, Polak-Ribiere-Polyak, Dai-Yuan, and Conjugate Descent)
are implemented in PETSc and may be chosen 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/SNES/SNESNCGSetType.html#SNESNCGSetType">SNESNCGSetType</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span> <span class="n">SNESNCGType</span> <span class="n">btype</span><span class="p">);</span>
</pre></div>
</div>
<p><strong>Anderson Mixing and Nonlinear GMRES Methods.</strong> Nonlinear GMRES and
Anderson Mixing methods combine the last <span class="math">\(m\)</span> iterates, plus a new
fixed-point iteration iterate, into a residual-minimizing new iterate.</p>
</div>
<div class="section" id="quasi-newton-methods">
<h3>Quasi-Newton Methods<a class="headerlink" href="#quasi-newton-methods" title="Permalink to this headline">¶</a></h3>
<p>Quasi-Newton methods store iterative rank-one updates to the Jacobian
instead of computing it directly. Three limited-memory quasi-Newton
methods are provided, L-BFGS, which are described in
Table <a class="reference internal" href="#tab-qndefaults"><span class="std std-ref">PETSc quasi-Newton solvers</span></a>. These all are encapsulated under
<code class="docutils literal notranslate"><span class="pre">-snes_type</span> <span class="pre">qn</span></code> and may be changed with <code class="docutils literal notranslate"><span class="pre">snes_qn_type</span></code>. The default
is L-BFGS, which provides symmetric updates to an approximate Jacobian.
This iteration is similar to the line search Newton methods.</p>
<table class="docutils align-default" id="tab-qndefaults">
<caption><span class="caption-number">Table 9 </span><span class="caption-text">PETSc quasi-Newton solvers</span><a class="headerlink" href="#tab-qndefaults" title="Permalink to this table">¶</a></caption>
<colgroup>
<col style="width: 25%" />
<col style="width: 25%" />
<col style="width: 25%" />
<col style="width: 25%" />
</colgroup>
<thead>
<tr class="row-odd"><th class="head"><p>QN Method</p></th>
<th class="head"><p><code class="docutils literal notranslate"><span class="pre">SNESQNType</span></code></p></th>
<th class="head"><p>Options Name</p></th>
<th class="head"><p>Default Line Search</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>L-BFGS</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">SNES_QN_LBFGS</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">lbfgs</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHCP.html#SNESLINESEARCHCP">SNESLINESEARCHCP</a></span></code></p></td>
</tr>
<tr class="row-odd"><td><p>“Good” Broyden</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">SNES_QN_BROYDEN</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">broyden</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHBASIC.html#SNESLINESEARCHBASIC">SNESLINESEARCHBASIC</a></span></code></p></td>
</tr>
<tr class="row-even"><td><p>“Bad” Broyden</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">SNES_QN_BADBROYEN</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">badbroyden</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLINESEARCHL2.html#SNESLINESEARCHL2">SNESLINESEARCHL2</a></span></code></p></td>
</tr>
</tbody>
</table>
<p>One may also control the form of the initial Jacobian approximation 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/SNES/SNESQNSetScaleType.html#SNESQNSetScaleType">SNESQNSetScaleType</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span> <span class="n">SNESQNScaleType</span> <span class="n">stype</span><span class="p">);</span>
</pre></div>
</div>
<p>and the restart type 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/SNES/SNESQNSetRestartType.html#SNESQNSetRestartType">SNESQNSetRestartType</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span> <span class="n">SNESQNRestartType</span> <span class="n">rtype</span><span class="p">);</span>
</pre></div>
</div>
</div>
<div class="section" id="the-full-approximation-scheme">
<h3>The Full Approximation Scheme<a class="headerlink" href="#the-full-approximation-scheme" title="Permalink to this headline">¶</a></h3>
<p>The Full Approximation Scheme is a nonlinear multigrid correction. At
each level, there is a recursive cycle control <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> instance, and
either one or two nonlinear solvers as smoothers (up and down). Problems
set up using the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> <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> interface are automatically
coarsened. FAS differs slightly from <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html#PCMG">PCMG</a></span></code>, in that the hierarchy is
constructed recursively. However, much of the interface is a one-to-one
map. We describe the “get” operations here, and it can be assumed that
each has a corresponding “set” operation. For instance, the number of
levels in the hierarchy may be retrieved 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/SNESFAS/SNESFASGetLevels.html#SNESFASGetLevels">SNESFASGetLevels</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">levels</span><span class="p">);</span>
</pre></div>
</div>
<p>There are four <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNESFAS/SNESFAS.html#SNESFAS">SNESFAS</a></span></code> cycle types, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESFASType.html#SNESFASType">SNES_FAS_MULTIPLICATIVE</a></span></code>,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESFASType.html#SNESFASType">SNES_FAS_ADDITIVE</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESFASType.html#SNESFASType">SNES_FAS_FULL</a></span></code>, and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESFASType.html#SNESFASType">SNES_FAS_KASKADE</a></span></code>. The
type may be set 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/SNESFAS/SNESFASSetType.html#SNESFASSetType">SNESFASSetType</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESFASType.html#SNESFASType">SNESFASType</a></span> <span class="n">fastype</span><span class="p">);.</span>
</pre></div>
</div>
<p>and the cycle type, 1 for V, 2 for W, may be set 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/SNESFAS/SNESFASSetCycles.html#SNESFASSetCycles">SNESFASSetCycles</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">cycles</span><span class="p">);.</span>
</pre></div>
</div>
<p>Much like the interface to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html#PCMG">PCMG</a></span></code> described in <a class="reference internal" href="ksp.html#sec-mg"><span class="std std-ref">Multigrid Preconditioners</span></a>, there are interfaces to recover the
various levels’ cycles and smoothers. The level smoothers may 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/SNESFAS/SNESFASGetSmoother.html#SNESFASGetSmoother">SNESFASGetSmoother</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">level</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">smooth</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNESFAS/SNESFASGetSmootherUp.html#SNESFASGetSmootherUp">SNESFASGetSmootherUp</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">level</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">smooth</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNESFAS/SNESFASGetSmootherDown.html#SNESFASGetSmootherDown">SNESFASGetSmootherDown</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">level</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">smooth</span><span class="p">);</span>
</pre></div>
</div>
<p>and the level cycles 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/SNESFAS/SNESFASGetCycleSNES.html#SNESFASGetCycleSNES">SNESFASGetCycleSNES</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">level</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">lsnes</span><span class="p">);.</span>
</pre></div>
</div>
<p>Also akin to <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html#PCMG">PCMG</a></span></code>, the restriction and prolongation at a level may
be acquired 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/SNESFAS/SNESFASGetInterpolation.html#SNESFASGetInterpolation">SNESFASGetInterpolation</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">level</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="o">*</span><span class="n">mat</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNESFAS/SNESFASGetRestriction.html#SNESFASGetRestriction">SNESFASGetRestriction</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">level</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="o">*</span><span class="n">mat</span><span class="p">);</span>
</pre></div>
</div>
<p>In addition, FAS requires special restriction for solution-like
variables, called injection. This may be set 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/SNESFAS/SNESFASGetInjection.html#SNESFASGetInjection">SNESFASGetInjection</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">level</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="o">*</span><span class="n">mat</span><span class="p">);.</span>
</pre></div>
</div>
<p>The coarse solve context may be acquired 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/SNESFAS/SNESFASGetCoarseSolve.html#SNESFASGetCoarseSolve">SNESFASGetCoarseSolve</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">smooth</span><span class="p">);</span>
</pre></div>
</div>
</div>
<div class="section" id="nonlinear-additive-schwarz">
<h3>Nonlinear Additive Schwarz<a class="headerlink" href="#nonlinear-additive-schwarz" title="Permalink to this headline">¶</a></h3>
<p>Nonlinear Additive Schwarz methods (NASM) take a number of local
nonlinear subproblems, solves them independently in parallel, and
combines those solutions into a new approximate solution.</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/SNES/SNESNASMSetSubdomains.html#SNESNASMSetSubdomains">SNESNASMSetSubdomains</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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/SNES/SNES.html#SNES">SNES</a></span> <span class="n">subsnes</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">iscatter</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">oscatter</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">gscatter</span><span class="p">[]);</span>
</pre></div>
</div>
<p>allows for the user to create these local subdomains. Problems set up
using the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> <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> interface are automatically decomposed. To
begin, the type of subdomain updates to the whole solution are limited
to two types borrowed from <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCASM.html#PCASM">PCASM</a></span></code>: <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCASMType.html#PCASMType">PC_ASM_BASIC</a></span></code>, in which the
overlapping updates added. <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCASMType.html#PCASMType">PC_ASM_RESTRICT</a></span></code> updates in a
nonoverlapping fashion. This may be set 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/SNES/SNESNASMSetType.html#SNESNASMSetType">SNESNASMSetType</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCASMType.html#PCASMType">PCASMType</a></span> <span class="n">type</span><span class="p">);.</span>
</pre></div>
</div>
<p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESASPIN.html#SNESASPIN">SNESASPIN</a></span></code> is a helper <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> type that sets up a nonlinearly
preconditioned Newton’s method using NASM as the preconditioner.</p>
</div>
</div>
<div class="section" id="general-options">
<h2>General Options<a class="headerlink" href="#general-options" title="Permalink to this headline">¶</a></h2>
<p>This section discusses options and routines that apply to all <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>
solvers and problem classes. In particular, we focus on convergence
tests, monitoring routines, and tools for checking derivative
computations.</p>
<div class="section" id="convergence-tests">
<span id="sec-snesconvergence"></span><h3>Convergence Tests<a class="headerlink" href="#convergence-tests" title="Permalink to this headline">¶</a></h3>
<p>Convergence of the nonlinear solvers can be detected in a variety of
ways; the user can even specify a customized test, as discussed below.
Most of the nonlinear solvers use <code class="docutils literal notranslate"><span class="pre">SNESConvergenceTestDefault()</span></code>,
however, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNEWTONTR.html#SNESNEWTONTR">SNESNEWTONTR</a></span></code> uses a method-specific additional convergence
test as well. The convergence tests involves several parameters, which
are set by default to values that should be reasonable for a wide range
of problems. The user can customize the parameters to the problem at
hand by using some of the following routines and options.</p>
<p>One method of convergence testing is to declare convergence when the
norm of the change in the solution between successive iterations is less
than some tolerance, <code class="docutils literal notranslate"><span class="pre">stol</span></code>. Convergence can also be determined based
on the norm of the function. Such a test can use either the absolute
size of the norm, <code class="docutils literal notranslate"><span class="pre">atol</span></code>, or its relative decrease, <code class="docutils literal notranslate"><span class="pre">rtol</span></code>, from an
initial guess. The following routine sets these parameters, which are
used in many of the default <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> convergence tests:</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/SNES/SNESSetTolerances.html#SNESSetTolerances">SNESSetTolerances</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">atol</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">rtol</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">stol</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">its</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">fcts</span><span class="p">);</span>
</pre></div>
</div>
<p>This routine also sets the maximum numbers of allowable nonlinear
iterations, <code class="docutils literal notranslate"><span class="pre">its</span></code>, and function evaluations, <code class="docutils literal notranslate"><span class="pre">fcts</span></code>. The
corresponding options database commands for setting these parameters are:</p>
<ul class="simple">
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_atol</span> <span class="pre"><atol></span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_rtol</span> <span class="pre"><rtol></span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_stol</span> <span class="pre"><stol></span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_max_it</span> <span class="pre"><its></span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-snes_max_funcs</span> <span class="pre"><fcts></span></code></p></li>
</ul>
<p>A related routine is <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetTolerances.html#SNESGetTolerances">SNESGetTolerances</a>()</span></code>.</p>
<p>Convergence tests for trust regions methods often use an additional
parameter that indicates the minimum allowable trust region radius. The
user can set this parameter with the option <code class="docutils literal notranslate"><span class="pre">-snes_trtol</span> <span class="pre"><trtol></span></code> or
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/SNES/SNESSetTrustRegionTolerance.html#SNESSetTrustRegionTolerance">SNESSetTrustRegionTolerance</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">trtol</span><span class="p">);</span>
</pre></div>
</div>
<p>Users can set their own customized convergence tests in <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> by
using the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetConvergenceTest.html#SNESSetConvergenceTest">SNESSetConvergenceTest</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">test</span><span class="p">)(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">it</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">xnorm</span><span class="p">,</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">gnorm</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">f</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESConvergedReason.html#SNESConvergedReason">SNESConvergedReason</a></span> <span class="n">reason</span><span class="p">,</span> <span class="kt">void</span> <span class="o">*</span><span class="n">cctx</span><span class="p">),</span><span class="kt">void</span> <span class="o">*</span><span class="n">cctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">destroy</span><span class="p">)(</span><span class="kt">void</span> <span class="o">*</span><span class="n">cctx</span><span class="p">));</span>
</pre></div>
</div>
<p>The final argument of the convergence test routine, <code class="docutils literal notranslate"><span class="pre">cctx</span></code>, denotes an
optional user-defined context for private data. When solving systems of
nonlinear equations, the arguments <code class="docutils literal notranslate"><span class="pre">xnorm</span></code>, <code class="docutils literal notranslate"><span class="pre">gnorm</span></code>, and <code class="docutils literal notranslate"><span class="pre">f</span></code> are
the current iterate norm, current step norm, and function norm,
respectively. <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESConvergedReason.html#SNESConvergedReason">SNESConvergedReason</a></span></code> should be set positive for
convergence and negative for divergence. See <code class="docutils literal notranslate"><span class="pre">include/petscsnes.h</span></code> for
a list of values for <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESConvergedReason.html#SNESConvergedReason">SNESConvergedReason</a></span></code>.</p>
</div>
<div class="section" id="convergence-monitoring">
<span id="sec-snesmonitor"></span><h3>Convergence Monitoring<a class="headerlink" href="#convergence-monitoring" title="Permalink to this headline">¶</a></h3>
<p>By default the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> solvers run silently without displaying
information about the iterations. The user can initiate monitoring 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/SNES/SNESMonitorSet.html#SNESMonitorSet">SNESMonitorSet</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">mon</span><span class="p">)(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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">its</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">norm</span><span class="p">,</span><span class="kt">void</span><span class="o">*</span> <span class="n">mctx</span><span class="p">),</span><span class="kt">void</span> <span class="o">*</span><span class="n">mctx</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="n">monitordestroy</span><span class="p">)(</span><span class="kt">void</span><span class="o">**</span><span class="p">));</span>
</pre></div>
</div>
<p>The routine, <code class="docutils literal notranslate"><span class="pre">mon</span></code>, indicates a user-defined monitoring routine, where
<code class="docutils literal notranslate"><span class="pre">its</span></code> and <code class="docutils literal notranslate"><span class="pre">mctx</span></code> respectively denote the iteration number and an
optional user-defined context for private data for the monitor routine.
The argument <code class="docutils literal notranslate"><span class="pre">norm</span></code> is the function norm.</p>
<p>The routine set by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESMonitorSet.html#SNESMonitorSet">SNESMonitorSet</a>()</span></code> is called once after every
successful step computation within the nonlinear solver. Hence, the user
can employ this routine for any application-specific computations that
should be done after the solution update. The option <code class="docutils literal notranslate"><span class="pre">-snes_monitor</span></code>
activates the default <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> monitor routine,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESMonitorDefault.html#SNESMonitorDefault">SNESMonitorDefault</a>()</span></code>, while <code class="docutils literal notranslate"><span class="pre">-snes_monitor_lg_residualnorm</span></code> draws
a simple line graph of the residual norm’s convergence.</p>
<p>One can cancel hardwired monitoring routines for <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> at runtime
with <code class="docutils literal notranslate"><span class="pre">-snes_monitor_cancel</span></code>.</p>
<p>As the Newton method converges so that the residual norm is small, say
<span class="math">\(10^{-10}\)</span>, many of the final digits printed with the
<code class="docutils literal notranslate"><span class="pre">-snes_monitor</span></code> option are meaningless. Worse, they are different on
different machines; due to different round-off rules used by, say, the
IBM RS6000 and the Sun SPARC. This makes testing between different
machines difficult. The option <code class="docutils literal notranslate"><span class="pre">-snes_monitor_short</span></code> causes PETSc to
print fewer of the digits of the residual norm as it gets smaller; thus
on most of the machines it will always print the same numbers making
cross-process testing easier.</p>
<p>The routines</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetSolution.html#SNESGetSolution">SNESGetSolution</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetFunction.html#SNESGetFunction">SNESGetFunction</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">r</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">,</span><span class="kt">int</span><span class="p">(</span><span class="o">**</span><span class="n">func</span><span class="p">)(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="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="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">));</span>
</pre></div>
</div>
<p>return the solution vector and function vector from a <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> context.
These routines are useful, for instance, if the convergence test
requires some property of the solution or function other than those
passed with routine arguments.</p>
</div>
<div class="section" id="checking-accuracy-of-derivatives">
<span id="sec-snesderivs"></span><h3>Checking Accuracy of Derivatives<a class="headerlink" href="#checking-accuracy-of-derivatives" title="Permalink to this headline">¶</a></h3>
<p>Since hand-coding routines for Jacobian matrix evaluation can be error
prone, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> provides easy-to-use support for checking these matrices
against finite difference versions. In the simplest form of comparison,
users can employ the option <code class="docutils literal notranslate"><span class="pre">-snes_test_jacobian</span></code> to compare the
matrices at several points. Although not exhaustive, this test will
generally catch obvious problems. One can compare the elements of the
two matrices by using the option <code class="docutils literal notranslate"><span class="pre">-snes_test_jacobian_view</span></code> , which
causes the two matrices to be printed to the screen.</p>
<p>Another means for verifying the correctness of a code for Jacobian
computation is running the problem with either the finite difference or
matrix-free variant, <code class="docutils literal notranslate"><span class="pre">-snes_fd</span></code> or <code class="docutils literal notranslate"><span class="pre">-snes_mf</span></code>; see <a class="reference internal" href="#sec-fdmatrix"><span class="std std-ref">Finite Difference Jacobian Approximations</span></a> or <a class="reference internal" href="#sec-nlmatrixfree"><span class="std std-ref">Matrix-Free Methods</span></a>.
If a
problem converges well with these matrix approximations but not with a
user-provided routine, the problem probably lies with the hand-coded
matrix. See the note in <a class="reference internal" href="#sec-snesjacobian"><span class="std std-ref">Jacobian Evaluation</span></a> about
assembling your Jabobian in the “preconditioner” slot <code class="docutils literal notranslate"><span class="pre">Pmat</span></code>.</p>
<p>The correctness of user provided <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATSHELL.html#MATSHELL">MATSHELL</a></span></code> Jacobians in general can be
checked with <code class="docutils literal notranslate"><span class="pre">MatShellTestMultTranspose()</span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatShellTestMult.html#MatShellTestMult">MatShellTestMult</a>()</span></code>.</p>
<p>The correctness of user provided <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATSHELL.html#MATSHELL">MATSHELL</a></span></code> Jacobians via <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetRHSJacobian.html#TSSetRHSJacobian">TSSetRHSJacobian</a>()</span></code>
can be checked with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSRHSJacobianTestTranspose.html#TSRHSJacobianTestTranspose">TSRHSJacobianTestTranspose</a>()</span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSRHSJacobianTest.html#TSRHSJacobianTest">TSRHSJacobianTest</a>()</span></code>
that check the correction of the matrix-transpose vector product and the
matrix-product. From the command line, these can be checked with</p>
<ul class="simple">
<li><p><code class="docutils literal notranslate"><span class="pre">-ts_rhs_jacobian_test_mult_transpose</span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-mat_shell_test_mult_transpose_view</span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-ts_rhs_jacobian_test_mult</span></code></p></li>
<li><p><code class="docutils literal notranslate"><span class="pre">-mat_shell_test_mult_view</span></code></p></li>
</ul>
</div>
</div>
<div class="section" id="inexact-newton-like-methods">
<h2>Inexact Newton-like Methods<a class="headerlink" href="#inexact-newton-like-methods" title="Permalink to this headline">¶</a></h2>
<p>Since exact solution of the linear Newton systems within <a class="reference internal" href="#equation-newton">(4)</a>
at each iteration can be costly, modifications
are often introduced that significantly reduce these expenses and yet
retain the rapid convergence of Newton’s method. Inexact or truncated
Newton techniques approximately solve the linear systems using an
iterative scheme. In comparison with using direct methods for solving
the Newton systems, iterative methods have the virtue of requiring
little space for matrix storage and potentially saving significant
computational work. Within the class of inexact Newton methods, of
particular interest are Newton-Krylov methods, where the subsidiary
iterative technique for solving the Newton system is chosen from the
class of Krylov subspace projection methods. Note that at runtime the
user can set any of the linear solver options discussed in <a class="reference internal" href="ksp.html#chapter-ksp"><span class="std std-ref">KSP: Linear System Solvers</span></a>,
such as <code class="docutils literal notranslate"><span class="pre">-ksp_type</span> <span class="pre"><ksp_method></span></code> and
<code class="docutils literal notranslate"><span class="pre">-pc_type</span> <span class="pre"><pc_method></span></code>, to set the Krylov subspace and preconditioner
methods.</p>
<p>Two levels of iterations occur for the inexact techniques, where during
each global or outer Newton iteration a sequence of subsidiary inner
iterations of a linear solver is performed. Appropriate control of the
accuracy to which the subsidiary iterative method solves the Newton
system at each global iteration is critical, since these inner
iterations determine the asymptotic convergence rate for inexact Newton
techniques. While the Newton systems must be solved well enough to
retain fast local convergence of the Newton’s iterates, use of excessive
inner iterations, particularly when <span class="math">\(\| \mathbf{x}_k - \mathbf{x}_* \|\)</span> is large,
is neither necessary nor economical. Thus, the number of required inner
iterations typically increases as the Newton process progresses, so that
the truncated iterates approach the true Newton iterates.</p>
<p>A sequence of nonnegative numbers <span class="math">\(\{\eta_k\}\)</span> can be used to
indicate the variable convergence criterion. In this case, when solving
a system of nonlinear equations, the update step of the Newton process
remains unchanged, and direct solution of the linear system is replaced
by iteration on the system until the residuals</p>
<div class="math">
\[\mathbf{r}_k^{(i)} = \mathbf{F}'(\mathbf{x}_k) \Delta \mathbf{x}_k + \mathbf{F}(\mathbf{x}_k)
\]</div>
<p>satisfy</p>
<div class="math">
\[\frac{ \| \mathbf{r}_k^{(i)} \| }{ \| \mathbf{F}(\mathbf{x}_k) \| } \leq \eta_k \leq \eta < 1.
\]</div>
<p>Here <span class="math">\(\mathbf{x}_0\)</span> is an initial approximation of the solution, and
<span class="math">\(\| \cdot \|\)</span> denotes an arbitrary norm in <span class="math">\(\Re^n\)</span> .</p>
<p>By default a constant relative convergence tolerance is used for solving
the subsidiary linear systems within the Newton-like methods of
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>. When solving a system of nonlinear equations, one can instead
employ the techniques of Eisenstat and Walker <span id="id3">[<a class="reference internal" href="#id333"><span>EW96</span></a>]</span>
to compute <span class="math">\(\eta_k\)</span> at each step of the nonlinear solver by using
the option <code class="docutils literal notranslate"><span class="pre">-snes_ksp_ew</span></code> . In addition, by adding one’s own
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span></code> convergence test (see <a class="reference internal" href="ksp.html#sec-convergencetests"><span class="std std-ref">Convergence Tests</span></a>), one can easily create one’s own,
problem-dependent, inner convergence tests.</p>
</div>
<div class="section" id="matrix-free-methods">
<span id="sec-nlmatrixfree"></span><h2>Matrix-Free Methods<a class="headerlink" href="#matrix-free-methods" title="Permalink to this headline">¶</a></h2>
<p>The <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> class fully supports matrix-free methods. The matrices
specified in the Jacobian evaluation routine need not be conventional
matrices; instead, they can point to the data required to implement a
particular matrix-free method. The matrix-free variant is allowed <em>only</em>
when the linear systems are solved by an iterative method in combination
with no preconditioning (<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCNONE.html#PCNONE">PCNONE</a></span></code> or <code class="docutils literal notranslate"><span class="pre">-pc_type</span></code> <code class="docutils literal notranslate"><span class="pre">none</span></code>), a
user-provided preconditioner matrix, or a user-provided preconditioner
shell (<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSHELL.html#PCSHELL">PCSHELL</a></span></code>, discussed in <a class="reference internal" href="ksp.html#sec-pc"><span class="std std-ref">Preconditioners</span></a>); that
is, obviously matrix-free methods cannot be used with a direct solver,
approximate factorization, or other preconditioner which requires access
to explicit matrix entries.</p>
<p>The user can create a matrix-free context for use within <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> 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/SNES/MatCreateSNESMF.html#MatCreateSNESMF">MatCreateSNESMF</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="o">*</span><span class="n">mat</span><span class="p">);</span>
</pre></div>
</div>
<p>This routine creates the data structures needed for the matrix-vector
products that arise within Krylov space iterative
methods <span id="id4">[<a class="reference internal" href="#id301"><span>BS90</span></a>]</span> by employing the matrix type
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATSHELL.html#MATSHELL">MATSHELL</a></span></code>, discussed in <a class="reference internal" href="mat.html#sec-matrixfree"><span class="std std-ref">Matrix-Free Matrices</span></a>.
The default <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>
matrix-free approximations can also be invoked with the command
<code class="docutils literal notranslate"><span class="pre">-snes_mf</span></code>. Or, one can retain the user-provided Jacobian
preconditioner, but replace the user-provided Jacobian matrix with the
default matrix free variant with the option <code class="docutils literal notranslate"><span class="pre">-snes_mf_operator</span></code>.</p>
<p>See also</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateMFFD.html#MatCreateMFFD">MatCreateMFFD</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/Mat/Mat.html#Mat">Mat</a></span> <span class="o">*</span><span class="n">mat</span><span class="p">);</span>
</pre></div>
</div>
<p>for users who need a matrix-free matrix but are not using <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>.</p>
<p>The user can set one parameter to control the Jacobian-vector product
approximation 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/Mat/MatMFFDSetFunctionError.html#MatMFFDSetFunctionError">MatMFFDSetFunctionError</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">mat</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">rerror</span><span class="p">);</span>
</pre></div>
</div>
<p>The parameter <code class="docutils literal notranslate"><span class="pre">rerror</span></code> should be set to the square root of the
relative error in the function evaluations, <span class="math">\(e_{rel}\)</span>; the default
is the square root of machine epsilon (about <span class="math">\(10^{-8}\)</span> in double
precision), which assumes that the functions are evaluated to full
floating-point precision accuracy. This parameter can also be set from
the options database with <code class="docutils literal notranslate"><span class="pre">-snes_mf_err</span> <span class="pre"><err></span></code></p>
<p>In addition, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> provides a way to register new routines to compute
the differencing parameter (<span class="math">\(h\)</span>); see the manual page for
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMFFDSetType.html#MatMFFDSetType">MatMFFDSetType</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/MatMFFDRegister.html#MatMFFDRegister">MatMFFDRegister</a>()</span></code>. We currently provide two
default routines accessible via <code class="docutils literal notranslate"><span class="pre">-snes_mf_type</span> <span class="pre"><default</span> <span class="pre">or</span> <span class="pre">wp></span></code>. For
the default approach there is one “tuning” parameter, set 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/Mat/MatMFFDDSSetUmin.html#MatMFFDDSSetUmin">MatMFFDDSSetUmin</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">mat</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">umin</span><span class="p">);</span>
</pre></div>
</div>
<p>This parameter, <code class="docutils literal notranslate"><span class="pre">umin</span></code> (or <span class="math">\(u_{min}\)</span>), is a bit involved; its
default is <span class="math">\(10^{-6}\)</span> . The Jacobian-vector product is approximated
via the formula</p>
<div class="math">
\[F'(u) a \approx \frac{F(u + h*a) - F(u)}{h}
\]</div>
<p>where <span class="math">\(h\)</span> is computed via</p>
<div class="math">
\[h = e_{\text{rel}} \cdot \begin{cases}
u^{T}a/\lVert a \rVert^2_2 & \text{if $|u^T a| > u_{\min} \lVert a \rVert_{1}$} \\
u_{\min} \operatorname{sign}(u^{T}a) \lVert a \rVert_{1}/\lVert a\rVert^2_2 & \text{otherwise}.
\end{cases}\]</div>
<p>This approach is taken from Brown and Saad
<span id="id5">[<a class="reference internal" href="#id301"><span>BS90</span></a>]</span>. The parameter can also be set from the
options database with <code class="docutils literal notranslate"><span class="pre">-snes_mf_umin</span> <span class="pre"><umin></span></code></p>
<p>The second approach, taken from Walker and Pernice,
<span id="id6">[<a class="reference internal" href="#id240"><span>PW98</span></a>]</span>, computes <span class="math">\(h\)</span> via</p>
<div class="math">
\[\begin{aligned}
h = \frac{\sqrt{1 + ||u||}e_{rel}}{||a||}\end{aligned}\]</div>
<p>This has no tunable parameters, but note that inside the nonlinear solve
for the entire <em>linear</em> iterative process <span class="math">\(u\)</span> does not change
hence <span class="math">\(\sqrt{1 + ||u||}\)</span> need be computed only once. This
information may be set with the options</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMFFDWPSetComputeNormU.html#MatMFFDWPSetComputeNormU">MatMFFDWPSetComputeNormU</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">mat</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="p">);</span>
</pre></div>
</div>
<p>or <code class="docutils literal notranslate"><span class="pre">-mat_mffd_compute_normu</span> <span class="pre"><true</span> <span class="pre">or</span> <span class="pre">false></span></code>. This information is used
to eliminate the redundant computation of these parameters, therefore
reducing the number of collective operations and improving the
efficiency of the application code.</p>
<p>It is also possible to monitor the differencing parameters h that are
computed via 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/Mat/MatMFFDSetHHistory.html#MatMFFDSetHHistory">MatMFFDSetHHistory</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="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="p">,</span><span class="kt">int</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMFFDResetHHistory.html#MatMFFDResetHHistory">MatMFFDResetHHistory</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="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="p">,</span><span class="kt">int</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMFFDGetH.html#MatMFFDGetH">MatMFFDGetH</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="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="p">);</span>
</pre></div>
</div>
<p>We include an explicit example of using matrix-free methods in <a class="reference external" href="#snes-ex3">ex3.c</a>.
Note that by using the option <code class="docutils literal notranslate"><span class="pre">-snes_mf</span></code> one can
easily convert any <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> code to use a matrix-free Newton-Krylov
method without a preconditioner. As shown in this example,
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFromOptions.html#SNESSetFromOptions">SNESSetFromOptions</a>()</span></code> must be called <em>after</em> <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a>()</span></code> to
enable runtime switching between the user-specified Jacobian and the
default <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> matrix-free form.</p>
<div class="admonition-listing-src-snes-tutorials-ex3-c admonition" id="snes-ex3">
<p class="admonition-title">Listing: <code class="docutils literal notranslate"><span class="pre">src/snes/tutorials/ex3.c</span></code></p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="k">static</span> <span class="kt">char</span> <span class="n">help</span><span class="p">[]</span> <span class="o">=</span> <span class="s">"Newton methods to solve u'' + u^{2} = f in parallel.</span><span class="se">\n</span><span class="s">\</span>
<span class="s">This example employs a user-defined monitoring routine and optionally a user-defined</span><span class="se">\n</span><span class="s">\</span>
<span class="s">routine to check candidate iterates produced by line search routines.</span><span class="se">\n</span><span class="s">\</span>
<span class="s">The command line options include:</span><span class="se">\n</span><span class="s">\</span>
<span class="s"> -pre_check_iterates : activate checking of iterates</span><span class="se">\n</span><span class="s">\</span>
<span class="s"> -post_check_iterates : activate checking of iterates</span><span class="se">\n</span><span class="s">\</span>
<span class="s"> -check_tol <tol>: set tolerance for iterate checking</span><span class="se">\n</span><span class="s">\</span>
<span class="s"> -user_precond : activate a (trivial) user-defined preconditioner</span><span class="se">\n\n</span><span class="s">"</span><span class="p">;</span>
<span class="cm">/*T</span>
<span class="cm"> Concepts: <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a>^basic parallel example</span>
<span class="cm"> Concepts: <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a>^setting a user-defined monitoring routine</span>
<span class="cm"> Processors: n</span>
<span class="cm">T*/</span>
<span class="cm">/*</span>
<span class="cm"> Include "petscdm.h" so that we can use data management objects (DMs)</span>
<span class="cm"> Include "petscdmda.h" so that we can use distributed arrays (DMDAs).</span>
<span class="cm"> Include "petscsnes.h" so that we can use <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> solvers. Note that this</span>
<span class="cm"> file automatically includes:</span>
<span class="cm"> petscsys.h - base PETSc routines</span>
<span class="cm"> petscvec.h - vectors</span>
<span class="cm"> petscmat.h - matrices</span>
<span class="cm"> petscis.h - index sets</span>
<span class="cm"> petscksp.h - Krylov subspace methods</span>
<span class="cm"> petscviewer.h - viewers</span>
<span class="cm"> petscpc.h - preconditioners</span>
<span class="cm"> petscksp.h - linear solvers</span>
<span class="cm">*/</span>
<span class="cp">#include</span> <span class="cpf"><petscdm.h></span><span class="cp"></span>
<span class="cp">#include</span> <span class="cpf"><petscdmda.h></span><span class="cp"></span>
<span class="cp">#include</span> <span class="cpf"><petscsnes.h></span><span class="cp"></span>
<span class="cm">/*</span>
<span class="cm"> User-defined routines.</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormJacobian</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="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="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="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormFunction</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="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="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormInitialGuess</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="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">Monitor</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span><span class="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">PreCheck</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span><span class="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="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="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span><span class="o">*</span><span class="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">PostCheck</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span><span class="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="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="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="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span><span class="o">*</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span><span class="o">*</span><span class="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">PostSetSubKSP</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span><span class="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="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="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="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span><span class="o">*</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span><span class="o">*</span><span class="p">,</span><span class="kt">void</span><span class="o">*</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">MatrixFreePreconditioner</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span><span class="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="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="p">);</span>
<span class="cm">/*</span>
<span class="cm"> User-defined application context</span>
<span class="cm">*/</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/DM/DM.html#DM">DM</a></span> <span class="n">da</span><span class="p">;</span> <span class="cm">/* distributed array */</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">F</span><span class="p">;</span> <span class="cm">/* right-hand-side of PDE */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscMPIInt.html#PetscMPIInt">PetscMPIInt</a></span> <span class="n">rank</span><span class="p">;</span> <span class="cm">/* rank of processor */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscMPIInt.html#PetscMPIInt">PetscMPIInt</a></span> <span class="n">size</span><span class="p">;</span> <span class="cm">/* size of communicator */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">h</span><span class="p">;</span> <span class="cm">/* mesh spacing */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">sjerr</span><span class="p">;</span> <span class="cm">/* if or not to test jacobian domain error */</span>
<span class="p">}</span> <span class="n">ApplicationCtx</span><span class="p">;</span>
<span class="cm">/*</span>
<span class="cm"> User-defined context for monitoring</span>
<span class="cm">*/</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/Viewer/PetscViewer.html#PetscViewer">PetscViewer</a></span> <span class="n">viewer</span><span class="p">;</span>
<span class="p">}</span> <span class="n">MonitorCtx</span><span class="p">;</span>
<span class="cm">/*</span>
<span class="cm"> User-defined context for checking candidate iterates that are</span>
<span class="cm"> determined by line search methods</span>
<span class="cm">*/</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/Vec/Vec.html#Vec">Vec</a></span> <span class="n">last_step</span><span class="p">;</span> <span class="cm">/* previous iterate */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">tolerance</span><span class="p">;</span> <span class="cm">/* tolerance for changes between successive iterates */</span>
<span class="n">ApplicationCtx</span> <span class="o">*</span><span class="n">user</span><span class="p">;</span>
<span class="p">}</span> <span class="n">StepCheckCtx</span><span class="p">;</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/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">its0</span><span class="p">;</span> <span class="cm">/* num of prevous outer <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a> iterations */</span>
<span class="p">}</span> <span class="n">SetSubKSPCtx</span><span class="p">;</span>
<span class="kt">int</span> <span class="nf">main</span><span class="p">(</span><span class="kt">int</span> <span class="n">argc</span><span class="p">,</span><span class="kt">char</span> <span class="o">**</span><span class="n">argv</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">;</span> <span class="cm">/* <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">linesearch</span><span class="p">;</span> <span class="cm">/* <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a> context */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="n">J</span><span class="p">;</span> <span class="cm">/* Jacobian matrix */</span>
<span class="n">ApplicationCtx</span> <span class="n">ctx</span><span class="p">;</span> <span class="cm">/* user-defined context */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/Vec.html#Vec">Vec</a></span> <span class="n">x</span><span class="p">,</span><span class="n">r</span><span class="p">,</span><span class="n">U</span><span class="p">,</span><span class="n">F</span><span class="p">;</span> <span class="cm">/* vectors */</span>
<span class="n">MonitorCtx</span> <span class="n">monP</span><span class="p">;</span> <span class="cm">/* monitoring context */</span>
<span class="n">StepCheckCtx</span> <span class="n">checkP</span><span class="p">;</span> <span class="cm">/* step-checking context */</span>
<span class="n">SetSubKSPCtx</span> <span class="n">checkP1</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">pre_check</span><span class="p">,</span><span class="n">post_check</span><span class="p">,</span><span class="n">post_setsubksp</span><span class="p">;</span> <span class="cm">/* flag indicating whether we're checking candidate iterates */</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">xp</span><span class="p">,</span><span class="o">*</span><span class="n">FF</span><span class="p">,</span><span class="o">*</span><span class="n">UU</span><span class="p">,</span><span class="n">none</span> <span class="o">=</span> <span class="mf">-1.0</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">its</span><span class="p">,</span><span class="n">N</span> <span class="o">=</span> <span class="mi">5</span><span class="p">,</span><span class="n">i</span><span class="p">,</span><span class="n">maxit</span><span class="p">,</span><span class="n">maxf</span><span class="p">,</span><span class="n">xs</span><span class="p">,</span><span class="n">xm</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">abstol</span><span class="p">,</span><span class="n">rtol</span><span class="p">,</span><span class="n">stol</span><span class="p">,</span><span class="n">norm</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="n">flg</span><span class="p">,</span><span class="n">viewinitial</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</a></span><span class="p">(</span><span class="o">&</span><span class="n">argc</span><span class="p">,</span><span class="o">&</span><span class="n">argv</span><span class="p">,(</span><span class="kt">char</span><span class="o">*</span><span class="p">)</span><span class="mi">0</span><span class="p">,</span><span class="n">help</span><span class="p">);</span><span class="k">if</span> <span class="p">(</span><span class="n">ierr</span><span class="p">)</span> <span class="k">return</span> <span class="n">ierr</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Comm_rank.html#MPI_Comm_rank">MPI_Comm_rank</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="o">&</span><span class="n">ctx</span><span class="p">.</span><span class="n">rank</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Comm_size.html#MPI_Comm_size">MPI_Comm_size</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="o">&</span><span class="n">ctx</span><span class="p">.</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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-n"</span><span class="p">,</span><span class="o">&</span><span class="n">N</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ctx</span><span class="p">.</span><span class="n">h</span> <span class="o">=</span> <span class="mf">1.0</span><span class="o">/</span><span class="p">(</span><span class="n">N</span><span class="mi">-1</span><span class="p">);</span>
<span class="n">ctx</span><span class="p">.</span><span class="n">sjerr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetBool.html#PetscOptionsGetBool">PetscOptionsGetBool</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-test_jacobian_domain_error"</span><span class="p">,</span><span class="o">&</span><span class="n">ctx</span><span class="p">.</span><span class="n">sjerr</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetBool.html#PetscOptionsGetBool">PetscOptionsGetBool</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-view_initial"</span><span class="p">,</span><span class="o">&</span><span class="n">viewinitial</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Create nonlinear solver context</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESCreate.html#SNESCreate">SNESCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="o">&</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Create vector data structures; set function evaluation routine</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="cm">/*</span>
<span class="cm"> Create distributed array (<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a>) to manage parallel grid and vectors</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/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/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMBoundaryType.html#DMBoundaryType">DM_BOUNDARY_NONE</a></span><span class="p">,</span><span class="n">N</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">ctx</span><span class="p">.</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMSetFromOptions.html#DMSetFromOptions">DMSetFromOptions</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMSetUp.html#DMSetUp">DMSetUp</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Extract global and local vectors from <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a>; then duplicate for remaining</span>
<span class="cm"> vectors that are the same types</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateGlobalVector.html#DMCreateGlobalVector">DMCreateGlobalVector</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</span><span class="n">da</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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">r</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">F</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span> <span class="n">ctx</span><span class="p">.</span><span class="n">F</span> <span class="o">=</span> <span class="n">F</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n">x</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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set function evaluation routine and vector. Whenever the nonlinear</span>
<span class="cm"> solver needs to compute the nonlinear function, it will call this</span>
<span class="cm"> routine.</span>
<span class="cm"> - Note that the final routine argument is the user-defined</span>
<span class="cm"> context that provides application-specific data for the</span>
<span class="cm"> function evaluation routine.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html#SNESSetFunction">SNESSetFunction</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">r</span><span class="p">,</span><span class="n">FormFunction</span><span class="p">,</span><span class="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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Create matrix data structure; set Jacobian evaluation routine</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="o">&</span><span class="n">J</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetSizes.html#MatSetSizes">MatSetSizes</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</a></span><span class="p">,</span><span class="n">N</span><span class="p">,</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetFromOptions.html#MatSetFromOptions">MatSetFromOptions</a></span><span class="p">(</span><span class="n">J</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSeqAIJSetPreallocation.html#MatSeqAIJSetPreallocation">MatSeqAIJSetPreallocation</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="mi">3</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMPIAIJSetPreallocation.html#MatMPIAIJSetPreallocation">MatMPIAIJSetPreallocation</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="mi">3</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="mi">3</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set Jacobian matrix data structure and default Jacobian evaluation</span>
<span class="cm"> routine. Whenever the nonlinear solver needs to compute the</span>
<span class="cm"> Jacobian matrix, it will call this routine.</span>
<span class="cm"> - Note that the final routine argument is the user-defined</span>
<span class="cm"> context that provides application-specific data for the</span>
<span class="cm"> Jacobian evaluation routine.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">FormJacobian</span><span class="p">,</span><span class="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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Optionally allow user-provided preconditioner</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-user_precond"</span><span class="p">,</span><span class="o">&</span><span class="n">flg</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">flg</span><span class="p">)</span> <span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span> <span class="n">ksp</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span> <span class="n">pc</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetKSP.html#SNESGetKSP">SNESGetKSP</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">ksp</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGetPC.html#KSPGetPC">KSPGetPC</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="o">&</span><span class="n">pc</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSetType.html#PCSetType">PCSetType</a></span><span class="p">(</span><span class="n">pc</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSHELL.html#PCSHELL">PCSHELL</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCShellSetApply.html#PCShellSetApply">PCShellSetApply</a></span><span class="p">(</span><span class="n">pc</span><span class="p">,</span><span class="n">MatrixFreePreconditioner</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Customize nonlinear solver; set runtime options</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="cm">/*</span>
<span class="cm"> Set an optional user-defined monitoring routine</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerDrawOpen.html#PetscViewerDrawOpen">PetscViewerDrawOpen</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="mi">400</span><span class="p">,</span><span class="mi">400</span><span class="p">,</span><span class="o">&</span><span class="n">monP</span><span class="p">.</span><span class="n">viewer</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESMonitorSet.html#SNESMonitorSet">SNESMonitorSet</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">Monitor</span><span class="p">,</span><span class="o">&</span><span class="n">monP</span><span class="p">,</span><span class="mi">0</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set names for some vectors to facilitate monitoring (optional)</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSetName.html#PetscObjectSetName">PetscObjectSetName</a></span><span class="p">((</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObject.html#PetscObject">PetscObject</a></span><span class="p">)</span><span class="n">x</span><span class="p">,</span><span class="s">"Approximate Solution"</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObjectSetName.html#PetscObjectSetName">PetscObjectSetName</a></span><span class="p">((</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscObject.html#PetscObject">PetscObject</a></span><span class="p">)</span><span class="n">U</span><span class="p">,</span><span class="s">"Exact Solution"</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a>/<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> runtime options, e.g.,</span>
<span class="cm"> -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc></span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFromOptions.html#SNESSetFromOptions">SNESSetFromOptions</a></span><span class="p">(</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set an optional user-defined routine to check the validity of candidate</span>
<span class="cm"> iterates that are determined by line search methods</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetLineSearch.html#SNESGetLineSearch">SNESGetLineSearch</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span> <span class="o">&</span><span class="n">linesearch</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-post_check_iterates"</span><span class="p">,</span><span class="o">&</span><span class="n">post_check</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">post_check</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"Activating post step checking routine</span><span class="se">\n</span><span class="s">"</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetPostCheck.html#SNESLineSearchSetPostCheck">SNESLineSearchSetPostCheck</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span><span class="n">PostCheck</span><span class="p">,</span><span class="o">&</span><span class="n">checkP</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="p">(</span><span class="n">checkP</span><span class="p">.</span><span class="n">last_step</span><span class="p">));</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">checkP</span><span class="p">.</span><span class="n">tolerance</span> <span class="o">=</span> <span class="mf">1.0</span><span class="p">;</span>
<span class="n">checkP</span><span class="p">.</span><span class="n">user</span> <span class="o">=</span> <span class="o">&</span><span class="n">ctx</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsGetReal.html#PetscOptionsGetReal">PetscOptionsGetReal</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-check_tol"</span><span class="p">,</span><span class="o">&</span><span class="n">checkP</span><span class="p">.</span><span class="n">tolerance</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-post_setsubksp"</span><span class="p">,</span><span class="o">&</span><span class="n">post_setsubksp</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">post_setsubksp</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"Activating post setsubksp</span><span class="se">\n</span><span class="s">"</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetPostCheck.html#SNESLineSearchSetPostCheck">SNESLineSearchSetPostCheck</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span><span class="n">PostSetSubKSP</span><span class="p">,</span><span class="o">&</span><span class="n">checkP1</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</a></span><span class="p">(</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="s">"-pre_check_iterates"</span><span class="p">,</span><span class="o">&</span><span class="n">pre_check</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">pre_check</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"Activating pre step checking routine</span><span class="se">\n</span><span class="s">"</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetPreCheck.html#SNESLineSearchSetPreCheck">SNESLineSearchSetPreCheck</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span><span class="n">PreCheck</span><span class="p">,</span><span class="o">&</span><span class="n">checkP</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Print parameters used for convergence testing (optional) ... just</span>
<span class="cm"> to demonstrate this routine; this information is also printed with</span>
<span class="cm"> the option -snes_view</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetTolerances.html#SNESGetTolerances">SNESGetTolerances</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">abstol</span><span class="p">,</span><span class="o">&</span><span class="n">rtol</span><span class="p">,</span><span class="o">&</span><span class="n">stol</span><span class="p">,</span><span class="o">&</span><span class="n">maxit</span><span class="p">,</span><span class="o">&</span><span class="n">maxf</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D</span><span class="se">\n</span><span class="s">"</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">abstol</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">rtol</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">stol</span><span class="p">,</span><span class="n">maxit</span><span class="p">,</span><span class="n">maxf</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Initialize application:</span>
<span class="cm"> Store right-hand-side of PDE and exact solution</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="cm">/*</span>
<span class="cm"> Get local grid boundaries (for 1-dimensional <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a>):</span>
<span class="cm"> xs, xm - starting grid index, width of local grid (no ghost points)</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetCorners.html#DMDAGetCorners">DMDAGetCorners</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</span><span class="n">da</span><span class="p">,</span><span class="o">&</span><span class="n">xs</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">xm</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Get pointers to vector data</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</span><span class="n">da</span><span class="p">,</span><span class="n">F</span><span class="p">,</span><span class="o">&</span><span class="n">FF</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</span><span class="n">da</span><span class="p">,</span><span class="n">U</span><span class="p">,</span><span class="o">&</span><span class="n">UU</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Compute local vector entries</span>
<span class="cm"> */</span>
<span class="n">xp</span> <span class="o">=</span> <span class="n">ctx</span><span class="p">.</span><span class="n">h</span><span class="o">*</span><span class="n">xs</span><span class="p">;</span>
<span class="k">for</span> <span class="p">(</span><span class="n">i</span><span class="o">=</span><span class="n">xs</span><span class="p">;</span> <span class="n">i</span><span class="o"><</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span><span class="p">;</span> <span class="n">i</span><span class="o">++</span><span class="p">)</span> <span class="p">{</span>
<span class="n">FF</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="mf">6.0</span><span class="o">*</span><span class="n">xp</span> <span class="o">+</span> <span class="n">PetscPowScalar</span><span class="p">(</span><span class="n">xp</span><span class="o">+</span><span class="mf">1.</span><span class="n">e</span><span class="mi">-12</span><span class="p">,</span><span class="mf">6.0</span><span class="p">);</span> <span class="cm">/* +1.e-12 is to prevent 0^6 */</span>
<span class="n">UU</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">xp</span><span class="o">*</span><span class="n">xp</span><span class="o">*</span><span class="n">xp</span><span class="p">;</span>
<span class="n">xp</span> <span class="o">+=</span> <span class="n">ctx</span><span class="p">.</span><span class="n">h</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Restore vectors</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</span><span class="n">da</span><span class="p">,</span><span class="n">F</span><span class="p">,</span><span class="o">&</span><span class="n">FF</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</span><span class="n">ctx</span><span class="p">.</span><span class="n">da</span><span class="p">,</span><span class="n">U</span><span class="p">,</span><span class="o">&</span><span class="n">UU</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">viewinitial</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a></span><span class="p">(</span><span class="n">U</span><span class="p">,</span><span class="mi">0</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</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">F</span><span class="p">,</span><span class="mi">0</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Evaluate initial guess; then solve nonlinear system</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="cm">/*</span>
<span class="cm"> Note: The user should initialize the vector, x, with the initial guess</span>
<span class="cm"> for the nonlinear solver prior to calling <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSolve.html#SNESSolve">SNESSolve</a>(). In particular,</span>
<span class="cm"> to employ an initial guess of zero, the user should explicitly set</span>
<span class="cm"> this vector to zero by calling <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a>().</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n">FormInitialGuess</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSolve.html#SNESSolve">SNESSolve</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="nb">NULL</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetIterationNumber.html#SNESGetIterationNumber">SNESGetIterationNumber</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">its</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"Number of <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> iterations = %D</span><span class="se">\n</span><span class="s">"</span><span class="p">,</span><span class="n">its</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -</span>
<span class="cm"> Check solution and clean up</span>
<span class="cm"> - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */</span>
<span class="cm">/*</span>
<span class="cm"> Check the error</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecAXPY.html#VecAXPY">VecAXPY</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">none</span><span class="p">,</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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecNorm.html#VecNorm">VecNorm</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/NORM_2.html#NORM_2">NORM_2</a></span><span class="p">,</span><span class="o">&</span><span class="n">norm</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"Norm of error %g Iterations %D</span><span class="se">\n</span><span class="s">"</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">norm</span><span class="p">,</span><span class="n">its</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">ctx</span><span class="p">.</span><span class="n">sjerr</span><span class="p">)</span> <span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESType.html#SNESType">SNESType</a></span> <span class="n">snestype</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetType.html#SNESGetType">SNESGetType</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">snestype</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> Type %s</span><span class="se">\n</span><span class="s">"</span><span class="p">,</span><span class="n">snestype</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Free work space. All PETSc objects should be destroyed when they</span>
<span class="cm"> are no longer needed.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerDestroy.html#PetscViewerDestroy">PetscViewerDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">monP</span><span class="p">.</span><span class="n">viewer</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">post_check</span><span class="p">)</span> <span class="p">{</span><span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">checkP</span><span class="p">.</span><span class="n">last_step</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">r</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">F</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatDestroy.html#MatDestroy">MatDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">J</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESDestroy.html#SNESDestroy">SNESDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMDestroy.html#DMDestroy">DMDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">ctx</span><span class="p">.</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/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFinalize.html#PetscFinalize">PetscFinalize</a></span><span class="p">();</span>
<span class="k">return</span> <span class="n">ierr</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> FormInitialGuess - Computes initial guess.</span>
<span class="cm"> Input/Output Parameter:</span>
<span class="cm">. x - the solution vector</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormInitialGuess</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="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="n">pfive</span> <span class="o">=</span> <span class="mf">.50</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSet.html#VecSet">VecSet</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">pfive</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> FormFunction - Evaluates nonlinear function, F(x).</span>
<span class="cm"> Input Parameters:</span>
<span class="cm">. snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="cm">. x - input vector</span>
<span class="cm">. ctx - optional user-defined context, as set by <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html#SNESSetFunction">SNESSetFunction</a>()</span>
<span class="cm"> Output Parameter:</span>
<span class="cm">. f - function vector</span>
<span class="cm"> Note:</span>
<span class="cm"> The user-defined context can contain any application-specific</span>
<span class="cm"> data needed for the function evaluation.</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormFunction</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">f</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">)</span>
<span class="p">{</span>
<span class="n">ApplicationCtx</span> <span class="o">*</span><span class="n">user</span> <span class="o">=</span> <span class="p">(</span><span class="n">ApplicationCtx</span><span class="o">*</span><span class="p">)</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/DM/DM.html#DM">DM</a></span> <span class="n">da</span> <span class="o">=</span> <span class="n">user</span><span class="o">-></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/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="o">*</span><span class="n">xx</span><span class="p">,</span><span class="o">*</span><span class="n">ff</span><span class="p">,</span><span class="o">*</span><span class="n">FF</span><span class="p">,</span><span class="n">d</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">i</span><span class="p">,</span><span class="n">M</span><span class="p">,</span><span class="n">xs</span><span class="p">,</span><span class="n">xm</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">xlocal</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</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">da</span><span class="p">,</span><span class="o">&</span><span class="n">xlocal</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Scatter ghost points to local vector, using the 2-step process</span>
<span class="cm"> <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGlobalToLocalBegin.html#DMGlobalToLocalBegin">DMGlobalToLocalBegin</a>(), <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGlobalToLocalEnd.html#DMGlobalToLocalEnd">DMGlobalToLocalEnd</a>().</span>
<span class="cm"> By placing code between these two statements, computations can</span>
<span class="cm"> be done while messages are in transition.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGlobalToLocalBegin.html#DMGlobalToLocalBegin">DMGlobalToLocalBegin</a></span><span class="p">(</span><span class="n">da</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">xlocal</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMGlobalToLocalEnd.html#DMGlobalToLocalEnd">DMGlobalToLocalEnd</a></span><span class="p">(</span><span class="n">da</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">xlocal</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Get pointers to vector data.</span>
<span class="cm"> - The vector xlocal includes ghost point; the vectors x and f do</span>
<span class="cm"> NOT include ghost points.</span>
<span class="cm"> - Using <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a>() allows accessing the values using global ordering</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">xlocal</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">f</span><span class="p">,</span><span class="o">&</span><span class="n">ff</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">user</span><span class="o">-></span><span class="n">F</span><span class="p">,</span><span class="o">&</span><span class="n">FF</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Get local grid boundaries (for 1-dimensional <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDA.html#DMDA">DMDA</a>):</span>
<span class="cm"> xs, xm - starting grid index, width of local grid (no ghost points)</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetCorners.html#DMDAGetCorners">DMDAGetCorners</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="o">&</span><span class="n">xs</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">xm</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetInfo.html#DMDAGetInfo">DMDAGetInfo</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">M</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Set function values for boundary points; define local interior grid point range:</span>
<span class="cm"> xsi - starting interior grid index</span>
<span class="cm"> xei - ending interior grid index</span>
<span class="cm"> */</span>
<span class="k">if</span> <span class="p">(</span><span class="n">xs</span> <span class="o">==</span> <span class="mi">0</span><span class="p">)</span> <span class="p">{</span> <span class="cm">/* left boundary */</span>
<span class="n">ff</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="mi">0</span><span class="p">];</span>
<span class="n">xs</span><span class="o">++</span><span class="p">;</span><span class="n">xm</span><span class="o">--</span><span class="p">;</span>
<span class="p">}</span>
<span class="k">if</span> <span class="p">(</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span> <span class="o">==</span> <span class="n">M</span><span class="p">)</span> <span class="p">{</span> <span class="cm">/* right boundary */</span>
<span class="n">ff</span><span class="p">[</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span><span class="mi">-1</span><span class="p">]</span> <span class="o">=</span> <span class="n">xx</span><span class="p">[</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span><span class="mi">-1</span><span class="p">]</span> <span class="o">-</span> <span class="mf">1.0</span><span class="p">;</span>
<span class="n">xm</span><span class="o">--</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Compute function over locally owned part of the grid (interior points only)</span>
<span class="cm"> */</span>
<span class="n">d</span> <span class="o">=</span> <span class="mf">1.0</span><span class="o">/</span><span class="p">(</span><span class="n">user</span><span class="o">-></span><span class="n">h</span><span class="o">*</span><span class="n">user</span><span class="o">-></span><span class="n">h</span><span class="p">);</span>
<span class="k">for</span> <span class="p">(</span><span class="n">i</span><span class="o">=</span><span class="n">xs</span><span class="p">;</span> <span class="n">i</span><span class="o"><</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span><span class="p">;</span> <span class="n">i</span><span class="o">++</span><span class="p">)</span> <span class="n">ff</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">d</span><span class="o">*</span><span class="p">(</span><span class="n">xx</span><span class="p">[</span><span class="n">i</span><span class="mi">-1</span><span class="p">]</span> <span class="o">-</span> <span class="mf">2.0</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">+</span> <span class="n">xx</span><span class="p">[</span><span class="n">i</span><span class="o">+</span><span class="mi">1</span><span class="p">])</span> <span class="o">+</span> <span class="n">xx</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">-</span> <span class="n">FF</span><span class="p">[</span><span class="n">i</span><span class="p">];</span>
<span class="cm">/*</span>
<span class="cm"> Restore vectors</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">xlocal</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">f</span><span class="p">,</span><span class="o">&</span><span class="n">ff</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">user</span><span class="o">-></span><span class="n">F</span><span class="p">,</span><span class="o">&</span><span class="n">FF</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMRestoreLocalVector.html#DMRestoreLocalVector">DMRestoreLocalVector</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="o">&</span><span class="n">xlocal</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> FormJacobian - Evaluates Jacobian matrix.</span>
<span class="cm"> Input Parameters:</span>
<span class="cm">. snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="cm">. x - input vector</span>
<span class="cm">. dummy - optional user-defined context (not used here)</span>
<span class="cm"> Output Parameters:</span>
<span class="cm">. jac - Jacobian matrix</span>
<span class="cm">. B - optionally different preconditioning matrix</span>
<span class="cm">. flag - flag indicating matrix structure</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">FormJacobian</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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/Mat/Mat.html#Mat">Mat</a></span> <span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat">Mat</a></span> <span class="n">B</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">)</span>
<span class="p">{</span>
<span class="n">ApplicationCtx</span> <span class="o">*</span><span class="n">user</span> <span class="o">=</span> <span class="p">(</span><span class="n">ApplicationCtx</span><span class="o">*</span><span class="p">)</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/PetscScalar.html#PetscScalar">PetscScalar</a></span> <span class="o">*</span><span class="n">xx</span><span class="p">,</span><span class="n">d</span><span class="p">,</span><span class="n">A</span><span class="p">[</span><span class="mi">3</span><span class="p">];</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">i</span><span class="p">,</span><span class="n">j</span><span class="p">[</span><span class="mi">3</span><span class="p">],</span><span class="n">M</span><span class="p">,</span><span class="n">xs</span><span class="p">,</span><span class="n">xm</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="o">=</span> <span class="n">user</span><span class="o">-></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/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="k">if</span> <span class="p">(</span><span class="n">user</span><span class="o">-></span><span class="n">sjerr</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobianDomainError.html#SNESSetJacobianDomainError">SNESSetJacobianDomainError</a></span><span class="p">(</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Get pointer to vector data</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArrayRead.html#DMDAVecGetArrayRead">DMDAVecGetArrayRead</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetCorners.html#DMDAGetCorners">DMDAGetCorners</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="o">&</span><span class="n">xs</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">xm</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Get range of locally owned matrix</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetInfo.html#DMDAGetInfo">DMDAGetInfo</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">M</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Determine starting and ending local indices for interior grid points.</span>
<span class="cm"> Set Jacobian entries for boundary points.</span>
<span class="cm"> */</span>
<span class="k">if</span> <span class="p">(</span><span class="n">xs</span> <span class="o">==</span> <span class="mi">0</span><span class="p">)</span> <span class="p">{</span> <span class="cm">/* left boundary */</span>
<span class="n">i</span> <span class="o">=</span> <span class="mi">0</span><span class="p">;</span> <span class="n">A</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="mf">1.0</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="o">&</span><span class="n">i</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="o">&</span><span class="n">i</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">xs</span><span class="o">++</span><span class="p">;</span><span class="n">xm</span><span class="o">--</span><span class="p">;</span>
<span class="p">}</span>
<span class="k">if</span> <span class="p">(</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span> <span class="o">==</span> <span class="n">M</span><span class="p">)</span> <span class="p">{</span> <span class="cm">/* right boundary */</span>
<span class="n">i</span> <span class="o">=</span> <span class="n">M</span><span class="mi">-1</span><span class="p">;</span>
<span class="n">A</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="mf">1.0</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="o">&</span><span class="n">i</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="o">&</span><span class="n">i</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">xm</span><span class="o">--</span><span class="p">;</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Interior grid points</span>
<span class="cm"> - Note that in this case we set all elements for a particular</span>
<span class="cm"> row at once.</span>
<span class="cm"> */</span>
<span class="n">d</span> <span class="o">=</span> <span class="mf">1.0</span><span class="o">/</span><span class="p">(</span><span class="n">user</span><span class="o">-></span><span class="n">h</span><span class="o">*</span><span class="n">user</span><span class="o">-></span><span class="n">h</span><span class="p">);</span>
<span class="k">for</span> <span class="p">(</span><span class="n">i</span><span class="o">=</span><span class="n">xs</span><span class="p">;</span> <span class="n">i</span><span class="o"><</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span><span class="p">;</span> <span class="n">i</span><span class="o">++</span><span class="p">)</span> <span class="p">{</span>
<span class="n">j</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span> <span class="o">-</span> <span class="mi">1</span><span class="p">;</span> <span class="n">j</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span><span class="p">;</span> <span class="n">j</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span> <span class="o">+</span> <span class="mi">1</span><span class="p">;</span>
<span class="n">A</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="n">A</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span> <span class="o">=</span> <span class="n">d</span><span class="p">;</span> <span class="n">A</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="mf">-2.0</span><span class="o">*</span><span class="n">d</span> <span class="o">+</span> <span class="mf">2.0</span><span class="o">*</span><span class="n">xx</span><span class="p">[</span><span class="n">i</span><span class="p">];</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="o">&</span><span class="n">i</span><span class="p">,</span><span class="mi">3</span><span class="p">,</span><span class="n">j</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/*</span>
<span class="cm"> Assemble matrix, using the 2-step process:</span>
<span class="cm"> <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a>(), <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a>().</span>
<span class="cm"> By placing code between these two statements, computations can be</span>
<span class="cm"> done while messages are in transition.</span>
<span class="cm"> Also, restore vector.</span>
<span class="cm"> */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArrayRead.html#DMDAVecRestoreArrayRead">DMDAVecRestoreArrayRead</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xx</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</a></span><span class="p">(</span><span class="n">jac</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAssemblyType.html#MatAssemblyType">MAT_FINAL_ASSEMBLY</a></span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> Monitor - Optional user-defined monitoring routine that views the</span>
<span class="cm"> current iterate with an x-window plot. Set by <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESMonitorSet.html#SNESMonitorSet">SNESMonitorSet</a>().</span>
<span class="cm"> Input Parameters:</span>
<span class="cm"> snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="cm"> its - iteration number</span>
<span class="cm"> norm - 2-norm function value (may be estimated)</span>
<span class="cm"> ctx - optional user-defined context for private data for the</span>
<span class="cm"> monitor routine, as set by <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESMonitorSet.html#SNESMonitorSet">SNESMonitorSet</a>()</span>
<span class="cm"> Note:</span>
<span class="cm"> See the manpage for <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Viewer/PetscViewerDrawOpen.html#PetscViewerDrawOpen">PetscViewerDrawOpen</a>() for useful runtime options,</span>
<span class="cm"> such as -nox to deactivate all x-window output.</span>
<span class="cm"> */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">Monitor</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="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">its</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">fnorm</span><span class="p">,</span><span class="kt">void</span> <span class="o">*</span><span class="n">ctx</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n">MonitorCtx</span> <span class="o">*</span><span class="n">monP</span> <span class="o">=</span> <span class="p">(</span><span class="n">MonitorCtx</span><span class="o">*</span><span class="p">)</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/Sys/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"iter = %D,<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> Function norm %g</span><span class="se">\n</span><span class="s">"</span><span class="p">,</span><span class="n">its</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">fnorm</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetSolution.html#SNESGetSolution">SNESGetSolution</a></span><span class="p">(</span><span class="n">snes</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/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecView.html#VecView">VecView</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">monP</span><span class="o">-></span><span class="n">viewer</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> PreCheck - Optional user-defined routine that checks the validity of</span>
<span class="cm"> candidate steps of a line search method. Set by <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetPreCheck.html#SNESLineSearchSetPreCheck">SNESLineSearchSetPreCheck</a>().</span>
<span class="cm"> Input Parameters:</span>
<span class="cm"> snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="cm"> xcurrent - current solution</span>
<span class="cm"> y - search direction and length</span>
<span class="cm"> Output Parameters:</span>
<span class="cm"> y - proposed step (search direction and length) (possibly changed)</span>
<span class="cm"> changed_y - tells if the step has changed or not</span>
<span class="cm"> */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">PreCheck</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">linesearch</span><span class="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">xcurrent</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/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">changed_y</span><span class="p">,</span> <span class="kt">void</span> <span class="o">*</span> <span class="n">ctx</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="o">*</span><span class="n">changed_y</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</a></span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> PostCheck - Optional user-defined routine that checks the validity of</span>
<span class="cm"> candidate steps of a line search method. Set by <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetPostCheck.html#SNESLineSearchSetPostCheck">SNESLineSearchSetPostCheck</a>().</span>
<span class="cm"> Input Parameters:</span>
<span class="cm"> snes - the <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> context</span>
<span class="cm"> ctx - optional user-defined context for private data for the</span>
<span class="cm"> monitor routine, as set by <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetPostCheck.html#SNESLineSearchSetPostCheck">SNESLineSearchSetPostCheck</a>()</span>
<span class="cm"> xcurrent - current solution</span>
<span class="cm"> y - search direction and length</span>
<span class="cm"> x - the new candidate iterate</span>
<span class="cm"> Output Parameters:</span>
<span class="cm"> y - proposed step (search direction and length) (possibly changed)</span>
<span class="cm"> x - current iterate (possibly modified)</span>
<span class="cm"> */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">PostCheck</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">linesearch</span><span class="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">xcurrent</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/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">changed_y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">changed_x</span><span class="p">,</span> <span class="kt">void</span> <span class="o">*</span> <span class="n">ctx</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span> <span class="n">i</span><span class="p">,</span><span class="n">iter</span><span class="p">,</span><span class="n">xs</span><span class="p">,</span><span class="n">xm</span><span class="p">;</span>
<span class="n">StepCheckCtx</span> <span class="o">*</span><span class="n">check</span><span class="p">;</span>
<span class="n">ApplicationCtx</span> <span class="o">*</span><span class="n">user</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">xa</span><span class="p">,</span><span class="o">*</span><span class="n">xa_last</span><span class="p">,</span><span class="n">tmp</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">rdiff</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/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="o">*</span><span class="n">changed_x</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</a></span><span class="p">;</span>
<span class="o">*</span><span class="n">changed_y</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchGetSNES.html#SNESLineSearchGetSNES">SNESLineSearchGetSNES</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="o">&</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">check</span> <span class="o">=</span> <span class="p">(</span><span class="n">StepCheckCtx</span><span class="o">*</span><span class="p">)</span><span class="n">ctx</span><span class="p">;</span>
<span class="n">user</span> <span class="o">=</span> <span class="n">check</span><span class="o">-></span><span class="n">user</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetIterationNumber.html#SNESGetIterationNumber">SNESGetIterationNumber</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">iter</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* iteration 1 indicates we are working on the second iteration */</span>
<span class="k">if</span> <span class="p">(</span><span class="n">iter</span> <span class="o">></span> <span class="mi">0</span><span class="p">)</span> <span class="p">{</span>
<span class="n">da</span> <span class="o">=</span> <span class="n">user</span><span class="o">-></span><span class="n">da</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">"Checking candidate step at iteration %D with tolerance %g</span><span class="se">\n</span><span class="s">"</span><span class="p">,</span><span class="n">iter</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">check</span><span class="o">-></span><span class="n">tolerance</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/* Access local array data */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">check</span><span class="o">-></span><span class="n">last_step</span><span class="p">,</span><span class="o">&</span><span class="n">xa_last</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecGetArray.html#DMDAVecGetArray">DMDAVecGetArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xa</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAGetCorners.html#DMDAGetCorners">DMDAGetCorners</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="o">&</span><span class="n">xs</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">xm</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> If we fail the user-defined check for validity of the candidate iterate,</span>
<span class="cm"> then modify the iterate as we like. (Note that the particular modification</span>
<span class="cm"> below is intended simply to demonstrate how to manipulate this data, not</span>
<span class="cm"> as a meaningful or appropriate choice.)</span>
<span class="cm"> */</span>
<span class="k">for</span> <span class="p">(</span><span class="n">i</span><span class="o">=</span><span class="n">xs</span><span class="p">;</span> <span class="n">i</span><span class="o"><</span><span class="n">xs</span><span class="o">+</span><span class="n">xm</span><span class="p">;</span> <span class="n">i</span><span class="o">++</span><span class="p">)</span> <span class="p">{</span>
<span class="k">if</span> <span class="p">(</span><span class="o">!</span><span class="n">PetscAbsScalar</span><span class="p">(</span><span class="n">xa</span><span class="p">[</span><span class="n">i</span><span class="p">]))</span> <span class="n">rdiff</span> <span class="o">=</span> <span class="mi">2</span><span class="o">*</span><span class="n">check</span><span class="o">-></span><span class="n">tolerance</span><span class="p">;</span>
<span class="k">else</span> <span class="n">rdiff</span> <span class="o">=</span> <span class="n">PetscAbsScalar</span><span class="p">((</span><span class="n">xa</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">-</span> <span class="n">xa_last</span><span class="p">[</span><span class="n">i</span><span class="p">])</span><span class="o">/</span><span class="n">xa</span><span class="p">[</span><span class="n">i</span><span class="p">]);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">rdiff</span> <span class="o">></span> <span class="n">check</span><span class="o">-></span><span class="n">tolerance</span><span class="p">)</span> <span class="p">{</span>
<span class="n">tmp</span> <span class="o">=</span> <span class="n">xa</span><span class="p">[</span><span class="n">i</span><span class="p">];</span>
<span class="n">xa</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="mf">.5</span><span class="o">*</span><span class="p">(</span><span class="n">xa</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">+</span> <span class="n">xa_last</span><span class="p">[</span><span class="n">i</span><span class="p">]);</span>
<span class="o">*</span><span class="n">changed_x</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_TRUE.html#PETSC_TRUE">PETSC_TRUE</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">" Altering entry %D: x=%g, x_last=%g, diff=%g, x_new=%g</span><span class="se">\n</span><span class="s">"</span><span class="p">,</span><span class="n">i</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">PetscAbsScalar</span><span class="p">(</span><span class="n">tmp</span><span class="p">),(</span><span class="kt">double</span><span class="p">)</span><span class="n">PetscAbsScalar</span><span class="p">(</span><span class="n">xa_last</span><span class="p">[</span><span class="n">i</span><span class="p">]),(</span><span class="kt">double</span><span class="p">)</span><span class="n">rdiff</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">PetscAbsScalar</span><span class="p">(</span><span class="n">xa</span><span class="p">[</span><span class="n">i</span><span class="p">]));</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="p">}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">check</span><span class="o">-></span><span class="n">last_step</span><span class="p">,</span><span class="o">&</span><span class="n">xa_last</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMDA/DMDAVecRestoreArray.html#DMDAVecRestoreArray">DMDAVecRestoreArray</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">x</span><span class="p">,</span><span class="o">&</span><span class="n">xa</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCopy.html#VecCopy">VecCopy</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">check</span><span class="o">-></span><span class="n">last_step</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a> is used</span>
<span class="cm"> e.g,</span>
<span class="cm"> mpiexec -n 8 ./ex3 -nox -n 10000 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 4 -sub_ksp_type gmres -sub_ksp_max_it 3 -post_setsubksp -sub_ksp_rtol 1.e-16</span>
<span class="cm"> Set by <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetPostCheck.html#SNESLineSearchSetPostCheck">SNESLineSearchSetPostCheck</a>().</span>
<span class="cm"> Input Parameters:</span>
<span class="cm"> linesearch - the LineSearch context</span>
<span class="cm"> xcurrent - current solution</span>
<span class="cm"> y - search direction and length</span>
<span class="cm"> x - the new candidate iterate</span>
<span class="cm"> Output Parameters:</span>
<span class="cm"> y - proposed step (search direction and length) (possibly changed)</span>
<span class="cm"> x - current iterate (possibly modified)</span>
<span class="cm"> */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">PostSetSubKSP</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearch.html#SNESLineSearch">SNESLineSearch</a></span> <span class="n">linesearch</span><span class="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">xcurrent</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/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">changed_y</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</a></span> <span class="o">*</span><span class="n">changed_x</span><span class="p">,</span> <span class="kt">void</span> <span class="o">*</span> <span class="n">ctx</span><span class="p">)</span>
<span class="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n">SetSubKSPCtx</span> <span class="o">*</span><span class="n">check</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">iter</span><span class="p">,</span><span class="n">its</span><span class="p">,</span><span class="n">sub_its</span><span class="p">,</span><span class="n">maxit</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a></span> <span class="n">ksp</span><span class="p">,</span><span class="n">sub_ksp</span><span class="p">,</span><span class="o">*</span><span class="n">sub_ksps</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span> <span class="n">pc</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">ksp_ratio</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionBeginUser.html#PetscFunctionBeginUser">PetscFunctionBeginUser</a></span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchGetSNES.html#SNESLineSearchGetSNES">SNESLineSearchGetSNES</a></span><span class="p">(</span><span class="n">linesearch</span><span class="p">,</span> <span class="o">&</span><span class="n">snes</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">check</span> <span class="o">=</span> <span class="p">(</span><span class="n">SetSubKSPCtx</span><span class="o">*</span><span class="p">)</span><span class="n">ctx</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetIterationNumber.html#SNESGetIterationNumber">SNESGetIterationNumber</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">iter</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetKSP.html#SNESGetKSP">SNESGetKSP</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="o">&</span><span class="n">ksp</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGetPC.html#KSPGetPC">KSPGetPC</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="o">&</span><span class="n">pc</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCBJacobiGetSubKSP.html#PCBJacobiGetSubKSP">PCBJacobiGetSubKSP</a></span><span class="p">(</span><span class="n">pc</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="o">&</span><span class="n">sub_ksps</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">sub_ksp</span> <span class="o">=</span> <span class="n">sub_ksps</span><span class="p">[</span><span class="mi">0</span><span class="p">];</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGetIterationNumber.html#KSPGetIterationNumber">KSPGetIterationNumber</a></span><span class="p">(</span><span class="n">ksp</span><span class="p">,</span><span class="o">&</span><span class="n">its</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span> <span class="cm">/* outer <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a> iteration number */</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGetIterationNumber.html#KSPGetIterationNumber">KSPGetIterationNumber</a></span><span class="p">(</span><span class="n">sub_ksp</span><span class="p">,</span><span class="o">&</span><span class="n">sub_its</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span> <span class="cm">/* inner <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a> iteration number */</span>
<span class="k">if</span> <span class="p">(</span><span class="n">iter</span><span class="p">)</span> <span class="p">{</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">" ...PostCheck snes iteration %D, ksp_it %D %D, subksp_it %D</span><span class="se">\n</span><span class="s">"</span><span class="p">,</span><span class="n">iter</span><span class="p">,</span><span class="n">check</span><span class="o">-></span><span class="n">its0</span><span class="p">,</span><span class="n">its</span><span class="p">,</span><span class="n">sub_its</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ksp_ratio</span> <span class="o">=</span> <span class="p">((</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span><span class="p">)(</span><span class="n">its</span><span class="p">))</span><span class="o">/</span><span class="n">check</span><span class="o">-></span><span class="n">its0</span><span class="p">;</span>
<span class="n">maxit</span> <span class="o">=</span> <span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</a></span><span class="p">)(</span><span class="n">ksp_ratio</span><span class="o">*</span><span class="n">sub_its</span> <span class="o">+</span> <span class="mf">0.5</span><span class="p">);</span>
<span class="k">if</span> <span class="p">(</span><span class="n">maxit</span> <span class="o"><</span> <span class="mi">2</span><span class="p">)</span> <span class="n">maxit</span> <span class="o">=</span> <span class="mi">2</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetTolerances.html#KSPSetTolerances">KSPSetTolerances</a></span><span class="p">(</span><span class="n">sub_ksp</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DEFAULT.html#PETSC_DEFAULT">PETSC_DEFAULT</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DEFAULT.html#PETSC_DEFAULT">PETSC_DEFAULT</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_DEFAULT.html#PETSC_DEFAULT">PETSC_DEFAULT</a></span><span class="p">,</span><span class="n">maxit</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscPrintf.html#PetscPrintf">PetscPrintf</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="s">" ...ksp_ratio %g, new maxit %D</span><span class="se">\n\n</span><span class="s">"</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">ksp_ratio</span><span class="p">,</span><span class="n">maxit</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="p">}</span>
<span class="n">check</span><span class="o">-></span><span class="n">its0</span> <span class="o">=</span> <span class="n">its</span><span class="p">;</span> <span class="cm">/* save current outer <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSP.html#KSP">KSP</a> iteration number */</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</a></span><span class="p">(</span><span class="mi">0</span><span class="p">);</span>
<span class="p">}</span>
<span class="cm">/* ------------------------------------------------------------------- */</span>
<span class="cm">/*</span>
<span class="cm"> MatrixFreePreconditioner - This routine demonstrates the use of a</span>
<span class="cm"> user-provided preconditioner. This code implements just the null</span>
<span class="cm"> preconditioner, which of course is not recommended for general use.</span>
<span class="cm"> Input Parameters:</span>
<span class="cm">+ pc - preconditioner</span>
<span class="cm">- x - input vector</span>
<span class="cm"> Output Parameter:</span>
<span class="cm">. y - preconditioned vector</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="nf">MatrixFreePreconditioner</span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span> <span class="n">pc</span><span class="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="p">{</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="n">ierr</span><span class="p">;</span>
<span class="n">ierr</span> <span class="o">=</span> <span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecCopy.html#VecCopy">VecCopy</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">y</span><span class="p">);</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</a></span><span class="p">(</span><span class="n">ierr</span><span class="p">);</span>
<span class="k">return</span> <span class="mi">0</span><span class="p">;</span>
<span class="p">}</span>
</pre></div>
</div>
</div>
<p>Table <a class="reference internal" href="#tab-jacobians"><span class="std std-ref">Jacobian Options</span></a> summarizes the various matrix situations
that <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> supports. In particular, different linear system matrices
and preconditioning matrices are allowed, as well as both matrix-free
and application-provided preconditioners. If <a class="reference external" href="#snes-ex3">ex3.c</a> is run with
the options <code class="docutils literal notranslate"><span class="pre">-snes_mf</span></code> and <code class="docutils literal notranslate"><span class="pre">-user_precond</span></code> then it uses a
matrix-free application of the matrix-vector multiple and a user
provided matrix free Jacobian.</p>
<table class="docutils align-default" id="tab-jacobians">
<caption><span class="caption-number">Table 10 </span><span class="caption-text">Jacobian Options</span><a class="headerlink" href="#tab-jacobians" title="Permalink to this table">¶</a></caption>
<colgroup>
<col style="width: 33%" />
<col style="width: 33%" />
<col style="width: 33%" />
</colgroup>
<tbody>
<tr class="row-odd"><td><p>Matrix Use</p></td>
<td><p>Conventional Matrix Formats</p></td>
<td><p>Matrix-free versions</p></td>
</tr>
<tr class="row-even"><td><p>Jacobian Matrix</p></td>
<td><p>Create matrix with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</a>()</span></code><span class="math">\(^*\)</span>. Assemble matrix with user-defined routine <span class="math">\(^\dagger\)</span></p></td>
<td><p>Create matrix with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateShell.html#MatCreateShell">MatCreateShell</a>()</span></code>. Use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatShellSetOperation.html#MatShellSetOperation">MatShellSetOperation</a>()</span></code> to set various matrix actions, or use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateMFFD.html#MatCreateMFFD">MatCreateMFFD</a>()</span></code> or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/MatCreateSNESMF.html#MatCreateSNESMF">MatCreateSNESMF</a>()</span></code>.</p></td>
</tr>
<tr class="row-odd"><td><p>Preconditioning Matrix</p></td>
<td><p>Create matrix with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</a>()</span></code><span class="math">\(^*\)</span>. Assemble matrix with user-defined routine <span class="math">\(^\dagger\)</span></p></td>
<td><p>Use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESGetKSP.html#SNESGetKSP">SNESGetKSP</a>()</span></code> and <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGetPC.html#KSPGetPC">KSPGetPC</a>()</span></code> to access the <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PC.html#PC">PC</a></span></code>, then use <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSetType.html#PCSetType">PCSetType</a>(pc,</span> <span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSHELL.html#PCSHELL">PCSHELL</a>)</span></code> followed by <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCShellSetApply.html#PCShellSetApply">PCShellSetApply</a>()</span></code>.</p></td>
</tr>
</tbody>
</table>
<div class="line-block">
<div class="line"><span class="math">\(^*\)</span> Use either the generic <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</a>()</span></code> or a format-specific variant such as <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateAIJ.html#MatCreateAIJ">MatCreateAIJ</a>()</span></code>.</div>
<div class="line"><span class="math">\(^\dagger\)</span> Set user-defined matrix formation routine with <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a>()</span></code> or with a <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> variant such as <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/DMDASNESSetJacobianLocal.html#DMDASNESSetJacobianLocal">DMDASNESSetJacobianLocal</a>()</span></code></div>
</div>
</div>
<div class="section" id="finite-difference-jacobian-approximations">
<span id="sec-fdmatrix"></span><h2>Finite Difference Jacobian Approximations<a class="headerlink" href="#finite-difference-jacobian-approximations" title="Permalink to this headline">¶</a></h2>
<p>PETSc provides some tools to help approximate the Jacobian matrices
efficiently via finite differences. These tools are intended for use in
certain situations where one is unable to compute Jacobian matrices
analytically, and matrix-free methods do not work well without a
preconditioner, due to very poor conditioning. The approximation
requires several steps:</p>
<ul class="simple">
<li><p>First, one colors the columns of the (not yet built) Jacobian matrix,
so that columns of the same color do not share any common rows.</p></li>
<li><p>Next, one creates a <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatFDColoring.html#MatFDColoring">MatFDColoring</a></span></code> data structure that will be
used later in actually computing the Jacobian.</p></li>
<li><p>Finally, one tells the nonlinear solvers of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> to use the
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESComputeJacobianDefaultColor.html#SNESComputeJacobianDefaultColor">SNESComputeJacobianDefaultColor</a>()</span></code> routine to compute the
Jacobians.</p></li>
</ul>
<p>A code fragment that demonstrates this process is given below.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISColoring.html#ISColoring">ISColoring</a></span> <span class="n">iscoloring</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatFDColoring.html#MatFDColoring">MatFDColoring</a></span> <span class="n">fdcoloring</span><span class="p">;</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatColoring.html#MatColoring">MatColoring</a></span> <span class="n">coloring</span><span class="p">;</span>
<span class="cm">/*</span>
<span class="cm"> This initializes the nonzero structure of the Jacobian. This is artificial</span>
<span class="cm"> because clearly if we had a routine to compute the Jacobian we wouldn't</span>
<span class="cm"> need to use finite differences.</span>
<span class="cm">*/</span>
<span class="n">FormJacobian</span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">x</span><span class="p">,</span> <span class="o">&</span><span class="n">J</span><span class="p">,</span> <span class="o">&</span><span class="n">J</span><span class="p">,</span> <span class="o">&</span><span class="n">user</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Color the matrix, i.e. determine groups of columns that share no common</span>
<span class="cm"> rows. These columns in the Jacobian can all be computed simultaneously.</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MatColoringCreate.html#MatColoringCreate">MatColoringCreate</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span> <span class="o">&</span><span class="n">coloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MatColoringSetType.html#MatColoringSetType">MatColoringSetType</a></span><span class="p">(</span><span class="n">coloring</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGSL.html#MATCOLORINGSL">MATCOLORINGSL</a></span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MatColoringSetFromOptions.html#MatColoringSetFromOptions">MatColoringSetFromOptions</a></span><span class="p">(</span><span class="n">coloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MatColoringApply.html#MatColoringApply">MatColoringApply</a></span><span class="p">(</span><span class="n">coloring</span><span class="p">,</span> <span class="o">&</span><span class="n">iscoloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MatColoringDestroy.html#MatColoringDestroy">MatColoringDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">coloring</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Create the data structure that <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESComputeJacobianDefaultColor.html#SNESComputeJacobianDefaultColor">SNESComputeJacobianDefaultColor</a>() uses</span>
<span class="cm"> to compute the actual Jacobians via finite differences.</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringCreate.html#MatFDColoringCreate">MatFDColoringCreate</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n">iscoloring</span><span class="p">,</span> <span class="o">&</span><span class="n">fdcoloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISColoringDestroy.html#ISColoringDestroy">ISColoringDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">iscoloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringSetFunction.html#MatFDColoringSetFunction">MatFDColoringSetFunction</a></span><span class="p">(</span><span class="n">fdcoloring</span><span class="p">,(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</a></span> <span class="p">(</span><span class="o">*</span><span class="p">)(</span><span class="kt">void</span><span class="p">))</span><span class="n">FormFunction</span><span class="p">,</span> <span class="o">&</span><span class="n">user</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringSetFromOptions.html#MatFDColoringSetFromOptions">MatFDColoringSetFromOptions</a></span><span class="p">(</span><span class="n">fdcoloring</span><span class="p">);</span>
<span class="cm">/*</span>
<span class="cm"> Tell <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a> to use the routine <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESComputeJacobianDefaultColor.html#SNESComputeJacobianDefaultColor">SNESComputeJacobianDefaultColor</a>()</span>
<span class="cm"> to compute Jacobians.</span>
<span class="cm">*/</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</a></span><span class="p">(</span><span class="n">snes</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n">J</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESComputeJacobianDefaultColor.html#SNESComputeJacobianDefaultColor">SNESComputeJacobianDefaultColor</a></span><span class="p">,</span><span class="n">fdcoloring</span><span class="p">);</span>
</pre></div>
</div>
<p>Of course, we are cheating a bit. If we do not have an analytic formula
for computing the Jacobian, then how do we know what its nonzero
structure is so that it may be colored? Determining the structure is
problem dependent, but fortunately, for most structured grid problems
(the class of problems for which PETSc was originally designed) if one
knows the stencil used for the nonlinear function one can usually fairly
easily obtain an estimate of the location of nonzeros in the matrix.
This is harder in the unstructured case, but one typically knows where the nonzero entries are from the mesh topology and distribution of degrees of freedom.
If using <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>) for unstructured meshes, the nonzero locations will be identified in <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateMatrix.html#DMCreateMatrix">DMCreateMatrix</a>()</span></code> and the procedure above can be used.
Most external packages for unstructured meshes have similar functionality.</p>
<p>One need not necessarily use a <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatColoring.html#MatColoring">MatColoring</a></span></code> object to determine a
coloring. For example, if a grid can be colored directly (without using
the associated matrix), then that coloring can be provided to
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringCreate.html#MatFDColoringCreate">MatFDColoringCreate</a>()</span></code>. Note that the user must always preset the
nonzero structure in the matrix regardless of which coloring routine is
used.</p>
<p>PETSc provides the following coloring algorithms, which can be selected using <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MatColoringSetType.html#MatColoringSetType">MatColoringSetType</a>()</span></code> or via the command line argument <code class="docutils literal notranslate"><span class="pre">-mat_coloring_type</span></code>.</p>
<table class="docutils align-default">
<colgroup>
<col style="width: 25%" />
<col style="width: 25%" />
<col style="width: 25%" />
<col style="width: 25%" />
</colgroup>
<thead>
<tr class="row-odd"><th class="head"><p>Algorithm</p></th>
<th class="head"><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatColoringType.html#MatColoringType">MatColoringType</a></span></code></p></th>
<th class="head"><p><code class="docutils literal notranslate"><span class="pre">-mat_coloring_type</span></code></p></th>
<th class="head"><p>Parallel</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>smallest-last <span id="id8">[<a class="reference internal" href="#id291"><span>MoreSGH84</span></a>]</span></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGSL.html#MATCOLORINGSL">MATCOLORINGSL</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">sl</span></code></p></td>
<td><p>No</p></td>
</tr>
<tr class="row-odd"><td><p>largest-first <span id="id9">[<a class="reference internal" href="#id291"><span>MoreSGH84</span></a>]</span></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGLF.html#MATCOLORINGLF">MATCOLORINGLF</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">lf</span></code></p></td>
<td><p>No</p></td>
</tr>
<tr class="row-even"><td><p>incidence-degree <span id="id10">[<a class="reference internal" href="#id291"><span>MoreSGH84</span></a>]</span></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGID.html#MATCOLORINGID">MATCOLORINGID</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">id</span></code></p></td>
<td><p>No</p></td>
</tr>
<tr class="row-odd"><td><p>Jones-Plassmann <span id="id11">[<a class="reference internal" href="#id288"><span>JP93</span></a>]</span></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGJP.html#MATCOLORINGJP">MATCOLORINGJP</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">jp</span></code></p></td>
<td><p>Yes</p></td>
</tr>
<tr class="row-even"><td><p>Greedy</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGGREEDY.html#MATCOLORINGGREEDY">MATCOLORINGGREEDY</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">greedy</span></code></p></td>
<td><p>Yes</p></td>
</tr>
<tr class="row-odd"><td><p>Natural (1 color per column)</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">MATCOLORINGNATURAL</span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">natural</span></code></p></td>
<td><p>Yes</p></td>
</tr>
<tr class="row-even"><td><p>Power (<span class="math">\(A^k\)</span> followed by 1-coloring)</p></td>
<td><p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGPOWER.html#MATCOLORINGPOWER">MATCOLORINGPOWER</a></span></code></p></td>
<td><p><code class="docutils literal notranslate"><span class="pre">power</span></code></p></td>
<td><p>Yes</p></td>
</tr>
</tbody>
</table>
<p>As for the matrix-free computation of Jacobians (<a class="reference internal" href="#sec-nlmatrixfree"><span class="std std-ref">Matrix-Free Methods</span></a>), two parameters affect the accuracy of the
finite difference Jacobian approximation. These are set with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringSetParameters.html#MatFDColoringSetParameters">MatFDColoringSetParameters</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatFDColoring.html#MatFDColoring">MatFDColoring</a></span> <span class="n">fdcoloring</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">rerror</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</a></span> <span class="n">umin</span><span class="p">);</span>
</pre></div>
</div>
<p>The parameter <code class="docutils literal notranslate"><span class="pre">rerror</span></code> is the square root of the relative error in the
function evaluations, <span class="math">\(e_{rel}\)</span>; the default is the square root of
machine epsilon (about <span class="math">\(10^{-8}\)</span> in double precision), which
assumes that the functions are evaluated approximately to floating-point
precision accuracy. The second parameter, <code class="docutils literal notranslate"><span class="pre">umin</span></code>, is a bit more
involved; its default is <span class="math">\(10e^{-6}\)</span> . Column <span class="math">\(i\)</span> of the
Jacobian matrix (denoted by <span class="math">\(F_{:i}\)</span>) is approximated by the
formula</p>
<div class="math">
\[F'_{:i} \approx \frac{F(u + h*dx_{i}) - F(u)}{h}
\]</div>
<p>where <span class="math">\(h\)</span> is computed via:</p>
<div class="math">
\[h = e_{\text{rel}} \cdot \begin{cases}
u_{i} & \text{if $|u_{i}| > u_{\min}$} \\
u_{\min} \cdot \operatorname{sign}(u_{i}) & \text{otherwise}.
\end{cases}\]</div>
<p>for <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATMFFD_DS.html#MATMFFD_DS">MATMFFD_DS</a></span></code> or:</p>
<div class="math">
\[h = e_{\text{rel}} \sqrt(\|u\|)\]</div>
<p>for <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATMFFD_WP.html#MATMFFD_WP">MATMFFD_WP</a></span></code> (default). These parameters may be set from the options
database with</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="o">-</span><span class="n">mat_fd_coloring_err</span> <span class="o"><</span><span class="n">err</span><span class="o">></span>
<span class="o">-</span><span class="n">mat_fd_coloring_umin</span> <span class="o"><</span><span class="n">umin</span><span class="o">></span>
<span class="o">-</span><span class="n">mat_fd_type</span> <span class="o"><</span><span class="n">htype</span><span class="o">></span>
</pre></div>
</div>
<p>Note that <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatColoring.html#MatColoring">MatColoring</a></span></code> type <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGSL.html#MATCOLORINGSL">MATCOLORINGSL</a></span></code>, <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGLF.html#MATCOLORINGLF">MATCOLORINGLF</a></span></code>, and
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGID.html#MATCOLORINGID">MATCOLORINGID</a></span></code> are sequential algorithms. <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGJP.html#MATCOLORINGJP">MATCOLORINGJP</a></span></code> and
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatOrderings/MATCOLORINGGREEDY.html#MATCOLORINGGREEDY">MATCOLORINGGREEDY</a></span></code> are parallel algorithms, although in practice they
may create more colors than the sequential algorithms. If one computes
the coloring <code class="docutils literal notranslate"><span class="pre">iscoloring</span></code> reasonably with a parallel algorithm or by
knowledge of the discretization, the routine <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringCreate.html#MatFDColoringCreate">MatFDColoringCreate</a>()</span></code>
is scalable. An example of this for 2D distributed arrays is given below
that uses the utility routine <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateColoring.html#DMCreateColoring">DMCreateColoring</a>()</span></code>.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateColoring.html#DMCreateColoring">DMCreateColoring</a></span><span class="p">(</span><span class="n">da</span><span class="p">,</span><span class="n">IS_COLORING_GHOSTED</span><span class="p">,</span> <span class="o">&</span><span class="n">iscoloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringCreate.html#MatFDColoringCreate">MatFDColoringCreate</a></span><span class="p">(</span><span class="n">J</span><span class="p">,</span><span class="n">iscoloring</span><span class="p">,</span> <span class="o">&</span><span class="n">fdcoloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringSetFromOptions.html#MatFDColoringSetFromOptions">MatFDColoringSetFromOptions</a></span><span class="p">(</span><span class="n">fdcoloring</span><span class="p">);</span>
<span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/ISColoringDestroy.html#ISColoringDestroy">ISColoringDestroy</a></span><span class="p">(</span> <span class="o">&</span><span class="n">iscoloring</span><span class="p">);</span>
</pre></div>
</div>
<p>Note that the routine <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/MatFD/MatFDColoringCreate.html#MatFDColoringCreate">MatFDColoringCreate</a>()</span></code> currently is only
supported for the AIJ and BAIJ matrix formats.</p>
</div>
<div class="section" id="variational-inequalities">
<span id="sec-vi"></span><h2>Variational Inequalities<a class="headerlink" href="#variational-inequalities" title="Permalink to this headline">¶</a></h2>
<p><code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> can also solve variational inequalities with box constraints.
These are nonlinear algebraic systems with additional inequality
constraints on some or all of the variables:
<span class="math">\(Lu_i \le u_i \le Hu_i\)</span>. Some or all of the lower bounds may be
negative infinity (indicated to PETSc with <code class="docutils literal notranslate"><span class="pre">SNES_VI_NINF</span></code>) and some or
all of the upper bounds may be infinity (indicated by <code class="docutils literal notranslate"><span class="pre">SNES_VI_INF</span></code>).
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/SNES/SNESVISetVariableBounds.html#SNESVISetVariableBounds">SNESVISetVariableBounds</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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">Lu</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">Hu</span><span class="p">);</span>
</pre></div>
</div>
<p>is used to indicate that one is solving a variational inequality. The
option <code class="docutils literal notranslate"><span class="pre">-snes_vi_monitor</span></code> turns on extra monitoring of the active set
associated with the bounds and <code class="docutils literal notranslate"><span class="pre">-snes_vi_type</span></code> allows selecting from
several VI solvers, the default is preferred.</p>
</div>
<div class="section" id="nonlinear-preconditioning">
<h2>Nonlinear Preconditioning<a class="headerlink" href="#nonlinear-preconditioning" title="Permalink to this headline">¶</a></h2>
<p>The mathematical framework of nonlinear preconditioning is explained in detail in <span id="id12">[<a class="reference internal" href="#id1963"><span>BKST15</span></a>]</span>.
Nonlinear preconditioning in PETSc involves the use of an inner <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>
instance to define the step for an outer <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> instance. The inner
instance may be extracted 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/SNES/SNESGetNPC.html#SNESGetNPC">SNESGetNPC</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="n">snes</span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="n">npc</span><span class="p">);</span>
</pre></div>
</div>
<p>and passed run-time options using the <code class="docutils literal notranslate"><span class="pre">-npc_</span></code> prefix. Nonlinear
preconditioning comes in two flavors: left and right. The side may be
changed using <code class="docutils literal notranslate"><span class="pre">-snes_npc_side</span></code> or <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetNPCSide.html#SNESSetNPCSide">SNESSetNPCSide</a>()</span></code>. Left nonlinear
preconditioning redefines the nonlinear function as the action of the
nonlinear preconditioner <span class="math">\(\mathbf{M}\)</span>;</p>
<div class="math">
\[\mathbf{F}_{M}(x) = \mathbf{M}(\mathbf{x},\mathbf{b}) - \mathbf{x}.
\]</div>
<p>Right nonlinear preconditioning redefines the nonlinear function as the
function on the action of the nonlinear preconditioner;</p>
<div class="math">
\[\mathbf{F}(\mathbf{M}(\mathbf{x},\mathbf{b})) = \mathbf{b},
\]</div>
<p>which can be interpreted as putting the preconditioner into “striking
distance” of the solution by outer acceleration.</p>
<p>In addition, basic patterns of solver composition are available with the
<code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESType.html#SNESType">SNESType</a></span></code> <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESCOMPOSITE.html#SNESCOMPOSITE">SNESCOMPOSITE</a></span></code>. This allows for two or more <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code>
instances to be combined additively or multiplicatively. By command
line, a set of <code class="docutils literal notranslate"><span class="pre"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span></code> types may be given by comma separated list
argument to <code class="docutils literal notranslate"><span class="pre">-snes_composite_sneses</span></code>. There are additive
(<code class="docutils literal notranslate"><span class="pre">SNES_COMPOSITE_ADDITIVE</span></code>), additive with optimal damping
(<code class="docutils literal notranslate"><span class="pre">SNES_COMPOSITE_ADDITIVEOPTIMAL</span></code>), and multiplicative
(<code class="docutils literal notranslate"><span class="pre">SNES_COMPOSITE_MULTIPLICATIVE</span></code>) variants which may be set 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/SNES/SNESCompositeSetType.html#SNESCompositeSetType">SNESCompositeSetType</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="p">,</span><span class="n">SNESCompositeType</span><span class="p">);</span>
</pre></div>
</div>
<p>New subsolvers may be added to the composite solver 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/SNES/SNESCompositeAddSNES.html#SNESCompositeAddSNES">SNESCompositeAddSNES</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESType.html#SNESType">SNESType</a></span><span class="p">);</span>
</pre></div>
</div>
<p>and 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/SNES/SNESCompositeGetSNES.html#SNESCompositeGetSNES">SNESCompositeGetSNES</a></span><span class="p">(</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span><span class="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="p">,</span><span class="n"><a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNES.html#SNES">SNES</a></span> <span class="o">*</span><span class="p">);</span>
</pre></div>
</div>
<hr><p id="id13"><dl class="citation">
<dt class="label" id="id301"><span class="brackets">BS90</span><span class="fn-backref">(<a href="#id4">1</a>,<a href="#id5">2</a>)</span></dt>
<dd><p>Peter N. Brown and Youcef Saad. Hybrid Krylov methods for nonlinear systems of equations. <em>SIAM J. Sci. Stat. Comput.</em>, 11:450–481, 1990.</p>
</dd>
<dt class="label" id="id299"><span class="brackets"><a class="fn-backref" href="#id1">DS83</a></span></dt>
<dd><p>J. E. Dennis, Jr. and Robert B. Schnabel. <em>Numerical Methods for Unconstrained Optimization and Nonlinear Equations</em>. Prentice-Hall, Englewood Cliffs, NJ, 1983.</p>
</dd>
<dt class="label" id="id333"><span class="brackets"><a class="fn-backref" href="#id3">EW96</a></span></dt>
<dd><p>S. C. Eisenstat and H. F. Walker. Choosing the forcing terms in an inexact Newton method. <em>SIAM J. Scientific Computing</em>, 17:16–32, 1996.</p>
</dd>
<dt class="label" id="id288"><span class="brackets"><a class="fn-backref" href="#id11">JP93</a></span></dt>
<dd><p>Mark T. Jones and Paul E. Plassmann. A parallel graph coloring heuristic. <em>SIAM J. Sci. Comput.</em>, 14(3):654–669, 1993.</p>
</dd>
<dt class="label" id="id291"><span class="brackets">MoreSGH84</span><span class="fn-backref">(<a href="#id2">1</a>,<a href="#id8">2</a>,<a href="#id9">3</a>,<a href="#id10">4</a>)</span></dt>
<dd><p>Jorge J. Moré, Danny C. Sorenson, Burton S. Garbow, and Kenneth E. Hillstrom. The MINPACK project. In Wayne R. Cowell, editor, <em>Sources and Development of Mathematical Software</em>, 88–111. 1984.</p>
</dd>
<dt class="label" id="id240"><span class="brackets"><a class="fn-backref" href="#id6">PW98</a></span></dt>
<dd><p>M. Pernice and H. F. Walker. NITSOL: a Newton iterative solver for nonlinear systems. <em>SIAM J. Sci. Stat. Comput.</em>, 19:302–318, 1998.</p>
</dd>
</dl>
</p>
<p id="id1263"><dl class="citation">
<dt class="label" id="id1963"><span class="brackets"><a class="fn-backref" href="#id12">BKST15</a></span></dt>
<dd><p>Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu. Composing scalable nonlinear algebraic solvers. <em>SIAM Review</em>, 57(4):535–565, 2015. <span><a class="reference external" href="#"></a></span>http://www.mcs.anl.gov/papers/P2010-0112.pdf. URL: <a class="reference external" href="http://www.mcs.anl.gov/papers/P2010-0112.pdf">http://www.mcs.anl.gov/papers/P2010-0112.pdf</a>, <a class="reference external" href="https://doi.org/10.1137/130936725">doi:10.1137/130936725</a>.</p>
</dd>
</dl>
</p>
</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="ts.html" title="TS: Scalable ODE and DAE Solvers"
>next</a> |</li>
<li class="right" >
<a href="ksp.html" title="KSP: Linear System Solvers"
>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>
|