Actual source code: ex8.c
petsc-3.6.4 2016-04-12
1: #include <petscsnes.h>
2: #include <petscdm.h>
3: #include <petscdmda.h>
5: static char help[] = "Parallel version of the minimum surface area problem using DMs.\n\
6: See ex10.c for the serial version. It solves a system of nonlinear equations in mixed\n\
7: complementarity form using semismooth newton algorithm.This example is based on a\n\
8: problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and\n\
9: boundary values along the edges of the domain, the objective is to find the\n\
10: surface with the minimal area that satisfies the boundary conditions.\n\
11: This application solves this problem using complimentarity -- We are actually\n\
12: solving the system (grad f)_i >= 0, if x_i == l_i \n\
13: (grad f)_i = 0, if l_i < x_i < u_i \n\
14: (grad f)_i <= 0, if x_i == u_i \n\
15: where f is the function to be minimized. \n\
16: \n\
17: The command line options are:\n\
18: -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
19: -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
20: -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise\n\
21: -lb <value>, lower bound on the variables\n\
22: -ub <value>, upper bound on the variables\n\n";
24: /*
25: User-defined application context - contains data needed by the
26: application-provided call-back routines, FormJacobian() and
27: FormFunction().
28: */
30: typedef struct {
31: DM da;
32: PetscScalar *bottom, *top, *left, *right;
33: PetscInt mx,my;
34: } AppCtx;
37: /* -------- User-defined Routines --------- */
39: extern PetscErrorCode MSA_BoundaryConditions(AppCtx*);
40: extern PetscErrorCode MSA_InitialPoint(AppCtx*, Vec);
41: extern PetscErrorCode FormGradient(SNES, Vec, Vec, void*);
42: extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void*);
46: int main(int argc, char **argv)
47: {
49: Vec x,r; /* solution and residual vectors */
50: Vec xl,xu; /* Bounds on the variables */
51: PetscBool flg_l,flg_u; /* flags to check if the bounds are set */
52: SNES snes; /* nonlinear solver context */
53: Mat J; /* Jacobian matrix */
54: PetscInt N; /* Number of elements in vector */
55: PetscScalar lb = .05;
56: PetscScalar ub = PETSC_INFINITY;
57: AppCtx user; /* user-defined work context */
58: PetscBool flg;
60: /* Initialize PETSc */
61: PetscInitialize(&argc, &argv, (char*)0, help);
63: #if defined(PETSC_USE_COMPLEX)
64: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example does not work for scalar type complex\n");
65: #endif
67: /* Check if lower and upper bounds are set */
68: PetscOptionsGetScalar(NULL, "-lb", &lb, &flg_l);
69: PetscOptionsGetScalar(NULL, "-ub", &ub, &flg_u);
71: /* Create distributed array to manage the 2d grid */
72: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.da);
73: DMDAGetIerr(user.da,PETSC_IGNORE,&user.mx,&user.my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
74: /* Extract global vectors from DMDA; */
75: DMCreateGlobalVector(user.da,&x);
76: VecDuplicate(x, &r);
78: N = user.mx*user.my;
79: DMSetMatType(user.da,MATAIJ);
80: DMCreateMatrix(user.da,&J);
82: /* Create nonlinear solver context */
83: SNESCreate(PETSC_COMM_WORLD,&snes);
85: /* Set function evaluation and Jacobian evaluation routines */
86: SNESSetFunction(snes,r,FormGradient,&user);
87: SNESSetJacobian(snes,J,J,FormJacobian,&user);
89: /* Set the boundary conditions */
90: MSA_BoundaryConditions(&user);
92: /* Set initial solution guess */
93: MSA_InitialPoint(&user, x);
96: /* Set Bounds on variables */
97: VecDuplicate(x, &xl);
98: VecDuplicate(x, &xu);
99: VecSet(xl, lb);
100: VecSet(xu, ub);
102: SNESVISetVariableBounds(snes,xl,xu);
104: SNESSetFromOptions(snes);
106: /* Solve the application */
107: SNESSolve(snes,NULL,x);
109: PetscOptionsHasName(NULL,"-view_sol",&flg);
110: if (flg) { VecView(x,PETSC_VIEWER_STDOUT_WORLD); }
112: /* Free memory */
113: VecDestroy(&x);
114: VecDestroy(&xl);
115: VecDestroy(&xu);
116: VecDestroy(&r);
117: MatDestroy(&J);
118: SNESDestroy(&snes);
120: /* Free user-created data structures */
121: DMDestroy(&user.da);
122: PetscFree(user.bottom);
123: PetscFree(user.top);
124: PetscFree(user.left);
125: PetscFree(user.right);
127: PetscFinalize();
129: return 0;
130: }
132: /* -------------------------------------------------------------------- */
136: /* FormGradient - Evaluates gradient of f.
138: Input Parameters:
139: . snes - the SNES context
140: . X - input vector
141: . ptr - optional user-defined context, as set by SNESSetFunction()
143: Output Parameters:
144: . G - vector containing the newly evaluated gradient
145: */
146: PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
147: {
148: AppCtx *user = (AppCtx*) ptr;
149: int ierr;
150: PetscInt i,j;
151: PetscInt mx=user->mx, my=user->my;
152: PetscScalar hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
153: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
154: PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
155: PetscScalar **g, **x;
156: PetscInt xs,xm,ys,ym;
157: Vec localX;
160: /* Initialize vector to zero */
161: VecSet(G,0.0);
163: /* Get local vector */
164: DMGetLocalVector(user->da,&localX);
165: /* Get ghost points */
166: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
167: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
168: /* Get pointer to local vector data */
169: DMDAVecGetArrayRead(user->da,localX, &x);
170: DMDAVecGetArray(user->da,G, &g);
172: DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
173: /* Compute function over the locally owned part of the mesh */
174: for (j=ys; j < ys+ym; j++) {
175: for (i=xs; i< xs+xm; i++) {
177: xc = x[j][i];
178: xlt=xrb=xl=xr=xb=xt=xc;
180: if (i==0) { /* left side */
181: xl = user->left[j+1];
182: xlt = user->left[j+2];
183: } else xl = x[j][i-1];
185: if (j==0) { /* bottom side */
186: xb = user->bottom[i+1];
187: xrb = user->bottom[i+2];
188: } else xb = x[j-1][i];
190: if (i+1 == mx) { /* right side */
191: xr = user->right[j+1];
192: xrb = user->right[j];
193: } else xr = x[j][i+1];
195: if (j+1==0+my) { /* top side */
196: xt = user->top[i+1];
197: xlt = user->top[i];
198: } else xt = x[j+1][i];
200: if (i>0 && j+1<my) xlt = x[j+1][i-1]; /* left top side */
201: if (j>0 && i+1<mx) xrb = x[j-1][i+1]; /* right bottom */
203: d1 = (xc-xl);
204: d2 = (xc-xr);
205: d3 = (xc-xt);
206: d4 = (xc-xb);
207: d5 = (xr-xrb);
208: d6 = (xrb-xb);
209: d7 = (xlt-xl);
210: d8 = (xt-xlt);
212: df1dxc = d1*hydhx;
213: df2dxc = (d1*hydhx + d4*hxdhy);
214: df3dxc = d3*hxdhy;
215: df4dxc = (d2*hydhx + d3*hxdhy);
216: df5dxc = d2*hydhx;
217: df6dxc = d4*hxdhy;
219: d1 /= hx;
220: d2 /= hx;
221: d3 /= hy;
222: d4 /= hy;
223: d5 /= hy;
224: d6 /= hx;
225: d7 /= hy;
226: d8 /= hx;
228: f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
229: f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
230: f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
231: f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
232: f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
233: f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
235: df1dxc /= f1;
236: df2dxc /= f2;
237: df3dxc /= f3;
238: df4dxc /= f4;
239: df5dxc /= f5;
240: df6dxc /= f6;
242: g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
244: }
245: }
247: /* Restore vectors */
248: DMDAVecRestoreArrayRead(user->da,localX, &x);
249: DMDAVecRestoreArray(user->da,G, &g);
250: DMRestoreLocalVector(user->da,&localX);
251: PetscLogFlops(67*mx*my);
252: return(0);
253: }
255: /* ------------------------------------------------------------------- */
258: /*
259: FormJacobian - Evaluates Jacobian matrix.
261: Input Parameters:
262: . snes - SNES context
263: . X - input vector
264: . ptr - optional user-defined context, as set by SNESSetJacobian()
266: Output Parameters:
267: . tH - Jacobian matrix
269: */
270: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *ptr)
271: {
272: AppCtx *user = (AppCtx*) ptr;
274: PetscInt i,j,k;
275: PetscInt mx=user->mx, my=user->my;
276: MatStencil row,col[7];
277: PetscScalar hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
278: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
279: PetscScalar hl,hr,ht,hb,hc,htl,hbr;
280: PetscScalar **x, v[7];
281: PetscBool assembled;
282: PetscInt xs,xm,ys,ym;
283: Vec localX;
286: /* Set various matrix options */
287: MatAssembled(H,&assembled);
288: if (assembled) {MatZeroEntries(H);}
290: /* Get local vector */
291: DMGetLocalVector(user->da,&localX);
292: /* Get ghost points */
293: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
294: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
296: /* Get pointers to vector data */
297: DMDAVecGetArrayRead(user->da,localX, &x);
299: DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
300: /* Compute Jacobian over the locally owned part of the mesh */
301: for (j=ys; j< ys+ym; j++) {
302: for (i=xs; i< xs+xm; i++) {
303: xc = x[j][i];
304: xlt=xrb=xl=xr=xb=xt=xc;
306: /* Left */
307: if (i==0) {
308: xl = user->left[j+1];
309: xlt = user->left[j+2];
310: } else xl = x[j][i-1];
312: /* Bottom */
313: if (j==0) {
314: xb = user->bottom[i+1];
315: xrb = user->bottom[i+2];
316: } else xb = x[j-1][i];
318: /* Right */
319: if (i+1 == mx) {
320: xr = user->right[j+1];
321: xrb = user->right[j];
322: } else xr = x[j][i+1];
324: /* Top */
325: if (j+1==my) {
326: xt = user->top[i+1];
327: xlt = user->top[i];
328: } else xt = x[j+1][i];
330: /* Top left */
331: if (i>0 && j+1<my) xlt = x[j+1][i-1];
333: /* Bottom right */
334: if (j>0 && i+1<mx) xrb = x[j-1][i+1];
336: d1 = (xc-xl)/hx;
337: d2 = (xc-xr)/hx;
338: d3 = (xc-xt)/hy;
339: d4 = (xc-xb)/hy;
340: d5 = (xrb-xr)/hy;
341: d6 = (xrb-xb)/hx;
342: d7 = (xlt-xl)/hy;
343: d8 = (xlt-xt)/hx;
345: f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
346: f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
347: f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
348: f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
349: f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
350: f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
353: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
354: (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
355: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
356: (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
357: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
358: (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
359: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
360: (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
362: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
363: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
365: hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
366: hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
367: (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
368: (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
370: hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0;
372: k =0;
373: row.i = i;row.j= j;
374: /* Bottom */
375: if (j>0) {
376: v[k] =hb;
377: col[k].i = i; col[k].j=j-1; k++;
378: }
380: /* Bottom right */
381: if (j>0 && i < mx -1) {
382: v[k] =hbr;
383: col[k].i = i+1; col[k].j = j-1; k++;
384: }
386: /* left */
387: if (i>0) {
388: v[k] = hl;
389: col[k].i = i-1; col[k].j = j; k++;
390: }
392: /* Centre */
393: v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++;
395: /* Right */
396: if (i < mx-1) {
397: v[k] = hr;
398: col[k].i= i+1; col[k].j = j;k++;
399: }
401: /* Top left */
402: if (i>0 && j < my-1) {
403: v[k] = htl;
404: col[k].i = i-1;col[k].j = j+1; k++;
405: }
407: /* Top */
408: if (j < my-1) {
409: v[k] = ht;
410: col[k].i = i; col[k].j = j+1; k++;
411: }
413: MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES);
414: }
415: }
417: /* Assemble the matrix */
418: MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
419: DMDAVecRestoreArrayRead(user->da,localX,&x);
420: MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
421: DMRestoreLocalVector(user->da,&localX);
423: PetscLogFlops(199*mx*my);
424: return(0);
425: }
427: /* ------------------------------------------------------------------- */
430: /*
431: MSA_BoundaryConditions - Calculates the boundary conditions for
432: the region.
434: Input Parameter:
435: . user - user-defined application context
437: Output Parameter:
438: . user - user-defined application context
439: */
440: PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
441: {
443: PetscInt i,j,k,limit=0,maxits=5;
444: PetscInt mx =user->mx,my=user->my;
445: PetscInt bsize=0, lsize=0, tsize=0, rsize=0;
446: PetscScalar one =1.0, two=2.0, three=3.0, tol=1e-10;
447: PetscScalar fnorm,det,hx,hy,xt=0,yt=0;
448: PetscScalar u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
449: PetscScalar b=-0.5, t=0.5, l=-0.5, r=0.5;
450: PetscScalar *boundary;
453: bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
455: PetscMalloc1(bsize, &user->bottom);
456: PetscMalloc1(tsize, &user->top);
457: PetscMalloc1(lsize, &user->left);
458: PetscMalloc1(rsize, &user->right);
460: hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
462: for (j=0; j<4; j++) {
463: if (j==0) {
464: yt = b;
465: xt = l;
466: limit = bsize;
467: boundary = user->bottom;
468: } else if (j==1) {
469: yt = t;
470: xt = l;
471: limit = tsize;
472: boundary = user->top;
473: } else if (j==2) {
474: yt = b;
475: xt = l;
476: limit = lsize;
477: boundary = user->left;
478: } else { /* if (j==3) */
479: yt = b;
480: xt = r;
481: limit = rsize;
482: boundary = user->right;
483: }
485: for (i=0; i<limit; i++) {
486: u1=xt;
487: u2=-yt;
488: for (k=0; k<maxits; k++) {
489: nf1 = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
490: nf2 = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
491: fnorm = PetscSqrtReal(nf1*nf1+nf2*nf2);
492: if (fnorm <= tol) break;
493: njac11 = one+u2*u2-u1*u1;
494: njac12 = two*u1*u2;
495: njac21 = -two*u1*u2;
496: njac22 = -one - u1*u1 + u2*u2;
497: det = njac11*njac22-njac21*njac12;
498: u1 = u1-(njac22*nf1-njac12*nf2)/det;
499: u2 = u2-(njac11*nf2-njac21*nf1)/det;
500: }
502: boundary[i]=u1*u1-u2*u2;
503: if (j==0 || j==1) xt=xt+hx;
504: else yt=yt+hy; /* if (j==2 || j==3) */
505: }
506: }
507: return(0);
508: }
510: /* ------------------------------------------------------------------- */
513: /*
514: MSA_InitialPoint - Calculates the initial guess in one of three ways.
516: Input Parameters:
517: . user - user-defined application context
518: . X - vector for initial guess
520: Output Parameters:
521: . X - newly computed initial guess
522: */
523: PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
524: {
526: PetscInt start=-1,i,j;
527: PetscScalar zero =0.0;
528: PetscBool flg;
531: PetscOptionsGetInt(NULL,"-start",&start,&flg);
533: if (flg && start==0) { /* The zero vector is reasonable */
535: VecSet(X, zero);
536: /* PLogIerr(user,"Min. Surface Area Problem: Start with 0 vector \n"); */
539: } else { /* Take an average of the boundary conditions */
540: PetscInt mx=user->mx,my=user->my;
541: PetscScalar **x;
542: PetscInt xs,xm,ys,ym;
544: /* Get pointers to vector data */
545: DMDAVecGetArray(user->da,X,&x);
546: DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
548: /* Perform local computations */
549: for (j=ys; j<ys+ym; j++) {
550: for (i=xs; i< xs+xm; i++) {
551: x[j][i] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+
552: ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0;
553: }
554: }
556: /* Restore vectors */
557: DMDAVecRestoreArray(user->da,X,&x);
559: }
560: return(0);
561: }