Actual source code: ex7.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  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,"-n",&n,NULL);
 37:   PetscOptionsHasName(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:   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++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
115:   ff[n-1] = xx[n-1] - 1.0;
116:   VecRestoreArray(x,&xx);
117:   VecRestoreArray(f,&ff);
118:   VecRestoreArray((Vec)dummy,&FF);
119:   return 0;
120: }

122: /*

124:    Example function that when differenced produces the same matrix free Jacobian as FormFunction()
125:    this is provided to show how a user can provide a different function
126: */
127: PetscErrorCode  OtherFunctionForDifferencing(void *dummy,Vec x,Vec f)
128: {

131:   FormFunction(NULL,x,f,dummy);
132:   VecShift(f,1.0);
133:   return 0;
134: }

136: /* --------------------  Form initial approximation ----------------- */

140: PetscErrorCode  FormInitialGuess(SNES snes,Vec x)
141: {
143:   PetscScalar    pfive = .50;
144:   VecSet(x,pfive);
145:   return 0;
146: }
149: /* --------------------  Evaluate Jacobian F'(x) -------------------- */
150: /*  Evaluates a matrix that is used to precondition the matrix-free
151:     jacobian. In this case, the explict preconditioner matrix is
152:     also EXACTLY the Jacobian. In general, it would be some lower
153:     order, simplified apprioximation */

155: PetscErrorCode  FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
156: {
157:   PetscScalar    *xx,A[3],d;
158:   PetscInt       i,n,j[3];
160:   AppCtx         *user = (AppCtx*) dummy;

162:   VecGetArray(x,&xx);
163:   VecGetSize(x,&n);
164:   d    = (PetscReal)(n - 1); d = d*d;

166:   i    = 0; A[0] = 1.0;
167:   MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES);
168:   for (i=1; i<n-1; i++) {
169:     j[0] = i - 1; j[1] = i;                   j[2] = i + 1;
170:     A[0] = d;     A[1] = -2.0*d + 2.0*xx[i];  A[2] = d;
171:     MatSetValues(B,1,&i,3,j,A,INSERT_VALUES);
172:   }
173:   i     = n-1; A[0] = 1.0;
174:   MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES);
175:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
176:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
177:   VecRestoreArray(x,&xx);

179:   if (user->variant) {
180:     MatMFFDSetBase(jac,x,NULL);
181:   }
182:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
183:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
184:   return 0;
185: }
186: /* --------------------  User-defined monitor ----------------------- */

190: PetscErrorCode  Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *dummy)
191: {
193:   MonitorCtx     *monP = (MonitorCtx*) dummy;
194:   Vec            x;
195:   MPI_Comm       comm;

197:   PetscObjectGetComm((PetscObject)snes,&comm);
198:   PetscFPrintf(comm,stdout,"iter = %D, SNES Function norm %g \n",its,(double)fnorm);
199:   SNESGetSolution(snes,&x);
200:   VecView(x,monP->viewer);
201:   return 0;
202: }