Actual source code: ex2.c
2: static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
3: This example employs a user-defined monitoring routine.\n\n";
5: /*T
6: Concepts: SNES^basic uniprocessor example
7: Concepts: SNES^setting a user-defined monitoring routine
8: Processors: 1
9: T*/
11: /*
12: Include "petscdraw.h" so that we can use PETSc drawing routines.
13: Include "petscsnes.h" so that we can use SNES solvers. Note that this
14: file automatically includes:
15: petsc.h - base PETSc routines petscvec.h - vectors
16: petscsys.h - system routines petscmat.h - matrices
17: petscis.h - index sets petscksp.h - Krylov subspace methods
18: petscviewer.h - viewers petscpc.h - preconditioners
19: petscksp.h - linear solvers
20: */
22: #include petscsnes.h
24: /*
25: User-defined routines
26: */
32: /*
33: User-defined context for monitoring
34: */
35: typedef struct {
36: PetscViewer viewer;
37: } MonitorCtx;
41: int main(int argc,char **argv)
42: {
43: SNES snes; /* SNES context */
44: Vec x,r,F,U; /* vectors */
45: Mat J; /* Jacobian matrix */
46: MonitorCtx monP; /* monitoring context */
48: PetscInt its,n = 5,i,maxit,maxf;
49: PetscMPIInt size;
50: PetscScalar h,xp,v,none = -1.0;
51: PetscReal abstol,rtol,stol,norm;
53: PetscInitialize(&argc,&argv,(char *)0,help);
54: MPI_Comm_size(PETSC_COMM_WORLD,&size);
55: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
56: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
57: h = 1.0/(n-1);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create nonlinear solver context
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: SNESCreate(PETSC_COMM_WORLD,&snes);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Create vector data structures; set function evaluation routine
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: /*
69: Note that we form 1 vector from scratch and then duplicate as needed.
70: */
71: VecCreate(PETSC_COMM_WORLD,&x);
72: VecSetSizes(x,PETSC_DECIDE,n);
73: VecSetFromOptions(x);
74: VecDuplicate(x,&r);
75: VecDuplicate(x,&F);
76: VecDuplicate(x,&U);
78: /*
79: Set function evaluation routine and vector
80: */
81: SNESSetFunction(snes,r,FormFunction,(void*)F);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create matrix data structure; set Jacobian evaluation routine
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&J);
89: MatSetFromOptions(J);
91: /*
92: Set Jacobian matrix data structure and default Jacobian evaluation
93: routine. User can override with:
94: -snes_fd : default finite differencing approximation of Jacobian
95: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
96: (unless user explicitly sets preconditioner)
97: -snes_mf_operator : form preconditioning matrix as set by the user,
98: but use matrix-free approx for Jacobian-vector
99: products within Newton-Krylov method
100: */
102: SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Customize nonlinear solver; set runtime options
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: /*
109: Set an optional user-defined monitoring routine
110: */
111: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer);
112: SNESSetMonitor(snes,Monitor,&monP,0);
114: /*
115: Set names for some vectors to facilitate monitoring (optional)
116: */
117: PetscObjectSetName((PetscObject)x,"Approximate Solution");
118: PetscObjectSetName((PetscObject)U,"Exact Solution");
120: /*
121: Set SNES/KSP/KSP/PC runtime options, e.g.,
122: -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
123: */
124: SNESSetFromOptions(snes);
126: /*
127: Print parameters used for convergence testing (optional) ... just
128: to demonstrate this routine; this information is also printed with
129: the option -snes_view
130: */
131: SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);
132: PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",abstol,rtol,stol,maxit,maxf);
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Initialize application:
136: Store right-hand-side of PDE and exact solution
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: xp = 0.0;
140: for (i=0; i<n; i++) {
141: v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
142: VecSetValues(F,1,&i,&v,INSERT_VALUES);
143: v= xp*xp*xp;
144: VecSetValues(U,1,&i,&v,INSERT_VALUES);
145: xp += h;
146: }
148: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: Evaluate initial guess; then solve nonlinear system
150: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: /*
152: Note: The user should initialize the vector, x, with the initial guess
153: for the nonlinear solver prior to calling SNESSolve(). In particular,
154: to employ an initial guess of zero, the user should explicitly set
155: this vector to zero by calling VecSet().
156: */
157: FormInitialGuess(x);
158: SNESSolve(snes,x);
159: SNESGetIterationNumber(snes,&its);
160: PetscPrintf(PETSC_COMM_WORLD,"number of Newton iterations = %D\n\n",its);
162: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163: Check solution and clean up
164: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166: /*
167: Check the error
168: */
169: VecAXPY(&none,U,x);
170: VecNorm(x,NORM_2,&norm);
171: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %D\n",norm,its);
174: /*
175: Free work space. All PETSc objects should be destroyed when they
176: are no longer needed.
177: */
178: VecDestroy(x); VecDestroy(r);
179: VecDestroy(U); VecDestroy(F);
180: MatDestroy(J); SNESDestroy(snes);
181: PetscViewerDestroy(monP.viewer);
182: PetscFinalize();
184: return 0;
185: }
186: /* ------------------------------------------------------------------- */
189: /*
190: FormInitialGuess - Computes initial guess.
192: Input/Output Parameter:
193: . x - the solution vector
194: */
195: PetscErrorCode FormInitialGuess(Vec x)
196: {
198: PetscScalar pfive = .50;
199: VecSet(&pfive,x);
200: return 0;
201: }
202: /* ------------------------------------------------------------------- */
205: /*
206: FormFunction - Evaluates nonlinear function, F(x).
208: Input Parameters:
209: . snes - the SNES context
210: . x - input vector
211: . ctx - optional user-defined context, as set by SNESSetFunction()
213: Output Parameter:
214: . f - function vector
216: Note:
217: The user-defined context can contain any application-specific data
218: needed for the function evaluation (such as various parameters, work
219: vectors, and grid information). In this program the context is just
220: a vector containing the right-hand-side of the discretized PDE.
221: */
223: PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
224: {
225: Vec g = (Vec)ctx;
226: PetscScalar *xx,*ff,*gg,d;
228: PetscInt i,n;
230: /*
231: Get pointers to vector data.
232: - For default PETSc vectors, VecGetArray() returns a pointer to
233: the data array. Otherwise, the routine is implementation dependent.
234: - You MUST call VecRestoreArray() when you no longer need access to
235: the array.
236: */
237: VecGetArray(x,&xx);
238: VecGetArray(f,&ff);
239: VecGetArray(g,&gg);
241: /*
242: Compute function
243: */
244: VecGetSize(x,&n);
245: d = (PetscReal)(n - 1); d = d*d;
246: ff[0] = xx[0];
247: for (i=1; i<n-1; i++) {
248: ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - gg[i];
249: }
250: ff[n-1] = xx[n-1] - 1.0;
252: /*
253: Restore vectors
254: */
255: VecRestoreArray(x,&xx);
256: VecRestoreArray(f,&ff);
257: VecRestoreArray(g,&gg);
258: return 0;
259: }
260: /* ------------------------------------------------------------------- */
263: /*
264: FormJacobian - Evaluates Jacobian matrix.
266: Input Parameters:
267: . snes - the SNES context
268: . x - input vector
269: . dummy - optional user-defined context (not used here)
271: Output Parameters:
272: . jac - Jacobian matrix
273: . B - optionally different preconditioning matrix
274: . flag - flag indicating matrix structure
275: */
277: PetscErrorCode FormJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure*flag,void *dummy)
278: {
279: PetscScalar *xx,A[3],d;
281: PetscInt i,n,j[3];
283: /*
284: Get pointer to vector data
285: */
286: VecGetArray(x,&xx);
288: /*
289: Compute Jacobian entries and insert into matrix.
290: - Note that in this case we set all elements for a particular
291: row at once.
292: */
293: VecGetSize(x,&n);
294: d = (PetscReal)(n - 1); d = d*d;
296: /*
297: Interior grid points
298: */
299: for (i=1; i<n-1; i++) {
300: j[0] = i - 1; j[1] = i; j[2] = i + 1;
301: A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
302: MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);
303: }
305: /*
306: Boundary points
307: */
308: i = 0; A[0] = 1.0;
309: MatSetValues(*jac,1,&i,1,&i,A,INSERT_VALUES);
310: i = n-1; A[0] = 1.0;
311: MatSetValues(*jac,1,&i,1,&i,A,INSERT_VALUES);
313: /*
314: Restore vector
315: */
316: VecRestoreArray(x,&xx);
318: /*
319: Assemble matrix
320: */
321: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
322: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
324: return 0;
325: }
326: /* ------------------------------------------------------------------- */
329: /*
330: Monitor - User-defined monitoring routine that views the
331: current iterate with an x-window plot.
333: Input Parameters:
334: snes - the SNES context
335: its - iteration number
336: norm - 2-norm function value (may be estimated)
337: ctx - optional user-defined context for private data for the
338: monitor routine, as set by SNESSetMonitor()
340: Note:
341: See the manpage for PetscViewerDrawOpen() for useful runtime options,
342: such as -nox to deactivate all x-window output.
343: */
344: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
345: {
347: MonitorCtx *monP = (MonitorCtx*) ctx;
348: Vec x;
350: PetscPrintf(PETSC_COMM_WORLD,"iter = %D, SNES Function norm %g\n",its,fnorm);
351: SNESGetSolution(snes,&x);
352: VecView(x,monP->viewer);
353: return 0;
354: }