Actual source code: ex2.c

  1: /*$Id: ex2.c,v 1.83 2001/08/07 21:31:17 bsmith Exp $*/

  3: static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.n
  4: This example employs a user-defined monitoring routine.nn";

  6: /*T
  7:    Concepts: SNES^basic uniprocessor example
  8:    Concepts: SNES^setting a user-defined monitoring routine
  9:    Processors: 1
 10: T*/

 12: /* 
 13:    Include "petscdraw.h" so that we can use PETSc drawing routines.
 14:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 15:    file automatically includes:
 16:      petsc.h       - base PETSc routines   petscvec.h - vectors
 17:      petscsys.h    - system routines       petscmat.h - matrices
 18:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 19:      petscviewer.h - viewers               petscpc.h  - preconditioners
 20:      petscsles.h   - linear solvers
 21: */

 23:  #include petscsnes.h

 25: /* 
 26:    User-defined routines
 27: */
 28: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 29: extern int FormFunction(SNES,Vec,Vec,void*);
 30: extern int FormInitialGuess(Vec);
 31: extern int Monitor(SNES,int,PetscReal,void *);

 33: /*
 34:    User-defined context for monitoring
 35: */
 36: typedef struct {
 37:    PetscViewer viewer;
 38: } MonitorCtx;

 40: int main(int argc,char **argv)
 41: {
 42:   SNES         snes;                   /* SNES context */
 43:   Vec          x,r,F,U;             /* vectors */
 44:   Mat          J;                      /* Jacobian matrix */
 45:   MonitorCtx   monP;                   /* monitoring context */
 46:   int          ierr,its,n = 5,i,maxit,maxf,size;
 47:   PetscScalar  h,xp,v,none = -1.0;
 48:   PetscReal    atol,rtol,stol,norm;

 50:   PetscInitialize(&argc,&argv,(char *)0,help);
 51:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 52:   if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
 53:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 54:   h = 1.0/(n-1);

 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 57:      Create nonlinear solver context
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 60:   SNESCreate(PETSC_COMM_WORLD,&snes);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:      Create vector data structures; set function evaluation routine
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 65:   /*
 66:      Note that we form 1 vector from scratch and then duplicate as needed.
 67:   */
 68:   VecCreate(PETSC_COMM_WORLD,&x);
 69:   VecSetSizes(x,PETSC_DECIDE,n);
 70:   VecSetFromOptions(x);
 71:   VecDuplicate(x,&r);
 72:   VecDuplicate(x,&F);
 73:   VecDuplicate(x,&U);

 75:   /* 
 76:      Set function evaluation routine and vector
 77:   */
 78:   SNESSetFunction(snes,r,FormFunction,(void*)F);


 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:      Create matrix data structure; set Jacobian evaluation routine
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&J);
 86:   MatSetFromOptions(J);

 88:   /* 
 89:      Set Jacobian matrix data structure and default Jacobian evaluation
 90:      routine. User can override with:
 91:      -snes_fd : default finite differencing approximation of Jacobian
 92:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
 93:                 (unless user explicitly sets preconditioner) 
 94:      -snes_mf_operator : form preconditioning matrix as set by the user,
 95:                          but use matrix-free approx for Jacobian-vector
 96:                          products within Newton-Krylov method
 97:   */

 99:   SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL);

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:      Customize nonlinear solver; set runtime options
103:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   /* 
106:      Set an optional user-defined monitoring routine
107:   */
108:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer);
109:   SNESSetMonitor(snes,Monitor,&monP,0);

111:   /*
112:      Set names for some vectors to facilitate monitoring (optional)
113:   */
114:   PetscObjectSetName((PetscObject)x,"Approximate Solution");
115:   PetscObjectSetName((PetscObject)U,"Exact Solution");

117:   /* 
118:      Set SNES/SLES/KSP/PC runtime options, e.g.,
119:          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
120:   */
121:   SNESSetFromOptions(snes);

123:   /* 
124:      Print parameters used for convergence testing (optional) ... just
125:      to demonstrate this routine; this information is also printed with
126:      the option -snes_view
127:   */
128:   SNESGetTolerances(snes,&atol,&rtol,&stol,&maxit,&maxf);
129:   PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%d, maxf=%dn",atol,rtol,stol,maxit,maxf);

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:      Initialize application:
133:      Store right-hand-side of PDE and exact solution
134:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

136:   xp = 0.0;
137:   for (i=0; i<n; i++) {
138:     v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
139:     VecSetValues(F,1,&i,&v,INSERT_VALUES);
140:     v= xp*xp*xp;
141:     VecSetValues(U,1,&i,&v,INSERT_VALUES);
142:     xp += h;
143:   }

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:      Evaluate initial guess; then solve nonlinear system
147:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:   /*
149:      Note: The user should initialize the vector, x, with the initial guess
150:      for the nonlinear solver prior to calling SNESSolve().  In particular,
151:      to employ an initial guess of zero, the user should explicitly set
152:      this vector to zero by calling VecSet().
153:   */
154:   FormInitialGuess(x);
155:   SNESSolve(snes,x,&its);
156:   PetscPrintf(PETSC_COMM_WORLD,"number of Newton iterations = %dnn",its);

158:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159:      Check solution and clean up
160:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

162:   /* 
163:      Check the error
164:   */
165:   VecAXPY(&none,U,x);
166:   ierr  = VecNorm(x,NORM_2,&norm);
167:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %dn",norm,its);


170:   /*
171:      Free work space.  All PETSc objects should be destroyed when they
172:      are no longer needed.
173:   */
174:   VecDestroy(x);  VecDestroy(r);
175:   VecDestroy(U);  VecDestroy(F);
176:   MatDestroy(J);  SNESDestroy(snes);
177:   PetscViewerDestroy(monP.viewer);
178:   PetscFinalize();

180:   return 0;
181: }
182: /* ------------------------------------------------------------------- */
183: /*
184:    FormInitialGuess - Computes initial guess.

186:    Input/Output Parameter:
187: .  x - the solution vector
188: */
189: int FormInitialGuess(Vec x)
190: {
191:    int    ierr;
192:    PetscScalar pfive = .50;
193:    VecSet(&pfive,x);
194:    return 0;
195: }
196: /* ------------------------------------------------------------------- */
197: /* 
198:    FormFunction - Evaluates nonlinear function, F(x).

200:    Input Parameters:
201: .  snes - the SNES context
202: .  x - input vector
203: .  ctx - optional user-defined context, as set by SNESSetFunction()

205:    Output Parameter:
206: .  f - function vector

208:    Note:
209:    The user-defined context can contain any application-specific data
210:    needed for the function evaluation (such as various parameters, work
211:    vectors, and grid information).  In this program the context is just
212:    a vector containing the right-hand-side of the discretized PDE.
213:  */

215: int FormFunction(SNES snes,Vec x,Vec f,void *ctx)
216: {
217:    Vec    g = (Vec)ctx;
218:    PetscScalar *xx,*ff,*gg,d;
219:    int    i,ierr,n;

221:   /*
222:      Get pointers to vector data.
223:        - For default PETSc vectors, VecGetArray() returns a pointer to
224:          the data array.  Otherwise, the routine is implementation dependent.
225:        - You MUST call VecRestoreArray() when you no longer need access to
226:          the array.
227:   */
228:    VecGetArray(x,&xx);
229:    VecGetArray(f,&ff);
230:    VecGetArray(g,&gg);

232:   /*
233:      Compute function
234:   */
235:    VecGetSize(x,&n);
236:    d = (PetscReal)(n - 1); d = d*d;
237:    ff[0]   = xx[0];
238:    for (i=1; i<n-1; i++) {
239:      ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - gg[i];
240:    }
241:    ff[n-1] = xx[n-1] - 1.0;

243:   /*
244:      Restore vectors
245:   */
246:   VecRestoreArray(x,&xx);
247:   VecRestoreArray(f,&ff);
248:   VecRestoreArray(g,&gg);
249:   return 0;
250: }
251: /* ------------------------------------------------------------------- */
252: /*
253:    FormJacobian - Evaluates Jacobian matrix.

255:    Input Parameters:
256: .  snes - the SNES context
257: .  x - input vector
258: .  dummy - optional user-defined context (not used here)

260:    Output Parameters:
261: .  jac - Jacobian matrix
262: .  B - optionally different preconditioning matrix
263: .  flag - flag indicating matrix structure
264: */

266: int FormJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure*flag,void *dummy)
267: {
268:   PetscScalar *xx,A[3],d;
269:   int    i,n,j[3],ierr;

271:   /*
272:      Get pointer to vector data
273:   */
274:   VecGetArray(x,&xx);

276:   /*
277:      Compute Jacobian entries and insert into matrix.
278:       - Note that in this case we set all elements for a particular
279:         row at once.
280:   */
281:   VecGetSize(x,&n);
282:   d = (PetscReal)(n - 1); d = d*d;

284:   /*
285:      Interior grid points
286:   */
287:   for (i=1; i<n-1; i++) {
288:     j[0] = i - 1; j[1] = i; j[2] = i + 1;
289:     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
290:     MatSetValues(*jac,1,&i,3,j,A,INSERT_VALUES);
291:   }

293:   /*
294:      Boundary points
295:   */
296:   i = 0;   A[0] = 1.0;
297:   MatSetValues(*jac,1,&i,1,&i,A,INSERT_VALUES);
298:   i = n-1; A[0] = 1.0;
299:   MatSetValues(*jac,1,&i,1,&i,A,INSERT_VALUES);

301:   /*
302:      Restore vector
303:   */
304:   VecRestoreArray(x,&xx);

306:   /*
307:      Assemble matrix
308:   */
309:   MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
310:   MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);

312:   return 0;
313: }
314: /* ------------------------------------------------------------------- */
315: /*
316:    Monitor - User-defined monitoring routine that views the
317:    current iterate with an x-window plot.

319:    Input Parameters:
320:    snes - the SNES context
321:    its - iteration number
322:    norm - 2-norm function value (may be estimated)
323:    ctx - optional user-defined context for private data for the 
324:          monitor routine, as set by SNESSetMonitor()

326:    Note:
327:    See the manpage for PetscViewerDrawOpen() for useful runtime options,
328:    such as -nox to deactivate all x-window output.
329:  */
330: int Monitor(SNES snes,int its,PetscReal fnorm,void *ctx)
331: {
332:   int        ierr;
333:   MonitorCtx *monP = (MonitorCtx*) ctx;
334:   Vec        x;

336:   PetscPrintf(PETSC_COMM_WORLD,"iter = %d, SNES Function norm %gn",its,fnorm);
337:   SNESGetSolution(snes,&x);
338:   VecView(x,monP->viewer);
339:   return 0;
340: }