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