Actual source code: ex3.c
2: static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel.\n\
3: This example employs a user-defined monitoring routine and optionally a user-defined\n\
4: routine to check candidate iterates produced by line search routines. This code also\n\
5: demonstrates use of the macro __FUNCT__ to define routine names for use in error handling.\n\
6: The command line options include:\n\
7: -check_iterates : activate checking of iterates\n\
8: -check_tol <tol>: set tolerance for iterate checking\n\n";
10: /*T
11: Concepts: SNES^basic parallel example
12: Concepts: SNES^setting a user-defined monitoring routine
13: Concepts: error handling^using the macro __FUNCT__ to define routine names;
14: Processors: n
15: T*/
17: /*
18: Include "petscdraw.h" so that we can use distributed arrays (DAs).
19: Include "petscdraw.h" so that we can use PETSc drawing routines.
20: Include "petscsnes.h" so that we can use SNES solvers. Note that this
21: file automatically includes:
22: petsc.h - base PETSc routines petscvec.h - vectors
23: petscsys.h - system routines petscmat.h - matrices
24: petscis.h - index sets petscksp.h - Krylov subspace methods
25: petscviewer.h - viewers petscpc.h - preconditioners
26: petscksp.h - linear solvers
27: */
29: #include petscda.h
30: #include petscsnes.h
32: /*
33: User-defined routines. Note that immediately before each routine below,
34: we define the macro __FUNCT__ to be a string containing the routine name.
35: If defined, this macro is used in the PETSc error handlers to provide a
36: complete traceback of routine names. All PETSc library routines use this
37: macro, and users can optionally employ it as well in their application
38: codes. Note that users can get a traceback of PETSc errors regardless of
39: whether they define __FUNCT__ in application codes; this macro merely
40: provides the added traceback detail of the application routine names.
41: */
42: PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
43: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
44: PetscErrorCode FormInitialGuess(Vec);
45: PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void *);
46: PetscErrorCode StepCheck(SNES,void *,Vec,PetscTruth *);
48: /*
49: User-defined application context
50: */
51: typedef struct {
52: DA da; /* distributed array */
53: Vec F; /* right-hand-side of PDE */
54: PetscMPIInt rank; /* rank of processor */
55: PetscMPIInt size; /* size of communicator */
56: PetscReal h; /* mesh spacing */
57: } ApplicationCtx;
59: /*
60: User-defined context for monitoring
61: */
62: typedef struct {
63: PetscViewer viewer;
64: } MonitorCtx;
66: /*
67: User-defined context for checking candidate iterates that are
68: determined by line search methods
69: */
70: typedef struct {
71: Vec last_step; /* previous iterate */
72: PetscReal tolerance; /* tolerance for changes between successive iterates */
73: } StepCheckCtx;
77: int main(int argc,char **argv)
78: {
79: SNES snes; /* SNES context */
80: Mat J; /* Jacobian matrix */
81: ApplicationCtx ctx; /* user-defined context */
82: Vec x,r,U,F; /* vectors */
83: MonitorCtx monP; /* monitoring context */
84: StepCheckCtx checkP; /* step-checking context */
85: PetscTruth step_check; /* flag indicating whether we're checking
86: candidate iterates */
87: PetscScalar xp,*FF,*UU,none = -1.0;
89: PetscInt its,N = 5,i,maxit,maxf,xs,xm;
90: PetscReal abstol,rtol,stol,norm;
92: PetscInitialize(&argc,&argv,(char *)0,help);
93: MPI_Comm_rank(PETSC_COMM_WORLD,&ctx.rank);
94: MPI_Comm_size(PETSC_COMM_WORLD,&ctx.size);
95: PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);
96: ctx.h = 1.0/(N-1);
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Create nonlinear solver context
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: SNESCreate(PETSC_COMM_WORLD,&snes);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Create vector data structures; set function evaluation routine
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: /*
109: Create distributed array (DA) to manage parallel grid and vectors
110: */
111: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,N,1,1,PETSC_NULL,&ctx.da);
113: /*
114: Extract global and local vectors from DA; then duplicate for remaining
115: vectors that are the same types
116: */
117: DACreateGlobalVector(ctx.da,&x);
118: VecDuplicate(x,&r);
119: VecDuplicate(x,&F); ctx.F = F;
120: VecDuplicate(x,&U);
122: /*
123: Set function evaluation routine and vector. Whenever the nonlinear
124: solver needs to compute the nonlinear function, it will call this
125: routine.
126: - Note that the final routine argument is the user-defined
127: context that provides application-specific data for the
128: function evaluation routine.
129: */
130: SNESSetFunction(snes,r,FormFunction,&ctx);
132: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: Create matrix data structure; set Jacobian evaluation routine
134: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&J);
137: MatSetFromOptions(J);
139: /*
140: Set Jacobian matrix data structure and default Jacobian evaluation
141: routine. Whenever the nonlinear solver needs to compute the
142: Jacobian matrix, it will call this routine.
143: - Note that the final routine argument is the user-defined
144: context that provides application-specific data for the
145: Jacobian evaluation routine.
146: */
147: SNESSetJacobian(snes,J,J,FormJacobian,&ctx);
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Customize nonlinear solver; set runtime options
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: /*
154: Set an optional user-defined monitoring routine
155: */
156: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer);
157: SNESSetMonitor(snes,Monitor,&monP,0);
159: /*
160: Set names for some vectors to facilitate monitoring (optional)
161: */
162: PetscObjectSetName((PetscObject)x,"Approximate Solution");
163: PetscObjectSetName((PetscObject)U,"Exact Solution");
165: /*
166: Set SNES/KSP/KSP/PC runtime options, e.g.,
167: -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
168: */
169: SNESSetFromOptions(snes);
171: /*
172: Set an optional user-defined routine to check the validity of candidate
173: iterates that are determined by line search methods
174: */
175: PetscOptionsHasName(PETSC_NULL,"-check_iterates",&step_check);
176: if (step_check) {
177: PetscPrintf(PETSC_COMM_WORLD,"Activating step checking routine\n");
178: SNESSetLineSearchCheck(snes,StepCheck,&checkP);
179: VecDuplicate(x,&(checkP.last_step));
180: checkP.tolerance = 1.0;
181: PetscOptionsGetReal(PETSC_NULL,"-check_tol",&checkP.tolerance,PETSC_NULL);
182: }
185: /*
186: Print parameters used for convergence testing (optional) ... just
187: to demonstrate this routine; this information is also printed with
188: the option -snes_view
189: */
190: SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);
191: PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",abstol,rtol,stol,maxit,maxf);
193: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194: Initialize application:
195: Store right-hand-side of PDE and exact solution
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198: /*
199: Get local grid boundaries (for 1-dimensional DA):
200: xs, xm - starting grid index, width of local grid (no ghost points)
201: */
202: DAGetCorners(ctx.da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
204: /*
205: Get pointers to vector data
206: */
207: DAVecGetArray(ctx.da,F,&FF);
208: DAVecGetArray(ctx.da,U,&UU);
210: /*
211: Compute local vector entries
212: */
213: xp = ctx.h*xs;
214: for (i=xs; i<xs+xm; i++) {
215: FF[i] = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
216: UU[i] = xp*xp*xp;
217: xp += ctx.h;
218: }
220: /*
221: Restore vectors
222: */
223: DAVecRestoreArray(ctx.da,F,&FF);
224: DAVecRestoreArray(ctx.da,U,&UU);
226: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227: Evaluate initial guess; then solve nonlinear system
228: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230: /*
231: Note: The user should initialize the vector, x, with the initial guess
232: for the nonlinear solver prior to calling SNESSolve(). In particular,
233: to employ an initial guess of zero, the user should explicitly set
234: this vector to zero by calling VecSet().
235: */
236: FormInitialGuess(x);
237: SNESSolve(snes,x);
238: SNESGetIterationNumber(snes,&its);
239: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n\n",its);
241: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242: Check solution and clean up
243: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
245: /*
246: Check the error
247: */
248: VecAXPY(&none,U,x);
249: VecNorm(x,NORM_2,&norm);
250: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %D\n",norm,its);
252: /*
253: Free work space. All PETSc objects should be destroyed when they
254: are no longer needed.
255: */
256: PetscViewerDestroy(monP.viewer);
257: if (step_check) {VecDestroy(checkP.last_step);}
258: VecDestroy(x);
259: VecDestroy(r);
260: VecDestroy(U);
261: VecDestroy(F);
262: MatDestroy(J);
263: SNESDestroy(snes);
264: DADestroy(ctx.da);
265: PetscFinalize();
267: return(0);
268: }
269: /* ------------------------------------------------------------------- */
272: /*
273: FormInitialGuess - Computes initial guess.
275: Input/Output Parameter:
276: . x - the solution vector
277: */
278: PetscErrorCode FormInitialGuess(Vec x)
279: {
281: PetscScalar pfive = .50;
284: VecSet(&pfive,x);
285: return(0);
286: }
287: /* ------------------------------------------------------------------- */
290: /*
291: FormFunction - Evaluates nonlinear function, F(x).
293: Input Parameters:
294: . snes - the SNES context
295: . x - input vector
296: . ctx - optional user-defined context, as set by SNESSetFunction()
298: Output Parameter:
299: . f - function vector
301: Note:
302: The user-defined context can contain any application-specific
303: data needed for the function evaluation.
304: */
305: PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
306: {
307: ApplicationCtx *user = (ApplicationCtx*) ctx;
308: DA da = user->da;
309: PetscScalar *xx,*ff,*FF,d;
311: PetscInt i,M,xs,xm;
312: Vec xlocal;
315: DAGetLocalVector(da,&xlocal);
316: /*
317: Scatter ghost points to local vector, using the 2-step process
318: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
319: By placing code between these two statements, computations can
320: be done while messages are in transition.
321: */
322: DAGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);
323: DAGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);
325: /*
326: Get pointers to vector data.
327: - The vector xlocal includes ghost point; the vectors x and f do
328: NOT include ghost points.
329: - Using DAVecGetArray() allows accessing the values using global ordering
330: */
331: DAVecGetArray(da,xlocal,&xx);
332: DAVecGetArray(da,f,&ff);
333: DAVecGetArray(da,user->F,&FF);
335: /*
336: Get local grid boundaries (for 1-dimensional DA):
337: xs, xm - starting grid index, width of local grid (no ghost points)
338: */
339: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
340: DAGetInfo(da,PETSC_NULL,&M,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
341: PETSC_NULL,PETSC_NULL,PETSC_NULL);
343: /*
344: Set function values for boundary points; define local interior grid point range:
345: xsi - starting interior grid index
346: xei - ending interior grid index
347: */
348: if (xs == 0) { /* left boundary */
349: ff[0] = xx[0];
350: xs++;xm--;
351: }
352: if (xs+xm == M) { /* right boundary */
353: ff[xs+xm-1] = xx[xs+xm-1] - 1.0;
354: xm--;
355: }
357: /*
358: Compute function over locally owned part of the grid (interior points only)
359: */
360: d = 1.0/(user->h*user->h);
361: for (i=xs; i<xs+xm; i++) {
362: ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
363: }
365: /*
366: Restore vectors
367: */
368: DAVecRestoreArray(da,xlocal,&xx);
369: DAVecRestoreArray(da,f,&ff);
370: DAVecRestoreArray(da,user->F,&FF);
371: DARestoreLocalVector(da,&xlocal);
372: return(0);
373: }
374: /* ------------------------------------------------------------------- */
377: /*
378: FormJacobian - Evaluates Jacobian matrix.
380: Input Parameters:
381: . snes - the SNES context
382: . x - input vector
383: . dummy - optional user-defined context (not used here)
385: Output Parameters:
386: . jac - Jacobian matrix
387: . B - optionally different preconditioning matrix
388: . flag - flag indicating matrix structure
389: */
390: PetscErrorCode FormJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure*flag,void *ctx)
391: {
392: ApplicationCtx *user = (ApplicationCtx*) ctx;
393: PetscScalar *xx,d,A[3];
395: PetscInt i,j[3],M,xs,xm;
396: DA da = user->da;
399: /*
400: Get pointer to vector data
401: */
402: DAVecGetArray(da,x,&xx);
403: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
405: /*
406: Get range of locally owned matrix
407: */
408: DAGetInfo(da,PETSC_NULL,&M,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
409: PETSC_NULL,PETSC_NULL,PETSC_NULL);
411: /*
412: Determine starting and ending local indices for interior grid points.
413: Set Jacobian entries for boundary points.
414: */
416: if (xs == 0) { /* left boundary */
417: i = 0; A[0] = 1.0;
418: MatSetValues(*jac,1,&i,1,&i,A,INSERT_VALUES);
419: xs++;xm--;
420: }
421: if (xs+xm == M) { /* right boundary */
422: i = M-1; A[0] = 1.0;
423: MatSetValues(*jac,1,&i,1,&i,A,INSERT_VALUES);
424: xm--;
425: }
427: /*
428: Interior grid points
429: - Note that in this case we set all elements for a particular
430: row at once.
431: */
432: d = 1.0/(user->h*user->h);
433: for (i=xs; i<xs+xm; i++) {
434: j[0] = i - 1; j[1] = i; j[2] = i + 1;
435: A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
436: MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);
437: }
439: /*
440: Assemble matrix, using the 2-step process:
441: MatAssemblyBegin(), MatAssemblyEnd().
442: By placing code between these two statements, computations can be
443: done while messages are in transition.
445: Also, restore vector.
446: */
448: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
449: DAVecRestoreArray(da,x,&xx);
450: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
452: *flag = SAME_NONZERO_PATTERN;
453: return(0);
454: }
455: /* ------------------------------------------------------------------- */
458: /*
459: Monitor - Optional user-defined monitoring routine that views the
460: current iterate with an x-window plot. Set by SNESSetMonitor().
462: Input Parameters:
463: snes - the SNES context
464: its - iteration number
465: norm - 2-norm function value (may be estimated)
466: ctx - optional user-defined context for private data for the
467: monitor routine, as set by SNESSetMonitor()
469: Note:
470: See the manpage for PetscViewerDrawOpen() for useful runtime options,
471: such as -nox to deactivate all x-window output.
472: */
473: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
474: {
476: MonitorCtx *monP = (MonitorCtx*) ctx;
477: Vec x;
480: PetscPrintf(PETSC_COMM_WORLD,"iter = %D,SNES Function norm %g\n",its,fnorm);
481: SNESGetSolution(snes,&x);
482: VecView(x,monP->viewer);
483: return(0);
484: }
485: /* ------------------------------------------------------------------- */
488: /*
489: StepCheck - Optional user-defined routine that checks the validity of
490: candidate steps of a line search method. Set by SNESSetLineSearchCheck().
492: Input Parameters:
493: snes - the SNES context
494: ctx - optional user-defined context for private data for the
495: monitor routine, as set by SNESSetLineSearchCheck()
496: x - the new candidate iterate
498: Output Parameters:
499: x - current iterate (possibly modified)
500: flg - flag indicating whether x has been modified (either
501: PETSC_TRUE of PETSC_FALSE)
502: */
503: PetscErrorCode StepCheck(SNES snes,void *ctx,Vec x,PetscTruth *flg)
504: {
506: PetscInt i,iter,xs,xm;
507: ApplicationCtx *user;
508: StepCheckCtx *check = (StepCheckCtx*) ctx;
509: PetscScalar *xa,*xa_last,tmp;
510: PetscReal rdiff;
511: DA da;
514: *flg = PETSC_FALSE;
515: SNESGetIterationNumber(snes,&iter);
517: if (iter > 1) {
518: SNESGetApplicationContext(snes,(void**)&user);
519: da = user->da;
520: PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %D with tolerance %g\n",
521: iter,check->tolerance);
523: /* Access local array data */
524: DAVecGetArray(da,check->last_step,&xa_last);
525: DAVecGetArray(da,x,&xa);
526: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
528: /*
529: If we fail the user-defined check for validity of the candidate iterate,
530: then modify the iterate as we like. (Note that the particular modification
531: below is intended simply to demonstrate how to manipulate this data, not
532: as a meaningful or appropriate choice.)
533: */
534: for (i=xs; i<xs+xm; i++) {
535: rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]);
536: if (rdiff > check->tolerance) {
537: tmp = xa[i];
538: xa[i] = (xa[i] + xa_last[i])/2.0;
539: *flg = PETSC_TRUE;
540: PetscPrintf(PETSC_COMM_WORLD," Altering entry %D: x=%g, x_last=%g, diff=%g, x_new=%g\n",
541: i,PetscAbsScalar(tmp),PetscAbsScalar(xa_last[i]),rdiff,PetscAbsScalar(xa[i]));
542: }
543: }
544: DAVecRestoreArray(da,check->last_step,&xa_last);
545: DAVecRestoreArray(da,x,&xa);
546: }
547: VecCopy(x,check->last_step);
549: return(0);
550: }