xref: /petsc/src/tao/complementarity/tutorials/minsurf1.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 #include <petsctao.h>
2 
3 static char  help[] =
4 "This example demonstrates use of the TAO package to\n\
5 solve an unconstrained system of equations.  This example is based on a\n\
6 problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and\n\
7 boundary values along the edges of the domain, the objective is to find the\n\
8 surface with the minimal area that satisfies the boundary conditions.\n\
9 This application solves this problem using complimentarity -- We are actually\n\
10 solving the system  (grad f)_i >= 0, if x_i == l_i \n\
11                     (grad f)_i = 0, if l_i < x_i < u_i \n\
12                     (grad f)_i <= 0, if x_i == u_i  \n\
13 where f is the function to be minimized. \n\
14 \n\
15 The command line options are:\n\
16   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
17   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
18   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";
19 
20 /*T
21    Concepts: TAO^Solving a complementarity problem
22    Routines: TaoCreate(); TaoDestroy();
23 
24    Processors: 1
25 T*/
26 
27 /*
28    User-defined application context - contains data needed by the
29    application-provided call-back routines, FormFunctionGradient(),
30    FormHessian().
31 */
32 typedef struct {
33   PetscInt  mx, my;
34   PetscReal *bottom, *top, *left, *right;
35 } AppCtx;
36 
37 /* -------- User-defined Routines --------- */
38 
39 static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
40 static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
41 PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
42 PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *);
43 
44 int main(int argc, char **argv)
45 {
46   Vec            x;                 /* solution vector */
47   Vec            c;                 /* Constraints function vector */
48   Vec            xl,xu;             /* Bounds on the variables */
49   PetscBool      flg;               /* A return variable when checking for user options */
50   Tao            tao;               /* TAO solver context */
51   Mat            J;                 /* Jacobian matrix */
52   PetscInt       N;                 /* Number of elements in vector */
53   PetscScalar    lb =  PETSC_NINFINITY;      /* lower bound constant */
54   PetscScalar    ub =  PETSC_INFINITY;      /* upper bound constant */
55   AppCtx         user;                    /* user-defined work context */
56 
57   /* Initialize PETSc, TAO */
58   CHKERRQ(PetscInitialize(&argc, &argv, (char *)0, help));
59 
60   /* Specify default dimension of the problem */
61   user.mx = 4; user.my = 4;
62 
63   /* Check for any command line arguments that override defaults */
64   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-mx", &user.mx, &flg));
65   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-my", &user.my, &flg));
66 
67   /* Calculate any derived values from parameters */
68   N = user.mx*user.my;
69 
70   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n"));
71   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"mx:%D, my:%D\n", user.mx,user.my));
72 
73   /* Create appropriate vectors and matrices */
74   CHKERRQ(VecCreateSeq(MPI_COMM_SELF, N, &x));
75   CHKERRQ(VecDuplicate(x, &c));
76   CHKERRQ(MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J));
77 
78   /* The TAO code begins here */
79 
80   /* Create TAO solver and set desired solution method */
81   CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao));
82   CHKERRQ(TaoSetType(tao,TAOSSILS));
83 
84   /* Set data structure */
85   CHKERRQ(TaoSetSolution(tao, x));
86 
87   /*  Set routines for constraints function and Jacobian evaluation */
88   CHKERRQ(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user));
89   CHKERRQ(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user));
90 
91   /* Set the variable bounds */
92   CHKERRQ(MSA_BoundaryConditions(&user));
93 
94   /* Set initial solution guess */
95   CHKERRQ(MSA_InitialPoint(&user, x));
96 
97   /* Set Bounds on variables */
98   CHKERRQ(VecDuplicate(x, &xl));
99   CHKERRQ(VecDuplicate(x, &xu));
100   CHKERRQ(VecSet(xl, lb));
101   CHKERRQ(VecSet(xu, ub));
102   CHKERRQ(TaoSetVariableBounds(tao,xl,xu));
103 
104   /* Check for any tao command line options */
105   CHKERRQ(TaoSetFromOptions(tao));
106 
107   /* Solve the application */
108   CHKERRQ(TaoSolve(tao));
109 
110   /* Free Tao data structures */
111   CHKERRQ(TaoDestroy(&tao));
112 
113   /* Free PETSc data structures */
114   CHKERRQ(VecDestroy(&x));
115   CHKERRQ(VecDestroy(&xl));
116   CHKERRQ(VecDestroy(&xu));
117   CHKERRQ(VecDestroy(&c));
118   CHKERRQ(MatDestroy(&J));
119 
120   /* Free user-created data structures */
121   CHKERRQ(PetscFree(user.bottom));
122   CHKERRQ(PetscFree(user.top));
123   CHKERRQ(PetscFree(user.left));
124   CHKERRQ(PetscFree(user.right));
125 
126   CHKERRQ(PetscFinalize());
127   return 0;
128 }
129 
130 /* -------------------------------------------------------------------- */
131 
132 /*  FormConstraints - Evaluates gradient of f.
133 
134     Input Parameters:
135 .   tao  - the TAO_APPLICATION context
136 .   X    - input vector
137 .   ptr  - optional user-defined context, as set by TaoSetConstraintsRoutine()
138 
139     Output Parameters:
140 .   G - vector containing the newly evaluated gradient
141 */
142 PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr)
143 {
144   AppCtx         *user = (AppCtx *) ptr;
145   PetscInt       i,j,row;
146   PetscInt       mx=user->mx, my=user->my;
147   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
148   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
149   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
150   PetscScalar    zero=0.0;
151   PetscScalar    *g, *x;
152 
153   PetscFunctionBegin;
154   /* Initialize vector to zero */
155   CHKERRQ(VecSet(G, zero));
156 
157   /* Get pointers to vector data */
158   CHKERRQ(VecGetArray(X, &x));
159   CHKERRQ(VecGetArray(G, &g));
160 
161   /* Compute function over the locally owned part of the mesh */
162   for (j=0; j<my; j++) {
163     for (i=0; i< mx; i++) {
164       row= j*mx + i;
165 
166       xc = x[row];
167       xlt=xrb=xl=xr=xb=xt=xc;
168 
169       if (i==0) { /* left side */
170         xl= user->left[j+1];
171         xlt = user->left[j+2];
172       } else {
173         xl = x[row-1];
174       }
175 
176       if (j==0) { /* bottom side */
177         xb=user->bottom[i+1];
178         xrb = user->bottom[i+2];
179       } else {
180         xb = x[row-mx];
181       }
182 
183       if (i+1 == mx) { /* right side */
184         xr=user->right[j+1];
185         xrb = user->right[j];
186       } else {
187         xr = x[row+1];
188       }
189 
190       if (j+1==0+my) { /* top side */
191         xt=user->top[i+1];
192         xlt = user->top[i];
193       }else {
194         xt = x[row+mx];
195       }
196 
197       if (i>0 && j+1<my) {
198         xlt = x[row-1+mx];
199       }
200       if (j>0 && i+1<mx) {
201         xrb = x[row+1-mx];
202       }
203 
204       d1 = (xc-xl);
205       d2 = (xc-xr);
206       d3 = (xc-xt);
207       d4 = (xc-xb);
208       d5 = (xr-xrb);
209       d6 = (xrb-xb);
210       d7 = (xlt-xl);
211       d8 = (xt-xlt);
212 
213       df1dxc = d1*hydhx;
214       df2dxc = (d1*hydhx + d4*hxdhy);
215       df3dxc = d3*hxdhy;
216       df4dxc = (d2*hydhx + d3*hxdhy);
217       df5dxc = d2*hydhx;
218       df6dxc = d4*hxdhy;
219 
220       d1 /= hx;
221       d2 /= hx;
222       d3 /= hy;
223       d4 /= hy;
224       d5 /= hy;
225       d6 /= hx;
226       d7 /= hy;
227       d8 /= hx;
228 
229       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
230       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
231       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
232       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
233       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
234       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
235 
236       df1dxc /= f1;
237       df2dxc /= f2;
238       df3dxc /= f3;
239       df4dxc /= f4;
240       df5dxc /= f5;
241       df6dxc /= f6;
242 
243       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
244     }
245   }
246 
247   /* Restore vectors */
248   CHKERRQ(VecRestoreArray(X, &x));
249   CHKERRQ(VecRestoreArray(G, &g));
250   CHKERRQ(PetscLogFlops(67*mx*my));
251   PetscFunctionReturn(0);
252 }
253 
254 /* ------------------------------------------------------------------- */
255 /*
256    FormJacobian - Evaluates Jacobian matrix.
257 
258    Input Parameters:
259 .  tao  - the TAO_APPLICATION context
260 .  X    - input vector
261 .  ptr  - optional user-defined context, as set by TaoSetJacobian()
262 
263    Output Parameters:
264 .  tH    - Jacobian matrix
265 
266 */
267 PetscErrorCode FormJacobian(Tao tao, Vec X, Mat H, Mat tHPre, void *ptr)
268 {
269   AppCtx            *user = (AppCtx *) ptr;
270   PetscInt          i,j,k,row;
271   PetscInt          mx=user->mx, my=user->my;
272   PetscInt          col[7];
273   PetscReal         hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
274   PetscReal         f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
275   PetscReal         hl,hr,ht,hb,hc,htl,hbr;
276   const PetscScalar *x;
277   PetscScalar       v[7];
278   PetscBool         assembled;
279 
280   /* Set various matrix options */
281   PetscFunctionBegin;
282   CHKERRQ(MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
283   CHKERRQ(MatAssembled(H,&assembled));
284   if (assembled) CHKERRQ(MatZeroEntries(H));
285 
286   /* Get pointers to vector data */
287   CHKERRQ(VecGetArrayRead(X, &x));
288 
289   /* Compute Jacobian over the locally owned part of the mesh */
290   for (i=0; i< mx; i++) {
291     for (j=0; j<my; j++) {
292       row= j*mx + i;
293 
294       xc = x[row];
295       xlt=xrb=xl=xr=xb=xt=xc;
296 
297       /* Left side */
298       if (i==0) {
299         xl  = user->left[j+1];
300         xlt = user->left[j+2];
301       } else {
302         xl = x[row-1];
303       }
304 
305       if (j==0) {
306         xb  = user->bottom[i+1];
307         xrb = user->bottom[i+2];
308       } else {
309         xb = x[row-mx];
310       }
311 
312       if (i+1 == mx) {
313         xr  = user->right[j+1];
314         xrb = user->right[j];
315       } else {
316         xr = x[row+1];
317       }
318 
319       if (j+1==my) {
320         xt  = user->top[i+1];
321         xlt = user->top[i];
322       }else {
323         xt = x[row+mx];
324       }
325 
326       if (i>0 && j+1<my) {
327         xlt = x[row-1+mx];
328       }
329       if (j>0 && i+1<mx) {
330         xrb = x[row+1-mx];
331       }
332 
333       d1 = (xc-xl)/hx;
334       d2 = (xc-xr)/hx;
335       d3 = (xc-xt)/hy;
336       d4 = (xc-xb)/hy;
337       d5 = (xrb-xr)/hy;
338       d6 = (xrb-xb)/hx;
339       d7 = (xlt-xl)/hy;
340       d8 = (xlt-xt)/hx;
341 
342       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
343       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
344       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
345       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
346       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
347       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
348 
349       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
350       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
351       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
352       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
353 
354       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
355       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
356 
357       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
358            (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) + (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
359 
360       hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0;  hc/=2.0;
361 
362       k=0;
363       if (j>0) {
364         v[k]=hb; col[k]=row - mx; k++;
365       }
366 
367       if (j>0 && i < mx -1) {
368         v[k]=hbr; col[k]=row - mx+1; k++;
369       }
370 
371       if (i>0) {
372         v[k]= hl; col[k]=row - 1; k++;
373       }
374 
375       v[k]= hc; col[k]=row; k++;
376 
377       if (i < mx-1) {
378         v[k]= hr; col[k]=row+1; k++;
379       }
380 
381       if (i>0 && j < my-1) {
382         v[k]= htl; col[k] = row+mx-1; k++;
383       }
384 
385       if (j < my-1) {
386         v[k]= ht; col[k] = row+mx; k++;
387       }
388 
389       /*
390          Set matrix values using local numbering, which was defined
391          earlier, in the main routine.
392       */
393       CHKERRQ(MatSetValues(H,1,&row,k,col,v,INSERT_VALUES));
394     }
395   }
396 
397   /* Restore vectors */
398   CHKERRQ(VecRestoreArrayRead(X,&x));
399 
400   /* Assemble the matrix */
401   CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
402   CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
403   CHKERRQ(PetscLogFlops(199*mx*my));
404   PetscFunctionReturn(0);
405 }
406 
407 /* ------------------------------------------------------------------- */
408 /*
409    MSA_BoundaryConditions -  Calculates the boundary conditions for
410    the region.
411 
412    Input Parameter:
413 .  user - user-defined application context
414 
415    Output Parameter:
416 .  user - user-defined application context
417 */
418 static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
419 {
420   PetscInt        i,j,k,limit=0,maxits=5;
421   PetscInt        mx=user->mx,my=user->my;
422   PetscInt        bsize=0, lsize=0, tsize=0, rsize=0;
423   PetscReal       one=1.0, two=2.0, three=3.0, tol=1e-10;
424   PetscReal       fnorm,det,hx,hy,xt=0,yt=0;
425   PetscReal       u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
426   PetscReal       b=-0.5, t=0.5, l=-0.5, r=0.5;
427   PetscReal       *boundary;
428 
429   PetscFunctionBegin;
430   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
431 
432   CHKERRQ(PetscMalloc1(bsize, &user->bottom));
433   CHKERRQ(PetscMalloc1(tsize, &user->top));
434   CHKERRQ(PetscMalloc1(lsize, &user->left));
435   CHKERRQ(PetscMalloc1(rsize, &user->right));
436 
437   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
438 
439   for (j=0; j<4; j++) {
440     if (j==0) {
441       yt=b;
442       xt=l;
443       limit=bsize;
444       boundary=user->bottom;
445     } else if (j==1) {
446       yt=t;
447       xt=l;
448       limit=tsize;
449       boundary=user->top;
450     } else if (j==2) {
451       yt=b;
452       xt=l;
453       limit=lsize;
454       boundary=user->left;
455     } else { /* if  (j==3) */
456       yt=b;
457       xt=r;
458       limit=rsize;
459       boundary=user->right;
460     }
461 
462     for (i=0; i<limit; i++) {
463       u1=xt;
464       u2=-yt;
465       for (k=0; k<maxits; k++) {
466         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
467         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
468         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
469         if (fnorm <= tol) break;
470         njac11=one+u2*u2-u1*u1;
471         njac12=two*u1*u2;
472         njac21=-two*u1*u2;
473         njac22=-one - u1*u1 + u2*u2;
474         det = njac11*njac22-njac21*njac12;
475         u1 = u1-(njac22*nf1-njac12*nf2)/det;
476         u2 = u2-(njac11*nf2-njac21*nf1)/det;
477       }
478 
479       boundary[i]=u1*u1-u2*u2;
480       if (j==0 || j==1) {
481         xt=xt+hx;
482       } else { /* if (j==2 || j==3) */
483         yt=yt+hy;
484       }
485     }
486   }
487   PetscFunctionReturn(0);
488 }
489 
490 /* ------------------------------------------------------------------- */
491 /*
492    MSA_InitialPoint - Calculates the initial guess in one of three ways.
493 
494    Input Parameters:
495 .  user - user-defined application context
496 .  X - vector for initial guess
497 
498    Output Parameters:
499 .  X - newly computed initial guess
500 */
501 static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
502 {
503   PetscInt       start=-1,i,j;
504   PetscScalar    zero=0.0;
505   PetscBool      flg;
506 
507   PetscFunctionBegin;
508   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
509 
510   if (flg && start==0) { /* The zero vector is reasonable */
511     CHKERRQ(VecSet(X, zero));
512   } else { /* Take an average of the boundary conditions */
513     PetscInt    row;
514     PetscInt    mx=user->mx,my=user->my;
515     PetscScalar *x;
516 
517     /* Get pointers to vector data */
518     CHKERRQ(VecGetArray(X,&x));
519 
520     /* Perform local computations */
521     for (j=0; j<my; j++) {
522       for (i=0; i< mx; i++) {
523         row=(j)*mx + (i);
524         x[row] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+ ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0;
525       }
526     }
527 
528     /* Restore vectors */
529     CHKERRQ(VecRestoreArray(X,&x));
530   }
531   PetscFunctionReturn(0);
532 }
533 
534 /*TEST
535 
536    build:
537       requires: !complex
538 
539    test:
540       args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5
541       requires: !single
542 
543    test:
544       suffix: 2
545       args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5
546 
547 TEST*/
548