Actual source code: ex2.c
1: /*$Id: ex2.c,v 1.69 2001/04/10 19:37:04 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: double 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: Scalar *work; /* work space */
24: Vec s,y,xvec; /* work space for computing Hessian */
25: Scalar 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,double*,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,double*,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,double*,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: double 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: PetscOptionsGetDouble(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(Scalar),&user.work);
79: } else {
80: user.work = 0;
81: }
83: /* Allocate vectors */
84: VecCreate(PETSC_COMM_SELF,PETSC_DECIDE,user.ndim,&user.y);
85: VecSetFromOptions(user.y);
86: VecDuplicate(user.y,&user.s);
87: VecDuplicate(user.y,&g);
88: VecDuplicate(user.y,&x);
89: VecGetLocalSize(x,&ldim);
91: /* Create nonlinear solver */
92: SNESCreate(PETSC_COMM_SELF,SNES_UNCONSTRAINED_MINIMIZATION,&snes);
93: SNESSetType(snes,type);
95: /* Set various routines */
96: SNESSetMinimizationFunction(snes,FormMinimizationFunction,&user);
97: SNESSetGradient(snes,g,FormGradient,&user);
99: /* Form Hessian matrix approx, using one of three methods:
100: (default) : explicitly form Hessian approximation
101: -snes_mf : employ default PETSc matrix-free code
102: -my_snes_mf : employ user-defined matrix-free code (since we just happen to
103: have a routine for matrix-vector products in this example)
104: */
105: PetscOptionsHasName(PETSC_NULL,"-my_snes_mf",&flg);
106: if (flg) {
107: MatCreateShell(PETSC_COMM_SELF,ldim,user.ndim,user.ndim,user.ndim,(void*)&user,&H);
108: if (user.problem == 1) {
109: MatShellSetOperation(H,MATOP_MULT,(void(*)())HessianProductMat1);
110: } else if (user.problem == 2) {
111: MatShellSetOperation(H,MATOP_MULT,(void(*)())HessianProductMat2);
112: }
113: SNESSetHessian(snes,H,H,MatrixFreeHessian,(void *)&user);
115: /* Set null preconditioner. Alternatively, set user-provided
116: preconditioner or explicitly form preconditioning matrix */
117: SNESGetSLES(snes,&sles);
118: SLESGetPC(sles,&pc);
119: PCSetType(pc,PCNONE);
120: } else {
121: MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,user.ndim,user.ndim,&H);
122: MatSetFromOptions(H);
123: MatSetOption(H,MAT_SYMMETRIC);
124: SNESSetHessian(snes,H,H,FormHessian,(void *)&user);
125: }
127: /* Set options; then solve minimization problem */
128: SNESSetFromOptions(snes);
129: if (user.problem == 1) {
130: FormInitialGuess1(&user,x);
131: } else if (user.problem == 2) {
132: FormInitialGuess2(&user,x);
133: }
134: SNESSolve(snes,x,&its);
135: SNESGetNumberUnsuccessfulSteps(snes,&nfails);
136: SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);
137: PetscPrintf(PETSC_COMM_SELF,"number of Newton iterations = %d, ",its);
138: PetscPrintf(PETSC_COMM_SELF,"number of unsuccessful steps = %dnn",nfails);
140: /* Free data structures */
141: if (user.work) {PetscFree(user.work);}
142: VecDestroy(user.s);
143: VecDestroy(user.y);
144: VecDestroy(x);
145: VecDestroy(g);
146: MatDestroy(H);
147: SNESDestroy(snes);
149: PetscFinalize();
150: return 0;
151: }
152: /* -------------------------------------------------------------------- */
153: /*
154: FormMinimizationFunction - Evaluates function f(x).
155: */
156: int FormMinimizationFunction(SNES snes,Vec x,double *f,void *ptr)
157: {
158: AppCtx *user = (AppCtx*)ptr;
161: if (user->problem == 1) {
162: EvalFunctionGradient1(snes,x,f,NULL,FunctionEval,user);
163: } else if (user->problem == 2) {
164: EvalFunctionGradient2(snes,x,f,NULL,FunctionEval,user);
165: } else SETERRQ(1,"FormMinimizationFunction: Invalid problem number.");
166: return 0;
167: }
168: /* -------------------------------------------------------------------- */
169: /*
170: FormGradient - Evaluates gradient g(x).
171: */
172: int FormGradient(SNES snes,Vec x,Vec g,void *ptr)
173: {
174: AppCtx *user = (AppCtx*)ptr;
177: if (user->problem == 1) {
178: EvalFunctionGradient1(snes,x,NULL,g,GradientEval,user);
179: } else if (user->problem == 2) {
180: EvalFunctionGradient2(snes,x,NULL,g,GradientEval,user);
181: } else SETERRQ(1,"FormGradient: Invalid problem number.");
182: return 0;
183: }
184: /* -------------------------------------------------------------------- */
185: /*
186: FormHessian - Forms Hessian matrix by computing a column at a time.
187: */
188: int FormHessian(SNES snes,Vec X,Mat *H,Mat *PrecH,MatStructure *flag,void *ptr)
189: {
190: AppCtx *user = (AppCtx*)ptr;
191: int i,j,ierr,ndim;
192: Scalar *y,zero = 0.0,one = 1.0;
194: ndim = user->ndim;
195: VecSet(&zero,user->s);
196: user->xvec = X; /* Set location of vector */
198: MatZeroEntries(*H);
199: for (j=0; j<ndim; j++) { /* loop over columns */
201: VecSetValues(user->s,1,&j,&one,INSERT_VALUES);
202: VecAssemblyBegin(user->s);
203: VecAssemblyEnd(user->s);
205: if (user->problem == 1) {
206: HessianProduct1(ptr,user->s,user->y);
207: } else if (user->problem == 2) {
208: HessianProduct2(ptr,user->s,user->y);
209: }
211: VecSetValues(user->s,1,&j,&zero,INSERT_VALUES);
212: VecAssemblyBegin(user->s);
213: VecAssemblyEnd(user->s);
215: VecGetArray(user->y,&y);
216: for (i=0; i<ndim; i++) {
217: if (y[i] != zero) {
218: MatSetValues(*H,1,&i,1,&j,&y[i],ADD_VALUES);
219: }
220: }
221: VecRestoreArray(user->y,&y);
222: }
223: MatAssemblyBegin(*H,MAT_FINAL_ASSEMBLY);
224: MatAssemblyEnd(*H,MAT_FINAL_ASSEMBLY);
226: return 0;
227: }
228: /* -------------------------------------------------------------------- */
229: /*
230: MatrixFreeHessian
231: */
232: int MatrixFreeHessian(SNES snes,Vec X,Mat *H,Mat *PrecH,MatStructure *flag,void *ptr)
233: {
234: AppCtx *user = (AppCtx*)ptr;
235: int ierr;
236: /* Sets location of vector for use in computing matrix-vector products
237: of the form H(X)*y */
238: user->xvec = X;
239: MatAssemblyBegin(*H,MAT_FINAL_ASSEMBLY);
240: MatAssemblyEnd(*H,MAT_FINAL_ASSEMBLY);
241: MatSNESMFSetBase(*H,X);
242: return 0;
243: }
245: /* ------------------------------------------------------------------ */
246: /* Elastic-Plastic Torsion Test Problem */
247: /* ------------------------------------------------------------------ */
249: /* -------------------- Form initial approximation ----------------- */
251: int FormInitialGuess1(AppCtx *user,Vec X)
252: {
253: int ierr,i,j,k,nx = user->mx,ny = user->my;
254: Scalar hx = user->hx,hy = user->hy,temp;
255: Scalar *x;
257: VecGetArray(X,&x);
258: for (j=0; j<ny; j++) {
259: temp = ((double)PetscMin(j+1,ny-j))*hy;
260: for (i=0; i<nx; i++) {
261: k = nx*j + i;
262: #if !defined(PETSC_USE_COMPLEX)
263: x[k] = PetscMin((PetscMin(i+1,nx-i))*hx,temp);
264: #else
265: x[k] = PetscMin(PetscRealPart(((double)PetscMin(i+1,nx-i))*hx),PetscRealPart(temp));
266: #endif
267: }
268: }
269: VecRestoreArray(X,&x);
270: return 0;
271: }
272: /* ---------- Evaluate function f(x) and/or gradient g(x) ----------- */
274: int EvalFunctionGradient1(SNES snes,Vec X,double *f,Vec gvec,FctGradFlag fg,AppCtx *user)
275: {
276: int ierr,nx = user->mx,ny = user->my,ind,i,j,k;
277: Scalar hx = user->hx,hy = user->hy,area,three = 3.0,p5 = 0.5,cdiv3;
278: Scalar zero = 0.0,vb,vl,vr,vt,dvdx,dvdy,flin = 0.0,fquad = 0.0;
279: Scalar val,v,*x;
281: cdiv3 = user->param/three;
282: VecGetArray(X,&x);
284: if (fg & GradientEval) {
285: VecSet(&zero,gvec);
286: }
288: /* Compute function and gradient over the lower triangular elements */
289: for (j=-1; j<ny; j++) {
290: for (i=-1; i<nx; i++) {
291: k = nx*j + i;
292: v = zero;
293: vr = zero;
294: vt = zero;
295: if (i >= 0 && j >= 0) v = x[k];
296: if (i < nx-1 && j > -1) vr = x[k+1];
297: if (i > -1 && j < ny-1) vt = x[k+nx];
298: dvdx = (vr-v)/hx;
299: dvdy = (vt-v)/hy;
300: if (fg & FunctionEval) {
301: fquad += dvdx*dvdx + dvdy*dvdy;
302: flin -= cdiv3*(v+vr+vt);
303: }
304: if (fg & GradientEval) {
305: if (i != -1 && j != -1) {
306: ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
307: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
308: }
309: if (i != nx-1 && j != -1) {
310: ind = k+1; val = dvdx/hx - cdiv3;
311: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
312: }
313: if (i != -1 && j != ny-1) {
314: ind = k+nx; val = dvdy/hy - cdiv3;
315: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
316: }
317: }
318: }
319: }
321: /* Compute function and gradient over the upper triangular elements */
322: for (j=0; j<=ny; j++) {
323: for (i=0; i<=nx; i++) {
324: k = nx*j + i;
325: vb = zero;
326: vl = zero;
327: v = zero;
328: if (i < nx && j > 0) vb = x[k-nx];
329: if (i > 0 && j < ny) vl = x[k-1];
330: if (i < nx && j < ny) v = x[k];
331: dvdx = (v-vl)/hx;
332: dvdy = (v-vb)/hy;
333: if (fg & FunctionEval) {
334: fquad = fquad + dvdx*dvdx + dvdy*dvdy;
335: flin = flin - cdiv3*(vb+vl+v);
336: } if (fg & GradientEval) {
337: if (i != nx && j != 0) {
338: ind = k-nx; val = - dvdy/hy - cdiv3;
339: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
340: }
341: if (i != 0 && j != ny) {
342: ind = k-1; val = - dvdx/hx - cdiv3;
343: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
344: }
345: if (i != nx && j != ny) {
346: ind = k; val = dvdx/hx + dvdy/hy - cdiv3;
347: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
348: }
349: }
350: }
351: }
352: VecRestoreArray(X,&x);
353: area = p5*hx*hy;
354: if (fg & FunctionEval) { /* Scale the function */
355: #if !defined(PETSC_USE_COMPLEX)
356: *f = area*(p5*fquad+flin);
357: #else
358: *f = PetscRealPart(area*(p5*fquad+flin));
359: #endif
360: } if (fg & GradientEval) { /* Scale the gradient */
361: VecAssemblyBegin(gvec);
362: VecAssemblyEnd(gvec);
363: VecScale((Scalar*)&area,gvec);
364: }
365: return 0;
366: }
368: int HessianProductMat1(Mat mat,Vec svec,Vec y)
369: {
370: void *ptr;
371: MatShellGetContext(mat,&ptr);
372: HessianProduct1(ptr,svec,y);
373: return 0;
374: }
375:
376: /* --------------------------------------------------------------------- */
377: /*
378: HessianProduct - Computes the matrix-vector product: y = f''(x)*s
379: */
380: int HessianProduct1(void *ptr,Vec svec,Vec y)
381: {
382: AppCtx *user = (AppCtx *)ptr;
383: int nx,ny,i,j,k,ierr,ind;
384: Scalar p5 = 0.5,one = 1.0,hx,hy;
385: Scalar v,vb,vl,vr,vt,hxhx,hyhy,zero = 0.0;
386: Scalar val,area,*x,*s,szero = 0.0;
388: nx = user->mx;
389: ny = user->my;
390: hx = user->hx;
391: hy = user->hy;
393: hxhx = one/(hx*hx);
394: hyhy = one/(hy*hy);
396: VecGetArray(user->xvec,&x);
397: VecGetArray(svec,&s);
398: VecSet(&szero,y);
400: /* Compute f''(x)*s over the lower triangular elements */
401: for (j=-1; j<ny; j++) {
402: for (i=-1; i<nx; i++) {
403: k = nx*j + i;
404: v = zero;
405: vr = zero;
406: vt = zero;
407: if (i != -1 && j != -1) v = s[k];
408: if (i != nx-1 && j != -1) {
409: vr = s[k+1];
410: ind = k+1; val = hxhx*(vr-v);
411: VecSetValues(y,1,&ind,&val,ADD_VALUES);
412: }
413: if (i != -1 && j != ny-1) {
414: vt = s[k+nx];
415: ind = k+nx; val = hyhy*(vt-v);
416: VecSetValues(y,1,&ind,&val,ADD_VALUES);
417: }
418: if (i != -1 && j != -1) {
419: ind = k; val = hxhx*(v-vr) + hyhy*(v-vt);
420: VecSetValues(y,1,&ind,&val,ADD_VALUES);
421: }
422: }
423: }
425: /* Compute f''(x)*s over the upper triangular elements */
426: for (j=0; j<=ny; j++) {
427: for (i=0; i<=nx; i++) {
428: k = nx*j + i;
429: v = zero;
430: vl = zero;
431: vb = zero;
432: if (i != nx && j != ny) v = s[k];
433: if (i != nx && j != 0) {
434: vb = s[k-nx];
435: ind = k-nx; val = hyhy*(vb-v);
436: VecSetValues(y,1,&ind,&val,ADD_VALUES);
437: }
438: if (i != 0 && j != ny) {
439: vl = s[k-1];
440: ind = k-1; val = hxhx*(vl-v);
441: VecSetValues(y,1,&ind,&val,ADD_VALUES);
442: }
443: if (i != nx && j != ny) {
444: ind = k; val = hxhx*(v-vl) + hyhy*(v-vb);
445: VecSetValues(y,1,&ind,&val,ADD_VALUES);
446: }
447: }
448: }
449: VecRestoreArray(svec,&s);
450: VecRestoreArray(user->xvec,&x);
451: VecAssemblyBegin(y);
452: VecAssemblyEnd(y);
454: /* Scale result by area */
455: area = p5*hx*hy;
456: VecScale(&area,y);
457: return 0;
458: }
460: /* ------------------------------------------------------------------ */
461: /* Minimal Surface Area Test Problem */
462: /* ------------------------------------------------------------------ */
464: /* -------------------- Form initial approximation ----------------- */
466: int FormInitialGuess2(AppCtx *user,Vec X)
467: {
468: int ierr,i,j,k,nx = user->mx,ny = user->my;
469: Scalar one = 1.0,p5 = 0.5,alphaj,betai;
470: Scalar hx = user->hx,hy = user->hy,*x;
471: Scalar *bottom,*top,*left,*right,xline,yline;
473: bottom = user->work;
474: top = &user->work[nx+2];
475: left = &user->work[2*nx+4];
476: right = &user->work[2*nx+ny+6];
478: /* Compute the boundary values once only */
479: BoundaryValues(user);
480: VecGetArray(X,&x);
481: for (j=0; j<ny; j++) {
482: alphaj = ((double)(j+1))*hy;
483: for (i=0; i<nx; i++) {
484: betai = ((double)(i+1))*hx;
485: yline = alphaj*top[i+1] + (one-alphaj)*bottom[i+1];
486: xline = betai*right[j+1] + (one-betai)*left[j+1];
487: k = nx*j + i;
488: x[k] = (yline+xline)*p5;
489: }
490: }
491: VecRestoreArray(X,&x);
492: return 0;
493: }
495: /* ---------- Evaluate function f(x) and/or gradient g(x) ----------- */
497: int EvalFunctionGradient2(SNES snes,Vec X,double *f,Vec gvec,FctGradFlag fg,AppCtx *user)
498: {
499: int ierr,nx = user->mx,ny = user->my,ind,i,j,k;
500: Scalar one = 1.0,p5 = 0.5,hx = user->hx,hy = user->hy,fl,fu,area;
501: Scalar *bottom,*top,*left,*right;
502: Scalar v=0.0,vb=0.0,vl=0.0,vr=0.0,vt=0.0,dvdx,dvdy;
503: Scalar zero = 0.0,val,*x;
505: bottom = user->work;
506: top = &user->work[nx+2];
507: left = &user->work[2*nx+4];
508: right = &user->work[2*nx+ny+6];
510: VecGetArray(X,&x);
511: if (fg & FunctionEval) {
512: *f = 0.0;
513: }
514: if (fg & GradientEval) {
515: VecSet(&zero,gvec);
516: }
518: /* Compute function and gradient over the lower triangular elements */
519: for (j=-1; j<ny; j++) {
520: for (i=-1; i<nx; i++) {
521: k = nx*j + i;
522: if (i >= 0 && j >= 0) {
523: v = x[k];
524: } else {
525: if (j == -1) v = bottom[i+1];
526: if (i == -1) v = left[j+1];
527: }
528: if (i<nx-1 && j>-1) {
529: vr = x[k+1];
530: } else {
531: if (i == nx-1) vr = right[j+1];
532: if (j == -1) vr = bottom[i+2];
533: }
534: if (i>-1 && j<ny-1) {
535: vt = x[k+nx];
536: } else {
537: if (i == -1) vt = left[j+2];
538: if (j == ny-1) vt = top[i+1];
539: }
540: dvdx = (vr-v)/hx;
541: dvdy = (vt-v)/hy;
542: fl = PetscSqrtScalar(one + dvdx*dvdx + dvdy*dvdy);
543: if (fg & FunctionEval) {
544: #if !defined(PETSC_USE_COMPLEX)
545: *f += fl;
546: #else
547: *f += PetscRealPart(fl);
548: #endif
549: }
550: if (fg & GradientEval) {
551: if (i>-1 && j>-1) {
552: ind = k; val = -(dvdx/hx+dvdy/hy)/fl;
553: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
554: }
555: if (i<nx-1 && j>-1) {
556: ind = k+1; val = (dvdx/hx)/fl;
557: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
558: }
559: if (i>-1 && j<ny-1) {
560: ind = k+nx; val = (dvdy/hy)/fl;
561: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
562: }
563: }
564: }
565: }
567: /* Compute function and gradient over the upper triangular elements */
568: for (j=0; j<=ny; j++) {
569: for (i=0; i<=nx; i++) {
570: k = nx*j + i;
571: if (i<nx && j>0) {
572: vb = x[k-nx];
573: } else {
574: if (j == 0) vb = bottom[i+1];
575: if (i == nx) vb = right[j];
576: }
577: if (i>0 && j<ny) {
578: vl = x[k-1];
579: } else {
580: if (j == ny) vl = top[i];
581: if (i == 0) vl = left[j+1];
582: }
583: if (i<nx && j<ny) {
584: v = x[k];
585: } else {
586: if (i == nx) v = right[j+1];
587: if (j == ny) v = top[i+1];
588: }
589: dvdx = (v-vl)/hx;
590: dvdy = (v-vb)/hy;
591: fu = PetscSqrtScalar(one + dvdx*dvdx + dvdy*dvdy);
592: if (fg & FunctionEval) {
593: #if !defined(PETSC_USE_COMPLEX)
594: *f += fu;
595: #else
596: *f += PetscRealPart(fu);
597: #endif
598: } if (fg & GradientEval) {
599: if (i<nx && j>0) {
600: ind = k-nx; val = -(dvdy/hy)/fu;
601: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
602: }
603: if (i>0 && j<ny) {
604: ind = k-1; val = -(dvdx/hx)/fu;
605: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
606: }
607: if (i<nx && j<ny) {
608: ind = k; val = (dvdx/hx+dvdy/hy)/fu;
609: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
610: }
611: }
612: }
613: }
614: VecRestoreArray(X,&x);
615: area = p5*hx*hy;
616: if (fg & FunctionEval) { /* Scale the function */
617: #if !defined(PETSC_USE_COMPLEX)
618: *f *= area;
619: #else
620: *f *= PetscRealPart(area);
621: #endif
622: } if (fg & GradientEval) { /* Scale the gradient */
623: VecAssemblyBegin(gvec);
624: VecAssemblyEnd(gvec);
625: VecScale(&area,gvec);
626: }
627: return 0;
628: }
629: /* --------------------------------------------------------------------- */
630: int HessianProductMat2(Mat mat,Vec svec,Vec y)
631: {
632: void *ptr;
633: MatShellGetContext(mat,&ptr);
634: HessianProduct2(ptr,svec,y);
635: return 0;
636: }
638: /*
639: HessianProduct2 - Computes the matrix-vector product: y = f''(x)*s
640: */
641: int HessianProduct2(void *ptr,Vec svec,Vec y)
642: {
643: AppCtx *user = (AppCtx*)ptr;
644: int nx,ny,i,j,k,ierr,ind;
645: Scalar one = 1.0,p5 = 0.5,hx,hy;
646: Scalar dzdy,dzdyhy,fl,fl3,fu,fu3,tl,tu,z,zb,zl,zr,zt;
647: Scalar *bottom,*top,*left,*right;
648: Scalar dvdx,dvdxhx,dvdy,dvdyhy,dzdx,dzdxhx;
649: Scalar v=0.0,vb=0.0,vl=0.0,vr=0.0,vt=0.0,zerod = 0.0;
650: Scalar val,area,zero = 0.0,*s,*x;
652: nx = user->mx;
653: ny = user->my;
654: hx = user->hx;
655: hy = user->hy;
657: bottom = user->work;
658: top = &user->work[nx+2];
659: left = &user->work[2*nx+4];
660: right = &user->work[2*nx+ny+6];
662: VecGetArray(user->xvec,&x);
663: VecGetArray(svec,&s);
664: VecSet(&zero,y);
666: /* Compute f''(x)*s over the lower triangular elements */
667: for (j=-1; j<ny; j++) {
668: for (i=-1; i<nx; i++) {
669: k = nx*j + i;
670: if (i != -1 && j != -1) {
671: v = x[k];
672: z = s[k];
673: } else {
674: if (j == -1) v = bottom[i+1];
675: if (i == -1) v = left[j+1];
676: z = zerod;
677: }
678: if (i != nx-1 && j != -1) {
679: vr = x[k+1];
680: zr = s[k+1];
681: } else {
682: if (i == nx-1) vr = right[j+1];
683: if (j == -1) vr = bottom[i+2];
684: zr = zerod;
685: }
686: if (i != -1 && j != ny-1) {
687: vt = x[k+nx];
688: zt = s[k+nx];
689: } else {
690: if (i == -1) vt = left[j+2];
691: if (j == ny-1) vt = top[i+1];
692: zt = zerod;
693: }
694: dvdx = (vr-v)/hx;
695: dvdy = (vt-v)/hy;
696: dzdx = (zr-z)/hx;
697: dzdy = (zt-z)/hy;
698: dvdxhx = dvdx/hx;
699: dvdyhy = dvdy/hy;
700: dzdxhx = dzdx/hx;
701: dzdyhy = dzdy/hy;
702: tl = one + dvdx*dvdx + dvdy*dvdy;
703: fl = PetscSqrtScalar(tl);
704: fl3 = fl*tl;
705: if (i != -1 && j != -1) {
706: ind = k;
707: val = (dvdx*dzdx+dvdy*dzdy)*(dvdxhx+dvdyhy)/fl3 - (dzdxhx+dzdyhy)/fl;
708: VecSetValues(y,1,&ind,&val,ADD_VALUES);
709: }
710: if (i != nx-1 && j != -1) {
711: ind = k+1;
712: val = dzdxhx/fl - (dvdx*dzdx+dvdy*dzdy)*dvdxhx/fl3;
713: VecSetValues(y,1,&ind,&val,ADD_VALUES);
714: }
715: if (i != -1 && j != ny-1) {
716: ind = k+nx;
717: val = dzdyhy/fl - (dvdx*dzdx+dvdy*dzdy)*dvdyhy/fl3;
718: VecSetValues(y,1,&ind,&val,ADD_VALUES);
719: }
720: }
721: }
723: /* Compute f''(x)*s over the upper triangular elements */
724: for (j=0; j<=ny; j++) {
725: for (i=0; i<=nx; i++) {
726: k = nx*j + i;
727: if (i != nx && j != 0) {
728: vb = x[k-nx];
729: zb = s[k-nx];
730: } else {
731: if (j == 0) vb = bottom[i+1];
732: if (i == nx) vb = right[j];
733: zb = zerod;
734: }
735: if (i != 0 && j != ny) {
736: vl = x[k-1];
737: zl = s[k-1];
738: } else {
739: if (j == ny) vl = top[i];
740: if (i == 0) vl = left[j+1];
741: zl = zerod;
742: }
743: if (i != nx && j != ny) {
744: v = x[k];
745: z = s[k];
746: } else {
747: if (i == nx) v = right[j+1];
748: if (j == ny) v = top[i+1];
749: z = zerod;
750: }
751: dvdx = (v-vl)/hx;
752: dvdy = (v-vb)/hy;
753: dzdx = (z-zl)/hx;
754: dzdy = (z-zb)/hy;
755: dvdxhx = dvdx/hx;
756: dvdyhy = dvdy/hy;
757: dzdxhx = dzdx/hx;
758: dzdyhy = dzdy/hy;
759: tu = one + dvdx*dvdx + dvdy*dvdy;
760: fu = PetscSqrtScalar(tu);
761: fu3 = fu*tu;
762: if (i != nx && j != ny) {
763: ind = k;
764: val = (dzdxhx+dzdyhy)/fu - (dvdx*dzdx+dvdy*dzdy)*(dvdxhx+dvdyhy)/fu3;
765: VecSetValues(y,1,&ind,&val,ADD_VALUES);
766: }
767: if (i != 0 && j != ny) {
768: ind = k-1;
769: val = (dvdx*dzdx+dvdy*dzdy)*dvdxhx/fu3 - dzdxhx/fu;
770: VecSetValues(y,1,&ind,&val,ADD_VALUES);
771: }
772: if (i != nx && j != 0) {
773: ind = k-nx;
774: val = (dvdx*dzdx+dvdy*dzdy)*dvdyhy/fu3 - dzdyhy/fu;
775: VecSetValues(y,1,&ind,&val,ADD_VALUES);
776: }
777: }
778: }
779: VecRestoreArray(svec,&s);
780: VecRestoreArray(user->xvec,&x);
781: VecAssemblyBegin(y);
782: VecAssemblyEnd(y);
784: /* Scale result by area */
785: area = p5*hx*hy;
786: VecScale(&area,y);
787: return 0;
788: }
789: /* ------------------------------------------------------------------- */
790: /*
791: BoundaryValues - For Minimal Surface Area problem. Computes Enneper's
792: boundary conditions (bottom, top, left, right) which are obtained by
793: defining:
794: bv(x,y) = u**2 - v**2, where u and v are the unique solutions of
795: x = u + u*(v**2) - (u**3)/3, y = -v - (u**2)*v + (v**3)/3.
796: */
797: int BoundaryValues(AppCtx *user)
798: {
799: int maxit=5,i,j,k,limit=0,nx = user->mx,ny = user->my;
800: double three=3.0,tol=1.0e-10;
801: Scalar one=1.0,two=2.0;
802: Scalar b=-.50,t=.50,l=-.50,r=.50,det,fnorm,xt=0.0,yt=0.0;
803: Scalar nf[2],njac[2][2],u[2],hx = user->hx,hy = user->hy;
804: Scalar *bottom,*top,*left,*right;
806: bottom = user->work;
807: top = &user->work[nx+2];
808: left = &user->work[2*nx+4];
809: right = &user->work[2*nx+ny+6];
811: for (j=0; j<4; j++) {
812: switch (j) {
813: case 0:
814: yt = b; xt = l; limit = nx + 2; break;
815: case 1:
816: yt = t; xt = l; limit = nx + 2; break;
817: case 2:
818: yt = b; xt = l; limit = ny + 2; break;
819: case 3:
820: yt = b; xt = r; limit = ny + 2; break;
821: default:
822: SETERRQ(1,"BoundaryValues: Only cases 0,1,2,3 are valid");
823: }
824: /* Use Newton's method to solve xt = u + u*(v**2) - (u**3)/3,
825: yt = -v - (u**2)*v + (v**3)/3. */
827: for (i=0; i<limit; i++) {
828: u[0] = xt;
829: u[1] = -yt;
830: for (k=0; k<maxit; k++) {
831: nf[0] = u[0] + u[0]*u[1]*u[1] - PetscPowScalar(u[0],three)/three - xt;
832: nf[1] = -u[1] - u[0]*u[0]*u[1] + PetscPowScalar(u[1],three)/three - yt;
833: fnorm = PetscSqrtScalar(nf[0]*nf[0]+nf[1]*nf[1]);
834: #if !defined(PETSC_USE_COMPLEX)
835: if (fnorm <= tol) break;
836: #else
837: if (PetscRealPart(fnorm) <= tol) break;
838: #endif
839: njac[0][0] = one + u[1]*u[1] - u[0]*u[0];
840: njac[0][1] = two*u[0]*u[1];
841: njac[1][0] = -two*u[0]*u[1];
842: njac[1][1] = -one - u[0]*u[0] + u[1]*u[1];
843: det = njac[0][0]*njac[1][1] - njac[0][1]*njac[1][0];
844: u[0] -= (njac[1][1]*nf[0]-njac[0][1]*nf[1])/det;
845: u[1] -= (njac[0][0]*nf[1]-njac[1][0]*nf[0])/det;
846: }
847: switch (j) {
848: case 0:
849: bottom[i] = u[0]*u[0] - u[1]*u[1]; xt += hx; break;
850: case 1:
851: top[i] = u[0]*u[0] - u[1]*u[1]; xt += hx; break;
852: case 2:
853: left[i] = u[0]*u[0] - u[1]*u[1]; yt += hy; break;
854: case 3:
855: right[i] = u[0]*u[0] - u[1]*u[1]; yt += hy; break;
856: default:
857: SETERRQ(1,"Only cases 0,1,2,3 are valid");
858: }
859: }
860: }
861: return 0;
862: }