Actual source code: jbearing2.c

petsc-master 2016-12-09
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 );if (ierr) return ierr;

 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,user.nx,user.ny,Nx,Ny,1,1,NULL,NULL,&user.dm);
104:   DMSetFromOptions(user.dm);
105:   DMSetUp(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);

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


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

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

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


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

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

228:   return 0;
229: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

316:      }

318:   }

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

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

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


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

334: }


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

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

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

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

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

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

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

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

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

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

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

412:     }

414:   }

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

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

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

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

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

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

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

470: }