xref: /petsc/src/snes/tutorials/ex15.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 /*
35d75802c7SJacob Faibussowitsch       mpiexec -n 2 ./ex15 -snes_monitor -ksp_monitor log_view
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 
main(int argc,char ** argv)100d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
101d71ae5a4SJacob Faibussowitsch {
102c4762a1bSJed Brown   SNES                snes;       /* nonlinear solver */
103c4762a1bSJed Brown   Vec                 x, r, b;    /* solution, residual, rhs vectors */
104c4762a1bSJed Brown   AppCtx              user;       /* user-defined work context */
105c4762a1bSJed Brown   PetscInt            its;        /* iterations for convergence */
106c4762a1bSJed Brown   SNESConvergedReason reason;     /* Check convergence */
107c4762a1bSJed Brown   PetscBool           alloc_star; /* Only allocate for the STAR stencil  */
108c4762a1bSJed Brown   PetscBool           write_output;
109c4762a1bSJed Brown   char                filename[PETSC_MAX_PATH_LEN] = "ex15.vts";
110c4762a1bSJed Brown   PetscReal           bratu_lambda_max = 6.81, bratu_lambda_min = 0.;
11185b95db3SBarry Smith   DM                  da;              /* distributed array data structure */
112c4762a1bSJed Brown   PreCheck            precheck = NULL; /* precheck context for version in this file */
113c4762a1bSJed Brown   PetscInt            use_precheck;    /* 0=none, 1=version in this file, 2=SNES-provided version */
114c4762a1bSJed Brown   PetscReal           precheck_angle;  /* When manually setting the SNES-provided precheck function */
115c4762a1bSJed Brown   SNESLineSearch      linesearch;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118c4762a1bSJed Brown      Initialize program
119c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120c4762a1bSJed Brown 
121327415f7SBarry Smith   PetscFunctionBeginUser;
122c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125c4762a1bSJed Brown      Initialize problem parameters
126c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1279371c9d4SSatish Balay   user.lambda    = 0.0;
1289371c9d4SSatish Balay   user.p         = 2.0;
1299371c9d4SSatish Balay   user.epsilon   = 1e-5;
1309371c9d4SSatish Balay   user.source    = 0.1;
1319371c9d4SSatish Balay   user.jtype     = JAC_NEWTON;
1329371c9d4SSatish Balay   user.initial   = -1;
1339371c9d4SSatish Balay   user.blocks[0] = 1;
1349371c9d4SSatish Balay   user.blocks[1] = 1;
1359371c9d4SSatish Balay   user.kappa     = 1e-3;
136c4762a1bSJed Brown   alloc_star     = PETSC_FALSE;
1379371c9d4SSatish Balay   use_precheck   = 0;
1389371c9d4SSatish Balay   precheck_angle = 10.;
139c4762a1bSJed Brown   user.picard    = PETSC_FALSE;
140d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "p-Bratu options", __FILE__);
141c4762a1bSJed Brown   {
142c4762a1bSJed Brown     PetscInt two = 2;
1439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda", "Bratu parameter", "", user.lambda, &user.lambda, NULL));
1449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-p", "Exponent `p' in p-Laplacian", "", user.p, &user.p, NULL));
1459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-epsilon", "Strain-regularization in p-Laplacian", "", user.epsilon, &user.epsilon, NULL));
1469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-source", "Constant source term", "", user.source, &user.source, NULL));
1479566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-jtype", "Jacobian approximation to assemble", "", JacTypes, (PetscEnum)user.jtype, (PetscEnum *)&user.jtype, NULL));
1489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-picard", "Solve with defect-correction Picard iteration", "", &user.picard));
14985b95db3SBarry Smith     if (user.picard) {
15085b95db3SBarry Smith       user.jtype = JAC_PICARD;
151e00437b9SBarry Smith       PetscCheck(user.p == 3, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Picard iteration is only supported for p == 3");
15285b95db3SBarry Smith       /* the Picard linearization only requires a 5 point stencil, while the Newton linearization requires a 9 point stencil */
15385b95db3SBarry Smith       /* hence allocating the 5 point stencil gives the same convergence as the 9 point stencil since the extra stencil points are not used */
1549566063dSJacob Faibussowitsch       PetscCall(PetscOptionsBool("-alloc_star", "Allocate for STAR stencil (5-point)", "", alloc_star, &alloc_star, NULL));
15585b95db3SBarry Smith     }
1569566063dSJacob 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));
1579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-initial", "Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)", "", user.initial, &user.initial, NULL));
158c4762a1bSJed Brown     if (use_precheck == 2) { /* Using library version, get the angle */
1599566063dSJacob Faibussowitsch       PetscCall(PetscOptionsReal("-precheck_angle", "Angle in degrees between successive search directions necessary to activate step correction", "", precheck_angle, &precheck_angle, NULL));
160c4762a1bSJed Brown     }
1619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-blocks", "number of coefficient interfaces in x and y direction", "", user.blocks, &two, NULL));
1629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-kappa", "diffusivity in odd regions", "", user.kappa, &user.kappa, NULL));
1639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-o", "Output solution in vts format", "", filename, filename, sizeof(filename), &write_output));
164c4762a1bSJed Brown   }
165d0609cedSBarry Smith   PetscOptionsEnd();
16648a46eb9SPierre 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));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169c4762a1bSJed Brown      Create nonlinear solver context
170c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1719566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
175c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1769566063dSJacob 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));
1779566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1789566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown      Extract global vectors from DM; then duplicate for remaining
182c4762a1bSJed Brown      vectors that are the same types
183c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1849566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
1859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
1869566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &b));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18985b95db3SBarry Smith      User can override with:
190c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
191c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
1927addb90fSBarry Smith      -snes_mf_operator : form matrix used to construct the preconditioner as set by the user,
193c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
194c4762a1bSJed Brown                          products within Newton-Krylov method
195c4762a1bSJed Brown 
196c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199c4762a1bSJed Brown      Set local function evaluation routine
200c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2019566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
2029566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
203c4762a1bSJed Brown   if (user.picard) {
204c4762a1bSJed Brown     /*
205c4762a1bSJed Brown         This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access
206c4762a1bSJed Brown         the SNES to set it
207c4762a1bSJed Brown     */
2089371c9d4SSatish Balay     PetscCall(DMDASNESSetPicardLocal(da, INSERT_VALUES, (PetscErrorCode (*)(DMDALocalInfo *, void *, void *, void *))FormFunctionPicardLocal, (PetscErrorCode (*)(DMDALocalInfo *, void *, Mat, Mat, void *))FormJacobianLocal, &user));
2099566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, SNESPicardComputeFunction, &user));
2109566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESPicardComputeJacobian, &user));
211c4762a1bSJed Brown   } else {
2129566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode (*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user));
2139566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetJacobianLocal(da, (PetscErrorCode (*)(DMDALocalInfo *, void *, Mat, Mat, void *))FormJacobianLocal, &user));
214c4762a1bSJed Brown   }
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
218c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2199566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
2209566063dSJacob Faibussowitsch   PetscCall(SNESSetNGS(snes, NonlinearGS, &user));
2219566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
222c4762a1bSJed Brown   /* Set up the precheck context if requested */
223c4762a1bSJed Brown   if (use_precheck == 1) { /* Use the precheck routines in this file */
2249566063dSJacob Faibussowitsch     PetscCall(PreCheckCreate(PETSC_COMM_WORLD, &precheck));
2259566063dSJacob Faibussowitsch     PetscCall(PreCheckSetFromOptions(precheck));
2269566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch, PreCheckFunction, precheck));
227c4762a1bSJed Brown   } else if (use_precheck == 2) { /* Use the version provided by the library */
2289566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &precheck_angle));
229c4762a1bSJed Brown   }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232c4762a1bSJed Brown      Evaluate initial guess
233c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
234c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
235c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
236c4762a1bSJed Brown      this vector to zero by calling VecSet().
237c4762a1bSJed Brown   */
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user, da, x));
2409566063dSJacob Faibussowitsch   PetscCall(FormRHS(&user, da, b));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
243c4762a1bSJed Brown      Solve nonlinear system
244c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2459566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, b, x));
2469566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
2479566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes, &reason));
248c4762a1bSJed Brown 
24963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s Number of nonlinear iterations = %" PetscInt_FMT "\n", SNESConvergedReasons[reason], its));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   if (write_output) {
252c4762a1bSJed Brown     PetscViewer viewer;
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD, filename, FILE_MODE_WRITE, &viewer));
2549566063dSJacob Faibussowitsch     PetscCall(VecView(x, viewer));
2559566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
256c4762a1bSJed Brown   }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
259c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
260c4762a1bSJed Brown      are no longer needed.
261c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262c4762a1bSJed Brown 
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
2669566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2679566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2689566063dSJacob Faibussowitsch   PetscCall(PreCheckDestroy(&precheck));
2699566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
270b122ec5aSJacob Faibussowitsch   return 0;
271c4762a1bSJed Brown }
272c4762a1bSJed Brown 
273c4762a1bSJed Brown /* ------------------------------------------------------------------- */
274c4762a1bSJed Brown /*
275c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    Input Parameters:
278c4762a1bSJed Brown    user - user-defined application context
279c4762a1bSJed Brown    X - vector
280c4762a1bSJed Brown 
281c4762a1bSJed Brown    Output Parameter:
282c4762a1bSJed Brown    X - vector
283c4762a1bSJed Brown  */
FormInitialGuess(AppCtx * user,DM da,Vec X)284d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormInitialGuess(AppCtx *user, DM da, Vec X)
285d71ae5a4SJacob Faibussowitsch {
286c4762a1bSJed Brown   PetscInt      i, j, Mx, My, xs, ys, xm, ym;
287c4762a1bSJed Brown   PetscReal     temp1, temp, hx, hy;
288c4762a1bSJed Brown   PetscScalar **x;
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   PetscFunctionBeginUser;
2919566063dSJacob 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));
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(Mx - 1);
294c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(My - 1);
295c4762a1bSJed Brown   temp1 = user->lambda / (user->lambda + 1.);
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /*
298c4762a1bSJed Brown      Get a pointer to vector data.
299c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
300c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
301c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
302c4762a1bSJed Brown          the array.
303c4762a1bSJed Brown   */
3049566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   /*
307c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DA):
308c4762a1bSJed Brown        xs, ys   - starting grid indices (no ghost points)
309c4762a1bSJed Brown        xm, ym   - widths of local grid (no ghost points)
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   */
3129566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   /*
315c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
316c4762a1bSJed Brown   */
317c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
318c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j, My - j - 1)) * hy;
319c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
320c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
321c4762a1bSJed Brown         /* boundary conditions are all zero Dirichlet */
322c4762a1bSJed Brown         x[j][i] = 0.0;
323c4762a1bSJed Brown       } else {
324c4762a1bSJed Brown         if (user->initial == -1) {
325c4762a1bSJed Brown           if (user->lambda != 0) {
326c4762a1bSJed Brown             x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp));
327c4762a1bSJed Brown           } else {
328c4762a1bSJed Brown             /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting
329c4762a1bSJed Brown              * with an exact solution. */
3309371c9d4SSatish Balay             const PetscReal xx = 2 * (PetscReal)i / (Mx - 1) - 1, yy = 2 * (PetscReal)j / (My - 1) - 1;
331c4762a1bSJed Brown             x[j][i] = (1 - xx * xx) * (1 - yy * yy) * xx * yy;
332c4762a1bSJed Brown           }
333c4762a1bSJed Brown         } else if (user->initial == 0) {
334c4762a1bSJed Brown           x[j][i] = 0.;
335c4762a1bSJed Brown         } else if (user->initial == 1) {
3369371c9d4SSatish Balay           const PetscReal xx = 2 * (PetscReal)i / (Mx - 1) - 1, yy = 2 * (PetscReal)j / (My - 1) - 1;
337c4762a1bSJed Brown           x[j][i] = (1 - xx * xx) * (1 - yy * yy) * xx * yy;
338c4762a1bSJed Brown         } else {
339c4762a1bSJed Brown           if (user->lambda != 0) {
340c4762a1bSJed Brown             x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp));
341c4762a1bSJed Brown           } else {
342c4762a1bSJed Brown             x[j][i] = 0.5 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp));
343c4762a1bSJed Brown           }
344c4762a1bSJed Brown         }
345c4762a1bSJed Brown       }
346c4762a1bSJed Brown     }
347c4762a1bSJed Brown   }
348c4762a1bSJed Brown   /*
349c4762a1bSJed Brown      Restore vector
350c4762a1bSJed Brown   */
3519566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
353c4762a1bSJed Brown }
354c4762a1bSJed Brown 
355c4762a1bSJed Brown /*
356c4762a1bSJed Brown    FormRHS - Forms constant RHS for the problem.
357c4762a1bSJed Brown 
358c4762a1bSJed Brown    Input Parameters:
359c4762a1bSJed Brown    user - user-defined application context
360c4762a1bSJed Brown    B - RHS vector
361c4762a1bSJed Brown 
362c4762a1bSJed Brown    Output Parameter:
363c4762a1bSJed Brown    B - vector
364c4762a1bSJed Brown  */
FormRHS(AppCtx * user,DM da,Vec B)365d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormRHS(AppCtx *user, DM da, Vec B)
366d71ae5a4SJacob Faibussowitsch {
367c4762a1bSJed Brown   PetscInt      i, j, Mx, My, xs, ys, xm, ym;
368c4762a1bSJed Brown   PetscReal     hx, hy;
369c4762a1bSJed Brown   PetscScalar **b;
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   PetscFunctionBeginUser;
3729566063dSJacob 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));
373c4762a1bSJed Brown 
374c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(Mx - 1);
375c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(My - 1);
3769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, B, &b));
3779566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
378c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
379c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
380c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
381c4762a1bSJed Brown         b[j][i] = 0.0;
382c4762a1bSJed Brown       } else {
383c4762a1bSJed Brown         b[j][i] = hx * hy * user->source;
384c4762a1bSJed Brown       }
385c4762a1bSJed Brown     }
386c4762a1bSJed Brown   }
3879566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, B, &b));
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
389c4762a1bSJed Brown }
390c4762a1bSJed Brown 
kappa(const AppCtx * ctx,PetscReal x,PetscReal y)391d71ae5a4SJacob Faibussowitsch static inline PetscReal kappa(const AppCtx *ctx, PetscReal x, PetscReal y)
392d71ae5a4SJacob Faibussowitsch {
393c4762a1bSJed Brown   return (((PetscInt)(x * ctx->blocks[0])) + ((PetscInt)(y * ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0;
394c4762a1bSJed Brown }
395c4762a1bSJed Brown /* p-Laplacian diffusivity */
eta(const AppCtx * ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)396d71ae5a4SJacob Faibussowitsch static inline PetscScalar eta(const AppCtx *ctx, PetscReal x, PetscReal y, PetscScalar ux, PetscScalar uy)
397d71ae5a4SJacob Faibussowitsch {
398c4762a1bSJed Brown   return kappa(ctx, x, y) * PetscPowScalar(PetscSqr(ctx->epsilon) + 0.5 * (ux * ux + uy * uy), 0.5 * (ctx->p - 2.));
399c4762a1bSJed Brown }
deta(const AppCtx * ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)400d71ae5a4SJacob Faibussowitsch static inline PetscScalar deta(const AppCtx *ctx, PetscReal x, PetscReal y, PetscScalar ux, PetscScalar uy)
401d71ae5a4SJacob Faibussowitsch {
4029371c9d4SSatish 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.);
403c4762a1bSJed Brown }
404c4762a1bSJed Brown 
405c4762a1bSJed Brown /* ------------------------------------------------------------------- */
406c4762a1bSJed Brown /*
407c4762a1bSJed Brown    FormFunctionLocal - Evaluates nonlinear function, F(x).
408c4762a1bSJed Brown  */
FormFunctionLocal(DMDALocalInfo * info,PetscScalar ** x,PetscScalar ** f,AppCtx * user)409d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user)
410d71ae5a4SJacob Faibussowitsch {
411c4762a1bSJed Brown   PetscReal   hx, hy, dhx, dhy, sc;
412c4762a1bSJed Brown   PetscInt    i, j;
413c4762a1bSJed Brown   PetscScalar eu;
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   PetscFunctionBeginUser;
416c4762a1bSJed Brown   hx  = 1.0 / (PetscReal)(info->mx - 1);
417c4762a1bSJed Brown   hy  = 1.0 / (PetscReal)(info->my - 1);
418c4762a1bSJed Brown   sc  = hx * hy * user->lambda;
419c4762a1bSJed Brown   dhx = 1 / hx;
420c4762a1bSJed Brown   dhy = 1 / hy;
421c4762a1bSJed Brown   /*
422c4762a1bSJed Brown      Compute function over the locally owned part of the grid
423c4762a1bSJed Brown   */
424c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
425c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
426c4762a1bSJed Brown       PetscReal xx = i * hx, yy = j * hy;
427c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
428c4762a1bSJed Brown         f[j][i] = x[j][i];
429c4762a1bSJed Brown       } else {
4309371c9d4SSatish 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);
431c4762a1bSJed Brown         if (sc) eu = PetscExpScalar(u);
432c4762a1bSJed Brown         else eu = 0.;
4330e3d61c9SBarry Smith         /* For p=2, these terms decay to:
4340e3d61c9SBarry Smith          uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx
4350e3d61c9SBarry Smith          uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy
4360e3d61c9SBarry Smith         */
437c4762a1bSJed Brown         f[j][i] = uxx + uyy - sc * eu;
438c4762a1bSJed Brown       }
439c4762a1bSJed Brown     }
440c4762a1bSJed Brown   }
4419566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm * info->ym * (72.0)));
4423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
443c4762a1bSJed Brown }
444c4762a1bSJed Brown 
445c4762a1bSJed Brown /*
446c4762a1bSJed Brown     This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation,
447c4762a1bSJed 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)
448c4762a1bSJed Brown 
449c4762a1bSJed Brown */
FormFunctionPicardLocal(DMDALocalInfo * info,PetscScalar ** x,PetscScalar ** f,AppCtx * user)450d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user)
451d71ae5a4SJacob Faibussowitsch {
452c4762a1bSJed Brown   PetscReal hx, hy, sc;
453c4762a1bSJed Brown   PetscInt  i, j;
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   PetscFunctionBeginUser;
456c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info->mx - 1);
457c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(info->my - 1);
458c4762a1bSJed Brown   sc = hx * hy * user->lambda;
459c4762a1bSJed Brown   /*
460c4762a1bSJed Brown      Compute function over the locally owned part of the grid
461c4762a1bSJed Brown   */
462c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
463c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
464c4762a1bSJed Brown       if (!(i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1)) {
465c4762a1bSJed Brown         const PetscScalar u = x[j][i];
466c4762a1bSJed Brown         f[j][i]             = sc * PetscExpScalar(u);
4670df40c35SBarry Smith       } else {
4680df40c35SBarry Smith         f[j][i] = 0.0; /* this is zero because the A(x) x term forces the x to be zero on the boundary */
469c4762a1bSJed Brown       }
470c4762a1bSJed Brown     }
471c4762a1bSJed Brown   }
4729566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm * info->ym * 2.0));
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
474c4762a1bSJed Brown }
475c4762a1bSJed Brown 
476c4762a1bSJed Brown /*
477c4762a1bSJed Brown    FormJacobianLocal - Evaluates Jacobian matrix.
478c4762a1bSJed Brown */
FormJacobianLocal(DMDALocalInfo * info,PetscScalar ** x,Mat J,Mat B,AppCtx * user)479d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat J, Mat B, AppCtx *user)
480d71ae5a4SJacob Faibussowitsch {
481c4762a1bSJed Brown   PetscInt    i, j;
482c4762a1bSJed Brown   MatStencil  col[9], row;
483c4762a1bSJed Brown   PetscScalar v[9];
484c4762a1bSJed Brown   PetscReal   hx, hy, hxdhy, hydhx, dhx, dhy, sc;
485c4762a1bSJed Brown 
486c4762a1bSJed Brown   PetscFunctionBeginUser;
4879566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
488c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(info->mx - 1);
489c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(info->my - 1);
490c4762a1bSJed Brown   sc    = hx * hy * user->lambda;
491c4762a1bSJed Brown   dhx   = 1 / hx;
492c4762a1bSJed Brown   dhy   = 1 / hy;
493c4762a1bSJed Brown   hxdhy = hx / hy;
494c4762a1bSJed Brown   hydhx = hy / hx;
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   /*
497c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
498c4762a1bSJed Brown       - PETSc parallel matrix formats are partitioned by
499c4762a1bSJed Brown         contiguous chunks of rows across the processors.
500c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
501c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
502c4762a1bSJed Brown         appropriate processor during matrix assembly).
503c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
504c4762a1bSJed Brown   */
505c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
506c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
507c4762a1bSJed Brown       PetscReal xx = i * hx, yy = j * hy;
5089371c9d4SSatish Balay       row.j = j;
5099371c9d4SSatish Balay       row.i = i;
510c4762a1bSJed Brown       /* boundary points */
511c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
512c4762a1bSJed Brown         v[0] = 1.0;
5139566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B, 1, &row, 1, &row, v, INSERT_VALUES));
514c4762a1bSJed Brown       } else {
515c4762a1bSJed Brown         /* interior grid points */
5169371c9d4SSatish 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);
517c4762a1bSJed Brown         /* interior grid points */
518c4762a1bSJed Brown         switch (user->jtype) {
519c4762a1bSJed Brown         case JAC_BRATU:
520c4762a1bSJed Brown           /* Jacobian from p=2 */
5219371c9d4SSatish Balay           v[0]     = -hxdhy;
5229371c9d4SSatish Balay           col[0].j = j - 1;
5239371c9d4SSatish Balay           col[0].i = i;
5249371c9d4SSatish Balay           v[1]     = -hydhx;
5259371c9d4SSatish Balay           col[1].j = j;
5269371c9d4SSatish Balay           col[1].i = i - 1;
5279371c9d4SSatish Balay           v[2]     = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(u);
5289371c9d4SSatish Balay           col[2].j = row.j;
5299371c9d4SSatish Balay           col[2].i = row.i;
5309371c9d4SSatish Balay           v[3]     = -hydhx;
5319371c9d4SSatish Balay           col[3].j = j;
5329371c9d4SSatish Balay           col[3].i = i + 1;
5339371c9d4SSatish Balay           v[4]     = -hxdhy;
5349371c9d4SSatish Balay           col[4].j = j + 1;
5359371c9d4SSatish Balay           col[4].i = i;
5369566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES));
537c4762a1bSJed Brown           break;
538c4762a1bSJed Brown         case JAC_PICARD:
539c4762a1bSJed Brown           /* Jacobian arising from Picard linearization */
5409371c9d4SSatish Balay           v[0]     = -hxdhy * e_S;
5419371c9d4SSatish Balay           col[0].j = j - 1;
5429371c9d4SSatish Balay           col[0].i = i;
5439371c9d4SSatish Balay           v[1]     = -hydhx * e_W;
5449371c9d4SSatish Balay           col[1].j = j;
5459371c9d4SSatish Balay           col[1].i = i - 1;
5469371c9d4SSatish Balay           v[2]     = (e_W + e_E) * hydhx + (e_S + e_N) * hxdhy;
5479371c9d4SSatish Balay           col[2].j = row.j;
5489371c9d4SSatish Balay           col[2].i = row.i;
5499371c9d4SSatish Balay           v[3]     = -hydhx * e_E;
5509371c9d4SSatish Balay           col[3].j = j;
5519371c9d4SSatish Balay           col[3].i = i + 1;
5529371c9d4SSatish Balay           v[4]     = -hxdhy * e_N;
5539371c9d4SSatish Balay           col[4].j = j + 1;
5549371c9d4SSatish Balay           col[4].i = i;
5559566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES));
556c4762a1bSJed Brown           break;
557c4762a1bSJed Brown         case JAC_STAR:
558c4762a1bSJed Brown           /* Full Jacobian, but only a star stencil */
5599371c9d4SSatish Balay           col[0].j = j - 1;
5609371c9d4SSatish Balay           col[0].i = i;
5619371c9d4SSatish Balay           col[1].j = j;
5629371c9d4SSatish Balay           col[1].i = i - 1;
5639371c9d4SSatish Balay           col[2].j = j;
5649371c9d4SSatish Balay           col[2].i = i;
5659371c9d4SSatish Balay           col[3].j = j;
5669371c9d4SSatish Balay           col[3].i = i + 1;
5679371c9d4SSatish Balay           col[4].j = j + 1;
5689371c9d4SSatish Balay           col[4].i = i;
569c4762a1bSJed Brown           v[0]     = -hxdhy * newt_S + cross_EW;
570c4762a1bSJed Brown           v[1]     = -hydhx * newt_W + cross_NS;
571c4762a1bSJed Brown           v[2]     = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * PetscExpScalar(u);
572c4762a1bSJed Brown           v[3]     = -hydhx * newt_E - cross_NS;
573c4762a1bSJed Brown           v[4]     = -hxdhy * newt_N - cross_EW;
5749566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES));
575c4762a1bSJed Brown           break;
576c4762a1bSJed Brown         case JAC_NEWTON:
5770e3d61c9SBarry Smith           /* The Jacobian is
5780e3d61c9SBarry Smith 
5790e3d61c9SBarry Smith            -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u
5800e3d61c9SBarry Smith 
5810e3d61c9SBarry Smith           */
5829371c9d4SSatish Balay           col[0].j = j - 1;
5839371c9d4SSatish Balay           col[0].i = i - 1;
5849371c9d4SSatish Balay           col[1].j = j - 1;
5859371c9d4SSatish Balay           col[1].i = i;
5869371c9d4SSatish Balay           col[2].j = j - 1;
5879371c9d4SSatish Balay           col[2].i = i + 1;
5889371c9d4SSatish Balay           col[3].j = j;
5899371c9d4SSatish Balay           col[3].i = i - 1;
5909371c9d4SSatish Balay           col[4].j = j;
5919371c9d4SSatish Balay           col[4].i = i;
5929371c9d4SSatish Balay           col[5].j = j;
5939371c9d4SSatish Balay           col[5].i = i + 1;
5949371c9d4SSatish Balay           col[6].j = j + 1;
5959371c9d4SSatish Balay           col[6].i = i - 1;
5969371c9d4SSatish Balay           col[7].j = j + 1;
5979371c9d4SSatish Balay           col[7].i = i;
5989371c9d4SSatish Balay           col[8].j = j + 1;
5999371c9d4SSatish Balay           col[8].i = i + 1;
600c4762a1bSJed Brown           v[0]     = -0.25 * (skew_S + skew_W);
601c4762a1bSJed Brown           v[1]     = -hxdhy * newt_S + cross_EW;
602c4762a1bSJed Brown           v[2]     = 0.25 * (skew_S + skew_E);
603c4762a1bSJed Brown           v[3]     = -hydhx * newt_W + cross_NS;
604c4762a1bSJed Brown           v[4]     = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * PetscExpScalar(u);
605c4762a1bSJed Brown           v[5]     = -hydhx * newt_E - cross_NS;
606c4762a1bSJed Brown           v[6]     = 0.25 * (skew_N + skew_W);
607c4762a1bSJed Brown           v[7]     = -hxdhy * newt_N - cross_EW;
608c4762a1bSJed Brown           v[8]     = -0.25 * (skew_N + skew_E);
6099566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B, 1, &row, 9, col, v, INSERT_VALUES));
610c4762a1bSJed Brown           break;
611d71ae5a4SJacob Faibussowitsch         default:
612d71ae5a4SJacob Faibussowitsch           SETERRQ(PetscObjectComm((PetscObject)info->da), PETSC_ERR_SUP, "Jacobian type %d not implemented", user->jtype);
613c4762a1bSJed Brown         }
614c4762a1bSJed Brown       }
615c4762a1bSJed Brown     }
616c4762a1bSJed Brown   }
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   /*
619c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
620c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
621c4762a1bSJed Brown   */
6229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
6239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
624c4762a1bSJed Brown 
625c4762a1bSJed Brown   if (J != B) {
6269566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
6279566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
628c4762a1bSJed Brown   }
629c4762a1bSJed Brown 
630c4762a1bSJed Brown   /*
631c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
632c4762a1bSJed Brown      matrix. If we do, it will generate an error.
633c4762a1bSJed Brown   */
63448a46eb9SPierre Jolivet   if (user->jtype == JAC_NEWTON) PetscCall(PetscLogFlops(info->xm * info->ym * (131.0)));
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
636c4762a1bSJed Brown }
637c4762a1bSJed Brown 
638c4762a1bSJed Brown /***********************************************************
639c4762a1bSJed Brown  * PreCheck implementation
640c4762a1bSJed Brown  ***********************************************************/
PreCheckSetFromOptions(PreCheck precheck)641d71ae5a4SJacob Faibussowitsch PetscErrorCode PreCheckSetFromOptions(PreCheck precheck)
642d71ae5a4SJacob Faibussowitsch {
643c4762a1bSJed Brown   PetscBool flg;
644c4762a1bSJed Brown 
645c4762a1bSJed Brown   PetscFunctionBeginUser;
646d0609cedSBarry Smith   PetscOptionsBegin(precheck->comm, NULL, "PreCheck Options", "none");
6479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-precheck_angle", "Angle in degrees between successive search directions necessary to activate step correction", "", precheck->angle, &precheck->angle, NULL));
648c4762a1bSJed Brown   flg = PETSC_FALSE;
6499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-precheck_monitor", "Monitor choices made by precheck routine", "", flg, &flg, NULL));
65048a46eb9SPierre Jolivet   if (flg) PetscCall(PetscViewerASCIIOpen(precheck->comm, "stdout", &precheck->monitor));
651d0609cedSBarry Smith   PetscOptionsEnd();
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
653c4762a1bSJed Brown }
654c4762a1bSJed Brown 
655c4762a1bSJed Brown /*
656c4762a1bSJed Brown   Compare the direction of the current and previous step, modify the current step accordingly
657c4762a1bSJed Brown */
PreCheckFunction(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool * changed,PetscCtx ctx)658*2a8381b2SBarry Smith PetscErrorCode PreCheckFunction(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, PetscCtx ctx)
659d71ae5a4SJacob Faibussowitsch {
660c4762a1bSJed Brown   PreCheck    precheck;
661c4762a1bSJed Brown   Vec         Ylast;
662c4762a1bSJed Brown   PetscScalar dot;
663c4762a1bSJed Brown   PetscInt    iter;
664c4762a1bSJed Brown   PetscReal   ynorm, ylastnorm, theta, angle_radians;
665c4762a1bSJed Brown   SNES        snes;
666c4762a1bSJed Brown 
667c4762a1bSJed Brown   PetscFunctionBeginUser;
6689566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
669c4762a1bSJed Brown   precheck = (PreCheck)ctx;
6709566063dSJacob Faibussowitsch   if (!precheck->Ylast) PetscCall(VecDuplicate(Y, &precheck->Ylast));
671c4762a1bSJed Brown   Ylast = precheck->Ylast;
6729566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
673c4762a1bSJed Brown   if (iter < 1) {
6749566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
675c4762a1bSJed Brown     *changed = PETSC_FALSE;
6763ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
677c4762a1bSJed Brown   }
678c4762a1bSJed Brown 
6799566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
6809566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
6819566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
682c4762a1bSJed Brown   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
683c4762a1bSJed Brown   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
684c4762a1bSJed Brown   angle_radians = precheck->angle * PETSC_PI / 180.;
685c4762a1bSJed Brown   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
686c4762a1bSJed Brown     /* Modify the step Y */
687c4762a1bSJed Brown     PetscReal alpha, ydiffnorm;
6889566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
6899566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
690c4762a1bSJed Brown     alpha = ylastnorm / ydiffnorm;
6919566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
6929566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
69348a46eb9SPierre 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));
694c4762a1bSJed Brown   } else {
6959566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
696c4762a1bSJed Brown     *changed = PETSC_FALSE;
69748a46eb9SPierre 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));
698c4762a1bSJed Brown   }
6993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
700c4762a1bSJed Brown }
701c4762a1bSJed Brown 
PreCheckDestroy(PreCheck * precheck)702d71ae5a4SJacob Faibussowitsch PetscErrorCode PreCheckDestroy(PreCheck *precheck)
703d71ae5a4SJacob Faibussowitsch {
704c4762a1bSJed Brown   PetscFunctionBeginUser;
7053ba16761SJacob Faibussowitsch   if (!*precheck) PetscFunctionReturn(PETSC_SUCCESS);
7069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*precheck)->Ylast));
7079566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*precheck)->monitor));
7089566063dSJacob Faibussowitsch   PetscCall(PetscFree(*precheck));
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
710c4762a1bSJed Brown }
711c4762a1bSJed Brown 
PreCheckCreate(MPI_Comm comm,PreCheck * precheck)712d71ae5a4SJacob Faibussowitsch PetscErrorCode PreCheckCreate(MPI_Comm comm, PreCheck *precheck)
713d71ae5a4SJacob Faibussowitsch {
714c4762a1bSJed Brown   PetscFunctionBeginUser;
7159566063dSJacob Faibussowitsch   PetscCall(PetscNew(precheck));
716c4762a1bSJed Brown 
717c4762a1bSJed Brown   (*precheck)->comm  = comm;
718c4762a1bSJed Brown   (*precheck)->angle = 10.; /* only active if angle is less than 10 degrees */
7193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
720c4762a1bSJed Brown }
721c4762a1bSJed Brown 
722c4762a1bSJed Brown /*
723c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
724c4762a1bSJed Brown 
725c4762a1bSJed Brown  */
NonlinearGS(SNES snes,Vec X,Vec B,PetscCtx ctx)726*2a8381b2SBarry Smith PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, PetscCtx ctx)
727d71ae5a4SJacob Faibussowitsch {
728c4762a1bSJed Brown   PetscInt      i, j, k, xs, ys, xm, ym, its, tot_its, sweeps, l, m;
729c4762a1bSJed Brown   PetscReal     hx, hy, hxdhy, hydhx, dhx, dhy, sc;
730c4762a1bSJed Brown   PetscScalar **x, **b, bij, F, F0 = 0, J, y, u, eu;
731c4762a1bSJed Brown   PetscReal     atol, rtol, stol;
732c4762a1bSJed Brown   DM            da;
733c4762a1bSJed Brown   AppCtx       *user = (AppCtx *)ctx;
734c4762a1bSJed Brown   Vec           localX, localB;
735c4762a1bSJed Brown   DMDALocalInfo info;
736c4762a1bSJed Brown 
737c4762a1bSJed Brown   PetscFunctionBeginUser;
7389566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
7399566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
740c4762a1bSJed Brown 
741c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(info.mx - 1);
742c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(info.my - 1);
743c4762a1bSJed Brown   sc    = hx * hy * user->lambda;
744c4762a1bSJed Brown   dhx   = 1 / hx;
745c4762a1bSJed Brown   dhy   = 1 / hy;
746c4762a1bSJed Brown   hxdhy = hx / hy;
747c4762a1bSJed Brown   hydhx = hy / hx;
748c4762a1bSJed Brown 
749c4762a1bSJed Brown   tot_its = 0;
7509566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetSweeps(snes, &sweeps));
7519566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its));
7529566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
75348a46eb9SPierre Jolivet   if (B) PetscCall(DMGetLocalVector(da, &localB));
754c4762a1bSJed Brown   if (B) {
7559566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, localB));
7569566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, localB));
757c4762a1bSJed Brown   }
7589566063dSJacob Faibussowitsch   if (B) PetscCall(DMDAVecGetArrayRead(da, localB, &b));
7599566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
7609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
7619566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
762c4762a1bSJed Brown   for (l = 0; l < sweeps; l++) {
763c4762a1bSJed Brown     /*
764c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
765c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
766c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
767c4762a1bSJed Brown      */
7689566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
769c4762a1bSJed Brown     for (m = 0; m < 2; m++) {
770c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
771c4762a1bSJed Brown         for (i = xs + (m + j) % 2; i < xs + xm; i += 2) {
772c4762a1bSJed Brown           PetscReal xx = i * hx, yy = j * hy;
773c4762a1bSJed Brown           if (B) bij = b[j][i];
774c4762a1bSJed Brown           else bij = 0.;
775c4762a1bSJed Brown 
776c4762a1bSJed Brown           if (i == 0 || j == 0 || i == info.mx - 1 || j == info.my - 1) {
777c4762a1bSJed Brown             /* boundary conditions are all zero Dirichlet */
778c4762a1bSJed Brown             x[j][i] = 0.0 + bij;
779c4762a1bSJed Brown           } else {
7809371c9d4SSatish 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];
7819371c9d4SSatish 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]);
782c4762a1bSJed Brown             u = x[j][i];
783c4762a1bSJed Brown             for (k = 0; k < its; k++) {
7849371c9d4SSatish 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);
785c4762a1bSJed Brown 
786c4762a1bSJed Brown               if (sc) eu = PetscExpScalar(u);
787c4762a1bSJed Brown               else eu = 0;
788c4762a1bSJed Brown 
789c4762a1bSJed Brown               F = uxx + uyy - sc * eu - bij;
790c4762a1bSJed Brown               if (k == 0) F0 = F;
791c4762a1bSJed Brown               J = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * eu;
792c4762a1bSJed Brown               y = F / J;
793c4762a1bSJed Brown               u -= y;
794c4762a1bSJed Brown               tot_its++;
795ad540459SPierre Jolivet               if (atol > PetscAbsReal(PetscRealPart(F)) || rtol * PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || stol * PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) break;
796c4762a1bSJed Brown             }
797c4762a1bSJed Brown             x[j][i] = u;
798c4762a1bSJed Brown           }
799c4762a1bSJed Brown         }
800c4762a1bSJed Brown       }
801c4762a1bSJed Brown     }
802c4762a1bSJed Brown     /*
803c4762a1bSJed Brown x     Restore vector
804c4762a1bSJed Brown      */
805c4762a1bSJed Brown   }
8069566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
8079566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X));
8089566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X));
8099566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(tot_its * (118.0)));
8109566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
811c4762a1bSJed Brown   if (B) {
8129566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, localB, &b));
8139566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da, &localB));
814c4762a1bSJed Brown   }
8153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
816c4762a1bSJed Brown }
817c4762a1bSJed Brown 
818c4762a1bSJed Brown /*TEST
819c4762a1bSJed Brown 
820c4762a1bSJed Brown    test:
821c4762a1bSJed Brown       nsize: 2
822c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype NEWTON
823c4762a1bSJed Brown       requires: !single
824c4762a1bSJed Brown 
825c4762a1bSJed Brown    test:
826c4762a1bSJed Brown       suffix: 2
827c4762a1bSJed Brown       nsize: 2
828c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -precheck 1
829c4762a1bSJed Brown       requires: !single
830c4762a1bSJed Brown 
831c4762a1bSJed Brown    test:
832c4762a1bSJed Brown       suffix: 3
833c4762a1bSJed Brown       nsize: 2
83485b95db3SBarry 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
83585b95db3SBarry Smith       requires: !single
83685b95db3SBarry Smith 
83785b95db3SBarry Smith    test:
83885b95db3SBarry Smith       suffix: 3_star
83985b95db3SBarry Smith       nsize: 2
84085b95db3SBarry 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
84185b95db3SBarry Smith       output_file: output/ex15_3.out
842c4762a1bSJed Brown       requires: !single
843c4762a1bSJed Brown 
844c4762a1bSJed Brown    test:
845c4762a1bSJed Brown       suffix: 4
846c4762a1bSJed 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
847c4762a1bSJed Brown       requires: !single
848c4762a1bSJed Brown 
849c4762a1bSJed Brown    test:
8504b0a5c37SStefano Zampini       suffix: 5
8514b0a5c37SStefano Zampini       args: -snes_monitor_short -snes_type newtontr -npc_snes_type ngs -snes_npc_side right -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -ksp_converged_reason -pc_type ilu
8524b0a5c37SStefano Zampini       requires: !single
8534b0a5c37SStefano Zampini 
8544b0a5c37SStefano Zampini    test:
855c4762a1bSJed Brown       suffix: lag_jac
856c4762a1bSJed Brown       nsize: 4
857c4762a1bSJed 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
858c4762a1bSJed Brown       requires: !single
859c4762a1bSJed Brown 
860c4762a1bSJed Brown    test:
861c4762a1bSJed Brown       suffix: lag_pc
862c4762a1bSJed Brown       nsize: 4
863c4762a1bSJed 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
864c4762a1bSJed Brown       requires: !single
865c4762a1bSJed Brown 
866c4762a1bSJed Brown    test:
867c4762a1bSJed Brown       suffix: nleqerr
868c4762a1bSJed 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
869c4762a1bSJed Brown       requires: !single
870c4762a1bSJed Brown 
871bbc1464cSBarry Smith    test:
872bbc1464cSBarry Smith       suffix: mf
873bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_mf_operator
874bbc1464cSBarry Smith       requires: !single
875bbc1464cSBarry Smith 
876bbc1464cSBarry Smith    test:
877bbc1464cSBarry Smith       suffix: mf_picard
878bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_mf_operator -picard
879bbc1464cSBarry Smith       requires: !single
880bbc1464cSBarry Smith       output_file: output/ex15_mf.out
881bbc1464cSBarry Smith 
88285b95db3SBarry Smith    test:
88385b95db3SBarry Smith       suffix: fd_picard
88485b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 2 -p 3 -ksp_rtol 1.e-12 -snes_fd -picard
88585b95db3SBarry Smith       requires: !single
88685b95db3SBarry Smith 
88785b95db3SBarry Smith    test:
88885b95db3SBarry Smith       suffix: fd_color_picard
88985b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_fd_color -picard
89085b95db3SBarry Smith       requires: !single
89185b95db3SBarry Smith       output_file: output/ex15_mf.out
89285b95db3SBarry Smith 
893c4762a1bSJed Brown TEST*/
894