xref: /petsc/src/snes/tutorials/ex15.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
116*327415f7SBarry Smith   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120c4762a1bSJed Brown      Initialize problem parameters
121c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122c4762a1bSJed Brown   user.lambda    = 0.0; user.p = 2.0; user.epsilon = 1e-5; user.source = 0.1; user.jtype = JAC_NEWTON;user.initial=-1;
123c4762a1bSJed Brown   user.blocks[0] = 1; user.blocks[1] = 1; user.kappa = 1e-3;
124c4762a1bSJed Brown   alloc_star     = PETSC_FALSE;
125c4762a1bSJed Brown   use_precheck   = 0; precheck_angle = 10.;
126c4762a1bSJed Brown   user.picard    = PETSC_FALSE;
127d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"p-Bratu options",__FILE__);
128c4762a1bSJed Brown   {
129c4762a1bSJed Brown     PetscInt two=2;
1309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda","Bratu parameter","",user.lambda,&user.lambda,NULL));
1319566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-p","Exponent `p' in p-Laplacian","",user.p,&user.p,NULL));
1329566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-epsilon","Strain-regularization in p-Laplacian","",user.epsilon,&user.epsilon,NULL));
1339566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-source","Constant source term","",user.source,&user.source,NULL));
1349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-jtype","Jacobian approximation to assemble","",JacTypes,(PetscEnum)user.jtype,(PetscEnum*)&user.jtype,NULL));
1359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-picard","Solve with defect-correction Picard iteration","",&user.picard));
13685b95db3SBarry Smith     if (user.picard) {
13785b95db3SBarry Smith       user.jtype = JAC_PICARD;
138e00437b9SBarry Smith       PetscCheck(user.p == 3,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Picard iteration is only supported for p == 3");
13985b95db3SBarry Smith       /* the Picard linearization only requires a 5 point stencil, while the Newton linearization requires a 9 point stencil */
14085b95db3SBarry Smith       /* hence allocating the 5 point stencil gives the same convergence as the 9 point stencil since the extra stencil points are not used */
1419566063dSJacob Faibussowitsch       PetscCall(PetscOptionsBool("-alloc_star","Allocate for STAR stencil (5-point)","",alloc_star,&alloc_star,NULL));
14285b95db3SBarry Smith     }
1439566063dSJacob 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));
1449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-initial","Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)","",user.initial,&user.initial,NULL));
145c4762a1bSJed Brown     if (use_precheck == 2) {    /* Using library version, get the angle */
1469566063dSJacob Faibussowitsch       PetscCall(PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck_angle,&precheck_angle,NULL));
147c4762a1bSJed Brown     }
1489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-blocks","number of coefficient interfaces in x and y direction","",user.blocks,&two,NULL));
1499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-kappa","diffusivity in odd regions","",user.kappa,&user.kappa,NULL));
1509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-o","Output solution in vts format","",filename,filename,sizeof(filename),&write_output));
151c4762a1bSJed Brown   }
152d0609cedSBarry Smith   PetscOptionsEnd();
153c4762a1bSJed Brown   if (user.lambda > bratu_lambda_max || user.lambda < bratu_lambda_min) {
1549566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING: lambda %g out of range for p=2\n",(double)user.lambda));
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158c4762a1bSJed Brown      Create nonlinear solver context
159c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1609566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
164c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1659566063dSJacob 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));
1669566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1679566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170c4762a1bSJed Brown      Extract global vectors from DM; then duplicate for remaining
171c4762a1bSJed Brown      vectors that are the same types
172c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1739566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&x));
1749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
1759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&b));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
17885b95db3SBarry Smith      User can override with:
179c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
180c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
181c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
182c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
183c4762a1bSJed Brown                          products within Newton-Krylov method
184c4762a1bSJed Brown 
185c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown      Set local function evaluation routine
189c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1909566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
1919566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes,da));
192c4762a1bSJed Brown   if (user.picard) {
193c4762a1bSJed Brown     /*
194c4762a1bSJed Brown         This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access
195c4762a1bSJed Brown         the SNES to set it
196c4762a1bSJed Brown     */
197d0609cedSBarry Smith     PetscCall(DMDASNESSetPicardLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionPicardLocal,
198d0609cedSBarry Smith                                      (PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user));
1999566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user));
2009566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes,NULL,NULL,SNESPicardComputeJacobian,&user));
201c4762a1bSJed Brown   } else {
2029566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user));
2039566063dSJacob Faibussowitsch     PetscCall(DMDASNESSetJacobianLocal(da,(PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user));
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
208c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2099566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
2109566063dSJacob Faibussowitsch   PetscCall(SNESSetNGS(snes,NonlinearGS,&user));
2119566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
212c4762a1bSJed Brown   /* Set up the precheck context if requested */
213c4762a1bSJed Brown   if (use_precheck == 1) {      /* Use the precheck routines in this file */
2149566063dSJacob Faibussowitsch     PetscCall(PreCheckCreate(PETSC_COMM_WORLD,&precheck));
2159566063dSJacob Faibussowitsch     PetscCall(PreCheckSetFromOptions(precheck));
2169566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch,PreCheckFunction,precheck));
217c4762a1bSJed Brown   } else if (use_precheck == 2) { /* Use the version provided by the library */
2189566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&precheck_angle));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222c4762a1bSJed Brown      Evaluate initial guess
223c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
224c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
225c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
226c4762a1bSJed Brown      this vector to zero by calling VecSet().
227c4762a1bSJed Brown   */
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user,da,x));
2309566063dSJacob Faibussowitsch   PetscCall(FormRHS(&user,da,b));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233c4762a1bSJed Brown      Solve nonlinear system
234c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2359566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,b,x));
2369566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
2379566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes,&reason));
238c4762a1bSJed Brown 
23963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s Number of nonlinear iterations = %" PetscInt_FMT "\n",SNESConvergedReasons[reason],its));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   if (write_output) {
242c4762a1bSJed Brown     PetscViewer viewer;
2439566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer));
2449566063dSJacob Faibussowitsch     PetscCall(VecView(x,viewer));
2459566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
250c4762a1bSJed Brown      are no longer needed.
251c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
252c4762a1bSJed Brown 
2539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
2569566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2589566063dSJacob Faibussowitsch   PetscCall(PreCheckDestroy(&precheck));
2599566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
260b122ec5aSJacob Faibussowitsch   return 0;
261c4762a1bSJed Brown }
262c4762a1bSJed Brown 
263c4762a1bSJed Brown /* ------------------------------------------------------------------- */
264c4762a1bSJed Brown /*
265c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
266c4762a1bSJed Brown 
267c4762a1bSJed Brown    Input Parameters:
268c4762a1bSJed Brown    user - user-defined application context
269c4762a1bSJed Brown    X - vector
270c4762a1bSJed Brown 
271c4762a1bSJed Brown    Output Parameter:
272c4762a1bSJed Brown    X - vector
273c4762a1bSJed Brown  */
274c4762a1bSJed Brown static PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
275c4762a1bSJed Brown {
276c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
277c4762a1bSJed Brown   PetscReal      temp1,temp,hx,hy;
278c4762a1bSJed Brown   PetscScalar    **x;
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   PetscFunctionBeginUser;
2819566063dSJacob 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));
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(Mx-1);
284c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(My-1);
285c4762a1bSJed Brown   temp1 = user->lambda / (user->lambda + 1.);
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /*
288c4762a1bSJed Brown      Get a pointer to vector data.
289c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
290c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
291c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
292c4762a1bSJed Brown          the array.
293c4762a1bSJed Brown   */
2949566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   /*
297c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DA):
298c4762a1bSJed Brown        xs, ys   - starting grid indices (no ghost points)
299c4762a1bSJed Brown        xm, ym   - widths of local grid (no ghost points)
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   */
3029566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /*
305c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
306c4762a1bSJed Brown   */
307c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
308c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
309c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
310c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
311c4762a1bSJed Brown         /* boundary conditions are all zero Dirichlet */
312c4762a1bSJed Brown         x[j][i] = 0.0;
313c4762a1bSJed Brown       } else {
314c4762a1bSJed Brown         if (user->initial == -1) {
315c4762a1bSJed Brown           if (user->lambda != 0) {
316c4762a1bSJed Brown             x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
317c4762a1bSJed Brown           } else {
318c4762a1bSJed Brown             /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting
319c4762a1bSJed Brown              * with an exact solution. */
320c4762a1bSJed Brown             const PetscReal
321c4762a1bSJed Brown               xx = 2*(PetscReal)i/(Mx-1) - 1,
322c4762a1bSJed Brown               yy = 2*(PetscReal)j/(My-1) - 1;
323c4762a1bSJed Brown             x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
324c4762a1bSJed Brown           }
325c4762a1bSJed Brown         } else if (user->initial == 0) {
326c4762a1bSJed Brown           x[j][i] = 0.;
327c4762a1bSJed Brown         } else if (user->initial == 1) {
328c4762a1bSJed Brown           const PetscReal
329c4762a1bSJed Brown             xx = 2*(PetscReal)i/(Mx-1) - 1,
330c4762a1bSJed Brown             yy = 2*(PetscReal)j/(My-1) - 1;
331c4762a1bSJed Brown           x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
332c4762a1bSJed Brown         } else {
333c4762a1bSJed Brown           if (user->lambda != 0) {
334c4762a1bSJed Brown             x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
335c4762a1bSJed Brown           } else {
336c4762a1bSJed Brown             x[j][i] = 0.5*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
337c4762a1bSJed Brown           }
338c4762a1bSJed Brown         }
339c4762a1bSJed Brown       }
340c4762a1bSJed Brown     }
341c4762a1bSJed Brown   }
342c4762a1bSJed Brown   /*
343c4762a1bSJed Brown      Restore vector
344c4762a1bSJed Brown   */
3459566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,X,&x));
346c4762a1bSJed Brown   PetscFunctionReturn(0);
347c4762a1bSJed Brown }
348c4762a1bSJed Brown 
349c4762a1bSJed Brown /*
350c4762a1bSJed Brown    FormRHS - Forms constant RHS for the problem.
351c4762a1bSJed Brown 
352c4762a1bSJed Brown    Input Parameters:
353c4762a1bSJed Brown    user - user-defined application context
354c4762a1bSJed Brown    B - RHS vector
355c4762a1bSJed Brown 
356c4762a1bSJed Brown    Output Parameter:
357c4762a1bSJed Brown    B - vector
358c4762a1bSJed Brown  */
359c4762a1bSJed Brown static PetscErrorCode FormRHS(AppCtx *user,DM da,Vec B)
360c4762a1bSJed Brown {
361c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
362c4762a1bSJed Brown   PetscReal      hx,hy;
363c4762a1bSJed Brown   PetscScalar    **b;
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   PetscFunctionBeginUser;
3669566063dSJacob 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));
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(Mx-1);
369c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(My-1);
3709566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,B,&b));
3719566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
372c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
373c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
374c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
375c4762a1bSJed Brown         b[j][i] = 0.0;
376c4762a1bSJed Brown       } else {
377c4762a1bSJed Brown         b[j][i] = hx*hy*user->source;
378c4762a1bSJed Brown       }
379c4762a1bSJed Brown     }
380c4762a1bSJed Brown   }
3819566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,B,&b));
382c4762a1bSJed Brown   PetscFunctionReturn(0);
383c4762a1bSJed Brown }
384c4762a1bSJed Brown 
3859fbee547SJacob Faibussowitsch static inline PetscReal kappa(const AppCtx *ctx,PetscReal x,PetscReal y)
386c4762a1bSJed Brown {
387c4762a1bSJed Brown   return (((PetscInt)(x*ctx->blocks[0])) + ((PetscInt)(y*ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0;
388c4762a1bSJed Brown }
389c4762a1bSJed Brown /* p-Laplacian diffusivity */
3909fbee547SJacob Faibussowitsch static inline PetscScalar eta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
391c4762a1bSJed Brown {
392c4762a1bSJed Brown   return kappa(ctx,x,y) * PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-2.));
393c4762a1bSJed Brown }
3949fbee547SJacob Faibussowitsch static inline PetscScalar deta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
395c4762a1bSJed Brown {
396c4762a1bSJed Brown   return (ctx->p == 2)
397c4762a1bSJed Brown          ? 0
398c4762a1bSJed 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.);
399c4762a1bSJed Brown }
400c4762a1bSJed Brown 
401c4762a1bSJed Brown /* ------------------------------------------------------------------- */
402c4762a1bSJed Brown /*
403c4762a1bSJed Brown    FormFunctionLocal - Evaluates nonlinear function, F(x).
404c4762a1bSJed Brown  */
405c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
406c4762a1bSJed Brown {
407c4762a1bSJed Brown   PetscReal      hx,hy,dhx,dhy,sc;
408c4762a1bSJed Brown   PetscInt       i,j;
409c4762a1bSJed Brown   PetscScalar    eu;
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   PetscFunctionBeginUser;
412c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info->mx-1);
413c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info->my-1);
414c4762a1bSJed Brown   sc     = hx*hy*user->lambda;
415c4762a1bSJed Brown   dhx    = 1/hx;
416c4762a1bSJed Brown   dhy    = 1/hy;
417c4762a1bSJed Brown   /*
418c4762a1bSJed Brown      Compute function over the locally owned part of the grid
419c4762a1bSJed Brown   */
420c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
421c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
422c4762a1bSJed Brown       PetscReal xx = i*hx,yy = j*hy;
423c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
424c4762a1bSJed Brown         f[j][i] = x[j][i];
425c4762a1bSJed Brown       } else {
426c4762a1bSJed Brown         const PetscScalar
427c4762a1bSJed Brown           u    = x[j][i],
428c4762a1bSJed Brown           ux_E = dhx*(x[j][i+1]-x[j][i]),
429c4762a1bSJed Brown           uy_E = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
430c4762a1bSJed Brown           ux_W = dhx*(x[j][i]-x[j][i-1]),
431c4762a1bSJed Brown           uy_W = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
432c4762a1bSJed Brown           ux_N = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
433c4762a1bSJed Brown           uy_N = dhy*(x[j+1][i]-x[j][i]),
434c4762a1bSJed Brown           ux_S = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
435c4762a1bSJed Brown           uy_S = dhy*(x[j][i]-x[j-1][i]),
436c4762a1bSJed Brown           e_E  = eta(user,xx,yy,ux_E,uy_E),
437c4762a1bSJed Brown           e_W  = eta(user,xx,yy,ux_W,uy_W),
438c4762a1bSJed Brown           e_N  = eta(user,xx,yy,ux_N,uy_N),
439c4762a1bSJed Brown           e_S  = eta(user,xx,yy,ux_S,uy_S),
440c4762a1bSJed Brown           uxx  = -hy * (e_E*ux_E - e_W*ux_W),
441c4762a1bSJed Brown           uyy  = -hx * (e_N*uy_N - e_S*uy_S);
442c4762a1bSJed Brown         if (sc) eu = PetscExpScalar(u);
443c4762a1bSJed Brown         else    eu = 0.;
4440e3d61c9SBarry Smith         /* For p=2, these terms decay to:
4450e3d61c9SBarry Smith          uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx
4460e3d61c9SBarry Smith          uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy
4470e3d61c9SBarry Smith         */
448c4762a1bSJed Brown         f[j][i] = uxx + uyy - sc*eu;
449c4762a1bSJed Brown       }
450c4762a1bSJed Brown     }
451c4762a1bSJed Brown   }
4529566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm*info->ym*(72.0)));
453c4762a1bSJed Brown   PetscFunctionReturn(0);
454c4762a1bSJed Brown }
455c4762a1bSJed Brown 
456c4762a1bSJed Brown /*
457c4762a1bSJed Brown     This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation,
458c4762a1bSJed 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)
459c4762a1bSJed Brown 
460c4762a1bSJed Brown */
461c4762a1bSJed Brown static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
462c4762a1bSJed Brown {
463c4762a1bSJed Brown   PetscReal hx,hy,sc;
464c4762a1bSJed Brown   PetscInt  i,j;
465c4762a1bSJed Brown 
466c4762a1bSJed Brown   PetscFunctionBeginUser;
467c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info->mx-1);
468c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info->my-1);
469c4762a1bSJed Brown   sc     = hx*hy*user->lambda;
470c4762a1bSJed Brown   /*
471c4762a1bSJed Brown      Compute function over the locally owned part of the grid
472c4762a1bSJed Brown   */
473c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
474c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
475c4762a1bSJed Brown       if (!(i == 0 || j == 0 || i == info->mx-1 || j == info->my-1)) {
476c4762a1bSJed Brown         const PetscScalar u = x[j][i];
477c4762a1bSJed Brown         f[j][i] = sc*PetscExpScalar(u);
4780df40c35SBarry Smith       } else {
4790df40c35SBarry Smith         f[j][i] = 0.0; /* this is zero because the A(x) x term forces the x to be zero on the boundary */
480c4762a1bSJed Brown       }
481c4762a1bSJed Brown     }
482c4762a1bSJed Brown   }
4839566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(info->xm*info->ym*2.0));
484c4762a1bSJed Brown   PetscFunctionReturn(0);
485c4762a1bSJed Brown }
486c4762a1bSJed Brown 
487c4762a1bSJed Brown /*
488c4762a1bSJed Brown    FormJacobianLocal - Evaluates Jacobian matrix.
489c4762a1bSJed Brown */
490c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat J,Mat B,AppCtx *user)
491c4762a1bSJed Brown {
492c4762a1bSJed Brown   PetscInt       i,j;
493c4762a1bSJed Brown   MatStencil     col[9],row;
494c4762a1bSJed Brown   PetscScalar    v[9];
495c4762a1bSJed Brown   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   PetscFunctionBeginUser;
4989566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
499c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(info->mx-1);
500c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(info->my-1);
501c4762a1bSJed Brown   sc    = hx*hy*user->lambda;
502c4762a1bSJed Brown   dhx   = 1/hx;
503c4762a1bSJed Brown   dhy   = 1/hy;
504c4762a1bSJed Brown   hxdhy = hx/hy;
505c4762a1bSJed Brown   hydhx = hy/hx;
506c4762a1bSJed Brown 
507c4762a1bSJed Brown   /*
508c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
509c4762a1bSJed Brown       - PETSc parallel matrix formats are partitioned by
510c4762a1bSJed Brown         contiguous chunks of rows across the processors.
511c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
512c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
513c4762a1bSJed Brown         appropriate processor during matrix assembly).
514c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
515c4762a1bSJed Brown   */
516c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
517c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
518c4762a1bSJed Brown       PetscReal xx = i*hx,yy = j*hy;
519c4762a1bSJed Brown       row.j = j; row.i = i;
520c4762a1bSJed Brown       /* boundary points */
521c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
522c4762a1bSJed Brown         v[0] = 1.0;
5239566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES));
524c4762a1bSJed Brown       } else {
525c4762a1bSJed Brown         /* interior grid points */
526c4762a1bSJed Brown         const PetscScalar
527c4762a1bSJed Brown           ux_E     = dhx*(x[j][i+1]-x[j][i]),
528c4762a1bSJed Brown           uy_E     = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
529c4762a1bSJed Brown           ux_W     = dhx*(x[j][i]-x[j][i-1]),
530c4762a1bSJed Brown           uy_W     = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
531c4762a1bSJed Brown           ux_N     = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
532c4762a1bSJed Brown           uy_N     = dhy*(x[j+1][i]-x[j][i]),
533c4762a1bSJed Brown           ux_S     = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
534c4762a1bSJed Brown           uy_S     = dhy*(x[j][i]-x[j-1][i]),
535c4762a1bSJed Brown           u        = x[j][i],
536c4762a1bSJed Brown           e_E      = eta(user,xx,yy,ux_E,uy_E),
537c4762a1bSJed Brown           e_W      = eta(user,xx,yy,ux_W,uy_W),
538c4762a1bSJed Brown           e_N      = eta(user,xx,yy,ux_N,uy_N),
539c4762a1bSJed Brown           e_S      = eta(user,xx,yy,ux_S,uy_S),
540c4762a1bSJed Brown           de_E     = deta(user,xx,yy,ux_E,uy_E),
541c4762a1bSJed Brown           de_W     = deta(user,xx,yy,ux_W,uy_W),
542c4762a1bSJed Brown           de_N     = deta(user,xx,yy,ux_N,uy_N),
543c4762a1bSJed Brown           de_S     = deta(user,xx,yy,ux_S,uy_S),
544c4762a1bSJed Brown           skew_E   = de_E*ux_E*uy_E,
545c4762a1bSJed Brown           skew_W   = de_W*ux_W*uy_W,
546c4762a1bSJed Brown           skew_N   = de_N*ux_N*uy_N,
547c4762a1bSJed Brown           skew_S   = de_S*ux_S*uy_S,
548c4762a1bSJed Brown           cross_EW = 0.25*(skew_E - skew_W),
549c4762a1bSJed Brown           cross_NS = 0.25*(skew_N - skew_S),
550c4762a1bSJed Brown           newt_E   = e_E+de_E*PetscSqr(ux_E),
551c4762a1bSJed Brown           newt_W   = e_W+de_W*PetscSqr(ux_W),
552c4762a1bSJed Brown           newt_N   = e_N+de_N*PetscSqr(uy_N),
553c4762a1bSJed Brown           newt_S   = e_S+de_S*PetscSqr(uy_S);
554c4762a1bSJed Brown         /* interior grid points */
555c4762a1bSJed Brown         switch (user->jtype) {
556c4762a1bSJed Brown         case JAC_BRATU:
557c4762a1bSJed Brown           /* Jacobian from p=2 */
558c4762a1bSJed Brown           v[0] = -hxdhy;                                           col[0].j = j-1;   col[0].i = i;
559c4762a1bSJed Brown           v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
560c4762a1bSJed Brown           v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(u);       col[2].j = row.j; col[2].i = row.i;
561c4762a1bSJed Brown           v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
562c4762a1bSJed Brown           v[4] = -hxdhy;                                           col[4].j = j+1;   col[4].i = i;
5639566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
564c4762a1bSJed Brown           break;
565c4762a1bSJed Brown         case JAC_PICARD:
566c4762a1bSJed Brown           /* Jacobian arising from Picard linearization */
567c4762a1bSJed Brown           v[0] = -hxdhy*e_S;                                           col[0].j = j-1;   col[0].i = i;
568c4762a1bSJed Brown           v[1] = -hydhx*e_W;                                           col[1].j = j;     col[1].i = i-1;
569c4762a1bSJed Brown           v[2] = (e_W+e_E)*hydhx + (e_S+e_N)*hxdhy;                    col[2].j = row.j; col[2].i = row.i;
570c4762a1bSJed Brown           v[3] = -hydhx*e_E;                                           col[3].j = j;     col[3].i = i+1;
571c4762a1bSJed Brown           v[4] = -hxdhy*e_N;                                           col[4].j = j+1;   col[4].i = i;
5729566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
573c4762a1bSJed Brown           break;
574c4762a1bSJed Brown         case JAC_STAR:
575c4762a1bSJed Brown           /* Full Jacobian, but only a star stencil */
576c4762a1bSJed Brown           col[0].j = j-1; col[0].i = i;
577c4762a1bSJed Brown           col[1].j = j;   col[1].i = i-1;
578c4762a1bSJed Brown           col[2].j = j;   col[2].i = i;
579c4762a1bSJed Brown           col[3].j = j;   col[3].i = i+1;
580c4762a1bSJed Brown           col[4].j = j+1; col[4].i = i;
581c4762a1bSJed Brown           v[0]     = -hxdhy*newt_S + cross_EW;
582c4762a1bSJed Brown           v[1]     = -hydhx*newt_W + cross_NS;
583c4762a1bSJed Brown           v[2]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
584c4762a1bSJed Brown           v[3]     = -hydhx*newt_E - cross_NS;
585c4762a1bSJed Brown           v[4]     = -hxdhy*newt_N - cross_EW;
5869566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
587c4762a1bSJed Brown           break;
588c4762a1bSJed Brown         case JAC_NEWTON:
5890e3d61c9SBarry Smith           /* The Jacobian is
5900e3d61c9SBarry Smith 
5910e3d61c9SBarry Smith            -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u
5920e3d61c9SBarry Smith 
5930e3d61c9SBarry Smith           */
594c4762a1bSJed Brown           col[0].j = j-1; col[0].i = i-1;
595c4762a1bSJed Brown           col[1].j = j-1; col[1].i = i;
596c4762a1bSJed Brown           col[2].j = j-1; col[2].i = i+1;
597c4762a1bSJed Brown           col[3].j = j;   col[3].i = i-1;
598c4762a1bSJed Brown           col[4].j = j;   col[4].i = i;
599c4762a1bSJed Brown           col[5].j = j;   col[5].i = i+1;
600c4762a1bSJed Brown           col[6].j = j+1; col[6].i = i-1;
601c4762a1bSJed Brown           col[7].j = j+1; col[7].i = i;
602c4762a1bSJed Brown           col[8].j = j+1; col[8].i = i+1;
603c4762a1bSJed Brown           v[0]     = -0.25*(skew_S + skew_W);
604c4762a1bSJed Brown           v[1]     = -hxdhy*newt_S + cross_EW;
605c4762a1bSJed Brown           v[2]     =  0.25*(skew_S + skew_E);
606c4762a1bSJed Brown           v[3]     = -hydhx*newt_W + cross_NS;
607c4762a1bSJed Brown           v[4]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
608c4762a1bSJed Brown           v[5]     = -hydhx*newt_E - cross_NS;
609c4762a1bSJed Brown           v[6]     =  0.25*(skew_N + skew_W);
610c4762a1bSJed Brown           v[7]     = -hxdhy*newt_N - cross_EW;
611c4762a1bSJed Brown           v[8]     = -0.25*(skew_N + skew_E);
6129566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,9,col,v,INSERT_VALUES));
613c4762a1bSJed Brown           break;
614c4762a1bSJed Brown         default:
61598921bdaSJacob Faibussowitsch           SETERRQ(PetscObjectComm((PetscObject)info->da),PETSC_ERR_SUP,"Jacobian type %d not implemented",user->jtype);
616c4762a1bSJed Brown         }
617c4762a1bSJed Brown       }
618c4762a1bSJed Brown     }
619c4762a1bSJed Brown   }
620c4762a1bSJed Brown 
621c4762a1bSJed Brown   /*
622c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
623c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
624c4762a1bSJed Brown   */
6259566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
6269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
627c4762a1bSJed Brown 
628c4762a1bSJed Brown   if (J != B) {
6299566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
6309566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
631c4762a1bSJed Brown   }
632c4762a1bSJed Brown 
633c4762a1bSJed Brown   /*
634c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
635c4762a1bSJed Brown      matrix. If we do, it will generate an error.
636c4762a1bSJed Brown   */
637c4762a1bSJed Brown   if (user->jtype == JAC_NEWTON) {
6389566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(info->xm*info->ym*(131.0)));
639c4762a1bSJed Brown   }
640c4762a1bSJed Brown   PetscFunctionReturn(0);
641c4762a1bSJed Brown }
642c4762a1bSJed Brown 
643c4762a1bSJed Brown /***********************************************************
644c4762a1bSJed Brown  * PreCheck implementation
645c4762a1bSJed Brown  ***********************************************************/
646c4762a1bSJed Brown PetscErrorCode PreCheckSetFromOptions(PreCheck precheck)
647c4762a1bSJed Brown {
648c4762a1bSJed Brown   PetscBool      flg;
649c4762a1bSJed Brown 
650c4762a1bSJed Brown   PetscFunctionBeginUser;
651d0609cedSBarry Smith   PetscOptionsBegin(precheck->comm,NULL,"PreCheck Options","none");
6529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck->angle,&precheck->angle,NULL));
653c4762a1bSJed Brown   flg  = PETSC_FALSE;
6549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-precheck_monitor","Monitor choices made by precheck routine","",flg,&flg,NULL));
655c4762a1bSJed Brown   if (flg) {
6569566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(precheck->comm,"stdout",&precheck->monitor));
657c4762a1bSJed Brown   }
658d0609cedSBarry Smith   PetscOptionsEnd();
659c4762a1bSJed Brown   PetscFunctionReturn(0);
660c4762a1bSJed Brown }
661c4762a1bSJed Brown 
662c4762a1bSJed Brown /*
663c4762a1bSJed Brown   Compare the direction of the current and previous step, modify the current step accordingly
664c4762a1bSJed Brown */
665c4762a1bSJed Brown PetscErrorCode PreCheckFunction(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed, void *ctx)
666c4762a1bSJed Brown {
667c4762a1bSJed Brown   PreCheck       precheck;
668c4762a1bSJed Brown   Vec            Ylast;
669c4762a1bSJed Brown   PetscScalar    dot;
670c4762a1bSJed Brown   PetscInt       iter;
671c4762a1bSJed Brown   PetscReal      ynorm,ylastnorm,theta,angle_radians;
672c4762a1bSJed Brown   SNES           snes;
673c4762a1bSJed Brown 
674c4762a1bSJed Brown   PetscFunctionBeginUser;
6759566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
676c4762a1bSJed Brown   precheck = (PreCheck)ctx;
6779566063dSJacob Faibussowitsch   if (!precheck->Ylast) PetscCall(VecDuplicate(Y,&precheck->Ylast));
678c4762a1bSJed Brown   Ylast = precheck->Ylast;
6799566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&iter));
680c4762a1bSJed Brown   if (iter < 1) {
6819566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
682c4762a1bSJed Brown     *changed = PETSC_FALSE;
683c4762a1bSJed Brown     PetscFunctionReturn(0);
684c4762a1bSJed Brown   }
685c4762a1bSJed Brown 
6869566063dSJacob Faibussowitsch   PetscCall(VecDot(Y,Ylast,&dot));
6879566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y,NORM_2,&ynorm));
6889566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast,NORM_2,&ylastnorm));
689c4762a1bSJed Brown   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
690c4762a1bSJed Brown   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
691c4762a1bSJed Brown   angle_radians = precheck->angle * PETSC_PI / 180.;
692c4762a1bSJed Brown   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
693c4762a1bSJed Brown     /* Modify the step Y */
694c4762a1bSJed Brown     PetscReal alpha,ydiffnorm;
6959566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast,-1.0,Y));
6969566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast,NORM_2,&ydiffnorm));
697c4762a1bSJed Brown     alpha = ylastnorm / ydiffnorm;
6989566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
6999566063dSJacob Faibussowitsch     PetscCall(VecScale(Y,alpha));
700c4762a1bSJed Brown     if (precheck->monitor) {
70163a3b9bcSJacob Faibussowitsch       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));
702c4762a1bSJed Brown     }
703c4762a1bSJed Brown   } else {
7049566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
705c4762a1bSJed Brown     *changed = PETSC_FALSE;
706c4762a1bSJed Brown     if (precheck->monitor) {
70763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(precheck->monitor,"Angle %g degrees exceeds threshold %g, no correction applied\n",(double)(theta*180./PETSC_PI),(double)precheck->angle));
708c4762a1bSJed Brown     }
709c4762a1bSJed Brown   }
710c4762a1bSJed Brown   PetscFunctionReturn(0);
711c4762a1bSJed Brown }
712c4762a1bSJed Brown 
713c4762a1bSJed Brown PetscErrorCode PreCheckDestroy(PreCheck *precheck)
714c4762a1bSJed Brown {
715c4762a1bSJed Brown   PetscFunctionBeginUser;
716c4762a1bSJed Brown   if (!*precheck) PetscFunctionReturn(0);
7179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*precheck)->Ylast));
7189566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*precheck)->monitor));
7199566063dSJacob Faibussowitsch   PetscCall(PetscFree(*precheck));
720c4762a1bSJed Brown   PetscFunctionReturn(0);
721c4762a1bSJed Brown }
722c4762a1bSJed Brown 
723c4762a1bSJed Brown PetscErrorCode PreCheckCreate(MPI_Comm comm,PreCheck *precheck)
724c4762a1bSJed Brown {
725c4762a1bSJed Brown   PetscFunctionBeginUser;
7269566063dSJacob Faibussowitsch   PetscCall(PetscNew(precheck));
727c4762a1bSJed Brown 
728c4762a1bSJed Brown   (*precheck)->comm  = comm;
729c4762a1bSJed Brown   (*precheck)->angle = 10.;     /* only active if angle is less than 10 degrees */
730c4762a1bSJed Brown   PetscFunctionReturn(0);
731c4762a1bSJed Brown }
732c4762a1bSJed Brown 
733c4762a1bSJed Brown /*
734c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
735c4762a1bSJed Brown 
736c4762a1bSJed Brown  */
737c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx)
738c4762a1bSJed Brown {
739c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,xm,ym,its,tot_its,sweeps,l,m;
740c4762a1bSJed Brown   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;
741c4762a1bSJed Brown   PetscScalar    **x,**b,bij,F,F0=0,J,y,u,eu;
742c4762a1bSJed Brown   PetscReal      atol,rtol,stol;
743c4762a1bSJed Brown   DM             da;
744c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
745c4762a1bSJed Brown   Vec            localX,localB;
746c4762a1bSJed Brown   DMDALocalInfo  info;
747c4762a1bSJed Brown 
748c4762a1bSJed Brown   PetscFunctionBeginUser;
7499566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
7509566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
751c4762a1bSJed Brown 
752c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info.mx-1);
753c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info.my-1);
754c4762a1bSJed Brown   sc     = hx*hy*user->lambda;
755c4762a1bSJed Brown   dhx    = 1/hx;
756c4762a1bSJed Brown   dhy    = 1/hy;
757c4762a1bSJed Brown   hxdhy  = hx/hy;
758c4762a1bSJed Brown   hydhx  = hy/hx;
759c4762a1bSJed Brown 
760c4762a1bSJed Brown   tot_its = 0;
7619566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetSweeps(snes,&sweeps));
7629566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its));
7639566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
764c4762a1bSJed Brown   if (B) {
7659566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(da,&localB));
766c4762a1bSJed Brown   }
767c4762a1bSJed Brown   if (B) {
7689566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB));
7699566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB));
770c4762a1bSJed Brown   }
7719566063dSJacob Faibussowitsch   if (B) PetscCall(DMDAVecGetArrayRead(da,localB,&b));
7729566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
7739566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
7749566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localX,&x));
775c4762a1bSJed Brown   for (l=0; l<sweeps; l++) {
776c4762a1bSJed Brown     /*
777c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
778c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
779c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
780c4762a1bSJed Brown      */
7819566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
782c4762a1bSJed Brown     for (m=0; m<2; m++) {
783c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
784c4762a1bSJed Brown         for (i=xs+(m+j)%2; i<xs+xm; i+=2) {
785c4762a1bSJed Brown           PetscReal xx = i*hx,yy = j*hy;
786c4762a1bSJed Brown           if (B) bij = b[j][i];
787c4762a1bSJed Brown           else   bij = 0.;
788c4762a1bSJed Brown 
789c4762a1bSJed Brown           if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
790c4762a1bSJed Brown             /* boundary conditions are all zero Dirichlet */
791c4762a1bSJed Brown             x[j][i] = 0.0 + bij;
792c4762a1bSJed Brown           } else {
793c4762a1bSJed Brown             const PetscScalar
794c4762a1bSJed Brown               u_E = x[j][i+1],
795c4762a1bSJed Brown               u_W = x[j][i-1],
796c4762a1bSJed Brown               u_N = x[j+1][i],
797c4762a1bSJed Brown               u_S = x[j-1][i];
798c4762a1bSJed Brown             const PetscScalar
799c4762a1bSJed Brown               uy_E   = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
800c4762a1bSJed Brown               uy_W   = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
801c4762a1bSJed Brown               ux_N   = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
802c4762a1bSJed Brown               ux_S   = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]);
803c4762a1bSJed Brown             u = x[j][i];
804c4762a1bSJed Brown             for (k=0; k<its; k++) {
805c4762a1bSJed Brown               const PetscScalar
806c4762a1bSJed Brown                 ux_E   = dhx*(u_E-u),
807c4762a1bSJed Brown                 ux_W   = dhx*(u-u_W),
808c4762a1bSJed Brown                 uy_N   = dhy*(u_N-u),
809c4762a1bSJed Brown                 uy_S   = dhy*(u-u_S),
810c4762a1bSJed Brown                 e_E    = eta(user,xx,yy,ux_E,uy_E),
811c4762a1bSJed Brown                 e_W    = eta(user,xx,yy,ux_W,uy_W),
812c4762a1bSJed Brown                 e_N    = eta(user,xx,yy,ux_N,uy_N),
813c4762a1bSJed Brown                 e_S    = eta(user,xx,yy,ux_S,uy_S),
814c4762a1bSJed Brown                 de_E   = deta(user,xx,yy,ux_E,uy_E),
815c4762a1bSJed Brown                 de_W   = deta(user,xx,yy,ux_W,uy_W),
816c4762a1bSJed Brown                 de_N   = deta(user,xx,yy,ux_N,uy_N),
817c4762a1bSJed Brown                 de_S   = deta(user,xx,yy,ux_S,uy_S),
818c4762a1bSJed Brown                 newt_E = e_E+de_E*PetscSqr(ux_E),
819c4762a1bSJed Brown                 newt_W = e_W+de_W*PetscSqr(ux_W),
820c4762a1bSJed Brown                 newt_N = e_N+de_N*PetscSqr(uy_N),
821c4762a1bSJed Brown                 newt_S = e_S+de_S*PetscSqr(uy_S),
822c4762a1bSJed Brown                 uxx    = -hy * (e_E*ux_E - e_W*ux_W),
823c4762a1bSJed Brown                 uyy    = -hx * (e_N*uy_N - e_S*uy_S);
824c4762a1bSJed Brown 
825c4762a1bSJed Brown               if (sc) eu = PetscExpScalar(u);
826c4762a1bSJed Brown               else    eu = 0;
827c4762a1bSJed Brown 
828c4762a1bSJed Brown               F = uxx + uyy - sc*eu - bij;
829c4762a1bSJed Brown               if (k == 0) F0 = F;
830c4762a1bSJed Brown               J  = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*eu;
831c4762a1bSJed Brown               y  = F/J;
832c4762a1bSJed Brown               u -= y;
833c4762a1bSJed Brown               tot_its++;
834c4762a1bSJed Brown               if (atol > PetscAbsReal(PetscRealPart(F)) ||
835c4762a1bSJed Brown                   rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) ||
836c4762a1bSJed Brown                   stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) {
837c4762a1bSJed Brown                 break;
838c4762a1bSJed Brown               }
839c4762a1bSJed Brown             }
840c4762a1bSJed Brown             x[j][i] = u;
841c4762a1bSJed Brown           }
842c4762a1bSJed Brown         }
843c4762a1bSJed Brown       }
844c4762a1bSJed Brown     }
845c4762a1bSJed Brown     /*
846c4762a1bSJed Brown x     Restore vector
847c4762a1bSJed Brown      */
848c4762a1bSJed Brown   }
8499566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localX,&x));
8509566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X));
8519566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X));
8529566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(tot_its*(118.0)));
8539566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
854c4762a1bSJed Brown   if (B) {
8559566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da,localB,&b));
8569566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da,&localB));
857c4762a1bSJed Brown   }
858c4762a1bSJed Brown   PetscFunctionReturn(0);
859c4762a1bSJed Brown }
860c4762a1bSJed Brown 
861c4762a1bSJed Brown /*TEST
862c4762a1bSJed Brown 
863c4762a1bSJed Brown    test:
864c4762a1bSJed Brown       nsize: 2
865c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype NEWTON
866c4762a1bSJed Brown       requires: !single
867c4762a1bSJed Brown 
868c4762a1bSJed Brown    test:
869c4762a1bSJed Brown       suffix: 2
870c4762a1bSJed Brown       nsize: 2
871c4762a1bSJed Brown       args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -precheck 1
872c4762a1bSJed Brown       requires: !single
873c4762a1bSJed Brown 
874c4762a1bSJed Brown    test:
875c4762a1bSJed Brown       suffix: 3
876c4762a1bSJed Brown       nsize: 2
87785b95db3SBarry 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
87885b95db3SBarry Smith       requires: !single
87985b95db3SBarry Smith 
88085b95db3SBarry Smith    test:
88185b95db3SBarry Smith       suffix: 3_star
88285b95db3SBarry Smith       nsize: 2
88385b95db3SBarry 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
88485b95db3SBarry Smith       output_file: output/ex15_3.out
885c4762a1bSJed Brown       requires: !single
886c4762a1bSJed Brown 
887c4762a1bSJed Brown    test:
888c4762a1bSJed Brown       suffix: 4
889c4762a1bSJed 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
890c4762a1bSJed Brown       requires: !single
891c4762a1bSJed Brown 
892c4762a1bSJed Brown    test:
893c4762a1bSJed Brown       suffix: lag_jac
894c4762a1bSJed Brown       nsize: 4
895c4762a1bSJed 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
896c4762a1bSJed Brown       requires: !single
897c4762a1bSJed Brown 
898c4762a1bSJed Brown    test:
899c4762a1bSJed Brown       suffix: lag_pc
900c4762a1bSJed Brown       nsize: 4
901c4762a1bSJed 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
902c4762a1bSJed Brown       requires: !single
903c4762a1bSJed Brown 
904c4762a1bSJed Brown    test:
905c4762a1bSJed Brown       suffix: nleqerr
906c4762a1bSJed 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
907c4762a1bSJed Brown       requires: !single
908c4762a1bSJed Brown 
909bbc1464cSBarry Smith    test:
910bbc1464cSBarry Smith       suffix: mf
911bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator
912bbc1464cSBarry Smith       requires: !single
913bbc1464cSBarry Smith 
914bbc1464cSBarry Smith    test:
915bbc1464cSBarry Smith       suffix: mf_picard
916bbc1464cSBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_mf_operator -picard
917bbc1464cSBarry Smith       requires: !single
918bbc1464cSBarry Smith       output_file: output/ex15_mf.out
919bbc1464cSBarry Smith 
92085b95db3SBarry Smith    test:
92185b95db3SBarry Smith       suffix: fd_picard
92285b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 2  -p 3 -ksp_rtol 1.e-12  -snes_fd -picard
92385b95db3SBarry Smith       requires: !single
92485b95db3SBarry Smith 
92585b95db3SBarry Smith    test:
92685b95db3SBarry Smith       suffix: fd_color_picard
92785b95db3SBarry Smith       args: -snes_monitor_short -pc_type lu -da_refine 4  -p 3 -ksp_rtol 1.e-12  -snes_fd_color -picard
92885b95db3SBarry Smith       requires: !single
92985b95db3SBarry Smith       output_file: output/ex15_mf.out
93085b95db3SBarry Smith 
931c4762a1bSJed Brown TEST*/
932