Actual source code: jbearing2.c

petsc-3.8.4 2018-03-24
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*);

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

 77:   /* Initialize PETSC and TAO */
 78:   PetscInitialize( &argc, &argv,(char *)0,help );if (ierr) return ierr;

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

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

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

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

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

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


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

120:   if (testgetdiag) {
121:     MatShellSetOperation(user.A,MATOP_GET_DIAGONAL,NULL);
122:   }

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

127:   /* The TAO code begins here */

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


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

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

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

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

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

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

167:   /* Check for any tao command line options */
168:   TaoSetFromOptions(tao);

170:   /* Solve the bound constrained problem */
171:   TaoSolve(tao);

173:   /* Free PETSc data structures */
174:   VecDestroy(&x);
175:   VecDestroy(&xl);
176:   VecDestroy(&xu);
177:   MatDestroy(&user.A);
178:   VecDestroy(&user.B);

180:   /* Free TAO data structures */
181:   TaoDestroy(&tao);
182:   DMDestroy(&user.dm);
183:   PetscFinalize();
184:   return ierr;
185: }


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

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

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


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

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

229:   return 0;
230: }

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

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

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

257:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
258:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

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

267:   VecGetArray(localX,&x);
268:   VecGetArray(G,&g);

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

279:     vdown=-(trule5+trule2)*hyhy;
280:     vleft=-hxhx*(trule2+trule4);
281:     vright= -hxhx*(trule1+trule3);
282:     vup=-hyhy*(trule1+trule6);
283:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);

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

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

290:        k=0;
291:        if (j>gys){
292:          v[k]=vdown; col[k]=row - gxm; k++;
293:        }

295:        if (i>gxs){
296:          v[k]= vleft; col[k]=row - 1; k++;
297:        }

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

301:        if (i+1 < gxs+gxm){
302:          v[k]= vright; col[k]=row+1; k++;
303:        }

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

315:      }

317:   }

319:   VecRestoreArray(localX,&x);
320:   VecRestoreArray(G,&g);

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

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


330:   PetscLogFlops((91 + 10*ym) * xm);
331:   return 0;

333: }


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

354:   nx=user->nx;
355:   ny=user->ny;
356:   hx=two*pi/(nx+1.0);
357:   hy=two*user->b/(ny+1.0);
358:   hxhy=hx*hy;
359:   hxhx=one/(hx*hx);
360:   hyhy=one/(hy*hy);

362:   /*
363:     Get local grid boundaries
364:   */
365:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
366:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);
367:   MatAssembled(hes,&assembled);
368:   if (assembled){MatZeroEntries(hes);}

370:   for (i=xs; i< xs+xm; i++){
371:     xi=(i+1)*hx;
372:     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
373:     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
374:     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
375:     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
376:     trule5=trule1; /* L(i,j-1) */
377:     trule6=trule2; /* U(i,j+1) */

379:     vdown=-(trule5+trule2)*hyhy;
380:     vleft=-hxhx*(trule2+trule4);
381:     vright= -hxhx*(trule1+trule3);
382:     vup=-hyhy*(trule1+trule6);
383:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
384:     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;

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

389:       k=0;
390:       if (j>gys){
391:         v[k]=vdown; col[k]=row - gxm; k++;
392:       }

394:       if (i>gxs){
395:         v[k]= vleft; col[k]=row - 1; k++;
396:       }

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

400:       if (i+1 < gxs+gxm){
401:         v[k]= vright; col[k]=row+1; k++;
402:       }

404:       if (j+1 <gys+gym){
405:         v[k]= vup; col[k] = row+gxm; k++;
406:       }
407:       MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES);

409:     }

411:   }

413:   /*
414:      Assemble matrix, using the 2-step process:
415:      MatAssemblyBegin(), MatAssemblyEnd().
416:      By placing code between these two statements, computations can be
417:      done while messages are in transition.
418:   */
419:   MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY);
420:   MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY);

422:   /*
423:     Tell the matrix we will never add a new nonzero location to the
424:     matrix. If we do it will generate an error.
425:   */
426:   MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
427:   MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE);

429:   PetscLogFlops(9*xm*ym+49*xm);
430:   MatNorm(hes,NORM_1,&hx);
431:   return 0;
432: }

434: PetscErrorCode Monitor(Tao tao, void *ctx)
435: {
436:   PetscErrorCode     ierr;
437:   PetscInt           its;
438:   PetscReal          f,gnorm,cnorm,xdiff;
439:   TaoConvergedReason reason;

442:   TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
443:   if (!(its%5)) {
444:     PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f);
445:   }
446:   return(0);
447: }

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

457:   TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
458:   if (its == 100) {
459:     TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS);
460:   }
461:   return(0);

463: }