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