Actual source code: ex8.c
petsc-3.3-p7 2013-05-11
1: #include <petscsnes.h>
2: #include <petscdmda.h>
4: static char help[] = "Parallel version of the minimum surface area problem using DMs.\n\
5: See ex10.c for the serial version. It solves a system of nonlinear equations in mixed\n\
6: complementarity form using semismooth newton algorithm.This example is based on a\n\
7: problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and\n\
8: boundary values along the edges of the domain, the objective is to find the\n\
9: surface with the minimal area that satisfies the boundary conditions.\n\
10: This application solves this problem using complimentarity -- We are actually\n\
11: solving the system (grad f)_i >= 0, if x_i == l_i \n\
12: (grad f)_i = 0, if l_i < x_i < u_i \n\
13: (grad f)_i <= 0, if x_i == u_i \n\
14: where f is the function to be minimized. \n\
15: \n\
16: The command line options are:\n\
17: -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
18: -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
19: -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise\n\
20: -lb <value>, lower bound on the variables\n\
21: -ub <value>, upper bound on the variables\n\n";
23: /*
24: User-defined application context - contains data needed by the
25: application-provided call-back routines, FormJacobian() and
26: FormFunction().
27: */
29: typedef struct {
30: DM da;
31: PetscScalar *bottom, *top, *left, *right;
32: PetscInt mx,my;
33: } AppCtx;
36: /* -------- User-defined Routines --------- */
38: extern PetscErrorCode MSA_BoundaryConditions(AppCtx *);
39: extern PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
40: extern PetscErrorCode FormGradient(SNES, Vec, Vec, void *);
41: extern PetscErrorCode FormJacobian(SNES, Vec, Mat *, Mat*, MatStructure*,void *);
45: int main(int argc, char **argv)
46: {
47: PetscErrorCode info; /* used to check for functions returning nonzeros */
48: Vec x,r; /* solution and residual vectors */
49: Vec xl,xu; /* Bounds on the variables */
50: PetscBool flg_l,flg_u; /* flags to check if the bounds are set */
51: SNES snes; /* nonlinear solver context */
52: Mat J; /* Jacobian matrix */
53: PetscInt N; /* Number of elements in vector */
54: PetscScalar lb = .05;
55: PetscScalar ub = SNES_VI_INF;
56: AppCtx user; /* user-defined work context */
57: PetscBool flg;
59: /* Initialize PETSc */
60: PetscInitialize(&argc, &argv, (char *)0, help );
62: #if defined(PETSC_USE_COMPLEX)
63: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example does not work for scalar type complex\n");
64: #endif
66: /* Check if lower and upper bounds are set */
67: info = PetscOptionsGetScalar(PETSC_NULL, "-lb", &lb, &flg_l);CHKERRQ(info);
68: info = PetscOptionsGetScalar(PETSC_NULL, "-ub", &ub, &flg_u);CHKERRQ(info);
70: /* Create distributed array to manage the 2d grid */
71: info = DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,&user.da);CHKERRQ(info);
72: info = DMDAGetInfo(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);CHKERRQ(info);
73: /* Extract global vectors from DMDA; */
74: info = DMCreateGlobalVector(user.da,&x);CHKERRQ(info);
75: info = VecDuplicate(x, &r); CHKERRQ(info);
77: N = user.mx*user.my;
78: info = DMCreateMatrix(user.da,MATAIJ,&J);CHKERRQ(info);
80: /* Create nonlinear solver context */
81: info = SNESCreate(PETSC_COMM_WORLD,&snes); CHKERRQ(info);
83: /* Set function evaluation and Jacobian evaluation routines */
84: info = SNESSetFunction(snes,r,FormGradient,&user);CHKERRQ(info);
85: info = SNESSetJacobian(snes,J,J,FormJacobian,&user);CHKERRQ(info);
87: /* Set the boundary conditions */
88: info = MSA_BoundaryConditions(&user); CHKERRQ(info);
90: /* Set initial solution guess */
91: info = MSA_InitialPoint(&user, x); CHKERRQ(info);
94: /* Set Bounds on variables */
95: info = VecDuplicate(x, &xl); CHKERRQ(info);
96: info = VecDuplicate(x, &xu); CHKERRQ(info);
97: info = VecSet(xl, lb); CHKERRQ(info);
98: info = VecSet(xu, ub); CHKERRQ(info);
100: info = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(info);
102: info = SNESSetFromOptions(snes);CHKERRQ(info);
104: /* Solve the application */
105: info = SNESSolve(snes,PETSC_NULL,x);CHKERRQ(info);
107: info = PetscOptionsHasName(PETSC_NULL,"-view_sol",&flg);CHKERRQ(info);
108: if (flg) { info = VecView(x,PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(info); }
110: /* Free memory */
111: info = VecDestroy(&x); CHKERRQ(info);
112: info = VecDestroy(&xl); CHKERRQ(info);
113: info = VecDestroy(&xu); CHKERRQ(info);
114: info = VecDestroy(&r); CHKERRQ(info);
115: info = MatDestroy(&J); CHKERRQ(info);
116: info = SNESDestroy(&snes); CHKERRQ(info);
118: /* Free user-created data structures */
119: info = DMDestroy(&user.da);CHKERRQ(info);
120: info = PetscFree(user.bottom); CHKERRQ(info);
121: info = PetscFree(user.top); CHKERRQ(info);
122: info = PetscFree(user.left); CHKERRQ(info);
123: info = PetscFree(user.right); CHKERRQ(info);
125: info = PetscFinalize();
127: return 0;
128: }
130: /* -------------------------------------------------------------------- */
134: /* FormGradient - Evaluates gradient of f.
136: Input Parameters:
137: . snes - the SNES context
138: . X - input vector
139: . ptr - optional user-defined context, as set by SNESSetFunction()
140:
141: Output Parameters:
142: . G - vector containing the newly evaluated gradient
143: */
144: PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr){
145: AppCtx *user = (AppCtx *) ptr;
146: int info;
147: PetscInt i,j;
148: PetscInt mx=user->mx, my=user->my;
149: PetscScalar hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
150: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
151: PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
152: PetscScalar **g, **x;
153: PetscInt xs,xm,ys,ym;
154: Vec localX;
157: /* Initialize vector to zero */
158: info = VecSet(G,0.0);CHKERRQ(info);
160: /* Get local vector */
161: info = DMGetLocalVector(user->da,&localX);CHKERRQ(info);
162: /* Get ghost points */
163: info = DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
164: info = DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
165: /* Get pointer to local vector data */
166: info = DMDAVecGetArray(user->da,localX, &x); CHKERRQ(info);
167: info = DMDAVecGetArray(user->da,G, &g); CHKERRQ(info);
169: info = DMDAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);CHKERRQ(info);
170: /* Compute function over the locally owned part of the mesh */
171: for (j=ys; j < ys+ym; j++){
172: for (i=xs; i< xs+xm; i++){
173:
174: xc = x[j][i];
175: xlt=xrb=xl=xr=xb=xt=xc;
176:
177: if (i==0){ /* left side */
178: xl= user->left[j+1];
179: xlt = user->left[j+2];
180: } else {
181: xl = x[j][i-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[j-1][i];
189: }
190:
191: if (i+1 == mx){ /* right side */
192: xr=user->right[j+1];
193: xrb = user->right[j];
194: } else {
195: xr = x[j][i+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[j+1][i];
203: }
205: if (i>0 && j+1<my){ /* left top side */
206: xlt = x[j+1][i-1];
207: }
208: if (j>0 && i+1<mx){ /* right bottom */
209: xrb = x[j-1][i+1];
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);
220:
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 = PetscSqrtReal( 1.0 + d1*d1 + d7*d7);
238: f2 = PetscSqrtReal( 1.0 + d1*d1 + d4*d4);
239: f3 = PetscSqrtReal( 1.0 + d3*d3 + d8*d8);
240: f4 = PetscSqrtReal( 1.0 + d3*d3 + d2*d2);
241: f5 = PetscSqrtReal( 1.0 + d2*d2 + d5*d5);
242: f6 = PetscSqrtReal( 1.0 + d4*d4 + d6*d6);
243:
244: df1dxc /= f1;
245: df2dxc /= f2;
246: df3dxc /= f3;
247: df4dxc /= f4;
248: df5dxc /= f5;
249: df6dxc /= f6;
251: g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc )/2.0;
252:
253: }
254: }
255:
256: /* Restore vectors */
257: info = DMDAVecRestoreArray(user->da,localX, &x); CHKERRQ(info);
258: info = DMDAVecRestoreArray(user->da,G, &g); CHKERRQ(info);
259: info = DMRestoreLocalVector(user->da,&localX);CHKERRQ(info);
260: info = PetscLogFlops(67*mx*my); CHKERRQ(info);
261: return(0);
262: }
264: /* ------------------------------------------------------------------- */
267: /*
268: FormJacobian - Evaluates Jacobian matrix.
270: Input Parameters:
271: . snes - SNES context
272: . X - input vector
273: . ptr - optional user-defined context, as set by SNESSetJacobian()
275: Output Parameters:
276: . tH - Jacobian matrix
278: */
279: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat *tH, Mat* tHPre, MatStructure* flag, void *ptr)
280: {
281: AppCtx *user = (AppCtx *) ptr;
282: Mat H = *tH;
283: PetscErrorCode info;
284: PetscInt i,j,k;
285: PetscInt mx=user->mx, my=user->my;
286: MatStencil row,col[7];
287: PetscScalar hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
288: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
289: PetscScalar hl,hr,ht,hb,hc,htl,hbr;
290: PetscScalar **x, v[7];
291: PetscBool assembled;
292: PetscInt xs,xm,ys,ym;
293: Vec localX;
296: /* Set various matrix options */
297: info = MatAssembled(H,&assembled); CHKERRQ(info);
298: if (assembled){info = MatZeroEntries(H); CHKERRQ(info);}
299: *flag=SAME_NONZERO_PATTERN;
301: /* Get local vector */
302: info = DMGetLocalVector(user->da,&localX);CHKERRQ(info);
303: /* Get ghost points */
304: info = DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
305: info = DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
306:
307: /* Get pointers to vector data */
308: info = DMDAVecGetArray(user->da,localX, &x); CHKERRQ(info);
310: info = DMDAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);CHKERRQ(info);
311: /* Compute Jacobian over the locally owned part of the mesh */
312: for (j=ys; j< ys+ym; j++){
313: for (i=xs; i< xs+xm; i++){
314: xc = x[j][i];
315: xlt=xrb=xl=xr=xb=xt=xc;
317: /* Left */
318: if (i==0){
319: xl= user->left[j+1];
320: xlt = user->left[j+2];
321: } else {
322: xl = x[j][i-1];
323: }
324:
325: /* Bottom */
326: if (j==0){
327: xb=user->bottom[i+1];
328: xrb = user->bottom[i+2];
329: } else {
330: xb = x[j-1][i];
331: }
332:
333: /* Right */
334: if (i+1 == mx){
335: xr=user->right[j+1];
336: xrb = user->right[j];
337: } else {
338: xr = x[j][i+1];
339: }
341: /* Top */
342: if (j+1==my){
343: xt=user->top[i+1];
344: xlt = user->top[i];
345: }else {
346: xt = x[j+1][i];
347: }
349: /* Top left */
350: if (i>0 && j+1<my){
351: xlt = x[j+1][i-1];
352: }
354: /* Bottom right */
355: if (j>0 && i+1<mx){
356: xrb = x[j-1][i+1];
357: }
359: d1 = (xc-xl)/hx;
360: d2 = (xc-xr)/hx;
361: d3 = (xc-xt)/hy;
362: d4 = (xc-xb)/hy;
363: d5 = (xrb-xr)/hy;
364: d6 = (xrb-xb)/hx;
365: d7 = (xlt-xl)/hy;
366: d8 = (xlt-xt)/hx;
367:
368: f1 = PetscSqrtReal( 1.0 + d1*d1 + d7*d7);
369: f2 = PetscSqrtReal( 1.0 + d1*d1 + d4*d4);
370: f3 = PetscSqrtReal( 1.0 + d3*d3 + d8*d8);
371: f4 = PetscSqrtReal( 1.0 + d3*d3 + d2*d2);
372: f5 = PetscSqrtReal( 1.0 + d2*d2 + d5*d5);
373: f6 = PetscSqrtReal( 1.0 + d4*d4 + d6*d6);
376: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
377: (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
378: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
379: (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
380: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
381: (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
382: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
383: (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
385: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
386: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
388: hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
389: hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
390: (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
391: (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
393: hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0;
395: k=0;
396: row.i = i;row.j= j;
397: /* Bottom */
398: if (j>0){
399: v[k]=hb;
400: col[k].i = i; col[k].j=j-1; k++;
401: }
402:
403: /* Bottom right */
404: if (j>0 && i < mx -1){
405: v[k]=hbr;
406: col[k].i = i+1; col[k].j = j-1; k++;
407: }
408:
409: /* left */
410: if (i>0){
411: v[k]= hl;
412: col[k].i = i-1; col[k].j = j; k++;
413: }
414:
415: /* Centre */
416: v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++;
417:
418: /* Right */
419: if (i < mx-1 ){
420: v[k]= hr;
421: col[k].i= i+1; col[k].j = j;k++;
422: }
423:
424: /* Top left */
425: if (i>0 && j < my-1 ){
426: v[k]= htl;
427: col[k].i = i-1;col[k].j = j+1; k++;
428: }
429:
430: /* Top */
431: if (j < my-1 ){
432: v[k]= ht;
433: col[k].i = i; col[k].j = j+1; k++;
434: }
435:
436: info = MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES);
437: CHKERRQ(info);
438: }
439: }
441: /* Assemble the matrix */
442: info = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
443: info = DMDAVecRestoreArray(user->da,localX,&x);CHKERRQ(info);
444: info = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
445: info = DMRestoreLocalVector(user->da,&localX);CHKERRQ(info);
447: info = PetscLogFlops(199*mx*my); CHKERRQ(info);
448: return(0);
449: }
451: /* ------------------------------------------------------------------- */
454: /*
455: MSA_BoundaryConditions - Calculates the boundary conditions for
456: the region.
458: Input Parameter:
459: . user - user-defined application context
461: Output Parameter:
462: . user - user-defined application context
463: */
464: PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
465: {
466: PetscErrorCode info;
467: PetscInt i,j,k,limit=0,maxits=5;
468: PetscInt mx=user->mx,my=user->my;
469: PetscInt bsize=0, lsize=0, tsize=0, rsize=0;
470: PetscScalar one=1.0, two=2.0, three=3.0, tol=1e-10;
471: PetscScalar fnorm,det,hx,hy,xt=0,yt=0;
472: PetscScalar u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
473: PetscScalar b=-0.5, t=0.5, l=-0.5, r=0.5;
474: PetscScalar *boundary;
477: bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
479: info = PetscMalloc(bsize*sizeof(PetscScalar), &user->bottom);CHKERRQ(info);
480: info = PetscMalloc(tsize*sizeof(PetscScalar), &user->top);CHKERRQ(info);
481: info = PetscMalloc(lsize*sizeof(PetscScalar), &user->left);CHKERRQ(info);
482: info = PetscMalloc(rsize*sizeof(PetscScalar), &user->right);CHKERRQ(info);
484: hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
486: for (j=0; j<4; j++){
487: if (j==0){
488: yt=b;
489: xt=l;
490: limit=bsize;
491: boundary=user->bottom;
492: } else if (j==1){
493: yt=t;
494: xt=l;
495: limit=tsize;
496: boundary=user->top;
497: } else if (j==2){
498: yt=b;
499: xt=l;
500: limit=lsize;
501: boundary=user->left;
502: } else { // if (j==3)
503: yt=b;
504: xt=r;
505: limit=rsize;
506: boundary=user->right;
507: }
509: for (i=0; i<limit; i++){
510: u1=xt;
511: u2=-yt;
512: for (k=0; k<maxits; k++){
513: nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
514: nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
515: fnorm=PetscSqrtReal(nf1*nf1+nf2*nf2);
516: if (fnorm <= tol) break;
517: njac11=one+u2*u2-u1*u1;
518: njac12=two*u1*u2;
519: njac21=-two*u1*u2;
520: njac22=-one - u1*u1 + u2*u2;
521: det = njac11*njac22-njac21*njac12;
522: u1 = u1-(njac22*nf1-njac12*nf2)/det;
523: u2 = u2-(njac11*nf2-njac21*nf1)/det;
524: }
526: boundary[i]=u1*u1-u2*u2;
527: if (j==0 || j==1) {
528: xt=xt+hx;
529: } else { // if (j==2 || j==3)
530: yt=yt+hy;
531: }
532: }
533: }
535: return(0);
536: }
538: /* ------------------------------------------------------------------- */
541: /*
542: MSA_InitialPoint - Calculates the initial guess in one of three ways.
544: Input Parameters:
545: . user - user-defined application context
546: . X - vector for initial guess
548: Output Parameters:
549: . X - newly computed initial guess
550: */
551: PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
552: {
553: PetscErrorCode info;
554: PetscInt start=-1,i,j;
555: PetscScalar zero=0.0;
556: PetscBool flg;
559: info = PetscOptionsGetInt(PETSC_NULL,"-start",&start,&flg); CHKERRQ(info);
561: if (flg && start==0){ /* The zero vector is reasonable */
562:
563: info = VecSet(X, zero); CHKERRQ(info);
564: /* PLogInfo(user,"Min. Surface Area Problem: Start with 0 vector \n"); */
567: } else { /* Take an average of the boundary conditions */
568: PetscInt mx=user->mx,my=user->my;
569: PetscScalar **x;
570: PetscInt xs,xm,ys,ym;
571:
572: /* Get pointers to vector data */
573: info = DMDAVecGetArray(user->da,X,&x); CHKERRQ(info);
574: info = DMDAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);CHKERRQ(info);
576: /* Perform local computations */
577: for (j=ys; j<ys+ym; j++){
578: for (i=xs; i< xs+xm; i++){
579: x[j][i] = ( ((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+
580: ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0;
581: }
582: }
583:
584: /* Restore vectors */
585: info = DMDAVecRestoreArray(user->da,X,&x); CHKERRQ(info);
586:
587: }
588: return(0);
589: }