xref: /petsc/src/tao/bound/tutorials/jbearing2.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3) !
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 */
7 
8 #include <petsctao.h>
9 #include <petscdmda.h>
10 
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";
23 
24 /*
25    User-defined application context - contains data needed by the
26    application-provided call-back routines, FormFunctionGradient(),
27    FormHessian().
28 */
29 typedef struct {
30   /* problem parameters */
31   PetscReal      ecc;          /* test problem parameter */
32   PetscReal      b;            /* A dimension of journal bearing */
33   PetscInt       nx,ny;        /* discretization in x, y directions */
34 
35   /* Working space */
36   DM          dm;           /* distributed array data structure */
37   Mat         A;            /* Quadratic Objective term */
38   Vec         B;            /* Linear Objective term */
39 } AppCtx;
40 
41 /* User-defined routines */
42 static PetscReal p(PetscReal xi, PetscReal ecc);
43 static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *,Vec,void *);
44 static PetscErrorCode FormHessian(Tao,Vec,Mat, Mat, void *);
45 static PetscErrorCode ComputeB(AppCtx*);
46 static PetscErrorCode Monitor(Tao, void*);
47 static PetscErrorCode ConvergenceTest(Tao, void*);
48 
49 int main(int argc, char **argv)
50 {
51   PetscInt           Nx, Ny;          /* number of processors in x- and y- directions */
52   PetscInt           m;               /* number of local elements in vectors */
53   Vec                x;               /* variables vector */
54   Vec                xl,xu;           /* bounds vectors */
55   PetscReal          d1000 = 1000;
56   PetscBool          flg,testgetdiag; /* A return variable when checking for user options */
57   Tao                tao;             /* Tao solver context */
58   KSP                ksp;
59   AppCtx             user;            /* user-defined work context */
60   PetscReal          zero = 0.0;      /* lower bound on all variables */
61 
62   /* Initialize PETSC and TAO */
63   PetscFunctionBeginUser;
64   PetscCall(PetscInitialize(&argc, &argv,(char *)0,help));
65 
66   /* Set the default values for the problem parameters */
67   user.nx = 50; user.ny = 50; user.ecc = 0.1; user.b = 10.0;
68   testgetdiag = PETSC_FALSE;
69 
70   /* Check for any command line arguments that override defaults */
71   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.nx,&flg));
72   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.ny,&flg));
73   PetscCall(PetscOptionsGetReal(NULL,NULL,"-ecc",&user.ecc,&flg));
74   PetscCall(PetscOptionsGetReal(NULL,NULL,"-b",&user.b,&flg));
75   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_getdiagonal",&testgetdiag,NULL));
76 
77   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n"));
78   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mx: %" PetscInt_FMT ",  my: %" PetscInt_FMT ",  ecc: %g \n\n",user.nx,user.ny,(double)user.ecc));
79 
80   /* Let Petsc determine the grid division */
81   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
82 
83   /*
84      A two dimensional distributed array will help define this problem,
85      which derives from an elliptic PDE on two dimensional domain.  From
86      the distributed array, Create the vectors.
87   */
88   PetscCall(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));
89   PetscCall(DMSetFromOptions(user.dm));
90   PetscCall(DMSetUp(user.dm));
91 
92   /*
93      Extract global and local vectors from DM; the vector user.B is
94      used solely as work space for the evaluation of the function,
95      gradient, and Hessian.  Duplicate for remaining vectors that are
96      the same types.
97   */
98   PetscCall(DMCreateGlobalVector(user.dm,&x)); /* Solution */
99   PetscCall(VecDuplicate(x,&user.B)); /* Linear objective */
100 
101   /*  Create matrix user.A to store quadratic, Create a local ordering scheme. */
102   PetscCall(VecGetLocalSize(x,&m));
103   PetscCall(DMCreateMatrix(user.dm,&user.A));
104 
105   if (testgetdiag) PetscCall(MatSetOperation(user.A,MATOP_GET_DIAGONAL,NULL));
106 
107   /* User defined function -- compute linear term of quadratic */
108   PetscCall(ComputeB(&user));
109 
110   /* The TAO code begins here */
111 
112   /*
113      Create the optimization solver
114      Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
115   */
116   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
117   PetscCall(TaoSetType(tao,TAOBLMVM));
118 
119   /* Set the initial vector */
120   PetscCall(VecSet(x, zero));
121   PetscCall(TaoSetSolution(tao,x));
122 
123   /* Set the user function, gradient, hessian evaluation routines and data structures */
124   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user));
125 
126   PetscCall(TaoSetHessian(tao,user.A,user.A,FormHessian,(void*)&user));
127 
128   /* Set a routine that defines the bounds */
129   PetscCall(VecDuplicate(x,&xl));
130   PetscCall(VecDuplicate(x,&xu));
131   PetscCall(VecSet(xl, zero));
132   PetscCall(VecSet(xu, d1000));
133   PetscCall(TaoSetVariableBounds(tao,xl,xu));
134 
135   PetscCall(TaoGetKSP(tao,&ksp));
136   if (ksp) PetscCall(KSPSetType(ksp,KSPCG));
137 
138   PetscCall(PetscOptionsHasName(NULL,NULL,"-testmonitor",&flg));
139   if (flg) {
140     PetscCall(TaoSetMonitor(tao,Monitor,&user,NULL));
141   }
142   PetscCall(PetscOptionsHasName(NULL,NULL,"-testconvergence",&flg));
143   if (flg) {
144     PetscCall(TaoSetConvergenceTest(tao,ConvergenceTest,&user));
145   }
146 
147   /* Check for any tao command line options */
148   PetscCall(TaoSetFromOptions(tao));
149 
150   /* Solve the bound constrained problem */
151   PetscCall(TaoSolve(tao));
152 
153   /* Free PETSc data structures */
154   PetscCall(VecDestroy(&x));
155   PetscCall(VecDestroy(&xl));
156   PetscCall(VecDestroy(&xu));
157   PetscCall(MatDestroy(&user.A));
158   PetscCall(VecDestroy(&user.B));
159 
160   /* Free TAO data structures */
161   PetscCall(TaoDestroy(&tao));
162   PetscCall(DMDestroy(&user.dm));
163   PetscCall(PetscFinalize());
164   return 0;
165 }
166 
167 static PetscReal p(PetscReal xi, PetscReal ecc)
168 {
169   PetscReal t=1.0+ecc*PetscCosScalar(xi);
170   return (t*t*t);
171 }
172 
173 PetscErrorCode ComputeB(AppCtx* user)
174 {
175   PetscInt       i,j,k;
176   PetscInt       nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
177   PetscReal      two=2.0, pi=4.0*atan(1.0);
178   PetscReal      hx,hy,ehxhy;
179   PetscReal      temp,*b;
180   PetscReal      ecc=user->ecc;
181 
182   PetscFunctionBegin;
183   nx=user->nx;
184   ny=user->ny;
185   hx=two*pi/(nx+1.0);
186   hy=two*user->b/(ny+1.0);
187   ehxhy = ecc*hx*hy;
188 
189   /*
190      Get local grid boundaries
191   */
192   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
193   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
194 
195   /* Compute the linear term in the objective function */
196   PetscCall(VecGetArray(user->B,&b));
197   for (i=xs; i<xs+xm; i++) {
198     temp=PetscSinScalar((i+1)*hx);
199     for (j=ys; j<ys+ym; j++) {
200       k=xm*(j-ys)+(i-xs);
201       b[k]=  - ehxhy*temp;
202     }
203   }
204   PetscCall(VecRestoreArray(user->B,&b));
205   PetscCall(PetscLogFlops(5.0*xm*ym+3.0*xm));
206   PetscFunctionReturn(0);
207 }
208 
209 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr)
210 {
211   AppCtx*        user=(AppCtx*)ptr;
212   PetscInt       i,j,k,kk;
213   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
214   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
215   PetscReal      hx,hy,hxhy,hxhx,hyhy;
216   PetscReal      xi,v[5];
217   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
218   PetscReal      vmiddle, vup, vdown, vleft, vright;
219   PetscReal      tt,f1,f2;
220   PetscReal      *x,*g,zero=0.0;
221   Vec            localX;
222 
223   PetscFunctionBegin;
224   nx=user->nx;
225   ny=user->ny;
226   hx=two*pi/(nx+1.0);
227   hy=two*user->b/(ny+1.0);
228   hxhy=hx*hy;
229   hxhx=one/(hx*hx);
230   hyhy=one/(hy*hy);
231 
232   PetscCall(DMGetLocalVector(user->dm,&localX));
233 
234   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
235   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
236 
237   PetscCall(VecSet(G, zero));
238   /*
239     Get local grid boundaries
240   */
241   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
242   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
243 
244   PetscCall(VecGetArray(localX,&x));
245   PetscCall(VecGetArray(G,&g));
246 
247   for (i=xs; i< xs+xm; i++) {
248     xi=(i+1)*hx;
249     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
250     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
251     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
252     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
253     trule5=trule1; /* L(i,j-1) */
254     trule6=trule2; /* U(i,j+1) */
255 
256     vdown=-(trule5+trule2)*hyhy;
257     vleft=-hxhx*(trule2+trule4);
258     vright= -hxhx*(trule1+trule3);
259     vup=-hyhy*(trule1+trule6);
260     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
261 
262     for (j=ys; j<ys+ym; j++) {
263 
264       row=(j-gys)*gxm + (i-gxs);
265        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
266 
267        k=0;
268        if (j>gys) {
269          v[k]=vdown; col[k]=row - gxm; k++;
270        }
271 
272        if (i>gxs) {
273          v[k]= vleft; col[k]=row - 1; k++;
274        }
275 
276        v[k]= vmiddle; col[k]=row; k++;
277 
278        if (i+1 < gxs+gxm) {
279          v[k]= vright; col[k]=row+1; k++;
280        }
281 
282        if (j+1 <gys+gym) {
283          v[k]= vup; col[k] = row+gxm; k++;
284        }
285        tt=0;
286        for (kk=0;kk<k;kk++) {
287          tt+=v[kk]*x[col[kk]];
288        }
289        row=(j-ys)*xm + (i-xs);
290        g[row]=tt;
291 
292      }
293 
294   }
295 
296   PetscCall(VecRestoreArray(localX,&x));
297   PetscCall(VecRestoreArray(G,&g));
298 
299   PetscCall(DMRestoreLocalVector(user->dm,&localX));
300 
301   PetscCall(VecDot(X,G,&f1));
302   PetscCall(VecDot(user->B,X,&f2));
303   PetscCall(VecAXPY(G, one, user->B));
304   *fcn = f1/2.0 + f2;
305 
306   PetscCall(PetscLogFlops((91 + 10.0*ym) * xm));
307   PetscFunctionReturn(0);
308 
309 }
310 
311 /*
312    FormHessian computes the quadratic term in the quadratic objective function
313    Notice that the objective function in this problem is quadratic (therefore a constant
314    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
315 */
316 PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
317 {
318   AppCtx*        user=(AppCtx*)ptr;
319   PetscInt       i,j,k;
320   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
321   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
322   PetscReal      hx,hy,hxhy,hxhx,hyhy;
323   PetscReal      xi,v[5];
324   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
325   PetscReal      vmiddle, vup, vdown, vleft, vright;
326   PetscBool      assembled;
327 
328   PetscFunctionBegin;
329   nx=user->nx;
330   ny=user->ny;
331   hx=two*pi/(nx+1.0);
332   hy=two*user->b/(ny+1.0);
333   hxhy=hx*hy;
334   hxhx=one/(hx*hx);
335   hyhy=one/(hy*hy);
336 
337   /*
338     Get local grid boundaries
339   */
340   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
341   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
342   PetscCall(MatAssembled(hes,&assembled));
343   if (assembled) PetscCall(MatZeroEntries(hes));
344 
345   for (i=xs; i< xs+xm; i++) {
346     xi=(i+1)*hx;
347     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
348     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
349     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
350     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
351     trule5=trule1; /* L(i,j-1) */
352     trule6=trule2; /* U(i,j+1) */
353 
354     vdown=-(trule5+trule2)*hyhy;
355     vleft=-hxhx*(trule2+trule4);
356     vright= -hxhx*(trule1+trule3);
357     vup=-hyhy*(trule1+trule6);
358     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
359     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
360 
361     for (j=ys; j<ys+ym; j++) {
362       row=(j-gys)*gxm + (i-gxs);
363 
364       k=0;
365       if (j>gys) {
366         v[k]=vdown; col[k]=row - gxm; k++;
367       }
368 
369       if (i>gxs) {
370         v[k]= vleft; col[k]=row - 1; k++;
371       }
372 
373       v[k]= vmiddle; col[k]=row; k++;
374 
375       if (i+1 < gxs+gxm) {
376         v[k]= vright; col[k]=row+1; k++;
377       }
378 
379       if (j+1 <gys+gym) {
380         v[k]= vup; col[k] = row+gxm; k++;
381       }
382       PetscCall(MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES));
383 
384     }
385 
386   }
387 
388   /*
389      Assemble matrix, using the 2-step process:
390      MatAssemblyBegin(), MatAssemblyEnd().
391      By placing code between these two statements, computations can be
392      done while messages are in transition.
393   */
394   PetscCall(MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY));
395   PetscCall(MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY));
396 
397   /*
398     Tell the matrix we will never add a new nonzero location to the
399     matrix. If we do it will generate an error.
400   */
401   PetscCall(MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
402   PetscCall(MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE));
403 
404   PetscCall(PetscLogFlops(9.0*xm*ym+49.0*xm));
405   PetscFunctionReturn(0);
406 }
407 
408 PetscErrorCode Monitor(Tao tao, void *ctx)
409 {
410   PetscInt           its;
411   PetscReal          f,gnorm,cnorm,xdiff;
412   TaoConvergedReason reason;
413 
414   PetscFunctionBegin;
415   PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
416   if (!(its%5)) {
417     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iteration=%" PetscInt_FMT "\tf=%g\n",its,(double)f));
418   }
419   PetscFunctionReturn(0);
420 }
421 
422 PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
423 {
424   PetscInt           its;
425   PetscReal          f,gnorm,cnorm,xdiff;
426   TaoConvergedReason reason;
427 
428   PetscFunctionBegin;
429   PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
430   if (its == 100) {
431     PetscCall(TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS));
432   }
433   PetscFunctionReturn(0);
434 
435 }
436 
437 /*TEST
438 
439    build:
440       requires: !complex
441 
442    test:
443       args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5
444       requires: !single
445 
446    test:
447       suffix: 2
448       nsize: 2
449       args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5
450       requires: !single
451 
452    test:
453       suffix: 3
454       nsize: 2
455       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4
456       requires: !single
457 
458    test:
459       suffix: 4
460       nsize: 2
461       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal
462       output_file: output/jbearing2_3.out
463       requires: !single
464 
465    test:
466       suffix: 5
467       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
468       requires: !single
469 
470    test:
471       suffix: 6
472       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4
473       requires: !single
474 
475    test:
476       suffix: 7
477       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5
478       requires: !single
479 
480    test:
481       suffix: 8
482       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5
483       requires: !single
484 
485    test:
486       suffix: 9
487       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5
488       requires: !single
489 
490    test:
491       suffix: 10
492       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
493       requires: !single
494 
495    test:
496       suffix: 11
497       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
498       requires: !single
499 
500    test:
501       suffix: 12
502       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
503       requires: !single
504 
505    test:
506      suffix: 13
507      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls
508      requires: !single
509 
510    test:
511      suffix: 14
512      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm
513      requires: !single
514 
515    test:
516      suffix: 15
517      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
518      requires: !single
519 
520    test:
521      suffix: 16
522      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
523      requires: !single
524 
525    test:
526      suffix: 17
527      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view
528      requires: !single
529 
530    test:
531      suffix: 18
532      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view
533      requires: !single
534 
535    test:
536      suffix: 19
537      args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
538      requires: !single
539 
540    test:
541       suffix: 20
542       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
543       requires: !single
544 
545    test:
546       suffix: 21
547       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
548       requires: !single
549 TEST*/
550