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