Actual source code: ex8.c

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

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

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

 29: typedef struct {
 30:   DM          da;
 31:   PetscScalar *bottom, *top, *left, *right;
 32:   PetscInt    mx,my;
 33: } AppCtx;


 36: /* -------- User-defined Routines --------- */

 38: extern PetscErrorCode MSA_BoundaryConditions(AppCtx*);
 39: extern PetscErrorCode MSA_InitialPoint(AppCtx*, Vec);
 40: extern PetscErrorCode FormGradient(SNES, Vec, Vec, void*);
 41: extern PetscErrorCode FormJacobian(SNES, Vec, Mat*, Mat*, MatStructure*,void*);

 45: int main(int argc, char **argv)
 46: {
 48:   Vec            x,r;               /* solution and residual vectors */
 49:   Vec            xl,xu;             /* Bounds on the variables */
 50:   PetscBool      flg_l,flg_u;       /* flags to check if the bounds are set */
 51:   SNES           snes;              /* nonlinear solver context */
 52:   Mat            J;                 /* Jacobian matrix */
 53:   PetscInt       N;                 /* Number of elements in vector */
 54:   PetscScalar    lb = .05;
 55:   PetscScalar    ub = PETSC_INFINITY;
 56:   AppCtx         user;              /* user-defined work context */
 57:   PetscBool      flg;

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

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

 66:   /* Check if lower and upper bounds are set */
 67:   PetscOptionsGetScalar(NULL, "-lb", &lb, &flg_l);
 68:   PetscOptionsGetScalar(NULL, "-ub", &ub, &flg_u);

 70:   /* Create distributed array to manage the 2d grid */
 71:   DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.da);
 72:   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);
 73:   /* Extract global vectors from DMDA; */
 74:   DMCreateGlobalVector(user.da,&x);
 75:   VecDuplicate(x, &r);

 77:   N    = user.mx*user.my;
 78:   DMSetMatType(user.da,MATAIJ);
 79:   DMCreateMatrix(user.da,&J);

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

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

 88:   /* Set the boundary conditions */
 89:   MSA_BoundaryConditions(&user);

 91:   /* Set initial solution guess */
 92:   MSA_InitialPoint(&user, x);


 95:   /* Set Bounds on variables */
 96:   VecDuplicate(x, &xl);
 97:   VecDuplicate(x, &xu);
 98:   VecSet(xl, lb);
 99:   VecSet(xu, ub);

101:   SNESVISetVariableBounds(snes,xl,xu);

103:   SNESSetFromOptions(snes);

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:   PetscFree(user.bottom);
122:   PetscFree(user.top);
123:   PetscFree(user.left);
124:   PetscFree(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;

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

162:   /* Get local vector */
163:   DMGetLocalVector(user->da,&localX);
164:   /* Get ghost points */
165:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
166:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
167:   /* Get pointer to local vector data */
168:   DMDAVecGetArray(user->da,localX, &x);
169:   DMDAVecGetArray(user->da,G, &g);

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

176:       xc = x[j][i];
177:       xlt=xrb=xl=xr=xb=xt=xc;

179:       if (i==0) { /* left side */
180:         xl  = user->left[j+1];
181:         xlt = user->left[j+2];
182:       } else xl = x[j][i-1];

184:       if (j==0) { /* bottom side */
185:         xb  = user->bottom[i+1];
186:         xrb = user->bottom[i+2];
187:       } else xb = x[j-1][i];

189:       if (i+1 == mx) { /* right side */
190:         xr  = user->right[j+1];
191:         xrb = user->right[j];
192:       } else xr = x[j][i+1];

194:       if (j+1==0+my) { /* top side */
195:         xt  = user->top[i+1];
196:         xlt = user->top[i];
197:       } else xt = x[j+1][i];

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

202:       d1 = (xc-xl);
203:       d2 = (xc-xr);
204:       d3 = (xc-xt);
205:       d4 = (xc-xb);
206:       d5 = (xr-xrb);
207:       d6 = (xrb-xb);
208:       d7 = (xlt-xl);
209:       d8 = (xt-xlt);

211:       df1dxc = d1*hydhx;
212:       df2dxc = (d1*hydhx + d4*hxdhy);
213:       df3dxc = d3*hxdhy;
214:       df4dxc = (d2*hydhx + d3*hxdhy);
215:       df5dxc = d2*hydhx;
216:       df6dxc = d4*hxdhy;

218:       d1 /= hx;
219:       d2 /= hx;
220:       d3 /= hy;
221:       d4 /= hy;
222:       d5 /= hy;
223:       d6 /= hx;
224:       d7 /= hy;
225:       d8 /= hx;

227:       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
228:       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
229:       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
230:       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
231:       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
232:       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);

234:       df1dxc /= f1;
235:       df2dxc /= f2;
236:       df3dxc /= f3;
237:       df4dxc /= f4;
238:       df5dxc /= f5;
239:       df6dxc /= f6;

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

243:     }
244:   }

246:   /* Restore vectors */
247:   DMDAVecRestoreArray(user->da,localX, &x);
248:   DMDAVecRestoreArray(user->da,G, &g);
249:   DMRestoreLocalVector(user->da,&localX);
250:   PetscLogFlops(67*mx*my);
251:   return(0);
252: }

254: /* ------------------------------------------------------------------- */
257: /*
258:    FormJacobian - Evaluates Jacobian matrix.

260:    Input Parameters:
261: .  snes - SNES context
262: .  X    - input vector
263: .  ptr  - optional user-defined context, as set by SNESSetJacobian()

265:    Output Parameters:
266: .  tH    - Jacobian matrix

268: */
269: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat *tH, Mat *tHPre, MatStructure *flag, void *ptr)
270: {
271:   AppCtx         *user = (AppCtx*) ptr;
272:   Mat            H     = *tH;
274:   PetscInt       i,j,k;
275:   PetscInt       mx=user->mx, my=user->my;
276:   MatStencil     row,col[7];
277:   PetscScalar    hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
278:   PetscScalar    f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
279:   PetscScalar    hl,hr,ht,hb,hc,htl,hbr;
280:   PetscScalar    **x, v[7];
281:   PetscBool      assembled;
282:   PetscInt       xs,xm,ys,ym;
283:   Vec            localX;

286:   /* Set various matrix options */
287:   MatAssembled(H,&assembled);
288:   if (assembled) {MatZeroEntries(H);}
289:   *flag=SAME_NONZERO_PATTERN;

291:   /* Get local vector */
292:   DMGetLocalVector(user->da,&localX);
293:   /* Get ghost points */
294:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
295:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

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

300:   DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
301:   /* Compute Jacobian over the locally owned part of the mesh */
302:   for (j=ys; j< ys+ym; j++) {
303:     for (i=xs; i< xs+xm; i++) {
304:       xc = x[j][i];
305:       xlt=xrb=xl=xr=xb=xt=xc;

307:       /* Left */
308:       if (i==0) {
309:         xl  = user->left[j+1];
310:         xlt = user->left[j+2];
311:       } else xl = x[j][i-1];

313:       /* Bottom */
314:       if (j==0) {
315:         xb  = user->bottom[i+1];
316:         xrb = user->bottom[i+2];
317:       } else xb = x[j-1][i];

319:       /* Right */
320:       if (i+1 == mx) {
321:         xr  = user->right[j+1];
322:         xrb = user->right[j];
323:       } else xr = x[j][i+1];

325:       /* Top */
326:       if (j+1==my) {
327:         xt  = user->top[i+1];
328:         xlt = user->top[i];
329:       } else xt = x[j+1][i];

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

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

337:       d1 = (xc-xl)/hx;
338:       d2 = (xc-xr)/hx;
339:       d3 = (xc-xt)/hy;
340:       d4 = (xc-xb)/hy;
341:       d5 = (xrb-xr)/hy;
342:       d6 = (xrb-xb)/hx;
343:       d7 = (xlt-xl)/hy;
344:       d8 = (xlt-xt)/hx;

346:       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
347:       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
348:       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
349:       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
350:       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
351:       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);


354:       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
355:            (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
356:       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
357:            (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
358:       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
359:            (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
360:       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
361:            (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);

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

366:       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
367:            hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
368:            (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
369:            (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);

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

373:       k     =0;
374:       row.i = i;row.j= j;
375:       /* Bottom */
376:       if (j>0) {
377:         v[k]     =hb;
378:         col[k].i = i; col[k].j=j-1; k++;
379:       }

381:       /* Bottom right */
382:       if (j>0 && i < mx -1) {
383:         v[k]     =hbr;
384:         col[k].i = i+1; col[k].j = j-1; k++;
385:       }

387:       /* left */
388:       if (i>0) {
389:         v[k]     = hl;
390:         col[k].i = i-1; col[k].j = j; k++;
391:       }

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

396:       /* Right */
397:       if (i < mx-1) {
398:         v[k]    = hr;
399:         col[k].i= i+1; col[k].j = j;k++;
400:       }

402:       /* Top left */
403:       if (i>0 && j < my-1) {
404:         v[k]     = htl;
405:         col[k].i = i-1;col[k].j = j+1; k++;
406:       }

408:       /* Top */
409:       if (j < my-1) {
410:         v[k]     = ht;
411:         col[k].i = i; col[k].j = j+1; k++;
412:       }

414:       MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES);
415:     }
416:   }

418:   /* Assemble the matrix */
419:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
420:   DMDAVecRestoreArray(user->da,localX,&x);
421:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
422:   DMRestoreLocalVector(user->da,&localX);

424:   PetscLogFlops(199*mx*my);
425:   return(0);
426: }

428: /* ------------------------------------------------------------------- */
431: /*
432:    MSA_BoundaryConditions -  Calculates the boundary conditions for
433:    the region.

435:    Input Parameter:
436: .  user - user-defined application context

438:    Output Parameter:
439: .  user - user-defined application context
440: */
441: PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
442: {
444:   PetscInt       i,j,k,limit=0,maxits=5;
445:   PetscInt       mx   =user->mx,my=user->my;
446:   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
447:   PetscScalar    one  =1.0, two=2.0, three=3.0, tol=1e-10;
448:   PetscScalar    fnorm,det,hx,hy,xt=0,yt=0;
449:   PetscScalar    u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
450:   PetscScalar    b=-0.5, t=0.5, l=-0.5, r=0.5;
451:   PetscScalar    *boundary;

454:   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;

456:   PetscMalloc1(bsize, &user->bottom);
457:   PetscMalloc1(tsize, &user->top);
458:   PetscMalloc1(lsize, &user->left);
459:   PetscMalloc1(rsize, &user->right);

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

463:   for (j=0; j<4; j++) {
464:     if (j==0) {
465:       yt       = b;
466:       xt       = l;
467:       limit    = bsize;
468:       boundary = user->bottom;
469:     } else if (j==1) {
470:       yt       = t;
471:       xt       = l;
472:       limit    = tsize;
473:       boundary = user->top;
474:     } else if (j==2) {
475:       yt       = b;
476:       xt       = l;
477:       limit    = lsize;
478:       boundary = user->left;
479:     } else { /* if  (j==3) */
480:       yt       = b;
481:       xt       = r;
482:       limit    = rsize;
483:       boundary = user->right;
484:     }

486:     for (i=0; i<limit; i++) {
487:       u1=xt;
488:       u2=-yt;
489:       for (k=0; k<maxits; k++) {
490:         nf1   = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
491:         nf2   = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
492:         fnorm = PetscSqrtReal(nf1*nf1+nf2*nf2);
493:         if (fnorm <= tol) break;
494:         njac11 = one+u2*u2-u1*u1;
495:         njac12 = two*u1*u2;
496:         njac21 = -two*u1*u2;
497:         njac22 = -one - u1*u1 + u2*u2;
498:         det    = njac11*njac22-njac21*njac12;
499:         u1     = u1-(njac22*nf1-njac12*nf2)/det;
500:         u2     = u2-(njac11*nf2-njac21*nf1)/det;
501:       }

503:       boundary[i]=u1*u1-u2*u2;
504:       if (j==0 || j==1) xt=xt+hx;
505:       else yt=yt+hy; /* if (j==2 || j==3) */
506:     }
507:   }
508:   return(0);
509: }

511: /* ------------------------------------------------------------------- */
514: /*
515:    MSA_InitialPoint - Calculates the initial guess in one of three ways.

517:    Input Parameters:
518: .  user - user-defined application context
519: .  X - vector for initial guess

521:    Output Parameters:
522: .  X - newly computed initial guess
523: */
524: PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
525: {
527:   PetscInt       start=-1,i,j;
528:   PetscScalar    zero =0.0;
529:   PetscBool      flg;

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

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

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


540:   } else { /* Take an average of the boundary conditions */
541:     PetscInt    mx=user->mx,my=user->my;
542:     PetscScalar **x;
543:     PetscInt    xs,xm,ys,ym;

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

549:     /* Perform local computations */
550:     for (j=ys; j<ys+ym; j++) {
551:       for (i=xs; i< xs+xm; i++) {
552:         x[j][i] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+
553:                    ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0;
554:       }
555:     }

557:     /* Restore vectors */
558:     DMDAVecRestoreArray(user->da,X,&x);

560:   }
561:   return(0);
562: }