Actual source code: minsurf1.c

  1: #include <petsctao.h>

  3: static char  help[] =
  4: "This example demonstrates use of the TAO package to\n\
  5: solve an unconstrained system of equations.  This example is based on a\n\
  6: problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and\n\
  7: boundary values along the edges of the domain, the objective is to find the\n\
  8: surface with the minimal area that satisfies the boundary conditions.\n\
  9: This application solves this problem using complimentarity -- We are actually\n\
 10: solving the system  (grad f)_i >= 0, if x_i == l_i \n\
 11:                     (grad f)_i = 0, if l_i < x_i < u_i \n\
 12:                     (grad f)_i <= 0, if x_i == u_i  \n\
 13: where f is the function to be minimized. \n\
 14: \n\
 15: The command line options are:\n\
 16:   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
 17:   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
 18:   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";

 20: /*T
 21:    Concepts: TAO^Solving a complementarity problem
 22:    Routines: TaoCreate(); TaoDestroy();

 24:    Processors: 1
 25: T*/

 27: /*
 28:    User-defined application context - contains data needed by the
 29:    application-provided call-back routines, FormFunctionGradient(),
 30:    FormHessian().
 31: */
 32: typedef struct {
 33:   PetscInt  mx, my;
 34:   PetscReal *bottom, *top, *left, *right;
 35: } AppCtx;

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

 39: static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
 40: static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
 41: PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
 42: PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *);

 44: int main(int argc, char **argv)
 45: {
 46:   Vec            x;                 /* solution vector */
 47:   Vec            c;                 /* Constraints function vector */
 48:   Vec            xl,xu;             /* Bounds on the variables */
 49:   PetscBool      flg;               /* A return variable when checking for user options */
 50:   Tao            tao;               /* TAO solver context */
 51:   Mat            J;                 /* Jacobian matrix */
 52:   PetscInt       N;                 /* Number of elements in vector */
 53:   PetscScalar    lb =  PETSC_NINFINITY;      /* lower bound constant */
 54:   PetscScalar    ub =  PETSC_INFINITY;      /* upper bound constant */
 55:   AppCtx         user;                    /* user-defined work context */

 57:   /* Initialize PETSc, TAO */
 58:   PetscInitialize(&argc, &argv, (char *)0, help);

 60:   /* Specify default dimension of the problem */
 61:   user.mx = 4; user.my = 4;

 63:   /* Check for any command line arguments that override defaults */
 64:   PetscOptionsGetInt(NULL,NULL, "-mx", &user.mx, &flg);
 65:   PetscOptionsGetInt(NULL,NULL, "-my", &user.my, &flg);

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

 70:   PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n");
 71:   PetscPrintf(PETSC_COMM_SELF,"mx:%D, my:%D\n", user.mx,user.my);

 73:   /* Create appropriate vectors and matrices */
 74:   VecCreateSeq(MPI_COMM_SELF, N, &x);
 75:   VecDuplicate(x, &c);
 76:   MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J);

 78:   /* The TAO code begins here */

 80:   /* Create TAO solver and set desired solution method */
 81:   TaoCreate(PETSC_COMM_SELF,&tao);
 82:   TaoSetType(tao,TAOSSILS);

 84:   /* Set data structure */
 85:   TaoSetSolution(tao, x);

 87:   /*  Set routines for constraints function and Jacobian evaluation */
 88:   TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user);
 89:   TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user);

 91:   /* Set the variable bounds */
 92:   MSA_BoundaryConditions(&user);

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

 97:   /* Set Bounds on variables */
 98:   VecDuplicate(x, &xl);
 99:   VecDuplicate(x, &xu);
100:   VecSet(xl, lb);
101:   VecSet(xu, ub);
102:   TaoSetVariableBounds(tao,xl,xu);

104:   /* Check for any tao command line options */
105:   TaoSetFromOptions(tao);

107:   /* Solve the application */
108:   TaoSolve(tao);

110:   /* Free Tao data structures */
111:   TaoDestroy(&tao);

113:   /* Free PETSc data structures */
114:   VecDestroy(&x);
115:   VecDestroy(&xl);
116:   VecDestroy(&xu);
117:   VecDestroy(&c);
118:   MatDestroy(&J);

120:   /* Free user-created data structures */
121:   PetscFree(user.bottom);
122:   PetscFree(user.top);
123:   PetscFree(user.left);
124:   PetscFree(user.right);

126:   PetscFinalize();
127:   return 0;
128: }

130: /* -------------------------------------------------------------------- */

132: /*  FormConstraints - Evaluates gradient of f.

134:     Input Parameters:
135: .   tao  - the TAO_APPLICATION context
136: .   X    - input vector
137: .   ptr  - optional user-defined context, as set by TaoSetConstraintsRoutine()

139:     Output Parameters:
140: .   G - vector containing the newly evaluated gradient
141: */
142: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr)
143: {
144:   AppCtx         *user = (AppCtx *) ptr;
145:   PetscInt       i,j,row;
146:   PetscInt       mx=user->mx, my=user->my;
147:   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
148:   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
149:   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
150:   PetscScalar    zero=0.0;
151:   PetscScalar    *g, *x;

153:   /* Initialize vector to zero */
154:   VecSet(G, zero);

156:   /* Get pointers to vector data */
157:   VecGetArray(X, &x);
158:   VecGetArray(G, &g);

160:   /* Compute function over the locally owned part of the mesh */
161:   for (j=0; j<my; j++) {
162:     for (i=0; i< mx; i++) {
163:       row= j*mx + i;

165:       xc = x[row];
166:       xlt=xrb=xl=xr=xb=xt=xc;

168:       if (i==0) { /* left side */
169:         xl= user->left[j+1];
170:         xlt = user->left[j+2];
171:       } else {
172:         xl = x[row-1];
173:       }

175:       if (j==0) { /* bottom side */
176:         xb=user->bottom[i+1];
177:         xrb = user->bottom[i+2];
178:       } else {
179:         xb = x[row-mx];
180:       }

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

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

196:       if (i>0 && j+1<my) {
197:         xlt = x[row-1+mx];
198:       }
199:       if (j>0 && i+1<mx) {
200:         xrb = x[row+1-mx];
201:       }

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

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

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

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

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

242:       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
243:     }
244:   }

246:   /* Restore vectors */
247:   VecRestoreArray(X, &x);
248:   VecRestoreArray(G, &g);
249:   PetscLogFlops(67*mx*my);
250:   return 0;
251: }

253: /* ------------------------------------------------------------------- */
254: /*
255:    FormJacobian - Evaluates Jacobian matrix.

257:    Input Parameters:
258: .  tao  - the TAO_APPLICATION context
259: .  X    - input vector
260: .  ptr  - optional user-defined context, as set by TaoSetJacobian()

262:    Output Parameters:
263: .  tH    - Jacobian matrix

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

279:   /* Set various matrix options */
280:   MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
281:   MatAssembled(H,&assembled);
282:   if (assembled) MatZeroEntries(H);

284:   /* Get pointers to vector data */
285:   VecGetArrayRead(X, &x);

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 {
300:         xl = x[row-1];
301:       }

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

310:       if (i+1 == mx) {
311:         xr  = user->right[j+1];
312:         xrb = user->right[j];
313:       } else {
314:         xr = x[row+1];
315:       }

317:       if (j+1==my) {
318:         xt  = user->top[i+1];
319:         xlt = user->top[i];
320:       }else {
321:         xt = x[row+mx];
322:       }

324:       if (i>0 && j+1<my) {
325:         xlt = x[row-1+mx];
326:       }
327:       if (j>0 && i+1<mx) {
328:         xrb = x[row+1-mx];
329:       }

331:       d1 = (xc-xl)/hx;
332:       d2 = (xc-xr)/hx;
333:       d3 = (xc-xt)/hy;
334:       d4 = (xc-xb)/hy;
335:       d5 = (xrb-xr)/hy;
336:       d6 = (xrb-xb)/hx;
337:       d7 = (xlt-xl)/hy;
338:       d8 = (xlt-xt)/hx;

340:       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
341:       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
342:       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
343:       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
344:       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
345:       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);

347:       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
348:       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
349:       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
350:       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);

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

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

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

360:       k=0;
361:       if (j>0) {
362:         v[k]=hb; col[k]=row - mx; k++;
363:       }

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

369:       if (i>0) {
370:         v[k]= hl; col[k]=row - 1; k++;
371:       }

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

375:       if (i < mx-1) {
376:         v[k]= hr; col[k]=row+1; k++;
377:       }

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

383:       if (j < my-1) {
384:         v[k]= ht; col[k] = row+mx; k++;
385:       }

387:       /*
388:          Set matrix values using local numbering, which was defined
389:          earlier, in the main routine.
390:       */
391:       MatSetValues(H,1,&row,k,col,v,INSERT_VALUES);
392:     }
393:   }

395:   /* Restore vectors */
396:   VecRestoreArrayRead(X,&x);

398:   /* Assemble the matrix */
399:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
400:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
401:   PetscLogFlops(199*mx*my);
402:   return 0;
403: }

405: /* ------------------------------------------------------------------- */
406: /*
407:    MSA_BoundaryConditions -  Calculates the boundary conditions for
408:    the region.

410:    Input Parameter:
411: .  user - user-defined application context

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

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

429:   PetscMalloc1(bsize, &user->bottom);
430:   PetscMalloc1(tsize, &user->top);
431:   PetscMalloc1(lsize, &user->left);
432:   PetscMalloc1(rsize, &user->right);

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

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

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

476:       boundary[i]=u1*u1-u2*u2;
477:       if (j==0 || j==1) {
478:         xt=xt+hx;
479:       } else { /* if (j==2 || j==3) */
480:         yt=yt+hy;
481:       }
482:     }
483:   }
484:   return 0;
485: }

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

491:    Input Parameters:
492: .  user - user-defined application context
493: .  X - vector for initial guess

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

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

506:   if (flg && start==0) { /* The zero vector is reasonable */
507:     VecSet(X, zero);
508:   } else { /* Take an average of the boundary conditions */
509:     PetscInt    row;
510:     PetscInt    mx=user->mx,my=user->my;
511:     PetscScalar *x;

513:     /* Get pointers to vector data */
514:     VecGetArray(X,&x);

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

524:     /* Restore vectors */
525:     VecRestoreArray(X,&x);
526:   }
527:   return 0;
528: }

530: /*TEST

532:    build:
533:       requires: !complex

535:    test:
536:       args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5
537:       requires: !single

539:    test:
540:       suffix: 2
541:       args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5

543: TEST*/