xref: /petsc/src/snes/tutorials/ex35.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown static const char help[] = "-Laplacian u = b as a nonlinear problem.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown 
5c4762a1bSJed Brown     The linear and nonlinear versions of these should give almost identical results on this problem
6c4762a1bSJed Brown 
7c4762a1bSJed Brown     Richardson
8c4762a1bSJed Brown       Nonlinear:
9c4762a1bSJed Brown         -snes_rtol 1.e-12 -snes_monitor -snes_type nrichardson -snes_linesearch_monitor
10c4762a1bSJed Brown 
11c4762a1bSJed Brown       Linear:
12c4762a1bSJed Brown         -snes_rtol 1.e-12 -snes_monitor -ksp_rtol 1.e-12  -ksp_monitor -ksp_type richardson -pc_type none -ksp_richardson_self_scale -info
13c4762a1bSJed Brown 
14c4762a1bSJed Brown     GMRES
15c4762a1bSJed Brown       Nonlinear:
16c4762a1bSJed Brown        -snes_rtol 1.e-12 -snes_monitor  -snes_type ngmres
17c4762a1bSJed Brown 
18c4762a1bSJed Brown       Linear:
19c4762a1bSJed Brown        -snes_rtol 1.e-12 -snes_monitor  -ksp_type gmres -ksp_monitor -ksp_rtol 1.e-12 -pc_type none
20c4762a1bSJed Brown 
21c4762a1bSJed Brown     CG
22c4762a1bSJed Brown        Nonlinear:
23c4762a1bSJed Brown             -snes_rtol 1.e-12 -snes_monitor  -snes_type ncg -snes_linesearch_monitor
24c4762a1bSJed Brown 
25c4762a1bSJed Brown        Linear:
26c4762a1bSJed Brown              -snes_rtol 1.e-12 -snes_monitor  -ksp_type cg -ksp_monitor -ksp_rtol 1.e-12 -pc_type none
27c4762a1bSJed Brown 
28c4762a1bSJed Brown     Multigrid
29c4762a1bSJed Brown        Linear:
30c4762a1bSJed Brown           1 level:
31c4762a1bSJed Brown             -snes_rtol 1.e-12 -snes_monitor  -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type none -mg_levels_ksp_monitor
32c4762a1bSJed Brown             -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor -ksp_rtol 1.e-12  -ksp_monitor_true_residual
33c4762a1bSJed Brown 
34c4762a1bSJed Brown           n levels:
35c4762a1bSJed Brown             -da_refine n
36c4762a1bSJed Brown 
37c4762a1bSJed Brown        Nonlinear:
38c4762a1bSJed Brown          1 level:
39c4762a1bSJed Brown            -snes_rtol 1.e-12 -snes_monitor  -snes_type fas -fas_levels_snes_monitor
40c4762a1bSJed Brown 
41c4762a1bSJed Brown           n levels:
42c4762a1bSJed Brown             -da_refine n  -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly
43c4762a1bSJed Brown 
44c4762a1bSJed Brown */
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*
47c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
48c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
49c4762a1bSJed Brown */
50c4762a1bSJed Brown #include <petscdm.h>
51c4762a1bSJed Brown #include <petscdmda.h>
52c4762a1bSJed Brown #include <petscsnes.h>
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*
55c4762a1bSJed Brown    User-defined routines
56c4762a1bSJed Brown */
57c4762a1bSJed Brown extern PetscErrorCode FormMatrix(DM, Mat);
58c4762a1bSJed Brown extern PetscErrorCode MyComputeFunction(SNES, Vec, Vec, void *);
59c4762a1bSJed Brown extern PetscErrorCode MyComputeJacobian(SNES, Vec, Mat, Mat, void *);
60c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES, Vec);
61c4762a1bSJed Brown 
62d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
63d71ae5a4SJacob Faibussowitsch {
64c4762a1bSJed Brown   SNES      snes;  /* nonlinear solver */
65c4762a1bSJed Brown   SNES      psnes; /* nonlinear Gauss-Seidel approximate solver */
66c4762a1bSJed Brown   Vec       x, b;  /* solution vector */
67c4762a1bSJed Brown   PetscInt  its;   /* iterations for convergence */
68c4762a1bSJed Brown   DM        da;
69c4762a1bSJed Brown   PetscBool use_ngs_as_npc = PETSC_FALSE; /* use the nonlinear Gauss-Seidel approximate solver */
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72c4762a1bSJed Brown      Initialize program
73c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74c4762a1bSJed Brown 
75327415f7SBarry Smith   PetscFunctionBeginUser;
769566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79c4762a1bSJed Brown      Create nonlinear solver context
80c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
819566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_ngs_as_npc", &use_ngs_as_npc, 0));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   if (use_ngs_as_npc) {
869566063dSJacob Faibussowitsch     PetscCall(SNESGetNPC(snes, &psnes));
879566063dSJacob Faibussowitsch     PetscCall(SNESSetType(psnes, SNESSHELL));
889566063dSJacob Faibussowitsch     PetscCall(SNESShellSetSolve(psnes, NonlinearGS));
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
93c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
949566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
959566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
969566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
979566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
989566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
991baa6e33SBarry Smith   if (use_ngs_as_npc) PetscCall(SNESShellSetContext(psnes, da));
100c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
102c4762a1bSJed Brown      vectors that are the same types
103c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1049566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
1059566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &b));
1069566063dSJacob Faibussowitsch   PetscCall(VecSet(b, 1.0));
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, MyComputeFunction, NULL));
1099566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, MyComputeJacobian, NULL));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
113c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1149566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Solve nonlinear system
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1199566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, b, x));
1209566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
127c4762a1bSJed Brown      are no longer needed.
128c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1319566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1329566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1339566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
134b122ec5aSJacob Faibussowitsch   return 0;
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /* ------------------------------------------------------------------- */
138d71ae5a4SJacob Faibussowitsch PetscErrorCode MyComputeFunction(SNES snes, Vec x, Vec F, void *ctx)
139d71ae5a4SJacob Faibussowitsch {
140c4762a1bSJed Brown   Mat J;
141c4762a1bSJed Brown   DM  dm;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   PetscFunctionBeginUser;
1449566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1459566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &J));
146c4762a1bSJed Brown   if (!J) {
1479566063dSJacob Faibussowitsch     PetscCall(DMSetMatType(dm, MATAIJ));
1489566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &J));
1499566063dSJacob Faibussowitsch     PetscCall(MatSetDM(J, NULL));
1509566063dSJacob Faibussowitsch     PetscCall(FormMatrix(dm, J));
1519566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContext(dm, J));
1529566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContextDestroy(dm, (PetscErrorCode(*)(void **))MatDestroy));
153c4762a1bSJed Brown   }
1549566063dSJacob Faibussowitsch   PetscCall(MatMult(J, x, F));
155*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156c4762a1bSJed Brown }
157c4762a1bSJed Brown 
158d71ae5a4SJacob Faibussowitsch PetscErrorCode MyComputeJacobian(SNES snes, Vec x, Mat J, Mat Jp, void *ctx)
159d71ae5a4SJacob Faibussowitsch {
160c4762a1bSJed Brown   DM dm;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   PetscFunctionBeginUser;
1639566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1649566063dSJacob Faibussowitsch   PetscCall(FormMatrix(dm, Jp));
165*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168d71ae5a4SJacob Faibussowitsch PetscErrorCode FormMatrix(DM da, Mat jac)
169d71ae5a4SJacob Faibussowitsch {
170c4762a1bSJed Brown   PetscInt      i, j, nrows = 0;
171c4762a1bSJed Brown   MatStencil    col[5], row, *rows;
172c4762a1bSJed Brown   PetscScalar   v[5], hx, hy, hxdhy, hydhx;
173c4762a1bSJed Brown   DMDALocalInfo info;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   PetscFunctionBeginUser;
1769566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
177c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(info.mx - 1);
178c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(info.my - 1);
179c4762a1bSJed Brown   hxdhy = hx / hy;
180c4762a1bSJed Brown   hydhx = hy / hx;
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(info.ym * info.xm, &rows));
183c4762a1bSJed Brown   /*
184c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
185c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
186c4762a1bSJed Brown         contiguous chunks of rows across the processors.
187c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
188c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
189c4762a1bSJed Brown         appropriate processor during matrix assembly).
190c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
191c4762a1bSJed Brown       - We can set matrix entries either using either
192c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues(), as discussed above.
193c4762a1bSJed Brown   */
194c4762a1bSJed Brown   for (j = info.ys; j < info.ys + info.ym; j++) {
195c4762a1bSJed Brown     for (i = info.xs; i < info.xs + info.xm; i++) {
1969371c9d4SSatish Balay       row.j = j;
1979371c9d4SSatish Balay       row.i = i;
198c4762a1bSJed Brown       /* boundary points */
199c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info.mx - 1 || j == info.my - 1) {
200c4762a1bSJed Brown         v[0] = 2.0 * (hydhx + hxdhy);
2019566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES));
202c4762a1bSJed Brown         rows[nrows].i   = i;
203c4762a1bSJed Brown         rows[nrows++].j = j;
204c4762a1bSJed Brown       } else {
205c4762a1bSJed Brown         /* interior grid points */
2069371c9d4SSatish Balay         v[0]     = -hxdhy;
2079371c9d4SSatish Balay         col[0].j = j - 1;
2089371c9d4SSatish Balay         col[0].i = i;
2099371c9d4SSatish Balay         v[1]     = -hydhx;
2109371c9d4SSatish Balay         col[1].j = j;
2119371c9d4SSatish Balay         col[1].i = i - 1;
2129371c9d4SSatish Balay         v[2]     = 2.0 * (hydhx + hxdhy);
2139371c9d4SSatish Balay         col[2].j = row.j;
2149371c9d4SSatish Balay         col[2].i = row.i;
2159371c9d4SSatish Balay         v[3]     = -hydhx;
2169371c9d4SSatish Balay         col[3].j = j;
2179371c9d4SSatish Balay         col[3].i = i + 1;
2189371c9d4SSatish Balay         v[4]     = -hxdhy;
2199371c9d4SSatish Balay         col[4].j = j + 1;
2209371c9d4SSatish Balay         col[4].i = i;
2219566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jac, 1, &row, 5, col, v, INSERT_VALUES));
222c4762a1bSJed Brown       }
223c4762a1bSJed Brown     }
224c4762a1bSJed Brown   }
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /*
227c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
228c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
229c4762a1bSJed Brown   */
2309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
2329566063dSJacob Faibussowitsch   PetscCall(MatZeroRowsColumnsStencil(jac, nrows, rows, 2.0 * (hydhx + hxdhy), NULL, NULL));
2339566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
234c4762a1bSJed Brown   /*
235c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
236c4762a1bSJed Brown      matrix. If we do, it will generate an error.
237c4762a1bSJed Brown   */
2389566063dSJacob Faibussowitsch   PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
239*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240c4762a1bSJed Brown }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown /* ------------------------------------------------------------------- */
243c4762a1bSJed Brown /*
244c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
245c4762a1bSJed Brown 
246c4762a1bSJed Brown  */
247d71ae5a4SJacob Faibussowitsch PetscErrorCode NonlinearGS(SNES snes, Vec X)
248d71ae5a4SJacob Faibussowitsch {
249c4762a1bSJed Brown   PetscInt      i, j, Mx, My, xs, ys, xm, ym, its, l;
250c4762a1bSJed Brown   PetscReal     hx, hy, hxdhy, hydhx;
251c4762a1bSJed Brown   PetscScalar **x, F, J, u, uxx, uyy;
252c4762a1bSJed Brown   DM            da;
253c4762a1bSJed Brown   Vec           localX;
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   PetscFunctionBeginUser;
2569566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, NULL, NULL, NULL, &its, NULL));
2579566063dSJacob Faibussowitsch   PetscCall(SNESShellGetContext(snes, &da));
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(Mx - 1);
262c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(My - 1);
263c4762a1bSJed Brown   hxdhy = hx / hy;
264c4762a1bSJed Brown   hydhx = hy / hx;
265c4762a1bSJed Brown 
2669566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   for (l = 0; l < its; l++) {
2699566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
2709566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
271c4762a1bSJed Brown     /*
272c4762a1bSJed Brown      Get a pointer to vector data.
273c4762a1bSJed Brown      - For default PETSc vectors, VecGetArray() returns a pointer to
274c4762a1bSJed Brown      the data array.  Otherwise, the routine is implementation dependent.
275c4762a1bSJed Brown      - You MUST call VecRestoreArray() when you no longer need access to
276c4762a1bSJed Brown      the array.
277c4762a1bSJed Brown      */
2789566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, localX, &x));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown     /*
281c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
282c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
283c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
284c4762a1bSJed Brown 
285c4762a1bSJed Brown      */
2869566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
287c4762a1bSJed Brown 
288c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
289c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
290c4762a1bSJed Brown         if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
291c4762a1bSJed Brown           /* boundary conditions are all zero Dirichlet */
292c4762a1bSJed Brown           x[j][i] = 0.0;
293c4762a1bSJed Brown         } else {
294c4762a1bSJed Brown           u   = x[j][i];
295c4762a1bSJed Brown           uxx = (2.0 * u - x[j][i - 1] - x[j][i + 1]) * hydhx;
296c4762a1bSJed Brown           uyy = (2.0 * u - x[j - 1][i] - x[j + 1][i]) * hxdhy;
297c4762a1bSJed Brown           F   = uxx + uyy;
298c4762a1bSJed Brown           J   = 2.0 * (hydhx + hxdhy);
299c4762a1bSJed Brown           u   = u - F / J;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown           x[j][i] = u;
302c4762a1bSJed Brown         }
303c4762a1bSJed Brown       }
304c4762a1bSJed Brown     }
305c4762a1bSJed Brown 
306c4762a1bSJed Brown     /*
307c4762a1bSJed Brown      Restore vector
308c4762a1bSJed Brown      */
3099566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, localX, &x));
3109566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X));
3119566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X));
312c4762a1bSJed Brown   }
3139566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
314*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
315c4762a1bSJed Brown }
316c4762a1bSJed Brown 
317c4762a1bSJed Brown /*TEST
318c4762a1bSJed Brown 
319c4762a1bSJed Brown    test:
320c4762a1bSJed Brown       args: -snes_monitor_short -snes_type nrichardson
321c4762a1bSJed Brown       requires: !single
322c4762a1bSJed Brown 
323c4762a1bSJed Brown    test:
324c4762a1bSJed Brown       suffix: 2
325c4762a1bSJed Brown       args: -snes_monitor_short -ksp_monitor_short -ksp_type richardson -pc_type none -ksp_richardson_self_scale
326c4762a1bSJed Brown       requires: !single
327c4762a1bSJed Brown 
328c4762a1bSJed Brown    test:
329c4762a1bSJed Brown       suffix: 3
330c4762a1bSJed Brown       args: -snes_monitor_short -snes_type ngmres
331c4762a1bSJed Brown 
332c4762a1bSJed Brown    test:
333c4762a1bSJed Brown       suffix: 4
334c4762a1bSJed Brown       args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none
335c4762a1bSJed Brown 
336c4762a1bSJed Brown    test:
337c4762a1bSJed Brown       suffix: 5
338c4762a1bSJed Brown       args: -snes_monitor_short -snes_type ncg
339c4762a1bSJed Brown 
340c4762a1bSJed Brown    test:
341c4762a1bSJed Brown       suffix: 6
342c4762a1bSJed Brown       args: -snes_monitor_short -ksp_type cg -ksp_monitor_short -pc_type none
343c4762a1bSJed Brown 
344c4762a1bSJed Brown    test:
345c4762a1bSJed Brown       suffix: 7
346c4762a1bSJed Brown       args: -da_refine 2 -snes_monitor_short -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type none -mg_levels_ksp_monitor_short -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor_short
347c4762a1bSJed Brown       requires: !single
348c4762a1bSJed Brown 
349c4762a1bSJed Brown    test:
350c4762a1bSJed Brown       suffix: 8
351c4762a1bSJed Brown       args: -da_refine 2 -snes_monitor_short -snes_type fas -fas_levels_snes_monitor_short -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_type fas -snes_rtol 1.e-5
352c4762a1bSJed Brown 
35341ba4c6cSHeeho Park    test:
35441ba4c6cSHeeho Park       suffix: 9
35541ba4c6cSHeeho Park       args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc
35641ba4c6cSHeeho Park 
35741ba4c6cSHeeho Park    test:
35841ba4c6cSHeeho Park       suffix: 10
35941ba4c6cSHeeho Park       args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc -snes_trdc_use_cauchy false
36041ba4c6cSHeeho Park 
361c4762a1bSJed Brown TEST*/
362