Actual source code: ex58.c
petsc-3.8.4 2018-03-24
1: #include <petscsnes.h>
2: #include <petscdm.h>
3: #include <petscdmda.h>
5: static const char help[] = "Parallel version of the minimum surface area problem in 2D using DMDA.\n\
6: It solves a system of nonlinear equations in mixed\n\
7: complementarity form.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: /*
31: This is a new version of the ../tests/ex8.c code
33: 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 pmat -ksp_type fgmres
35: Or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
36: multigrid levels, it will be determined automatically based on the number of refinements done)
38: ./ex58 -pc_type mg -ksp_monitor -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3
39: -mg_levels_ksp_monitor -snes_vi_monitor -mg_levels_pc_type sor -pc_mg_type full
42: */
44: typedef struct {
45: PetscScalar *bottom, *top, *left, *right;
46: PetscScalar lb,ub;
47: } AppCtx;
50: /* -------- User-defined Routines --------- */
52: extern PetscErrorCode FormBoundaryConditions(SNES,AppCtx**);
53: extern PetscErrorCode DestroyBoundaryConditions(AppCtx**);
54: extern PetscErrorCode ComputeInitialGuess(SNES, Vec,void*);
55: extern PetscErrorCode FormGradient(SNES, Vec, Vec, void*);
56: extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void*);
57: extern PetscErrorCode FormBounds(SNES,Vec,Vec);
59: int main(int argc, char **argv)
60: {
62: Vec x,r; /* solution and residual vectors */
63: SNES snes; /* nonlinear solver context */
64: Mat J; /* Jacobian matrix */
65: DM da;
67: PetscInitialize(&argc, &argv, (char*)0, help);if (ierr) return ierr;
69: /* Create distributed array to manage the 2d grid */
70: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
71: DMSetFromOptions(da);
72: DMSetUp(da);
74: /* Extract global vectors from DMDA; */
75: DMCreateGlobalVector(da,&x);
76: VecDuplicate(x, &r);
78: DMSetMatType(da,MATAIJ);
79: DMCreateMatrix(da,&J);
81: /* Create nonlinear solver context */
82: SNESCreate(PETSC_COMM_WORLD,&snes);
83: SNESSetDM(snes,da);
85: /* Set function evaluation and Jacobian evaluation routines */
86: SNESSetFunction(snes,r,FormGradient,NULL);
87: SNESSetJacobian(snes,J,J,FormJacobian,NULL);
89: SNESSetComputeApplicationContext(snes,(PetscErrorCode (*)(SNES,void**))FormBoundaryConditions,(PetscErrorCode (*)(void**))DestroyBoundaryConditions);
91: SNESSetComputeInitialGuess(snes,ComputeInitialGuess,NULL);
93: SNESVISetComputeVariableBounds(snes,FormBounds);
95: SNESSetFromOptions(snes);
97: /* Solve the application */
98: SNESSolve(snes,NULL,x);
100: /* Free memory */
101: VecDestroy(&x);
102: VecDestroy(&r);
103: MatDestroy(&J);
104: SNESDestroy(&snes);
106: /* Free user-created data structures */
107: DMDestroy(&da);
109: PetscFinalize();
110: return ierr;
111: }
113: /* -------------------------------------------------------------------- */
115: /* FormBounds - sets the upper and lower bounds
117: Input Parameters:
118: . snes - the SNES context
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: /* -------------------------------------------------------------------- */
138: /* FormGradient - Evaluates gradient of f.
140: Input Parameters:
141: . snes - the SNES context
142: . X - input vector
143: . ptr - optional user-defined context, as set by SNESSetFunction()
145: Output Parameters:
146: . G - vector containing the newly evaluated gradient
147: */
148: PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
149: {
150: AppCtx *user;
151: int ierr;
152: PetscInt i,j;
153: PetscInt mx, my;
154: PetscScalar hx,hy, hydhx, hxdhy;
155: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
156: PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
157: PetscScalar **g, **x;
158: PetscInt xs,xm,ys,ym;
159: Vec localX;
160: DM da;
163: SNESGetDM(snes,&da);
164: SNESGetApplicationContext(snes,(void**)&user);
165: 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);
166: hx = 1.0/(mx+1);hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;
168: VecSet(G,0.0);
170: /* Get local vector */
171: DMGetLocalVector(da,&localX);
172: /* Get ghost points */
173: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
174: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
175: /* Get pointer to local vector data */
176: DMDAVecGetArray(da,localX, &x);
177: DMDAVecGetArray(da,G, &g);
179: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
180: /* Compute function over the locally owned part of the mesh */
181: for (j=ys; j < ys+ym; j++) {
182: for (i=xs; i< xs+xm; i++) {
184: xc = x[j][i];
185: xlt=xrb=xl=xr=xb=xt=xc;
187: if (i==0) { /* left side */
188: xl = user->left[j+1];
189: xlt = user->left[j+2];
190: } else xl = x[j][i-1];
192: if (j==0) { /* bottom side */
193: xb = user->bottom[i+1];
194: xrb = user->bottom[i+2];
195: } else xb = x[j-1][i];
197: if (i+1 == mx) { /* right side */
198: xr = user->right[j+1];
199: xrb = user->right[j];
200: } else xr = x[j][i+1];
202: if (j+1==0+my) { /* top side */
203: xt = user->top[i+1];
204: xlt = user->top[i];
205: } else xt = x[j+1][i];
207: if (i>0 && j+1<my) xlt = x[j+1][i-1]; /* left top side */
208: if (j>0 && i+1<mx) xrb = x[j-1][i+1]; /* right bottom */
210: d1 = (xc-xl);
211: d2 = (xc-xr);
212: d3 = (xc-xt);
213: d4 = (xc-xb);
214: d5 = (xr-xrb);
215: d6 = (xrb-xb);
216: d7 = (xlt-xl);
217: d8 = (xt-xlt);
219: df1dxc = d1*hydhx;
220: df2dxc = (d1*hydhx + d4*hxdhy);
221: df3dxc = d3*hxdhy;
222: df4dxc = (d2*hydhx + d3*hxdhy);
223: df5dxc = d2*hydhx;
224: df6dxc = d4*hxdhy;
226: d1 /= hx;
227: d2 /= hx;
228: d3 /= hy;
229: d4 /= hy;
230: d5 /= hy;
231: d6 /= hx;
232: d7 /= hy;
233: d8 /= hx;
235: f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
236: f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
237: f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
238: f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
239: f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
240: f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
242: df1dxc /= f1;
243: df2dxc /= f2;
244: df3dxc /= f3;
245: df4dxc /= f4;
246: df5dxc /= f5;
247: df6dxc /= f6;
249: g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
251: }
252: }
254: /* Restore vectors */
255: DMDAVecRestoreArray(da,localX, &x);
256: DMDAVecRestoreArray(da,G, &g);
257: DMRestoreLocalVector(da,&localX);
258: PetscLogFlops(67*mx*my);
259: return(0);
260: }
262: /* ------------------------------------------------------------------- */
263: /*
264: FormJacobian - Evaluates Jacobian matrix.
266: Input Parameters:
267: . snes - SNES context
268: . X - input vector
269: . ptr - optional user-defined context, as set by SNESSetJacobian()
271: Output Parameters:
272: . tH - Jacobian matrix
274: */
275: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *ptr)
276: {
277: AppCtx *user;
279: PetscInt i,j,k;
280: PetscInt mx, my;
281: MatStencil row,col[7];
282: PetscScalar hx, hy, hydhx, hxdhy;
283: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
284: PetscScalar hl,hr,ht,hb,hc,htl,hbr;
285: PetscScalar **x, v[7];
286: PetscBool assembled;
287: PetscInt xs,xm,ys,ym;
288: Vec localX;
289: DM da;
292: SNESGetDM(snes,&da);
293: SNESGetApplicationContext(snes,(void**)&user);
294: 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);
295: hx = 1.0/(mx+1); hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;
297: /* Set various matrix options */
298: MatAssembled(H,&assembled);
299: if (assembled) {MatZeroEntries(H);}
301: /* Get local vector */
302: DMGetLocalVector(da,&localX);
303: /* Get ghost points */
304: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
305: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
307: /* Get pointers to vector data */
308: DMDAVecGetArray(da,localX, &x);
310: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
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 xl = x[j][i-1];
323: /* Bottom */
324: if (j==0) {
325: xb =user->bottom[i+1];
326: xrb = user->bottom[i+2];
327: } else xb = x[j-1][i];
329: /* Right */
330: if (i+1 == mx) {
331: xr =user->right[j+1];
332: xrb = user->right[j];
333: } else xr = x[j][i+1];
335: /* Top */
336: if (j+1==my) {
337: xt =user->top[i+1];
338: xlt = user->top[i];
339: } else xt = x[j+1][i];
341: /* Top left */
342: if (i>0 && j+1<my) xlt = x[j+1][i-1];
344: /* Bottom right */
345: if (j>0 && i+1<mx) xrb = x[j-1][i+1];
347: d1 = (xc-xl)/hx;
348: d2 = (xc-xr)/hx;
349: d3 = (xc-xt)/hy;
350: d4 = (xc-xb)/hy;
351: d5 = (xrb-xr)/hy;
352: d6 = (xrb-xb)/hx;
353: d7 = (xlt-xl)/hy;
354: d8 = (xlt-xt)/hx;
356: f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
357: f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
358: f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
359: f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
360: f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
361: f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
364: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
365: (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
366: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
367: (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
368: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
369: (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
370: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
371: (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
373: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
374: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
376: hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
377: hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
378: (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2.0*d1*d4)/(f2*f2*f2) +
379: (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2.0*d2*d3)/(f4*f4*f4);
381: hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0;
383: k =0;
384: row.i = i;row.j= j;
385: /* Bottom */
386: if (j>0) {
387: v[k] =hb;
388: col[k].i = i; col[k].j=j-1; k++;
389: }
391: /* Bottom right */
392: if (j>0 && i < mx -1) {
393: v[k] =hbr;
394: col[k].i = i+1; col[k].j = j-1; k++;
395: }
397: /* left */
398: if (i>0) {
399: v[k] = hl;
400: col[k].i = i-1; col[k].j = j; k++;
401: }
403: /* Centre */
404: v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++;
406: /* Right */
407: if (i < mx-1) {
408: v[k] = hr;
409: col[k].i= i+1; col[k].j = j;k++;
410: }
412: /* Top left */
413: if (i>0 && j < my-1) {
414: v[k] = htl;
415: col[k].i = i-1;col[k].j = j+1; k++;
416: }
418: /* Top */
419: if (j < my-1) {
420: v[k] = ht;
421: col[k].i = i; col[k].j = j+1; k++;
422: }
424: MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES);
425: }
426: }
428: /* Assemble the matrix */
429: MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
430: DMDAVecRestoreArray(da,localX,&x);
431: MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
432: DMRestoreLocalVector(da,&localX);
434: PetscLogFlops(199*mx*my);
435: return(0);
436: }
438: /* ------------------------------------------------------------------- */
439: /*
440: FormBoundaryConditions - Calculates the boundary conditions for
441: the region.
443: Input Parameter:
444: . user - user-defined application context
446: Output Parameter:
447: . user - user-defined application context
448: */
449: PetscErrorCode FormBoundaryConditions(SNES snes,AppCtx **ouser)
450: {
452: PetscInt i,j,k,limit=0,maxits=5;
453: PetscInt mx,my;
454: PetscInt bsize=0, lsize=0, tsize=0, rsize=0;
455: PetscScalar one =1.0, two=2.0, three=3.0;
456: PetscScalar det,hx,hy,xt=0,yt=0;
457: PetscReal fnorm, tol=1e-10;
458: PetscScalar u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
459: PetscScalar b=-0.5, t=0.5, l=-0.5, r=0.5;
460: PetscScalar *boundary;
461: AppCtx *user;
462: DM da;
465: SNESGetDM(snes,&da);
466: PetscNew(&user);
467: *ouser = user;
468: user->lb = .05;
469: user->ub = PETSC_INFINITY;
470: 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);
472: /* Check if lower and upper bounds are set */
473: PetscOptionsGetScalar(NULL,NULL, "-lb", &user->lb, 0);
474: PetscOptionsGetScalar(NULL,NULL, "-ub", &user->ub, 0);
475: bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
477: PetscMalloc1(bsize, &user->bottom);
478: PetscMalloc1(tsize, &user->top);
479: PetscMalloc1(lsize, &user->left);
480: PetscMalloc1(rsize, &user->right);
482: hx= (r-l)/(mx+1.0); hy=(t-b)/(my+1.0);
484: for (j=0; j<4; j++) {
485: if (j==0) {
486: yt = b;
487: xt = l;
488: limit = bsize;
489: boundary = user->bottom;
490: } else if (j==1) {
491: yt = t;
492: xt = l;
493: limit = tsize;
494: boundary = user->top;
495: } else if (j==2) {
496: yt = b;
497: xt = l;
498: limit = lsize;
499: boundary = user->left;
500: } else { /* if (j==3) */
501: yt = b;
502: xt = r;
503: limit = rsize;
504: boundary = user->right;
505: }
507: for (i=0; i<limit; i++) {
508: u1=xt;
509: u2=-yt;
510: for (k=0; k<maxits; k++) {
511: nf1 = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
512: nf2 = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
513: fnorm = PetscRealPart(PetscSqrtScalar(nf1*nf1+nf2*nf2));
514: if (fnorm <= tol) break;
515: njac11=one+u2*u2-u1*u1;
516: njac12=two*u1*u2;
517: njac21=-two*u1*u2;
518: njac22=-one - u1*u1 + u2*u2;
519: det = njac11*njac22-njac21*njac12;
520: u1 = u1-(njac22*nf1-njac12*nf2)/det;
521: u2 = u2-(njac11*nf2-njac21*nf1)/det;
522: }
524: boundary[i]=u1*u1-u2*u2;
525: if (j==0 || j==1) xt=xt+hx;
526: else yt=yt+hy; /* if (j==2 || j==3) */
527: }
528: }
529: return(0);
530: }
532: PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser)
533: {
535: AppCtx *user = *ouser;
538: PetscFree(user->bottom);
539: PetscFree(user->top);
540: PetscFree(user->left);
541: PetscFree(user->right);
542: PetscFree(*ouser);
543: return(0);
544: }
547: /* ------------------------------------------------------------------- */
548: /*
549: ComputeInitialGuess - Calculates the initial guess
551: Input Parameters:
552: . user - user-defined application context
553: . X - vector for initial guess
555: Output Parameters:
556: . X - newly computed initial guess
557: */
558: PetscErrorCode ComputeInitialGuess(SNES snes, Vec X,void *dummy)
559: {
561: PetscInt i,j,mx,my;
562: DM da;
563: AppCtx *user;
564: PetscScalar **x;
565: PetscInt xs,xm,ys,ym;
568: SNESGetDM(snes,&da);
569: SNESGetApplicationContext(snes,(void**)&user);
571: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
572: 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);
574: /* Get pointers to vector data */
575: DMDAVecGetArray(da,X,&x);
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.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;
580: }
581: }
582: /* Restore vectors */
583: DMDAVecRestoreArray(da,X,&x);
584: return(0);
585: }