xref: /petsc/src/snes/tutorials/ex9.c (revision 8434afd195968570cfdb5bc7b9cfc0a316d974ae)
1c4762a1bSJed Brown static const char help[] = "Solves obstacle problem in 2D as a variational inequality\n\
2c4762a1bSJed Brown or nonlinear complementarity problem.  This is a form of the Laplace equation in\n\
3c4762a1bSJed Brown which the solution u is constrained to be above a given function psi.  In the\n\
4c4762a1bSJed Brown problem here an exact solution is known.\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /*  On a square S = {-2<x<2,-2<y<2}, the PDE
7c4762a1bSJed Brown     u_{xx} + u_{yy} = 0
8c4762a1bSJed Brown is solved on the set where membrane is above obstacle (u(x,y) >= psi(x,y)).
9c4762a1bSJed Brown Here psi is the upper hemisphere of the unit ball.  On the boundary of S
10c4762a1bSJed Brown we have Dirichlet boundary conditions from the exact solution.  Uses centered
11c4762a1bSJed Brown FD scheme.  This example contributed by Ed Bueler.
12c4762a1bSJed Brown 
13c4762a1bSJed Brown Example usage:
14c4762a1bSJed Brown   * get help:
15c4762a1bSJed Brown     ./ex9 -help
16c4762a1bSJed Brown   * monitor run:
17c4762a1bSJed Brown     ./ex9 -da_refine 2 -snes_vi_monitor
18c4762a1bSJed Brown   * use other SNESVI type (default is SNESVINEWTONRSLS):
19c4762a1bSJed Brown     ./ex9 -da_refine 2 -snes_vi_monitor -snes_type vinewtonssls
20c4762a1bSJed Brown   * use FD evaluation of Jacobian by coloring, instead of analytical:
21c4762a1bSJed Brown     ./ex9 -da_refine 2 -snes_fd_color
22c4762a1bSJed Brown   * X windows visualizations:
23c4762a1bSJed Brown     ./ex9 -snes_monitor_solution draw -draw_pause 1 -da_refine 4
24c4762a1bSJed Brown     ./ex9 -snes_vi_monitor_residual -draw_pause 1 -da_refine 4
25c4762a1bSJed Brown   * full-cycle multigrid:
26c4762a1bSJed Brown     ./ex9 -snes_converged_reason -snes_grid_sequence 4 -pc_type mg
27c4762a1bSJed Brown   * serial convergence evidence:
28c4762a1bSJed Brown     for M in 3 4 5 6 7; do ./ex9 -snes_grid_sequence $M -pc_type mg; done
29c4762a1bSJed Brown   * FIXME sporadic parallel bug:
30c4762a1bSJed Brown     mpiexec -n 4 ./ex9 -snes_converged_reason -snes_grid_sequence 4 -pc_type mg
31c4762a1bSJed Brown */
32c4762a1bSJed Brown 
33c4762a1bSJed Brown #include <petsc.h>
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /* z = psi(x,y) is the hemispherical obstacle, but made C^1 with "skirt" at r=r0 */
36d71ae5a4SJacob Faibussowitsch PetscReal psi(PetscReal x, PetscReal y)
37d71ae5a4SJacob Faibussowitsch {
38c4762a1bSJed Brown   const PetscReal r = x * x + y * y, r0 = 0.9, psi0 = PetscSqrtReal(1.0 - r0 * r0), dpsi0 = -r0 / psi0;
39c4762a1bSJed Brown   if (r <= r0) {
40c4762a1bSJed Brown     return PetscSqrtReal(1.0 - r);
41c4762a1bSJed Brown   } else {
42c4762a1bSJed Brown     return psi0 + dpsi0 * (r - r0);
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*  This exact solution solves a 1D radial free-boundary problem for the
47c4762a1bSJed Brown Laplace equation, on the interval 0 < r < 2, with above obstacle psi(x,y).
48c4762a1bSJed Brown The Laplace equation applies where u(r) > psi(r),
49c4762a1bSJed Brown     u''(r) + r^-1 u'(r) = 0
50c4762a1bSJed Brown with boundary conditions including free b.c.s at an unknown location r = a:
51c4762a1bSJed Brown     u(a) = psi(a),  u'(a) = psi'(a),  u(2) = 0
52c4762a1bSJed Brown The solution is  u(r) = - A log(r) + B   on  r > a.  The boundary conditions
53c4762a1bSJed Brown can then be reduced to a root-finding problem for a:
54c4762a1bSJed Brown     a^2 (log(2) - log(a)) = 1 - a^2
55c4762a1bSJed Brown The solution is a = 0.697965148223374 (giving residual 1.5e-15).  Then
56c4762a1bSJed Brown A = a^2*(1-a^2)^(-0.5) and B = A*log(2) are as given below in the code.  */
57d71ae5a4SJacob Faibussowitsch PetscReal u_exact(PetscReal x, PetscReal y)
58d71ae5a4SJacob Faibussowitsch {
599371c9d4SSatish Balay   const PetscReal afree = 0.697965148223374, A = 0.680259411891719, B = 0.471519893402112;
60c4762a1bSJed Brown   PetscReal       r;
61c4762a1bSJed Brown   r = PetscSqrtReal(x * x + y * y);
62c4762a1bSJed Brown   return (r <= afree) ? psi(x, y)                 /* active set; on the obstacle */
63c4762a1bSJed Brown                       : -A * PetscLogReal(r) + B; /* solves laplace eqn */
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown extern PetscErrorCode FormExactSolution(DMDALocalInfo *, Vec);
67c4762a1bSJed Brown extern PetscErrorCode FormBounds(SNES, Vec, Vec);
68c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo *, PetscReal **, PetscReal **, void *);
69c4762a1bSJed Brown extern PetscErrorCode FormJacobianLocal(DMDALocalInfo *, PetscReal **, Mat, Mat, void *);
70c4762a1bSJed Brown 
71d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
72d71ae5a4SJacob Faibussowitsch {
73c4762a1bSJed Brown   SNES          snes;
74c4762a1bSJed Brown   DM            da, da_after;
75c4762a1bSJed Brown   Vec           u, u_exact;
76c4762a1bSJed Brown   DMDALocalInfo info;
77c4762a1bSJed Brown   PetscReal     error1, errorinf;
78c4762a1bSJed Brown 
79327415f7SBarry Smith   PetscFunctionBeginUser;
809566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
81c4762a1bSJed Brown 
82d0609cedSBarry Smith   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 5, 5, /* 5x5 coarse grid; override with -da_grid_x,_y */
83d0609cedSBarry Smith                          PETSC_DECIDE, PETSC_DECIDE, 1, 1,                                              /* dof=1 and s = 1 (stencil extends out one cell) */
84d0609cedSBarry Smith                          NULL, NULL, &da));
859566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
869566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
879566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, -2.0, 2.0, -2.0, 2.0, 0.0, 1.0));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &u));
909566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
91c4762a1bSJed Brown 
929566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
939566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
949566063dSJacob Faibussowitsch   PetscCall(SNESSetType(snes, SNESVINEWTONRSLS));
959566063dSJacob Faibussowitsch   PetscCall(SNESVISetComputeVariableBounds(snes, &FormBounds));
96*8434afd1SBarry Smith   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunctionFn *)FormFunctionLocal, NULL));
97*8434afd1SBarry Smith   PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobianFn *)FormJacobianLocal, NULL));
989566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* solve nonlinear system */
1019566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
1029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
104c4762a1bSJed Brown   /* DMDA after solve may be different, e.g. with -snes_grid_sequence */
1059566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da_after));
1069566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u)); /* do not destroy u */
1079566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da_after, &info));
1089566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &u_exact));
1099566063dSJacob Faibussowitsch   PetscCall(FormExactSolution(&info, u_exact));
1109566063dSJacob Faibussowitsch   PetscCall(VecAXPY(u, -1.0, u_exact)); /* u <-- u - u_exact */
1119566063dSJacob Faibussowitsch   PetscCall(VecNorm(u, NORM_1, &error1));
112c4762a1bSJed Brown   error1 /= (PetscReal)info.mx * (PetscReal)info.my; /* average error */
1139566063dSJacob Faibussowitsch   PetscCall(VecNorm(u, NORM_INFINITY, &errorinf));
11463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "errors on %" PetscInt_FMT " x %" PetscInt_FMT " grid:  av |u-uexact|  = %.3e,  |u-uexact|_inf = %.3e\n", info.mx, info.my, (double)error1, (double)errorinf));
1159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u_exact));
1169566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1179566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1189566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
119b122ec5aSJacob Faibussowitsch   return 0;
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122d71ae5a4SJacob Faibussowitsch PetscErrorCode FormExactSolution(DMDALocalInfo *info, Vec u)
123d71ae5a4SJacob Faibussowitsch {
124c4762a1bSJed Brown   PetscInt    i, j;
125c4762a1bSJed Brown   PetscReal **au, dx, dy, x, y;
1263ba16761SJacob Faibussowitsch 
1273ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
128c4762a1bSJed Brown   dx = 4.0 / (PetscReal)(info->mx - 1);
129c4762a1bSJed Brown   dy = 4.0 / (PetscReal)(info->my - 1);
1309566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(info->da, u, &au));
131c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
132c4762a1bSJed Brown     y = -2.0 + j * dy;
133c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
134c4762a1bSJed Brown       x        = -2.0 + i * dx;
135c4762a1bSJed Brown       au[j][i] = u_exact(x, y);
136c4762a1bSJed Brown     }
137c4762a1bSJed Brown   }
1389566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(info->da, u, &au));
1393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
140c4762a1bSJed Brown }
141c4762a1bSJed Brown 
142d71ae5a4SJacob Faibussowitsch PetscErrorCode FormBounds(SNES snes, Vec Xl, Vec Xu)
143d71ae5a4SJacob Faibussowitsch {
144c4762a1bSJed Brown   DM            da;
145c4762a1bSJed Brown   DMDALocalInfo info;
146c4762a1bSJed Brown   PetscInt      i, j;
147c4762a1bSJed Brown   PetscReal   **aXl, dx, dy, x, y;
148c4762a1bSJed Brown 
1493ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1509566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
1519566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
152c4762a1bSJed Brown   dx = 4.0 / (PetscReal)(info.mx - 1);
153c4762a1bSJed Brown   dy = 4.0 / (PetscReal)(info.my - 1);
1549566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xl, &aXl));
155c4762a1bSJed Brown   for (j = info.ys; j < info.ys + info.ym; j++) {
156c4762a1bSJed Brown     y = -2.0 + j * dy;
157c4762a1bSJed Brown     for (i = info.xs; i < info.xs + info.xm; i++) {
158c4762a1bSJed Brown       x         = -2.0 + i * dx;
159c4762a1bSJed Brown       aXl[j][i] = psi(x, y);
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown   }
1629566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xl, &aXl));
1639566063dSJacob Faibussowitsch   PetscCall(VecSet(Xu, PETSC_INFINITY));
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
167d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **au, PetscScalar **af, void *user)
168d71ae5a4SJacob Faibussowitsch {
169c4762a1bSJed Brown   PetscInt  i, j;
170c4762a1bSJed Brown   PetscReal dx, dy, x, y, ue, un, us, uw;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   PetscFunctionBeginUser;
173c4762a1bSJed Brown   dx = 4.0 / (PetscReal)(info->mx - 1);
174c4762a1bSJed Brown   dy = 4.0 / (PetscReal)(info->my - 1);
175c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
176c4762a1bSJed Brown     y = -2.0 + j * dy;
177c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
178c4762a1bSJed Brown       x = -2.0 + i * dx;
179c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
180c4762a1bSJed Brown         af[j][i] = 4.0 * (au[j][i] - u_exact(x, y));
181c4762a1bSJed Brown       } else {
182c4762a1bSJed Brown         uw       = (i - 1 == 0) ? u_exact(x - dx, y) : au[j][i - 1];
183c4762a1bSJed Brown         ue       = (i + 1 == info->mx - 1) ? u_exact(x + dx, y) : au[j][i + 1];
184c4762a1bSJed Brown         us       = (j - 1 == 0) ? u_exact(x, y - dy) : au[j - 1][i];
185c4762a1bSJed Brown         un       = (j + 1 == info->my - 1) ? u_exact(x, y + dy) : au[j + 1][i];
186c4762a1bSJed Brown         af[j][i] = -(dy / dx) * (uw - 2.0 * au[j][i] + ue) - (dx / dy) * (us - 2.0 * au[j][i] + un);
187c4762a1bSJed Brown       }
188c4762a1bSJed Brown     }
189c4762a1bSJed Brown   }
1909566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(12.0 * info->ym * info->xm));
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **au, Mat A, Mat jac, void *user)
195d71ae5a4SJacob Faibussowitsch {
196c4762a1bSJed Brown   PetscInt   i, j, n;
197c4762a1bSJed Brown   MatStencil col[5], row;
198c4762a1bSJed Brown   PetscReal  v[5], dx, dy, oxx, oyy;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   PetscFunctionBeginUser;
201c4762a1bSJed Brown   dx  = 4.0 / (PetscReal)(info->mx - 1);
202c4762a1bSJed Brown   dy  = 4.0 / (PetscReal)(info->my - 1);
203c4762a1bSJed Brown   oxx = dy / dx;
204c4762a1bSJed Brown   oyy = dx / dy;
205c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
206c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
2079371c9d4SSatish Balay       row.j = j;
2089371c9d4SSatish Balay       row.i = i;
209c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { /* boundary */
210c4762a1bSJed Brown         v[0] = 4.0;
2119566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES));
212c4762a1bSJed Brown       } else { /* interior grid points */
2139371c9d4SSatish Balay         v[0]     = 2.0 * (oxx + oyy);
2149371c9d4SSatish Balay         col[0].j = j;
2159371c9d4SSatish Balay         col[0].i = i;
216c4762a1bSJed Brown         n        = 1;
217c4762a1bSJed Brown         if (i - 1 > 0) {
2189371c9d4SSatish Balay           v[n]       = -oxx;
2199371c9d4SSatish Balay           col[n].j   = j;
2209371c9d4SSatish Balay           col[n++].i = i - 1;
221c4762a1bSJed Brown         }
222c4762a1bSJed Brown         if (i + 1 < info->mx - 1) {
2239371c9d4SSatish Balay           v[n]       = -oxx;
2249371c9d4SSatish Balay           col[n].j   = j;
2259371c9d4SSatish Balay           col[n++].i = i + 1;
226c4762a1bSJed Brown         }
227c4762a1bSJed Brown         if (j - 1 > 0) {
2289371c9d4SSatish Balay           v[n]       = -oyy;
2299371c9d4SSatish Balay           col[n].j   = j - 1;
2309371c9d4SSatish Balay           col[n++].i = i;
231c4762a1bSJed Brown         }
232c4762a1bSJed Brown         if (j + 1 < info->my - 1) {
2339371c9d4SSatish Balay           v[n]       = -oyy;
2349371c9d4SSatish Balay           col[n].j   = j + 1;
2359371c9d4SSatish Balay           col[n++].i = i;
236c4762a1bSJed Brown         }
2379566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jac, 1, &row, n, col, v, INSERT_VALUES));
238c4762a1bSJed Brown       }
239c4762a1bSJed Brown     }
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Assemble matrix, using the 2-step process: */
2439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
245c4762a1bSJed Brown   if (A != jac) {
2469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
248c4762a1bSJed Brown   }
2499566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * info->ym * info->xm));
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251c4762a1bSJed Brown }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown /*TEST
254c4762a1bSJed Brown 
255c4762a1bSJed Brown    build:
256c4762a1bSJed Brown       requires: !complex
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    test:
259c4762a1bSJed Brown       suffix: 1
260c4762a1bSJed Brown       requires: !single
261c4762a1bSJed Brown       nsize: 1
262c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonrsls
263c4762a1bSJed Brown 
264c4762a1bSJed Brown    test:
265c4762a1bSJed Brown       suffix: 2
266c4762a1bSJed Brown       requires: !single
267c4762a1bSJed Brown       nsize: 2
268c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonssls
269c4762a1bSJed Brown 
270c4762a1bSJed Brown    test:
271c4762a1bSJed Brown       suffix: 3
272c4762a1bSJed Brown       requires: !single
273c4762a1bSJed Brown       nsize: 2
274c4762a1bSJed Brown       args: -snes_grid_sequence 2 -snes_vi_monitor -snes_type vinewtonrsls
275c4762a1bSJed Brown 
276c4762a1bSJed Brown    test:
277c4762a1bSJed Brown       suffix: mg
278c4762a1bSJed Brown       requires: !single
279c4762a1bSJed Brown       nsize: 4
280c4762a1bSJed Brown       args: -snes_grid_sequence 3 -snes_converged_reason -pc_type mg
281c4762a1bSJed Brown 
282c4762a1bSJed Brown    test:
283c4762a1bSJed Brown       suffix: 4
284c4762a1bSJed Brown       nsize: 1
285c4762a1bSJed Brown       args: -mat_is_symmetric
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    test:
288c4762a1bSJed Brown       suffix: 5
289c4762a1bSJed Brown       nsize: 1
290c4762a1bSJed Brown       args: -ksp_converged_reason -snes_fd_color
291c4762a1bSJed Brown 
292c4762a1bSJed Brown    test:
293c4762a1bSJed Brown       suffix: 6
294c4762a1bSJed Brown       requires: !single
295c4762a1bSJed Brown       nsize: 2
296c4762a1bSJed Brown       args: -snes_grid_sequence 2 -pc_type mg -snes_monitor_short -ksp_converged_reason
297c4762a1bSJed Brown 
298c4762a1bSJed Brown    test:
299c4762a1bSJed Brown       suffix: 7
300c4762a1bSJed Brown       nsize: 2
301c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type multiplicative -snes_composite_sneses vinewtonrsls,vinewtonssls -sub_0_snes_vi_monitor -sub_1_snes_vi_monitor
302c4762a1bSJed Brown       TODO: fix nasty memory leak in SNESCOMPOSITE
303c4762a1bSJed Brown 
304c4762a1bSJed Brown    test:
305c4762a1bSJed Brown       suffix: 8
306c4762a1bSJed Brown       nsize: 2
307c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additive -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
308c4762a1bSJed Brown       TODO: fix nasty memory leak in SNESCOMPOSITE
309c4762a1bSJed Brown 
310c4762a1bSJed Brown    test:
311c4762a1bSJed Brown       suffix: 9
312c4762a1bSJed Brown       nsize: 2
313c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additiveoptimal -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
314c4762a1bSJed Brown       TODO: fix nasty memory leak in SNESCOMPOSITE
315c4762a1bSJed Brown 
316c4762a1bSJed Brown TEST*/
317