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