Actual source code: ex7.c
2: static char help[] = "Solves u`` + u^{2} = f with Newton-like methods. Using\n\
3: matrix-free techniques with user-provided explicit preconditioner matrix.\n\n";
5: #include petscsnes.h
12: typedef struct {
13: PetscViewer viewer;
14: } MonitorCtx;
16: typedef struct {
17: Mat precond;
18: PetscTruth variant;
19: } AppCtx;
23: int main(int argc,char **argv)
24: {
25: SNES snes; /* SNES context */
26: SNESType type = SNESLS; /* default nonlinear solution method */
27: Vec x,r,F,U,work; /* vectors */
28: Mat J,B; /* Jacobian matrix-free, explicit preconditioner */
29: MonitorCtx monP; /* monitoring context */
30: AppCtx user; /* user-defined work context */
31: PetscScalar h,xp = 0.0,v;
32: PetscInt its,n = 5,i;
35: PetscInitialize(&argc,&argv,(char *)0,help);
36: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
37: PetscOptionsHasName(PETSC_NULL,"-variant",&user.variant);
38: h = 1.0/(n-1);
40: /* Set up data structures */
41: PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&monP.viewer);
42: VecCreateSeq(PETSC_COMM_SELF,n,&x);
43: PetscObjectSetName((PetscObject)x,"Approximate Solution");
44: VecDuplicate(x,&r);
45: VecDuplicate(x,&F);
46: VecDuplicate(x,&U);
47: PetscObjectSetName((PetscObject)U,"Exact Solution");
49: /* create explict matrix preconditioner */
50: MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,PETSC_NULL,&B);
51: user.precond = B;
53: /* Store right-hand-side of PDE and exact solution */
54: for (i=0; i<n; i++) {
55: v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
56: VecSetValues(F,1,&i,&v,INSERT_VALUES);
57: v= xp*xp*xp;
58: VecSetValues(U,1,&i,&v,INSERT_VALUES);
59: xp += h;
60: }
62: /* Create nonlinear solver */
63: SNESCreate(PETSC_COMM_WORLD,&snes);
64: SNESSetType(snes,type);
66: if (user.variant) {
67: MatCreateMF(x,&J);
68: VecDuplicate(x,&work);
69: MatSNESMFSetFunction(J,work,FormFunction,F);
70: } else {
71: /* create matrix free matrix for Jacobian */
72: MatCreateSNESMF(snes,x,&J);
73: }
75: /* Set various routines and options */
76: SNESSetFunction(snes,r,FormFunction,F);
77: SNESSetJacobian(snes,J,B,FormJacobian,&user);
78: SNESMonitorSet(snes,Monitor,&monP,0);
79: SNESSetFromOptions(snes);
81: /* Solve nonlinear system */
82: FormInitialGuess(snes,x);
83: SNESSolve(snes,PETSC_NULL,x);
84: SNESGetIterationNumber(snes,&its);
85: PetscPrintf(PETSC_COMM_SELF,"number of Newton iterations = %D\n\n",its);
87: /* Free data structures */
88: if (user.variant) {
89: VecDestroy(work);
90: }
91: VecDestroy(x); VecDestroy(r);
92: VecDestroy(U); VecDestroy(F);
93: MatDestroy(J); MatDestroy(B);
94: SNESDestroy(snes);
95: PetscViewerDestroy(monP.viewer);
96: PetscFinalize();
98: return 0;
99: }
100: /* -------------------- Evaluate Function F(x) --------------------- */
102: PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *dummy)
103: {
104: PetscScalar *xx,*ff,*FF,d;
105: PetscInt i,n;
108: VecGetArray(x,&xx);
109: VecGetArray(f,&ff);
110: VecGetArray((Vec) dummy,&FF);
111: VecGetSize(x,&n);
112: d = (PetscReal)(n - 1); d = d*d;
113: ff[0] = xx[0];
114: for (i=1; i<n-1; i++) {
115: ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
116: }
117: ff[n-1] = xx[n-1] - 1.0;
118: VecRestoreArray(x,&xx);
119: VecRestoreArray(f,&ff);
120: VecRestoreArray((Vec)dummy,&FF);
121: return 0;
122: }
123: /* -------------------- Form initial approximation ----------------- */
127: PetscErrorCode FormInitialGuess(SNES snes,Vec x)
128: {
129: PetscErrorCode ierr;
130: PetscScalar pfive = .50;
131: VecSet(x,pfive);
132: return 0;
133: }
136: /* -------------------- Evaluate Jacobian F'(x) -------------------- */
137: /* Evaluates a matrix that is used to precondition the matrix-free
138: jacobian. In this case, the explict preconditioner matrix is
139: also EXACTLY the Jacobian. In general, it would be some lower
140: order, simplified apprioximation */
142: PetscErrorCode FormJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure*flag,void *dummy)
143: {
144: PetscScalar *xx,A[3],d;
145: PetscInt i,n,j[3],iter;
147: AppCtx *user = (AppCtx*) dummy;
149: SNESGetIterationNumber(snes,&iter);
151: if (iter%2 ==0) { /* Compute new preconditioner matrix */
152: PetscPrintf(PETSC_COMM_SELF,"iter=%D, computing new preconditioning matrix\n",iter+1);
153: *B = user->precond;
154: VecGetArray(x,&xx);
155: VecGetSize(x,&n);
156: d = (PetscReal)(n - 1); d = d*d;
158: i = 0; A[0] = 1.0;
159: MatSetValues(*B,1,&i,1,&i,&A[0],INSERT_VALUES);
160: for (i=1; i<n-1; i++) {
161: j[0] = i - 1; j[1] = i; j[2] = i + 1;
162: A[0] = d; A[1] = -2.0*d + 2.0*xx[i]; A[2] = d;
163: MatSetValues(*B,1,&i,3,j,A,INSERT_VALUES);
164: }
165: i = n-1; A[0] = 1.0;
166: MatSetValues(*B,1,&i,1,&i,&A[0],INSERT_VALUES);
167: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
168: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
169: VecRestoreArray(x,&xx);
170: *flag = SAME_NONZERO_PATTERN;
171: } else { /* reuse preconditioner from last iteration */
172: PetscPrintf(PETSC_COMM_SELF,"iter=%D, using old preconditioning matrix\n",iter+1);
173: *flag = SAME_PRECONDITIONER;
174: }
175: if (user->variant) {
176: MatSNESMFSetBase(*jac,x);
177: }
178: MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
179: MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
180: return 0;
181: }
182: /* -------------------- User-defined monitor ----------------------- */
186: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *dummy)
187: {
189: MonitorCtx *monP = (MonitorCtx*) dummy;
190: Vec x;
191: MPI_Comm comm;
193: PetscObjectGetComm((PetscObject)snes,&comm);
194: PetscFPrintf(comm,stdout,"iter = %D, SNES Function norm %G \n",its,fnorm);
195: SNESGetSolution(snes,&x);
196: VecView(x,monP->viewer);
197: return 0;
198: }