Actual source code: ex16.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <petscsnes.h>
  2: #include <petscdmda.h>

  4: static  char help[]=
  5: "This example is an implementation of minimal surface area with \n\
  6: a plate problem from the TAO package (examples/plate2.c) \n\
  7: This example is based on a problem from the MINPACK-2 test suite.\n\
  8: Given a rectangular 2-D domain, boundary values along the edges of \n\
  9: the domain, and a plate represented by lower boundary conditions, \n\
 10: the objective is to find the surface with the minimal area that \n\
 11: satisfies the boundary conditions.\n\
 12: The command line options are:\n\
 13:   -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\
 14:   -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\
 15:   -bheight <ht>, where <ht> = height of the plate\n\
 16:   -start <st>, where <st> =0 for zero vector, <st> != 0 \n\
 17:                for an average of the boundary conditions\n\n";

 19: /*
 20:    User-defined application context - contains data needed by the
 21:    application-provided call-back routines, FormJacobian() and
 22:    FormFunction().
 23: */

 25: typedef struct {
 26:   DM          da;
 27:   Vec         Bottom, Top, Left, Right;
 28:   PetscScalar bheight;
 29:   PetscInt    mx,my,bmx,bmy;
 30: } AppCtx;


 33: /* -------- User-defined Routines --------- */

 35: PetscErrorCode MSA_BoundaryConditions(AppCtx*);
 36: PetscErrorCode MSA_InitialPoint(AppCtx*, Vec);
 37: PetscErrorCode MSA_Plate(Vec,Vec,void*);
 38: PetscErrorCode FormGradient(SNES, Vec, Vec, void*);
 39: PetscErrorCode FormJacobian(SNES, Vec, Mat*, Mat*, MatStructure*,void*);

 43: int main(int argc, char **argv)
 44: {
 46:   Vec            x,r;               /* solution and residual vectors */
 47:   Vec            xl,xu;             /* Bounds on the variables */
 48:   SNES           snes;              /* nonlinear solver context */
 49:   Mat            J;                 /* Jacobian matrix */
 50:   AppCtx         user;              /* user-defined work context */
 51:   PetscBool      flg;

 53:   /* Initialize PETSc */
 54:   PetscInitialize(&argc, &argv, (char*)0, help);

 56: #if defined(PETSC_USE_COMPLEX)
 57:   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example does not work for scalar type complex\n");
 58: #endif

 60:   /* Create distributed array to manage the 2d grid */
 61:   DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,-10,-10,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.da);
 62:   DMDAGetIerr(user.da,PETSC_IGNORE,&user.mx,&user.my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

 64:   user.bheight = 0.1;
 65:   PetscOptionsGetScalar(NULL,"-bheight",&user.bheight,&flg);

 67:   user.bmx = user.mx/2; user.bmy = user.my/2;
 68:   PetscOptionsGetInt(NULL,"-bmx",&user.bmx,&flg);
 69:   PetscOptionsGetInt(NULL,"-bmy",&user.bmy,&flg);

 71:   PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n");
 72:   PetscPrintf(PETSC_COMM_WORLD,"mx:%d, my:%d, bmx:%d, bmy:%d, height:%4.2f\n",
 73:               user.mx,user.my,user.bmx,user.bmy,user.bheight);

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

 79:   DMSetMatType(user.da,MATAIJ);
 80:   DMCreateMatrix(user.da,&J);

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

 85:   /*  Set function evaluation and Jacobian evaluation  routines */
 86:   SNESSetDM(snes,user.da);
 87:   SNESSetFunction(snes,r,FormGradient,&user);
 88:   SNESSetJacobian(snes,J,J,FormJacobian,&user);

 90:   /* Set the boundary conditions */
 91:   MSA_BoundaryConditions(&user);

 93:   /* Set initial solution guess */
 94:   MSA_InitialPoint(&user, x);

 96:   SNESSetFromOptions(snes);

 98:   /* Set Bounds on variables */
 99:   VecDuplicate(x, &xl);
100:   VecDuplicate(x, &xu);
101:   MSA_Plate(xl,xu,&user);

103:   SNESVISetVariableBounds(snes,xl,xu);

105:   /* Solve the application */
106:   SNESSolve(snes,NULL,x);

108:   PetscOptionsHasName(NULL,"-view_sol",&flg);
109:   if (flg) { VecView(x,PETSC_VIEWER_STDOUT_WORLD); }

111:   /* Free memory */
112:   VecDestroy(&x);
113:   VecDestroy(&xl);
114:   VecDestroy(&xu);
115:   VecDestroy(&r);
116:   MatDestroy(&J);
117:   SNESDestroy(&snes);

119:   /* Free user-created data structures */
120:   DMDestroy(&user.da);
121:   VecDestroy(&user.Bottom);
122:   VecDestroy(&user.Top);
123:   VecDestroy(&user.Left);
124:   VecDestroy(&user.Right);

126:   PetscFinalize();

128:   return 0;
129: }

131: /* -------------------------------------------------------------------- */

135: /*  FormGradient - Evaluates gradient of f.

137:     Input Parameters:
138: .   snes  - the SNES context
139: .   X     - input vector
140: .   ptr   - optional user-defined context, as set by SNESSetFunction()

142:     Output Parameters:
143: .   G - vector containing the newly evaluated gradient
144: */
145: PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
146: {
147:   AppCtx      *user = (AppCtx*) ptr;
148:   int         ierr;
149:   PetscInt    i,j;
150:   PetscInt    mx=user->mx, my=user->my;
151:   PetscScalar hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
152:   PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
153:   PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
154:   PetscScalar **g, **x;
155:   PetscInt    xs,xm,ys,ym;
156:   Vec         localX;
157:   PetscScalar *top,*bottom,*left,*right;

160:   /* Initialize vector to zero */
161:   VecSet(G,0.0);

163:   /* Get local vector */
164:   DMGetLocalVector(user->da,&localX);
165:   VecGetArray(user->Top,&top);
166:   VecGetArray(user->Bottom,&bottom);
167:   VecGetArray(user->Left,&left);
168:   VecGetArray(user->Right,&right);

170:   /* Get ghost points */
171:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
172:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
173:   /* Get pointers to local vector data */
174:   DMDAVecGetArray(user->da,localX, &x);
175:   DMDAVecGetArray(user->da,G, &g);

177:   DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
178:   /* Compute function over the locally owned part of the mesh */
179:   for (j=ys; j < ys+ym; j++) {
180:     for (i=xs; i< xs+xm; i++) {

182:       xc = x[j][i];
183:       xlt=xrb=xl=xr=xb=xt=xc;

185:       if (i==0) { /* left side */
186:         xl  = left[j-ys+1];
187:         xlt = left[j-ys+2];
188:       } else xl = x[j][i-1];

190:       if (j==0) { /* bottom side */
191:         xb  = bottom[i-xs+1];
192:         xrb = bottom[i-xs+2];
193:       } else xb = x[j-1][i];

195:       if (i+1 == mx) { /* right side */
196:         xr  = right[j-ys+1];
197:         xrb = right[j-ys];
198:       } else xr = x[j][i+1];

200:       if (j+1==0+my) { /* top side */
201:         xt  = top[i-xs+1];
202:         xlt = top[i-xs];
203:       } else xt = x[j+1][i];

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

208:       d1 = (xc-xl);
209:       d2 = (xc-xr);
210:       d3 = (xc-xt);
211:       d4 = (xc-xb);
212:       d5 = (xr-xrb);
213:       d6 = (xrb-xb);
214:       d7 = (xlt-xl);
215:       d8 = (xt-xlt);

217:       df1dxc = d1*hydhx;
218:       df2dxc = (d1*hydhx + d4*hxdhy);
219:       df3dxc = d3*hxdhy;
220:       df4dxc = (d2*hydhx + d3*hxdhy);
221:       df5dxc = d2*hydhx;
222:       df6dxc = d4*hxdhy;

224:       d1 /= hx;
225:       d2 /= hx;
226:       d3 /= hy;
227:       d4 /= hy;
228:       d5 /= hy;
229:       d6 /= hx;
230:       d7 /= hy;
231:       d8 /= hx;

233:       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
234:       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
235:       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
236:       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
237:       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
238:       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);

240:       df1dxc /= f1;
241:       df2dxc /= f2;
242:       df3dxc /= f3;
243:       df4dxc /= f4;
244:       df5dxc /= f5;
245:       df6dxc /= f6;

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

249:     }
250:   }

252:   /* Restore vectors */
253:   DMDAVecRestoreArray(user->da,localX, &x);
254:   DMDAVecRestoreArray(user->da,G, &g);
255:   DMRestoreLocalVector(user->da,&localX);

257:   VecRestoreArray(user->Left,&left);
258:   VecRestoreArray(user->Top,&top);
259:   VecRestoreArray(user->Bottom,&bottom);
260:   VecRestoreArray(user->Right,&right);

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 *tH, Mat *tHPre, MatStructure *flag, void *ptr)
282: {
283:   AppCtx         *user = (AppCtx*) ptr;
284:   Mat            H     = *tH;
286:   PetscInt       i,j,k;
287:   PetscInt       mx=user->mx, my=user->my;
288:   MatStencil     row,col[7];
289:   PetscScalar    hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
290:   PetscScalar    f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
291:   PetscScalar    hl,hr,ht,hb,hc,htl,hbr;
292:   PetscScalar    **x, v[7];
293:   PetscBool      assembled;
294:   PetscInt       xs,xm,ys,ym;
295:   Vec            localX;
296:   PetscScalar    *top,*bottom,*left,*right;

299:   /* Set various matrix options */
300:   MatAssembled(H,&assembled);
301:   if (assembled) {MatZeroEntries(H);}
302:   *flag=SAME_NONZERO_PATTERN;

304:   /* Get local vectors */
305:   DMGetLocalVector(user->da,&localX);
306:   VecGetArray(user->Top,&top);
307:   VecGetArray(user->Bottom,&bottom);
308:   VecGetArray(user->Left,&left);
309:   VecGetArray(user->Right,&right);

311:   /* Get ghost points */
312:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
313:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

315:   /* Get pointers to vector data */
316:   DMDAVecGetArray(user->da,localX, &x);

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

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

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

337:       /* Right */
338:       if (i+1 == mx) {
339:         xr  = right[j+1];
340:         xrb = right[j];
341:       } else xr = x[j][i+1];

343:       /* Top */
344:       if (j+1==my) {
345:         xt  = top[i+1];
346:         xlt = top[i];
347:       } else xt = x[j+1][i];

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

436:   VecRestoreArray(user->Left,&left);
437:   VecRestoreArray(user->Top,&top);
438:   VecRestoreArray(user->Bottom,&bottom);
439:   VecRestoreArray(user->Right,&right);

441:   /* Assemble the matrix */
442:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
443:   DMDAVecRestoreArray(user->da,localX,&x);
444:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
445:   DMRestoreLocalVector(user->da,&localX);

447:   PetscLogFlops(199*mx*my);
448:   return(0);
449: }

451: /* ------------------------------------------------------------------- */
454: /*
455:    MSA_BoundaryConditions -  Calculates the boundary conditions for
456:    the region.

458:    Input Parameter:
459: .  user - user-defined application context

461:    Output Parameter:
462: .  user - user-defined application context
463: */
464: PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
465: {
467:   PetscInt       i,j,k,limit=0,maxits=5;
468:   PetscInt       mx=user->mx,my=user->my;
469:   PetscInt       xs,ys,xm,ym;
470:   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
471:   PetscScalar    one  =1.0, two=2.0, three=3.0, tol=1e-10;
472:   PetscScalar    fnorm,det,hx,hy,xt=0,yt=0;
473:   PetscScalar    u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
474:   PetscScalar    b=-0.5, t=0.5, l=-0.5, r=0.5;
475:   PetscScalar    *boundary;
476:   Vec            Bottom,Top,Right,Left;
477:   PetscScalar    scl=1.0;
478:   PetscBool      flg;

481:   DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);

483:   bsize=xm+2; lsize=ym+2; rsize=ym+2; tsize=xm+2;

485:   VecCreateMPI(PETSC_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom);
486:   VecCreateMPI(PETSC_COMM_WORLD,tsize,PETSC_DECIDE,&Top);
487:   VecCreateMPI(PETSC_COMM_WORLD,lsize,PETSC_DECIDE,&Left);
488:   VecCreateMPI(PETSC_COMM_WORLD,rsize,PETSC_DECIDE,&Right);

490:   user->Top    = Top;
491:   user->Left   = Left;
492:   user->Bottom = Bottom;
493:   user->Right  = Right;

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

497:   for (j=0; j<4; j++) {
498:     if (j==0) {
499:       yt    = b;
500:       xt    = l+hx*xs;
501:       limit = bsize;
502:       VecGetArray(Bottom,&boundary);
503:     } else if (j==1) {
504:       yt   = t;
505:       xt   = l+hx*xs;
506:       limit= tsize;
507:       VecGetArray(Top,&boundary);
508:     } else if (j==2) {
509:       yt   = b+hy*ys;
510:       xt   = l;
511:       limit= lsize;
512:       VecGetArray(Left,&boundary);
513:     } else { /* if  (j==3) */
514:       yt   = b+hy*ys;
515:       xt   = r;
516:       limit= rsize;
517:       VecGetArray(Right,&boundary);
518:     }

520:     for (i=0; i<limit; i++) {
521:       u1=xt;
522:       u2=-yt;
523:       for (k=0; k<maxits; k++) {
524:         nf1   = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
525:         nf2   = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
526:         fnorm = PetscSqrtReal(nf1*nf1+nf2*nf2);
527:         if (fnorm <= tol) break;
528:         njac11 = one+u2*u2-u1*u1;
529:         njac12 = two*u1*u2;
530:         njac21 = -two*u1*u2;
531:         njac22 = -one - u1*u1 + u2*u2;
532:         det    = njac11*njac22-njac21*njac12;
533:         u1     = u1-(njac22*nf1-njac12*nf2)/det;
534:         u2     = u2-(njac11*nf2-njac21*nf1)/det;
535:       }

537:       boundary[i]=u1*u1-u2*u2;
538:       if (j==0 || j==1) xt=xt+hx;
539:       else yt=yt+hy; /* if (j==2 || j==3) */
540:     }

542:     if (j==0) {
543:       VecRestoreArray(Bottom,&boundary);
544:     } else if (j==1) {
545:       VecRestoreArray(Top,&boundary);
546:     } else if (j==2) {
547:       VecRestoreArray(Left,&boundary);
548:     } else if (j==3) {
549:       VecRestoreArray(Right,&boundary);
550:     }

552:   }

554:   /* Scale the boundary if desired */

556:   PetscOptionsGetReal(NULL,"-bottom",&scl,&flg);
557:   if (flg) {
558:     VecScale(Bottom, scl);
559:   }

561:   PetscOptionsGetReal(NULL,"-top",&scl,&flg);
562:   if (flg) {
563:     VecScale(Top, scl);
564:   }

566:   PetscOptionsGetReal(NULL,"-right",&scl,&flg);
567:   if (flg) {
568:     VecScale(Right, scl);
569:   }

571:   PetscOptionsGetReal(NULL,"-left",&scl,&flg);
572:   if (flg) {
573:     VecScale(Left, scl);
574:   }
575:   return(0);
576: }

578: /* ------------------------------------------------------------------- */
581: /*
582:    MSA_InitialPoint - Calculates the initial guess in one of three ways.

584:    Input Parameters:
585: .  user - user-defined application context
586: .  X - vector for initial guess

588:    Output Parameters:
589: .  X - newly computed initial guess
590: */
591: PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
592: {
594:   PetscInt       start =-1,i,j;
595:   PetscScalar    zero  =0.0;
596:   PetscBool      flg;
597:   PetscScalar    *left,*right,*bottom,*top;

600:   PetscOptionsGetInt(NULL,"-start",&start,&flg);

602:   if (flg && start==0) { /* The zero vector is reasonable */

604:     VecSet(X, zero);
605:     /* PLogIerr(user,"Min. Surface Area Problem: Start with 0 vector \n"); */


608:   } else { /* Take an average of the boundary conditions */
609:     PetscInt    mx=user->mx,my=user->my;
610:     PetscScalar **x;
611:     PetscInt    xs,xm,ys,ym;

613:     VecGetArray(user->Top,&top);
614:     VecGetArray(user->Bottom,&bottom);
615:     VecGetArray(user->Left,&left);
616:     VecGetArray(user->Right,&right);

618:     /* Get pointers to vector data */
619:     DMDAVecGetArray(user->da,X,&x);
620:     DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);

622:     /* Perform local computations */
623:     for (j=ys; j<ys+ym; j++) {
624:       for (i=xs; i< xs+xm; i++) {
625:         x[j][i] = ((j+1)*bottom[i-xs+1]/my+(my-j+1)*top[i-xs+1]/(my+2)+
626:                    (i+1)*left[j-ys+1]/mx+(mx-i+1)*right[j-ys+1]/(mx+2))/2.0;
627:       }
628:     }

630:     /* Restore vectors */
631:     DMDAVecRestoreArray(user->da,X,&x);
632:     VecRestoreArray(user->Left,&left);
633:     VecRestoreArray(user->Top,&top);
634:     VecRestoreArray(user->Bottom,&bottom);
635:     VecRestoreArray(user->Right,&right);
636:   }
637:   return(0);
638: }

642: /*
643:    MSA_Plate -  Calculates an obstacle for surface to stretch over.
644: */
645: PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx)
646: {
647:   AppCtx         *user=(AppCtx*)ctx;
649:   PetscInt       i,j;
650:   PetscInt       xs,ys,xm,ym;
651:   PetscInt       mx=user->mx, my=user->my, bmy, bmx;
652:   PetscScalar    t1,t2,t3;
653:   PetscScalar    **xl;
654:   PetscScalar    lb=-PETSC_INFINITY, ub=PETSC_INFINITY;
655:   PetscBool      cylinder;

657:   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
658:   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);

660:   bmy=user->bmy, bmx=user->bmx;

662:   DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
663:   VecSet(XL, lb);
664:   DMDAVecGetArray(user->da,XL,&xl);
665:   VecSet(XU, ub);

667:   PetscOptionsHasName(NULL,"-cylinder",&cylinder);
668:   /* Compute the optional lower box */
669:   if (cylinder) {
670:     for (i=xs; i< xs+xm; i++) {
671:       for (j=ys; j<ys+ym; j++) {
672:         t1=(2.0*i-mx)*bmy;
673:         t2=(2.0*j-my)*bmx;
674:         t3=bmx*bmx*bmy*bmy;
675:         if (t1*t1 + t2*t2 <= t3) xl[j][i] = user->bheight;
676:       }
677:     }
678:   } else {
679:     /* Compute the optional lower box */
680:     for (i=xs; i< xs+xm; i++) {
681:       for (j=ys; j<ys+ym; j++) {
682:         if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
683:             j>=(my-bmy)/2 && j<my-(my-bmy)/2) {
684:           xl[j][i] = user->bheight;
685:         }
686:       }
687:     }
688:   }

690:   DMDAVecRestoreArray(user->da,XL,&xl);
691:   return(0);
692: }