Actual source code: ex58.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: #include <petscsnes.h>
  2: #include <petscdm.h>
  3: #include <petscdmda.h>

  5: static const char help[] = "Parallel version of the minimum surface area problem in 2D using DMDA.\n\
  6:  It solves a system of nonlinear equations in mixed\n\
  7: complementarity form.This example is based on a\n\
  8: problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and\n\
  9: boundary values along the edges of the domain, the objective is to find the\n\
 10: surface with the minimal area that satisfies the boundary conditions.\n\
 11: This application solves this problem using complimentarity -- We are actually\n\
 12: solving the system  (grad f)_i >= 0, if x_i == l_i \n\
 13:                     (grad f)_i = 0, if l_i < x_i < u_i \n\
 14:                     (grad f)_i <= 0, if x_i == u_i  \n\
 15: where f is the function to be minimized. \n\
 16: \n\
 17: The command line options are:\n\
 18:   -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
 19:   -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
 20:   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise\n\
 21:   -lb <value>, lower bound on the variables\n\
 22:   -ub <value>, upper bound on the variables\n\n";

 24: /*
 25:    User-defined application context - contains data needed by the
 26:    application-provided call-back routines, FormJacobian() and
 27:    FormFunction().
 28: */

 30: /*
 31:      This is a new version of the ../tests/ex8.c code

 33:      Run, for example, with the options ./ex58 -snes_vi_monitor -ksp_monitor -mg_levels_ksp_monitor -pc_type mg -pc_mg_levels 2 -pc_mg_galerkin -ksp_type fgmres

 35:      Or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
 36:          multigrid levels, it will be determined automatically based on the number of refinements done)

 38:       ./ex58 -pc_type mg -ksp_monitor  -snes_view -pc_mg_galerkin -snes_grid_sequence 3
 39:              -mg_levels_ksp_monitor -snes_vi_monitor -mg_levels_pc_type sor -pc_mg_type full


 42: */

 44: typedef struct {
 45:   PetscScalar *bottom, *top, *left, *right;
 46:   PetscScalar lb,ub;
 47: } AppCtx;


 50: /* -------- User-defined Routines --------- */

 52: extern PetscErrorCode FormBoundaryConditions(SNES,AppCtx**);
 53: extern PetscErrorCode DestroyBoundaryConditions(AppCtx**);
 54: extern PetscErrorCode ComputeInitialGuess(SNES, Vec,void*);
 55: extern PetscErrorCode FormGradient(SNES, Vec, Vec, void*);
 56: extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void*);
 57: extern PetscErrorCode FormBounds(SNES,Vec,Vec);

 61: int main(int argc, char **argv)
 62: {
 64:   Vec            x,r;               /* solution and residual vectors */
 65:   SNES           snes;              /* nonlinear solver context */
 66:   Mat            J;                 /* Jacobian matrix */
 67:   DM             da;

 69:   PetscInitialize(&argc, &argv, (char*)0, help);

 71:   /* Create distributed array to manage the 2d grid */
 72:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);

 74:   /* Extract global vectors from DMDA; */
 75:   DMCreateGlobalVector(da,&x);
 76:   VecDuplicate(x, &r);

 78:   DMSetMatType(da,MATAIJ);
 79:   DMCreateMatrix(da,&J);

 81:   /* Create nonlinear solver context */
 82:   SNESCreate(PETSC_COMM_WORLD,&snes);
 83:   SNESSetDM(snes,da);

 85:   /*  Set function evaluation and Jacobian evaluation  routines */
 86:   SNESSetFunction(snes,r,FormGradient,NULL);
 87:   SNESSetJacobian(snes,J,J,FormJacobian,NULL);

 89:   SNESSetComputeApplicationContext(snes,(PetscErrorCode (*)(SNES,void**))FormBoundaryConditions,(PetscErrorCode (*)(void**))DestroyBoundaryConditions);

 91:   SNESSetComputeInitialGuess(snes,ComputeInitialGuess,NULL);

 93:   SNESVISetComputeVariableBounds(snes,FormBounds);

 95:   SNESSetFromOptions(snes);

 97:   /* Solve the application */
 98:   SNESSolve(snes,NULL,x);

100:   /* Free memory */
101:   VecDestroy(&x);
102:   VecDestroy(&r);
103:   MatDestroy(&J);
104:   SNESDestroy(&snes);

106:   /* Free user-created data structures */
107:   DMDestroy(&da);

109:   PetscFinalize();
110:   return 0;
111: }

113: /* -------------------------------------------------------------------- */

117: /*  FormBounds - sets the upper and lower bounds

119:     Input Parameters:
120: .   snes  - the SNES context

122:     Output Parameters:
123: .   xl - lower bounds
124: .   xu - upper bounds
125: */
126: PetscErrorCode FormBounds(SNES snes, Vec xl, Vec xu)
127: {
129:   AppCtx         *ctx;

132:   SNESGetApplicationContext(snes,&ctx);
133:   VecSet(xl,ctx->lb);
134:   VecSet(xu,ctx->ub);
135:   return(0);
136: }

138: /* -------------------------------------------------------------------- */

142: /*  FormGradient - Evaluates gradient of f.

144:     Input Parameters:
145: .   snes  - the SNES context
146: .   X     - input vector
147: .   ptr   - optional user-defined context, as set by SNESSetFunction()

149:     Output Parameters:
150: .   G - vector containing the newly evaluated gradient
151: */
152: PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
153: {
154:   AppCtx      *user;
155:   int         ierr;
156:   PetscInt    i,j;
157:   PetscInt    mx, my;
158:   PetscScalar hx,hy, hydhx, hxdhy;
159:   PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
160:   PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
161:   PetscScalar **g, **x;
162:   PetscInt    xs,xm,ys,ym;
163:   Vec         localX;
164:   DM          da;

167:   SNESGetDM(snes,&da);
168:   SNESGetApplicationContext(snes,(void**)&user);
169:   DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
170:   hx   = 1.0/(mx+1);hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;

172:   VecSet(G,0.0);

174:   /* Get local vector */
175:   DMGetLocalVector(da,&localX);
176:   /* Get ghost points */
177:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
178:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
179:   /* Get pointer to local vector data */
180:   DMDAVecGetArray(da,localX, &x);
181:   DMDAVecGetArray(da,G, &g);

183:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
184:   /* Compute function over the locally owned part of the mesh */
185:   for (j=ys; j < ys+ym; j++) {
186:     for (i=xs; i< xs+xm; i++) {

188:       xc = x[j][i];
189:       xlt=xrb=xl=xr=xb=xt=xc;

191:       if (i==0) { /* left side */
192:         xl  = user->left[j+1];
193:         xlt = user->left[j+2];
194:       } else xl = x[j][i-1];

196:       if (j==0) { /* bottom side */
197:         xb  = user->bottom[i+1];
198:         xrb = user->bottom[i+2];
199:       } else xb = x[j-1][i];

201:       if (i+1 == mx) { /* right side */
202:         xr  = user->right[j+1];
203:         xrb = user->right[j];
204:       } else xr = x[j][i+1];

206:       if (j+1==0+my) { /* top side */
207:         xt  = user->top[i+1];
208:         xlt = user->top[i];
209:       } else xt = x[j+1][i];

211:       if (i>0 && j+1<my) xlt = x[j+1][i-1]; /* left top side */
212:       if (j>0 && i+1<mx) xrb = x[j-1][i+1]; /* right bottom */

214:       d1 = (xc-xl);
215:       d2 = (xc-xr);
216:       d3 = (xc-xt);
217:       d4 = (xc-xb);
218:       d5 = (xr-xrb);
219:       d6 = (xrb-xb);
220:       d7 = (xlt-xl);
221:       d8 = (xt-xlt);

223:       df1dxc = d1*hydhx;
224:       df2dxc = (d1*hydhx + d4*hxdhy);
225:       df3dxc = d3*hxdhy;
226:       df4dxc = (d2*hydhx + d3*hxdhy);
227:       df5dxc = d2*hydhx;
228:       df6dxc = d4*hxdhy;

230:       d1 /= hx;
231:       d2 /= hx;
232:       d3 /= hy;
233:       d4 /= hy;
234:       d5 /= hy;
235:       d6 /= hx;
236:       d7 /= hy;
237:       d8 /= hx;

239:       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
240:       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
241:       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
242:       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
243:       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
244:       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);

246:       df1dxc /= f1;
247:       df2dxc /= f2;
248:       df3dxc /= f3;
249:       df4dxc /= f4;
250:       df5dxc /= f5;
251:       df6dxc /= f6;

253:       g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;

255:     }
256:   }

258:   /* Restore vectors */
259:   DMDAVecRestoreArray(da,localX, &x);
260:   DMDAVecRestoreArray(da,G, &g);
261:   DMRestoreLocalVector(da,&localX);
262:   PetscLogFlops(67*mx*my);
263:   return(0);
264: }

266: /* ------------------------------------------------------------------- */
269: /*
270:    FormJacobian - Evaluates Jacobian matrix.

272:    Input Parameters:
273: .  snes - SNES context
274: .  X    - input vector
275: .  ptr  - optional user-defined context, as set by SNESSetJacobian()

277:    Output Parameters:
278: .  tH    - Jacobian matrix

280: */
281: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *ptr)
282: {
283:   AppCtx         *user;
285:   PetscInt       i,j,k;
286:   PetscInt       mx, my;
287:   MatStencil     row,col[7];
288:   PetscScalar    hx, hy, hydhx, hxdhy;
289:   PetscScalar    f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
290:   PetscScalar    hl,hr,ht,hb,hc,htl,hbr;
291:   PetscScalar    **x, v[7];
292:   PetscBool      assembled;
293:   PetscInt       xs,xm,ys,ym;
294:   Vec            localX;
295:   DM             da;

298:   SNESGetDM(snes,&da);
299:   SNESGetApplicationContext(snes,(void**)&user);
300:   DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
301:   hx   = 1.0/(mx+1); hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;

303: /* Set various matrix options */
304:   MatAssembled(H,&assembled);
305:   if (assembled) {MatZeroEntries(H);}

307:   /* Get local vector */
308:   DMGetLocalVector(da,&localX);
309:   /* Get ghost points */
310:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
311:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

313:   /* Get pointers to vector data */
314:   DMDAVecGetArray(da,localX, &x);

316:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
317:   /* Compute Jacobian over the locally owned part of the mesh */
318:   for (j=ys; j< ys+ym; j++) {
319:     for (i=xs; i< xs+xm; i++) {
320:       xc = x[j][i];
321:       xlt=xrb=xl=xr=xb=xt=xc;

323:       /* Left */
324:       if (i==0) {
325:         xl  = user->left[j+1];
326:         xlt = user->left[j+2];
327:       } else xl = x[j][i-1];

329:       /* Bottom */
330:       if (j==0) {
331:         xb  =user->bottom[i+1];
332:         xrb = user->bottom[i+2];
333:       } else xb = x[j-1][i];

335:       /* Right */
336:       if (i+1 == mx) {
337:         xr  =user->right[j+1];
338:         xrb = user->right[j];
339:       } else xr = x[j][i+1];

341:       /* Top */
342:       if (j+1==my) {
343:         xt  =user->top[i+1];
344:         xlt = user->top[i];
345:       } else xt = x[j+1][i];

347:       /* Top left */
348:       if (i>0 && j+1<my) xlt = x[j+1][i-1];

350:       /* Bottom right */
351:       if (j>0 && i+1<mx) xrb = x[j-1][i+1];

353:       d1 = (xc-xl)/hx;
354:       d2 = (xc-xr)/hx;
355:       d3 = (xc-xt)/hy;
356:       d4 = (xc-xb)/hy;
357:       d5 = (xrb-xr)/hy;
358:       d6 = (xrb-xb)/hx;
359:       d7 = (xlt-xl)/hy;
360:       d8 = (xlt-xt)/hx;

362:       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
363:       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
364:       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
365:       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
366:       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
367:       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);


370:       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
371:            (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
372:       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
373:            (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
374:       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
375:            (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
376:       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
377:            (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);

379:       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
380:       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);

382:       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
383:            hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
384:            (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2.0*d1*d4)/(f2*f2*f2) +
385:            (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2.0*d2*d3)/(f4*f4*f4);

387:       hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0;  hc/=2.0;

389:       k     =0;
390:       row.i = i;row.j= j;
391:       /* Bottom */
392:       if (j>0) {
393:         v[k]     =hb;
394:         col[k].i = i; col[k].j=j-1; k++;
395:       }

397:       /* Bottom right */
398:       if (j>0 && i < mx -1) {
399:         v[k]     =hbr;
400:         col[k].i = i+1; col[k].j = j-1; k++;
401:       }

403:       /* left */
404:       if (i>0) {
405:         v[k]     = hl;
406:         col[k].i = i-1; col[k].j = j; k++;
407:       }

409:       /* Centre */
410:       v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++;

412:       /* Right */
413:       if (i < mx-1) {
414:         v[k]    = hr;
415:         col[k].i= i+1; col[k].j = j;k++;
416:       }

418:       /* Top left */
419:       if (i>0 && j < my-1) {
420:         v[k]     = htl;
421:         col[k].i = i-1;col[k].j = j+1; k++;
422:       }

424:       /* Top */
425:       if (j < my-1) {
426:         v[k]     = ht;
427:         col[k].i = i; col[k].j = j+1; k++;
428:       }

430:       MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES);
431:     }
432:   }

434:   /* Assemble the matrix */
435:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
436:   DMDAVecRestoreArray(da,localX,&x);
437:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
438:   DMRestoreLocalVector(da,&localX);

440:   PetscLogFlops(199*mx*my);
441:   return(0);
442: }

444: /* ------------------------------------------------------------------- */
447: /*
448:    FormBoundaryConditions -  Calculates the boundary conditions for
449:    the region.

451:    Input Parameter:
452: .  user - user-defined application context

454:    Output Parameter:
455: .  user - user-defined application context
456: */
457: PetscErrorCode FormBoundaryConditions(SNES snes,AppCtx **ouser)
458: {
460:   PetscInt       i,j,k,limit=0,maxits=5;
461:   PetscInt       mx,my;
462:   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
463:   PetscScalar    one  =1.0, two=2.0, three=3.0;
464:   PetscScalar    det,hx,hy,xt=0,yt=0;
465:   PetscReal      fnorm, tol=1e-10;
466:   PetscScalar    u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
467:   PetscScalar    b=-0.5, t=0.5, l=-0.5, r=0.5;
468:   PetscScalar    *boundary;
469:   AppCtx         *user;
470:   DM             da;

473:   SNESGetDM(snes,&da);
474:   PetscNew(&user);
475:   *ouser   = user;
476:   user->lb = .05;
477:   user->ub = PETSC_INFINITY;
478:   DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

480:   /* Check if lower and upper bounds are set */
481:   PetscOptionsGetScalar(NULL,NULL, "-lb", &user->lb, 0);
482:   PetscOptionsGetScalar(NULL,NULL, "-ub", &user->ub, 0);
483:   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;

485:   PetscMalloc1(bsize, &user->bottom);
486:   PetscMalloc1(tsize, &user->top);
487:   PetscMalloc1(lsize, &user->left);
488:   PetscMalloc1(rsize, &user->right);

490:   hx= (r-l)/(mx+1.0); hy=(t-b)/(my+1.0);

492:   for (j=0; j<4; j++) {
493:     if (j==0) {
494:       yt       = b;
495:       xt       = l;
496:       limit    = bsize;
497:       boundary = user->bottom;
498:     } else if (j==1) {
499:       yt       = t;
500:       xt       = l;
501:       limit    = tsize;
502:       boundary = user->top;
503:     } else if (j==2) {
504:       yt       = b;
505:       xt       = l;
506:       limit    = lsize;
507:       boundary = user->left;
508:     } else { /* if  (j==3) */
509:       yt       = b;
510:       xt       = r;
511:       limit    = rsize;
512:       boundary = user->right;
513:     }

515:     for (i=0; i<limit; i++) {
516:       u1=xt;
517:       u2=-yt;
518:       for (k=0; k<maxits; k++) {
519:         nf1   = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
520:         nf2   = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
521:         fnorm = PetscRealPart(PetscSqrtScalar(nf1*nf1+nf2*nf2));
522:         if (fnorm <= tol) break;
523:         njac11=one+u2*u2-u1*u1;
524:         njac12=two*u1*u2;
525:         njac21=-two*u1*u2;
526:         njac22=-one - u1*u1 + u2*u2;
527:         det   = njac11*njac22-njac21*njac12;
528:         u1    = u1-(njac22*nf1-njac12*nf2)/det;
529:         u2    = u2-(njac11*nf2-njac21*nf1)/det;
530:       }

532:       boundary[i]=u1*u1-u2*u2;
533:       if (j==0 || j==1) xt=xt+hx;
534:       else yt=yt+hy; /* if (j==2 || j==3) */
535:     }
536:   }
537:   return(0);
538: }

542: PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser)
543: {
545:   AppCtx         *user = *ouser;

548:   PetscFree(user->bottom);
549:   PetscFree(user->top);
550:   PetscFree(user->left);
551:   PetscFree(user->right);
552:   PetscFree(*ouser);
553:   return(0);
554: }


557: /* ------------------------------------------------------------------- */
560: /*
561:    ComputeInitialGuess - Calculates the initial guess

563:    Input Parameters:
564: .  user - user-defined application context
565: .  X - vector for initial guess

567:    Output Parameters:
568: .  X - newly computed initial guess
569: */
570: PetscErrorCode ComputeInitialGuess(SNES snes, Vec X,void *dummy)
571: {
573:   PetscInt       i,j,mx,my;
574:   DM             da;
575:   AppCtx         *user;
576:   PetscScalar    **x;
577:   PetscInt       xs,xm,ys,ym;

580:   SNESGetDM(snes,&da);
581:   SNESGetApplicationContext(snes,(void**)&user);

583:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
584:   DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

586:   /* Get pointers to vector data */
587:   DMDAVecGetArray(da,X,&x);
588:   /* Perform local computations */
589:   for (j=ys; j<ys+ym; j++) {
590:     for (i=xs; i< xs+xm; i++) {
591:       x[j][i] = (((j+1.0)*user->bottom[i+1]+(my-j+1.0)*user->top[i+1])/(my+2.0)+((i+1.0)*user->left[j+1]+(mx-i+1.0)*user->right[j+1])/(mx+2.0))/2.0;
592:     }
593:   }
594:   /* Restore vectors */
595:   DMDAVecRestoreArray(da,X,&x);
596:   return(0);
597: }