xref: /petsc/src/snes/tutorials/ex15.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static const char help[] = "p-Bratu nonlinear PDE in 2d.\n\
2c4762a1bSJed Brown We solve the  p-Laplacian (nonlinear diffusion) combined with\n\
3c4762a1bSJed Brown the Bratu (solid fuel ignition) nonlinearity in a 2D rectangular\n\
4c4762a1bSJed Brown domain, using distributed arrays (DAs) to partition the parallel grid.\n\
5c4762a1bSJed Brown The command line options include:\n\
6c4762a1bSJed Brown   -p <2>: `p' in p-Laplacian term\n\
7c4762a1bSJed Brown   -epsilon <1e-05>: Strain-regularization in p-Laplacian\n\
8c4762a1bSJed Brown   -lambda <6>: Bratu parameter\n\
9c4762a1bSJed Brown   -blocks <bx,by>: number of coefficient interfaces in x and y direction\n\
10c4762a1bSJed Brown   -kappa <1e-3>: diffusivity in odd regions\n\
11c4762a1bSJed Brown \n";
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*F
14c4762a1bSJed Brown     The $p$-Bratu problem is a combination of the $p$-Laplacian (nonlinear diffusion) and the Brutu solid fuel ignition problem.
15c4762a1bSJed Brown     This problem is modeled by the partial differential equation
16c4762a1bSJed Brown 
17c4762a1bSJed Brown \begin{equation*}
18c4762a1bSJed Brown         -\nabla\cdot (\eta \nabla u) - \lambda \exp(u) = 0
19c4762a1bSJed Brown \end{equation*}
20c4762a1bSJed Brown 
21c4762a1bSJed Brown     on $\Omega = (-1,1)^2$ with closure
22c4762a1bSJed Brown 
23c4762a1bSJed Brown \begin{align*}
24c4762a1bSJed Brown         \eta(\gamma) &= (\epsilon^2 + \gamma)^{(p-2)/2} & \gamma &= \frac 1 2 |\nabla u|^2
25c4762a1bSJed Brown \end{align*}
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     and boundary conditions $u = 0$ for $(x,y) \in \partial \Omega$
28c4762a1bSJed Brown 
29c4762a1bSJed Brown     A 9-point finite difference stencil is used to discretize
30c4762a1bSJed Brown     the boundary value problem to obtain a nonlinear system of equations.
31c4762a1bSJed Brown     This would be a 5-point stencil if not for the $p$-Laplacian's nonlinearity.
32c4762a1bSJed Brown F*/
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /*
35c4762a1bSJed Brown       mpiexec -n 2 ./ex15 -snes_monitor -ksp_monitor log_summary
36c4762a1bSJed Brown */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown /*
39c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
40c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
41c4762a1bSJed Brown    file automatically includes:
42c4762a1bSJed Brown      petsc.h       - base PETSc routines   petscvec.h - vectors
43c4762a1bSJed Brown      petscsys.h    - system routines       petscmat.h - matrices
44c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
45c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
46c4762a1bSJed Brown      petscksp.h   - linear solvers
47c4762a1bSJed Brown */
48c4762a1bSJed Brown 
49c4762a1bSJed Brown #include <petscdm.h>
50c4762a1bSJed Brown #include <petscdmda.h>
51c4762a1bSJed Brown #include <petscsnes.h>
52c4762a1bSJed Brown 
539371c9d4SSatish Balay typedef enum {
549371c9d4SSatish Balay   JAC_BRATU,
559371c9d4SSatish Balay   JAC_PICARD,
569371c9d4SSatish Balay   JAC_STAR,
579371c9d4SSatish Balay   JAC_NEWTON
589371c9d4SSatish Balay } JacType;
59c4762a1bSJed Brown static const char *const JacTypes[] = {"BRATU", "PICARD", "STAR", "NEWTON", "JacType", "JAC_", 0};
60c4762a1bSJed Brown 
61c4762a1bSJed Brown /*
62c4762a1bSJed Brown    User-defined application context - contains data needed by the
63c4762a1bSJed Brown    application-provided call-back routines, FormJacobianLocal() and
64c4762a1bSJed Brown    FormFunctionLocal().
65c4762a1bSJed Brown */
66c4762a1bSJed Brown typedef struct {
67c4762a1bSJed Brown   PetscReal lambda;  /* Bratu parameter */
68c4762a1bSJed Brown   PetscReal p;       /* Exponent in p-Laplacian */
69c4762a1bSJed Brown   PetscReal epsilon; /* Regularization */
70c4762a1bSJed Brown   PetscReal source;  /* Source term */
71c4762a1bSJed Brown   JacType   jtype;   /* What type of Jacobian to assemble */
72c4762a1bSJed Brown   PetscBool picard;
73c4762a1bSJed Brown   PetscInt  blocks[2];
74c4762a1bSJed Brown   PetscReal kappa;
75c4762a1bSJed Brown   PetscInt  initial; /* initial conditions type */
76c4762a1bSJed Brown } AppCtx;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown /*
79c4762a1bSJed Brown    User-defined routines
80c4762a1bSJed Brown */
81c4762a1bSJed Brown static PetscErrorCode FormRHS(AppCtx *, DM, Vec);
82c4762a1bSJed Brown static PetscErrorCode FormInitialGuess(AppCtx *, DM, Vec);
83c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal(DMDALocalInfo *, PetscScalar **, PetscScalar **, AppCtx *);
84c4762a1bSJed Brown static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *, PetscScalar **, PetscScalar **, AppCtx *);
85c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal(DMDALocalInfo *, PetscScalar **, Mat, Mat, AppCtx *);
86c4762a1bSJed Brown static PetscErrorCode NonlinearGS(SNES, Vec, Vec, void *);
87c4762a1bSJed Brown 
88c4762a1bSJed Brown typedef struct _n_PreCheck *PreCheck;
89c4762a1bSJed Brown struct _n_PreCheck {
90c4762a1bSJed Brown   MPI_Comm    comm;
91c4762a1bSJed Brown   PetscReal   angle;
92c4762a1bSJed Brown   Vec         Ylast;
93c4762a1bSJed Brown   PetscViewer monitor;
94c4762a1bSJed Brown };
95c4762a1bSJed Brown PetscErrorCode PreCheckCreate(MPI_Comm, PreCheck *);
96c4762a1bSJed Brown PetscErrorCode PreCheckDestroy(PreCheck *);
97c4762a1bSJed Brown PetscErrorCode PreCheckFunction(SNESLineSearch, Vec, Vec, PetscBool *, void *);
98c4762a1bSJed Brown PetscErrorCode PreCheckSetFromOptions(PreCheck);
99c4762a1bSJed Brown 
1009371c9d4SSatish Balay int main(int argc, char **argv) {
101c4762a1bSJed Brown   SNES                snes;       /* nonlinear solver */
102c4762a1bSJed Brown   Vec                 x, r, b;    /* solution, residual, rhs vectors */
103c4762a1bSJed Brown   AppCtx              user;       /* user-defined work context */
104c4762a1bSJed Brown   PetscInt            its;        /* iterations for convergence */
105c4762a1bSJed Brown   SNESConvergedReason reason;     /* Check convergence */
106c4762a1bSJed Brown   PetscBool           alloc_star; /* Only allocate for the STAR stencil  */
107c4762a1bSJed Brown   PetscBool           write_output;
108c4762a1bSJed Brown   char                filename[PETSC_MAX_PATH_LEN] = "ex15.vts";
109c4762a1bSJed Brown   PetscReal           bratu_lambda_max = 6.81, bratu_lambda_min = 0.;
11085b95db3SBarry Smith   DM                  da;              /* distributed array data structure */
111c4762a1bSJed Brown   PreCheck            precheck = NULL; /* precheck context for version in this file */
112c4762a1bSJed Brown   PetscInt            use_precheck;    /* 0=none, 1=version in this file, 2=SNES-provided version */
113c4762a1bSJed Brown   PetscReal           precheck_angle;  /* When manually setting the SNES-provided precheck function */
114c4762a1bSJed Brown   SNESLineSearch      linesearch;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Initialize program
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4762a1bSJed Brown 
120327415f7SBarry Smith   PetscFunctionBeginUser;
1219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown      Initialize problem parameters
125c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1269371c9d4SSatish Balay   user.lambda    = 0.0;
1279371c9d4SSatish Balay   user.p         = 2.0;
1289371c9d4SSatish Balay   user.epsilon   = 1e-5;
1299371c9d4SSatish Balay   user.source    = 0.1;
1309371c9d4SSatish Balay   user.jtype     = JAC_NEWTON;
1319371c9d4SSatish Balay   user.initial   = -1;
1329371c9d4SSatish Balay   user.blocks[0] = 1;
1339371c9d4SSatish Balay   user.blocks[1] = 1;
1349371c9d4SSatish Balay   user.kappa     = 1e-3;
135c4762a1bSJed Brown   alloc_star     = PETSC_FALSE;
1369371c9d4SSatish Balay   use_precheck   = 0;
1379371c9d4SSatish Balay   precheck_angle = 10.;
138c4762a1bSJed Brown   user.picard    = PETSC_FALSE;
139d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "p-Bratu options", __FILE__);
140c4762a1bSJed Brown   {
141c4762a1bSJed Brown     PetscInt two = 2;
1429566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda", "Bratu parameter", "", user.lambda, &user.lambda, NULL));
1439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-p", "Exponent `p' in p-Laplacian", "", user.p, &user.p, NULL));
1449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-epsilon", "Strain-regularization in p-Laplacian", "", user.epsilon, &user.epsilon, NULL));
1459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-source", "Constant source term", "", user.source, &user.source, NULL));
1469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-jtype", "Jacobian approximation to assemble", "", JacTypes, (PetscEnum)user.jtype, (PetscEnum *)&user.jtype, NULL));
1479566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-picard", "Solve with defect-correction Picard iteration", "", &user.picard));
14885b95db3SBarry Smith     if (user.picard) {
14985b95db3SBarry Smith       user.jtype = JAC_PICARD;
150e00437b9SBarry Smith       PetscCheck(user.p == 3, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Picard iteration is only supported for p == 3");
15185b95db3SBarry Smith       /* the Picard linearization only requires a 5 point stencil, while the Newton linearization requires a 9 point stencil */
15285b95db3SBarry Smith       /* hence allocating the 5 point stencil gives the same convergence as the 9 point stencil since the extra stencil points are not used */
1539566063dSJacob Faibussowitsch       PetscCall(PetscOptionsBool("-alloc_star", "Allocate for STAR stencil (5-point)", "", alloc_star, &alloc_star, NULL));
15485b95db3SBarry Smith     }
1559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-precheck", "Use a pre-check correction intended for use with Picard iteration 1=this version, 2=library", "", use_precheck, &use_precheck, NULL));
1569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-initial", "Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)", "", user.initial, &user.initial, NULL));
157c4762a1bSJed Brown     if (use_precheck == 2) { /* Using library version, get the angle */
1589566063dSJacob Faibussowitsch       PetscCall(PetscOptionsReal("-precheck_angle", "Angle in degrees between successive search directions necessary to activate step correction", "", precheck_angle, &precheck_angle, NULL));
159c4762a1bSJed Brown     }
1609566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-blocks", "number of coefficient interfaces in x and y direction", "", user.blocks, &two, NULL));
1619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-kappa", "diffusivity in odd regions", "", user.kappa, &user.kappa, NULL));
1629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-o", "Output solution in vts format", "", filename, filename, sizeof(filename), &write_output));
163c4762a1bSJed Brown   }
164d0609cedSBarry Smith   PetscOptionsEnd();
165*48a46eb9SPierre Jolivet   if (user.lambda > bratu_lambda_max || user.lambda < bratu_lambda_min) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING: lambda %g out of range for p=2\n", (double)user.lambda));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168c4762a1bSJed Brown      Create nonlinear solver context
169c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1709566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
174c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1759566063dSJacob Faibussowitsch   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));
1769566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1779566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown      Extract global vectors from DM; then duplicate for remaining
181c4762a1bSJed Brown      vectors that are the same types
182c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1839566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
1849566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
1859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &b));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18885b95db3SBarry Smith      User can override with:
189c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
190c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
191c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
192c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
193c4762a1bSJed Brown                          products within Newton-Krylov method
194c4762a1bSJed Brown 
195c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown      Set local function evaluation routine
199c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2009566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
2019566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
202c4762a1bSJed Brown   if (user.picard) {
203c4762a1bSJed Brown     /*
204c4762a1bSJed Brown         This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access
205c4762a1bSJed Brown         the SNES to set it
206c4762a1bSJed Brown     */
2079371c9d4SSatish Balay     PetscCall(DMDASNESSetPicardLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionPicardLocal, (PetscErrorCode(*)(DMDALocalInfo *, void *, Mat, Mat, void *))FormJacobianLocal, &user));
2089566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, SNESPicardComputeFunction, &user));
2099566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESPicardComputeJacobian, &user));
210c4762a1bSJed Brown   } else {
2119566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user));
2129566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetJacobianLocal(da, (PetscErrorCode(*)(DMDALocalInfo *, void *, Mat, Mat, void *))FormJacobianLocal, &user));
213c4762a1bSJed Brown   }
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
217c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2189566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
2199566063dSJacob Faibussowitsch   PetscCall(SNESSetNGS(snes, NonlinearGS, &user));
2209566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
221c4762a1bSJed Brown   /* Set up the precheck context if requested */
222c4762a1bSJed Brown   if (use_precheck == 1) { /* Use the precheck routines in this file */
2239566063dSJacob Faibussowitsch     PetscCall(PreCheckCreate(PETSC_COMM_WORLD, &precheck));
2249566063dSJacob Faibussowitsch     PetscCall(PreCheckSetFromOptions(precheck));
2259566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch, PreCheckFunction, precheck));
226c4762a1bSJed Brown   } else if (use_precheck == 2) { /* Use the version provided by the library */
2279566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &precheck_angle));
228c4762a1bSJed Brown   }
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231c4762a1bSJed Brown      Evaluate initial guess
232c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
233c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
234c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
235c4762a1bSJed Brown      this vector to zero by calling VecSet().
236c4762a1bSJed Brown   */
237c4762a1bSJed Brown 
2389566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user, da, x));
2399566063dSJacob Faibussowitsch   PetscCall(FormRHS(&user, da, b));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242c4762a1bSJed Brown      Solve nonlinear system
243c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2449566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, b, x));
2459566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
2469566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes, &reason));
247c4762a1bSJed Brown 
24863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s Number of nonlinear iterations = %" PetscInt_FMT "\n", SNESConvergedReasons[reason], its));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   if (write_output) {
251c4762a1bSJed Brown     PetscViewer viewer;
2529566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD, filename, FILE_MODE_WRITE, &viewer));
2539566063dSJacob Faibussowitsch     PetscCall(VecView(x, viewer));
2549566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
259c4762a1bSJed Brown      are no longer needed.
260c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
261c4762a1bSJed Brown 
2629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
2659566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2669566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2679566063dSJacob Faibussowitsch   PetscCall(PreCheckDestroy(&precheck));
2689566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
269b122ec5aSJacob Faibussowitsch   return 0;
270c4762a1bSJed Brown }
271c4762a1bSJed Brown 
272c4762a1bSJed Brown /* ------------------------------------------------------------------- */
273c4762a1bSJed Brown /*
274c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
275c4762a1bSJed Brown 
276c4762a1bSJed Brown    Input Parameters:
277c4762a1bSJed Brown    user - user-defined application context
278c4762a1bSJed Brown    X - vector
279c4762a1bSJed Brown 
280c4762a1bSJed Brown    Output Parameter:
281c4762a1bSJed Brown    X - vector
282c4762a1bSJed Brown  */
2839371c9d4SSatish Balay static PetscErrorCode FormInitialGuess(AppCtx *user, DM da, Vec X) {
284c4762a1bSJed Brown   PetscInt      i, j, Mx, My, xs, ys, xm, ym;
285c4762a1bSJed Brown   PetscReal     temp1, temp, hx, hy;
286c4762a1bSJed Brown   PetscScalar **x;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   PetscFunctionBeginUser;
2899566063dSJacob 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));
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(Mx - 1);
292c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(My - 1);
293c4762a1bSJed Brown   temp1 = user->lambda / (user->lambda + 1.);
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /*
296c4762a1bSJed Brown      Get a pointer to vector data.
297c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
298c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
299c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
300c4762a1bSJed Brown          the array.
301c4762a1bSJed Brown   */
3029566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /*
305c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DA):
306c4762a1bSJed Brown        xs, ys   - starting grid indices (no ghost points)
307c4762a1bSJed Brown        xm, ym   - widths of local grid (no ghost points)
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   */
3109566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /*
313c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
314c4762a1bSJed Brown   */
315c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
316c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j, My - j - 1)) * hy;
317c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
318c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
319c4762a1bSJed Brown         /* boundary conditions are all zero Dirichlet */
320c4762a1bSJed Brown         x[j][i] = 0.0;
321c4762a1bSJed Brown       } else {
322c4762a1bSJed Brown         if (user->initial == -1) {
323c4762a1bSJed Brown           if (user->lambda != 0) {
324c4762a1bSJed Brown             x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp));
325c4762a1bSJed Brown           } else {
326c4762a1bSJed Brown             /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting
327c4762a1bSJed Brown              * with an exact solution. */
3289371c9d4SSatish Balay             const PetscReal xx = 2 * (PetscReal)i / (Mx - 1) - 1, yy = 2 * (PetscReal)j / (My - 1) - 1;
329c4762a1bSJed Brown             x[j][i] = (1 - xx * xx) * (1 - yy * yy) * xx * yy;
330c4762a1bSJed Brown           }
331c4762a1bSJed Brown         } else if (user->initial == 0) {
332c4762a1bSJed Brown           x[j][i] = 0.;
333c4762a1bSJed Brown         } else if (user->initial == 1) {
3349371c9d4SSatish Balay           const PetscReal xx = 2 * (PetscReal)i / (Mx - 1) - 1, yy = 2 * (PetscReal)j / (My - 1) - 1;
335c4762a1bSJed Brown           x[j][i] = (1 - xx * xx) * (1 - yy * yy) * xx * yy;
336c4762a1bSJed Brown         } else {
337c4762a1bSJed Brown           if (user->lambda != 0) {
338c4762a1bSJed Brown             x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp));
339c4762a1bSJed Brown           } else {
340c4762a1bSJed Brown             x[j][i] = 0.5 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp));
341c4762a1bSJed Brown           }
342c4762a1bSJed Brown         }
343c4762a1bSJed Brown       }
344c4762a1bSJed Brown     }
345c4762a1bSJed Brown   }
346c4762a1bSJed Brown   /*
347c4762a1bSJed Brown      Restore vector
348c4762a1bSJed Brown   */
3499566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
350c4762a1bSJed Brown   PetscFunctionReturn(0);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown 
353c4762a1bSJed Brown /*
354c4762a1bSJed Brown    FormRHS - Forms constant RHS for the problem.
355c4762a1bSJed Brown 
356c4762a1bSJed Brown    Input Parameters:
357c4762a1bSJed Brown    user - user-defined application context
358c4762a1bSJed Brown    B - RHS vector
359c4762a1bSJed Brown 
360c4762a1bSJed Brown    Output Parameter:
361c4762a1bSJed Brown    B - vector
362c4762a1bSJed Brown  */
3639371c9d4SSatish Balay static PetscErrorCode FormRHS(AppCtx *user, DM da, Vec B) {
364c4762a1bSJed Brown   PetscInt      i, j, Mx, My, xs, ys, xm, ym;
365c4762a1bSJed Brown   PetscReal     hx, hy;
366c4762a1bSJed Brown   PetscScalar **b;
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   PetscFunctionBeginUser;
3699566063dSJacob 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));
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(Mx - 1);
372c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(My - 1);
3739566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, B, &b));
3749566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
375c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
376c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
377c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
378c4762a1bSJed Brown         b[j][i] = 0.0;
379c4762a1bSJed Brown       } else {
380c4762a1bSJed Brown         b[j][i] = hx * hy * user->source;
381c4762a1bSJed Brown       }
382c4762a1bSJed Brown     }
383c4762a1bSJed Brown   }
3849566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, B, &b));
385c4762a1bSJed Brown   PetscFunctionReturn(0);
386c4762a1bSJed Brown }
387c4762a1bSJed Brown 
3889371c9d4SSatish Balay static inline PetscReal kappa(const AppCtx *ctx, PetscReal x, PetscReal y) {
389c4762a1bSJed Brown   return (((PetscInt)(x * ctx->blocks[0])) + ((PetscInt)(y * ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0;
390c4762a1bSJed Brown }
391c4762a1bSJed Brown /* p-Laplacian diffusivity */
3929371c9d4SSatish Balay static inline PetscScalar eta(const AppCtx *ctx, PetscReal x, PetscReal y, PetscScalar ux, PetscScalar uy) {
393c4762a1bSJed Brown   return kappa(ctx, x, y) * PetscPowScalar(PetscSqr(ctx->epsilon) + 0.5 * (ux * ux + uy * uy), 0.5 * (ctx->p - 2.));
394c4762a1bSJed Brown }
3959371c9d4SSatish Balay static inline PetscScalar deta(const AppCtx *ctx, PetscReal x, PetscReal y, PetscScalar ux, PetscScalar uy) {
3969371c9d4SSatish Balay   return (ctx->p == 2) ? 0 : kappa(ctx, x, y) * PetscPowScalar(PetscSqr(ctx->epsilon) + 0.5 * (ux * ux + uy * uy), 0.5 * (ctx->p - 4)) * 0.5 * (ctx->p - 2.);
397c4762a1bSJed Brown }
398c4762a1bSJed Brown 
399c4762a1bSJed Brown /* ------------------------------------------------------------------- */
400c4762a1bSJed Brown /*
401c4762a1bSJed Brown    FormFunctionLocal - Evaluates nonlinear function, F(x).
402c4762a1bSJed Brown  */
4039371c9d4SSatish Balay static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user) {
404c4762a1bSJed Brown   PetscReal   hx, hy, dhx, dhy, sc;
405c4762a1bSJed Brown   PetscInt    i, j;
406c4762a1bSJed Brown   PetscScalar eu;
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   PetscFunctionBeginUser;
409c4762a1bSJed Brown   hx  = 1.0 / (PetscReal)(info->mx - 1);
410c4762a1bSJed Brown   hy  = 1.0 / (PetscReal)(info->my - 1);
411c4762a1bSJed Brown   sc  = hx * hy * user->lambda;
412c4762a1bSJed Brown   dhx = 1 / hx;
413c4762a1bSJed Brown   dhy = 1 / hy;
414c4762a1bSJed Brown   /*
415c4762a1bSJed Brown      Compute function over the locally owned part of the grid
416c4762a1bSJed Brown   */
417c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
418c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
419c4762a1bSJed Brown       PetscReal xx = i * hx, yy = j * hy;
420c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
421c4762a1bSJed Brown         f[j][i] = x[j][i];
422c4762a1bSJed Brown       } else {
4239371c9d4SSatish Balay         const PetscScalar u = x[j][i], ux_E = dhx * (x[j][i + 1] - x[j][i]), uy_E = 0.25 * dhy * (x[j + 1][i] + x[j + 1][i + 1] - x[j - 1][i] - x[j - 1][i + 1]), ux_W = dhx * (x[j][i] - x[j][i - 1]), uy_W = 0.25 * dhy * (x[j + 1][i - 1] + x[j + 1][i] - x[j - 1][i - 1] - x[j - 1][i]), ux_N = 0.25 * dhx * (x[j][i + 1] + x[j + 1][i + 1] - x[j][i - 1] - x[j + 1][i - 1]), uy_N = dhy * (x[j + 1][i] - x[j][i]), ux_S = 0.25 * dhx * (x[j - 1][i + 1] + x[j][i + 1] - x[j - 1][i - 1] - x[j][i - 1]), uy_S = dhy * (x[j][i] - x[j - 1][i]), e_E = eta(user, xx, yy, ux_E, uy_E), e_W = eta(user, xx, yy, ux_W, uy_W), e_N = eta(user, xx, yy, ux_N, uy_N), e_S = eta(user, xx, yy, ux_S, uy_S), uxx = -hy * (e_E * ux_E - e_W * ux_W), uyy = -hx * (e_N * uy_N - e_S * uy_S);
424c4762a1bSJed Brown         if (sc) eu = PetscExpScalar(u);
425c4762a1bSJed Brown         else eu = 0.;
4260e3d61c9SBarry Smith         /* For p=2, these terms decay to:
4270e3d61c9SBarry Smith          uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx
4280e3d61c9SBarry Smith          uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy
4290e3d61c9SBarry Smith         */
430c4762a1bSJed Brown         f[j][i] = uxx + uyy - sc * eu;
431c4762a1bSJed Brown       }
432c4762a1bSJed Brown     }
433c4762a1bSJed Brown   }
4349566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm * info->ym * (72.0)));
435c4762a1bSJed Brown   PetscFunctionReturn(0);
436c4762a1bSJed Brown }
437c4762a1bSJed Brown 
438c4762a1bSJed Brown /*
439c4762a1bSJed Brown     This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation,
440c4762a1bSJed Brown     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)
441c4762a1bSJed Brown 
442c4762a1bSJed Brown */
4439371c9d4SSatish Balay static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user) {
444c4762a1bSJed Brown   PetscReal hx, hy, sc;
445c4762a1bSJed Brown   PetscInt  i, j;
446c4762a1bSJed Brown 
447c4762a1bSJed Brown   PetscFunctionBeginUser;
448c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info->mx - 1);
449c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(info->my - 1);
450c4762a1bSJed Brown   sc = hx * hy * user->lambda;
451c4762a1bSJed Brown   /*
452c4762a1bSJed Brown      Compute function over the locally owned part of the grid
453c4762a1bSJed Brown   */
454c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
455c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
456c4762a1bSJed Brown       if (!(i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1)) {
457c4762a1bSJed Brown         const PetscScalar u = x[j][i];
458c4762a1bSJed Brown         f[j][i]             = sc * PetscExpScalar(u);
4590df40c35SBarry Smith       } else {
4600df40c35SBarry Smith         f[j][i] = 0.0; /* this is zero because the A(x) x term forces the x to be zero on the boundary */
461c4762a1bSJed Brown       }
462c4762a1bSJed Brown     }
463c4762a1bSJed Brown   }
4649566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm * info->ym * 2.0));
465c4762a1bSJed Brown   PetscFunctionReturn(0);
466c4762a1bSJed Brown }
467c4762a1bSJed Brown 
468c4762a1bSJed Brown /*
469c4762a1bSJed Brown    FormJacobianLocal - Evaluates Jacobian matrix.
470c4762a1bSJed Brown */
4719371c9d4SSatish Balay static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat J, Mat B, AppCtx *user) {
472c4762a1bSJed Brown   PetscInt    i, j;
473c4762a1bSJed Brown   MatStencil  col[9], row;
474c4762a1bSJed Brown   PetscScalar v[9];
475c4762a1bSJed Brown   PetscReal   hx, hy, hxdhy, hydhx, dhx, dhy, sc;
476c4762a1bSJed Brown 
477c4762a1bSJed Brown   PetscFunctionBeginUser;
4789566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
479c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(info->mx - 1);
480c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(info->my - 1);
481c4762a1bSJed Brown   sc    = hx * hy * user->lambda;
482c4762a1bSJed Brown   dhx   = 1 / hx;
483c4762a1bSJed Brown   dhy   = 1 / hy;
484c4762a1bSJed Brown   hxdhy = hx / hy;
485c4762a1bSJed Brown   hydhx = hy / hx;
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   /*
488c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
489c4762a1bSJed Brown       - PETSc parallel matrix formats are partitioned by
490c4762a1bSJed Brown         contiguous chunks of rows across the processors.
491c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
492c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
493c4762a1bSJed Brown         appropriate processor during matrix assembly).
494c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
495c4762a1bSJed Brown   */
496c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
497c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
498c4762a1bSJed Brown       PetscReal xx = i * hx, yy = j * hy;
4999371c9d4SSatish Balay       row.j = j;
5009371c9d4SSatish Balay       row.i = i;
501c4762a1bSJed Brown       /* boundary points */
502c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
503c4762a1bSJed Brown         v[0] = 1.0;
5049566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B, 1, &row, 1, &row, v, INSERT_VALUES));
505c4762a1bSJed Brown       } else {
506c4762a1bSJed Brown         /* interior grid points */
5079371c9d4SSatish Balay         const PetscScalar ux_E = dhx * (x[j][i + 1] - x[j][i]), uy_E = 0.25 * dhy * (x[j + 1][i] + x[j + 1][i + 1] - x[j - 1][i] - x[j - 1][i + 1]), ux_W = dhx * (x[j][i] - x[j][i - 1]), uy_W = 0.25 * dhy * (x[j + 1][i - 1] + x[j + 1][i] - x[j - 1][i - 1] - x[j - 1][i]), ux_N = 0.25 * dhx * (x[j][i + 1] + x[j + 1][i + 1] - x[j][i - 1] - x[j + 1][i - 1]), uy_N = dhy * (x[j + 1][i] - x[j][i]), ux_S = 0.25 * dhx * (x[j - 1][i + 1] + x[j][i + 1] - x[j - 1][i - 1] - x[j][i - 1]), uy_S = dhy * (x[j][i] - x[j - 1][i]), u = x[j][i], e_E = eta(user, xx, yy, ux_E, uy_E), e_W = eta(user, xx, yy, ux_W, uy_W), e_N = eta(user, xx, yy, ux_N, uy_N), e_S = eta(user, xx, yy, ux_S, uy_S), de_E = deta(user, xx, yy, ux_E, uy_E), de_W = deta(user, xx, yy, ux_W, uy_W), de_N = deta(user, xx, yy, ux_N, uy_N), de_S = deta(user, xx, yy, ux_S, uy_S), skew_E = de_E * ux_E * uy_E, skew_W = de_W * ux_W * uy_W, skew_N = de_N * ux_N * uy_N, skew_S = de_S * ux_S * uy_S, cross_EW = 0.25 * (skew_E - skew_W), cross_NS = 0.25 * (skew_N - skew_S), newt_E = e_E + de_E * PetscSqr(ux_E), newt_W = e_W + de_W * PetscSqr(ux_W), newt_N = e_N + de_N * PetscSqr(uy_N), newt_S = e_S + de_S * PetscSqr(uy_S);
508c4762a1bSJed Brown         /* interior grid points */
509c4762a1bSJed Brown         switch (user->jtype) {
510c4762a1bSJed Brown         case JAC_BRATU:
511c4762a1bSJed Brown           /* Jacobian from p=2 */
5129371c9d4SSatish Balay           v[0]     = -hxdhy;
5139371c9d4SSatish Balay           col[0].j = j - 1;
5149371c9d4SSatish Balay           col[0].i = i;
5159371c9d4SSatish Balay           v[1]     = -hydhx;
5169371c9d4SSatish Balay           col[1].j = j;
5179371c9d4SSatish Balay           col[1].i = i - 1;
5189371c9d4SSatish Balay           v[2]     = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(u);
5199371c9d4SSatish Balay           col[2].j = row.j;
5209371c9d4SSatish Balay           col[2].i = row.i;
5219371c9d4SSatish Balay           v[3]     = -hydhx;
5229371c9d4SSatish Balay           col[3].j = j;
5239371c9d4SSatish Balay           col[3].i = i + 1;
5249371c9d4SSatish Balay           v[4]     = -hxdhy;
5259371c9d4SSatish Balay           col[4].j = j + 1;
5269371c9d4SSatish Balay           col[4].i = i;
5279566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES));
528c4762a1bSJed Brown           break;
529c4762a1bSJed Brown         case JAC_PICARD:
530c4762a1bSJed Brown           /* Jacobian arising from Picard linearization */
5319371c9d4SSatish Balay           v[0]     = -hxdhy * e_S;
5329371c9d4SSatish Balay           col[0].j = j - 1;
5339371c9d4SSatish Balay           col[0].i = i;
5349371c9d4SSatish Balay           v[1]     = -hydhx * e_W;
5359371c9d4SSatish Balay           col[1].j = j;
5369371c9d4SSatish Balay           col[1].i = i - 1;
5379371c9d4SSatish Balay           v[2]     = (e_W + e_E) * hydhx + (e_S + e_N) * hxdhy;
5389371c9d4SSatish Balay           col[2].j = row.j;
5399371c9d4SSatish Balay           col[2].i = row.i;
5409371c9d4SSatish Balay           v[3]     = -hydhx * e_E;
5419371c9d4SSatish Balay           col[3].j = j;
5429371c9d4SSatish Balay           col[3].i = i + 1;
5439371c9d4SSatish Balay           v[4]     = -hxdhy * e_N;
5449371c9d4SSatish Balay           col[4].j = j + 1;
5459371c9d4SSatish Balay           col[4].i = i;
5469566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES));
547c4762a1bSJed Brown           break;
548c4762a1bSJed Brown         case JAC_STAR:
549c4762a1bSJed Brown           /* Full Jacobian, but only a star stencil */
5509371c9d4SSatish Balay           col[0].j = j - 1;
5519371c9d4SSatish Balay           col[0].i = i;
5529371c9d4SSatish Balay           col[1].j = j;
5539371c9d4SSatish Balay           col[1].i = i - 1;
5549371c9d4SSatish Balay           col[2].j = j;
5559371c9d4SSatish Balay           col[2].i = i;
5569371c9d4SSatish Balay           col[3].j = j;
5579371c9d4SSatish Balay           col[3].i = i + 1;
5589371c9d4SSatish Balay           col[4].j = j + 1;
5599371c9d4SSatish Balay           col[4].i = i;
560c4762a1bSJed Brown           v[0]     = -hxdhy * newt_S + cross_EW;
561c4762a1bSJed Brown           v[1]     = -hydhx * newt_W + cross_NS;
562c4762a1bSJed Brown           v[2]     = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * PetscExpScalar(u);
563c4762a1bSJed Brown           v[3]     = -hydhx * newt_E - cross_NS;
564c4762a1bSJed Brown           v[4]     = -hxdhy * newt_N - cross_EW;
5659566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES));
566c4762a1bSJed Brown           break;
567c4762a1bSJed Brown         case JAC_NEWTON:
5680e3d61c9SBarry Smith           /* The Jacobian is
5690e3d61c9SBarry Smith 
5700e3d61c9SBarry Smith            -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u
5710e3d61c9SBarry Smith 
5720e3d61c9SBarry Smith           */
5739371c9d4SSatish Balay           col[0].j = j - 1;
5749371c9d4SSatish Balay           col[0].i = i - 1;
5759371c9d4SSatish Balay           col[1].j = j - 1;
5769371c9d4SSatish Balay           col[1].i = i;
5779371c9d4SSatish Balay           col[2].j = j - 1;
5789371c9d4SSatish Balay           col[2].i = i + 1;
5799371c9d4SSatish Balay           col[3].j = j;
5809371c9d4SSatish Balay           col[3].i = i - 1;
5819371c9d4SSatish Balay           col[4].j = j;
5829371c9d4SSatish Balay           col[4].i = i;
5839371c9d4SSatish Balay           col[5].j = j;
5849371c9d4SSatish Balay           col[5].i = i + 1;
5859371c9d4SSatish Balay           col[6].j = j + 1;
5869371c9d4SSatish Balay           col[6].i = i - 1;
5879371c9d4SSatish Balay           col[7].j = j + 1;
5889371c9d4SSatish Balay           col[7].i = i;
5899371c9d4SSatish Balay           col[8].j = j + 1;
5909371c9d4SSatish Balay           col[8].i = i + 1;
591c4762a1bSJed Brown           v[0]     = -0.25 * (skew_S + skew_W);
592c4762a1bSJed Brown           v[1]     = -hxdhy * newt_S + cross_EW;
593c4762a1bSJed Brown           v[2]     = 0.25 * (skew_S + skew_E);
594c4762a1bSJed Brown           v[3]     = -hydhx * newt_W + cross_NS;
595c4762a1bSJed Brown           v[4]     = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * PetscExpScalar(u);
596c4762a1bSJed Brown           v[5]     = -hydhx * newt_E - cross_NS;
597c4762a1bSJed Brown           v[6]     = 0.25 * (skew_N + skew_W);
598c4762a1bSJed Brown           v[7]     = -hxdhy * newt_N - cross_EW;
599c4762a1bSJed Brown           v[8]     = -0.25 * (skew_N + skew_E);
6009566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 9, col, v, INSERT_VALUES));
601c4762a1bSJed Brown           break;
6029371c9d4SSatish Balay         default: SETERRQ(PetscObjectComm((PetscObject)info->da), PETSC_ERR_SUP, "Jacobian type %d not implemented", user->jtype);
603c4762a1bSJed Brown         }
604c4762a1bSJed Brown       }
605c4762a1bSJed Brown     }
606c4762a1bSJed Brown   }
607c4762a1bSJed Brown 
608c4762a1bSJed Brown   /*
609c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
610c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
611c4762a1bSJed Brown   */
6129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
6139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   if (J != B) {
6169566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
6179566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
618c4762a1bSJed Brown   }
619c4762a1bSJed Brown 
620c4762a1bSJed Brown   /*
621c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
622c4762a1bSJed Brown      matrix. If we do, it will generate an error.
623c4762a1bSJed Brown   */
624*48a46eb9SPierre Jolivet   if (user->jtype == JAC_NEWTON) PetscCall(PetscLogFlops(info->xm * info->ym * (131.0)));
625c4762a1bSJed Brown   PetscFunctionReturn(0);
626c4762a1bSJed Brown }
627c4762a1bSJed Brown 
628c4762a1bSJed Brown /***********************************************************
629c4762a1bSJed Brown  * PreCheck implementation
630c4762a1bSJed Brown  ***********************************************************/
6319371c9d4SSatish Balay PetscErrorCode PreCheckSetFromOptions(PreCheck precheck) {
632c4762a1bSJed Brown   PetscBool flg;
633c4762a1bSJed Brown 
634c4762a1bSJed Brown   PetscFunctionBeginUser;
635d0609cedSBarry Smith   PetscOptionsBegin(precheck->comm, NULL, "PreCheck Options", "none");
6369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-precheck_angle", "Angle in degrees between successive search directions necessary to activate step correction", "", precheck->angle, &precheck->angle, NULL));
637c4762a1bSJed Brown   flg = PETSC_FALSE;
6389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-precheck_monitor", "Monitor choices made by precheck routine", "", flg, &flg, NULL));
639*48a46eb9SPierre Jolivet   if (flg) PetscCall(PetscViewerASCIIOpen(precheck->comm, "stdout", &precheck->monitor));
640d0609cedSBarry Smith   PetscOptionsEnd();
641c4762a1bSJed Brown   PetscFunctionReturn(0);
642c4762a1bSJed Brown }
643c4762a1bSJed Brown 
644c4762a1bSJed Brown /*
645c4762a1bSJed Brown   Compare the direction of the current and previous step, modify the current step accordingly
646c4762a1bSJed Brown */
6479371c9d4SSatish Balay PetscErrorCode PreCheckFunction(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) {
648c4762a1bSJed Brown   PreCheck    precheck;
649c4762a1bSJed Brown   Vec         Ylast;
650c4762a1bSJed Brown   PetscScalar dot;
651c4762a1bSJed Brown   PetscInt    iter;
652c4762a1bSJed Brown   PetscReal   ynorm, ylastnorm, theta, angle_radians;
653c4762a1bSJed Brown   SNES        snes;
654c4762a1bSJed Brown 
655c4762a1bSJed Brown   PetscFunctionBeginUser;
6569566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
657c4762a1bSJed Brown   precheck = (PreCheck)ctx;
6589566063dSJacob Faibussowitsch   if (!precheck->Ylast) PetscCall(VecDuplicate(Y, &precheck->Ylast));
659c4762a1bSJed Brown   Ylast = precheck->Ylast;
6609566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
661c4762a1bSJed Brown   if (iter < 1) {
6629566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
663c4762a1bSJed Brown     *changed = PETSC_FALSE;
664c4762a1bSJed Brown     PetscFunctionReturn(0);
665c4762a1bSJed Brown   }
666c4762a1bSJed Brown 
6679566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
6689566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
6699566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
670c4762a1bSJed Brown   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
671c4762a1bSJed Brown   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
672c4762a1bSJed Brown   angle_radians = precheck->angle * PETSC_PI / 180.;
673c4762a1bSJed Brown   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
674c4762a1bSJed Brown     /* Modify the step Y */
675c4762a1bSJed Brown     PetscReal alpha, ydiffnorm;
6769566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
6779566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
678c4762a1bSJed Brown     alpha = ylastnorm / ydiffnorm;
6799566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
6809566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
681*48a46eb9SPierre Jolivet     if (precheck->monitor) 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));
682c4762a1bSJed Brown   } else {
6839566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
684c4762a1bSJed Brown     *changed = PETSC_FALSE;
685*48a46eb9SPierre Jolivet     if (precheck->monitor) PetscCall(PetscViewerASCIIPrintf(precheck->monitor, "Angle %g degrees exceeds threshold %g, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)precheck->angle));
686c4762a1bSJed Brown   }
687c4762a1bSJed Brown   PetscFunctionReturn(0);
688c4762a1bSJed Brown }
689c4762a1bSJed Brown 
6909371c9d4SSatish Balay PetscErrorCode PreCheckDestroy(PreCheck *precheck) {
691c4762a1bSJed Brown   PetscFunctionBeginUser;
692c4762a1bSJed Brown   if (!*precheck) PetscFunctionReturn(0);
6939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*precheck)->Ylast));
6949566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*precheck)->monitor));
6959566063dSJacob Faibussowitsch   PetscCall(PetscFree(*precheck));
696c4762a1bSJed Brown   PetscFunctionReturn(0);
697c4762a1bSJed Brown }
698c4762a1bSJed Brown 
6999371c9d4SSatish Balay PetscErrorCode PreCheckCreate(MPI_Comm comm, PreCheck *precheck) {
700c4762a1bSJed Brown   PetscFunctionBeginUser;
7019566063dSJacob Faibussowitsch   PetscCall(PetscNew(precheck));
702c4762a1bSJed Brown 
703c4762a1bSJed Brown   (*precheck)->comm  = comm;
704c4762a1bSJed Brown   (*precheck)->angle = 10.; /* only active if angle is less than 10 degrees */
705c4762a1bSJed Brown   PetscFunctionReturn(0);
706c4762a1bSJed Brown }
707c4762a1bSJed Brown 
708c4762a1bSJed Brown /*
709c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
710c4762a1bSJed Brown 
711c4762a1bSJed Brown  */
7129371c9d4SSatish Balay PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx) {
713c4762a1bSJed Brown   PetscInt      i, j, k, xs, ys, xm, ym, its, tot_its, sweeps, l, m;
714c4762a1bSJed Brown   PetscReal     hx, hy, hxdhy, hydhx, dhx, dhy, sc;
715c4762a1bSJed Brown   PetscScalar **x, **b, bij, F, F0 = 0, J, y, u, eu;
716c4762a1bSJed Brown   PetscReal     atol, rtol, stol;
717c4762a1bSJed Brown   DM            da;
718c4762a1bSJed Brown   AppCtx       *user = (AppCtx *)ctx;
719c4762a1bSJed Brown   Vec           localX, localB;
720c4762a1bSJed Brown   DMDALocalInfo info;
721c4762a1bSJed Brown 
722c4762a1bSJed Brown   PetscFunctionBeginUser;
7239566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
7249566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
725c4762a1bSJed Brown 
726c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(info.mx - 1);
727c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(info.my - 1);
728c4762a1bSJed Brown   sc    = hx * hy * user->lambda;
729c4762a1bSJed Brown   dhx   = 1 / hx;
730c4762a1bSJed Brown   dhy   = 1 / hy;
731c4762a1bSJed Brown   hxdhy = hx / hy;
732c4762a1bSJed Brown   hydhx = hy / hx;
733c4762a1bSJed Brown 
734c4762a1bSJed Brown   tot_its = 0;
7359566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetSweeps(snes, &sweeps));
7369566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its));
7379566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
738*48a46eb9SPierre Jolivet   if (B) PetscCall(DMGetLocalVector(da, &localB));
739c4762a1bSJed Brown   if (B) {
7409566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, localB));
7419566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, localB));
742c4762a1bSJed Brown   }
7439566063dSJacob Faibussowitsch   if (B) PetscCall(DMDAVecGetArrayRead(da, localB, &b));
7449566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
7459566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
7469566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
747c4762a1bSJed Brown   for (l = 0; l < sweeps; l++) {
748c4762a1bSJed Brown     /*
749c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
750c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
751c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
752c4762a1bSJed Brown      */
7539566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
754c4762a1bSJed Brown     for (m = 0; m < 2; m++) {
755c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
756c4762a1bSJed Brown         for (i = xs + (m + j) % 2; i < xs + xm; i += 2) {
757c4762a1bSJed Brown           PetscReal xx = i * hx, yy = j * hy;
758c4762a1bSJed Brown           if (B) bij = b[j][i];
759c4762a1bSJed Brown           else bij = 0.;
760c4762a1bSJed Brown 
761c4762a1bSJed Brown           if (i == 0 || j == 0 || i == info.mx - 1 || j == info.my - 1) {
762c4762a1bSJed Brown             /* boundary conditions are all zero Dirichlet */
763c4762a1bSJed Brown             x[j][i] = 0.0 + bij;
764c4762a1bSJed Brown           } else {
7659371c9d4SSatish Balay             const PetscScalar u_E = x[j][i + 1], u_W = x[j][i - 1], u_N = x[j + 1][i], u_S = x[j - 1][i];
7669371c9d4SSatish Balay             const PetscScalar uy_E = 0.25 * dhy * (x[j + 1][i] + x[j + 1][i + 1] - x[j - 1][i] - x[j - 1][i + 1]), uy_W = 0.25 * dhy * (x[j + 1][i - 1] + x[j + 1][i] - x[j - 1][i - 1] - x[j - 1][i]), ux_N = 0.25 * dhx * (x[j][i + 1] + x[j + 1][i + 1] - x[j][i - 1] - x[j + 1][i - 1]), ux_S = 0.25 * dhx * (x[j - 1][i + 1] + x[j][i + 1] - x[j - 1][i - 1] - x[j][i - 1]);
767c4762a1bSJed Brown             u = x[j][i];
768c4762a1bSJed Brown             for (k = 0; k < its; k++) {
7699371c9d4SSatish Balay               const PetscScalar ux_E = dhx * (u_E - u), ux_W = dhx * (u - u_W), uy_N = dhy * (u_N - u), uy_S = dhy * (u - u_S), e_E = eta(user, xx, yy, ux_E, uy_E), e_W = eta(user, xx, yy, ux_W, uy_W), e_N = eta(user, xx, yy, ux_N, uy_N), e_S = eta(user, xx, yy, ux_S, uy_S), de_E = deta(user, xx, yy, ux_E, uy_E), de_W = deta(user, xx, yy, ux_W, uy_W), de_N = deta(user, xx, yy, ux_N, uy_N), de_S = deta(user, xx, yy, ux_S, uy_S), newt_E = e_E + de_E * PetscSqr(ux_E), newt_W = e_W + de_W * PetscSqr(ux_W), newt_N = e_N + de_N * PetscSqr(uy_N), newt_S = e_S + de_S * PetscSqr(uy_S), uxx = -hy * (e_E * ux_E - e_W * ux_W), uyy = -hx * (e_N * uy_N - e_S * uy_S);
770c4762a1bSJed Brown 
771c4762a1bSJed Brown               if (sc) eu = PetscExpScalar(u);
772c4762a1bSJed Brown               else eu = 0;
773c4762a1bSJed Brown 
774c4762a1bSJed Brown               F = uxx + uyy - sc * eu - bij;
775c4762a1bSJed Brown               if (k == 0) F0 = F;
776c4762a1bSJed Brown               J = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * eu;
777c4762a1bSJed Brown               y = F / J;
778c4762a1bSJed Brown               u -= y;
779c4762a1bSJed Brown               tot_its++;
7809371c9d4SSatish Balay               if (atol > PetscAbsReal(PetscRealPart(F)) || rtol * PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || stol * PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) { break; }
781c4762a1bSJed Brown             }
782c4762a1bSJed Brown             x[j][i] = u;
783c4762a1bSJed Brown           }
784c4762a1bSJed Brown         }
785c4762a1bSJed Brown       }
786c4762a1bSJed Brown     }
787c4762a1bSJed Brown     /*
788c4762a1bSJed Brown x     Restore vector
789c4762a1bSJed Brown      */
790c4762a1bSJed Brown   }
7919566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
7929566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X));
7939566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X));
7949566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(tot_its * (118.0)));
7959566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
796c4762a1bSJed Brown   if (B) {
7979566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, localB, &b));
7989566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da, &localB));
799c4762a1bSJed Brown   }
800c4762a1bSJed Brown   PetscFunctionReturn(0);
801c4762a1bSJed Brown }
802c4762a1bSJed Brown 
803c4762a1bSJed Brown /*TEST
804c4762a1bSJed Brown 
805c4762a1bSJed Brown    test:
806c4762a1bSJed Brown       nsize: 2
807c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype NEWTON
808c4762a1bSJed Brown       requires: !single
809c4762a1bSJed Brown 
810c4762a1bSJed Brown    test:
811c4762a1bSJed Brown       suffix: 2
812c4762a1bSJed Brown       nsize: 2
813c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -precheck 1
814c4762a1bSJed Brown       requires: !single
815c4762a1bSJed Brown 
816c4762a1bSJed Brown    test:
817c4762a1bSJed Brown       suffix: 3
818c4762a1bSJed Brown       nsize: 2
81985b95db3SBarry Smith       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -picard -precheck 1 -p 3
82085b95db3SBarry Smith       requires: !single
82185b95db3SBarry Smith 
82285b95db3SBarry Smith    test:
82385b95db3SBarry Smith       suffix: 3_star
82485b95db3SBarry Smith       nsize: 2
82585b95db3SBarry Smith       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
82685b95db3SBarry Smith       output_file: output/ex15_3.out
827c4762a1bSJed Brown       requires: !single
828c4762a1bSJed Brown 
829c4762a1bSJed Brown    test:
830c4762a1bSJed Brown       suffix: 4
831c4762a1bSJed Brown       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
832c4762a1bSJed Brown       requires: !single
833c4762a1bSJed Brown 
834c4762a1bSJed Brown    test:
835c4762a1bSJed Brown       suffix: lag_jac
836c4762a1bSJed Brown       nsize: 4
837c4762a1bSJed Brown       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
838c4762a1bSJed Brown       requires: !single
839c4762a1bSJed Brown 
840c4762a1bSJed Brown    test:
841c4762a1bSJed Brown       suffix: lag_pc
842c4762a1bSJed Brown       nsize: 4
843c4762a1bSJed Brown       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
844c4762a1bSJed Brown       requires: !single
845c4762a1bSJed Brown 
846c4762a1bSJed Brown    test:
847c4762a1bSJed Brown       suffix: nleqerr
848c4762a1bSJed Brown       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
849c4762a1bSJed Brown       requires: !single
850c4762a1bSJed Brown 
851bbc1464cSBarry Smith    test:
852bbc1464cSBarry Smith       suffix: mf
853bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator
854bbc1464cSBarry Smith       requires: !single
855bbc1464cSBarry Smith 
856bbc1464cSBarry Smith    test:
857bbc1464cSBarry Smith       suffix: mf_picard
858bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator -picard
859bbc1464cSBarry Smith       requires: !single
860bbc1464cSBarry Smith       output_file: output/ex15_mf.out
861bbc1464cSBarry Smith 
86285b95db3SBarry Smith    test:
86385b95db3SBarry Smith       suffix: fd_picard
86485b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 2  -p 3 -ksp_rtol 1.e-12  -snes_fd -picard
86585b95db3SBarry Smith       requires: !single
86685b95db3SBarry Smith 
86785b95db3SBarry Smith    test:
86885b95db3SBarry Smith       suffix: fd_color_picard
86985b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_fd_color -picard
87085b95db3SBarry Smith       requires: !single
87185b95db3SBarry Smith       output_file: output/ex15_mf.out
87285b95db3SBarry Smith 
873c4762a1bSJed Brown TEST*/
874