Actual source code: minsurf1.c
petsc-dev 2014-02-02
1: #include <petsctao.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: TaoCreate(); TaoDestroy();
24: Processors: 1
25: T*/
28: /*
29: User-defined application context - contains data needed by the
30: application-provided call-back routines, FormFunctionGradient(),
31: FormHessian().
32: */
33: typedef struct {
34: PetscInt mx, my;
35: PetscReal *bottom, *top, *left, *right;
36: } AppCtx;
39: /* -------- User-defined Routines --------- */
41: static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
42: static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
43: PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
44: PetscErrorCode FormJacobian(Tao, Vec, Mat *, Mat*, MatStructure*,void *);
48: int main(int argc, char **argv)
49: {
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 */
62: /* Initialize PETSc, TAO */
63: PetscInitialize(&argc, &argv, (char *)0, help );
65: /* Specify default dimension of the problem */
66: user.mx = 4; user.my = 4;
68: /* Check for any command line arguments that override defaults */
69: PetscOptionsGetInt(NULL, "-mx", &user.mx, &flg);
70: PetscOptionsGetInt(NULL, "-my", &user.my, &flg);
72: /* Calculate any derived values from parameters */
73: N = user.mx*user.my;
75: PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n");
76: PetscPrintf(PETSC_COMM_SELF,"mx:%d, my:%d\n", user.mx,user.my);
78: /* Create appropriate vectors and matrices */
79: VecCreateSeq(MPI_COMM_SELF, N, &x);
80: VecDuplicate(x, &c);
81: MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J);
83: /* The TAO code begins here */
85: /* Create TAO solver and set desired solution method */
86: TaoCreate(PETSC_COMM_SELF,&tao);
87: TaoSetType(tao,"tao_ssils");
89: /* Set data structure */
90: TaoSetInitialVector(tao, x);
92: /* Set routines for constraints function and Jacobian evaluation */
93: TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user);
94: TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user);
96: /* Set the variable bounds */
97: MSA_BoundaryConditions(&user);
99: /* Set initial solution guess */
100: MSA_InitialPoint(&user, x);
102: /* Set Bounds on variables */
103: VecDuplicate(x, &xl);
104: VecDuplicate(x, &xu);
105: VecSet(xl, lb);
106: VecSet(xu, ub);
107: TaoSetVariableBounds(tao,xl,xu);
109: /* Check for any tao command line options */
110: TaoSetFromOptions(tao);
112: /* Solve the application */
113: TaoSolve(tao);
115: /* Free Tao data structures */
116: TaoDestroy(&tao);
118: /* Free PETSc data structures */
119: VecDestroy(&x);
120: VecDestroy(&xl);
121: VecDestroy(&xu);
122: VecDestroy(&c);
123: MatDestroy(&J);
125: /* Free user-created data structures */
126: PetscFree(user.bottom);
127: PetscFree(user.top);
128: PetscFree(user.left);
129: PetscFree(user.right);
131: PetscFinalize();
132: return 0;
133: }
135: /* -------------------------------------------------------------------- */
139: /* FormConstraints - Evaluates gradient of f.
141: Input Parameters:
142: . tao - the TAO_APPLICATION context
143: . X - input vector
144: . ptr - optional user-defined context, as set by TaoSetConstraintsRoutine()
146: Output Parameters:
147: . G - vector containing the newly evaluated gradient
148: */
149: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr)
150: {
151: AppCtx *user = (AppCtx *) ptr;
153: PetscInt i,j,row;
154: PetscInt mx=user->mx, my=user->my;
155: PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
156: PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
157: PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
158: PetscScalar zero=0.0;
159: PetscScalar *g, *x;
162: /* Initialize vector to zero */
163: VecSet(G, zero);
165: /* Get pointers to vector data */
166: VecGetArray(X, &x);
167: VecGetArray(G, &g);
169: /* Compute function over the locally owned part of the mesh */
170: for (j=0; j<my; j++){
171: for (i=0; i< mx; i++){
172: row= j*mx + i;
174: xc = x[row];
175: xlt=xrb=xl=xr=xb=xt=xc;
177: if (i==0){ /* left side */
178: xl= user->left[j+1];
179: xlt = user->left[j+2];
180: } else {
181: xl = x[row-1];
182: }
184: if (j==0){ /* bottom side */
185: xb=user->bottom[i+1];
186: xrb = user->bottom[i+2];
187: } else {
188: xb = x[row-mx];
189: }
191: if (i+1 == mx){ /* right side */
192: xr=user->right[j+1];
193: xrb = user->right[j];
194: } else {
195: xr = x[row+1];
196: }
198: if (j+1==0+my){ /* top side */
199: xt=user->top[i+1];
200: xlt = user->top[i];
201: }else {
202: xt = x[row+mx];
203: }
205: if (i>0 && j+1<my){
206: xlt = x[row-1+mx];
207: }
208: if (j>0 && i+1<mx){
209: xrb = x[row+1-mx];
210: }
212: d1 = (xc-xl);
213: d2 = (xc-xr);
214: d3 = (xc-xt);
215: d4 = (xc-xb);
216: d5 = (xr-xrb);
217: d6 = (xrb-xb);
218: d7 = (xlt-xl);
219: d8 = (xt-xlt);
221: df1dxc = d1*hydhx;
222: df2dxc = ( d1*hydhx + d4*hxdhy );
223: df3dxc = d3*hxdhy;
224: df4dxc = ( d2*hydhx + d3*hxdhy );
225: df5dxc = d2*hydhx;
226: df6dxc = d4*hxdhy;
228: d1 /= hx;
229: d2 /= hx;
230: d3 /= hy;
231: d4 /= hy;
232: d5 /= hy;
233: d6 /= hx;
234: d7 /= hy;
235: d8 /= hx;
237: f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7);
238: f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
239: f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8);
240: f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
241: f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5);
242: f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6);
244: df1dxc /= f1;
245: df2dxc /= f2;
246: df3dxc /= f3;
247: df4dxc /= f4;
248: df5dxc /= f5;
249: df6dxc /= f6;
251: g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc )/2.0;
252: }
253: }
255: /* Restore vectors */
256: VecRestoreArray(X, &x);
257: VecRestoreArray(G, &g);
258: PetscLogFlops(67*mx*my);
259: return(0);
260: }
262: /* ------------------------------------------------------------------- */
265: /*
266: FormJacobian - Evaluates Jacobian matrix.
268: Input Parameters:
269: . tao - the TAO_APPLICATION context
270: . X - input vector
271: . ptr - optional user-defined context, as set by TaoSetJacobian()
273: Output Parameters:
274: . tH - Jacobian matrix
276: */
277: PetscErrorCode FormJacobian(Tao tao, Vec X, Mat *tH, Mat* tHPre, MatStructure* flag, void *ptr)
278: {
279: AppCtx *user = (AppCtx *) ptr;
280: Mat H = *tH;
282: PetscInt i,j,k,row;
283: PetscInt mx=user->mx, my=user->my;
284: PetscInt col[7];
285: PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
286: PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
287: PetscReal hl,hr,ht,hb,hc,htl,hbr;
288: PetscScalar *x, v[7];
289: PetscBool assembled;
291: /* Set various matrix options */
292: MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
293: MatAssembled(H,&assembled);
294: if (assembled){MatZeroEntries(H); }
295: *flag=SAME_NONZERO_PATTERN;
297: /* Get pointers to vector data */
298: VecGetArray(X, &x);
300: /* Compute Jacobian over the locally owned part of the mesh */
301: for (i=0; i< mx; i++){
302: for (j=0; j<my; j++){
303: row= j*mx + i;
305: xc = x[row];
306: xlt=xrb=xl=xr=xb=xt=xc;
308: /* Left side */
309: if (i==0){
310: xl= user->left[j+1];
311: xlt = user->left[j+2];
312: } else {
313: xl = x[row-1];
314: }
316: if (j==0){
317: xb=user->bottom[i+1];
318: xrb = user->bottom[i+2];
319: } else {
320: xb = x[row-mx];
321: }
323: if (i+1 == mx){
324: xr=user->right[j+1];
325: xrb = user->right[j];
326: } else {
327: xr = x[row+1];
328: }
330: if (j+1==my){
331: xt=user->top[i+1];
332: xlt = user->top[i];
333: }else {
334: xt = x[row+mx];
335: }
337: if (i>0 && j+1<my){
338: xlt = x[row-1+mx];
339: }
340: if (j>0 && i+1<mx){
341: xrb = x[row+1-mx];
342: }
345: d1 = (xc-xl)/hx;
346: d2 = (xc-xr)/hx;
347: d3 = (xc-xt)/hy;
348: d4 = (xc-xb)/hy;
349: d5 = (xrb-xr)/hy;
350: d6 = (xrb-xb)/hx;
351: d7 = (xlt-xl)/hy;
352: d8 = (xlt-xt)/hx;
354: f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7);
355: f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
356: f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8);
357: f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
358: f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5);
359: f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6);
362: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
363: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
364: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
365: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
367: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
368: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
370: 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) +
371: (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);
373: hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0;
375: k=0;
376: if (j>0){
377: v[k]=hb; col[k]=row - mx; k++;
378: }
380: if (j>0 && i < mx -1){
381: v[k]=hbr; col[k]=row - mx+1; k++;
382: }
384: if (i>0){
385: v[k]= hl; col[k]=row - 1; k++;
386: }
388: v[k]= hc; col[k]=row; k++;
390: if (i < mx-1 ){
391: v[k]= hr; col[k]=row+1; k++;
392: }
394: if (i>0 && j < my-1 ){
395: v[k]= htl; col[k] = row+mx-1; k++;
396: }
398: if (j < my-1 ){
399: v[k]= ht; col[k] = row+mx; k++;
400: }
402: /*
403: Set matrix values using local numbering, which was defined
404: earlier, in the main routine.
405: */
406: MatSetValues(H,1,&row,k,col,v,INSERT_VALUES);
407: }
408: }
410: /* Restore vectors */
411: VecRestoreArray(X,&x);
413: /* Assemble the matrix */
414: MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
415: MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
416: PetscLogFlops(199*mx*my);
417: return(0);
418: }
420: /* ------------------------------------------------------------------- */
423: /*
424: MSA_BoundaryConditions - Calculates the boundary conditions for
425: the region.
427: Input Parameter:
428: . user - user-defined application context
430: Output Parameter:
431: . user - user-defined application context
432: */
433: static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
434: {
435: PetscErrorCode ierr;
436: PetscInt i,j,k,limit=0,maxits=5;
437: PetscInt mx=user->mx,my=user->my;
438: PetscInt bsize=0, lsize=0, tsize=0, rsize=0;
439: PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10;
440: PetscReal fnorm,det,hx,hy,xt=0,yt=0;
441: PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
442: PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5;
443: PetscReal *boundary;
446: bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
448: PetscMalloc1(bsize, &user->bottom);
449: PetscMalloc1(tsize, &user->top);
450: PetscMalloc1(lsize, &user->left);
451: PetscMalloc1(rsize, &user->right);
453: hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
455: for (j=0; j<4; j++){
456: if (j==0){
457: yt=b;
458: xt=l;
459: limit=bsize;
460: boundary=user->bottom;
461: } else if (j==1){
462: yt=t;
463: xt=l;
464: limit=tsize;
465: boundary=user->top;
466: } else if (j==2){
467: yt=b;
468: xt=l;
469: limit=lsize;
470: boundary=user->left;
471: } else { /* if (j==3) */
472: yt=b;
473: xt=r;
474: limit=rsize;
475: boundary=user->right;
476: }
478: for (i=0; i<limit; i++){
479: u1=xt;
480: u2=-yt;
481: for (k=0; k<maxits; k++){
482: nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
483: nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
484: fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
485: if (fnorm <= tol) break;
486: njac11=one+u2*u2-u1*u1;
487: njac12=two*u1*u2;
488: njac21=-two*u1*u2;
489: njac22=-one - u1*u1 + u2*u2;
490: det = njac11*njac22-njac21*njac12;
491: u1 = u1-(njac22*nf1-njac12*nf2)/det;
492: u2 = u2-(njac11*nf2-njac21*nf1)/det;
493: }
495: boundary[i]=u1*u1-u2*u2;
496: if (j==0 || j==1) {
497: xt=xt+hx;
498: } else { /* if (j==2 || j==3) */
499: yt=yt+hy;
500: }
501: }
502: }
503: return(0);
504: }
506: /* ------------------------------------------------------------------- */
509: /*
510: MSA_InitialPoint - Calculates the initial guess in one of three ways.
512: Input Parameters:
513: . user - user-defined application context
514: . X - vector for initial guess
516: Output Parameters:
517: . X - newly computed initial guess
518: */
519: static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
520: {
522: PetscInt start=-1,i,j;
523: PetscScalar zero=0.0;
524: PetscBool flg;
527: PetscOptionsGetInt(NULL,"-start",&start,&flg);
529: if (flg && start==0){ /* The zero vector is reasonable */
530: VecSet(X, zero);
531: } else { /* Take an average of the boundary conditions */
532: PetscInt row;
533: PetscInt mx=user->mx,my=user->my;
534: PetscScalar *x;
536: /* Get pointers to vector data */
537: VecGetArray(X,&x);
539: /* Perform local computations */
540: for (j=0; j<my; j++){
541: for (i=0; i< mx; i++){
542: row=(j)*mx + (i);
543: 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;
544: }
545: }
547: /* Restore vectors */
548: VecRestoreArray(X,&x);
549: }
550: return(0);
551: }