Actual source code: ex58.c
petsc-3.3-p7 2013-05-11
1: #include <petscsnes.h>
2: #include <petscdmda.h>
4: static const char help[] = "Parallel version of the minimum surface area problem in 2D using DMDA.\n\
5: It solves a system of nonlinear equations in mixed\n\
6: complementarity form.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: /*
30: This is a new version of the ../tests/ex8.c code
32: Run, for example, with the options ./ex58 -snes_vi_monitor -ksp_monitor -mg_levels_ksp_monitor -pc_type mg -pc_mg_levels 2 -pc_mg_galerkin -ksp_type fgmres
34: Or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
35: multigrid levels, it will be determined automatically based on the number of refinements done)
37: ./ex58 -pc_type mg -ksp_monitor -snes_view -pc_mg_galerkin -snes_grid_sequence 3
38: -mg_levels_ksp_monitor -snes_vi_monitor -mg_levels_pc_type sor -pc_mg_type full
41: */
43: typedef struct {
44: PetscScalar *bottom, *top, *left, *right;
45: PetscScalar lb,ub;
46: } AppCtx;
49: /* -------- User-defined Routines --------- */
51: extern PetscErrorCode FormBoundaryConditions(SNES,AppCtx **);
52: extern PetscErrorCode DestroyBoundaryConditions(AppCtx **);
53: extern PetscErrorCode ComputeInitialGuess(SNES, Vec,void*);
54: extern PetscErrorCode FormGradient(SNES, Vec, Vec, void *);
55: extern PetscErrorCode FormJacobian(SNES, Vec, Mat *, Mat*, MatStructure*,void *);
56: extern PetscErrorCode FormBounds(SNES,Vec,Vec);
60: int main(int argc, char **argv)
61: {
62: PetscErrorCode ierr; /* used to check for functions returning nonzeros */
63: Vec x,r; /* solution and residual vectors */
64: SNES snes; /* nonlinear solver context */
65: Mat J; /* Jacobian matrix */
66: DM da;
68: PetscInitialize(&argc, &argv, (char *)0, help );
70: /* Create distributed array to manage the 2d grid */
71: 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,&da);
73: /* Extract global vectors from DMDA; */
74: DMCreateGlobalVector(da,&x);
75: VecDuplicate(x, &r);
77: DMCreateMatrix(da,MATAIJ,&J);
79: /* Create nonlinear solver context */
80: SNESCreate(PETSC_COMM_WORLD,&snes);
81: SNESSetDM(snes,da);
83: /* Set function evaluation and Jacobian evaluation routines */
84: SNESSetFunction(snes,r,FormGradient,PETSC_NULL);
85: SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL);
87: SNESSetComputeApplicationContext(snes,(PetscErrorCode (*)(SNES,void**))FormBoundaryConditions,(PetscErrorCode (*)(void**))DestroyBoundaryConditions);
89: SNESSetComputeInitialGuess(snes,ComputeInitialGuess,PETSC_NULL);
91: SNESVISetComputeVariableBounds(snes,FormBounds);
93: SNESSetFromOptions(snes);
95: /* Solve the application */
96: SNESSolve(snes,PETSC_NULL,x);
98: /* Free memory */
99: VecDestroy(&x);
100: VecDestroy(&r);
101: MatDestroy(&J);
102: SNESDestroy(&snes);
104: /* Free user-created data structures */
105: DMDestroy(&da);
107: PetscFinalize();
108: return 0;
109: }
111: /* -------------------------------------------------------------------- */
115: /* FormBounds - sets the upper and lower bounds
117: Input Parameters:
118: . snes - the SNES context
119:
120: Output Parameters:
121: . xl - lower bounds
122: . xu - upper bounds
123: */
124: PetscErrorCode FormBounds(SNES snes, Vec xl, Vec xu)
125: {
127: AppCtx *ctx;
130: SNESGetApplicationContext(snes,&ctx);
131: VecSet(xl,ctx->lb);
132: VecSet(xu,ctx->ub);
133: return(0);
134: }
136: /* -------------------------------------------------------------------- */
140: /* FormGradient - Evaluates gradient of f.
142: Input Parameters:
143: . snes - the SNES context
144: . X - input vector
145: . ptr - optional user-defined context, as set by SNESSetFunction()
146:
147: Output Parameters:
148: . G - vector containing the newly evaluated gradient
149: */
150: PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
151: {
152: AppCtx *user;
153: int ierr;
154: PetscInt i,j;
155: PetscInt mx, my;
156: PetscScalar hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
157: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
158: PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
159: PetscScalar **g, **x;
160: PetscInt xs,xm,ys,ym;
161: Vec localX;
162: DM da;
165: SNESGetDM(snes,&da);
166: SNESGetApplicationContext(snes,(void**)&user);
167: DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
168: hx=1.0/(mx+1);hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;
170: VecSet(G,0.0);
172: /* Get local vector */
173: DMGetLocalVector(da,&localX);
174: /* Get ghost points */
175: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
176: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
177: /* Get pointer to local vector data */
178: DMDAVecGetArray(da,localX, &x);
179: DMDAVecGetArray(da,G, &g);
181: DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
182: /* Compute function over the locally owned part of the mesh */
183: for (j=ys; j < ys+ym; j++){
184: for (i=xs; i< xs+xm; i++){
185:
186: xc = x[j][i];
187: xlt=xrb=xl=xr=xb=xt=xc;
188:
189: if (i==0){ /* left side */
190: xl= user->left[j+1];
191: xlt = user->left[j+2];
192: } else {
193: xl = x[j][i-1];
194: }
196: if (j==0){ /* bottom side */
197: xb=user->bottom[i+1];
198: xrb = user->bottom[i+2];
199: } else {
200: xb = x[j-1][i];
201: }
202:
203: if (i+1 == mx){ /* right side */
204: xr=user->right[j+1];
205: xrb = user->right[j];
206: } else {
207: xr = x[j][i+1];
208: }
210: if (j+1==0+my){ /* top side */
211: xt=user->top[i+1];
212: xlt = user->top[i];
213: }else {
214: xt = x[j+1][i];
215: }
217: if (i>0 && j+1<my){ /* left top side */
218: xlt = x[j+1][i-1];
219: }
220: if (j>0 && i+1<mx){ /* right bottom */
221: xrb = x[j-1][i+1];
222: }
224: d1 = (xc-xl);
225: d2 = (xc-xr);
226: d3 = (xc-xt);
227: d4 = (xc-xb);
228: d5 = (xr-xrb);
229: d6 = (xrb-xb);
230: d7 = (xlt-xl);
231: d8 = (xt-xlt);
232:
233: df1dxc = d1*hydhx;
234: df2dxc = ( d1*hydhx + d4*hxdhy );
235: df3dxc = d3*hxdhy;
236: df4dxc = ( d2*hydhx + d3*hxdhy );
237: df5dxc = d2*hydhx;
238: df6dxc = d4*hxdhy;
240: d1 /= hx;
241: d2 /= hx;
242: d3 /= hy;
243: d4 /= hy;
244: d5 /= hy;
245: d6 /= hx;
246: d7 /= hy;
247: d8 /= hx;
249: f1 = sqrt( 1.0 + d1*d1 + d7*d7);
250: f2 = sqrt( 1.0 + d1*d1 + d4*d4);
251: f3 = sqrt( 1.0 + d3*d3 + d8*d8);
252: f4 = sqrt( 1.0 + d3*d3 + d2*d2);
253: f5 = sqrt( 1.0 + d2*d2 + d5*d5);
254: f6 = sqrt( 1.0 + d4*d4 + d6*d6);
255:
256: df1dxc /= f1;
257: df2dxc /= f2;
258: df3dxc /= f3;
259: df4dxc /= f4;
260: df5dxc /= f5;
261: df6dxc /= f6;
263: g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc )/2.0;
264:
265: }
266: }
267:
268: /* Restore vectors */
269: DMDAVecRestoreArray(da,localX, &x);
270: DMDAVecRestoreArray(da,G, &g);
271: DMRestoreLocalVector(da,&localX);
272: PetscLogFlops(67*mx*my);
273: return(0);
274: }
276: /* ------------------------------------------------------------------- */
279: /*
280: FormJacobian - Evaluates Jacobian matrix.
282: Input Parameters:
283: . snes - SNES context
284: . X - input vector
285: . ptr - optional user-defined context, as set by SNESSetJacobian()
287: Output Parameters:
288: . tH - Jacobian matrix
290: */
291: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat *tH, Mat* tHPre, MatStructure* flag, void *ptr)
292: {
293: AppCtx *user;
294: Mat H = *tH;
295: PetscErrorCode ierr;
296: PetscInt i,j,k;
297: PetscInt mx, my;
298: MatStencil row,col[7];
299: PetscScalar hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
300: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
301: PetscScalar hl,hr,ht,hb,hc,htl,hbr;
302: PetscScalar **x, v[7];
303: PetscBool assembled;
304: PetscInt xs,xm,ys,ym;
305: Vec localX;
306: DM da;
309: SNESGetDM(snes,&da);
310: SNESGetApplicationContext(snes,(void**)&user);
311: DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
312: hx=1.0/(mx+1); hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;
314: /* Set various matrix options */
315: MatAssembled(H,&assembled);
316: if (assembled){MatZeroEntries(H); }
317: *flag=SAME_NONZERO_PATTERN;
319: /* Get local vector */
320: DMGetLocalVector(da,&localX);
321: /* Get ghost points */
322: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
323: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
324:
325: /* Get pointers to vector data */
326: DMDAVecGetArray(da,localX, &x);
328: DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
329: /* Compute Jacobian over the locally owned part of the mesh */
330: for (j=ys; j< ys+ym; j++){
331: for (i=xs; i< xs+xm; i++){
332: xc = x[j][i];
333: xlt=xrb=xl=xr=xb=xt=xc;
335: /* Left */
336: if (i==0){
337: xl= user->left[j+1];
338: xlt = user->left[j+2];
339: } else {
340: xl = x[j][i-1];
341: }
342:
343: /* Bottom */
344: if (j==0){
345: xb=user->bottom[i+1];
346: xrb = user->bottom[i+2];
347: } else {
348: xb = x[j-1][i];
349: }
350:
351: /* Right */
352: if (i+1 == mx){
353: xr=user->right[j+1];
354: xrb = user->right[j];
355: } else {
356: xr = x[j][i+1];
357: }
359: /* Top */
360: if (j+1==my){
361: xt=user->top[i+1];
362: xlt = user->top[i];
363: }else {
364: xt = x[j+1][i];
365: }
367: /* Top left */
368: if (i>0 && j+1<my){
369: xlt = x[j+1][i-1];
370: }
372: /* Bottom right */
373: if (j>0 && i+1<mx){
374: xrb = x[j-1][i+1];
375: }
377: d1 = (xc-xl)/hx;
378: d2 = (xc-xr)/hx;
379: d3 = (xc-xt)/hy;
380: d4 = (xc-xb)/hy;
381: d5 = (xrb-xr)/hy;
382: d6 = (xrb-xb)/hx;
383: d7 = (xlt-xl)/hy;
384: d8 = (xlt-xt)/hx;
385:
386: f1 = sqrt( 1.0 + d1*d1 + d7*d7);
387: f2 = sqrt( 1.0 + d1*d1 + d4*d4);
388: f3 = sqrt( 1.0 + d3*d3 + d8*d8);
389: f4 = sqrt( 1.0 + d3*d3 + d2*d2);
390: f5 = sqrt( 1.0 + d2*d2 + d5*d5);
391: f6 = sqrt( 1.0 + d4*d4 + d6*d6);
394: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
395: (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
396: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
397: (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
398: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
399: (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
400: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
401: (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
403: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
404: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
406: hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
407: hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
408: (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2.0*d1*d4)/(f2*f2*f2) +
409: (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2.0*d2*d3)/(f4*f4*f4);
411: hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0;
413: k=0;
414: row.i = i;row.j= j;
415: /* Bottom */
416: if (j>0){
417: v[k]=hb;
418: col[k].i = i; col[k].j=j-1; k++;
419: }
420:
421: /* Bottom right */
422: if (j>0 && i < mx -1){
423: v[k]=hbr;
424: col[k].i = i+1; col[k].j = j-1; k++;
425: }
426:
427: /* left */
428: if (i>0){
429: v[k]= hl;
430: col[k].i = i-1; col[k].j = j; k++;
431: }
432:
433: /* Centre */
434: v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++;
435:
436: /* Right */
437: if (i < mx-1 ){
438: v[k]= hr;
439: col[k].i= i+1; col[k].j = j;k++;
440: }
441:
442: /* Top left */
443: if (i>0 && j < my-1 ){
444: v[k]= htl;
445: col[k].i = i-1;col[k].j = j+1; k++;
446: }
447:
448: /* Top */
449: if (j < my-1 ){
450: v[k]= ht;
451: col[k].i = i; col[k].j = j+1; k++;
452: }
453:
454: MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES);
455:
456: }
457: }
459: /* Assemble the matrix */
460: MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
461: DMDAVecRestoreArray(da,localX,&x);
462: MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
463: DMRestoreLocalVector(da,&localX);
465: PetscLogFlops(199*mx*my);
466: return(0);
467: }
469: /* ------------------------------------------------------------------- */
472: /*
473: FormBoundaryConditions - Calculates the boundary conditions for
474: the region.
476: Input Parameter:
477: . user - user-defined application context
479: Output Parameter:
480: . user - user-defined application context
481: */
482: PetscErrorCode FormBoundaryConditions(SNES snes,AppCtx **ouser)
483: {
484: PetscErrorCode ierr;
485: PetscInt i,j,k,limit=0,maxits=5;
486: PetscInt mx,my;
487: PetscInt bsize=0, lsize=0, tsize=0, rsize=0;
488: PetscScalar one=1.0, two=2.0, three=3.0;
489: PetscScalar det,hx,hy,xt=0,yt=0;
490: PetscReal fnorm, tol=1e-10;
491: PetscScalar u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
492: PetscScalar b=-0.5, t=0.5, l=-0.5, r=0.5;
493: PetscScalar *boundary;
494: AppCtx *user;
495: DM da;
498: SNESGetDM(snes,&da);
499: PetscNew(AppCtx,&user);
500: *ouser = user;
501: user->lb = .05;
502: user->ub = SNES_VI_INF;
503: DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
505: /* Check if lower and upper bounds are set */
506: PetscOptionsGetScalar(PETSC_NULL, "-lb", &user->lb, 0);
507: PetscOptionsGetScalar(PETSC_NULL, "-ub", &user->ub, 0);
508: bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
510: PetscMalloc(bsize*sizeof(PetscScalar), &user->bottom);
511: PetscMalloc(tsize*sizeof(PetscScalar), &user->top);
512: PetscMalloc(lsize*sizeof(PetscScalar), &user->left);
513: PetscMalloc(rsize*sizeof(PetscScalar), &user->right);
515: hx= (r-l)/(mx+1.0); hy=(t-b)/(my+1.0);
517: for (j=0; j<4; j++){
518: if (j==0){
519: yt=b;
520: xt=l;
521: limit=bsize;
522: boundary=user->bottom;
523: } else if (j==1){
524: yt=t;
525: xt=l;
526: limit=tsize;
527: boundary=user->top;
528: } else if (j==2){
529: yt=b;
530: xt=l;
531: limit=lsize;
532: boundary=user->left;
533: } else { // if (j==3)
534: yt=b;
535: xt=r;
536: limit=rsize;
537: boundary=user->right;
538: }
540: for (i=0; i<limit; i++){
541: u1=xt;
542: u2=-yt;
543: for (k=0; k<maxits; k++){
544: nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
545: nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
546: fnorm=PetscRealPart(sqrt(nf1*nf1+nf2*nf2));
547: if (fnorm <= tol) break;
548: njac11=one+u2*u2-u1*u1;
549: njac12=two*u1*u2;
550: njac21=-two*u1*u2;
551: njac22=-one - u1*u1 + u2*u2;
552: det = njac11*njac22-njac21*njac12;
553: u1 = u1-(njac22*nf1-njac12*nf2)/det;
554: u2 = u2-(njac11*nf2-njac21*nf1)/det;
555: }
557: boundary[i]=u1*u1-u2*u2;
558: if (j==0 || j==1) {
559: xt=xt+hx;
560: } else { // if (j==2 || j==3)
561: yt=yt+hy;
562: }
563: }
564: }
565: return(0);
566: }
570: PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser)
571: {
572: PetscErrorCode ierr;
573: AppCtx *user = *ouser;
576: PetscFree(user->bottom);
577: PetscFree(user->top);
578: PetscFree(user->left);
579: PetscFree(user->right);
580: PetscFree(*ouser);
581: return(0);
582: }
585: /* ------------------------------------------------------------------- */
588: /*
589: ComputeInitialGuess - Calculates the initial guess
591: Input Parameters:
592: . user - user-defined application context
593: . X - vector for initial guess
595: Output Parameters:
596: . X - newly computed initial guess
597: */
598: PetscErrorCode ComputeInitialGuess(SNES snes, Vec X,void *dummy)
599: {
600: PetscErrorCode ierr;
601: PetscInt i,j,mx,my;
602: DM da;
603: AppCtx *user;
604: PetscScalar **x;
605: PetscInt xs,xm,ys,ym;
608: SNESGetDM(snes,&da);
609: SNESGetApplicationContext(snes,(void**)&user);
611: DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
612: DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
614: /* Get pointers to vector data */
615: DMDAVecGetArray(da,X,&x);
616: /* Perform local computations */
617: for (j=ys; j<ys+ym; j++){
618: for (i=xs; i< xs+xm; i++){
619: x[j][i] = ( ((j+1.0)*user->bottom[i+1]+(my-j+1.0)*user->top[i+1])/(my+2.0)+((i+1.0)*user->left[j+1]+(mx-i+1.0)*user->right[j+1])/(mx+2.0))/2.0;
620: }
621: }
622: /* Restore vectors */
623: DMDAVecRestoreArray(da,X,&x);
624: return(0);
625: }