Actual source code: jbearing2.c

  1: /*$Id$*/

  3: /*
  4:   Include "tao.h" so we can use TAO solvers with PETSc support.  
  5:   Include "petscda.h" so that we can use distributed arrays (DAs) for managing
  6:   the parallel mesh.
  7: */

  9: #include "petscda.h"
 10: #include "tao.h"
 11: #include <math.h>  /* for cos() sin(0), and atan() */

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

 26: /*T
 27:    Concepts: TAO - Solving a bound constrained minimization problem
 28:    Routines: TaoInitialize(); TaoFinalize();
 29:    Routines: TaoApplicationCreate();  TaoAppDestroy();
 30:    Routines: TaoAppSetObjectiveAndGradientRoutine(); TaoAppSetHessianRoutine();
 31:    Routines: TaoAppSetInitialSolutionVec(); TaoAppSetHessianMat();
 32:    Routines: TaoCreate(); TaoDestroy();
 33:    Routines: TaoSetOptions(); TaoApplyGetKSP()
 34:    Routines: TaoSolveApplication(); TaoGetTerminationReason();
 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:   int         nx,ny;        /* discretization in x, y directions */

 49:   /* Working space */
 50:   DA          da;           /* 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 int FormFunctionGradient(TAO_APPLICATION, Vec, double *,Vec,void *);
 58: static int FormHessian(TAO_APPLICATION,Vec,Mat *, Mat *, MatStructure *, void *);
 59: static int ComputeB(AppCtx*);
 60: static int ComputeBounds(TAO_APPLICATION, Vec, Vec, void *);

 64: int main( int argc, char **argv )
 65: {
 66:   int        info;               /* used to check for functions returning nonzeros */
 67:   int        Nx, Ny;             /* number of processors in x- and y- directions */
 68:   int        m, N;               /* number of local and global elements in vectors */
 69:   Vec        x;                  /* variables vector */
 70:   PetscTruth   flg;              /* A return variable when checking for user options */
 71:   TAO_SOLVER tao;                /* TAO_SOLVER solver context */
 72:   TaoMethod  method = "tao_gpcg";/* default minimization method */
 73:   TAO_APPLICATION jbearingapp;   /* The PETSc application */
 74:   ISLocalToGlobalMapping isltog; /* local-to-global mapping object */
 75:   int nloc;                      /* The number of local elements */
 76:   int *ltog;                     /* mapping of local elements to global elements */
 77:   TaoTerminateReason reason;
 78:   AppCtx     user;               /* user-defined work context */
 79:   PetscReal     zero=0.0;           /* lower bound on all variables */
 80:   KSP      ksp;

 82:   
 83:   /* Initialize PETSC and TAO */
 84:   PetscInitialize( &argc, &argv,(char *)0,help );
 85:   TaoInitialize( &argc, &argv,(char *)0,help );

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

 90:   /* Check for any command line arguments that override defaults */
 91:   info = PetscOptionsGetInt(TAO_NULL,"-mx",&user.nx,&flg); CHKERRQ(info);
 92:   info = PetscOptionsGetInt(TAO_NULL,"-my",&user.ny,&flg); CHKERRQ(info);
 93:   info = PetscOptionsGetReal(TAO_NULL,"-ecc",&user.ecc,&flg); CHKERRQ(info);
 94:   info = PetscOptionsGetReal(TAO_NULL,"-b",&user.b,&flg); CHKERRQ(info);


 97:   PetscPrintf(MPI_COMM_WORLD,"\n---- Journal Bearing Problem -----\n");
 98:   PetscPrintf(MPI_COMM_WORLD,"mx: %d,  my: %d,  ecc: %4.3f \n\n",
 99:               user.nx,user.ny,user.ecc);

101:   /* Calculate any derived values from parameters */
102:   N = user.nx*user.ny; 

104:   /* Let Petsc determine the grid division */
105:   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;

107:   /*
108:      A two dimensional distributed array will help define this problem,
109:      which derives from an elliptic PDE on two dimensional domain.  From
110:      the distributed array, Create the vectors.
111:   */
112:   info = DACreate2d(MPI_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.nx,
113:                     user.ny,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.da); CHKERRQ(info);

115:   /*
116:      Extract global and local vectors from DA; the vector user.B is
117:      used solely as work space for the evaluation of the function, 
118:      gradient, and Hessian.  Duplicate for remaining vectors that are 
119:      the same types.
120:   */
121:   info = DACreateGlobalVector(user.da,&x); CHKERRQ(info); /* Solution */
122:   info = VecDuplicate(x,&user.B); CHKERRQ(info); /* Linear objective */


125:   /*  Create matrix user.A to store quadratic, Create a local ordering scheme. */
126:   info = VecGetLocalSize(x,&m); CHKERRQ(info);
127:   info = MatCreateMPIAIJ(MPI_COMM_WORLD,m,m,N,N,5,TAO_NULL,3,TAO_NULL,&user.A); CHKERRQ(info);
128:   
129:   info = DAGetGlobalIndices(user.da,&nloc,&ltog); CHKERRQ(info);
130:   info = ISLocalToGlobalMappingCreate(MPI_COMM_SELF,nloc,ltog,&isltog); 
131:   CHKERRQ(info);
132:   info = MatSetLocalToGlobalMapping(user.A,isltog); CHKERRQ(info);
133:   info = ISLocalToGlobalMappingDestroy(isltog); CHKERRQ(info);


136:   /* User defined function -- compute linear term of quadratic */
137:   info = ComputeB(&user); CHKERRQ(info);


140:   /* The TAO code begins here */

142:   /* 
143:      Create the optimization solver, Petsc application 
144:      Suitable methods: "tao_gpcg","tao_bqpip","tao_tron","tao_blmvm" 
145:   */
146:   info = TaoCreate(MPI_COMM_WORLD,method,&tao); CHKERRQ(info);
147:   info = TaoApplicationCreate(MPI_COMM_WORLD,&jbearingapp); CHKERRQ(info);

149:   /* Set the initial vector */
150:   info = VecSet(x, zero); CHKERRQ(info);
151:   info = TaoAppSetInitialSolutionVec(jbearingapp,x); CHKERRQ(info);

153:   /* Set the user function, gradient, hessian evaluation routines and data structures */
154:   info = TaoAppSetObjectiveAndGradientRoutine(jbearingapp,FormFunctionGradient,(void*) &user); 
155:   CHKERRQ(info);
156:   
157:   info = TaoAppSetHessianMat(jbearingapp,user.A,user.A); CHKERRQ(info);
158:   info = TaoAppSetHessianRoutine(jbearingapp,FormHessian,(void*)&user); 
159:   CHKERRQ(info);

161:   /* Set a routine that defines the bounds */
162:   info = TaoAppSetVariableBoundsRoutine(jbearingapp,ComputeBounds,(void*)&user); CHKERRQ(info);

164:   info = TaoAppGetKSP(jbearingapp,&ksp); CHKERRQ(info);
165:   if (ksp) {                                         
166:     info = KSPSetType(ksp,KSPCG); CHKERRQ(info);
167:   }



171:   /* Check for any tao command line options */
172:   info = TaoSetOptions(jbearingapp,tao); CHKERRQ(info);

174:   /* Solve the bound constrained problem */
175:   info = TaoSolveApplication(jbearingapp,tao); CHKERRQ(info);

177:   info = TaoGetTerminationReason(tao,&reason); CHKERRQ(info);
178:   if (reason <= 0)
179:     PetscPrintf(MPI_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");

181:   /* Free TAO data structures */
182:   info = TaoDestroy(tao); CHKERRQ(info);
183:   info = TaoAppDestroy(jbearingapp); CHKERRQ(info);

185:   /* Free PETSc data structures */
186:   info = VecDestroy(x); CHKERRQ(info); 
187:   info = MatDestroy(user.A); CHKERRQ(info);
188:   info = VecDestroy(user.B); CHKERRQ(info); 
189:   info = DADestroy(user.da); CHKERRQ(info);

191:   TaoFinalize();
192:   PetscFinalize();

194:   return 0;
195: }

199: static int ComputeBounds(TAO_APPLICATION taoapp, Vec xl, Vec xu, void *ctx){
200:   int info;
201:   PetscReal zero=0.0, d1000=1000;
202:   /* Set the upper and lower bounds */
203:   info = VecSet(xl, zero); CHKERRQ(info);
204:   info = VecSet(xu, d1000); CHKERRQ(info);
205:   return 0;
206: }

208: static PetscReal p(PetscReal xi, PetscReal ecc)
209: { 
210:   PetscReal t=1.0+ecc*cos(xi); 
211:   return (t*t*t); 
212: }

216: int ComputeB(AppCtx* user)
217: {
218:   int i,j,k,info;
219:   int nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
220:   PetscReal two=2.0, pi=4.0*atan(1.0);
221:   PetscReal hx,hy,ehxhy;
222:   PetscReal temp,*b;
223:   PetscReal ecc=user->ecc;

225:   nx=user->nx;
226:   ny=user->ny;
227:   hx=two*pi/(nx+1.0);
228:   hy=two*user->b/(ny+1.0);
229:   ehxhy = ecc*hx*hy;


232:   /*
233:      Get local grid boundaries
234:   */
235:   info = DAGetCorners(user->da,&xs,&ys,TAO_NULL,&xm,&ym,TAO_NULL); CHKERRQ(info);
236:   info = DAGetGhostCorners(user->da,&gxs,&gys,TAO_NULL,&gxm,&gym,TAO_NULL); CHKERRQ(info);
237:   

239:   /* Compute the linear term in the objective function */  
240:   info = VecGetArray(user->B,&b); CHKERRQ(info);
241:   for (i=xs; i<xs+xm; i++){
242:     temp=sin((i+1)*hx);
243:     for (j=ys; j<ys+ym; j++){
244:       k=xm*(j-ys)+(i-xs);
245:       b[k]=  - ehxhy*temp;
246:     }
247:   }
248:   info = VecRestoreArray(user->B,&b); CHKERRQ(info);
249:   info = PetscLogFlops(5*xm*ym+3*xm); CHKERRQ(info);

251:   return 0;
252: }

256: int FormFunctionGradient(TAO_APPLICATION taoapp, Vec X, double *fcn,Vec G,void *ptr)
257: {
258:   AppCtx* user=(AppCtx*)ptr;
259:   int i,j,k,kk,info;
260:   int col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
261:   PetscReal one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
262:   PetscReal hx,hy,hxhy,hxhx,hyhy;
263:   PetscReal xi,v[5];
264:   PetscReal ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
265:   PetscReal vmiddle, vup, vdown, vleft, vright;
266:   PetscReal tt,f1,f2;
267:   PetscReal *x,*g,zero=0.0;
268:   Vec localX;

270:   nx=user->nx;
271:   ny=user->ny;
272:   hx=two*pi/(nx+1.0);
273:   hy=two*user->b/(ny+1.0);
274:   hxhy=hx*hy;
275:   hxhx=one/(hx*hx);
276:   hyhy=one/(hy*hy);

278:   info = DAGetLocalVector(user->da,&localX);CHKERRQ(info);

280:   info = DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX); CHKERRQ(info);
281:   info = DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX); CHKERRQ(info);

283:   info = VecSet(G, zero); CHKERRQ(info);
284:   /*
285:     Get local grid boundaries
286:   */
287:   info = DAGetCorners(user->da,&xs,&ys,TAO_NULL,&xm,&ym,TAO_NULL); CHKERRQ(info);
288:   info = DAGetGhostCorners(user->da,&gxs,&gys,TAO_NULL,&gxm,&gym,TAO_NULL); CHKERRQ(info);
289:   
290:   info = VecGetArray(localX,&x); CHKERRQ(info);
291:   info = VecGetArray(G,&g); CHKERRQ(info);

293:   for (i=xs; i< xs+xm; i++){
294:     xi=(i+1)*hx;
295:     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
296:     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
297:     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
298:     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
299:     trule5=trule1; /* L(i,j-1) */
300:     trule6=trule2; /* U(i,j+1) */

302:     vdown=-(trule5+trule2)*hyhy;
303:     vleft=-hxhx*(trule2+trule4);
304:     vright= -hxhx*(trule1+trule3);
305:     vup=-hyhy*(trule1+trule6);
306:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);

308:     for (j=ys; j<ys+ym; j++){
309:       
310:       row=(j-gys)*gxm + (i-gxs);
311:        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
312:        
313:        k=0;
314:        if (j>gys){ 
315:          v[k]=vdown; col[k]=row - gxm; k++;
316:        }
317:        
318:        if (i>gxs){
319:          v[k]= vleft; col[k]=row - 1; k++;
320:        }

322:        v[k]= vmiddle; col[k]=row; k++;
323:        
324:        if (i+1 < gxs+gxm){
325:          v[k]= vright; col[k]=row+1; k++;
326:        }
327:        
328:        if (j+1 <gys+gym){
329:          v[k]= vup; col[k] = row+gxm; k++;
330:        }
331:        tt=0;
332:        for (kk=0;kk<k;kk++){
333:          tt+=v[kk]*x[col[kk]];
334:        }
335:        row=(j-ys)*xm + (i-xs);
336:        g[row]=tt;

338:      }

340:   }

342:   info = VecRestoreArray(localX,&x); CHKERRQ(info);
343:   info = VecRestoreArray(G,&g); CHKERRQ(info);

345:   info = DARestoreLocalVector(user->da,&localX); CHKERRQ(info);

347:   info = VecDot(X,G,&f1); CHKERRQ(info);
348:   info = VecDot(user->B,X,&f2); CHKERRQ(info);
349:   info = VecAXPY(G, one, user->B); CHKERRQ(info);
350:   *fcn = f1/2.0 + f2;

352:   info = PetscLogFlops((91 + 10*ym) * xm); CHKERRQ(info);
353:   return 0;

355: }



361: /* 
362:    FormHessian computes the quadratic term in the quadratic objective function 
363:    Notice that the objective function in this problem is quadratic (therefore a constant
364:    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
365: */
366: int FormHessian(TAO_APPLICATION taoapp,Vec X,Mat *H, Mat *Hpre, MatStructure *flg, void *ptr)
367: {
368:   AppCtx* user=(AppCtx*)ptr;
369:   int i,j,k,info;
370:   int col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
371:   PetscReal one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
372:   PetscReal hx,hy,hxhy,hxhx,hyhy;
373:   PetscReal xi,v[5];
374:   PetscReal ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
375:   PetscReal vmiddle, vup, vdown, vleft, vright;
376:   Mat hes=*H;
377:   PetscTruth assembled;

379:   nx=user->nx;
380:   ny=user->ny;
381:   hx=two*pi/(nx+1.0);
382:   hy=two*user->b/(ny+1.0);
383:   hxhy=hx*hy;
384:   hxhx=one/(hx*hx);
385:   hyhy=one/(hy*hy);

387:   *flg=SAME_NONZERO_PATTERN;
388:   /*
389:     Get local grid boundaries
390:   */
391:   info = DAGetCorners(user->da,&xs,&ys,TAO_NULL,&xm,&ym,TAO_NULL); CHKERRQ(info);
392:   info = DAGetGhostCorners(user->da,&gxs,&gys,TAO_NULL,&gxm,&gym,TAO_NULL); CHKERRQ(info);
393:   
394:   info = MatAssembled(hes,&assembled); CHKERRQ(info);
395:   if (assembled){info = MatZeroEntries(hes);  CHKERRQ(info);}

397:   for (i=xs; i< xs+xm; i++){
398:     xi=(i+1)*hx;
399:     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
400:     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
401:     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
402:     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
403:     trule5=trule1; /* L(i,j-1) */
404:     trule6=trule2; /* U(i,j+1) */

406:     vdown=-(trule5+trule2)*hyhy;
407:     vleft=-hxhx*(trule2+trule4);
408:     vright= -hxhx*(trule1+trule3);
409:     vup=-hyhy*(trule1+trule6);
410:     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
411:     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;

413:     for (j=ys; j<ys+ym; j++){
414:       row=(j-gys)*gxm + (i-gxs);
415:        
416:       k=0;
417:       if (j>gys){ 
418:         v[k]=vdown; col[k]=row - gxm; k++;
419:       }
420:        
421:       if (i>gxs){
422:         v[k]= vleft; col[k]=row - 1; k++;
423:       }

425:       v[k]= vmiddle; col[k]=row; k++;
426:        
427:       if (i+1 < gxs+gxm){
428:         v[k]= vright; col[k]=row+1; k++;
429:       }
430:        
431:       if (j+1 <gys+gym){
432:         v[k]= vup; col[k] = row+gxm; k++;
433:       }
434:       info = MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES); CHKERRQ(info);
435:        
436:     }

438:   }

440:   /* 
441:      Assemble matrix, using the 2-step process:
442:      MatAssemblyBegin(), MatAssemblyEnd().
443:      By placing code between these two statements, computations can be
444:      done while messages are in transition.
445:   */
446:   info = MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
447:   info = MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY); CHKERRQ(info);

449:   /*
450:     Tell the matrix we will never add a new nonzero location to the
451:     matrix. If we do it will generate an error.
452:   */
453:   info = MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR); CHKERRQ(info);
454:   info = MatSetOption(hes,MAT_SYMMETRIC); CHKERRQ(info);

456:   info = PetscLogFlops(9*xm*ym+49*xm); CHKERRQ(info);

458:   return 0;
459: }