Actual source code: jbearing2.c
petsc-3.8.0 2017-09-26
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; /* 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;
83: /* Check for any command line arguments that override defaults */
84: PetscOptionsGetInt(NULL,NULL,"-mx",&user.nx,&flg);
85: PetscOptionsGetInt(NULL,NULL,"-my",&user.ny,&flg);
86: PetscOptionsGetReal(NULL,NULL,"-ecc",&user.ecc,&flg);
87: PetscOptionsGetReal(NULL,NULL,"-b",&user.b,&flg);
90: PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n");
91: PetscPrintf(PETSC_COMM_WORLD,"mx: %D, my: %D, ecc: %g \n\n",user.nx,user.ny,(double)user.ecc);
93: /* Let Petsc determine the grid division */
94: Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
96: /*
97: A two dimensional distributed array will help define this problem,
98: which derives from an elliptic PDE on two dimensional domain. From
99: the distributed array, Create the vectors.
100: */
101: 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);
102: DMSetFromOptions(user.dm);
103: DMSetUp(user.dm);
105: /*
106: Extract global and local vectors from DM; the vector user.B is
107: used solely as work space for the evaluation of the function,
108: gradient, and Hessian. Duplicate for remaining vectors that are
109: the same types.
110: */
111: DMCreateGlobalVector(user.dm,&x); /* Solution */
112: VecDuplicate(x,&user.B); /* Linear objective */
115: /* Create matrix user.A to store quadratic, Create a local ordering scheme. */
116: VecGetLocalSize(x,&m);
117: DMCreateMatrix(user.dm,&user.A);
119: /* User defined function -- compute linear term of quadratic */
120: ComputeB(&user);
122: /* The TAO code begins here */
124: /*
125: Create the optimization solver
126: Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
127: */
128: TaoCreate(PETSC_COMM_WORLD,&tao);
129: TaoSetType(tao,TAOBLMVM);
132: /* Set the initial vector */
133: VecSet(x, zero);
134: TaoSetInitialVector(tao,x);
136: /* Set the user function, gradient, hessian evaluation routines and data structures */
137: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*) &user);
139: TaoSetHessianRoutine(tao,user.A,user.A,FormHessian,(void*)&user);
141: /* Set a routine that defines the bounds */
142: VecDuplicate(x,&xl);
143: VecDuplicate(x,&xu);
144: VecSet(xl, zero);
145: VecSet(xu, d1000);
146: TaoSetVariableBounds(tao,xl,xu);
148: TaoGetKSP(tao,&ksp);
149: if (ksp) {
150: KSPSetType(ksp,KSPCG);
151: }
153: PetscOptionsHasName(NULL,NULL,"-testmonitor",&flg);
154: if (flg) {
155: TaoSetMonitor(tao,Monitor,&user,NULL);
156: }
157: PetscOptionsHasName(NULL,NULL,"-testconvergence",&flg);
158: if (flg) {
159: TaoSetConvergenceTest(tao,ConvergenceTest,&user);
160: }
162: /* Check for any tao command line options */
163: TaoSetFromOptions(tao);
165: /* Solve the bound constrained problem */
166: TaoSolve(tao);
168: /* Free PETSc data structures */
169: VecDestroy(&x);
170: VecDestroy(&xl);
171: VecDestroy(&xu);
172: MatDestroy(&user.A);
173: VecDestroy(&user.B);
175: /* Free TAO data structures */
176: TaoDestroy(&tao);
177: DMDestroy(&user.dm);
178: PetscFinalize();
179: return ierr;
180: }
183: static PetscReal p(PetscReal xi, PetscReal ecc)
184: {
185: PetscReal t=1.0+ecc*PetscCosScalar(xi);
186: return (t*t*t);
187: }
189: PetscErrorCode ComputeB(AppCtx* user)
190: {
192: PetscInt i,j,k;
193: PetscInt nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
194: PetscReal two=2.0, pi=4.0*atan(1.0);
195: PetscReal hx,hy,ehxhy;
196: PetscReal temp,*b;
197: PetscReal ecc=user->ecc;
199: nx=user->nx;
200: ny=user->ny;
201: hx=two*pi/(nx+1.0);
202: hy=two*user->b/(ny+1.0);
203: ehxhy = ecc*hx*hy;
206: /*
207: Get local grid boundaries
208: */
209: DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
210: DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);
212: /* Compute the linear term in the objective function */
213: VecGetArray(user->B,&b);
214: for (i=xs; i<xs+xm; i++){
215: temp=PetscSinScalar((i+1)*hx);
216: for (j=ys; j<ys+ym; j++){
217: k=xm*(j-ys)+(i-xs);
218: b[k]= - ehxhy*temp;
219: }
220: }
221: VecRestoreArray(user->B,&b);
222: PetscLogFlops(5*xm*ym+3*xm);
224: return 0;
225: }
227: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr)
228: {
229: AppCtx* user=(AppCtx*)ptr;
231: PetscInt i,j,k,kk;
232: PetscInt col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
233: PetscReal one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
234: PetscReal hx,hy,hxhy,hxhx,hyhy;
235: PetscReal xi,v[5];
236: PetscReal ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
237: PetscReal vmiddle, vup, vdown, vleft, vright;
238: PetscReal tt,f1,f2;
239: PetscReal *x,*g,zero=0.0;
240: Vec localX;
242: nx=user->nx;
243: ny=user->ny;
244: hx=two*pi/(nx+1.0);
245: hy=two*user->b/(ny+1.0);
246: hxhy=hx*hy;
247: hxhx=one/(hx*hx);
248: hyhy=one/(hy*hy);
250: DMGetLocalVector(user->dm,&localX);
252: DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
253: DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);
255: VecSet(G, zero);
256: /*
257: Get local grid boundaries
258: */
259: DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
260: DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);
262: VecGetArray(localX,&x);
263: VecGetArray(G,&g);
265: for (i=xs; i< xs+xm; i++){
266: xi=(i+1)*hx;
267: trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
268: trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
269: trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
270: trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
271: trule5=trule1; /* L(i,j-1) */
272: trule6=trule2; /* U(i,j+1) */
274: vdown=-(trule5+trule2)*hyhy;
275: vleft=-hxhx*(trule2+trule4);
276: vright= -hxhx*(trule1+trule3);
277: vup=-hyhy*(trule1+trule6);
278: vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
280: for (j=ys; j<ys+ym; j++){
282: row=(j-gys)*gxm + (i-gxs);
283: v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
285: k=0;
286: if (j>gys){
287: v[k]=vdown; col[k]=row - gxm; k++;
288: }
290: if (i>gxs){
291: v[k]= vleft; col[k]=row - 1; k++;
292: }
294: v[k]= vmiddle; col[k]=row; k++;
296: if (i+1 < gxs+gxm){
297: v[k]= vright; col[k]=row+1; k++;
298: }
300: if (j+1 <gys+gym){
301: v[k]= vup; col[k] = row+gxm; k++;
302: }
303: tt=0;
304: for (kk=0;kk<k;kk++){
305: tt+=v[kk]*x[col[kk]];
306: }
307: row=(j-ys)*xm + (i-xs);
308: g[row]=tt;
310: }
312: }
314: VecRestoreArray(localX,&x);
315: VecRestoreArray(G,&g);
317: DMRestoreLocalVector(user->dm,&localX);
319: VecDot(X,G,&f1);
320: VecDot(user->B,X,&f2);
321: VecAXPY(G, one, user->B);
322: *fcn = f1/2.0 + f2;
325: PetscLogFlops((91 + 10*ym) * xm);
326: return 0;
328: }
331: /*
332: FormHessian computes the quadratic term in the quadratic objective function
333: Notice that the objective function in this problem is quadratic (therefore a constant
334: hessian). If using a nonquadratic solver, then you might want to reconsider this function
335: */
336: PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
337: {
338: AppCtx* user=(AppCtx*)ptr;
340: PetscInt i,j,k;
341: PetscInt col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
342: PetscReal one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
343: PetscReal hx,hy,hxhy,hxhx,hyhy;
344: PetscReal xi,v[5];
345: PetscReal ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
346: PetscReal vmiddle, vup, vdown, vleft, vright;
347: PetscBool assembled;
349: nx=user->nx;
350: ny=user->ny;
351: hx=two*pi/(nx+1.0);
352: hy=two*user->b/(ny+1.0);
353: hxhy=hx*hy;
354: hxhx=one/(hx*hx);
355: hyhy=one/(hy*hy);
357: /*
358: Get local grid boundaries
359: */
360: DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
361: DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);
362: MatAssembled(hes,&assembled);
363: if (assembled){MatZeroEntries(hes);}
365: for (i=xs; i< xs+xm; i++){
366: xi=(i+1)*hx;
367: trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
368: trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
369: trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
370: trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
371: trule5=trule1; /* L(i,j-1) */
372: trule6=trule2; /* U(i,j+1) */
374: vdown=-(trule5+trule2)*hyhy;
375: vleft=-hxhx*(trule2+trule4);
376: vright= -hxhx*(trule1+trule3);
377: vup=-hyhy*(trule1+trule6);
378: vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
379: v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
381: for (j=ys; j<ys+ym; j++){
382: row=(j-gys)*gxm + (i-gxs);
384: k=0;
385: if (j>gys){
386: v[k]=vdown; col[k]=row - gxm; k++;
387: }
389: if (i>gxs){
390: v[k]= vleft; col[k]=row - 1; k++;
391: }
393: v[k]= vmiddle; col[k]=row; k++;
395: if (i+1 < gxs+gxm){
396: v[k]= vright; col[k]=row+1; k++;
397: }
399: if (j+1 <gys+gym){
400: v[k]= vup; col[k] = row+gxm; k++;
401: }
402: MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES);
404: }
406: }
408: /*
409: Assemble matrix, using the 2-step process:
410: MatAssemblyBegin(), MatAssemblyEnd().
411: By placing code between these two statements, computations can be
412: done while messages are in transition.
413: */
414: MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY);
415: MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY);
417: /*
418: Tell the matrix we will never add a new nonzero location to the
419: matrix. If we do it will generate an error.
420: */
421: MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
422: MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE);
424: PetscLogFlops(9*xm*ym+49*xm);
425: MatNorm(hes,NORM_1,&hx);
426: return 0;
427: }
429: PetscErrorCode Monitor(Tao tao, void *ctx)
430: {
431: PetscErrorCode ierr;
432: PetscInt its;
433: PetscReal f,gnorm,cnorm,xdiff;
434: TaoConvergedReason reason;
437: TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
438: if (!(its%5)) {
439: PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f);
440: }
441: return(0);
442: }
444: PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
445: {
446: PetscErrorCode ierr;
447: PetscInt its;
448: PetscReal f,gnorm,cnorm,xdiff;
449: TaoConvergedReason reason;
452: TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);
453: if (its == 100) {
454: TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS);
455: }
456: return(0);
458: }