Actual source code: jbearing2.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: /*
  2:   Include "petsctao.h" so we can use TAO solvers
  3:   Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing
  4:   Include "petscksp.h" so we can set KSP type
  5:   the parallel mesh.
  6: */

  8: #include <petsctao.h>
  9: #include <petscdmda.h>

 11: static  char help[]=
 12: "This example demonstrates use of the TAO package to \n\
 13: solve a bound constrained minimization problem.  This example is based on \n\
 14: the problem DPJB from the MINPACK-2 test suite.  This pressure journal \n\
 15: bearing problem is an example of elliptic variational problem defined over \n\
 16: a two dimensional rectangle.  By discretizing the domain into triangular \n\
 17: elements, the pressure surrounding the journal bearing is defined as the \n\
 18: minimum of a quadratic function whose variables are bounded below by zero.\n\
 19: The command line options are:\n\
 20:   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
 21:   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
 22:  \n";

 24: /*T
 25:    Concepts: TAO^Solving a bound constrained minimization problem
 26:    Routines: TaoCreate();
 27:    Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
 28:    Routines: TaoSetHessianRoutine();
 29:    Routines: TaoSetVariableBounds();
 30:    Routines: TaoSetMonitor(); TaoSetConvergenceTest();
 31:    Routines: TaoSetInitialVector();
 32:    Routines: TaoSetFromOptions();
 33:    Routines: TaoSolve();
 34:    Routines: TaoGetConvergedReason(); TaoDestroy();
 35:    Processors: n
 36: T*/

 38: /*
 39:    User-defined application context - contains data needed by the
 40:    application-provided call-back routines, FormFunctionGradient(),
 41:    FormHessian().
 42: */
 43: typedef struct {
 44:   /* problem parameters */
 45:   PetscReal      ecc;          /* test problem parameter */
 46:   PetscReal      b;            /* A dimension of journal bearing */
 47:   PetscInt       nx,ny;        /* discretization in x, y directions */

 49:   /* Working space */
 50:   DM          dm;           /* distributed array data structure */
 51:   Mat         A;            /* Quadratic Objective term */
 52:   Vec         B;            /* Linear Objective term */
 53: } AppCtx;

 55: /* User-defined routines */
 56: static PetscReal p(PetscReal xi, PetscReal ecc);
 57: static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *,Vec,void *);
 58: static PetscErrorCode FormHessian(Tao,Vec,Mat, Mat, void *);
 59: static PetscErrorCode ComputeB(AppCtx*);
 60: static PetscErrorCode Monitor(Tao, void*);
 61: static PetscErrorCode ConvergenceTest(Tao, void*);

 65: int main( int argc, char **argv )
 66: {
 67:   PetscErrorCode     ierr;               /* used to check for functions returning nonzeros */
 68:   PetscInt           Nx, Ny;             /* number of processors in x- and y- directions */
 69:   PetscInt           m;               /* number of local elements in vectors */
 70:   Vec                x;                  /* variables vector */
 71:   Vec                xl,xu;                  /* bounds vectors */
 72:   PetscReal          d1000 = 1000;
 73:   PetscBool          flg;              /* A return variable when checking for user options */
 74:   Tao                tao;                /* Tao solver context */
 75:   TaoConvergedReason reason;
 76:   KSP                ksp;
 77:   AppCtx             user;               /* user-defined work context */
 78:   PetscReal          zero=0.0;           /* lower bound on all variables */

 80:   /* Initialize PETSC and TAO */
 81:   PetscInitialize( &argc, &argv,(char *)0,help );

 83:   /* Set the default values for the problem parameters */
 84:   user.nx = 50; user.ny = 50; user.ecc = 0.1; user.b = 10.0;

 86:   /* Check for any command line arguments that override defaults */
 87:   PetscOptionsGetInt(NULL,"-mx",&user.nx,&flg);
 88:   PetscOptionsGetInt(NULL,"-my",&user.ny,&flg);
 89:   PetscOptionsGetReal(NULL,"-ecc",&user.ecc,&flg);
 90:   PetscOptionsGetReal(NULL,"-b",&user.b,&flg);


 93:   PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n");
 94:   PetscPrintf(PETSC_COMM_WORLD,"mx: %D,  my: %D,  ecc: %g \n\n",user.nx,user.ny,(double)user.ecc);

 96:   /* Let Petsc determine the grid division */
 97:   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;

 99:   /*
100:      A two dimensional distributed array will help define this problem,
101:      which derives from an elliptic PDE on two dimensional domain.  From
102:      the distributed array, Create the vectors.
103:   */
104:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,
105:                       user.nx,user.ny,Nx,Ny,1,1,NULL,NULL,
106:                       &user.dm);

108:   /*
109:      Extract global and local vectors from DM; the vector user.B is
110:      used solely as work space for the evaluation of the function,
111:      gradient, and Hessian.  Duplicate for remaining vectors that are
112:      the same types.
113:   */
114:   DMCreateGlobalVector(user.dm,&x); /* Solution */
115:   VecDuplicate(x,&user.B); /* Linear objective */


118:   /*  Create matrix user.A to store quadratic, Create a local ordering scheme. */
119:   VecGetLocalSize(x,&m);
120:   DMCreateMatrix(user.dm,&user.A);

122:   /* User defined function -- compute linear term of quadratic */
123:   ComputeB(&user);

125:   /* The TAO code begins here */

127:   /*
128:      Create the optimization solver
129:      Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
130:   */
131:   TaoCreate(PETSC_COMM_WORLD,&tao);
132:   TaoSetType(tao,TAOBLMVM);


135:   /* Set the initial vector */
136:   VecSet(x, zero);
137:   TaoSetInitialVector(tao,x);

139:   /* Set the user function, gradient, hessian evaluation routines and data structures */
140:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*) &user);

142:   TaoSetHessianRoutine(tao,user.A,user.A,FormHessian,(void*)&user);

144:   /* Set a routine that defines the bounds */
145:   VecDuplicate(x,&xl);
146:   VecDuplicate(x,&xu);
147:   VecSet(xl, zero);
148:   VecSet(xu, d1000);
149:   TaoSetVariableBounds(tao,xl,xu);

151:   TaoGetKSP(tao,&ksp);
152:   if (ksp) {
153:     KSPSetType(ksp,KSPCG);
154:   }

156:   PetscOptionsHasName(NULL,"-testmonitor",&flg);
157:   if (flg) {
158:     TaoSetMonitor(tao,Monitor,&user,NULL);
159:   }
160:   PetscOptionsHasName(NULL,"-testconvergence",&flg);
161:   if (flg) {
162:     TaoSetConvergenceTest(tao,ConvergenceTest,&user);
163:   }

165:   /* Check for any tao command line options */
166:   TaoSetFromOptions(tao);

168:   /* Solve the bound constrained problem */
169:   TaoSolve(tao);

171:   TaoGetConvergedReason(tao,&reason);
172:   if (reason <= 0) {
173:     PetscPrintf(PETSC_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");
174:   }

176:   /* Free PETSc data structures */
177:   VecDestroy(&x);
178:   VecDestroy(&xl);
179:   VecDestroy(&xu);
180:   MatDestroy(&user.A);
181:   VecDestroy(&user.B);
182:   /* Free TAO data structures */
183:   TaoDestroy(&tao);

185:   DMDestroy(&user.dm);

187:   PetscFinalize();

189:   return 0;
190: }


193: static PetscReal p(PetscReal xi, PetscReal ecc)
194: {
195:   PetscReal t=1.0+ecc*PetscCosScalar(xi);
196:   return (t*t*t);
197: }

201: PetscErrorCode ComputeB(AppCtx* user)
202: {
204:   PetscInt       i,j,k;
205:   PetscInt       nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
206:   PetscReal      two=2.0, pi=4.0*atan(1.0);
207:   PetscReal      hx,hy,ehxhy;
208:   PetscReal      temp,*b;
209:   PetscReal      ecc=user->ecc;

211:   nx=user->nx;
212:   ny=user->ny;
213:   hx=two*pi/(nx+1.0);
214:   hy=two*user->b/(ny+1.0);
215:   ehxhy = ecc*hx*hy;


218:   /*
219:      Get local grid boundaries
220:   */
221:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
222:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

224:   /* Compute the linear term in the objective function */
225:   VecGetArray(user->B,&b);
226:   for (i=xs; i<xs+xm; i++){
227:     temp=PetscSinScalar((i+1)*hx);
228:     for (j=ys; j<ys+ym; j++){
229:       k=xm*(j-ys)+(i-xs);
230:       b[k]=  - ehxhy*temp;
231:     }
232:   }
233:   VecRestoreArray(user->B,&b);
234:   PetscLogFlops(5*xm*ym+3*xm);

236:   return 0;
237: }

241: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr)
242: {
243:   AppCtx*        user=(AppCtx*)ptr;
245:   PetscInt       i,j,k,kk;
246:   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
247:   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
248:   PetscReal      hx,hy,hxhy,hxhx,hyhy;
249:   PetscReal      xi,v[5];
250:   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
251:   PetscReal      vmiddle, vup, vdown, vleft, vright;
252:   PetscReal      tt,f1,f2;
253:   PetscReal      *x,*g,zero=0.0;
254:   Vec            localX;

256:   nx=user->nx;
257:   ny=user->ny;
258:   hx=two*pi/(nx+1.0);
259:   hy=two*user->b/(ny+1.0);
260:   hxhy=hx*hy;
261:   hxhx=one/(hx*hx);
262:   hyhy=one/(hy*hy);

264:   DMGetLocalVector(user->dm,&localX);

266:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
267:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

269:   VecSet(G, zero);
270:   /*
271:     Get local grid boundaries
272:   */
273:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
274:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

276:   VecGetArray(localX,&x);
277:   VecGetArray(G,&g);

279:   for (i=xs; i< xs+xm; i++){
280:     xi=(i+1)*hx;
281:     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
282:     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
283:     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
284:     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
285:     trule5=trule1; /* L(i,j-1) */
286:     trule6=trule2; /* U(i,j+1) */

288:     vdown=-(trule5+trule2)*hyhy;
289:     vleft=-hxhx*(trule2+trule4);
290:     vright= -hxhx*(trule1+trule3);
291:     vup=-hyhy*(trule1+trule6);
292:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);

294:     for (j=ys; j<ys+ym; j++){

296:       row=(j-gys)*gxm + (i-gxs);
297:        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;

299:        k=0;
300:        if (j>gys){
301:          v[k]=vdown; col[k]=row - gxm; k++;
302:        }

304:        if (i>gxs){
305:          v[k]= vleft; col[k]=row - 1; k++;
306:        }

308:        v[k]= vmiddle; col[k]=row; k++;

310:        if (i+1 < gxs+gxm){
311:          v[k]= vright; col[k]=row+1; k++;
312:        }

314:        if (j+1 <gys+gym){
315:          v[k]= vup; col[k] = row+gxm; k++;
316:        }
317:        tt=0;
318:        for (kk=0;kk<k;kk++){
319:          tt+=v[kk]*x[col[kk]];
320:        }
321:        row=(j-ys)*xm + (i-xs);
322:        g[row]=tt;

324:      }

326:   }

328:   VecRestoreArray(localX,&x);
329:   VecRestoreArray(G,&g);

331:   DMRestoreLocalVector(user->dm,&localX);

333:   VecDot(X,G,&f1);
334:   VecDot(user->B,X,&f2);
335:   VecAXPY(G, one, user->B);
336:   *fcn = f1/2.0 + f2;


339:   PetscLogFlops((91 + 10*ym) * xm);
340:   return 0;

342: }


347: /*
348:    FormHessian computes the quadratic term in the quadratic objective function
349:    Notice that the objective function in this problem is quadratic (therefore a constant
350:    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
351: */
352: PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
353: {
354:   AppCtx*        user=(AppCtx*)ptr;
356:   PetscInt       i,j,k;
357:   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
358:   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
359:   PetscReal      hx,hy,hxhy,hxhx,hyhy;
360:   PetscReal      xi,v[5];
361:   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
362:   PetscReal      vmiddle, vup, vdown, vleft, vright;
363:   PetscBool      assembled;

365:   nx=user->nx;
366:   ny=user->ny;
367:   hx=two*pi/(nx+1.0);
368:   hy=two*user->b/(ny+1.0);
369:   hxhy=hx*hy;
370:   hxhx=one/(hx*hx);
371:   hyhy=one/(hy*hy);

373:   /*
374:     Get local grid boundaries
375:   */
376:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
377:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);
378:   MatAssembled(hes,&assembled);
379:   if (assembled){MatZeroEntries(hes);}

381:   for (i=xs; i< xs+xm; i++){
382:     xi=(i+1)*hx;
383:     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
384:     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
385:     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
386:     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
387:     trule5=trule1; /* L(i,j-1) */
388:     trule6=trule2; /* U(i,j+1) */

390:     vdown=-(trule5+trule2)*hyhy;
391:     vleft=-hxhx*(trule2+trule4);
392:     vright= -hxhx*(trule1+trule3);
393:     vup=-hyhy*(trule1+trule6);
394:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
395:     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;

397:     for (j=ys; j<ys+ym; j++){
398:       row=(j-gys)*gxm + (i-gxs);

400:       k=0;
401:       if (j>gys){
402:         v[k]=vdown; col[k]=row - gxm; k++;
403:       }

405:       if (i>gxs){
406:         v[k]= vleft; col[k]=row - 1; k++;
407:       }

409:       v[k]= vmiddle; col[k]=row; k++;

411:       if (i+1 < gxs+gxm){
412:         v[k]= vright; col[k]=row+1; k++;
413:       }

415:       if (j+1 <gys+gym){
416:         v[k]= vup; col[k] = row+gxm; k++;
417:       }
418:       MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES);

420:     }

422:   }

424:   /*
425:      Assemble matrix, using the 2-step process:
426:      MatAssemblyBegin(), MatAssemblyEnd().
427:      By placing code between these two statements, computations can be
428:      done while messages are in transition.
429:   */
430:   MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY);
431:   MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY);

433:   /*
434:     Tell the matrix we will never add a new nonzero location to the
435:     matrix. If we do it will generate an error.
436:   */
437:   MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
438:   MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE);

440:   PetscLogFlops(9*xm*ym+49*xm);
441:   MatNorm(hes,NORM_1,&hx);
442:   return 0;
443: }

447: PetscErrorCode Monitor(Tao tao, void *ctx)
448: {
449:   PetscErrorCode     ierr;
450:   PetscInt           its;
451:   PetscReal          f,gnorm,cnorm,xdiff;
452:   TaoConvergedReason reason;

455:   TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
456:   if (!(its%5)) {
457:     PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f);
458:   }
459:   return(0);
460: }

464: PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
465: {
466:   PetscErrorCode     ierr;
467:   PetscInt           its;
468:   PetscReal          f,gnorm,cnorm,xdiff;
469:   TaoConvergedReason reason;

472:   TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
473:   if (its == 100) {
474:     TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS);
475:   }
476:   return(0);

478: }