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