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

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

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


 34: /* -------- User-defined Routines --------- */

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

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

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

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

 61:   /* Create distributed array to manage the 2d grid */
 62:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,-10,-10,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.da);
 63:   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);

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

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

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

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

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

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

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

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

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

 97:   SNESSetFromOptions(snes);

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

104:   SNESVISetVariableBounds(snes,xl,xu);

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

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

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

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

127:   PetscFinalize();

129:   return 0;
130: }

132: /* -------------------------------------------------------------------- */

136: /*  FormGradient - Evaluates gradient of f.

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

250:     }
251:   }

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

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

263:   PetscLogFlops(67*mx*my);
264:   return(0);
265: }

267: /* ------------------------------------------------------------------- */
270: /*
271:    FormJacobian - Evaluates Jacobian matrix.

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

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

281: */
282: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *ptr)
283: {
284:   AppCtx         *user = (AppCtx*) ptr;
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);}

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

450: /* ------------------------------------------------------------------- */
453: /*
454:    MSA_BoundaryConditions -  Calculates the boundary conditions for
455:    the region.

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

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

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

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

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

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

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

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

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

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

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

551:   }

553:   /* Scale the boundary if desired */

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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