Actual source code: ex16.c
petsc-3.4.5 2014-06-29
1: #include <petscsnes.h>
2: #include <petscdmda.h>
4: static char help[]=
5: "This example is an implementation of minimal surface area with \n\
6: a plate problem from the TAO package (examples/plate2.c) \n\
7: This example is based on a problem from the MINPACK-2 test suite.\n\
8: Given a rectangular 2-D domain, boundary values along the edges of \n\
9: the domain, and a plate represented by lower boundary conditions, \n\
10: the objective is to find the surface with the minimal area that \n\
11: satisfies the boundary conditions.\n\
12: The command line options are:\n\
13: -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\
14: -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\
15: -bheight <ht>, where <ht> = height of the plate\n\
16: -start <st>, where <st> =0 for zero vector, <st> != 0 \n\
17: for an average of the boundary conditions\n\n";
19: /*
20: User-defined application context - contains data needed by the
21: application-provided call-back routines, FormJacobian() and
22: FormFunction().
23: */
25: typedef struct {
26: DM da;
27: Vec Bottom, Top, Left, Right;
28: PetscScalar bheight;
29: PetscInt mx,my,bmx,bmy;
30: } AppCtx;
33: /* -------- User-defined Routines --------- */
35: PetscErrorCode MSA_BoundaryConditions(AppCtx*);
36: PetscErrorCode MSA_InitialPoint(AppCtx*, Vec);
37: PetscErrorCode MSA_Plate(Vec,Vec,void*);
38: PetscErrorCode FormGradient(SNES, Vec, Vec, void*);
39: PetscErrorCode FormJacobian(SNES, Vec, Mat*, Mat*, MatStructure*,void*);
43: int main(int argc, char **argv)
44: {
45: PetscErrorCode info; /* used to check for functions returning nonzeros */
46: Vec x,r; /* solution and residual vectors */
47: Vec xl,xu; /* Bounds on the variables */
48: SNES snes; /* nonlinear solver context */
49: Mat J; /* Jacobian matrix */
50: AppCtx user; /* user-defined work context */
51: PetscBool flg;
53: /* Initialize PETSc */
54: PetscInitialize(&argc, &argv, (char*)0, help);
56: #if defined(PETSC_USE_COMPLEX)
57: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example does not work for scalar type complex\n");
58: #endif
60: /* Create distributed array to manage the 2d grid */
61: info = DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,-10,-10,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.da);CHKERRQ(info);
62: 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);
64: user.bheight = 0.1;
65: info = PetscOptionsGetScalar(NULL,"-bheight",&user.bheight,&flg);CHKERRQ(info);
67: user.bmx = user.mx/2; user.bmy = user.my/2;
68: info = PetscOptionsGetInt(NULL,"-bmx",&user.bmx,&flg);CHKERRQ(info);
69: info = PetscOptionsGetInt(NULL,"-bmy",&user.bmy,&flg);CHKERRQ(info);
71: PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n");
72: PetscPrintf(PETSC_COMM_WORLD,"mx:%d, my:%d, bmx:%d, bmy:%d, height:%4.2f\n",
73: user.mx,user.my,user.bmx,user.bmy,user.bheight);
75: /* Extract global vectors from DMDA; */
76: info = DMCreateGlobalVector(user.da,&x);CHKERRQ(info);
77: info = VecDuplicate(x, &r);CHKERRQ(info);
79: info = DMCreateMatrix(user.da,MATAIJ,&J);CHKERRQ(info);
81: /* Create nonlinear solver context */
82: info = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(info);
84: /* Set function evaluation and Jacobian evaluation routines */
85: info = SNESSetDM(snes,user.da);CHKERRQ(info);
86: info = SNESSetFunction(snes,r,FormGradient,&user);CHKERRQ(info);
87: info = SNESSetJacobian(snes,J,J,FormJacobian,&user);CHKERRQ(info);
89: /* Set the boundary conditions */
90: info = MSA_BoundaryConditions(&user);CHKERRQ(info);
92: /* Set initial solution guess */
93: info = MSA_InitialPoint(&user, x);CHKERRQ(info);
95: info = SNESSetFromOptions(snes);CHKERRQ(info);
97: /* Set Bounds on variables */
98: info = VecDuplicate(x, &xl);CHKERRQ(info);
99: info = VecDuplicate(x, &xu);CHKERRQ(info);
100: info = MSA_Plate(xl,xu,&user);CHKERRQ(info);
102: info = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(info);
104: /* Solve the application */
105: info = SNESSolve(snes,NULL,x);CHKERRQ(info);
107: info = PetscOptionsHasName(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 = VecDestroy(&user.Bottom);CHKERRQ(info);
121: info = VecDestroy(&user.Top);CHKERRQ(info);
122: info = VecDestroy(&user.Left);CHKERRQ(info);
123: info = VecDestroy(&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()
141: Output Parameters:
142: . G - vector containing the newly evaluated gradient
143: */
144: PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
145: {
146: AppCtx *user = (AppCtx*) ptr;
147: int info;
148: PetscInt i,j;
149: PetscInt mx=user->mx, my=user->my;
150: PetscScalar hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
151: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
152: PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
153: PetscScalar **g, **x;
154: PetscInt xs,xm,ys,ym;
155: Vec localX;
156: PetscScalar *top,*bottom,*left,*right;
159: /* Initialize vector to zero */
160: info = VecSet(G,0.0);CHKERRQ(info);
162: /* Get local vector */
163: info = DMGetLocalVector(user->da,&localX);CHKERRQ(info);
164: info = VecGetArray(user->Top,&top);CHKERRQ(info);
165: info = VecGetArray(user->Bottom,&bottom);CHKERRQ(info);
166: info = VecGetArray(user->Left,&left);CHKERRQ(info);
167: info = VecGetArray(user->Right,&right);CHKERRQ(info);
169: /* Get ghost points */
170: info = DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
171: info = DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
172: /* Get pointers to local vector data */
173: info = DMDAVecGetArray(user->da,localX, &x);CHKERRQ(info);
174: info = DMDAVecGetArray(user->da,G, &g);CHKERRQ(info);
176: info = DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(info);
177: /* Compute function over the locally owned part of the mesh */
178: for (j=ys; j < ys+ym; j++) {
179: for (i=xs; i< xs+xm; i++) {
181: xc = x[j][i];
182: xlt=xrb=xl=xr=xb=xt=xc;
184: if (i==0) { /* left side */
185: xl = left[j-ys+1];
186: xlt = left[j-ys+2];
187: } else xl = x[j][i-1];
189: if (j==0) { /* bottom side */
190: xb = bottom[i-xs+1];
191: xrb = bottom[i-xs+2];
192: } else xb = x[j-1][i];
194: if (i+1 == mx) { /* right side */
195: xr = right[j-ys+1];
196: xrb = right[j-ys];
197: } else xr = x[j][i+1];
199: if (j+1==0+my) { /* top side */
200: xt = top[i-xs+1];
201: xlt = top[i-xs];
202: } else xt = x[j+1][i];
204: if (i>0 && j+1<my) xlt = x[j+1][i-1]; /* left top side */
205: if (j>0 && i+1<mx) xrb = x[j-1][i+1]; /* right bottom */
207: d1 = (xc-xl);
208: d2 = (xc-xr);
209: d3 = (xc-xt);
210: d4 = (xc-xb);
211: d5 = (xr-xrb);
212: d6 = (xrb-xb);
213: d7 = (xlt-xl);
214: d8 = (xt-xlt);
216: df1dxc = d1*hydhx;
217: df2dxc = (d1*hydhx + d4*hxdhy);
218: df3dxc = d3*hxdhy;
219: df4dxc = (d2*hydhx + d3*hxdhy);
220: df5dxc = d2*hydhx;
221: df6dxc = d4*hxdhy;
223: d1 /= hx;
224: d2 /= hx;
225: d3 /= hy;
226: d4 /= hy;
227: d5 /= hy;
228: d6 /= hx;
229: d7 /= hy;
230: d8 /= hx;
232: f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
233: f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
234: f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
235: f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
236: f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
237: f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
239: df1dxc /= f1;
240: df2dxc /= f2;
241: df3dxc /= f3;
242: df4dxc /= f4;
243: df5dxc /= f5;
244: df6dxc /= f6;
246: g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
248: }
249: }
251: /* Restore vectors */
252: info = DMDAVecRestoreArray(user->da,localX, &x);CHKERRQ(info);
253: info = DMDAVecRestoreArray(user->da,G, &g);CHKERRQ(info);
254: info = DMRestoreLocalVector(user->da,&localX);CHKERRQ(info);
256: info = VecRestoreArray(user->Left,&left);CHKERRQ(info);
257: info = VecRestoreArray(user->Top,&top);CHKERRQ(info);
258: info = VecRestoreArray(user->Bottom,&bottom);CHKERRQ(info);
259: info = VecRestoreArray(user->Right,&right);CHKERRQ(info);
261: info = PetscLogFlops(67*mx*my);CHKERRQ(info);
262: return(0);
263: }
265: /* ------------------------------------------------------------------- */
268: /*
269: FormJacobian - Evaluates Jacobian matrix.
271: Input Parameters:
272: . snes - SNES context
273: . X - input vector
274: . ptr - optional user-defined context, as set by SNESSetJacobian()
276: Output Parameters:
277: . tH - Jacobian matrix
279: */
280: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat *tH, Mat *tHPre, MatStructure *flag, void *ptr)
281: {
282: AppCtx *user = (AppCtx*) ptr;
283: Mat H = *tH;
284: PetscErrorCode info;
285: PetscInt i,j,k;
286: PetscInt mx=user->mx, my=user->my;
287: MatStencil row,col[7];
288: PetscScalar hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
289: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
290: PetscScalar hl,hr,ht,hb,hc,htl,hbr;
291: PetscScalar **x, v[7];
292: PetscBool assembled;
293: PetscInt xs,xm,ys,ym;
294: Vec localX;
295: PetscScalar *top,*bottom,*left,*right;
298: /* Set various matrix options */
299: info = MatAssembled(H,&assembled);CHKERRQ(info);
300: if (assembled) {info = MatZeroEntries(H);CHKERRQ(info);}
301: *flag=SAME_NONZERO_PATTERN;
303: /* Get local vectors */
304: info = DMGetLocalVector(user->da,&localX);CHKERRQ(info);
305: info = VecGetArray(user->Top,&top);CHKERRQ(info);
306: info = VecGetArray(user->Bottom,&bottom);CHKERRQ(info);
307: info = VecGetArray(user->Left,&left);CHKERRQ(info);
308: info = VecGetArray(user->Right,&right);CHKERRQ(info);
310: /* Get ghost points */
311: info = DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
312: info = DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);CHKERRQ(info);
314: /* Get pointers to vector data */
315: info = DMDAVecGetArray(user->da,localX, &x);CHKERRQ(info);
317: info = DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(info);
318: /* Compute Jacobian over the locally owned part of the mesh */
319: for (j=ys; j< ys+ym; j++) {
320: for (i=xs; i< xs+xm; i++) {
321: xc = x[j][i];
322: xlt=xrb=xl=xr=xb=xt=xc;
324: /* Left */
325: if (i==0) {
326: xl = left[j+1];
327: xlt = left[j+2];
328: } else xl = x[j][i-1];
330: /* Bottom */
331: if (j==0) {
332: xb = bottom[i+1];
333: xrb = bottom[i+2];
334: } else xb = x[j-1][i];
336: /* Right */
337: if (i+1 == mx) {
338: xr = right[j+1];
339: xrb = right[j];
340: } else xr = x[j][i+1];
342: /* Top */
343: if (j+1==my) {
344: xt = top[i+1];
345: xlt = top[i];
346: } else xt = x[j+1][i];
348: /* Top left */
349: if (i>0 && j+1<my) xlt = x[j+1][i-1];
351: /* Bottom right */
352: if (j>0 && i+1<mx) xrb = x[j-1][i+1];
354: d1 = (xc-xl)/hx;
355: d2 = (xc-xr)/hx;
356: d3 = (xc-xt)/hy;
357: d4 = (xc-xb)/hy;
358: d5 = (xrb-xr)/hy;
359: d6 = (xrb-xb)/hx;
360: d7 = (xlt-xl)/hy;
361: d8 = (xlt-xt)/hx;
363: f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
364: f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
365: f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
366: f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
367: f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
368: f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
371: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
372: (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
373: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
374: (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
375: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
376: (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
377: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
378: (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
380: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
381: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
383: hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
384: hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
385: (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
386: (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
388: hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0;
390: k = 0;
391: row.i = i;row.j = j;
392: /* Bottom */
393: if (j>0) {
394: v[k] = hb;
395: col[k].i = i; col[k].j=j-1; k++;
396: }
398: /* Bottom right */
399: if (j>0 && i < mx -1) {
400: v[k] = hbr;
401: col[k].i = i+1; col[k].j = j-1; k++;
402: }
404: /* left */
405: if (i>0) {
406: v[k] = hl;
407: col[k].i = i-1; col[k].j = j; k++;
408: }
410: /* Centre */
411: v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++;
413: /* Right */
414: if (i < mx-1) {
415: v[k] = hr;
416: col[k].i= i+1; col[k].j = j;k++;
417: }
419: /* Top left */
420: if (i>0 && j < my-1) {
421: v[k] = htl;
422: col[k].i = i-1;col[k].j = j+1; k++;
423: }
425: /* Top */
426: if (j < my-1) {
427: v[k] = ht;
428: col[k].i = i; col[k].j = j+1; k++;
429: }
431: info = MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(info);
432: }
433: }
435: info = VecRestoreArray(user->Left,&left);CHKERRQ(info);
436: info = VecRestoreArray(user->Top,&top);CHKERRQ(info);
437: info = VecRestoreArray(user->Bottom,&bottom);CHKERRQ(info);
438: info = VecRestoreArray(user->Right,&right);CHKERRQ(info);
440: /* Assemble the matrix */
441: info = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(info);
442: info = DMDAVecRestoreArray(user->da,localX,&x);CHKERRQ(info);
443: info = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(info);
444: info = DMRestoreLocalVector(user->da,&localX);CHKERRQ(info);
446: info = PetscLogFlops(199*mx*my);CHKERRQ(info);
447: return(0);
448: }
450: /* ------------------------------------------------------------------- */
453: /*
454: MSA_BoundaryConditions - Calculates the boundary conditions for
455: the region.
457: Input Parameter:
458: . user - user-defined application context
460: Output Parameter:
461: . user - user-defined application context
462: */
463: PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
464: {
465: PetscErrorCode info;
466: PetscInt i,j,k,limit=0,maxits=5;
467: PetscInt mx=user->mx,my=user->my;
468: PetscInt xs,ys,xm,ym;
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;
475: Vec Bottom,Top,Right,Left;
476: PetscScalar scl=1.0;
477: PetscBool flg;
480: info = DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(info);
482: bsize=xm+2; lsize=ym+2; rsize=ym+2; tsize=xm+2;
484: info = VecCreateMPI(PETSC_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom);CHKERRQ(info);
485: info = VecCreateMPI(PETSC_COMM_WORLD,tsize,PETSC_DECIDE,&Top);CHKERRQ(info);
486: info = VecCreateMPI(PETSC_COMM_WORLD,lsize,PETSC_DECIDE,&Left);CHKERRQ(info);
487: info = VecCreateMPI(PETSC_COMM_WORLD,rsize,PETSC_DECIDE,&Right);CHKERRQ(info);
489: user->Top = Top;
490: user->Left = Left;
491: user->Bottom = Bottom;
492: user->Right = Right;
494: hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
496: for (j=0; j<4; j++) {
497: if (j==0) {
498: yt = b;
499: xt = l+hx*xs;
500: limit = bsize;
501: info = VecGetArray(Bottom,&boundary);CHKERRQ(info);
502: } else if (j==1) {
503: yt = t;
504: xt = l+hx*xs;
505: limit= tsize;
506: info = VecGetArray(Top,&boundary);CHKERRQ(info);
507: } else if (j==2) {
508: yt = b+hy*ys;
509: xt = l;
510: limit= lsize;
511: info = VecGetArray(Left,&boundary);CHKERRQ(info);
512: } else { /* if (j==3) */
513: yt = b+hy*ys;
514: xt = r;
515: limit= rsize;
516: info = VecGetArray(Right,&boundary);CHKERRQ(info);
517: }
519: for (i=0; i<limit; i++) {
520: u1=xt;
521: u2=-yt;
522: for (k=0; k<maxits; k++) {
523: nf1 = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
524: nf2 = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
525: fnorm = PetscSqrtReal(nf1*nf1+nf2*nf2);
526: if (fnorm <= tol) break;
527: njac11 = one+u2*u2-u1*u1;
528: njac12 = two*u1*u2;
529: njac21 = -two*u1*u2;
530: njac22 = -one - u1*u1 + u2*u2;
531: det = njac11*njac22-njac21*njac12;
532: u1 = u1-(njac22*nf1-njac12*nf2)/det;
533: u2 = u2-(njac11*nf2-njac21*nf1)/det;
534: }
536: boundary[i]=u1*u1-u2*u2;
537: if (j==0 || j==1) xt=xt+hx;
538: else yt=yt+hy; /* if (j==2 || j==3) */
539: }
541: if (j==0) {
542: info = VecRestoreArray(Bottom,&boundary);CHKERRQ(info);
543: } else if (j==1) {
544: info = VecRestoreArray(Top,&boundary);CHKERRQ(info);
545: } else if (j==2) {
546: info = VecRestoreArray(Left,&boundary);CHKERRQ(info);
547: } else if (j==3) {
548: info = VecRestoreArray(Right,&boundary);CHKERRQ(info);
549: }
551: }
553: /* Scale the boundary if desired */
555: info = PetscOptionsGetReal(NULL,"-bottom",&scl,&flg);CHKERRQ(info);
556: if (flg) {
557: info = VecScale(Bottom, scl);CHKERRQ(info);
558: }
560: info = PetscOptionsGetReal(NULL,"-top",&scl,&flg);CHKERRQ(info);
561: if (flg) {
562: info = VecScale(Top, scl);CHKERRQ(info);
563: }
565: info = PetscOptionsGetReal(NULL,"-right",&scl,&flg);CHKERRQ(info);
566: if (flg) {
567: info = VecScale(Right, scl);CHKERRQ(info);
568: }
570: info = PetscOptionsGetReal(NULL,"-left",&scl,&flg);CHKERRQ(info);
571: if (flg) {
572: info = VecScale(Left, scl);CHKERRQ(info);
573: }
574: return(0);
575: }
577: /* ------------------------------------------------------------------- */
580: /*
581: MSA_InitialPoint - Calculates the initial guess in one of three ways.
583: Input Parameters:
584: . user - user-defined application context
585: . X - vector for initial guess
587: Output Parameters:
588: . X - newly computed initial guess
589: */
590: PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
591: {
592: PetscErrorCode info;
593: PetscInt start =-1,i,j;
594: PetscScalar zero =0.0;
595: PetscBool flg;
596: PetscScalar *left,*right,*bottom,*top;
599: info = PetscOptionsGetInt(NULL,"-start",&start,&flg);CHKERRQ(info);
601: if (flg && start==0) { /* The zero vector is reasonable */
603: info = VecSet(X, zero);CHKERRQ(info);
604: /* PLogInfo(user,"Min. Surface Area Problem: Start with 0 vector \n"); */
607: } else { /* Take an average of the boundary conditions */
608: PetscInt mx=user->mx,my=user->my;
609: PetscScalar **x;
610: PetscInt xs,xm,ys,ym;
612: info = VecGetArray(user->Top,&top);CHKERRQ(info);
613: info = VecGetArray(user->Bottom,&bottom);CHKERRQ(info);
614: info = VecGetArray(user->Left,&left);CHKERRQ(info);
615: info = VecGetArray(user->Right,&right);CHKERRQ(info);
617: /* Get pointers to vector data */
618: info = DMDAVecGetArray(user->da,X,&x);CHKERRQ(info);
619: info = DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(info);
621: /* Perform local computations */
622: for (j=ys; j<ys+ym; j++) {
623: for (i=xs; i< xs+xm; i++) {
624: x[j][i] = ((j+1)*bottom[i-xs+1]/my+(my-j+1)*top[i-xs+1]/(my+2)+
625: (i+1)*left[j-ys+1]/mx+(mx-i+1)*right[j-ys+1]/(mx+2))/2.0;
626: }
627: }
629: /* Restore vectors */
630: info = DMDAVecRestoreArray(user->da,X,&x);CHKERRQ(info);
631: info = VecRestoreArray(user->Left,&left);CHKERRQ(info);
632: info = VecRestoreArray(user->Top,&top);CHKERRQ(info);
633: info = VecRestoreArray(user->Bottom,&bottom);CHKERRQ(info);
634: info = VecRestoreArray(user->Right,&right);CHKERRQ(info);
635: }
636: return(0);
637: }
641: /*
642: MSA_Plate - Calculates an obstacle for surface to stretch over.
643: */
644: PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx)
645: {
646: AppCtx *user=(AppCtx*)ctx;
647: PetscErrorCode info;
648: PetscInt i,j;
649: PetscInt xs,ys,xm,ym;
650: PetscInt mx=user->mx, my=user->my, bmy, bmx;
651: PetscScalar t1,t2,t3;
652: PetscScalar **xl;
653: PetscScalar lb=-SNES_VI_INF, ub=SNES_VI_INF;
654: PetscBool cylinder;
656: user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
657: user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
659: bmy=user->bmy, bmx=user->bmx;
661: info = DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(info);
662: info = VecSet(XL, lb);CHKERRQ(info);
663: info = DMDAVecGetArray(user->da,XL,&xl);CHKERRQ(info);
664: info = VecSet(XU, ub);CHKERRQ(info);
666: info = PetscOptionsHasName(NULL,"-cylinder",&cylinder);CHKERRQ(info);
667: /* Compute the optional lower box */
668: if (cylinder) {
669: for (i=xs; i< xs+xm; i++) {
670: for (j=ys; j<ys+ym; j++) {
671: t1=(2.0*i-mx)*bmy;
672: t2=(2.0*j-my)*bmx;
673: t3=bmx*bmx*bmy*bmy;
674: if (t1*t1 + t2*t2 <= t3) xl[j][i] = user->bheight;
675: }
676: }
677: } else {
678: /* Compute the optional lower box */
679: for (i=xs; i< xs+xm; i++) {
680: for (j=ys; j<ys+ym; j++) {
681: if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
682: j>=(my-bmy)/2 && j<my-(my-bmy)/2) {
683: xl[j][i] = user->bheight;
684: }
685: }
686: }
687: }
689: info = DMDAVecRestoreArray(user->da,XL,&xl);CHKERRQ(info);
690: return(0);
691: }