Actual source code: ex10.c

petsc-3.4.5 2014-06-29
  1: #include <petscsnes.h>
  2: #include <../src/snes/impls/vi/viimpl.h>

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

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

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


 37: /* -------- User-defined Routines --------- */

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

 46: int main(int argc, char **argv)
 47: {
 48:   PetscErrorCode info;              /* used to check for functions returning nonzeros */
 49:   Vec            x,r;               /* solution and residual vectors */
 50:   Vec            xl,xu;             /* Bounds on the variables */
 51:   PetscBool      flg;               /* A return variable when checking for user options */
 52:   SNES           snes;              /* nonlinear solver context */
 53:   Mat            J;                 /* Jacobian matrix */
 54:   PetscInt       N;                 /* Number of elements in vector */
 55:   PetscScalar    lb = SNES_VI_NINF; /* lower bound constant */
 56:   PetscScalar    ub = SNES_VI_INF;  /* upper bound constant */
 57:   AppCtx         user;              /* user-defined work context */

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

 62: #if defined(PETSC_USE_COMPLEX)
 63:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This example does not work for scalar type complex\n");
 64: #endif
 65:   /* Specify default dimension of the problem */
 66:   user.mx = 4; user.my = 4;

 68:   /* Check for any command line arguments that override defaults */
 69:   info = PetscOptionsGetInt(NULL, "-mx", &user.mx, &flg);CHKERRQ(info);
 70:   info = PetscOptionsGetInt(NULL, "-my", &user.my, &flg);CHKERRQ(info);
 71:   info = PetscOptionsGetScalar(NULL, "-lb", &lb, &flg);CHKERRQ(info);
 72:   info = PetscOptionsGetScalar(NULL, "-ub", &ub, &flg);CHKERRQ(info);


 75:   /* Calculate any derived values from parameters */
 76:   N = user.mx*user.my;

 78:   PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n");
 79:   PetscPrintf(PETSC_COMM_SELF,"mx:%d, my:%d\n", user.mx,user.my);

 81:   /* Create appropriate vectors and matrices */
 82:   info = VecCreateSeq(MPI_COMM_SELF, N, &x);
 83:   info = VecDuplicate(x, &r);CHKERRQ(info);
 84:   info = MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J);CHKERRQ(info);

 86:   /* Create nonlinear solver context */
 87:   info = SNESCreate(PETSC_COMM_SELF,&snes);CHKERRQ(info);


 90:   /*  Set function evaluation and Jacobian evaluation  routines */
 91:   info = SNESSetFunction(snes,r,FormGradient,&user);CHKERRQ(info);
 92:   info = SNESSetJacobian(snes,J,J,FormJacobian,&user);CHKERRQ(info);

 94:   /* Set the variable bounds */
 95:   info = MSA_BoundaryConditions(&user);CHKERRQ(info);

 97:   /* Set initial solution guess */
 98:   info = MSA_InitialPoint(&user, x);CHKERRQ(info);

100:   info = SNESSetFromOptions(snes);CHKERRQ(info);

102:   /* Set Bounds on variables */
103:   info = VecDuplicate(x, &xl);CHKERRQ(info);
104:   info = VecDuplicate(x, &xu);CHKERRQ(info);
105:   info = VecSet(xl, lb);CHKERRQ(info);
106:   info = VecSet(xu, ub);CHKERRQ(info);

108:   info = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(info);

110:   /* Solve the application */
111:   info = SNESSolve(snes,NULL,x);CHKERRQ(info);

113:   info = VecView(x,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(info);

115:   /* Free memory */
116:   info = VecDestroy(&x);CHKERRQ(info);
117:   info = VecDestroy(&xl);CHKERRQ(info);
118:   info = VecDestroy(&xu);CHKERRQ(info);
119:   info = VecDestroy(&r);CHKERRQ(info);
120:   info = MatDestroy(&J);CHKERRQ(info);
121:   info = SNESDestroy(&snes);CHKERRQ(info);

123:   /* Free user-created data structures */
124:   info = PetscFree(user.bottom);CHKERRQ(info);
125:   info = PetscFree(user.top);CHKERRQ(info);
126:   info = PetscFree(user.left);CHKERRQ(info);
127:   info = PetscFree(user.right);CHKERRQ(info);

129:   info = PetscFinalize();

131:   return 0;
132: }

134: /* -------------------------------------------------------------------- */

138: /*  FormGradient - Evaluates gradient of f.

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

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

160:   /* Initialize vector to zero */
161:   info = VecSet(G, zero);CHKERRQ(info);

163:   /* Get pointers to vector data */
164:   info = VecGetArray(X, &x);CHKERRQ(info);
165:   info = VecGetArray(G, &g);CHKERRQ(info);

167:   /* Compute function over the locally owned part of the mesh */
168:   for (j=0; j<my; j++) {
169:     for (i=0; i< mx; i++) {
170:       row= j*mx + i;

172:       xc = x[row];
173:       xlt=xrb=xl=xr=xb=xt=xc;

175:       if (i==0) { /* left side */
176:         xl  = user->left[j+1];
177:         xlt = user->left[j+2];
178:       } else xl = x[row-1];

180:       if (j==0) { /* bottom side */
181:         xb  = user->bottom[i+1];
182:         xrb = user->bottom[i+2];
183:       } else xb = x[row-mx];

185:       if (i+1 == mx) { /* right side */
186:         xr  = user->right[j+1];
187:         xrb = user->right[j];
188:       } else xr = x[row+1];

190:       if (j+1==0+my) { /* top side */
191:         xt  = user->top[i+1];
192:         xlt = user->top[i];
193:       } else xt = x[row+mx];

195:       if (i>0 && j+1<my) xlt = x[row-1+mx];
196:       if (j>0 && i+1<mx) xrb = x[row+1-mx];

198:       d1 = (xc-xl);
199:       d2 = (xc-xr);
200:       d3 = (xc-xt);
201:       d4 = (xc-xb);
202:       d5 = (xr-xrb);
203:       d6 = (xrb-xb);
204:       d7 = (xlt-xl);
205:       d8 = (xt-xlt);

207:       df1dxc = d1*hydhx;
208:       df2dxc = (d1*hydhx + d4*hxdhy);
209:       df3dxc = d3*hxdhy;
210:       df4dxc = (d2*hydhx + d3*hxdhy);
211:       df5dxc = d2*hydhx;
212:       df6dxc = d4*hxdhy;

214:       d1 /= hx;
215:       d2 /= hx;
216:       d3 /= hy;
217:       d4 /= hy;
218:       d5 /= hy;
219:       d6 /= hx;
220:       d7 /= hy;
221:       d8 /= hx;

223:       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
224:       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
225:       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
226:       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
227:       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
228:       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);

230:       df1dxc /= f1;
231:       df2dxc /= f2;
232:       df3dxc /= f3;
233:       df4dxc /= f4;
234:       df5dxc /= f5;
235:       df6dxc /= f6;

237:       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;

239:     }
240:   }

242:   /* Restore vectors */
243:   info = VecRestoreArray(X, &x);CHKERRQ(info);
244:   info = VecRestoreArray(G, &g);CHKERRQ(info);
245:   info = PetscLogFlops(67*mx*my);CHKERRQ(info);
246:   return 0;
247: }

249: /* ------------------------------------------------------------------- */
252: /*
253:    FormJacobian - Evaluates Jacobian matrix.

255:    Input Parameters:
256: .  snes - SNES context
257: .  X    - input vector
258: .  ptr  - optional user-defined context, as set by SNESSetJacobian()

260:    Output Parameters:
261: .  tH    - Jacobian matrix

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

278:   /* Set various matrix options */
279:   info = MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(info);
280:   info = MatAssembled(H,&assembled);CHKERRQ(info);
281:   if (assembled) {info = MatZeroEntries(H);CHKERRQ(info);}
282:   *flag=SAME_NONZERO_PATTERN;

284:   /* Get pointers to vector data */
285:   info = VecGetArray(X, &x);CHKERRQ(info);

287:   /* Compute Jacobian over the locally owned part of the mesh */
288:   for (i=0; i< mx; i++) {
289:     for (j=0; j<my; j++) {
290:       row= j*mx + i;

292:       xc = x[row];
293:       xlt=xrb=xl=xr=xb=xt=xc;

295:       /* Left side */
296:       if (i==0) {
297:         xl  = user->left[j+1];
298:         xlt = user->left[j+2];
299:       } else xl = x[row-1];

301:       if (j==0) {
302:         xb  = user->bottom[i+1];
303:         xrb = user->bottom[i+2];
304:       } else xb = x[row-mx];

306:       if (i+1 == mx) {
307:         xr  = user->right[j+1];
308:         xrb = user->right[j];
309:       } else xr = x[row+1];

311:       if (j+1==my) {
312:         xt  = user->top[i+1];
313:         xlt = user->top[i];
314:       } else xt = x[row+mx];

316:       if (i>0 && j+1<my) xlt = x[row-1+mx];
317:       if (j>0 && i+1<mx) xrb = x[row+1-mx];

319:       d1 = (xc-xl)/hx;
320:       d2 = (xc-xr)/hx;
321:       d3 = (xc-xt)/hy;
322:       d4 = (xc-xb)/hy;
323:       d5 = (xrb-xr)/hy;
324:       d6 = (xrb-xb)/hx;
325:       d7 = (xlt-xl)/hy;
326:       d8 = (xlt-xt)/hx;

328:       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
329:       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
330:       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
331:       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
332:       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
333:       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);


336:       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
337:            (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
338:       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
339:            (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
340:       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
341:            (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
342:       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
343:            (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);

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

348:       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
349:            hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
350:            (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
351:            (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);

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

355:       k=0;
356:       if (j>0) {
357:         v[k]=hb; col[k]=row - mx; k++;
358:       }

360:       if (j>0 && i < mx-1) {
361:         v[k]=hbr; col[k]=row - mx+1; k++;
362:       }

364:       if (i>0) {
365:         v[k]= hl; col[k]=row - 1; k++;
366:       }

368:       v[k]= hc; col[k]=row; k++;

370:       if (i < mx-1) {
371:         v[k]= hr; col[k]=row+1; k++;
372:       }

374:       if (i>0 && j < my-1) {
375:         v[k]= htl; col[k] = row+mx-1; k++;
376:       }

378:       if (j < my-1) {
379:         v[k]= ht; col[k] = row+mx; k++;
380:       }

382:       /*
383:          Set matrix values using local numbering, which was defined
384:          earlier, in the main routine.
385:       */
386:       info = MatSetValues(H,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(info);
387:     }
388:   }

390:   /* Restore vectors */
391:   info = VecRestoreArray(X,&x);CHKERRQ(info);

393:   /* Assemble the matrix */
394:   info = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(info);
395:   info = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(info);
396:   info = PetscLogFlops(199*mx*my);CHKERRQ(info);
397:   return 0;
398: }

400: /* ------------------------------------------------------------------- */
403: /*
404:    MSA_BoundaryConditions -  Calculates the boundary conditions for
405:    the region.

407:    Input Parameter:
408: .  user - user-defined application context

410:    Output Parameter:
411: .  user - user-defined application context
412: */
413: PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
414: {
415:   PetscErrorCode info;
416:   PetscInt       i,j,k,limit=0,maxits=5;
417:   PetscInt       mx   =user->mx,my=user->my;
418:   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
419:   PetscScalar    one  =1.0, two=2.0, three=3.0, tol=1e-10;
420:   PetscScalar    fnorm,det,hx,hy,xt=0,yt=0;
421:   PetscScalar    u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
422:   PetscScalar    b=-0.5, t=0.5, l=-0.5, r=0.5;
423:   PetscScalar    *boundary;

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

427:   info = PetscMalloc(bsize*sizeof(PetscScalar), &user->bottom);CHKERRQ(info);
428:   info = PetscMalloc(tsize*sizeof(PetscScalar), &user->top);CHKERRQ(info);
429:   info = PetscMalloc(lsize*sizeof(PetscScalar), &user->left);CHKERRQ(info);
430:   info = PetscMalloc(rsize*sizeof(PetscScalar), &user->right);CHKERRQ(info);

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

434:   for (j=0; j<4; j++) {
435:     if (j==0) {
436:       yt       = b;
437:       xt       = l;
438:       limit    = bsize;
439:       boundary = user->bottom;
440:     } else if (j==1) {
441:       yt       = t;
442:       xt       = l;
443:       limit    = tsize;
444:       boundary = user->top;
445:     } else if (j==2) {
446:       yt       = b;
447:       xt       = l;
448:       limit    = lsize;
449:       boundary = user->left;
450:     } else { /* if  (j==3) */
451:       yt       = b;
452:       xt       = r;
453:       limit    = rsize;
454:       boundary = user->right;
455:     }

457:     for (i=0; i<limit; i++) {
458:       u1=xt;
459:       u2=-yt;
460:       for (k=0; k<maxits; k++) {
461:         nf1   = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
462:         nf2   = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
463:         fnorm = PetscSqrtReal(nf1*nf1+nf2*nf2);
464:         if (fnorm <= tol) break;
465:         njac11 = one+u2*u2-u1*u1;
466:         njac12 = two*u1*u2;
467:         njac21 = -two*u1*u2;
468:         njac22 = -one - u1*u1 + u2*u2;
469:         det    = njac11*njac22-njac21*njac12;
470:         u1     = u1-(njac22*nf1-njac12*nf2)/det;
471:         u2     = u2-(njac11*nf2-njac21*nf1)/det;
472:       }

474:       boundary[i]=u1*u1-u2*u2;
475:       if (j==0 || j==1) xt=xt+hx;
476:       else yt=yt+hy; /* if (j==2 || j==3) */
477:     }
478:   }

480:   return 0;
481: }

483: /* ------------------------------------------------------------------- */
486: /*
487:    MSA_InitialPoint - Calculates the initial guess in one of three ways.

489:    Input Parameters:
490: .  user - user-defined application context
491: .  X - vector for initial guess

493:    Output Parameters:
494: .  X - newly computed initial guess
495: */
496: PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
497: {
498:   PetscErrorCode info;
499:   PetscInt       start=-1,i,j;
500:   PetscScalar    zero =0.0;
501:   PetscBool      flg;

503:   info = PetscOptionsGetInt(NULL,"-start",&start,&flg);CHKERRQ(info);

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

507:     info = VecSet(X, zero);CHKERRQ(info);
508:     /* PLogInfo(user,"Min. Surface Area Problem: Start with 0 vector \n"); */


511:   } else { /* Take an average of the boundary conditions */

513:     PetscInt    row;
514:     PetscInt    mx=user->mx,my=user->my;
515:     PetscScalar *x;

517:     /* Get pointers to vector data */
518:     info = VecGetArray(X,&x);CHKERRQ(info);

520:     /* Perform local computations */
521:     for (j=0; j<my; j++) {
522:       for (i=0; i< mx; i++) {
523:         row    =(j)*mx + (i);
524:         x[row] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+
525:                   ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0;
526:       }
527:     }

529:     /* Restore vectors */
530:     info = VecRestoreArray(X,&x);CHKERRQ(info);

532:   }
533:   return 0;
534: }