Actual source code: ex15.c

petsc-3.4.5 2014-06-29
  1: static const char help[] = "p-Bratu nonlinear PDE in 2d.\n\
  2: We solve the  p-Laplacian (nonlinear diffusion) combined with\n\
  3: the Bratu (solid fuel ignition) nonlinearity in a 2D rectangular\n\
  4: domain, using distributed arrays (DAs) to partition the parallel grid.\n\
  5: The command line options include:\n\
  6:   -p <2>: `p' in p-Laplacian term\n\
  7:   -epsilon <1e-05>: Strain-regularization in p-Laplacian\n\
  8:   -lambda <6>: Bratu parameter\n\
  9:   -blocks <bx,by>: number of coefficient interfaces in x and y direction\n\
 10:   -kappa <1e-3>: diffusivity in odd regions\n\
 11: \n";


The $p$-Bratu problem is a combination of the $p$-Laplacian (nonlinear diffusion) and the Brutu solid fuel ignition problem.
This problem is modeled by the partial differential equation

\begin{equation*}
-\nabla\cdot (\eta \nabla u) - \lambda \exp(u) = 0
\end{equation*}

on $\Omega = (-1,1)^2$ with closure

\begin{align*}
\eta(\gamma) &= (\epsilon^2 + \gamma)^{(p-2)/2} & \gamma &= \frac 1 2 |\nabla u|^2
\end{align*}

and boundary conditions $u = 0$ for $(x,y) \in \partial \Omega$

A 9-point finite difference stencil is used to discretize
the boundary value problem to obtain a nonlinear system of equations.
This would be a 5-point stencil if not for the $p$-Laplacian's nonlinearity.
 35: /*
 36:       mpiexec -n 2 ./ex15 -snes_monitor -ksp_monitor log_summary
 37: */

 39: /*
 40:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 41:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 42:    file automatically includes:
 43:      petsc.h       - base PETSc routines   petscvec.h - vectors
 44:      petscsys.h    - system routines       petscmat.h - matrices
 45:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 46:      petscviewer.h - viewers               petscpc.h  - preconditioners
 47:      petscksp.h   - linear solvers
 48: */
 49: #include <petscdmda.h>
 50: #include <petscsnes.h>

 52: /* These functions _should_ be internal, but currently have a reverse dependency so cannot be set with
 53:  * DMDASNESSetPicardLocal.  This hack needs to be fixed in PETSc. */
 54: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES,Vec,Vec,void*);
 55: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);

 57: typedef enum {JAC_BRATU,JAC_PICARD,JAC_STAR,JAC_NEWTON} JacType;
 58: static const char *const JacTypes[] = {"BRATU","PICARD","STAR","NEWTON","JacType","JAC_",0};

 60: /*
 61:    User-defined application context - contains data needed by the
 62:    application-provided call-back routines, FormJacobianLocal() and
 63:    FormFunctionLocal().
 64: */
 65: typedef struct {
 66:   PassiveReal lambda;         /* Bratu parameter */
 67:   PassiveReal p;              /* Exponent in p-Laplacian */
 68:   PassiveReal epsilon;        /* Regularization */
 69:   PassiveReal source;         /* Source term */
 70:   JacType     jtype;          /* What type of Jacobian to assemble */
 71:   PetscBool   picard;
 72:   PetscInt    blocks[2];
 73:   PetscReal   kappa;
 74:   PetscInt    initial;        /* initial conditions type */
 75: } AppCtx;

 77: /*
 78:    User-defined routines
 79: */
 80: static PetscErrorCode FormRHS(AppCtx*,DM,Vec);
 81: static PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
 82: static PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
 83: static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
 84: static PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,Mat,MatStructure*,AppCtx*);
 85: static PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);

 87: typedef struct _n_PreCheck *PreCheck;
 88: struct _n_PreCheck {
 89:   MPI_Comm    comm;
 90:   PetscReal   angle;
 91:   Vec         Ylast;
 92:   PetscViewer monitor;
 93: };
 94: PetscErrorCode PreCheckCreate(MPI_Comm,PreCheck*);
 95: PetscErrorCode PreCheckDestroy(PreCheck*);
 96: PetscErrorCode PreCheckFunction(SNESLineSearch,Vec,Vec,PetscBool*,void*);
 97: PetscErrorCode PreCheckSetFromOptions(PreCheck);

101: int main(int argc,char **argv)
102: {
103:   SNES                snes;                    /* nonlinear solver */
104:   Vec                 x,r,b;                   /* solution, residual, rhs vectors */
105:   Mat                 A,B;                     /* Jacobian and preconditioning matrices */
106:   AppCtx              user;                    /* user-defined work context */
107:   PetscInt            its;                     /* iterations for convergence */
108:   SNESConvergedReason reason;                  /* Check convergence */
109:   PetscBool           alloc_star;              /* Only allocate for the STAR stencil  */
110:   PetscBool           write_output;
111:   char                filename[PETSC_MAX_PATH_LEN] = "ex15.vts";
112:   PetscReal           bratu_lambda_max             = 6.81,bratu_lambda_min = 0.;
113:   DM                  da,dastar;               /* distributed array data structure */
114:   PreCheck            precheck = NULL;    /* precheck context for version in this file */
115:   PetscInt            use_precheck;      /* 0=none, 1=version in this file, 2=SNES-provided version */
116:   PetscReal           precheck_angle;    /* When manually setting the SNES-provided precheck function */
117:   PetscErrorCode      ierr;
118:   SNESLineSearch      linesearch;

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:      Initialize program
122:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

124:   PetscInitialize(&argc,&argv,(char*)0,help);

126:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127:      Initialize problem parameters
128:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129:   user.lambda    = 0.0; user.p = 2.0; user.epsilon = 1e-5; user.source = 0.1; user.jtype = JAC_NEWTON;user.initial=-1;
130:   user.blocks[0] = 1; user.blocks[1] = 1; user.kappa = 1e-3;
131:   alloc_star     = PETSC_FALSE;
132:   use_precheck   = 0; precheck_angle = 10.;
133:   user.picard    = PETSC_FALSE;
134:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"p-Bratu options",__FILE__);
135:   {
136:     PetscInt two=2;
137:     PetscOptionsReal("-lambda","Bratu parameter","",user.lambda,&user.lambda,NULL);
138:     PetscOptionsReal("-p","Exponent `p' in p-Laplacian","",user.p,&user.p,NULL);
139:     PetscOptionsReal("-epsilon","Strain-regularization in p-Laplacian","",user.epsilon,&user.epsilon,NULL);
140:     PetscOptionsReal("-source","Constant source term","",user.source,&user.source,NULL);
141:     PetscOptionsEnum("-jtype","Jacobian approximation to assemble","",JacTypes,(PetscEnum)user.jtype,(PetscEnum*)&user.jtype,NULL);
142:     PetscOptionsName("-picard","Solve with defect-correction Picard iteration","",&user.picard);
143:     if (user.picard) {user.jtype = JAC_PICARD; user.p = 3;}
144:     PetscOptionsBool("-alloc_star","Allocate for STAR stencil (5-point)","",alloc_star,&alloc_star,NULL);
145:     PetscOptionsInt("-precheck","Use a pre-check correction intended for use with Picard iteration 1=this version, 2=library","",use_precheck,&use_precheck,NULL);
146:     PetscOptionsInt("-initial","Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)","",user.initial,&user.initial,NULL);
147:     if (use_precheck == 2) {    /* Using library version, get the angle */
148:       PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck_angle,&precheck_angle,NULL);
149:     }
150:     PetscOptionsIntArray("-blocks","number of coefficient interfaces in x and y direction","",user.blocks,&two,NULL);
151:     PetscOptionsReal("-kappa","diffusivity in odd regions","",user.kappa,&user.kappa,NULL);
152:     PetscOptionsString("-o","Output solution in vts format","",filename,filename,sizeof(filename),&write_output);
153:   }
154:   PetscOptionsEnd();
155:   if (user.lambda > bratu_lambda_max || user.lambda < bratu_lambda_min) {
156:     PetscPrintf(PETSC_COMM_WORLD,"WARNING: lambda %g out of range for p=2\n",user.lambda);
157:   }

159:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160:      Create nonlinear solver context
161:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162:   SNESCreate(PETSC_COMM_WORLD,&snes);

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:      Create distributed array (DMDA) to manage parallel grid and vectors
166:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167:   DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE,
168:                       1,1,NULL,NULL,&da);
169:   DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,
170:                       1,1,NULL,NULL,&dastar);


173:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174:      Extract global vectors from DM; then duplicate for remaining
175:      vectors that are the same types
176:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177:   DMCreateGlobalVector(da,&x);
178:   VecDuplicate(x,&r);
179:   VecDuplicate(x,&b);

181:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182:      Create matrix data structure; set Jacobian evaluation routine

184:      Set Jacobian matrix data structure and default Jacobian evaluation
185:      routine. User can override with:
186:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
187:                 (unless user explicitly sets preconditioner)
188:      -snes_mf_operator : form preconditioning matrix as set by the user,
189:                          but use matrix-free approx for Jacobian-vector
190:                          products within Newton-Krylov method

192:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193:   /* B can be type of MATAIJ,MATBAIJ or MATSBAIJ */
194:   DMCreateMatrix(alloc_star ? dastar : da,MATAIJ,&B);
195:   A    = B;

197:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198:      Set local function evaluation routine
199:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200:   DMSetApplicationContext(da, &user);
201:   SNESSetDM(snes,da);
202:   if (user.picard) {
203:     /*
204:         This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access
205:         the SNES to set it
206:     */
207:     DMDASNESSetPicardLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionPicardLocal,
208:                                   (PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*))FormJacobianLocal,&user);
209:     SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);
210:     SNESSetJacobian(snes,NULL,NULL,SNESPicardComputeJacobian,&user);
211:   } else {
212:     DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);
213:     DMDASNESSetJacobianLocal(da,(PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*))FormJacobianLocal,&user);
214:   }


217:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218:      Customize nonlinear solver; set runtime options
219:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220:   SNESSetFromOptions(snes);
221:   SNESSetGS(snes,NonlinearGS,&user);
222:   SNESGetLineSearch(snes, &linesearch);
223:   /* Set up the precheck context if requested */
224:   if (use_precheck == 1) {      /* Use the precheck routines in this file */
225:     PreCheckCreate(PETSC_COMM_WORLD,&precheck);
226:     PreCheckSetFromOptions(precheck);
227:     SNESLineSearchSetPreCheck(linesearch,PreCheckFunction,precheck);
228:   } else if (use_precheck == 2) { /* Use the version provided by the library */
229:     SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&precheck_angle);
230:   }

232:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233:      Evaluate initial guess
234:      Note: The user should initialize the vector, x, with the initial guess
235:      for the nonlinear solver prior to calling SNESSolve().  In particular,
236:      to employ an initial guess of zero, the user should explicitly set
237:      this vector to zero by calling VecSet().
238:   */

240:   FormInitialGuess(&user,da,x);
241:   FormRHS(&user,da,b);

243:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244:      Solve nonlinear system
245:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246:   SNESSolve(snes,b,x);
247:   SNESGetIterationNumber(snes,&its);
248:   SNESGetConvergedReason(snes,&reason);

250:   PetscPrintf(PETSC_COMM_WORLD,"%s Number of nonlinear iterations = %D\n",SNESConvergedReasons[reason],its);

252:   if (write_output) {
253:     PetscViewer viewer;
254:     PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer);
255:     VecView(x,viewer);
256:     PetscViewerDestroy(&viewer);
257:   }

259:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260:      Free work space.  All PETSc objects should be destroyed when they
261:      are no longer needed.
262:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

264:   if (A != B) {
265:     MatDestroy(&A);
266:   }
267:   MatDestroy(&B);
268:   VecDestroy(&x);
269:   VecDestroy(&r);
270:   VecDestroy(&b);
271:   SNESDestroy(&snes);
272:   DMDestroy(&da);
273:   DMDestroy(&dastar);
274:   PreCheckDestroy(&precheck);
275:   PetscFinalize();
276:   return 0;
277: }

279: /* ------------------------------------------------------------------- */
282: /*
283:    FormInitialGuess - Forms initial approximation.

285:    Input Parameters:
286:    user - user-defined application context
287:    X - vector

289:    Output Parameter:
290:    X - vector
291:  */
292: static PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
293: {
294:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
296:   PetscReal      temp1,temp,hx,hy;
297:   PetscScalar    **x;

300:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
301:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

303:   hx    = 1.0/(PetscReal)(Mx-1);
304:   hy    = 1.0/(PetscReal)(My-1);
305:   temp1 = user->lambda / (user->lambda + 1.);

307:   /*
308:      Get a pointer to vector data.
309:        - For default PETSc vectors, VecGetArray() returns a pointer to
310:          the data array.  Otherwise, the routine is implementation dependent.
311:        - You MUST call VecRestoreArray() when you no longer need access to
312:          the array.
313:   */
314:   DMDAVecGetArray(da,X,&x);

316:   /*
317:      Get local grid boundaries (for 2-dimensional DA):
318:        xs, ys   - starting grid indices (no ghost points)
319:        xm, ym   - widths of local grid (no ghost points)

321:   */
322:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

324:   /*
325:      Compute initial guess over the locally owned part of the grid
326:   */
327:   for (j=ys; j<ys+ym; j++) {
328:     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
329:     for (i=xs; i<xs+xm; i++) {
330:       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
331:         /* boundary conditions are all zero Dirichlet */
332:         x[j][i] = 0.0;
333:       } else {
334:         if (user->initial == -1) {
335:           if (user->lambda != 0) {
336:             x[j][i] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
337:           } else {
338:             /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting
339:              * with an exact solution. */
340:             const PetscReal
341:               xx = 2*(PetscReal)i/(Mx-1) - 1,
342:               yy = 2*(PetscReal)j/(My-1) - 1;
343:             x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
344:           }
345:         } else if (user->initial == 0) {
346:           x[j][i] = 0.;
347:         } else if (user->initial == 1) {
348:           const PetscReal
349:             xx = 2*(PetscReal)i/(Mx-1) - 1,
350:             yy = 2*(PetscReal)j/(My-1) - 1;
351:           x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy;
352:         } else {
353:           if (user->lambda != 0) {
354:             x[j][i] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
355:           } else {
356:             x[j][i] = 0.5*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
357:           }
358:         }
359:       }
360:     }
361:   }
362:   /*
363:      Restore vector
364:   */
365:   DMDAVecRestoreArray(da,X,&x);
366:   return(0);
367: }

371: /*
372:    FormRHS - Forms constant RHS for the problem.

374:    Input Parameters:
375:    user - user-defined application context
376:    B - RHS vector

378:    Output Parameter:
379:    B - vector
380:  */
381: static PetscErrorCode FormRHS(AppCtx *user,DM da,Vec B)
382: {
383:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
385:   PetscReal      hx,hy;
386:   PetscScalar    **b;

389:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
390:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

392:   hx    = 1.0/(PetscReal)(Mx-1);
393:   hy    = 1.0/(PetscReal)(My-1);
394:   DMDAVecGetArray(da,B,&b);
395:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
396:   for (j=ys; j<ys+ym; j++) {
397:     for (i=xs; i<xs+xm; i++) {
398:       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
399:         b[j][i] = 0.0;
400:       } else {
401:         b[j][i] = hx*hy*user->source;
402:       }
403:     }
404:   }
405:   DMDAVecRestoreArray(da,B,&b);
406:   return(0);
407: }

409: PETSC_STATIC_INLINE PetscReal kappa(const AppCtx *ctx,PetscReal x,PetscReal y)
410: {
411:   return (((PetscInt)(x*ctx->blocks[0])) + ((PetscInt)(y*ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0;
412: }
413: /* p-Laplacian diffusivity */
414: PETSC_STATIC_INLINE PetscScalar eta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
415: {
416:   return kappa(ctx,x,y) * PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-2.));
417: }
418: PETSC_STATIC_INLINE PetscScalar deta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy)
419: {
420:   return (ctx->p == 2)
421:          ? 0
422:          : kappa(ctx,x,y)*PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-4)) * 0.5 * (ctx->p-2.);
423: }


426: /* ------------------------------------------------------------------- */
429: /*
430:    FormFunctionLocal - Evaluates nonlinear function, F(x).
431:  */
432: static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
433: {
434:   PetscReal      hx,hy,dhx,dhy,sc;
435:   PetscInt       i,j;
436:   PetscScalar    eu;


441:   hx     = 1.0/(PetscReal)(info->mx-1);
442:   hy     = 1.0/(PetscReal)(info->my-1);
443:   sc     = hx*hy*user->lambda;
444:   dhx    = 1/hx;
445:   dhy    = 1/hy;
446:   /*
447:      Compute function over the locally owned part of the grid
448:   */
449:   for (j=info->ys; j<info->ys+info->ym; j++) {
450:     for (i=info->xs; i<info->xs+info->xm; i++) {
451:       PetscReal xx = i*hx,yy = j*hy;
452:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
453:         f[j][i] = x[j][i];
454:       } else {
455:         const PetscScalar
456:           u    = x[j][i],
457:           ux_E = dhx*(x[j][i+1]-x[j][i]),
458:           uy_E = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
459:           ux_W = dhx*(x[j][i]-x[j][i-1]),
460:           uy_W = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
461:           ux_N = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
462:           uy_N = dhy*(x[j+1][i]-x[j][i]),
463:           ux_S = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
464:           uy_S = dhy*(x[j][i]-x[j-1][i]),
465:           e_E  = eta(user,xx,yy,ux_E,uy_E),
466:           e_W  = eta(user,xx,yy,ux_W,uy_W),
467:           e_N  = eta(user,xx,yy,ux_N,uy_N),
468:           e_S  = eta(user,xx,yy,ux_S,uy_S),
469:           uxx  = -hy * (e_E*ux_E - e_W*ux_W),
470:           uyy  = -hx * (e_N*uy_N - e_S*uy_S);
471:         if (sc) eu = PetscExpScalar(u);
472:         else    eu = 0.;
473:         /** For p=2, these terms decay to:
474:         * uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx
475:         * uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy
476:         **/
477:         f[j][i] = uxx + uyy - sc*eu;
478:       }
479:     }
480:   }
481:   PetscLogFlops(info->xm*info->ym*(72.0));
482:   return(0);
483: }

487: /*
488:     This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation,
489:     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)

491: */
492: static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
493: {
494:   PetscReal hx,hy,sc;
495:   PetscInt  i,j;

499:   hx     = 1.0/(PetscReal)(info->mx-1);
500:   hy     = 1.0/(PetscReal)(info->my-1);
501:   sc     = hx*hy*user->lambda;
502:   /*
503:      Compute function over the locally owned part of the grid
504:   */
505:   for (j=info->ys; j<info->ys+info->ym; j++) {
506:     for (i=info->xs; i<info->xs+info->xm; i++) {
507:       if (!(i == 0 || j == 0 || i == info->mx-1 || j == info->my-1)) {
508:         const PetscScalar u = x[j][i];
509:         f[j][i] = sc*PetscExpScalar(u);
510:       }
511:     }
512:   }
513:   PetscLogFlops(info->xm*info->ym*2.0);
514:   return(0);
515: }

519: /*
520:    FormJacobianLocal - Evaluates Jacobian matrix.
521: */
522: static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat J,Mat B,MatStructure *str,AppCtx *user)
523: {
525:   PetscInt       i,j;
526:   MatStencil     col[9],row;
527:   PetscScalar    v[9];
528:   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;

531:   hx    = 1.0/(PetscReal)(info->mx-1);
532:   hy    = 1.0/(PetscReal)(info->my-1);
533:   sc    = hx*hy*user->lambda;
534:   dhx   = 1/hx;
535:   dhy   = 1/hy;
536:   hxdhy = hx/hy;
537:   hydhx = hy/hx;

539:   /*
540:      Compute entries for the locally owned part of the Jacobian.
541:       - PETSc parallel matrix formats are partitioned by
542:         contiguous chunks of rows across the processors.
543:       - Each processor needs to insert only elements that it owns
544:         locally (but any non-local elements will be sent to the
545:         appropriate processor during matrix assembly).
546:       - Here, we set all entries for a particular row at once.
547:   */
548:   for (j=info->ys; j<info->ys+info->ym; j++) {
549:     for (i=info->xs; i<info->xs+info->xm; i++) {
550:       PetscReal xx = i*hx,yy = j*hy;
551:       row.j = j; row.i = i;
552:       /* boundary points */
553:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
554:         v[0] = 1.0;
555:         MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);
556:       } else {
557:         /* interior grid points */
558:         const PetscScalar
559:           ux_E     = dhx*(x[j][i+1]-x[j][i]),
560:           uy_E     = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
561:           ux_W     = dhx*(x[j][i]-x[j][i-1]),
562:           uy_W     = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
563:           ux_N     = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
564:           uy_N     = dhy*(x[j+1][i]-x[j][i]),
565:           ux_S     = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]),
566:           uy_S     = dhy*(x[j][i]-x[j-1][i]),
567:           u        = x[j][i],
568:           e_E      = eta(user,xx,yy,ux_E,uy_E),
569:           e_W      = eta(user,xx,yy,ux_W,uy_W),
570:           e_N      = eta(user,xx,yy,ux_N,uy_N),
571:           e_S      = eta(user,xx,yy,ux_S,uy_S),
572:           de_E     = deta(user,xx,yy,ux_E,uy_E),
573:           de_W     = deta(user,xx,yy,ux_W,uy_W),
574:           de_N     = deta(user,xx,yy,ux_N,uy_N),
575:           de_S     = deta(user,xx,yy,ux_S,uy_S),
576:           skew_E   = de_E*ux_E*uy_E,
577:           skew_W   = de_W*ux_W*uy_W,
578:           skew_N   = de_N*ux_N*uy_N,
579:           skew_S   = de_S*ux_S*uy_S,
580:           cross_EW = 0.25*(skew_E - skew_W),
581:           cross_NS = 0.25*(skew_N - skew_S),
582:           newt_E   = e_E+de_E*PetscSqr(ux_E),
583:           newt_W   = e_W+de_W*PetscSqr(ux_W),
584:           newt_N   = e_N+de_N*PetscSqr(uy_N),
585:           newt_S   = e_S+de_S*PetscSqr(uy_S);
586:         /* interior grid points */
587:         switch (user->jtype) {
588:         case JAC_BRATU:
589:           /* Jacobian from p=2 */
590:           v[0] = -hxdhy;                                           col[0].j = j-1;   col[0].i = i;
591:           v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
592:           v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(u);       col[2].j = row.j; col[2].i = row.i;
593:           v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
594:           v[4] = -hxdhy;                                           col[4].j = j+1;   col[4].i = i;
595:           MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES);
596:           break;
597:         case JAC_PICARD:
598:           /* Jacobian arising from Picard linearization */
599:           v[0] = -hxdhy*e_S;                                           col[0].j = j-1;   col[0].i = i;
600:           v[1] = -hydhx*e_W;                                           col[1].j = j;     col[1].i = i-1;
601:           v[2] = (e_W+e_E)*hydhx + (e_S+e_N)*hxdhy;                    col[2].j = row.j; col[2].i = row.i;
602:           v[3] = -hydhx*e_E;                                           col[3].j = j;     col[3].i = i+1;
603:           v[4] = -hxdhy*e_N;                                           col[4].j = j+1;   col[4].i = i;
604:           MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES);
605:           break;
606:         case JAC_STAR:
607:           /* Full Jacobian, but only a star stencil */
608:           col[0].j = j-1; col[0].i = i;
609:           col[1].j = j;   col[1].i = i-1;
610:           col[2].j = j;   col[2].i = i;
611:           col[3].j = j;   col[3].i = i+1;
612:           col[4].j = j+1; col[4].i = i;
613:           v[0]     = -hxdhy*newt_S + cross_EW;
614:           v[1]     = -hydhx*newt_W + cross_NS;
615:           v[2]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
616:           v[3]     = -hydhx*newt_E - cross_NS;
617:           v[4]     = -hxdhy*newt_N - cross_EW;
618:           MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES);
619:           break;
620:         case JAC_NEWTON:
621:           /** The Jacobian is
622:           *
623:           * -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u
624:           *
625:           **/
626:           col[0].j = j-1; col[0].i = i-1;
627:           col[1].j = j-1; col[1].i = i;
628:           col[2].j = j-1; col[2].i = i+1;
629:           col[3].j = j;   col[3].i = i-1;
630:           col[4].j = j;   col[4].i = i;
631:           col[5].j = j;   col[5].i = i+1;
632:           col[6].j = j+1; col[6].i = i-1;
633:           col[7].j = j+1; col[7].i = i;
634:           col[8].j = j+1; col[8].i = i+1;
635:           v[0]     = -0.25*(skew_S + skew_W);
636:           v[1]     = -hxdhy*newt_S + cross_EW;
637:           v[2]     =  0.25*(skew_S + skew_E);
638:           v[3]     = -hydhx*newt_W + cross_NS;
639:           v[4]     = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u);
640:           v[5]     = -hydhx*newt_E - cross_NS;
641:           v[6]     =  0.25*(skew_N + skew_W);
642:           v[7]     = -hxdhy*newt_N - cross_EW;
643:           v[8]     = -0.25*(skew_N + skew_E);
644:           MatSetValuesStencil(B,1,&row,9,col,v,INSERT_VALUES);
645:           break;
646:         default:
647:           SETERRQ1(PetscObjectComm((PetscObject)info->da),PETSC_ERR_SUP,"Jacobian type %d not implemented",user->jtype);
648:         }
649:       }
650:     }
651:   }

653:   /*
654:      Assemble matrix, using the 2-step process:
655:        MatAssemblyBegin(), MatAssemblyEnd().
656:   */
657:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
658:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

660:   if (J != B) {
661:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
662:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
663:   }
664:   *str = SAME_NONZERO_PATTERN;

666:   /*
667:      Tell the matrix we will never add a new nonzero location to the
668:      matrix. If we do, it will generate an error.
669:   */
670:   if (user->jtype == JAC_NEWTON) {
671:     PetscLogFlops(info->xm*info->ym*(131.0));
672:   }
673:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
674:   return(0);
675: }

677: /***********************************************************
678:  * PreCheck implementation
679:  ***********************************************************/
682: PetscErrorCode PreCheckSetFromOptions(PreCheck precheck)
683: {
685:   PetscBool      flg;

688:   PetscOptionsBegin(precheck->comm,NULL,"PreCheck Options","none");
689:   PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck->angle,&precheck->angle,NULL);
690:   flg  = PETSC_FALSE;
691:   PetscOptionsBool("-precheck_monitor","Monitor choices made by precheck routine","",flg,&flg,NULL);
692:   if (flg) {
693:     PetscViewerASCIIOpen(precheck->comm,"stdout",&precheck->monitor);
694:   }
695:   PetscOptionsEnd();
696:   return(0);
697: }

701: /*
702:   Compare the direction of the current and previous step, modify the current step accordingly
703: */
704: PetscErrorCode PreCheckFunction(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed, void *ctx)
705: {
707:   PreCheck       precheck;
708:   Vec            Ylast;
709:   PetscScalar    dot;
710:   PetscInt       iter;
711:   PetscReal      ynorm,ylastnorm,theta,angle_radians;
712:   SNES           snes;

715:   SNESLineSearchGetSNES(linesearch, &snes);
716:   precheck = (PreCheck)ctx;
717:   if (!precheck->Ylast) {VecDuplicate(Y,&precheck->Ylast);}
718:   Ylast = precheck->Ylast;
719:   SNESGetIterationNumber(snes,&iter);
720:   if (iter < 1) {
721:     VecCopy(Y,Ylast);
722:     *changed = PETSC_FALSE;
723:     return(0);
724:   }

726:   VecDot(Y,Ylast,&dot);
727:   VecNorm(Y,NORM_2,&ynorm);
728:   VecNorm(Ylast,NORM_2,&ylastnorm);
729:   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
730:   theta         = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
731:   angle_radians = precheck->angle * PETSC_PI / 180.;
732:   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
733:     /* Modify the step Y */
734:     PetscReal alpha,ydiffnorm;
735:     VecAXPY(Ylast,-1.0,Y);
736:     VecNorm(Ylast,NORM_2,&ydiffnorm);
737:     alpha = ylastnorm / ydiffnorm;
738:     VecCopy(Y,Ylast);
739:     VecScale(Y,alpha);
740:     if (precheck->monitor) {
741:       PetscViewerASCIIPrintf(precheck->monitor,"Angle %E degrees less than threshold %G, corrected step by alpha=%G\n",theta*180./PETSC_PI,precheck->angle,alpha);
742:     }
743:   } else {
744:     VecCopy(Y,Ylast);
745:     *changed = PETSC_FALSE;
746:     if (precheck->monitor) {
747:       PetscViewerASCIIPrintf(precheck->monitor,"Angle %E degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,precheck->angle);
748:     }
749:   }
750:   return(0);
751: }

755: PetscErrorCode PreCheckDestroy(PreCheck *precheck)
756: {

760:   if (!*precheck) return(0);
761:   VecDestroy(&(*precheck)->Ylast);
762:   PetscViewerDestroy(&(*precheck)->monitor);
763:   PetscFree(*precheck);
764:   return(0);
765: }

769: PetscErrorCode PreCheckCreate(MPI_Comm comm,PreCheck *precheck)
770: {

774:   PetscMalloc(sizeof(struct _n_PreCheck),precheck);
775:   PetscMemzero(*precheck,sizeof(struct _n_PreCheck));

777:   (*precheck)->comm  = comm;
778:   (*precheck)->angle = 10.;     /* only active if angle is less than 10 degrees */
779:   return(0);
780: }

784: /*
785:       Applies some sweeps on nonlinear Gauss-Seidel on each process

787:  */
788: PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx)
789: {
790:   PetscInt       i,j,k,xs,ys,xm,ym,its,tot_its,sweeps,l,m;
792:   PetscReal      hx,hy,hxdhy,hydhx,dhx,dhy,sc;
793:   PetscScalar    **x,**b,bij,F,F0=0,J,y,u,eu;
794:   PetscReal      atol,rtol,stol;
795:   DM             da;
796:   AppCtx         *user = (AppCtx*)ctx;
797:   Vec            localX,localB;
798:   DMDALocalInfo  info;

801:   SNESGetDM(snes,&da);
802:   DMDAGetLocalInfo(da,&info);

804:   hx     = 1.0/(PetscReal)(info.mx-1);
805:   hy     = 1.0/(PetscReal)(info.my-1);
806:   sc     = hx*hy*user->lambda;
807:   dhx    = 1/hx;
808:   dhy    = 1/hy;
809:   hxdhy  = hx/hy;
810:   hydhx  = hy/hx;

812:   tot_its = 0;
813:   SNESGSGetSweeps(snes,&sweeps);
814:   SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);
815:   DMGetLocalVector(da,&localX);
816:   if (B) {
817:     DMGetLocalVector(da,&localB);
818:   }
819:   if (B) {
820:     DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB);
821:     DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB);
822:   }
823:   if (B) DMDAVecGetArray(da,localB,&b);
824:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
825:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
826:   DMDAVecGetArray(da,localX,&x);
827:   for (l=0; l<sweeps; l++) {
828:     /*
829:      Get local grid boundaries (for 2-dimensional DMDA):
830:      xs, ys   - starting grid indices (no ghost points)
831:      xm, ym   - widths of local grid (no ghost points)
832:      */
833:     DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
834:     for (m=0; m<2; m++) {
835:       for (j=ys; j<ys+ym; j++) {
836:         for (i=xs+(m+j)%2; i<xs+xm; i+=2) {
837:           PetscReal xx = i*hx,yy = j*hy;
838:           if (B) bij = b[j][i];
839:           else   bij = 0.;

841:           if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
842:             /* boundary conditions are all zero Dirichlet */
843:             x[j][i] = 0.0 + bij;
844:           } else {
845:             const PetscScalar
846:               u_E = x[j][i+1],
847:               u_W = x[j][i-1],
848:               u_N = x[j+1][i],
849:               u_S = x[j-1][i];
850:             const PetscScalar
851:               uy_E   = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]),
852:               uy_W   = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]),
853:               ux_N   = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]),
854:               ux_S   = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]);
855:             u = x[j][i];
856:             for (k=0; k<its; k++) {
857:               const PetscScalar
858:                 ux_E   = dhx*(u_E-u),
859:                 ux_W   = dhx*(u-u_W),
860:                 uy_N   = dhy*(u_N-u),
861:                 uy_S   = dhy*(u-u_S),
862:                 e_E    = eta(user,xx,yy,ux_E,uy_E),
863:                 e_W    = eta(user,xx,yy,ux_W,uy_W),
864:                 e_N    = eta(user,xx,yy,ux_N,uy_N),
865:                 e_S    = eta(user,xx,yy,ux_S,uy_S),
866:                 de_E   = deta(user,xx,yy,ux_E,uy_E),
867:                 de_W   = deta(user,xx,yy,ux_W,uy_W),
868:                 de_N   = deta(user,xx,yy,ux_N,uy_N),
869:                 de_S   = deta(user,xx,yy,ux_S,uy_S),
870:                 newt_E = e_E+de_E*PetscSqr(ux_E),
871:                 newt_W = e_W+de_W*PetscSqr(ux_W),
872:                 newt_N = e_N+de_N*PetscSqr(uy_N),
873:                 newt_S = e_S+de_S*PetscSqr(uy_S),
874:                 uxx    = -hy * (e_E*ux_E - e_W*ux_W),
875:                 uyy    = -hx * (e_N*uy_N - e_S*uy_S);

877:               if (sc) eu = PetscExpScalar(u);
878:               else    eu = 0;

880:               F = uxx + uyy - sc*eu - bij;
881:               if (k == 0) F0 = F;
882:               J  = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*eu;
883:               y  = F/J;
884:               u -= y;
885:               tot_its++;
886:               if (atol > PetscAbsReal(PetscRealPart(F)) ||
887:                   rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) ||
888:                   stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) {
889:                 break;
890:               }
891:             }
892:             x[j][i] = u;
893:           }
894:         }
895:       }
896:     }
897:     /*
898: x     Restore vector
899:      */
900:   }
901:   DMDAVecRestoreArray(da,localX,&x);
902:   DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);
903:   DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);
904:   PetscLogFlops(tot_its*(118.0));
905:   DMRestoreLocalVector(da,&localX);
906:   if (B) {
907:     DMDAVecRestoreArray(da,localB,&b);
908:     DMRestoreLocalVector(da,&localB);
909:   }
910:   return(0);
911: }