xref: /petsc/src/snes/tutorials/ex9.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static const char help[] = "Solves obstacle problem in 2D as a variational inequality\n\
2*c4762a1bSJed Brown or nonlinear complementarity problem.  This is a form of the Laplace equation in\n\
3*c4762a1bSJed Brown which the solution u is constrained to be above a given function psi.  In the\n\
4*c4762a1bSJed Brown problem here an exact solution is known.\n";
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown /*  On a square S = {-2<x<2,-2<y<2}, the PDE
7*c4762a1bSJed Brown     u_{xx} + u_{yy} = 0
8*c4762a1bSJed Brown is solved on the set where membrane is above obstacle (u(x,y) >= psi(x,y)).
9*c4762a1bSJed Brown Here psi is the upper hemisphere of the unit ball.  On the boundary of S
10*c4762a1bSJed Brown we have Dirichlet boundary conditions from the exact solution.  Uses centered
11*c4762a1bSJed Brown FD scheme.  This example contributed by Ed Bueler.
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown Example usage:
14*c4762a1bSJed Brown   * get help:
15*c4762a1bSJed Brown     ./ex9 -help
16*c4762a1bSJed Brown   * monitor run:
17*c4762a1bSJed Brown     ./ex9 -da_refine 2 -snes_vi_monitor
18*c4762a1bSJed Brown   * use other SNESVI type (default is SNESVINEWTONRSLS):
19*c4762a1bSJed Brown     ./ex9 -da_refine 2 -snes_vi_monitor -snes_type vinewtonssls
20*c4762a1bSJed Brown   * use FD evaluation of Jacobian by coloring, instead of analytical:
21*c4762a1bSJed Brown     ./ex9 -da_refine 2 -snes_fd_color
22*c4762a1bSJed Brown   * X windows visualizations:
23*c4762a1bSJed Brown     ./ex9 -snes_monitor_solution draw -draw_pause 1 -da_refine 4
24*c4762a1bSJed Brown     ./ex9 -snes_vi_monitor_residual -draw_pause 1 -da_refine 4
25*c4762a1bSJed Brown   * full-cycle multigrid:
26*c4762a1bSJed Brown     ./ex9 -snes_converged_reason -snes_grid_sequence 4 -pc_type mg
27*c4762a1bSJed Brown   * serial convergence evidence:
28*c4762a1bSJed Brown     for M in 3 4 5 6 7; do ./ex9 -snes_grid_sequence $M -pc_type mg; done
29*c4762a1bSJed Brown   * FIXME sporadic parallel bug:
30*c4762a1bSJed Brown     mpiexec -n 4 ./ex9 -snes_converged_reason -snes_grid_sequence 4 -pc_type mg
31*c4762a1bSJed Brown */
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown #include <petsc.h>
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown /* z = psi(x,y) is the hemispherical obstacle, but made C^1 with "skirt" at r=r0 */
36*c4762a1bSJed Brown PetscReal psi(PetscReal x, PetscReal y)
37*c4762a1bSJed Brown {
38*c4762a1bSJed Brown     const PetscReal  r = x * x + y * y,r0 = 0.9,psi0 = PetscSqrtReal(1.0 - r0*r0),dpsi0 = - r0 / psi0;
39*c4762a1bSJed Brown     if (r <= r0) {
40*c4762a1bSJed Brown       return PetscSqrtReal(1.0 - r);
41*c4762a1bSJed Brown     } else {
42*c4762a1bSJed Brown       return psi0 + dpsi0 * (r - r0);
43*c4762a1bSJed Brown     }
44*c4762a1bSJed Brown }
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown /*  This exact solution solves a 1D radial free-boundary problem for the
47*c4762a1bSJed Brown Laplace equation, on the interval 0 < r < 2, with above obstacle psi(x,y).
48*c4762a1bSJed Brown The Laplace equation applies where u(r) > psi(r),
49*c4762a1bSJed Brown     u''(r) + r^-1 u'(r) = 0
50*c4762a1bSJed Brown with boundary conditions including free b.c.s at an unknown location r = a:
51*c4762a1bSJed Brown     u(a) = psi(a),  u'(a) = psi'(a),  u(2) = 0
52*c4762a1bSJed Brown The solution is  u(r) = - A log(r) + B   on  r > a.  The boundary conditions
53*c4762a1bSJed Brown can then be reduced to a root-finding problem for a:
54*c4762a1bSJed Brown     a^2 (log(2) - log(a)) = 1 - a^2
55*c4762a1bSJed Brown The solution is a = 0.697965148223374 (giving residual 1.5e-15).  Then
56*c4762a1bSJed Brown A = a^2*(1-a^2)^(-0.5) and B = A*log(2) are as given below in the code.  */
57*c4762a1bSJed Brown PetscReal u_exact(PetscReal x, PetscReal y)
58*c4762a1bSJed Brown {
59*c4762a1bSJed Brown     const PetscReal afree = 0.697965148223374,
60*c4762a1bSJed Brown                     A     = 0.680259411891719,
61*c4762a1bSJed Brown                     B     = 0.471519893402112;
62*c4762a1bSJed Brown     PetscReal  r;
63*c4762a1bSJed Brown     r = PetscSqrtReal(x * x + y * y);
64*c4762a1bSJed Brown     return (r <= afree) ? psi(x,y)  /* active set; on the obstacle */
65*c4762a1bSJed Brown                         : - A * PetscLogReal(r) + B; /* solves laplace eqn */
66*c4762a1bSJed Brown }
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown extern PetscErrorCode FormExactSolution(DMDALocalInfo*,Vec);
69*c4762a1bSJed Brown extern PetscErrorCode FormBounds(SNES,Vec,Vec);
70*c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscReal**,PetscReal**,void*);
71*c4762a1bSJed Brown extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscReal**,Mat,Mat,void*);
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown int main(int argc,char **argv)
74*c4762a1bSJed Brown {
75*c4762a1bSJed Brown   PetscErrorCode      ierr;
76*c4762a1bSJed Brown   SNES                snes;
77*c4762a1bSJed Brown   DM                  da, da_after;
78*c4762a1bSJed Brown   Vec                 u, u_exact;
79*c4762a1bSJed Brown   DMDALocalInfo       info;
80*c4762a1bSJed Brown   PetscReal           error1,errorinf;
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
83*c4762a1bSJed Brown 
84*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,
85*c4762a1bSJed Brown                       DMDA_STENCIL_STAR,5,5, /* 5x5 coarse grid; override with -da_grid_x,_y */
86*c4762a1bSJed Brown                       PETSC_DECIDE,PETSC_DECIDE,
87*c4762a1bSJed Brown                       1,1,  /* dof=1 and s = 1 (stencil extends out one cell) */
88*c4762a1bSJed Brown                       NULL,NULL,&da);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,-2.0,2.0,-2.0,2.0,0.0,1.0);CHKERRQ(ierr);
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&u);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = VecSet(u,0.0);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = SNESSetDM(snes,da);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = SNESVISetComputeVariableBounds(snes,&FormBounds);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,NULL);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,NULL);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown   /* solve nonlinear system */
105*c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,u);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = VecDestroy(&u); CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = DMDestroy(&da); CHKERRQ(ierr);
108*c4762a1bSJed Brown   /* DMDA after solve may be different, e.g. with -snes_grid_sequence */
109*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&da_after); CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = SNESGetSolution(snes,&u); CHKERRQ(ierr); /* do not destroy u */
111*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da_after,&info); CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = VecDuplicate(u,&u_exact); CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = FormExactSolution(&info,u_exact); CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = VecAXPY(u,-1.0,u_exact);CHKERRQ(ierr); /* u <-- u - u_exact */
115*c4762a1bSJed Brown   ierr = VecNorm(u,NORM_1,&error1);CHKERRQ(ierr);
116*c4762a1bSJed Brown   error1 /= (PetscReal)info.mx * (PetscReal)info.my; /* average error */
117*c4762a1bSJed Brown   ierr = VecNorm(u,NORM_INFINITY,&errorinf);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"errors on %D x %D grid:  av |u-uexact|  = %.3e,  |u-uexact|_inf = %.3e\n",info.mx,info.my,(double)error1,(double)errorinf);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = VecDestroy(&u_exact); CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscFinalize();
123*c4762a1bSJed Brown   return ierr;
124*c4762a1bSJed Brown }
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown PetscErrorCode FormExactSolution(DMDALocalInfo *info, Vec u)
127*c4762a1bSJed Brown {
128*c4762a1bSJed Brown   PetscErrorCode ierr;
129*c4762a1bSJed Brown   PetscInt       i,j;
130*c4762a1bSJed Brown   PetscReal      **au, dx, dy, x, y;
131*c4762a1bSJed Brown   dx = 4.0 / (PetscReal)(info->mx-1);
132*c4762a1bSJed Brown   dy = 4.0 / (PetscReal)(info->my-1);
133*c4762a1bSJed Brown   ierr = DMDAVecGetArray(info->da, u, &au);CHKERRQ(ierr);
134*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
135*c4762a1bSJed Brown     y = -2.0 + j * dy;
136*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
137*c4762a1bSJed Brown       x = -2.0 + i * dx;
138*c4762a1bSJed Brown       au[j][i] = u_exact(x,y);
139*c4762a1bSJed Brown     }
140*c4762a1bSJed Brown   }
141*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(info->da, u, &au);CHKERRQ(ierr);
142*c4762a1bSJed Brown   return 0;
143*c4762a1bSJed Brown }
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown PetscErrorCode FormBounds(SNES snes, Vec Xl, Vec Xu)
146*c4762a1bSJed Brown {
147*c4762a1bSJed Brown   PetscErrorCode ierr;
148*c4762a1bSJed Brown   DM             da;
149*c4762a1bSJed Brown   DMDALocalInfo  info;
150*c4762a1bSJed Brown   PetscInt       i, j;
151*c4762a1bSJed Brown   PetscReal      **aXl, dx, dy, x, y;
152*c4762a1bSJed Brown 
153*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&da);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info); CHKERRQ(ierr);
155*c4762a1bSJed Brown   dx = 4.0 / (PetscReal)(info.mx-1);
156*c4762a1bSJed Brown   dy = 4.0 / (PetscReal)(info.my-1);
157*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da, Xl, &aXl);CHKERRQ(ierr);
158*c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
159*c4762a1bSJed Brown     y = -2.0 + j * dy;
160*c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
161*c4762a1bSJed Brown       x = -2.0 + i * dx;
162*c4762a1bSJed Brown       aXl[j][i] = psi(x,y);
163*c4762a1bSJed Brown     }
164*c4762a1bSJed Brown   }
165*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da, Xl, &aXl);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = VecSet(Xu,PETSC_INFINITY);CHKERRQ(ierr);
167*c4762a1bSJed Brown   return 0;
168*c4762a1bSJed Brown }
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **au, PetscScalar **af, void *user)
171*c4762a1bSJed Brown {
172*c4762a1bSJed Brown   PetscErrorCode ierr;
173*c4762a1bSJed Brown   PetscInt       i,j;
174*c4762a1bSJed Brown   PetscReal      dx,dy,x,y,ue,un,us,uw;
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown   PetscFunctionBeginUser;
177*c4762a1bSJed Brown   dx = 4.0 / (PetscReal)(info->mx-1);
178*c4762a1bSJed Brown   dy = 4.0 / (PetscReal)(info->my-1);
179*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
180*c4762a1bSJed Brown     y = -2.0 + j * dy;
181*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
182*c4762a1bSJed Brown       x = -2.0 + i * dx;
183*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
184*c4762a1bSJed Brown         af[j][i] = 4.0 * (au[j][i] - u_exact(x,y));
185*c4762a1bSJed Brown       } else {
186*c4762a1bSJed Brown         uw = (i-1 == 0)          ? u_exact(x-dx,y) : au[j][i-1];
187*c4762a1bSJed Brown         ue = (i+1 == info->mx-1) ? u_exact(x+dx,y) : au[j][i+1];
188*c4762a1bSJed Brown         us = (j-1 == 0)          ? u_exact(x,y-dy) : au[j-1][i];
189*c4762a1bSJed Brown         un = (j+1 == info->my-1) ? u_exact(x,y+dy) : au[j+1][i];
190*c4762a1bSJed Brown         af[j][i] = - (dy/dx) * (uw - 2.0 * au[j][i] + ue) - (dx/dy) * (us - 2.0 * au[j][i] + un);
191*c4762a1bSJed Brown       }
192*c4762a1bSJed Brown     }
193*c4762a1bSJed Brown   }
194*c4762a1bSJed Brown   ierr = PetscLogFlops(12.0*info->ym*info->xm);CHKERRQ(ierr);
195*c4762a1bSJed Brown   PetscFunctionReturn(0);
196*c4762a1bSJed Brown }
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **au, Mat A, Mat jac, void *user)
199*c4762a1bSJed Brown {
200*c4762a1bSJed Brown   PetscErrorCode ierr;
201*c4762a1bSJed Brown   PetscInt       i,j,n;
202*c4762a1bSJed Brown   MatStencil     col[5],row;
203*c4762a1bSJed Brown   PetscReal      v[5],dx,dy,oxx,oyy;
204*c4762a1bSJed Brown 
205*c4762a1bSJed Brown   PetscFunctionBeginUser;
206*c4762a1bSJed Brown   dx  = 4.0 / (PetscReal)(info->mx-1);
207*c4762a1bSJed Brown   dy  = 4.0 / (PetscReal)(info->my-1);
208*c4762a1bSJed Brown   oxx = dy / dx;
209*c4762a1bSJed Brown   oyy = dx / dy;
210*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
211*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
212*c4762a1bSJed Brown       row.j = j; row.i = i;
213*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { /* boundary */
214*c4762a1bSJed Brown         v[0] = 4.0;
215*c4762a1bSJed Brown         ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
216*c4762a1bSJed Brown       } else { /* interior grid points */
217*c4762a1bSJed Brown         v[0] = 2.0 * (oxx + oyy);  col[0].j = j;  col[0].i = i;
218*c4762a1bSJed Brown         n = 1;
219*c4762a1bSJed Brown         if (i-1 > 0) {
220*c4762a1bSJed Brown           v[n] = -oxx;  col[n].j = j;  col[n++].i = i-1;
221*c4762a1bSJed Brown         }
222*c4762a1bSJed Brown         if (i+1 < info->mx-1) {
223*c4762a1bSJed Brown           v[n] = -oxx;  col[n].j = j;  col[n++].i = i+1;
224*c4762a1bSJed Brown         }
225*c4762a1bSJed Brown         if (j-1 > 0) {
226*c4762a1bSJed Brown           v[n] = -oyy;  col[n].j = j-1;  col[n++].i = i;
227*c4762a1bSJed Brown         }
228*c4762a1bSJed Brown         if (j+1 < info->my-1) {
229*c4762a1bSJed Brown           v[n] = -oyy;  col[n].j = j+1;  col[n++].i = i;
230*c4762a1bSJed Brown         }
231*c4762a1bSJed Brown         ierr = MatSetValuesStencil(jac,1,&row,n,col,v,INSERT_VALUES);CHKERRQ(ierr);
232*c4762a1bSJed Brown       }
233*c4762a1bSJed Brown     }
234*c4762a1bSJed Brown   }
235*c4762a1bSJed Brown 
236*c4762a1bSJed Brown   /* Assemble matrix, using the 2-step process: */
237*c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
239*c4762a1bSJed Brown   if (A != jac) {
240*c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
241*c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
242*c4762a1bSJed Brown   }
243*c4762a1bSJed Brown   ierr = PetscLogFlops(2.0*info->ym*info->xm);CHKERRQ(ierr);
244*c4762a1bSJed Brown   PetscFunctionReturn(0);
245*c4762a1bSJed Brown }
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown /*TEST
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown    build:
250*c4762a1bSJed Brown       requires: !complex
251*c4762a1bSJed Brown 
252*c4762a1bSJed Brown    test:
253*c4762a1bSJed Brown       suffix: 1
254*c4762a1bSJed Brown       requires: !single
255*c4762a1bSJed Brown       nsize: 1
256*c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonrsls
257*c4762a1bSJed Brown 
258*c4762a1bSJed Brown    test:
259*c4762a1bSJed Brown       suffix: 2
260*c4762a1bSJed Brown       requires: !single
261*c4762a1bSJed Brown       nsize: 2
262*c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonssls
263*c4762a1bSJed Brown 
264*c4762a1bSJed Brown    test:
265*c4762a1bSJed Brown       suffix: 3
266*c4762a1bSJed Brown       requires: !single
267*c4762a1bSJed Brown       nsize: 2
268*c4762a1bSJed Brown       args: -snes_grid_sequence 2 -snes_vi_monitor -snes_type vinewtonrsls
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown    test:
271*c4762a1bSJed Brown       suffix: mg
272*c4762a1bSJed Brown       requires: !single
273*c4762a1bSJed Brown       nsize: 4
274*c4762a1bSJed Brown       args: -snes_grid_sequence 3 -snes_converged_reason -pc_type mg
275*c4762a1bSJed Brown 
276*c4762a1bSJed Brown    test:
277*c4762a1bSJed Brown       suffix: 4
278*c4762a1bSJed Brown       nsize: 1
279*c4762a1bSJed Brown       args: -mat_is_symmetric
280*c4762a1bSJed Brown 
281*c4762a1bSJed Brown    test:
282*c4762a1bSJed Brown       suffix: 5
283*c4762a1bSJed Brown       nsize: 1
284*c4762a1bSJed Brown       args: -ksp_converged_reason -snes_fd_color
285*c4762a1bSJed Brown 
286*c4762a1bSJed Brown    test:
287*c4762a1bSJed Brown       suffix: 6
288*c4762a1bSJed Brown       requires: !single
289*c4762a1bSJed Brown       nsize: 2
290*c4762a1bSJed Brown       args: -snes_grid_sequence 2 -pc_type mg -snes_monitor_short -ksp_converged_reason
291*c4762a1bSJed Brown 
292*c4762a1bSJed Brown    test:
293*c4762a1bSJed Brown       suffix: 7
294*c4762a1bSJed Brown       nsize: 2
295*c4762a1bSJed 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
296*c4762a1bSJed Brown       TODO: fix nasty memory leak in SNESCOMPOSITE
297*c4762a1bSJed Brown 
298*c4762a1bSJed Brown    test:
299*c4762a1bSJed Brown       suffix: 8
300*c4762a1bSJed Brown       nsize: 2
301*c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additive -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
302*c4762a1bSJed Brown       TODO: fix nasty memory leak in SNESCOMPOSITE
303*c4762a1bSJed Brown 
304*c4762a1bSJed Brown    test:
305*c4762a1bSJed Brown       suffix: 9
306*c4762a1bSJed Brown       nsize: 2
307*c4762a1bSJed Brown       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additiveoptimal -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
308*c4762a1bSJed Brown       TODO: fix nasty memory leak in SNESCOMPOSITE
309*c4762a1bSJed Brown 
310*c4762a1bSJed Brown TEST*/
311