Actual source code: ex6.c
1: /*$Id: ex6.c,v 1.71 2001/08/07 03:04:16 balay Exp $*/
3: static char help[] = "u`` + u^{2} = f. Different matrices for the Jacobian and the preconditioner.n
4: Demonstrates the use of matrix-free Newton-Krylov methods in conjunctionn
5: with a user-provided preconditioner. Input arguments are:n
6: -snes_mf : Use matrix-free Newton methodsn
7: -user_precond : Employ a user-defined preconditioner. Used only withn
8: matrix-free methods in this example.nn";
10: /*T
11: Concepts: SNES^different matrices for the Jacobian and preconditioner;
12: Concepts: SNES^matrix-free methods
13: Concepts: SNES^user-provided preconditioner;
14: Concepts: matrix-free methods
15: Concepts: user-provided preconditioner;
16: Processors: 1
17: T*/
19: /*
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: petscsles.h - linear solvers
27: */
28: #include petscsnes.h
30: /*
31: User-defined routines
32: */
33: int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
34: int FormFunction(SNES,Vec,Vec,void*);
35: int MatrixFreePreconditioner(void*,Vec,Vec);
37: int main(int argc,char **argv)
38: {
39: SNES snes; /* SNES context */
40: SLES sles; /* SLES context */
41: PC pc; /* PC context */
42: KSP ksp; /* KSP context */
43: Vec x,r,F; /* vectors */
44: Mat J,JPrec; /* Jacobian,preconditioner matrices */
45: int ierr,it,n = 5,i,size;
46: int *Shistit = 0,Khistl = 200,Shistl = 10;
47: PetscReal h,xp = 0.0,*Khist = 0,*Shist = 0;
48: PetscScalar v,pfive = .5;
49: PetscTruth flg;
51: PetscInitialize(&argc,&argv,(char *)0,help);
52: MPI_Comm_size(PETSC_COMM_WORLD,&size);
53: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
54: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
55: h = 1.0/(n-1);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create nonlinear solver context
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: SNESCreate(PETSC_COMM_WORLD,&snes);
63: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: Create vector data structures; set function evaluation routine
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: VecCreate(PETSC_COMM_SELF,&x);
68: VecSetSizes(x,PETSC_DECIDE,n);
69: VecSetFromOptions(x);
70: VecDuplicate(x,&r);
71: VecDuplicate(x,&F);
73: SNESSetFunction(snes,r,FormFunction,(void*)F);
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Create matrix data structures; set Jacobian evaluation routine
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79: MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,PETSC_NULL,&J);
80: MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,1,PETSC_NULL,&JPrec);
82: /*
83: Note that in this case we create separate matrices for the Jacobian
84: and preconditioner matrix. Both of these are computed in the
85: routine FormJacobian()
86: */
87: SNESSetJacobian(snes,J,JPrec,FormJacobian,0);
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Customize nonlinear solver; set runtime options
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: /* Set preconditioner for matrix-free method */
94: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&flg);
95: if (flg) {
96: SNESGetSLES(snes,&sles);
97: SLESGetPC(sles,&pc);
98: PetscOptionsHasName(PETSC_NULL,"-user_precond",&flg);
99: if (flg) { /* user-defined precond */
100: PCSetType(pc,PCSHELL);
101: PCShellSetApply(pc,MatrixFreePreconditioner,PETSC_NULL);
102: } else {PCSetType(pc,PCNONE);}
103: }
105: SNESSetFromOptions(snes);
107: /*
108: Save all the linear residuals for all the Newton steps; this enables us
109: to retain complete convergence history for printing after the conclusion
110: of SNESSolve(). Alternatively, one could use the monitoring options
111: -snes_monitor -ksp_monitor
112: to see this information during the solver's execution; however, such
113: output during the run distorts performance evaluation data. So, the
114: following is a good option when monitoring code performance, for example
115: when using -log_summary.
116: */
117: PetscOptionsHasName(PETSC_NULL,"-rhistory",&flg);
118: if (flg) {
119: SNESGetSLES(snes,&sles);
120: SLESGetKSP(sles,&ksp);
121: PetscMalloc(Khistl*sizeof(PetscReal),&Khist);
122: KSPSetResidualHistory(ksp,Khist,Khistl,PETSC_FALSE);
123: PetscMalloc(Shistl*sizeof(PetscReal),&Shist);
124: PetscMalloc(Shistl*sizeof(int),&Shistit);
125: SNESSetConvergenceHistory(snes,Shist,Shistit,Shistl,PETSC_FALSE);
126: }
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Initialize application:
130: Store right-hand-side of PDE and exact solution
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: xp = 0.0;
134: for (i=0; i<n; i++) {
135: v = 6.0*xp + pow(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
136: VecSetValues(F,1,&i,&v,INSERT_VALUES);
137: xp += h;
138: }
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Evaluate initial guess; then solve nonlinear system
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: VecSet(&pfive,x);
145: SNESSolve(snes,x,&it);
146: PetscPrintf(PETSC_COMM_SELF,"Newton iterations = %dnn",it);
148: PetscOptionsHasName(PETSC_NULL,"-rhistory",&flg);
149: if (flg) {
150: KSPGetResidualHistory(ksp,PETSC_NULL,&Khistl);
151: PetscRealView(Khistl,Khist,PETSC_VIEWER_STDOUT_SELF);
152: PetscFree(Khist);
153: SNESGetConvergenceHistory(snes,PETSC_NULL,PETSC_NULL,&Shistl);
154: PetscRealView(Shistl,Shist,PETSC_VIEWER_STDOUT_SELF);
155: PetscIntView(Shistl,Shistit,PETSC_VIEWER_STDOUT_SELF);
156: PetscFree(Shist);
157: PetscFree(Shistit);
158: }
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Free work space. All PETSc objects should be destroyed when they
162: are no longer needed.
163: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165: VecDestroy(x); VecDestroy(r);
166: VecDestroy(F); MatDestroy(J);
167: MatDestroy(JPrec); SNESDestroy(snes);
168: PetscFinalize();
170: return 0;
171: }
172: /* ------------------------------------------------------------------- */
173: /*
174: FormInitialGuess - Forms initial approximation.
176: Input Parameters:
177: user - user-defined application context
178: X - vector
180: Output Parameter:
181: X - vector
182: */
183: int FormFunction(SNES snes,Vec x,Vec f,void *dummy)
184: {
185: PetscScalar *xx,*ff,*FF,d;
186: int i,ierr,n;
188: VecGetArray(x,&xx);
189: VecGetArray(f,&ff);
190: VecGetArray((Vec)dummy,&FF);
191: VecGetSize(x,&n);
192: d = (PetscReal)(n - 1); d = d*d;
193: ff[0] = xx[0];
194: for (i=1; i<n-1; i++) {
195: ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
196: }
197: ff[n-1] = xx[n-1] - 1.0;
198: VecRestoreArray(x,&xx);
199: VecRestoreArray(f,&ff);
200: VecRestoreArray((Vec)dummy,&FF);
201: return 0;
202: }
203: /* ------------------------------------------------------------------- */
204: /*
205: FormJacobian - This routine demonstrates the use of different
206: matrices for the Jacobian and preconditioner
208: Input Parameters:
209: . snes - the SNES context
210: . x - input vector
211: . ptr - optional user-defined context, as set by SNESSetJacobian()
213: Output Parameters:
214: . A - Jacobian matrix
215: . B - different preconditioning matrix
216: . flag - flag indicating matrix structure
217: */
218: int FormJacobian(SNES snes,Vec x,Mat *jac,Mat *prejac,MatStructure *flag,void *dummy)
219: {
220: PetscScalar *xx,A[3],d;
221: int i,n,j[3],ierr;
223: VecGetArray(x,&xx);
224: VecGetSize(x,&n);
225: d = (PetscReal)(n - 1); d = d*d;
227: /* Form Jacobian. Also form a different preconditioning matrix that
228: has only the diagonal elements. */
229: i = 0; A[0] = 1.0;
230: MatSetValues(*jac,1,&i,1,&i,&A[0],INSERT_VALUES);
231: MatSetValues(*prejac,1,&i,1,&i,&A[0],INSERT_VALUES);
232: for (i=1; i<n-1; i++) {
233: j[0] = i - 1; j[1] = i; j[2] = i + 1;
234: A[0] = d; A[1] = -2.0*d + 2.0*xx[i]; A[2] = d;
235: MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);
236: MatSetValues(*prejac,1,&i,1,&i,&A[1],INSERT_VALUES);
237: }
238: i = n-1; A[0] = 1.0;
239: MatSetValues(*jac,1,&i,1,&i,&A[0],INSERT_VALUES);
240: MatSetValues(*prejac,1,&i,1,&i,&A[0],INSERT_VALUES);
242: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
243: MatAssemblyBegin(*prejac,MAT_FINAL_ASSEMBLY);
244: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
245: MatAssemblyEnd(*prejac,MAT_FINAL_ASSEMBLY);
247: VecRestoreArray(x,&xx);
248: *flag = SAME_NONZERO_PATTERN;
249: return 0;
250: }
251: /* ------------------------------------------------------------------- */
252: /*
253: MatrixFreePreconditioner - This routine demonstrates the use of a
254: user-provided preconditioner. This code implements just the null
255: preconditioner, which of course is not recommended for general use.
257: Input Parameters:
258: . ctx - optional user-defined context, as set by PCShellSetApply()
259: . x - input vector
261: Output Parameter:
262: . y - preconditioned vector
263: */
264: int MatrixFreePreconditioner(void *ctx,Vec x,Vec y)
265: {
267: VecCopy(x,y);
268: return 0;
269: }