Actual source code: jbearing2.c

petsc-3.7.5 2017-01-01
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: 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:   KSP                ksp;
 76:   AppCtx             user;               /* user-defined work context */
 77:   PetscReal          zero=0.0;           /* lower bound on all variables */

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

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

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


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

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

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

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


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

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

124:   /* The TAO code begins here */

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


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

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

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

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

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

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

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

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

170:   /* Free PETSc data structures */
171:   VecDestroy(&x);
172:   VecDestroy(&xl);
173:   VecDestroy(&xu);
174:   MatDestroy(&user.A);
175:   VecDestroy(&user.B);
176:   /* Free TAO data structures */
177:   TaoDestroy(&tao);

179:   DMDestroy(&user.dm);

181:   PetscFinalize();

183:   return 0;
184: }


187: static PetscReal p(PetscReal xi, PetscReal ecc)
188: {
189:   PetscReal t=1.0+ecc*PetscCosScalar(xi);
190:   return (t*t*t);
191: }

195: PetscErrorCode ComputeB(AppCtx* user)
196: {
198:   PetscInt       i,j,k;
199:   PetscInt       nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
200:   PetscReal      two=2.0, pi=4.0*atan(1.0);
201:   PetscReal      hx,hy,ehxhy;
202:   PetscReal      temp,*b;
203:   PetscReal      ecc=user->ecc;

205:   nx=user->nx;
206:   ny=user->ny;
207:   hx=two*pi/(nx+1.0);
208:   hy=two*user->b/(ny+1.0);
209:   ehxhy = ecc*hx*hy;


212:   /*
213:      Get local grid boundaries
214:   */
215:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
216:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

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

230:   return 0;
231: }

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

250:   nx=user->nx;
251:   ny=user->ny;
252:   hx=two*pi/(nx+1.0);
253:   hy=two*user->b/(ny+1.0);
254:   hxhy=hx*hy;
255:   hxhx=one/(hx*hx);
256:   hyhy=one/(hy*hy);

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

260:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
261:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

263:   VecSet(G, zero);
264:   /*
265:     Get local grid boundaries
266:   */
267:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
268:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

270:   VecGetArray(localX,&x);
271:   VecGetArray(G,&g);

273:   for (i=xs; i< xs+xm; i++){
274:     xi=(i+1)*hx;
275:     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
276:     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
277:     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
278:     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
279:     trule5=trule1; /* L(i,j-1) */
280:     trule6=trule2; /* U(i,j+1) */

282:     vdown=-(trule5+trule2)*hyhy;
283:     vleft=-hxhx*(trule2+trule4);
284:     vright= -hxhx*(trule1+trule3);
285:     vup=-hyhy*(trule1+trule6);
286:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);

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

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

293:        k=0;
294:        if (j>gys){
295:          v[k]=vdown; col[k]=row - gxm; k++;
296:        }

298:        if (i>gxs){
299:          v[k]= vleft; col[k]=row - 1; k++;
300:        }

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

304:        if (i+1 < gxs+gxm){
305:          v[k]= vright; col[k]=row+1; k++;
306:        }

308:        if (j+1 <gys+gym){
309:          v[k]= vup; col[k] = row+gxm; k++;
310:        }
311:        tt=0;
312:        for (kk=0;kk<k;kk++){
313:          tt+=v[kk]*x[col[kk]];
314:        }
315:        row=(j-ys)*xm + (i-xs);
316:        g[row]=tt;

318:      }

320:   }

322:   VecRestoreArray(localX,&x);
323:   VecRestoreArray(G,&g);

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

327:   VecDot(X,G,&f1);
328:   VecDot(user->B,X,&f2);
329:   VecAXPY(G, one, user->B);
330:   *fcn = f1/2.0 + f2;


333:   PetscLogFlops((91 + 10*ym) * xm);
334:   return 0;

336: }


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

359:   nx=user->nx;
360:   ny=user->ny;
361:   hx=two*pi/(nx+1.0);
362:   hy=two*user->b/(ny+1.0);
363:   hxhy=hx*hy;
364:   hxhx=one/(hx*hx);
365:   hyhy=one/(hy*hy);

367:   /*
368:     Get local grid boundaries
369:   */
370:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
371:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);
372:   MatAssembled(hes,&assembled);
373:   if (assembled){MatZeroEntries(hes);}

375:   for (i=xs; i< xs+xm; i++){
376:     xi=(i+1)*hx;
377:     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
378:     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
379:     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
380:     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
381:     trule5=trule1; /* L(i,j-1) */
382:     trule6=trule2; /* U(i,j+1) */

384:     vdown=-(trule5+trule2)*hyhy;
385:     vleft=-hxhx*(trule2+trule4);
386:     vright= -hxhx*(trule1+trule3);
387:     vup=-hyhy*(trule1+trule6);
388:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
389:     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;

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

394:       k=0;
395:       if (j>gys){
396:         v[k]=vdown; col[k]=row - gxm; k++;
397:       }

399:       if (i>gxs){
400:         v[k]= vleft; col[k]=row - 1; k++;
401:       }

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

405:       if (i+1 < gxs+gxm){
406:         v[k]= vright; col[k]=row+1; k++;
407:       }

409:       if (j+1 <gys+gym){
410:         v[k]= vup; col[k] = row+gxm; k++;
411:       }
412:       MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES);

414:     }

416:   }

418:   /*
419:      Assemble matrix, using the 2-step process:
420:      MatAssemblyBegin(), MatAssemblyEnd().
421:      By placing code between these two statements, computations can be
422:      done while messages are in transition.
423:   */
424:   MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY);
425:   MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY);

427:   /*
428:     Tell the matrix we will never add a new nonzero location to the
429:     matrix. If we do it will generate an error.
430:   */
431:   MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
432:   MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE);

434:   PetscLogFlops(9*xm*ym+49*xm);
435:   MatNorm(hes,NORM_1,&hx);
436:   return 0;
437: }

441: PetscErrorCode Monitor(Tao tao, void *ctx)
442: {
443:   PetscErrorCode     ierr;
444:   PetscInt           its;
445:   PetscReal          f,gnorm,cnorm,xdiff;
446:   TaoConvergedReason reason;

449:   TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
450:   if (!(its%5)) {
451:     PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f);
452:   }
453:   return(0);
454: }

458: PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
459: {
460:   PetscErrorCode     ierr;
461:   PetscInt           its;
462:   PetscReal          f,gnorm,cnorm,xdiff;
463:   TaoConvergedReason reason;

466:   TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
467:   if (its == 100) {
468:     TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS);
469:   }
470:   return(0);

472: }