xref: /petsc/src/snes/tutorials/ex15.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
53c4762a1bSJed Brown typedef enum {JAC_BRATU,JAC_PICARD,JAC_STAR,JAC_NEWTON} JacType;
54c4762a1bSJed Brown static const char *const JacTypes[] = {"BRATU","PICARD","STAR","NEWTON","JacType","JAC_",0};
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*
57c4762a1bSJed Brown    User-defined application context - contains data needed by the
58c4762a1bSJed Brown    application-provided call-back routines, FormJacobianLocal() and
59c4762a1bSJed Brown    FormFunctionLocal().
60c4762a1bSJed Brown */
61c4762a1bSJed Brown typedef struct {
62c4762a1bSJed Brown   PetscReal   lambda;         /* Bratu parameter */
63c4762a1bSJed Brown   PetscReal   p;              /* Exponent in p-Laplacian */
64c4762a1bSJed Brown   PetscReal   epsilon;        /* Regularization */
65c4762a1bSJed Brown   PetscReal   source;         /* Source term */
66c4762a1bSJed Brown   JacType     jtype;          /* What type of Jacobian to assemble */
67c4762a1bSJed Brown   PetscBool   picard;
68c4762a1bSJed Brown   PetscInt    blocks[2];
69c4762a1bSJed Brown   PetscReal   kappa;
70c4762a1bSJed Brown   PetscInt    initial;        /* initial conditions type */
71c4762a1bSJed Brown } AppCtx;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*
74c4762a1bSJed Brown    User-defined routines
75c4762a1bSJed Brown */
76c4762a1bSJed Brown static PetscErrorCode FormRHS(AppCtx*,DM,Vec);
77c4762a1bSJed Brown static PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
78c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
79c4762a1bSJed Brown static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
80c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,Mat,AppCtx*);
81c4762a1bSJed Brown static PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);
82c4762a1bSJed Brown 
83c4762a1bSJed Brown typedef struct _n_PreCheck *PreCheck;
84c4762a1bSJed Brown struct _n_PreCheck {
85c4762a1bSJed Brown   MPI_Comm    comm;
86c4762a1bSJed Brown   PetscReal   angle;
87c4762a1bSJed Brown   Vec         Ylast;
88c4762a1bSJed Brown   PetscViewer monitor;
89c4762a1bSJed Brown };
90c4762a1bSJed Brown PetscErrorCode PreCheckCreate(MPI_Comm,PreCheck*);
91c4762a1bSJed Brown PetscErrorCode PreCheckDestroy(PreCheck*);
92c4762a1bSJed Brown PetscErrorCode PreCheckFunction(SNESLineSearch,Vec,Vec,PetscBool*,void*);
93c4762a1bSJed Brown PetscErrorCode PreCheckSetFromOptions(PreCheck);
94c4762a1bSJed Brown 
95c4762a1bSJed Brown int main(int argc,char **argv)
96c4762a1bSJed Brown {
97c4762a1bSJed Brown   SNES                snes;                    /* nonlinear solver */
98c4762a1bSJed Brown   Vec                 x,r,b;                   /* solution, residual, rhs vectors */
99c4762a1bSJed Brown   AppCtx              user;                    /* user-defined work context */
100c4762a1bSJed Brown   PetscInt            its;                     /* iterations for convergence */
101c4762a1bSJed Brown   SNESConvergedReason reason;                  /* Check convergence */
102c4762a1bSJed Brown   PetscBool           alloc_star;              /* Only allocate for the STAR stencil  */
103c4762a1bSJed Brown   PetscBool           write_output;
104c4762a1bSJed Brown   char                filename[PETSC_MAX_PATH_LEN] = "ex15.vts";
105c4762a1bSJed Brown   PetscReal           bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
10685b95db3SBarry Smith   DM                  da;                      /* distributed array data structure */
107c4762a1bSJed Brown   PreCheck            precheck = NULL;         /* precheck context for version in this file */
108c4762a1bSJed Brown   PetscInt            use_precheck;            /* 0=none, 1=version in this file, 2=SNES-provided version */
109c4762a1bSJed Brown   PetscReal           precheck_angle;          /* When manually setting the SNES-provided precheck function */
110c4762a1bSJed Brown   SNESLineSearch      linesearch;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown      Initialize program
114c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119c4762a1bSJed Brown      Initialize problem parameters
120c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121c4762a1bSJed Brown   user.lambda    = 0.0; user.p = 2.0; user.epsilon = 1e-5; user.source = 0.1; user.jtype = JAC_NEWTON;user.initial=-1;
122c4762a1bSJed Brown   user.blocks[0] = 1; user.blocks[1] = 1; user.kappa = 1e-3;
123c4762a1bSJed Brown   alloc_star     = PETSC_FALSE;
124c4762a1bSJed Brown   use_precheck   = 0; precheck_angle = 10.;
125c4762a1bSJed Brown   user.picard    = PETSC_FALSE;
126*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"p-Bratu options",__FILE__);
127c4762a1bSJed Brown   {
128c4762a1bSJed Brown     PetscInt two=2;
1299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda","Bratu parameter","",user.lambda,&user.lambda,NULL));
1309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-p","Exponent `p' in p-Laplacian","",user.p,&user.p,NULL));
1319566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-epsilon","Strain-regularization in p-Laplacian","",user.epsilon,&user.epsilon,NULL));
1329566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-source","Constant source term","",user.source,&user.source,NULL));
1339566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-jtype","Jacobian approximation to assemble","",JacTypes,(PetscEnum)user.jtype,(PetscEnum*)&user.jtype,NULL));
1349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-picard","Solve with defect-correction Picard iteration","",&user.picard));
13585b95db3SBarry Smith     if (user.picard) {
13685b95db3SBarry Smith       user.jtype = JAC_PICARD;
137e00437b9SBarry Smith       PetscCheck(user.p == 3,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Picard iteration is only supported for p == 3");
13885b95db3SBarry Smith       /* the Picard linearization only requires a 5 point stencil, while the Newton linearization requires a 9 point stencil */
13985b95db3SBarry Smith       /* hence allocating the 5 point stencil gives the same convergence as the 9 point stencil since the extra stencil points are not used */
1409566063dSJacob Faibussowitsch       PetscCall(PetscOptionsBool("-alloc_star","Allocate for STAR stencil (5-point)","",alloc_star,&alloc_star,NULL));
14185b95db3SBarry Smith     }
1429566063dSJacob 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));
1439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-initial","Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)","",user.initial,&user.initial,NULL));
144c4762a1bSJed Brown     if (use_precheck == 2) {    /* Using library version, get the angle */
1459566063dSJacob Faibussowitsch       PetscCall(PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck_angle,&precheck_angle,NULL));
146c4762a1bSJed Brown     }
1479566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-blocks","number of coefficient interfaces in x and y direction","",user.blocks,&two,NULL));
1489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-kappa","diffusivity in odd regions","",user.kappa,&user.kappa,NULL));
1499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-o","Output solution in vts format","",filename,filename,sizeof(filename),&write_output));
150c4762a1bSJed Brown   }
151*d0609cedSBarry Smith   PetscOptionsEnd();
152c4762a1bSJed Brown   if (user.lambda > bratu_lambda_max || user.lambda < bratu_lambda_min) {
1539566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING: lambda %g out of range for p=2\n",(double)user.lambda));
154c4762a1bSJed Brown   }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157c4762a1bSJed Brown      Create nonlinear solver context
158c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1599566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
163c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1649566063dSJacob 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));
1659566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1669566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169c4762a1bSJed Brown      Extract global vectors from DM; then duplicate for remaining
170c4762a1bSJed Brown      vectors that are the same types
171c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1729566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&x));
1739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
1749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&b));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
17785b95db3SBarry Smith      User can override with:
178c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
179c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
180c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
181c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
182c4762a1bSJed Brown                          products within Newton-Krylov method
183c4762a1bSJed Brown 
184c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown      Set local function evaluation routine
188c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1899566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
1909566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes,da));
191c4762a1bSJed Brown   if (user.picard) {
192c4762a1bSJed Brown     /*
193c4762a1bSJed Brown         This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access
194c4762a1bSJed Brown         the SNES to set it
195c4762a1bSJed Brown     */
196*d0609cedSBarry Smith     PetscCall(DMDASNESSetPicardLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionPicardLocal,
197*d0609cedSBarry Smith                                      (PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user));
1989566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user));
1999566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes,NULL,NULL,SNESPicardComputeJacobian,&user));
200c4762a1bSJed Brown   } else {
2019566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user));
2029566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetJacobianLocal(da,(PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user));
203c4762a1bSJed Brown   }
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
207c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2089566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
2099566063dSJacob Faibussowitsch   PetscCall(SNESSetNGS(snes,NonlinearGS,&user));
2109566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
211c4762a1bSJed Brown   /* Set up the precheck context if requested */
212c4762a1bSJed Brown   if (use_precheck == 1) {      /* Use the precheck routines in this file */
2139566063dSJacob Faibussowitsch     PetscCall(PreCheckCreate(PETSC_COMM_WORLD,&precheck));
2149566063dSJacob Faibussowitsch     PetscCall(PreCheckSetFromOptions(precheck));
2159566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch,PreCheckFunction,precheck));
216c4762a1bSJed Brown   } else if (use_precheck == 2) { /* Use the version provided by the library */
2179566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&precheck_angle));
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221c4762a1bSJed Brown      Evaluate initial guess
222c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
223c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
224c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
225c4762a1bSJed Brown      this vector to zero by calling VecSet().
226c4762a1bSJed Brown   */
227c4762a1bSJed Brown 
2289566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user,da,x));
2299566063dSJacob Faibussowitsch   PetscCall(FormRHS(&user,da,b));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232c4762a1bSJed Brown      Solve nonlinear system
233c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2349566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,b,x));
2359566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
2369566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes,&reason));
237c4762a1bSJed Brown 
2389566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s Number of nonlinear iterations = %D\n",SNESConvergedReasons[reason],its));
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   if (write_output) {
241c4762a1bSJed Brown     PetscViewer viewer;
2429566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer));
2439566063dSJacob Faibussowitsch     PetscCall(VecView(x,viewer));
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
245c4762a1bSJed Brown   }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
249c4762a1bSJed Brown      are no longer needed.
250c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
251c4762a1bSJed Brown 
2529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
2559566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2579566063dSJacob Faibussowitsch   PetscCall(PreCheckDestroy(&precheck));
2589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
259b122ec5aSJacob Faibussowitsch   return 0;
260c4762a1bSJed Brown }
261c4762a1bSJed Brown 
262c4762a1bSJed Brown /* ------------------------------------------------------------------- */
263c4762a1bSJed Brown /*
264c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
265c4762a1bSJed Brown 
266c4762a1bSJed Brown    Input Parameters:
267c4762a1bSJed Brown    user - user-defined application context
268c4762a1bSJed Brown    X - vector
269c4762a1bSJed Brown 
270c4762a1bSJed Brown    Output Parameter:
271c4762a1bSJed Brown    X - vector
272c4762a1bSJed Brown  */
273c4762a1bSJed Brown static PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
274c4762a1bSJed Brown {
275c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
276c4762a1bSJed Brown   PetscReal      temp1,temp,hx,hy;
277c4762a1bSJed Brown   PetscScalar    **x;
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   PetscFunctionBeginUser;
2809566063dSJacob 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));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(Mx-1);
283c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(My-1);
284c4762a1bSJed Brown   temp1 = user->lambda / (user->lambda + 1.);
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /*
287c4762a1bSJed Brown      Get a pointer to vector data.
288c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
289c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
290c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
291c4762a1bSJed Brown          the array.
292c4762a1bSJed Brown   */
2939566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /*
296c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DA):
297c4762a1bSJed Brown        xs, ys   - starting grid indices (no ghost points)
298c4762a1bSJed Brown        xm, ym   - widths of local grid (no ghost points)
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   */
3019566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   /*
304c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
305c4762a1bSJed Brown   */
306c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
307c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
308c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
309c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
310c4762a1bSJed Brown         /* boundary conditions are all zero Dirichlet */
311c4762a1bSJed Brown         x[j][i] = 0.0;
312c4762a1bSJed Brown       } else {
313c4762a1bSJed Brown         if (user->initial == -1) {
314c4762a1bSJed Brown           if (user->lambda != 0) {
315c4762a1bSJed Brown             x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
316c4762a1bSJed Brown           } else {
317c4762a1bSJed Brown             /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting
318c4762a1bSJed Brown              * with an exact solution. */
319c4762a1bSJed Brown             const PetscReal
320c4762a1bSJed Brown               xx = 2*(PetscReal)i/(Mx-1) - 1,
321c4762a1bSJed Brown               yy = 2*(PetscReal)j/(My-1) - 1;
322c4762a1bSJed Brown             x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
323c4762a1bSJed Brown           }
324c4762a1bSJed Brown         } else if (user->initial == 0) {
325c4762a1bSJed Brown           x[j][i] = 0.;
326c4762a1bSJed Brown         } else if (user->initial == 1) {
327c4762a1bSJed Brown           const PetscReal
328c4762a1bSJed Brown             xx = 2*(PetscReal)i/(Mx-1) - 1,
329c4762a1bSJed Brown             yy = 2*(PetscReal)j/(My-1) - 1;
330c4762a1bSJed Brown           x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
331c4762a1bSJed Brown         } else {
332c4762a1bSJed Brown           if (user->lambda != 0) {
333c4762a1bSJed Brown             x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
334c4762a1bSJed Brown           } else {
335c4762a1bSJed Brown             x[j][i] = 0.5*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
336c4762a1bSJed Brown           }
337c4762a1bSJed Brown         }
338c4762a1bSJed Brown       }
339c4762a1bSJed Brown     }
340c4762a1bSJed Brown   }
341c4762a1bSJed Brown   /*
342c4762a1bSJed Brown      Restore vector
343c4762a1bSJed Brown   */
3449566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,X,&x));
345c4762a1bSJed Brown   PetscFunctionReturn(0);
346c4762a1bSJed Brown }
347c4762a1bSJed Brown 
348c4762a1bSJed Brown /*
349c4762a1bSJed Brown    FormRHS - Forms constant RHS for the problem.
350c4762a1bSJed Brown 
351c4762a1bSJed Brown    Input Parameters:
352c4762a1bSJed Brown    user - user-defined application context
353c4762a1bSJed Brown    B - RHS vector
354c4762a1bSJed Brown 
355c4762a1bSJed Brown    Output Parameter:
356c4762a1bSJed Brown    B - vector
357c4762a1bSJed Brown  */
358c4762a1bSJed Brown static PetscErrorCode FormRHS(AppCtx *user,DM da,Vec B)
359c4762a1bSJed Brown {
360c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
361c4762a1bSJed Brown   PetscReal      hx,hy;
362c4762a1bSJed Brown   PetscScalar    **b;
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   PetscFunctionBeginUser;
3659566063dSJacob 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));
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(Mx-1);
368c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(My-1);
3699566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,B,&b));
3709566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
371c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
372c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
373c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
374c4762a1bSJed Brown         b[j][i] = 0.0;
375c4762a1bSJed Brown       } else {
376c4762a1bSJed Brown         b[j][i] = hx*hy*user->source;
377c4762a1bSJed Brown       }
378c4762a1bSJed Brown     }
379c4762a1bSJed Brown   }
3809566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,B,&b));
381c4762a1bSJed Brown   PetscFunctionReturn(0);
382c4762a1bSJed Brown }
383c4762a1bSJed Brown 
3849fbee547SJacob Faibussowitsch static inline PetscReal kappa(const AppCtx *ctx,PetscReal x,PetscReal y)
385c4762a1bSJed Brown {
386c4762a1bSJed Brown   return (((PetscInt)(x*ctx->blocks[0])) + ((PetscInt)(y*ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0;
387c4762a1bSJed Brown }
388c4762a1bSJed Brown /* p-Laplacian diffusivity */
3899fbee547SJacob Faibussowitsch static inline PetscScalar eta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
390c4762a1bSJed Brown {
391c4762a1bSJed Brown   return kappa(ctx,x,y) * PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-2.));
392c4762a1bSJed Brown }
3939fbee547SJacob Faibussowitsch static inline PetscScalar deta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
394c4762a1bSJed Brown {
395c4762a1bSJed Brown   return (ctx->p == 2)
396c4762a1bSJed Brown          ? 0
397c4762a1bSJed Brown          : kappa(ctx,x,y)*PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-4)) * 0.5 * (ctx->p-2.);
398c4762a1bSJed Brown }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown /* ------------------------------------------------------------------- */
401c4762a1bSJed Brown /*
402c4762a1bSJed Brown    FormFunctionLocal - Evaluates nonlinear function, F(x).
403c4762a1bSJed Brown  */
404c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
405c4762a1bSJed Brown {
406c4762a1bSJed Brown   PetscReal      hx,hy,dhx,dhy,sc;
407c4762a1bSJed Brown   PetscInt       i,j;
408c4762a1bSJed Brown   PetscScalar    eu;
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   PetscFunctionBeginUser;
411c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info->mx-1);
412c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info->my-1);
413c4762a1bSJed Brown   sc     = hx*hy*user->lambda;
414c4762a1bSJed Brown   dhx    = 1/hx;
415c4762a1bSJed Brown   dhy    = 1/hy;
416c4762a1bSJed Brown   /*
417c4762a1bSJed Brown      Compute function over the locally owned part of the grid
418c4762a1bSJed Brown   */
419c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
420c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
421c4762a1bSJed Brown       PetscReal xx = i*hx,yy = j*hy;
422c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
423c4762a1bSJed Brown         f[j][i] = x[j][i];
424c4762a1bSJed Brown       } else {
425c4762a1bSJed Brown         const PetscScalar
426c4762a1bSJed Brown           u    = x[j][i],
427c4762a1bSJed Brown           ux_E = dhx*(x[j][i+1]-x[j][i]),
428c4762a1bSJed Brown           uy_E = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
429c4762a1bSJed Brown           ux_W = dhx*(x[j][i]-x[j][i-1]),
430c4762a1bSJed Brown           uy_W = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
431c4762a1bSJed Brown           ux_N = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
432c4762a1bSJed Brown           uy_N = dhy*(x[j+1][i]-x[j][i]),
433c4762a1bSJed Brown           ux_S = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
434c4762a1bSJed Brown           uy_S = dhy*(x[j][i]-x[j-1][i]),
435c4762a1bSJed Brown           e_E  = eta(user,xx,yy,ux_E,uy_E),
436c4762a1bSJed Brown           e_W  = eta(user,xx,yy,ux_W,uy_W),
437c4762a1bSJed Brown           e_N  = eta(user,xx,yy,ux_N,uy_N),
438c4762a1bSJed Brown           e_S  = eta(user,xx,yy,ux_S,uy_S),
439c4762a1bSJed Brown           uxx  = -hy * (e_E*ux_E - e_W*ux_W),
440c4762a1bSJed Brown           uyy  = -hx * (e_N*uy_N - e_S*uy_S);
441c4762a1bSJed Brown         if (sc) eu = PetscExpScalar(u);
442c4762a1bSJed Brown         else    eu = 0.;
4430e3d61c9SBarry Smith         /* For p=2, these terms decay to:
4440e3d61c9SBarry Smith          uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx
4450e3d61c9SBarry Smith          uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy
4460e3d61c9SBarry Smith         */
447c4762a1bSJed Brown         f[j][i] = uxx + uyy - sc*eu;
448c4762a1bSJed Brown       }
449c4762a1bSJed Brown     }
450c4762a1bSJed Brown   }
4519566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm*info->ym*(72.0)));
452c4762a1bSJed Brown   PetscFunctionReturn(0);
453c4762a1bSJed Brown }
454c4762a1bSJed Brown 
455c4762a1bSJed Brown /*
456c4762a1bSJed Brown     This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation,
457c4762a1bSJed 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)
458c4762a1bSJed Brown 
459c4762a1bSJed Brown */
460c4762a1bSJed Brown static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
461c4762a1bSJed Brown {
462c4762a1bSJed Brown   PetscReal hx,hy,sc;
463c4762a1bSJed Brown   PetscInt  i,j;
464c4762a1bSJed Brown 
465c4762a1bSJed Brown   PetscFunctionBeginUser;
466c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info->mx-1);
467c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info->my-1);
468c4762a1bSJed Brown   sc     = hx*hy*user->lambda;
469c4762a1bSJed Brown   /*
470c4762a1bSJed Brown      Compute function over the locally owned part of the grid
471c4762a1bSJed Brown   */
472c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
473c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
474c4762a1bSJed Brown       if (!(i == 0 || j == 0 || i == info->mx-1 || j == info->my-1)) {
475c4762a1bSJed Brown         const PetscScalar u = x[j][i];
476c4762a1bSJed Brown         f[j][i] = sc*PetscExpScalar(u);
4770df40c35SBarry Smith       } else {
4780df40c35SBarry Smith         f[j][i] = 0.0; /* this is zero because the A(x) x term forces the x to be zero on the boundary */
479c4762a1bSJed Brown       }
480c4762a1bSJed Brown     }
481c4762a1bSJed Brown   }
4829566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm*info->ym*2.0));
483c4762a1bSJed Brown   PetscFunctionReturn(0);
484c4762a1bSJed Brown }
485c4762a1bSJed Brown 
486c4762a1bSJed Brown /*
487c4762a1bSJed Brown    FormJacobianLocal - Evaluates Jacobian matrix.
488c4762a1bSJed Brown */
489c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat J,Mat B,AppCtx *user)
490c4762a1bSJed Brown {
491c4762a1bSJed Brown   PetscInt       i,j;
492c4762a1bSJed Brown   MatStencil     col[9],row;
493c4762a1bSJed Brown   PetscScalar    v[9];
494c4762a1bSJed Brown   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   PetscFunctionBeginUser;
4979566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
498c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(info->mx-1);
499c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(info->my-1);
500c4762a1bSJed Brown   sc    = hx*hy*user->lambda;
501c4762a1bSJed Brown   dhx   = 1/hx;
502c4762a1bSJed Brown   dhy   = 1/hy;
503c4762a1bSJed Brown   hxdhy = hx/hy;
504c4762a1bSJed Brown   hydhx = hy/hx;
505c4762a1bSJed Brown 
506c4762a1bSJed Brown   /*
507c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
508c4762a1bSJed Brown       - PETSc parallel matrix formats are partitioned by
509c4762a1bSJed Brown         contiguous chunks of rows across the processors.
510c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
511c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
512c4762a1bSJed Brown         appropriate processor during matrix assembly).
513c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
514c4762a1bSJed Brown   */
515c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
516c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
517c4762a1bSJed Brown       PetscReal xx = i*hx,yy = j*hy;
518c4762a1bSJed Brown       row.j = j; row.i = i;
519c4762a1bSJed Brown       /* boundary points */
520c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
521c4762a1bSJed Brown         v[0] = 1.0;
5229566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES));
523c4762a1bSJed Brown       } else {
524c4762a1bSJed Brown         /* interior grid points */
525c4762a1bSJed Brown         const PetscScalar
526c4762a1bSJed Brown           ux_E     = dhx*(x[j][i+1]-x[j][i]),
527c4762a1bSJed Brown           uy_E     = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
528c4762a1bSJed Brown           ux_W     = dhx*(x[j][i]-x[j][i-1]),
529c4762a1bSJed Brown           uy_W     = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
530c4762a1bSJed Brown           ux_N     = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
531c4762a1bSJed Brown           uy_N     = dhy*(x[j+1][i]-x[j][i]),
532c4762a1bSJed Brown           ux_S     = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
533c4762a1bSJed Brown           uy_S     = dhy*(x[j][i]-x[j-1][i]),
534c4762a1bSJed Brown           u        = x[j][i],
535c4762a1bSJed Brown           e_E      = eta(user,xx,yy,ux_E,uy_E),
536c4762a1bSJed Brown           e_W      = eta(user,xx,yy,ux_W,uy_W),
537c4762a1bSJed Brown           e_N      = eta(user,xx,yy,ux_N,uy_N),
538c4762a1bSJed Brown           e_S      = eta(user,xx,yy,ux_S,uy_S),
539c4762a1bSJed Brown           de_E     = deta(user,xx,yy,ux_E,uy_E),
540c4762a1bSJed Brown           de_W     = deta(user,xx,yy,ux_W,uy_W),
541c4762a1bSJed Brown           de_N     = deta(user,xx,yy,ux_N,uy_N),
542c4762a1bSJed Brown           de_S     = deta(user,xx,yy,ux_S,uy_S),
543c4762a1bSJed Brown           skew_E   = de_E*ux_E*uy_E,
544c4762a1bSJed Brown           skew_W   = de_W*ux_W*uy_W,
545c4762a1bSJed Brown           skew_N   = de_N*ux_N*uy_N,
546c4762a1bSJed Brown           skew_S   = de_S*ux_S*uy_S,
547c4762a1bSJed Brown           cross_EW = 0.25*(skew_E - skew_W),
548c4762a1bSJed Brown           cross_NS = 0.25*(skew_N - skew_S),
549c4762a1bSJed Brown           newt_E   = e_E+de_E*PetscSqr(ux_E),
550c4762a1bSJed Brown           newt_W   = e_W+de_W*PetscSqr(ux_W),
551c4762a1bSJed Brown           newt_N   = e_N+de_N*PetscSqr(uy_N),
552c4762a1bSJed Brown           newt_S   = e_S+de_S*PetscSqr(uy_S);
553c4762a1bSJed Brown         /* interior grid points */
554c4762a1bSJed Brown         switch (user->jtype) {
555c4762a1bSJed Brown         case JAC_BRATU:
556c4762a1bSJed Brown           /* Jacobian from p=2 */
557c4762a1bSJed Brown           v[0] = -hxdhy;                                           col[0].j = j-1;   col[0].i = i;
558c4762a1bSJed Brown           v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
559c4762a1bSJed Brown           v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(u);       col[2].j = row.j; col[2].i = row.i;
560c4762a1bSJed Brown           v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
561c4762a1bSJed Brown           v[4] = -hxdhy;                                           col[4].j = j+1;   col[4].i = i;
5629566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
563c4762a1bSJed Brown           break;
564c4762a1bSJed Brown         case JAC_PICARD:
565c4762a1bSJed Brown           /* Jacobian arising from Picard linearization */
566c4762a1bSJed Brown           v[0] = -hxdhy*e_S;                                           col[0].j = j-1;   col[0].i = i;
567c4762a1bSJed Brown           v[1] = -hydhx*e_W;                                           col[1].j = j;     col[1].i = i-1;
568c4762a1bSJed Brown           v[2] = (e_W+e_E)*hydhx + (e_S+e_N)*hxdhy;                    col[2].j = row.j; col[2].i = row.i;
569c4762a1bSJed Brown           v[3] = -hydhx*e_E;                                           col[3].j = j;     col[3].i = i+1;
570c4762a1bSJed Brown           v[4] = -hxdhy*e_N;                                           col[4].j = j+1;   col[4].i = i;
5719566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
572c4762a1bSJed Brown           break;
573c4762a1bSJed Brown         case JAC_STAR:
574c4762a1bSJed Brown           /* Full Jacobian, but only a star stencil */
575c4762a1bSJed Brown           col[0].j = j-1; col[0].i = i;
576c4762a1bSJed Brown           col[1].j = j;   col[1].i = i-1;
577c4762a1bSJed Brown           col[2].j = j;   col[2].i = i;
578c4762a1bSJed Brown           col[3].j = j;   col[3].i = i+1;
579c4762a1bSJed Brown           col[4].j = j+1; col[4].i = i;
580c4762a1bSJed Brown           v[0]     = -hxdhy*newt_S + cross_EW;
581c4762a1bSJed Brown           v[1]     = -hydhx*newt_W + cross_NS;
582c4762a1bSJed Brown           v[2]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
583c4762a1bSJed Brown           v[3]     = -hydhx*newt_E - cross_NS;
584c4762a1bSJed Brown           v[4]     = -hxdhy*newt_N - cross_EW;
5859566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
586c4762a1bSJed Brown           break;
587c4762a1bSJed Brown         case JAC_NEWTON:
5880e3d61c9SBarry Smith           /* The Jacobian is
5890e3d61c9SBarry Smith 
5900e3d61c9SBarry Smith            -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u
5910e3d61c9SBarry Smith 
5920e3d61c9SBarry Smith           */
593c4762a1bSJed Brown           col[0].j = j-1; col[0].i = i-1;
594c4762a1bSJed Brown           col[1].j = j-1; col[1].i = i;
595c4762a1bSJed Brown           col[2].j = j-1; col[2].i = i+1;
596c4762a1bSJed Brown           col[3].j = j;   col[3].i = i-1;
597c4762a1bSJed Brown           col[4].j = j;   col[4].i = i;
598c4762a1bSJed Brown           col[5].j = j;   col[5].i = i+1;
599c4762a1bSJed Brown           col[6].j = j+1; col[6].i = i-1;
600c4762a1bSJed Brown           col[7].j = j+1; col[7].i = i;
601c4762a1bSJed Brown           col[8].j = j+1; col[8].i = i+1;
602c4762a1bSJed Brown           v[0]     = -0.25*(skew_S + skew_W);
603c4762a1bSJed Brown           v[1]     = -hxdhy*newt_S + cross_EW;
604c4762a1bSJed Brown           v[2]     =  0.25*(skew_S + skew_E);
605c4762a1bSJed Brown           v[3]     = -hydhx*newt_W + cross_NS;
606c4762a1bSJed Brown           v[4]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
607c4762a1bSJed Brown           v[5]     = -hydhx*newt_E - cross_NS;
608c4762a1bSJed Brown           v[6]     =  0.25*(skew_N + skew_W);
609c4762a1bSJed Brown           v[7]     = -hxdhy*newt_N - cross_EW;
610c4762a1bSJed Brown           v[8]     = -0.25*(skew_N + skew_E);
6119566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,9,col,v,INSERT_VALUES));
612c4762a1bSJed Brown           break;
613c4762a1bSJed Brown         default:
61498921bdaSJacob Faibussowitsch           SETERRQ(PetscObjectComm((PetscObject)info->da),PETSC_ERR_SUP,"Jacobian type %d not implemented",user->jtype);
615c4762a1bSJed Brown         }
616c4762a1bSJed Brown       }
617c4762a1bSJed Brown     }
618c4762a1bSJed Brown   }
619c4762a1bSJed Brown 
620c4762a1bSJed Brown   /*
621c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
622c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
623c4762a1bSJed Brown   */
6249566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
6259566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
626c4762a1bSJed Brown 
627c4762a1bSJed Brown   if (J != B) {
6289566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
6299566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
630c4762a1bSJed Brown   }
631c4762a1bSJed Brown 
632c4762a1bSJed Brown   /*
633c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
634c4762a1bSJed Brown      matrix. If we do, it will generate an error.
635c4762a1bSJed Brown   */
636c4762a1bSJed Brown   if (user->jtype == JAC_NEWTON) {
6379566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(info->xm*info->ym*(131.0)));
638c4762a1bSJed Brown   }
639c4762a1bSJed Brown   PetscFunctionReturn(0);
640c4762a1bSJed Brown }
641c4762a1bSJed Brown 
642c4762a1bSJed Brown /***********************************************************
643c4762a1bSJed Brown  * PreCheck implementation
644c4762a1bSJed Brown  ***********************************************************/
645c4762a1bSJed Brown PetscErrorCode PreCheckSetFromOptions(PreCheck precheck)
646c4762a1bSJed Brown {
647c4762a1bSJed Brown   PetscBool      flg;
648c4762a1bSJed Brown 
649c4762a1bSJed Brown   PetscFunctionBeginUser;
650*d0609cedSBarry Smith   PetscOptionsBegin(precheck->comm,NULL,"PreCheck Options","none");
6519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck->angle,&precheck->angle,NULL));
652c4762a1bSJed Brown   flg  = PETSC_FALSE;
6539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-precheck_monitor","Monitor choices made by precheck routine","",flg,&flg,NULL));
654c4762a1bSJed Brown   if (flg) {
6559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(precheck->comm,"stdout",&precheck->monitor));
656c4762a1bSJed Brown   }
657*d0609cedSBarry Smith   PetscOptionsEnd();
658c4762a1bSJed Brown   PetscFunctionReturn(0);
659c4762a1bSJed Brown }
660c4762a1bSJed Brown 
661c4762a1bSJed Brown /*
662c4762a1bSJed Brown   Compare the direction of the current and previous step, modify the current step accordingly
663c4762a1bSJed Brown */
664c4762a1bSJed Brown PetscErrorCode PreCheckFunction(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed, void *ctx)
665c4762a1bSJed Brown {
666c4762a1bSJed Brown   PreCheck       precheck;
667c4762a1bSJed Brown   Vec            Ylast;
668c4762a1bSJed Brown   PetscScalar    dot;
669c4762a1bSJed Brown   PetscInt       iter;
670c4762a1bSJed Brown   PetscReal      ynorm,ylastnorm,theta,angle_radians;
671c4762a1bSJed Brown   SNES           snes;
672c4762a1bSJed Brown 
673c4762a1bSJed Brown   PetscFunctionBeginUser;
6749566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
675c4762a1bSJed Brown   precheck = (PreCheck)ctx;
6769566063dSJacob Faibussowitsch   if (!precheck->Ylast) PetscCall(VecDuplicate(Y,&precheck->Ylast));
677c4762a1bSJed Brown   Ylast = precheck->Ylast;
6789566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&iter));
679c4762a1bSJed Brown   if (iter < 1) {
6809566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
681c4762a1bSJed Brown     *changed = PETSC_FALSE;
682c4762a1bSJed Brown     PetscFunctionReturn(0);
683c4762a1bSJed Brown   }
684c4762a1bSJed Brown 
6859566063dSJacob Faibussowitsch   PetscCall(VecDot(Y,Ylast,&dot));
6869566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y,NORM_2,&ynorm));
6879566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast,NORM_2,&ylastnorm));
688c4762a1bSJed Brown   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
689c4762a1bSJed Brown   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
690c4762a1bSJed Brown   angle_radians = precheck->angle * PETSC_PI / 180.;
691c4762a1bSJed Brown   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
692c4762a1bSJed Brown     /* Modify the step Y */
693c4762a1bSJed Brown     PetscReal alpha,ydiffnorm;
6949566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast,-1.0,Y));
6959566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast,NORM_2,&ydiffnorm));
696c4762a1bSJed Brown     alpha = ylastnorm / ydiffnorm;
6979566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
6989566063dSJacob Faibussowitsch     PetscCall(VecScale(Y,alpha));
699c4762a1bSJed Brown     if (precheck->monitor) {
7009566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(precheck->monitor,"Angle %E degrees less than threshold %g, corrected step by alpha=%g\n",(double)(theta*180./PETSC_PI),(double)precheck->angle,(double)alpha));
701c4762a1bSJed Brown     }
702c4762a1bSJed Brown   } else {
7039566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
704c4762a1bSJed Brown     *changed = PETSC_FALSE;
705c4762a1bSJed Brown     if (precheck->monitor) {
7069566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(precheck->monitor,"Angle %E degrees exceeds threshold %g, no correction applied\n",(double)(theta*180./PETSC_PI),(double)precheck->angle));
707c4762a1bSJed Brown     }
708c4762a1bSJed Brown   }
709c4762a1bSJed Brown   PetscFunctionReturn(0);
710c4762a1bSJed Brown }
711c4762a1bSJed Brown 
712c4762a1bSJed Brown PetscErrorCode PreCheckDestroy(PreCheck *precheck)
713c4762a1bSJed Brown {
714c4762a1bSJed Brown   PetscFunctionBeginUser;
715c4762a1bSJed Brown   if (!*precheck) PetscFunctionReturn(0);
7169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*precheck)->Ylast));
7179566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*precheck)->monitor));
7189566063dSJacob Faibussowitsch   PetscCall(PetscFree(*precheck));
719c4762a1bSJed Brown   PetscFunctionReturn(0);
720c4762a1bSJed Brown }
721c4762a1bSJed Brown 
722c4762a1bSJed Brown PetscErrorCode PreCheckCreate(MPI_Comm comm,PreCheck *precheck)
723c4762a1bSJed Brown {
724c4762a1bSJed Brown   PetscFunctionBeginUser;
7259566063dSJacob Faibussowitsch   PetscCall(PetscNew(precheck));
726c4762a1bSJed Brown 
727c4762a1bSJed Brown   (*precheck)->comm  = comm;
728c4762a1bSJed Brown   (*precheck)->angle = 10.;     /* only active if angle is less than 10 degrees */
729c4762a1bSJed Brown   PetscFunctionReturn(0);
730c4762a1bSJed Brown }
731c4762a1bSJed Brown 
732c4762a1bSJed Brown /*
733c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
734c4762a1bSJed Brown 
735c4762a1bSJed Brown  */
736c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx)
737c4762a1bSJed Brown {
738c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,xm,ym,its,tot_its,sweeps,l,m;
739c4762a1bSJed Brown   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;
740c4762a1bSJed Brown   PetscScalar    **x,**b,bij,F,F0=0,J,y,u,eu;
741c4762a1bSJed Brown   PetscReal      atol,rtol,stol;
742c4762a1bSJed Brown   DM             da;
743c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
744c4762a1bSJed Brown   Vec            localX,localB;
745c4762a1bSJed Brown   DMDALocalInfo  info;
746c4762a1bSJed Brown 
747c4762a1bSJed Brown   PetscFunctionBeginUser;
7489566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
7499566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
750c4762a1bSJed Brown 
751c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info.mx-1);
752c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info.my-1);
753c4762a1bSJed Brown   sc     = hx*hy*user->lambda;
754c4762a1bSJed Brown   dhx    = 1/hx;
755c4762a1bSJed Brown   dhy    = 1/hy;
756c4762a1bSJed Brown   hxdhy  = hx/hy;
757c4762a1bSJed Brown   hydhx  = hy/hx;
758c4762a1bSJed Brown 
759c4762a1bSJed Brown   tot_its = 0;
7609566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetSweeps(snes,&sweeps));
7619566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its));
7629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
763c4762a1bSJed Brown   if (B) {
7649566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(da,&localB));
765c4762a1bSJed Brown   }
766c4762a1bSJed Brown   if (B) {
7679566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB));
7689566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB));
769c4762a1bSJed Brown   }
7709566063dSJacob Faibussowitsch   if (B) PetscCall(DMDAVecGetArrayRead(da,localB,&b));
7719566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
7729566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
7739566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localX,&x));
774c4762a1bSJed Brown   for (l=0; l<sweeps; l++) {
775c4762a1bSJed Brown     /*
776c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
777c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
778c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
779c4762a1bSJed Brown      */
7809566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
781c4762a1bSJed Brown     for (m=0; m<2; m++) {
782c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
783c4762a1bSJed Brown         for (i=xs+(m+j)%2; i<xs+xm; i+=2) {
784c4762a1bSJed Brown           PetscReal xx = i*hx,yy = j*hy;
785c4762a1bSJed Brown           if (B) bij = b[j][i];
786c4762a1bSJed Brown           else   bij = 0.;
787c4762a1bSJed Brown 
788c4762a1bSJed Brown           if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
789c4762a1bSJed Brown             /* boundary conditions are all zero Dirichlet */
790c4762a1bSJed Brown             x[j][i] = 0.0 + bij;
791c4762a1bSJed Brown           } else {
792c4762a1bSJed Brown             const PetscScalar
793c4762a1bSJed Brown               u_E = x[j][i+1],
794c4762a1bSJed Brown               u_W = x[j][i-1],
795c4762a1bSJed Brown               u_N = x[j+1][i],
796c4762a1bSJed Brown               u_S = x[j-1][i];
797c4762a1bSJed Brown             const PetscScalar
798c4762a1bSJed Brown               uy_E   = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
799c4762a1bSJed Brown               uy_W   = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
800c4762a1bSJed Brown               ux_N   = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
801c4762a1bSJed Brown               ux_S   = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]);
802c4762a1bSJed Brown             u = x[j][i];
803c4762a1bSJed Brown             for (k=0; k<its; k++) {
804c4762a1bSJed Brown               const PetscScalar
805c4762a1bSJed Brown                 ux_E   = dhx*(u_E-u),
806c4762a1bSJed Brown                 ux_W   = dhx*(u-u_W),
807c4762a1bSJed Brown                 uy_N   = dhy*(u_N-u),
808c4762a1bSJed Brown                 uy_S   = dhy*(u-u_S),
809c4762a1bSJed Brown                 e_E    = eta(user,xx,yy,ux_E,uy_E),
810c4762a1bSJed Brown                 e_W    = eta(user,xx,yy,ux_W,uy_W),
811c4762a1bSJed Brown                 e_N    = eta(user,xx,yy,ux_N,uy_N),
812c4762a1bSJed Brown                 e_S    = eta(user,xx,yy,ux_S,uy_S),
813c4762a1bSJed Brown                 de_E   = deta(user,xx,yy,ux_E,uy_E),
814c4762a1bSJed Brown                 de_W   = deta(user,xx,yy,ux_W,uy_W),
815c4762a1bSJed Brown                 de_N   = deta(user,xx,yy,ux_N,uy_N),
816c4762a1bSJed Brown                 de_S   = deta(user,xx,yy,ux_S,uy_S),
817c4762a1bSJed Brown                 newt_E = e_E+de_E*PetscSqr(ux_E),
818c4762a1bSJed Brown                 newt_W = e_W+de_W*PetscSqr(ux_W),
819c4762a1bSJed Brown                 newt_N = e_N+de_N*PetscSqr(uy_N),
820c4762a1bSJed Brown                 newt_S = e_S+de_S*PetscSqr(uy_S),
821c4762a1bSJed Brown                 uxx    = -hy * (e_E*ux_E - e_W*ux_W),
822c4762a1bSJed Brown                 uyy    = -hx * (e_N*uy_N - e_S*uy_S);
823c4762a1bSJed Brown 
824c4762a1bSJed Brown               if (sc) eu = PetscExpScalar(u);
825c4762a1bSJed Brown               else    eu = 0;
826c4762a1bSJed Brown 
827c4762a1bSJed Brown               F = uxx + uyy - sc*eu - bij;
828c4762a1bSJed Brown               if (k == 0) F0 = F;
829c4762a1bSJed Brown               J  = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*eu;
830c4762a1bSJed Brown               y  = F/J;
831c4762a1bSJed Brown               u -= y;
832c4762a1bSJed Brown               tot_its++;
833c4762a1bSJed Brown               if (atol > PetscAbsReal(PetscRealPart(F)) ||
834c4762a1bSJed Brown                   rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) ||
835c4762a1bSJed Brown                   stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) {
836c4762a1bSJed Brown                 break;
837c4762a1bSJed Brown               }
838c4762a1bSJed Brown             }
839c4762a1bSJed Brown             x[j][i] = u;
840c4762a1bSJed Brown           }
841c4762a1bSJed Brown         }
842c4762a1bSJed Brown       }
843c4762a1bSJed Brown     }
844c4762a1bSJed Brown     /*
845c4762a1bSJed Brown x     Restore vector
846c4762a1bSJed Brown      */
847c4762a1bSJed Brown   }
8489566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localX,&x));
8499566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X));
8509566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X));
8519566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(tot_its*(118.0)));
8529566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
853c4762a1bSJed Brown   if (B) {
8549566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da,localB,&b));
8559566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da,&localB));
856c4762a1bSJed Brown   }
857c4762a1bSJed Brown   PetscFunctionReturn(0);
858c4762a1bSJed Brown }
859c4762a1bSJed Brown 
860c4762a1bSJed Brown /*TEST
861c4762a1bSJed Brown 
862c4762a1bSJed Brown    test:
863c4762a1bSJed Brown       nsize: 2
864c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype NEWTON
865c4762a1bSJed Brown       requires: !single
866c4762a1bSJed Brown 
867c4762a1bSJed Brown    test:
868c4762a1bSJed Brown       suffix: 2
869c4762a1bSJed Brown       nsize: 2
870c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -precheck 1
871c4762a1bSJed Brown       requires: !single
872c4762a1bSJed Brown 
873c4762a1bSJed Brown    test:
874c4762a1bSJed Brown       suffix: 3
875c4762a1bSJed Brown       nsize: 2
87685b95db3SBarry 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
87785b95db3SBarry Smith       requires: !single
87885b95db3SBarry Smith 
87985b95db3SBarry Smith    test:
88085b95db3SBarry Smith       suffix: 3_star
88185b95db3SBarry Smith       nsize: 2
88285b95db3SBarry 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
88385b95db3SBarry Smith       output_file: output/ex15_3.out
884c4762a1bSJed Brown       requires: !single
885c4762a1bSJed Brown 
886c4762a1bSJed Brown    test:
887c4762a1bSJed Brown       suffix: 4
888c4762a1bSJed 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
889c4762a1bSJed Brown       requires: !single
890c4762a1bSJed Brown 
891c4762a1bSJed Brown    test:
892c4762a1bSJed Brown       suffix: lag_jac
893c4762a1bSJed Brown       nsize: 4
894c4762a1bSJed 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
895c4762a1bSJed Brown       requires: !single
896c4762a1bSJed Brown 
897c4762a1bSJed Brown    test:
898c4762a1bSJed Brown       suffix: lag_pc
899c4762a1bSJed Brown       nsize: 4
900c4762a1bSJed 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
901c4762a1bSJed Brown       requires: !single
902c4762a1bSJed Brown 
903c4762a1bSJed Brown    test:
904c4762a1bSJed Brown       suffix: nleqerr
905c4762a1bSJed 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
906c4762a1bSJed Brown       requires: !single
907c4762a1bSJed Brown 
908bbc1464cSBarry Smith    test:
909bbc1464cSBarry Smith       suffix: mf
910bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator
911bbc1464cSBarry Smith       requires: !single
912bbc1464cSBarry Smith 
913bbc1464cSBarry Smith    test:
914bbc1464cSBarry Smith       suffix: mf_picard
915bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator -picard
916bbc1464cSBarry Smith       requires: !single
917bbc1464cSBarry Smith       output_file: output/ex15_mf.out
918bbc1464cSBarry Smith 
91985b95db3SBarry Smith    test:
92085b95db3SBarry Smith       suffix: fd_picard
92185b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 2  -p 3 -ksp_rtol 1.e-12  -snes_fd -picard
92285b95db3SBarry Smith       requires: !single
92385b95db3SBarry Smith 
92485b95db3SBarry Smith    test:
92585b95db3SBarry Smith       suffix: fd_color_picard
92685b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_fd_color -picard
92785b95db3SBarry Smith       requires: !single
92885b95db3SBarry Smith       output_file: output/ex15_mf.out
92985b95db3SBarry Smith 
930c4762a1bSJed Brown TEST*/
931