Actual source code: ex38.c
2: static char help[] = "u`` + u^{2} = f. Different matrices for the Jacobian and the preconditioner.\n\
3: Demonstrates the use of matrix-free Newton-Krylov methods in conjunction\n\
4: with a user-provided preconditioner. Input arguments are:\n\
5: -snes_mf : Use matrix-free Newton methods\n\
6: -user_precond : Employ a user-defined preconditioner. Used only with\n\
7: matrix-free methods in this example.\n\n";
8: /*
9: Modified from ex6.c by Mike McCourt <mccomic@iit.edu>
10: for testing SNESLineSearchSet()
11: */
13: /*T
14: Concepts: SNES^different matrices for the Jacobian and preconditioner;
15: Concepts: SNES^matrix-free methods
16: Concepts: SNES^user-provided preconditioner;
17: Concepts: matrix-free methods
18: Concepts: user-provided preconditioner;
19: Processors: 1
20: T*/
22: /*
23: Include "petscsnes.h" so that we can use SNES solvers. Note that this
24: file automatically includes:
25: petsc.h - base PETSc routines petscvec.h - vectors
26: petscsys.h - system routines petscmat.h - matrices
27: petscis.h - index sets petscksp.h - Krylov subspace methods
28: petscviewer.h - viewers petscpc.h - preconditioners
29: petscksp.h - linear solvers
30: */
31: #include <iostream>
32: using namespace std;
33: #include petscsnes.h
35: struct AppCtx{int testint;};
37: /*
38: User-defined routines
39: */
40: PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
41: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
42: PetscErrorCode MatrixFreePreconditioner(void*,Vec,Vec);
43: PetscErrorCode FormLineSearch(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*);
45: int main(int argc,char **argv)
46: {
47: SNES snes; /* SNES context */
48: KSP ksp; /* KSP context */
49: PC pc; /* PC context */
50: Vec x,r,F; /* vectors */
51: Mat J,JPrec; /* Jacobian,preconditioner matrices */
53: PetscInt it,n = 5,i;
54: PetscMPIInt size;
55: PetscInt *Shistit = 0,Khistl = 200,Shistl = 10;
56: PetscReal h,xp = 0.0,*Khist = 0,*Shist = 0;
57: PetscScalar v,pfive = .5;
58: PetscTruth flg;
59: AppCtx user;
61: PetscInitialize(&argc,&argv,(char *)0,help);
62: MPI_Comm_size(PETSC_COMM_WORLD,&size);
63: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
64: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
65: h = 1.0/(n-1);
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create nonlinear solver context
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: SNESCreate(PETSC_COMM_WORLD,&snes);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create vector data structures; set function evaluation routine
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: VecCreate(PETSC_COMM_SELF,&x);
78: VecSetSizes(x,PETSC_DECIDE,n);
79: VecSetFromOptions(x);
80: VecDuplicate(x,&r);
81: VecDuplicate(x,&F);
83: SNESSetFunction(snes,r,FormFunction,(void*)F);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Create matrix data structures; set Jacobian evaluation routine
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,PETSC_NULL,&J);
90: MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,1,PETSC_NULL,&JPrec);
92: /*
93: Note that in this case we create separate matrices for the Jacobian
94: and preconditioner matrix. Both of these are computed in the
95: routine FormJacobian()
96: */
97: SNESSetJacobian(snes,J,JPrec,FormJacobian,0);
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Customize nonlinear solver; set runtime options
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: /* Set preconditioner for matrix-free method */
104: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&flg);
105: if (flg) {
106: SNESGetKSP(snes,&ksp);
107: KSPGetPC(ksp,&pc);
108: PetscOptionsHasName(PETSC_NULL,"-user_precond",&flg);
109: if (flg) { /* user-defined precond */
110: PCSetType(pc,PCSHELL);
111: PCShellSetApply(pc,MatrixFreePreconditioner);
112: } else {PCSetType(pc,PCNONE);}
113: }
115: SNESSetFromOptions(snes);
116: user.testint = 0;
117: SNESLineSearchSet(snes,FormLineSearch,&user);
119: /*
120: Save all the linear residuals for all the Newton steps; this enables us
121: to retain complete convergence history for printing after the conclusion
122: of SNESSolve(). Alternatively, one could use the monitoring options
123: -snes_monitor -ksp_monitor
124: to see this information during the solver's execution; however, such
125: output during the run distorts performance evaluation data. So, the
126: following is a good option when monitoring code performance, for example
127: when using -log_summary.
128: */
129: PetscOptionsHasName(PETSC_NULL,"-rhistory",&flg);
130: if (flg) {
131: SNESGetKSP(snes,&ksp);
132: PetscMalloc(Khistl*sizeof(PetscReal),&Khist);
133: KSPSetResidualHistory(ksp,Khist,Khistl,PETSC_FALSE);
134: PetscMalloc(Shistl*sizeof(PetscReal),&Shist);
135: PetscMalloc(Shistl*sizeof(PetscInt),&Shistit);
136: SNESSetConvergenceHistory(snes,Shist,Shistit,Shistl,PETSC_FALSE);
137: }
139: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140: Initialize application:
141: Store right-hand-side of PDE and exact solution
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: xp = 0.0;
145: for (i=0; i<n; i++) {
146: v = 6.0*xp + pow(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
147: VecSetValues(F,1,&i,&v,INSERT_VALUES);
148: xp += h;
149: }
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Evaluate initial guess; then solve nonlinear system
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: VecSet(x,pfive);
156: SNESSolve(snes,PETSC_NULL,x);
157: SNESGetIterationNumber(snes,&it);
158: PetscPrintf(PETSC_COMM_SELF,"Newton iterations = %D\n\n",it);
160: PetscOptionsHasName(PETSC_NULL,"-rhistory",&flg);
161: if (flg) {
162: KSPGetResidualHistory(ksp,PETSC_NULL,&Khistl);
163: PetscRealView(Khistl,Khist,PETSC_VIEWER_STDOUT_SELF);
164: PetscFree(Khist);
165: SNESGetConvergenceHistory(snes,PETSC_NULL,PETSC_NULL,&Shistl);
166: PetscRealView(Shistl,Shist,PETSC_VIEWER_STDOUT_SELF);
167: PetscIntView(Shistl,Shistit,PETSC_VIEWER_STDOUT_SELF);
168: PetscFree(Shist);
169: PetscFree(Shistit);
170: }
172: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173: Free work space. All PETSc objects should be destroyed when they
174: are no longer needed.
175: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177: VecDestroy(x); VecDestroy(r);
178: VecDestroy(F); MatDestroy(J);
179: MatDestroy(JPrec); SNESDestroy(snes);
180: PetscFinalize();
182: return 0;
183: }
185: PetscErrorCode FormLineSearch(SNES snes,void* user,Vec X,Vec F,Vec G,Vec Y,Vec W,PetscReal fnorm,
186: PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
187: {
189: PetscScalar mone=-1.0;
190: AppCtx *myguy = (AppCtx*)user;
191: *flag=PETSC_TRUE;
193: cout << "Inside FormLineSearch \n user.testint=" << myguy->testint << endl;
194: ierr=VecNorm(Y,NORM_2,ynorm);
195: ierr=VecWAXPY(W,mone,Y,X); /* W = -Y + X */
196: ierr=SNESComputeFunction(snes,W,G);
197: ierr=VecNorm(G,NORM_2,gnorm);
198: myguy->testint++;
199: return ierr;
200: }
202: /* ------------------------------------------------------------------- */
203: /*
204: FormInitialGuess - Forms initial approximation.
206: Input Parameters:
207: user - user-defined application context
208: X - vector
210: Output Parameter:
211: X - vector
212: */
213: PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *dummy)
214: {
215: PetscScalar *xx,*ff,*FF,d;
217: PetscInt i,n;
219: VecGetArray(x,&xx);
220: VecGetArray(f,&ff);
221: VecGetArray((Vec)dummy,&FF);
222: VecGetSize(x,&n);
223: d = (PetscReal)(n - 1); d = d*d;
224: ff[0] = xx[0];
225: for (i=1; i<n-1; i++) {
226: ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
227: }
228: ff[n-1] = xx[n-1] - 1.0;
229: VecRestoreArray(x,&xx);
230: VecRestoreArray(f,&ff);
231: VecRestoreArray((Vec)dummy,&FF);
232: return 0;
233: }
234: /* ------------------------------------------------------------------- */
235: /*
236: FormJacobian - This routine demonstrates the use of different
237: matrices for the Jacobian and preconditioner
239: Input Parameters:
240: . snes - the SNES context
241: . x - input vector
242: . ptr - optional user-defined context, as set by SNESSetJacobian()
244: Output Parameters:
245: . A - Jacobian matrix
246: . B - different preconditioning matrix
247: . flag - flag indicating matrix structure
248: */
249: PetscErrorCode FormJacobian(SNES snes,Vec x,Mat *jac,Mat *prejac,MatStructure *flag,void *dummy)
250: {
251: PetscScalar *xx,A[3],d;
252: PetscInt i,n,j[3];
255: VecGetArray(x,&xx);
256: VecGetSize(x,&n);
257: d = (PetscReal)(n - 1); d = d*d;
259: /* Form Jacobian. Also form a different preconditioning matrix that
260: has only the diagonal elements. */
261: i = 0; A[0] = 1.0;
262: MatSetValues(*jac,1,&i,1,&i,&A[0],INSERT_VALUES);
263: MatSetValues(*prejac,1,&i,1,&i,&A[0],INSERT_VALUES);
264: for (i=1; i<n-1; i++) {
265: j[0] = i - 1; j[1] = i; j[2] = i + 1;
266: A[0] = d; A[1] = -2.0*d + 2.0*xx[i]; A[2] = d;
267: MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);
268: MatSetValues(*prejac,1,&i,1,&i,&A[1],INSERT_VALUES);
269: }
270: i = n-1; A[0] = 1.0;
271: MatSetValues(*jac,1,&i,1,&i,&A[0],INSERT_VALUES);
272: MatSetValues(*prejac,1,&i,1,&i,&A[0],INSERT_VALUES);
274: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
275: MatAssemblyBegin(*prejac,MAT_FINAL_ASSEMBLY);
276: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
277: MatAssemblyEnd(*prejac,MAT_FINAL_ASSEMBLY);
279: VecRestoreArray(x,&xx);
280: *flag = SAME_NONZERO_PATTERN;
281: return 0;
282: }
283: /* ------------------------------------------------------------------- */
284: /*
285: MatrixFreePreconditioner - This routine demonstrates the use of a
286: user-provided preconditioner. This code implements just the null
287: preconditioner, which of course is not recommended for general use.
289: Input Parameters:
290: . ctx - optional user-defined context, as set by PCShellSetContext()
291: . x - input vector
293: Output Parameter:
294: . y - preconditioned vector
295: */
296: PetscErrorCode MatrixFreePreconditioner(void *ctx,Vec x,Vec y)
297: {
299: VecCopy(x,y);
300: return 0;
301: }