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