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