Actual source code: ex2.c
1: /*$Id: ex2.c,v 1.72 2001/08/07 21:31:16 bsmith Exp $*/
3: static char help[] = "Demonstrates use of the SNES package to solve unconstrained minimization problems on a single processor. These examples are based onn
4: problems from the MINPACK-2 test suite. The command line options are:n
5: -mx xg, where xg = number of grid points in the 1st coordinate directionn
6: -my yg, where yg = number of grid points in the 2nd coordinate directionn
7: -p problem_number, where the possible problem numbers are:n
8: 1: Elastic-Plastic Torsion (dept)n
9: 2: Minimal Surface Area (dmsa)n
10: -snes_mf: use matrix-free methodsn
11: -par param, where param = angle of twist per unit length (problem 1 only)nn";
13: #include petscsnes.h
15: /* User-defined application context */
16: typedef struct {
17: int problem; /* test problem number */
18: PetscReal param; /* nonlinearity parameter */
19: int mx; /* discretization in x-direction */
20: int my; /* discretization in y-direction */
21: int ndim; /* problem dimension */
22: int number; /* test problem number */
23: PetscScalar *work; /* work space */
24: Vec s,y,xvec; /* work space for computing Hessian */
25: PetscScalar hx,hy;
26: } AppCtx;
28: /* Flag to indicate evaluation of function and/or gradient */
29: typedef enum {FunctionEval=1,GradientEval=2} FctGradFlag;
31: /* General routines */
32: extern int FormHessian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
33: extern int MatrixFreeHessian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
34: extern int FormMinimizationFunction(SNES,Vec,PetscReal*,void*);
35: extern int FormGradient(SNES,Vec,Vec,void*);
37: /* For Elastic-Plastic Torsion test problem */
38: extern int HessianProduct1(void *,Vec,Vec);
39: extern int HessianProductMat1(Mat,Vec,Vec);
40: extern int FormInitialGuess1(AppCtx*,Vec);
41: extern int EvalFunctionGradient1(SNES,Vec,PetscReal*,Vec,FctGradFlag,AppCtx*);
43: /* For Minimal Surface Area test problem */
44: extern int HessianProduct2(void *,Vec,Vec);
45: extern int HessianProductMat2(Mat,Vec,Vec);
46: extern int FormInitialGuess2(AppCtx*,Vec);
47: extern int EvalFunctionGradient2(SNES,Vec,PetscReal*,Vec,FctGradFlag,AppCtx*);
48: extern int BoundaryValues(AppCtx*);
50: int main(int argc,char **argv)
51: {
52: SNES snes; /* SNES context */
53: SNESType type = SNESUMTR; /* nonlinear solution method */
54: Vec x,g; /* solution, gradient vectors */
55: Mat H; /* Hessian matrix */
56: SLES sles; /* linear solver */
57: PC pc; /* preconditioner */
58: AppCtx user; /* application context */
59: int mx=10; /* discretization of problem in x-direction */
60: int my=10; /* discretization of problem in y-direction */
61: int ierr,its,nfails,ldim;
62: PetscReal one = 1.0;
63: PetscTruth flg;
65: PetscInitialize(&argc,&argv,(char *)0,help);
67: /* Set up user-defined work space */
68: user.problem = 1;
69: PetscOptionsGetInt(PETSC_NULL,"-p",&user.problem,PETSC_NULL);
70: user.param = 5.0;
71: PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
72: if (user.problem != 1 && user.problem != 2) SETERRQ(1,"Invalid problem number");
73: PetscOptionsGetInt(PETSC_NULL,"-my",&my,PETSC_NULL);
74: PetscOptionsGetInt(PETSC_NULL,"-mx",&mx,PETSC_NULL);
75: user.ndim = mx * my; user.mx = mx; user.my = my;
76: user.hx = one/(mx+1); user.hy = one/(my+1);
77: if (user.problem == 2) {
78: PetscMalloc(2*(mx+my+4)*sizeof(PetscScalar),&user.work);
79: } else {
80: user.work = 0;
81: }
83: /* Allocate vectors */
84: VecCreate(PETSC_COMM_SELF,&user.y);
85: VecSetSizes(user.y,PETSC_DECIDE,user.ndim);
86: VecSetFromOptions(user.y);
87: VecDuplicate(user.y,&user.s);
88: VecDuplicate(user.y,&g);
89: VecDuplicate(user.y,&x);
90: VecGetLocalSize(x,&ldim);
92: /* Create nonlinear solver */
93: SNESCreate(PETSC_COMM_SELF,SNES_UNCONSTRAINED_MINIMIZATION,&snes);
94: SNESSetType(snes,type);
96: /* Set various routines */
97: SNESSetMinimizationFunction(snes,FormMinimizationFunction,&user);
98: SNESSetGradient(snes,g,FormGradient,&user);
100: /* Form Hessian matrix approx, using one of three methods:
101: (default) : explicitly form Hessian approximation
102: -snes_mf : employ default PETSc matrix-free code
103: -my_snes_mf : employ user-defined matrix-free code (since we just happen to
104: have a routine for matrix-vector products in this example)
105: */
106: PetscOptionsHasName(PETSC_NULL,"-my_snes_mf",&flg);
107: if (flg) {
108: MatCreateShell(PETSC_COMM_SELF,ldim,user.ndim,user.ndim,user.ndim,(void*)&user,&H);
109: if (user.problem == 1) {
110: MatShellSetOperation(H,MATOP_MULT,(void(*)(void))HessianProductMat1);
111: } else if (user.problem == 2) {
112: MatShellSetOperation(H,MATOP_MULT,(void(*)(void))HessianProductMat2);
113: }
114: SNESSetHessian(snes,H,H,MatrixFreeHessian,(void *)&user);
116: /* Set null preconditioner. Alternatively, set user-provided
117: preconditioner or explicitly form preconditioning matrix */
118: SNESGetSLES(snes,&sles);
119: SLESGetPC(sles,&pc);
120: PCSetType(pc,PCNONE);
121: } else {
122: MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,user.ndim,user.ndim,&H);
123: MatSetFromOptions(H);
124: MatSetOption(H,MAT_SYMMETRIC);
125: SNESSetHessian(snes,H,H,FormHessian,(void *)&user);
126: }
128: /* Set options; then solve minimization problem */
129: SNESSetFromOptions(snes);
130: if (user.problem == 1) {
131: FormInitialGuess1(&user,x);
132: } else if (user.problem == 2) {
133: FormInitialGuess2(&user,x);
134: }
135: SNESSolve(snes,x,&its);
136: SNESGetNumberUnsuccessfulSteps(snes,&nfails);
137: SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);
138: PetscPrintf(PETSC_COMM_SELF,"number of Newton iterations = %d, ",its);
139: PetscPrintf(PETSC_COMM_SELF,"number of unsuccessful steps = %dnn",nfails);
141: /* Free data structures */
142: if (user.work) {PetscFree(user.work);}
143: VecDestroy(user.s);
144: VecDestroy(user.y);
145: VecDestroy(x);
146: VecDestroy(g);
147: MatDestroy(H);
148: SNESDestroy(snes);
150: PetscFinalize();
151: return 0;
152: }
153: /* -------------------------------------------------------------------- */
154: /*
155: FormMinimizationFunction - Evaluates function f(x).
156: */
157: int FormMinimizationFunction(SNES snes,Vec x,PetscReal *f,void *ptr)
158: {
159: AppCtx *user = (AppCtx*)ptr;
162: if (user->problem == 1) {
163: EvalFunctionGradient1(snes,x,f,NULL,FunctionEval,user);
164: } else if (user->problem == 2) {
165: EvalFunctionGradient2(snes,x,f,NULL,FunctionEval,user);
166: } else SETERRQ(1,"FormMinimizationFunction: Invalid problem number.");
167: return 0;
168: }
169: /* -------------------------------------------------------------------- */
170: /*
171: FormGradient - Evaluates gradient g(x).
172: */
173: int FormGradient(SNES snes,Vec x,Vec g,void *ptr)
174: {
175: AppCtx *user = (AppCtx*)ptr;
178: if (user->problem == 1) {
179: EvalFunctionGradient1(snes,x,NULL,g,GradientEval,user);
180: } else if (user->problem == 2) {
181: EvalFunctionGradient2(snes,x,NULL,g,GradientEval,user);
182: } else SETERRQ(1,"FormGradient: Invalid problem number.");
183: return 0;
184: }
185: /* -------------------------------------------------------------------- */
186: /*
187: FormHessian - Forms Hessian matrix by computing a column at a time.
188: */
189: int FormHessian(SNES snes,Vec X,Mat *H,Mat *PrecH,MatStructure *flag,void *ptr)
190: {
191: AppCtx *user = (AppCtx*)ptr;
192: int i,j,ierr,ndim;
193: PetscScalar *y,zero = 0.0,one = 1.0;
195: ndim = user->ndim;
196: VecSet(&zero,user->s);
197: user->xvec = X; /* Set location of vector */
199: MatZeroEntries(*H);
200: for (j=0; j<ndim; j++) { /* loop over columns */
202: VecSetValues(user->s,1,&j,&one,INSERT_VALUES);
203: VecAssemblyBegin(user->s);
204: VecAssemblyEnd(user->s);
206: if (user->problem == 1) {
207: HessianProduct1(ptr,user->s,user->y);
208: } else if (user->problem == 2) {
209: HessianProduct2(ptr,user->s,user->y);
210: }
212: VecSetValues(user->s,1,&j,&zero,INSERT_VALUES);
213: VecAssemblyBegin(user->s);
214: VecAssemblyEnd(user->s);
216: VecGetArray(user->y,&y);
217: for (i=0; i<ndim; i++) {
218: if (y[i] != zero) {
219: MatSetValues(*H,1,&i,1,&j,&y[i],ADD_VALUES);
220: }
221: }
222: VecRestoreArray(user->y,&y);
223: }
224: MatAssemblyBegin(*H,MAT_FINAL_ASSEMBLY);
225: MatAssemblyEnd(*H,MAT_FINAL_ASSEMBLY);
227: return 0;
228: }
229: /* -------------------------------------------------------------------- */
230: /*
231: MatrixFreeHessian
232: */
233: int MatrixFreeHessian(SNES snes,Vec X,Mat *H,Mat *PrecH,MatStructure *flag,void *ptr)
234: {
235: AppCtx *user = (AppCtx*)ptr;
236: int ierr;
237: /* Sets location of vector for use in computing matrix-vector products
238: of the form H(X)*y */
239: user->xvec = X;
240: MatAssemblyBegin(*H,MAT_FINAL_ASSEMBLY);
241: MatAssemblyEnd(*H,MAT_FINAL_ASSEMBLY);
242: MatSNESMFSetBase(*H,X);
243: return 0;
244: }
246: /* ------------------------------------------------------------------ */
247: /* Elastic-Plastic Torsion Test Problem */
248: /* ------------------------------------------------------------------ */
250: /* -------------------- Form initial approximation ----------------- */
252: int FormInitialGuess1(AppCtx *user,Vec X)
253: {
254: int ierr,i,j,k,nx = user->mx,ny = user->my;
255: PetscScalar hx = user->hx,hy = user->hy,temp;
256: PetscScalar *x;
258: VecGetArray(X,&x);
259: for (j=0; j<ny; j++) {
260: temp = ((PetscReal)PetscMin(j+1,ny-j))*hy;
261: for (i=0; i<nx; i++) {
262: k = nx*j + i;
263: #if !defined(PETSC_USE_COMPLEX)
264: x[k] = PetscMin((PetscMin(i+1,nx-i))*hx,temp);
265: #else
266: x[k] = PetscMin(PetscRealPart(((PetscReal)PetscMin(i+1,nx-i))*hx),PetscRealPart(temp));
267: #endif
268: }
269: }
270: VecRestoreArray(X,&x);
271: return 0;
272: }
273: /* ---------- Evaluate function f(x) and/or gradient g(x) ----------- */
275: int EvalFunctionGradient1(SNES snes,Vec X,PetscReal *f,Vec gvec,FctGradFlag fg,AppCtx *user)
276: {
277: int ierr,nx = user->mx,ny = user->my,ind,i,j,k;
278: PetscScalar hx = user->hx,hy = user->hy,area,three = 3.0,p5 = 0.5,cdiv3;
279: PetscScalar zero = 0.0,vb,vl,vr,vt,dvdx,dvdy,flin = 0.0,fquad = 0.0;
280: PetscScalar val,v,*x;
282: cdiv3 = user->param/three;
283: VecGetArray(X,&x);
285: if (fg & GradientEval) {
286: VecSet(&zero,gvec);
287: }
289: /* Compute function and gradient over the lower triangular elements */
290: for (j=-1; j<ny; j++) {
291: for (i=-1; i<nx; i++) {
292: k = nx*j + i;
293: v = zero;
294: vr = zero;
295: vt = zero;
296: if (i >= 0 && j >= 0) v = x[k];
297: if (i < nx-1 && j > -1) vr = x[k+1];
298: if (i > -1 && j < ny-1) vt = x[k+nx];
299: dvdx = (vr-v)/hx;
300: dvdy = (vt-v)/hy;
301: if (fg & FunctionEval) {
302: fquad += dvdx*dvdx + dvdy*dvdy;
303: flin -= cdiv3*(v+vr+vt);
304: }
305: if (fg & GradientEval) {
306: if (i != -1 && j != -1) {
307: ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
308: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
309: }
310: if (i != nx-1 && j != -1) {
311: ind = k+1; val = dvdx/hx - cdiv3;
312: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
313: }
314: if (i != -1 && j != ny-1) {
315: ind = k+nx; val = dvdy/hy - cdiv3;
316: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
317: }
318: }
319: }
320: }
322: /* Compute function and gradient over the upper triangular elements */
323: for (j=0; j<=ny; j++) {
324: for (i=0; i<=nx; i++) {
325: k = nx*j + i;
326: vb = zero;
327: vl = zero;
328: v = zero;
329: if (i < nx && j > 0) vb = x[k-nx];
330: if (i > 0 && j < ny) vl = x[k-1];
331: if (i < nx && j < ny) v = x[k];
332: dvdx = (v-vl)/hx;
333: dvdy = (v-vb)/hy;
334: if (fg & FunctionEval) {
335: fquad = fquad + dvdx*dvdx + dvdy*dvdy;
336: flin = flin - cdiv3*(vb+vl+v);
337: } if (fg & GradientEval) {
338: if (i != nx && j != 0) {
339: ind = k-nx; val = - dvdy/hy - cdiv3;
340: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
341: }
342: if (i != 0 && j != ny) {
343: ind = k-1; val = - dvdx/hx - cdiv3;
344: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
345: }
346: if (i != nx && j != ny) {
347: ind = k; val = dvdx/hx + dvdy/hy - cdiv3;
348: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
349: }
350: }
351: }
352: }
353: VecRestoreArray(X,&x);
354: area = p5*hx*hy;
355: if (fg & FunctionEval) { /* Scale the function */
356: #if !defined(PETSC_USE_COMPLEX)
357: *f = area*(p5*fquad+flin);
358: #else
359: *f = PetscRealPart(area*(p5*fquad+flin));
360: #endif
361: } if (fg & GradientEval) { /* Scale the gradient */
362: VecAssemblyBegin(gvec);
363: VecAssemblyEnd(gvec);
364: VecScale((PetscScalar*)&area,gvec);
365: }
366: return 0;
367: }
369: int HessianProductMat1(Mat mat,Vec svec,Vec y)
370: {
371: void *ptr;
372: MatShellGetContext(mat,&ptr);
373: HessianProduct1(ptr,svec,y);
374: return 0;
375: }
376:
377: /* --------------------------------------------------------------------- */
378: /*
379: HessianProduct - Computes the matrix-vector product: y = f''(x)*s
380: */
381: int HessianProduct1(void *ptr,Vec svec,Vec y)
382: {
383: AppCtx *user = (AppCtx *)ptr;
384: int nx,ny,i,j,k,ierr,ind;
385: PetscScalar p5 = 0.5,one = 1.0,hx,hy;
386: PetscScalar v,vb,vl,vr,vt,hxhx,hyhy,zero = 0.0;
387: PetscScalar val,area,*x,*s,szero = 0.0;
389: nx = user->mx;
390: ny = user->my;
391: hx = user->hx;
392: hy = user->hy;
394: hxhx = one/(hx*hx);
395: hyhy = one/(hy*hy);
397: VecGetArray(user->xvec,&x);
398: VecGetArray(svec,&s);
399: VecSet(&szero,y);
401: /* Compute f''(x)*s over the lower triangular elements */
402: for (j=-1; j<ny; j++) {
403: for (i=-1; i<nx; i++) {
404: k = nx*j + i;
405: v = zero;
406: vr = zero;
407: vt = zero;
408: if (i != -1 && j != -1) v = s[k];
409: if (i != nx-1 && j != -1) {
410: vr = s[k+1];
411: ind = k+1; val = hxhx*(vr-v);
412: VecSetValues(y,1,&ind,&val,ADD_VALUES);
413: }
414: if (i != -1 && j != ny-1) {
415: vt = s[k+nx];
416: ind = k+nx; val = hyhy*(vt-v);
417: VecSetValues(y,1,&ind,&val,ADD_VALUES);
418: }
419: if (i != -1 && j != -1) {
420: ind = k; val = hxhx*(v-vr) + hyhy*(v-vt);
421: VecSetValues(y,1,&ind,&val,ADD_VALUES);
422: }
423: }
424: }
426: /* Compute f''(x)*s over the upper triangular elements */
427: for (j=0; j<=ny; j++) {
428: for (i=0; i<=nx; i++) {
429: k = nx*j + i;
430: v = zero;
431: vl = zero;
432: vb = zero;
433: if (i != nx && j != ny) v = s[k];
434: if (i != nx && j != 0) {
435: vb = s[k-nx];
436: ind = k-nx; val = hyhy*(vb-v);
437: VecSetValues(y,1,&ind,&val,ADD_VALUES);
438: }
439: if (i != 0 && j != ny) {
440: vl = s[k-1];
441: ind = k-1; val = hxhx*(vl-v);
442: VecSetValues(y,1,&ind,&val,ADD_VALUES);
443: }
444: if (i != nx && j != ny) {
445: ind = k; val = hxhx*(v-vl) + hyhy*(v-vb);
446: VecSetValues(y,1,&ind,&val,ADD_VALUES);
447: }
448: }
449: }
450: VecRestoreArray(svec,&s);
451: VecRestoreArray(user->xvec,&x);
452: VecAssemblyBegin(y);
453: VecAssemblyEnd(y);
455: /* Scale result by area */
456: area = p5*hx*hy;
457: VecScale(&area,y);
458: return 0;
459: }
461: /* ------------------------------------------------------------------ */
462: /* Minimal Surface Area Test Problem */
463: /* ------------------------------------------------------------------ */
465: /* -------------------- Form initial approximation ----------------- */
467: int FormInitialGuess2(AppCtx *user,Vec X)
468: {
469: int ierr,i,j,k,nx = user->mx,ny = user->my;
470: PetscScalar one = 1.0,p5 = 0.5,alphaj,betai;
471: PetscScalar hx = user->hx,hy = user->hy,*x;
472: PetscScalar *bottom,*top,*left,*right,xline,yline;
474: bottom = user->work;
475: top = &user->work[nx+2];
476: left = &user->work[2*nx+4];
477: right = &user->work[2*nx+ny+6];
479: /* Compute the boundary values once only */
480: BoundaryValues(user);
481: VecGetArray(X,&x);
482: for (j=0; j<ny; j++) {
483: alphaj = ((PetscReal)(j+1))*hy;
484: for (i=0; i<nx; i++) {
485: betai = ((PetscReal)(i+1))*hx;
486: yline = alphaj*top[i+1] + (one-alphaj)*bottom[i+1];
487: xline = betai*right[j+1] + (one-betai)*left[j+1];
488: k = nx*j + i;
489: x[k] = (yline+xline)*p5;
490: }
491: }
492: VecRestoreArray(X,&x);
493: return 0;
494: }
496: /* ---------- Evaluate function f(x) and/or gradient g(x) ----------- */
498: int EvalFunctionGradient2(SNES snes,Vec X,PetscReal *f,Vec gvec,FctGradFlag fg,AppCtx *user)
499: {
500: int ierr,nx = user->mx,ny = user->my,ind,i,j,k;
501: PetscScalar one = 1.0,p5 = 0.5,hx = user->hx,hy = user->hy,fl,fu,area;
502: PetscScalar *bottom,*top,*left,*right;
503: PetscScalar v=0.0,vb=0.0,vl=0.0,vr=0.0,vt=0.0,dvdx,dvdy;
504: PetscScalar zero = 0.0,val,*x;
506: bottom = user->work;
507: top = &user->work[nx+2];
508: left = &user->work[2*nx+4];
509: right = &user->work[2*nx+ny+6];
511: VecGetArray(X,&x);
512: if (fg & FunctionEval) {
513: *f = 0.0;
514: }
515: if (fg & GradientEval) {
516: VecSet(&zero,gvec);
517: }
519: /* Compute function and gradient over the lower triangular elements */
520: for (j=-1; j<ny; j++) {
521: for (i=-1; i<nx; i++) {
522: k = nx*j + i;
523: if (i >= 0 && j >= 0) {
524: v = x[k];
525: } else {
526: if (j == -1) v = bottom[i+1];
527: if (i == -1) v = left[j+1];
528: }
529: if (i<nx-1 && j>-1) {
530: vr = x[k+1];
531: } else {
532: if (i == nx-1) vr = right[j+1];
533: if (j == -1) vr = bottom[i+2];
534: }
535: if (i>-1 && j<ny-1) {
536: vt = x[k+nx];
537: } else {
538: if (i == -1) vt = left[j+2];
539: if (j == ny-1) vt = top[i+1];
540: }
541: dvdx = (vr-v)/hx;
542: dvdy = (vt-v)/hy;
543: fl = PetscSqrtScalar(one + dvdx*dvdx + dvdy*dvdy);
544: if (fg & FunctionEval) {
545: #if !defined(PETSC_USE_COMPLEX)
546: *f += fl;
547: #else
548: *f += PetscRealPart(fl);
549: #endif
550: }
551: if (fg & GradientEval) {
552: if (i>-1 && j>-1) {
553: ind = k; val = -(dvdx/hx+dvdy/hy)/fl;
554: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
555: }
556: if (i<nx-1 && j>-1) {
557: ind = k+1; val = (dvdx/hx)/fl;
558: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
559: }
560: if (i>-1 && j<ny-1) {
561: ind = k+nx; val = (dvdy/hy)/fl;
562: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
563: }
564: }
565: }
566: }
568: /* Compute function and gradient over the upper triangular elements */
569: for (j=0; j<=ny; j++) {
570: for (i=0; i<=nx; i++) {
571: k = nx*j + i;
572: if (i<nx && j>0) {
573: vb = x[k-nx];
574: } else {
575: if (j == 0) vb = bottom[i+1];
576: if (i == nx) vb = right[j];
577: }
578: if (i>0 && j<ny) {
579: vl = x[k-1];
580: } else {
581: if (j == ny) vl = top[i];
582: if (i == 0) vl = left[j+1];
583: }
584: if (i<nx && j<ny) {
585: v = x[k];
586: } else {
587: if (i == nx) v = right[j+1];
588: if (j == ny) v = top[i+1];
589: }
590: dvdx = (v-vl)/hx;
591: dvdy = (v-vb)/hy;
592: fu = PetscSqrtScalar(one + dvdx*dvdx + dvdy*dvdy);
593: if (fg & FunctionEval) {
594: #if !defined(PETSC_USE_COMPLEX)
595: *f += fu;
596: #else
597: *f += PetscRealPart(fu);
598: #endif
599: } if (fg & GradientEval) {
600: if (i<nx && j>0) {
601: ind = k-nx; val = -(dvdy/hy)/fu;
602: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
603: }
604: if (i>0 && j<ny) {
605: ind = k-1; val = -(dvdx/hx)/fu;
606: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
607: }
608: if (i<nx && j<ny) {
609: ind = k; val = (dvdx/hx+dvdy/hy)/fu;
610: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
611: }
612: }
613: }
614: }
615: VecRestoreArray(X,&x);
616: area = p5*hx*hy;
617: if (fg & FunctionEval) { /* Scale the function */
618: #if !defined(PETSC_USE_COMPLEX)
619: *f *= area;
620: #else
621: *f *= PetscRealPart(area);
622: #endif
623: } if (fg & GradientEval) { /* Scale the gradient */
624: VecAssemblyBegin(gvec);
625: VecAssemblyEnd(gvec);
626: VecScale(&area,gvec);
627: }
628: return 0;
629: }
630: /* --------------------------------------------------------------------- */
631: int HessianProductMat2(Mat mat,Vec svec,Vec y)
632: {
633: void *ptr;
634: MatShellGetContext(mat,&ptr);
635: HessianProduct2(ptr,svec,y);
636: return 0;
637: }
639: /*
640: HessianProduct2 - Computes the matrix-vector product: y = f''(x)*s
641: */
642: int HessianProduct2(void *ptr,Vec svec,Vec y)
643: {
644: AppCtx *user = (AppCtx*)ptr;
645: int nx,ny,i,j,k,ierr,ind;
646: PetscScalar one = 1.0,p5 = 0.5,hx,hy;
647: PetscScalar dzdy,dzdyhy,fl,fl3,fu,fu3,tl,tu,z,zb,zl,zr,zt;
648: PetscScalar *bottom,*top,*left,*right;
649: PetscScalar dvdx,dvdxhx,dvdy,dvdyhy,dzdx,dzdxhx;
650: PetscScalar v=0.0,vb=0.0,vl=0.0,vr=0.0,vt=0.0,zerod = 0.0;
651: PetscScalar val,area,zero = 0.0,*s,*x;
653: nx = user->mx;
654: ny = user->my;
655: hx = user->hx;
656: hy = user->hy;
658: bottom = user->work;
659: top = &user->work[nx+2];
660: left = &user->work[2*nx+4];
661: right = &user->work[2*nx+ny+6];
663: VecGetArray(user->xvec,&x);
664: VecGetArray(svec,&s);
665: VecSet(&zero,y);
667: /* Compute f''(x)*s over the lower triangular elements */
668: for (j=-1; j<ny; j++) {
669: for (i=-1; i<nx; i++) {
670: k = nx*j + i;
671: if (i != -1 && j != -1) {
672: v = x[k];
673: z = s[k];
674: } else {
675: if (j == -1) v = bottom[i+1];
676: if (i == -1) v = left[j+1];
677: z = zerod;
678: }
679: if (i != nx-1 && j != -1) {
680: vr = x[k+1];
681: zr = s[k+1];
682: } else {
683: if (i == nx-1) vr = right[j+1];
684: if (j == -1) vr = bottom[i+2];
685: zr = zerod;
686: }
687: if (i != -1 && j != ny-1) {
688: vt = x[k+nx];
689: zt = s[k+nx];
690: } else {
691: if (i == -1) vt = left[j+2];
692: if (j == ny-1) vt = top[i+1];
693: zt = zerod;
694: }
695: dvdx = (vr-v)/hx;
696: dvdy = (vt-v)/hy;
697: dzdx = (zr-z)/hx;
698: dzdy = (zt-z)/hy;
699: dvdxhx = dvdx/hx;
700: dvdyhy = dvdy/hy;
701: dzdxhx = dzdx/hx;
702: dzdyhy = dzdy/hy;
703: tl = one + dvdx*dvdx + dvdy*dvdy;
704: fl = PetscSqrtScalar(tl);
705: fl3 = fl*tl;
706: if (i != -1 && j != -1) {
707: ind = k;
708: val = (dvdx*dzdx+dvdy*dzdy)*(dvdxhx+dvdyhy)/fl3 - (dzdxhx+dzdyhy)/fl;
709: VecSetValues(y,1,&ind,&val,ADD_VALUES);
710: }
711: if (i != nx-1 && j != -1) {
712: ind = k+1;
713: val = dzdxhx/fl - (dvdx*dzdx+dvdy*dzdy)*dvdxhx/fl3;
714: VecSetValues(y,1,&ind,&val,ADD_VALUES);
715: }
716: if (i != -1 && j != ny-1) {
717: ind = k+nx;
718: val = dzdyhy/fl - (dvdx*dzdx+dvdy*dzdy)*dvdyhy/fl3;
719: VecSetValues(y,1,&ind,&val,ADD_VALUES);
720: }
721: }
722: }
724: /* Compute f''(x)*s over the upper triangular elements */
725: for (j=0; j<=ny; j++) {
726: for (i=0; i<=nx; i++) {
727: k = nx*j + i;
728: if (i != nx && j != 0) {
729: vb = x[k-nx];
730: zb = s[k-nx];
731: } else {
732: if (j == 0) vb = bottom[i+1];
733: if (i == nx) vb = right[j];
734: zb = zerod;
735: }
736: if (i != 0 && j != ny) {
737: vl = x[k-1];
738: zl = s[k-1];
739: } else {
740: if (j == ny) vl = top[i];
741: if (i == 0) vl = left[j+1];
742: zl = zerod;
743: }
744: if (i != nx && j != ny) {
745: v = x[k];
746: z = s[k];
747: } else {
748: if (i == nx) v = right[j+1];
749: if (j == ny) v = top[i+1];
750: z = zerod;
751: }
752: dvdx = (v-vl)/hx;
753: dvdy = (v-vb)/hy;
754: dzdx = (z-zl)/hx;
755: dzdy = (z-zb)/hy;
756: dvdxhx = dvdx/hx;
757: dvdyhy = dvdy/hy;
758: dzdxhx = dzdx/hx;
759: dzdyhy = dzdy/hy;
760: tu = one + dvdx*dvdx + dvdy*dvdy;
761: fu = PetscSqrtScalar(tu);
762: fu3 = fu*tu;
763: if (i != nx && j != ny) {
764: ind = k;
765: val = (dzdxhx+dzdyhy)/fu - (dvdx*dzdx+dvdy*dzdy)*(dvdxhx+dvdyhy)/fu3;
766: VecSetValues(y,1,&ind,&val,ADD_VALUES);
767: }
768: if (i != 0 && j != ny) {
769: ind = k-1;
770: val = (dvdx*dzdx+dvdy*dzdy)*dvdxhx/fu3 - dzdxhx/fu;
771: VecSetValues(y,1,&ind,&val,ADD_VALUES);
772: }
773: if (i != nx && j != 0) {
774: ind = k-nx;
775: val = (dvdx*dzdx+dvdy*dzdy)*dvdyhy/fu3 - dzdyhy/fu;
776: VecSetValues(y,1,&ind,&val,ADD_VALUES);
777: }
778: }
779: }
780: VecRestoreArray(svec,&s);
781: VecRestoreArray(user->xvec,&x);
782: VecAssemblyBegin(y);
783: VecAssemblyEnd(y);
785: /* Scale result by area */
786: area = p5*hx*hy;
787: VecScale(&area,y);
788: return 0;
789: }
790: /* ------------------------------------------------------------------- */
791: /*
792: BoundaryValues - For Minimal Surface Area problem. Computes Enneper's
793: boundary conditions (bottom, top, left, right) which are obtained by
794: defining:
795: bv(x,y) = u**2 - v**2, where u and v are the unique solutions of
796: x = u + u*(v**2) - (u**3)/3, y = -v - (u**2)*v + (v**3)/3.
797: */
798: int BoundaryValues(AppCtx *user)
799: {
800: int maxit=5,i,j,k,limit=0,nx = user->mx,ny = user->my;
801: PetscReal three=3.0,tol=1.0e-10;
802: PetscScalar one=1.0,two=2.0;
803: PetscScalar b=-.50,t=.50,l=-.50,r=.50,det,fnorm,xt=0.0,yt=0.0;
804: PetscScalar nf[2],njac[2][2],u[2],hx = user->hx,hy = user->hy;
805: PetscScalar *bottom,*top,*left,*right;
807: bottom = user->work;
808: top = &user->work[nx+2];
809: left = &user->work[2*nx+4];
810: right = &user->work[2*nx+ny+6];
812: for (j=0; j<4; j++) {
813: switch (j) {
814: case 0:
815: yt = b; xt = l; limit = nx + 2; break;
816: case 1:
817: yt = t; xt = l; limit = nx + 2; break;
818: case 2:
819: yt = b; xt = l; limit = ny + 2; break;
820: case 3:
821: yt = b; xt = r; limit = ny + 2; break;
822: default:
823: SETERRQ(1,"BoundaryValues: Only cases 0,1,2,3 are valid");
824: }
825: /* Use Newton's method to solve xt = u + u*(v**2) - (u**3)/3,
826: yt = -v - (u**2)*v + (v**3)/3. */
828: for (i=0; i<limit; i++) {
829: u[0] = xt;
830: u[1] = -yt;
831: for (k=0; k<maxit; k++) {
832: nf[0] = u[0] + u[0]*u[1]*u[1] - PetscPowScalar(u[0],three)/three - xt;
833: nf[1] = -u[1] - u[0]*u[0]*u[1] + PetscPowScalar(u[1],three)/three - yt;
834: fnorm = PetscSqrtScalar(nf[0]*nf[0]+nf[1]*nf[1]);
835: #if !defined(PETSC_USE_COMPLEX)
836: if (fnorm <= tol) break;
837: #else
838: if (PetscRealPart(fnorm) <= tol) break;
839: #endif
840: njac[0][0] = one + u[1]*u[1] - u[0]*u[0];
841: njac[0][1] = two*u[0]*u[1];
842: njac[1][0] = -two*u[0]*u[1];
843: njac[1][1] = -one - u[0]*u[0] + u[1]*u[1];
844: det = njac[0][0]*njac[1][1] - njac[0][1]*njac[1][0];
845: u[0] -= (njac[1][1]*nf[0]-njac[0][1]*nf[1])/det;
846: u[1] -= (njac[0][0]*nf[1]-njac[1][0]*nf[0])/det;
847: }
848: switch (j) {
849: case 0:
850: bottom[i] = u[0]*u[0] - u[1]*u[1]; xt += hx; break;
851: case 1:
852: top[i] = u[0]*u[0] - u[1]*u[1]; xt += hx; break;
853: case 2:
854: left[i] = u[0]*u[0] - u[1]*u[1]; yt += hy; break;
855: case 3:
856: right[i] = u[0]*u[0] - u[1]*u[1]; yt += hy; break;
857: default:
858: SETERRQ(1,"Only cases 0,1,2,3 are valid");
859: }
860: }
861: }
862: return 0;
863: }