Actual source code: ex7.c
petsc-3.7.3 2016-08-01
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>
7: extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
8: extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
9: extern PetscErrorCode OtherFunctionForDifferencing(void*,Vec,Vec);
10: extern PetscErrorCode FormInitialGuess(SNES,Vec);
11: extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
13: typedef struct {
14: PetscViewer viewer;
15: } MonitorCtx;
17: typedef struct {
18: PetscBool variant;
19: } AppCtx;
23: int main(int argc,char **argv)
24: {
25: SNES snes; /* SNES context */
26: SNESType type = SNESNEWTONLS; /* default nonlinear solution method */
27: Vec x,r,F,U; /* 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(NULL,NULL,"-n",&n,NULL);
37: PetscOptionsHasName(NULL,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,NULL,&B);
52: /* Store right-hand-side of PDE and exact solution */
53: for (i=0; i<n; i++) {
54: v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
55: VecSetValues(F,1,&i,&v,INSERT_VALUES);
56: v = xp*xp*xp;
57: VecSetValues(U,1,&i,&v,INSERT_VALUES);
58: xp += h;
59: }
61: /* Create nonlinear solver */
62: SNESCreate(PETSC_COMM_WORLD,&snes);
63: SNESSetType(snes,type);
65: /* Set various routines and options */
66: SNESSetFunction(snes,r,FormFunction,F);
67: if (user.variant) {
68: /* this approach is not normally needed, one should use the MatCreateSNESMF() below usually */
69: MatCreateMFFD(PETSC_COMM_WORLD,n,n,n,n,&J);
70: MatMFFDSetFunction(J,(PetscErrorCode (*)(void*, Vec, Vec))SNESComputeFunction,snes);
71: } else {
72: /* create matrix free matrix for Jacobian */
73: MatCreateSNESMF(snes,&J);
74: /* demonstrates differencing a different function than FormFunction() to apply a matrix operator */
75: /* note we use the same context for this function as FormFunction, the F vector */
76: MatMFFDSetFunction(J,OtherFunctionForDifferencing,F);
77: }
79: /* Set various routines and options */
80: SNESSetJacobian(snes,J,B,FormJacobian,&user);
81: SNESMonitorSet(snes,Monitor,&monP,0);
82: SNESSetFromOptions(snes);
84: /* Solve nonlinear system */
85: FormInitialGuess(snes,x);
86: SNESSolve(snes,NULL,x);
87: SNESGetIterationNumber(snes,&its);
88: PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its);
90: /* Free data structures */
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: const PetscScalar *xx;
105: PetscScalar *ff,*FF,d;
106: PetscInt i,n;
107: PetscErrorCode ierr;
109: VecGetArrayRead(x,&xx);
110: VecGetArray(f,&ff);
111: VecGetArray((Vec) dummy,&FF);
112: VecGetSize(x,&n);
113: d = (PetscReal)(n - 1); d = d*d;
114: ff[0] = xx[0];
115: for (i=1; i<n-1; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
116: ff[n-1] = xx[n-1] - 1.0;
117: VecRestoreArrayRead(x,&xx);
118: VecRestoreArray(f,&ff);
119: VecRestoreArray((Vec)dummy,&FF);
120: return 0;
121: }
123: /*
125: Example function that when differenced produces the same matrix free Jacobian as FormFunction()
126: this is provided to show how a user can provide a different function
127: */
128: PetscErrorCode OtherFunctionForDifferencing(void *dummy,Vec x,Vec f)
129: {
132: FormFunction(NULL,x,f,dummy);
133: VecShift(f,1.0);
134: return 0;
135: }
137: /* -------------------- Form initial approximation ----------------- */
141: PetscErrorCode FormInitialGuess(SNES snes,Vec x)
142: {
144: PetscScalar pfive = .50;
145: VecSet(x,pfive);
146: return 0;
147: }
150: /* -------------------- Evaluate Jacobian F'(x) -------------------- */
151: /* Evaluates a matrix that is used to precondition the matrix-free
152: jacobian. In this case, the explict preconditioner matrix is
153: also EXACTLY the Jacobian. In general, it would be some lower
154: order, simplified apprioximation */
156: PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
157: {
158: const PetscScalar *xx;
159: PetscScalar A[3],d;
160: PetscInt i,n,j[3];
161: PetscErrorCode ierr;
162: AppCtx *user = (AppCtx*) dummy;
164: VecGetArrayRead(x,&xx);
165: VecGetSize(x,&n);
166: d = (PetscReal)(n - 1); d = d*d;
168: i = 0; A[0] = 1.0;
169: MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES);
170: for (i=1; i<n-1; i++) {
171: j[0] = i - 1; j[1] = i; j[2] = i + 1;
172: A[0] = d; A[1] = -2.0*d + 2.0*xx[i]; A[2] = d;
173: MatSetValues(B,1,&i,3,j,A,INSERT_VALUES);
174: }
175: i = n-1; A[0] = 1.0;
176: MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES);
177: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
178: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
179: VecRestoreArrayRead(x,&xx);
181: if (user->variant) {
182: MatMFFDSetBase(jac,x,NULL);
183: }
184: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
185: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
186: return 0;
187: }
188: /* -------------------- User-defined monitor ----------------------- */
192: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *dummy)
193: {
195: MonitorCtx *monP = (MonitorCtx*) dummy;
196: Vec x;
197: MPI_Comm comm;
199: PetscObjectGetComm((PetscObject)snes,&comm);
200: PetscFPrintf(comm,stdout,"iter = %D, SNES Function norm %g \n",its,(double)fnorm);
201: SNESGetSolution(snes,&x);
202: VecView(x,monP->viewer);
203: return 0;
204: }