xref: /petsc/src/snes/tutorials/ex15.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
1 static const char help[] = "p-Bratu nonlinear PDE in 2d.\n\
2 We solve the  p-Laplacian (nonlinear diffusion) combined with\n\
3 the Bratu (solid fuel ignition) nonlinearity in a 2D rectangular\n\
4 domain, using distributed arrays (DAs) to partition the parallel grid.\n\
5 The command line options include:\n\
6   -p <2>: `p' in p-Laplacian term\n\
7   -epsilon <1e-05>: Strain-regularization in p-Laplacian\n\
8   -lambda <6>: Bratu parameter\n\
9   -blocks <bx,by>: number of coefficient interfaces in x and y direction\n\
10   -kappa <1e-3>: diffusivity in odd regions\n\
11 \n";
12 
13 /*F
14     The $p$-Bratu problem is a combination of the $p$-Laplacian (nonlinear diffusion) and the Brutu solid fuel ignition problem.
15     This problem is modeled by the partial differential equation
16 
17 \begin{equation*}
18         -\nabla\cdot (\eta \nabla u) - \lambda \exp(u) = 0
19 \end{equation*}
20 
21     on $\Omega = (-1,1)^2$ with closure
22 
23 \begin{align*}
24         \eta(\gamma) &= (\epsilon^2 + \gamma)^{(p-2)/2} & \gamma &= \frac 1 2 |\nabla u|^2
25 \end{align*}
26 
27     and boundary conditions $u = 0$ for $(x,y) \in \partial \Omega$
28 
29     A 9-point finite difference stencil is used to discretize
30     the boundary value problem to obtain a nonlinear system of equations.
31     This would be a 5-point stencil if not for the $p$-Laplacian's nonlinearity.
32 F*/
33 
34 /*
35       mpiexec -n 2 ./ex15 -snes_monitor -ksp_monitor log_summary
36 */
37 
38 /*
39    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
40    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
41    file automatically includes:
42      petsc.h       - base PETSc routines   petscvec.h - vectors
43      petscsys.h    - system routines       petscmat.h - matrices
44      petscis.h     - index sets            petscksp.h - Krylov subspace methods
45      petscviewer.h - viewers               petscpc.h  - preconditioners
46      petscksp.h   - linear solvers
47 */
48 
49 #include <petscdm.h>
50 #include <petscdmda.h>
51 #include <petscsnes.h>
52 
53 typedef enum {JAC_BRATU,JAC_PICARD,JAC_STAR,JAC_NEWTON} JacType;
54 static const char *const JacTypes[] = {"BRATU","PICARD","STAR","NEWTON","JacType","JAC_",0};
55 
56 /*
57    User-defined application context - contains data needed by the
58    application-provided call-back routines, FormJacobianLocal() and
59    FormFunctionLocal().
60 */
61 typedef struct {
62   PetscReal   lambda;         /* Bratu parameter */
63   PetscReal   p;              /* Exponent in p-Laplacian */
64   PetscReal   epsilon;        /* Regularization */
65   PetscReal   source;         /* Source term */
66   JacType     jtype;          /* What type of Jacobian to assemble */
67   PetscBool   picard;
68   PetscInt    blocks[2];
69   PetscReal   kappa;
70   PetscInt    initial;        /* initial conditions type */
71 } AppCtx;
72 
73 /*
74    User-defined routines
75 */
76 static PetscErrorCode FormRHS(AppCtx*,DM,Vec);
77 static PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
78 static PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
79 static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
80 static PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,Mat,AppCtx*);
81 static PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);
82 
83 typedef struct _n_PreCheck *PreCheck;
84 struct _n_PreCheck {
85   MPI_Comm    comm;
86   PetscReal   angle;
87   Vec         Ylast;
88   PetscViewer monitor;
89 };
90 PetscErrorCode PreCheckCreate(MPI_Comm,PreCheck*);
91 PetscErrorCode PreCheckDestroy(PreCheck*);
92 PetscErrorCode PreCheckFunction(SNESLineSearch,Vec,Vec,PetscBool*,void*);
93 PetscErrorCode PreCheckSetFromOptions(PreCheck);
94 
95 int main(int argc,char **argv)
96 {
97   SNES                snes;                    /* nonlinear solver */
98   Vec                 x,r,b;                   /* solution, residual, rhs vectors */
99   AppCtx              user;                    /* user-defined work context */
100   PetscInt            its;                     /* iterations for convergence */
101   SNESConvergedReason reason;                  /* Check convergence */
102   PetscBool           alloc_star;              /* Only allocate for the STAR stencil  */
103   PetscBool           write_output;
104   char                filename[PETSC_MAX_PATH_LEN] = "ex15.vts";
105   PetscReal           bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
106   DM                  da;                      /* distributed array data structure */
107   PreCheck            precheck = NULL;         /* precheck context for version in this file */
108   PetscInt            use_precheck;            /* 0=none, 1=version in this file, 2=SNES-provided version */
109   PetscReal           precheck_angle;          /* When manually setting the SNES-provided precheck function */
110   SNESLineSearch      linesearch;
111 
112   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113      Initialize program
114      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115 
116   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
117 
118   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119      Initialize problem parameters
120   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121   user.lambda    = 0.0; user.p = 2.0; user.epsilon = 1e-5; user.source = 0.1; user.jtype = JAC_NEWTON;user.initial=-1;
122   user.blocks[0] = 1; user.blocks[1] = 1; user.kappa = 1e-3;
123   alloc_star     = PETSC_FALSE;
124   use_precheck   = 0; precheck_angle = 10.;
125   user.picard    = PETSC_FALSE;
126   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"p-Bratu options",__FILE__);
127   {
128     PetscInt two=2;
129     PetscCall(PetscOptionsReal("-lambda","Bratu parameter","",user.lambda,&user.lambda,NULL));
130     PetscCall(PetscOptionsReal("-p","Exponent `p' in p-Laplacian","",user.p,&user.p,NULL));
131     PetscCall(PetscOptionsReal("-epsilon","Strain-regularization in p-Laplacian","",user.epsilon,&user.epsilon,NULL));
132     PetscCall(PetscOptionsReal("-source","Constant source term","",user.source,&user.source,NULL));
133     PetscCall(PetscOptionsEnum("-jtype","Jacobian approximation to assemble","",JacTypes,(PetscEnum)user.jtype,(PetscEnum*)&user.jtype,NULL));
134     PetscCall(PetscOptionsName("-picard","Solve with defect-correction Picard iteration","",&user.picard));
135     if (user.picard) {
136       user.jtype = JAC_PICARD;
137       PetscCheck(user.p == 3,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Picard iteration is only supported for p == 3");
138       /* the Picard linearization only requires a 5 point stencil, while the Newton linearization requires a 9 point stencil */
139       /* hence allocating the 5 point stencil gives the same convergence as the 9 point stencil since the extra stencil points are not used */
140       PetscCall(PetscOptionsBool("-alloc_star","Allocate for STAR stencil (5-point)","",alloc_star,&alloc_star,NULL));
141     }
142     PetscCall(PetscOptionsInt("-precheck","Use a pre-check correction intended for use with Picard iteration 1=this version, 2=library","",use_precheck,&use_precheck,NULL));
143     PetscCall(PetscOptionsInt("-initial","Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)","",user.initial,&user.initial,NULL));
144     if (use_precheck == 2) {    /* Using library version, get the angle */
145       PetscCall(PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck_angle,&precheck_angle,NULL));
146     }
147     PetscCall(PetscOptionsIntArray("-blocks","number of coefficient interfaces in x and y direction","",user.blocks,&two,NULL));
148     PetscCall(PetscOptionsReal("-kappa","diffusivity in odd regions","",user.kappa,&user.kappa,NULL));
149     PetscCall(PetscOptionsString("-o","Output solution in vts format","",filename,filename,sizeof(filename),&write_output));
150   }
151   PetscOptionsEnd();
152   if (user.lambda > bratu_lambda_max || user.lambda < bratu_lambda_min) {
153     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING: lambda %g out of range for p=2\n",(double)user.lambda));
154   }
155 
156   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157      Create nonlinear solver context
158      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
160 
161   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162      Create distributed array (DMDA) to manage parallel grid and vectors
163   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,alloc_star ? DMDA_STENCIL_STAR : DMDA_STENCIL_BOX,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
165   PetscCall(DMSetFromOptions(da));
166   PetscCall(DMSetUp(da));
167 
168   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169      Extract global vectors from DM; then duplicate for remaining
170      vectors that are the same types
171    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172   PetscCall(DMCreateGlobalVector(da,&x));
173   PetscCall(VecDuplicate(x,&r));
174   PetscCall(VecDuplicate(x,&b));
175 
176   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177      User can override with:
178      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
179                 (unless user explicitly sets preconditioner)
180      -snes_mf_operator : form preconditioning matrix as set by the user,
181                          but use matrix-free approx for Jacobian-vector
182                          products within Newton-Krylov method
183 
184      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185 
186   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187      Set local function evaluation routine
188   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189   PetscCall(DMSetApplicationContext(da, &user));
190   PetscCall(SNESSetDM(snes,da));
191   if (user.picard) {
192     /*
193         This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access
194         the SNES to set it
195     */
196     PetscCall(DMDASNESSetPicardLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionPicardLocal,
197                                      (PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user));
198     PetscCall(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user));
199     PetscCall(SNESSetJacobian(snes,NULL,NULL,SNESPicardComputeJacobian,&user));
200   } else {
201     PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user));
202     PetscCall(DMDASNESSetJacobianLocal(da,(PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user));
203   }
204 
205   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206      Customize nonlinear solver; set runtime options
207    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208   PetscCall(SNESSetFromOptions(snes));
209   PetscCall(SNESSetNGS(snes,NonlinearGS,&user));
210   PetscCall(SNESGetLineSearch(snes, &linesearch));
211   /* Set up the precheck context if requested */
212   if (use_precheck == 1) {      /* Use the precheck routines in this file */
213     PetscCall(PreCheckCreate(PETSC_COMM_WORLD,&precheck));
214     PetscCall(PreCheckSetFromOptions(precheck));
215     PetscCall(SNESLineSearchSetPreCheck(linesearch,PreCheckFunction,precheck));
216   } else if (use_precheck == 2) { /* Use the version provided by the library */
217     PetscCall(SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&precheck_angle));
218   }
219 
220   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221      Evaluate initial guess
222      Note: The user should initialize the vector, x, with the initial guess
223      for the nonlinear solver prior to calling SNESSolve().  In particular,
224      to employ an initial guess of zero, the user should explicitly set
225      this vector to zero by calling VecSet().
226   */
227 
228   PetscCall(FormInitialGuess(&user,da,x));
229   PetscCall(FormRHS(&user,da,b));
230 
231   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232      Solve nonlinear system
233      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
234   PetscCall(SNESSolve(snes,b,x));
235   PetscCall(SNESGetIterationNumber(snes,&its));
236   PetscCall(SNESGetConvergedReason(snes,&reason));
237 
238   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s Number of nonlinear iterations = %" PetscInt_FMT "\n",SNESConvergedReasons[reason],its));
239 
240   if (write_output) {
241     PetscViewer viewer;
242     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer));
243     PetscCall(VecView(x,viewer));
244     PetscCall(PetscViewerDestroy(&viewer));
245   }
246 
247   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248      Free work space.  All PETSc objects should be destroyed when they
249      are no longer needed.
250    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
251 
252   PetscCall(VecDestroy(&x));
253   PetscCall(VecDestroy(&r));
254   PetscCall(VecDestroy(&b));
255   PetscCall(SNESDestroy(&snes));
256   PetscCall(DMDestroy(&da));
257   PetscCall(PreCheckDestroy(&precheck));
258   PetscCall(PetscFinalize());
259   return 0;
260 }
261 
262 /* ------------------------------------------------------------------- */
263 /*
264    FormInitialGuess - Forms initial approximation.
265 
266    Input Parameters:
267    user - user-defined application context
268    X - vector
269 
270    Output Parameter:
271    X - vector
272  */
273 static PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
274 {
275   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
276   PetscReal      temp1,temp,hx,hy;
277   PetscScalar    **x;
278 
279   PetscFunctionBeginUser;
280   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));
281 
282   hx    = 1.0/(PetscReal)(Mx-1);
283   hy    = 1.0/(PetscReal)(My-1);
284   temp1 = user->lambda / (user->lambda + 1.);
285 
286   /*
287      Get a pointer to vector data.
288        - For default PETSc vectors, VecGetArray() returns a pointer to
289          the data array.  Otherwise, the routine is implementation dependent.
290        - You MUST call VecRestoreArray() when you no longer need access to
291          the array.
292   */
293   PetscCall(DMDAVecGetArray(da,X,&x));
294 
295   /*
296      Get local grid boundaries (for 2-dimensional DA):
297        xs, ys   - starting grid indices (no ghost points)
298        xm, ym   - widths of local grid (no ghost points)
299 
300   */
301   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
302 
303   /*
304      Compute initial guess over the locally owned part of the grid
305   */
306   for (j=ys; j<ys+ym; j++) {
307     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
308     for (i=xs; i<xs+xm; i++) {
309       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
310         /* boundary conditions are all zero Dirichlet */
311         x[j][i] = 0.0;
312       } else {
313         if (user->initial == -1) {
314           if (user->lambda != 0) {
315             x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
316           } else {
317             /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting
318              * with an exact solution. */
319             const PetscReal
320               xx = 2*(PetscReal)i/(Mx-1) - 1,
321               yy = 2*(PetscReal)j/(My-1) - 1;
322             x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
323           }
324         } else if (user->initial == 0) {
325           x[j][i] = 0.;
326         } else if (user->initial == 1) {
327           const PetscReal
328             xx = 2*(PetscReal)i/(Mx-1) - 1,
329             yy = 2*(PetscReal)j/(My-1) - 1;
330           x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
331         } else {
332           if (user->lambda != 0) {
333             x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
334           } else {
335             x[j][i] = 0.5*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
336           }
337         }
338       }
339     }
340   }
341   /*
342      Restore vector
343   */
344   PetscCall(DMDAVecRestoreArray(da,X,&x));
345   PetscFunctionReturn(0);
346 }
347 
348 /*
349    FormRHS - Forms constant RHS for the problem.
350 
351    Input Parameters:
352    user - user-defined application context
353    B - RHS vector
354 
355    Output Parameter:
356    B - vector
357  */
358 static PetscErrorCode FormRHS(AppCtx *user,DM da,Vec B)
359 {
360   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
361   PetscReal      hx,hy;
362   PetscScalar    **b;
363 
364   PetscFunctionBeginUser;
365   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));
366 
367   hx    = 1.0/(PetscReal)(Mx-1);
368   hy    = 1.0/(PetscReal)(My-1);
369   PetscCall(DMDAVecGetArray(da,B,&b));
370   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
371   for (j=ys; j<ys+ym; j++) {
372     for (i=xs; i<xs+xm; i++) {
373       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
374         b[j][i] = 0.0;
375       } else {
376         b[j][i] = hx*hy*user->source;
377       }
378     }
379   }
380   PetscCall(DMDAVecRestoreArray(da,B,&b));
381   PetscFunctionReturn(0);
382 }
383 
384 static inline PetscReal kappa(const AppCtx *ctx,PetscReal x,PetscReal y)
385 {
386   return (((PetscInt)(x*ctx->blocks[0])) + ((PetscInt)(y*ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0;
387 }
388 /* p-Laplacian diffusivity */
389 static inline PetscScalar eta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
390 {
391   return kappa(ctx,x,y) * PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-2.));
392 }
393 static inline PetscScalar deta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
394 {
395   return (ctx->p == 2)
396          ? 0
397          : kappa(ctx,x,y)*PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-4)) * 0.5 * (ctx->p-2.);
398 }
399 
400 /* ------------------------------------------------------------------- */
401 /*
402    FormFunctionLocal - Evaluates nonlinear function, F(x).
403  */
404 static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
405 {
406   PetscReal      hx,hy,dhx,dhy,sc;
407   PetscInt       i,j;
408   PetscScalar    eu;
409 
410   PetscFunctionBeginUser;
411   hx     = 1.0/(PetscReal)(info->mx-1);
412   hy     = 1.0/(PetscReal)(info->my-1);
413   sc     = hx*hy*user->lambda;
414   dhx    = 1/hx;
415   dhy    = 1/hy;
416   /*
417      Compute function over the locally owned part of the grid
418   */
419   for (j=info->ys; j<info->ys+info->ym; j++) {
420     for (i=info->xs; i<info->xs+info->xm; i++) {
421       PetscReal xx = i*hx,yy = j*hy;
422       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
423         f[j][i] = x[j][i];
424       } else {
425         const PetscScalar
426           u    = x[j][i],
427           ux_E = dhx*(x[j][i+1]-x[j][i]),
428           uy_E = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
429           ux_W = dhx*(x[j][i]-x[j][i-1]),
430           uy_W = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
431           ux_N = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
432           uy_N = dhy*(x[j+1][i]-x[j][i]),
433           ux_S = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
434           uy_S = dhy*(x[j][i]-x[j-1][i]),
435           e_E  = eta(user,xx,yy,ux_E,uy_E),
436           e_W  = eta(user,xx,yy,ux_W,uy_W),
437           e_N  = eta(user,xx,yy,ux_N,uy_N),
438           e_S  = eta(user,xx,yy,ux_S,uy_S),
439           uxx  = -hy * (e_E*ux_E - e_W*ux_W),
440           uyy  = -hx * (e_N*uy_N - e_S*uy_S);
441         if (sc) eu = PetscExpScalar(u);
442         else    eu = 0.;
443         /* For p=2, these terms decay to:
444          uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx
445          uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy
446         */
447         f[j][i] = uxx + uyy - sc*eu;
448       }
449     }
450   }
451   PetscCall(PetscLogFlops(info->xm*info->ym*(72.0)));
452   PetscFunctionReturn(0);
453 }
454 
455 /*
456     This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation,
457     that is FormFunction applies A(x) x - b(x) while this applies b(x) because for Picard we think of it as solving A(x) x = b(x)
458 
459 */
460 static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
461 {
462   PetscReal hx,hy,sc;
463   PetscInt  i,j;
464 
465   PetscFunctionBeginUser;
466   hx     = 1.0/(PetscReal)(info->mx-1);
467   hy     = 1.0/(PetscReal)(info->my-1);
468   sc     = hx*hy*user->lambda;
469   /*
470      Compute function over the locally owned part of the grid
471   */
472   for (j=info->ys; j<info->ys+info->ym; j++) {
473     for (i=info->xs; i<info->xs+info->xm; i++) {
474       if (!(i == 0 || j == 0 || i == info->mx-1 || j == info->my-1)) {
475         const PetscScalar u = x[j][i];
476         f[j][i] = sc*PetscExpScalar(u);
477       } else {
478         f[j][i] = 0.0; /* this is zero because the A(x) x term forces the x to be zero on the boundary */
479       }
480     }
481   }
482   PetscCall(PetscLogFlops(info->xm*info->ym*2.0));
483   PetscFunctionReturn(0);
484 }
485 
486 /*
487    FormJacobianLocal - Evaluates Jacobian matrix.
488 */
489 static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat J,Mat B,AppCtx *user)
490 {
491   PetscInt       i,j;
492   MatStencil     col[9],row;
493   PetscScalar    v[9];
494   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;
495 
496   PetscFunctionBeginUser;
497   PetscCall(MatZeroEntries(B));
498   hx    = 1.0/(PetscReal)(info->mx-1);
499   hy    = 1.0/(PetscReal)(info->my-1);
500   sc    = hx*hy*user->lambda;
501   dhx   = 1/hx;
502   dhy   = 1/hy;
503   hxdhy = hx/hy;
504   hydhx = hy/hx;
505 
506   /*
507      Compute entries for the locally owned part of the Jacobian.
508       - PETSc parallel matrix formats are partitioned by
509         contiguous chunks of rows across the processors.
510       - Each processor needs to insert only elements that it owns
511         locally (but any non-local elements will be sent to the
512         appropriate processor during matrix assembly).
513       - Here, we set all entries for a particular row at once.
514   */
515   for (j=info->ys; j<info->ys+info->ym; j++) {
516     for (i=info->xs; i<info->xs+info->xm; i++) {
517       PetscReal xx = i*hx,yy = j*hy;
518       row.j = j; row.i = i;
519       /* boundary points */
520       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
521         v[0] = 1.0;
522         PetscCall(MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES));
523       } else {
524         /* interior grid points */
525         const PetscScalar
526           ux_E     = dhx*(x[j][i+1]-x[j][i]),
527           uy_E     = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
528           ux_W     = dhx*(x[j][i]-x[j][i-1]),
529           uy_W     = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
530           ux_N     = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
531           uy_N     = dhy*(x[j+1][i]-x[j][i]),
532           ux_S     = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
533           uy_S     = dhy*(x[j][i]-x[j-1][i]),
534           u        = x[j][i],
535           e_E      = eta(user,xx,yy,ux_E,uy_E),
536           e_W      = eta(user,xx,yy,ux_W,uy_W),
537           e_N      = eta(user,xx,yy,ux_N,uy_N),
538           e_S      = eta(user,xx,yy,ux_S,uy_S),
539           de_E     = deta(user,xx,yy,ux_E,uy_E),
540           de_W     = deta(user,xx,yy,ux_W,uy_W),
541           de_N     = deta(user,xx,yy,ux_N,uy_N),
542           de_S     = deta(user,xx,yy,ux_S,uy_S),
543           skew_E   = de_E*ux_E*uy_E,
544           skew_W   = de_W*ux_W*uy_W,
545           skew_N   = de_N*ux_N*uy_N,
546           skew_S   = de_S*ux_S*uy_S,
547           cross_EW = 0.25*(skew_E - skew_W),
548           cross_NS = 0.25*(skew_N - skew_S),
549           newt_E   = e_E+de_E*PetscSqr(ux_E),
550           newt_W   = e_W+de_W*PetscSqr(ux_W),
551           newt_N   = e_N+de_N*PetscSqr(uy_N),
552           newt_S   = e_S+de_S*PetscSqr(uy_S);
553         /* interior grid points */
554         switch (user->jtype) {
555         case JAC_BRATU:
556           /* Jacobian from p=2 */
557           v[0] = -hxdhy;                                           col[0].j = j-1;   col[0].i = i;
558           v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
559           v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(u);       col[2].j = row.j; col[2].i = row.i;
560           v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
561           v[4] = -hxdhy;                                           col[4].j = j+1;   col[4].i = i;
562           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
563           break;
564         case JAC_PICARD:
565           /* Jacobian arising from Picard linearization */
566           v[0] = -hxdhy*e_S;                                           col[0].j = j-1;   col[0].i = i;
567           v[1] = -hydhx*e_W;                                           col[1].j = j;     col[1].i = i-1;
568           v[2] = (e_W+e_E)*hydhx + (e_S+e_N)*hxdhy;                    col[2].j = row.j; col[2].i = row.i;
569           v[3] = -hydhx*e_E;                                           col[3].j = j;     col[3].i = i+1;
570           v[4] = -hxdhy*e_N;                                           col[4].j = j+1;   col[4].i = i;
571           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
572           break;
573         case JAC_STAR:
574           /* Full Jacobian, but only a star stencil */
575           col[0].j = j-1; col[0].i = i;
576           col[1].j = j;   col[1].i = i-1;
577           col[2].j = j;   col[2].i = i;
578           col[3].j = j;   col[3].i = i+1;
579           col[4].j = j+1; col[4].i = i;
580           v[0]     = -hxdhy*newt_S + cross_EW;
581           v[1]     = -hydhx*newt_W + cross_NS;
582           v[2]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
583           v[3]     = -hydhx*newt_E - cross_NS;
584           v[4]     = -hxdhy*newt_N - cross_EW;
585           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
586           break;
587         case JAC_NEWTON:
588           /* The Jacobian is
589 
590            -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u
591 
592           */
593           col[0].j = j-1; col[0].i = i-1;
594           col[1].j = j-1; col[1].i = i;
595           col[2].j = j-1; col[2].i = i+1;
596           col[3].j = j;   col[3].i = i-1;
597           col[4].j = j;   col[4].i = i;
598           col[5].j = j;   col[5].i = i+1;
599           col[6].j = j+1; col[6].i = i-1;
600           col[7].j = j+1; col[7].i = i;
601           col[8].j = j+1; col[8].i = i+1;
602           v[0]     = -0.25*(skew_S + skew_W);
603           v[1]     = -hxdhy*newt_S + cross_EW;
604           v[2]     =  0.25*(skew_S + skew_E);
605           v[3]     = -hydhx*newt_W + cross_NS;
606           v[4]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
607           v[5]     = -hydhx*newt_E - cross_NS;
608           v[6]     =  0.25*(skew_N + skew_W);
609           v[7]     = -hxdhy*newt_N - cross_EW;
610           v[8]     = -0.25*(skew_N + skew_E);
611           PetscCall(MatSetValuesStencil(B,1,&row,9,col,v,INSERT_VALUES));
612           break;
613         default:
614           SETERRQ(PetscObjectComm((PetscObject)info->da),PETSC_ERR_SUP,"Jacobian type %d not implemented",user->jtype);
615         }
616       }
617     }
618   }
619 
620   /*
621      Assemble matrix, using the 2-step process:
622        MatAssemblyBegin(), MatAssemblyEnd().
623   */
624   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
625   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
626 
627   if (J != B) {
628     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
629     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
630   }
631 
632   /*
633      Tell the matrix we will never add a new nonzero location to the
634      matrix. If we do, it will generate an error.
635   */
636   if (user->jtype == JAC_NEWTON) {
637     PetscCall(PetscLogFlops(info->xm*info->ym*(131.0)));
638   }
639   PetscFunctionReturn(0);
640 }
641 
642 /***********************************************************
643  * PreCheck implementation
644  ***********************************************************/
645 PetscErrorCode PreCheckSetFromOptions(PreCheck precheck)
646 {
647   PetscBool      flg;
648 
649   PetscFunctionBeginUser;
650   PetscOptionsBegin(precheck->comm,NULL,"PreCheck Options","none");
651   PetscCall(PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck->angle,&precheck->angle,NULL));
652   flg  = PETSC_FALSE;
653   PetscCall(PetscOptionsBool("-precheck_monitor","Monitor choices made by precheck routine","",flg,&flg,NULL));
654   if (flg) {
655     PetscCall(PetscViewerASCIIOpen(precheck->comm,"stdout",&precheck->monitor));
656   }
657   PetscOptionsEnd();
658   PetscFunctionReturn(0);
659 }
660 
661 /*
662   Compare the direction of the current and previous step, modify the current step accordingly
663 */
664 PetscErrorCode PreCheckFunction(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed, void *ctx)
665 {
666   PreCheck       precheck;
667   Vec            Ylast;
668   PetscScalar    dot;
669   PetscInt       iter;
670   PetscReal      ynorm,ylastnorm,theta,angle_radians;
671   SNES           snes;
672 
673   PetscFunctionBeginUser;
674   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
675   precheck = (PreCheck)ctx;
676   if (!precheck->Ylast) PetscCall(VecDuplicate(Y,&precheck->Ylast));
677   Ylast = precheck->Ylast;
678   PetscCall(SNESGetIterationNumber(snes,&iter));
679   if (iter < 1) {
680     PetscCall(VecCopy(Y,Ylast));
681     *changed = PETSC_FALSE;
682     PetscFunctionReturn(0);
683   }
684 
685   PetscCall(VecDot(Y,Ylast,&dot));
686   PetscCall(VecNorm(Y,NORM_2,&ynorm));
687   PetscCall(VecNorm(Ylast,NORM_2,&ylastnorm));
688   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
689   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
690   angle_radians = precheck->angle * PETSC_PI / 180.;
691   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
692     /* Modify the step Y */
693     PetscReal alpha,ydiffnorm;
694     PetscCall(VecAXPY(Ylast,-1.0,Y));
695     PetscCall(VecNorm(Ylast,NORM_2,&ydiffnorm));
696     alpha = ylastnorm / ydiffnorm;
697     PetscCall(VecCopy(Y,Ylast));
698     PetscCall(VecScale(Y,alpha));
699     if (precheck->monitor) {
700       PetscCall(PetscViewerASCIIPrintf(precheck->monitor,"Angle %g degrees less than threshold %g, corrected step by alpha=%g\n",(double)(theta*180./PETSC_PI),(double)precheck->angle,(double)alpha));
701     }
702   } else {
703     PetscCall(VecCopy(Y,Ylast));
704     *changed = PETSC_FALSE;
705     if (precheck->monitor) {
706       PetscCall(PetscViewerASCIIPrintf(precheck->monitor,"Angle %g degrees exceeds threshold %g, no correction applied\n",(double)(theta*180./PETSC_PI),(double)precheck->angle));
707     }
708   }
709   PetscFunctionReturn(0);
710 }
711 
712 PetscErrorCode PreCheckDestroy(PreCheck *precheck)
713 {
714   PetscFunctionBeginUser;
715   if (!*precheck) PetscFunctionReturn(0);
716   PetscCall(VecDestroy(&(*precheck)->Ylast));
717   PetscCall(PetscViewerDestroy(&(*precheck)->monitor));
718   PetscCall(PetscFree(*precheck));
719   PetscFunctionReturn(0);
720 }
721 
722 PetscErrorCode PreCheckCreate(MPI_Comm comm,PreCheck *precheck)
723 {
724   PetscFunctionBeginUser;
725   PetscCall(PetscNew(precheck));
726 
727   (*precheck)->comm  = comm;
728   (*precheck)->angle = 10.;     /* only active if angle is less than 10 degrees */
729   PetscFunctionReturn(0);
730 }
731 
732 /*
733       Applies some sweeps on nonlinear Gauss-Seidel on each process
734 
735  */
736 PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx)
737 {
738   PetscInt       i,j,k,xs,ys,xm,ym,its,tot_its,sweeps,l,m;
739   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;
740   PetscScalar    **x,**b,bij,F,F0=0,J,y,u,eu;
741   PetscReal      atol,rtol,stol;
742   DM             da;
743   AppCtx         *user = (AppCtx*)ctx;
744   Vec            localX,localB;
745   DMDALocalInfo  info;
746 
747   PetscFunctionBeginUser;
748   PetscCall(SNESGetDM(snes,&da));
749   PetscCall(DMDAGetLocalInfo(da,&info));
750 
751   hx     = 1.0/(PetscReal)(info.mx-1);
752   hy     = 1.0/(PetscReal)(info.my-1);
753   sc     = hx*hy*user->lambda;
754   dhx    = 1/hx;
755   dhy    = 1/hy;
756   hxdhy  = hx/hy;
757   hydhx  = hy/hx;
758 
759   tot_its = 0;
760   PetscCall(SNESNGSGetSweeps(snes,&sweeps));
761   PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its));
762   PetscCall(DMGetLocalVector(da,&localX));
763   if (B) {
764     PetscCall(DMGetLocalVector(da,&localB));
765   }
766   if (B) {
767     PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB));
768     PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB));
769   }
770   if (B) PetscCall(DMDAVecGetArrayRead(da,localB,&b));
771   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
772   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
773   PetscCall(DMDAVecGetArray(da,localX,&x));
774   for (l=0; l<sweeps; l++) {
775     /*
776      Get local grid boundaries (for 2-dimensional DMDA):
777      xs, ys   - starting grid indices (no ghost points)
778      xm, ym   - widths of local grid (no ghost points)
779      */
780     PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
781     for (m=0; m<2; m++) {
782       for (j=ys; j<ys+ym; j++) {
783         for (i=xs+(m+j)%2; i<xs+xm; i+=2) {
784           PetscReal xx = i*hx,yy = j*hy;
785           if (B) bij = b[j][i];
786           else   bij = 0.;
787 
788           if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
789             /* boundary conditions are all zero Dirichlet */
790             x[j][i] = 0.0 + bij;
791           } else {
792             const PetscScalar
793               u_E = x[j][i+1],
794               u_W = x[j][i-1],
795               u_N = x[j+1][i],
796               u_S = x[j-1][i];
797             const PetscScalar
798               uy_E   = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
799               uy_W   = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
800               ux_N   = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
801               ux_S   = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]);
802             u = x[j][i];
803             for (k=0; k<its; k++) {
804               const PetscScalar
805                 ux_E   = dhx*(u_E-u),
806                 ux_W   = dhx*(u-u_W),
807                 uy_N   = dhy*(u_N-u),
808                 uy_S   = dhy*(u-u_S),
809                 e_E    = eta(user,xx,yy,ux_E,uy_E),
810                 e_W    = eta(user,xx,yy,ux_W,uy_W),
811                 e_N    = eta(user,xx,yy,ux_N,uy_N),
812                 e_S    = eta(user,xx,yy,ux_S,uy_S),
813                 de_E   = deta(user,xx,yy,ux_E,uy_E),
814                 de_W   = deta(user,xx,yy,ux_W,uy_W),
815                 de_N   = deta(user,xx,yy,ux_N,uy_N),
816                 de_S   = deta(user,xx,yy,ux_S,uy_S),
817                 newt_E = e_E+de_E*PetscSqr(ux_E),
818                 newt_W = e_W+de_W*PetscSqr(ux_W),
819                 newt_N = e_N+de_N*PetscSqr(uy_N),
820                 newt_S = e_S+de_S*PetscSqr(uy_S),
821                 uxx    = -hy * (e_E*ux_E - e_W*ux_W),
822                 uyy    = -hx * (e_N*uy_N - e_S*uy_S);
823 
824               if (sc) eu = PetscExpScalar(u);
825               else    eu = 0;
826 
827               F = uxx + uyy - sc*eu - bij;
828               if (k == 0) F0 = F;
829               J  = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*eu;
830               y  = F/J;
831               u -= y;
832               tot_its++;
833               if (atol > PetscAbsReal(PetscRealPart(F)) ||
834                   rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) ||
835                   stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) {
836                 break;
837               }
838             }
839             x[j][i] = u;
840           }
841         }
842       }
843     }
844     /*
845 x     Restore vector
846      */
847   }
848   PetscCall(DMDAVecRestoreArray(da,localX,&x));
849   PetscCall(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X));
850   PetscCall(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X));
851   PetscCall(PetscLogFlops(tot_its*(118.0)));
852   PetscCall(DMRestoreLocalVector(da,&localX));
853   if (B) {
854     PetscCall(DMDAVecRestoreArrayRead(da,localB,&b));
855     PetscCall(DMRestoreLocalVector(da,&localB));
856   }
857   PetscFunctionReturn(0);
858 }
859 
860 /*TEST
861 
862    test:
863       nsize: 2
864       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype NEWTON
865       requires: !single
866 
867    test:
868       suffix: 2
869       nsize: 2
870       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -precheck 1
871       requires: !single
872 
873    test:
874       suffix: 3
875       nsize: 2
876       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -picard -precheck 1 -p 3
877       requires: !single
878 
879    test:
880       suffix: 3_star
881       nsize: 2
882       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -picard -precheck 1 -p 3 -alloc_star
883       output_file: output/ex15_3.out
884       requires: !single
885 
886    test:
887       suffix: 4
888       args: -snes_monitor_short -snes_type newtonls -npc_snes_type ngs -snes_npc_side left -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -ksp_monitor_short -pc_type none
889       requires: !single
890 
891    test:
892       suffix: lag_jac
893       nsize: 4
894       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 6.0 -lambda 0 -jtype NEWTON -snes_type ngmres -npc_snes_type newtonls -npc_snes_lag_jacobian 5 -npc_pc_type asm -npc_ksp_converged_reason -npc_snes_lag_jacobian_persists
895       requires: !single
896 
897    test:
898       suffix: lag_pc
899       nsize: 4
900       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 6.0 -lambda 0 -jtype NEWTON -snes_type ngmres -npc_snes_type newtonls -npc_snes_lag_preconditioner 5 -npc_pc_type asm -npc_ksp_converged_reason -npc_snes_lag_preconditioner_persists
901       requires: !single
902 
903    test:
904       suffix: nleqerr
905       args: -snes_monitor_short -snes_type newtonls -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -snes_linesearch_monitor -pc_type lu -snes_linesearch_type nleqerr
906       requires: !single
907 
908    test:
909       suffix: mf
910       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator
911       requires: !single
912 
913    test:
914       suffix: mf_picard
915       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator -picard
916       requires: !single
917       output_file: output/ex15_mf.out
918 
919    test:
920       suffix: fd_picard
921       args: -snes_monitor_short -pc_type lu -da_refine 2  -p 3 -ksp_rtol 1.e-12  -snes_fd -picard
922       requires: !single
923 
924    test:
925       suffix: fd_color_picard
926       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_fd_color -picard
927       requires: !single
928       output_file: output/ex15_mf.out
929 
930 TEST*/
931