Actual source code: ex3.c
1: /*$Id: ex3.c,v 1.68 2001/04/10 19:37:04 bsmith Exp $*/
3: static char help[] = "Demonstrates use of the SNES package to solve unconstrained minimization problems in parallel. This example is based on then
4: Elastic-Plastic Torsion (dept) problem from the MINPACK-2 test suite.n
5: The command line options are:n
6: -mx xg, where xg = number of grid points in the 1st coordinate directionn
7: -my yg, where yg = number of grid points in the 2nd coordinate directionn
8: -par param, where param = angle of twist per unit lengthn
9: -snes_mf: use matrix-free methodsn
10: -defaultH: use default finite difference approximation of Hessiannn";
12: #include "petscsnes.h"
13: #include "petscda.h"
15: /* User-defined application context */
16: typedef struct {
17: double param; /* nonlinearity parameter */
18: int mx; /* discretization in x-direction */
19: int my; /* discretization in y-direction */
20: int ndim; /* problem dimension */
21: int number; /* test problem number */
22: Vec s,y,xvec; /* work space for computing Hessian */
23: Scalar hx,hy;
24: Vec localX,localS; /* ghosted local vector */
25: DA da; /* distributed array data structure */
26: } AppCtx;
28: /* Flag to indicate evaluation of function and/or gradient */
29: typedef enum {FunctionEval=1,GradientEval=2} FctGradFlag;
31: extern int FormHessian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
32: extern int MatrixFreeHessian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
33: extern int FormMinimizationFunction(SNES,Vec,double*,void*);
34: extern int FormGradient(SNES,Vec,Vec,void*);
35: extern int HessianProduct(void *,Vec,Vec);
36: extern int HessianProductMat(Mat,Vec,Vec);
37: extern int FormInitialGuess(AppCtx*,Vec);
38: extern int EvalFunctionGradient(SNES,Vec,double*,Vec,FctGradFlag,AppCtx*);
40: int main(int argc,char **argv)
41: {
42: SNES snes; /* SNES context */
43: SNESType type = SNESUMTR; /* nonlinear solution method */
44: Vec x,g; /* solution, gradient vectors */
45: Mat H; /* Hessian matrix */
46: AppCtx user; /* application context */
47: int mx=10; /* discretization in x-direction */
48: int my=10; /* discretization in y-direction */
49: int Nx=PETSC_DECIDE; /* processors in x-direction */
50: int Ny=PETSC_DECIDE; /* processors in y-direction */
51: int ierr,its,ldim,nfails,size;
52: PetscTruth flg;
53: double one = 1.0;
54: SLES sles;
55: PC pc;
57: PetscInitialize(&argc,&argv,(char *)0,help);
58: MPI_Comm_size(PETSC_COMM_WORLD,&size);
59: PetscOptionsGetInt(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
60: PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
61: if (Nx*Ny != size && (Nx != PETSC_DECIDE && Ny != PETSC_DECIDE))
62: SETERRQ(1,"Incompatible number of processors: Nx * Ny != size");
64: /* Set up user-defined work space */
65: user.param = 5.0;
66: PetscOptionsGetDouble(PETSC_NULL,"-par",&user.param,PETSC_NULL);
67: PetscOptionsGetInt(PETSC_NULL,"-my",&my,PETSC_NULL);
68: PetscOptionsGetInt(PETSC_NULL,"-mx",&mx,PETSC_NULL);
69: user.ndim = mx * my; user.mx = mx; user.my = my;
70: user.hx = one/(mx+1); user.hy = one/(my+1);
72: /* Set up distributed array and vectors */
73: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,user.mx,
74: user.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.da);
75: DACreateGlobalVector(user.da,&x);
76: DACreateLocalVector(user.da,&user.localX);
77: VecDuplicate(x,&user.s);
78: VecDuplicate(user.localX,&user.localS);
79: VecDuplicate(x,&g);
80: VecDuplicate(g,&user.y);
82: /* Create nonlinear solver */
83: SNESCreate(PETSC_COMM_WORLD,SNES_UNCONSTRAINED_MINIMIZATION,&snes);
84: SNESSetType(snes,type);
86: /* Set various routines */
87: SNESSetMinimizationFunction(snes,FormMinimizationFunction,
88: (void *)&user);
89: SNESSetGradient(snes,g,FormGradient,(void *)&user);
91: /* Either explicitly form Hessian matrix approx or use matrix-free version */
92: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&flg);
93: if (flg) {
94: VecGetLocalSize(x,&ldim);
95: MatCreateShell(PETSC_COMM_WORLD,ldim,user.ndim,user.ndim,user.ndim,
96: (void*)&user,&H);
97: MatShellSetOperation(H,MATOP_MULT,(void(*)())HessianProductMat);
98: SNESSetHessian(snes,H,H,MatrixFreeHessian,(void *)&user);
100: /* Set null preconditioner. Alternatively, set user-provided
101: preconditioner or explicitly form preconditioning matrix */
102: SNESGetSLES(snes,&sles);
103: SLESGetPC(sles,&pc);
104: PCSetType(pc,PCNONE);
105: } else {
106: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,user.ndim,user.ndim,&H);
107: MatSetFromOptions(H);
108: MatSetOption(H,MAT_SYMMETRIC);
109: PetscOptionsHasName(PETSC_NULL,"-defaultH",&flg);
110: if (flg) SNESSetHessian(snes,H,H,SNESDefaultComputeHessian,(void *)&user);
111: else SNESSetHessian(snes,H,H,FormHessian,(void *)&user);
112: }
114: /* Set options; then solve minimization problem */
115: SNESSetFromOptions(snes);
116: FormInitialGuess(&user,x);
117: SNESSolve(snes,x,&its);
118: SNESGetNumberUnsuccessfulSteps(snes,&nfails);
119: SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);
120: PetscPrintf(PETSC_COMM_WORLD,"number of Newton iterations = %d, ",its);
121: PetscPrintf(PETSC_COMM_WORLD,"number of unsuccessful steps = %dnn",nfails);
123: /* Free data structures */
124: VecDestroy(user.s);
125: VecDestroy(user.localX);
126: VecDestroy(user.localS);
127: VecDestroy(user.y);
128: VecDestroy(x);
129: VecDestroy(g);
130: MatDestroy(H);
131: SNESDestroy(snes);
132: DADestroy(user.da);
134: PetscFinalize();
135: return 0;
136: }
137: /* -------------------------------------------------------------------- */
138: /*
139: FormMinimizationFunction - Evaluates function f(x).
140: */
141: int FormMinimizationFunction(SNES snes,Vec x,double *f,void *ptr)
142: {
143: AppCtx *user = (AppCtx*)ptr;
144: return EvalFunctionGradient(snes,x,f,NULL,FunctionEval,user);
145: }
146: /* -------------------------------------------------------------------- */
147: /*
148: FormGradient - Evaluates gradient g(x).
149: */
150: int FormGradient(SNES snes,Vec x,Vec g,void *ptr)
151: {
152: AppCtx *user = (AppCtx*)ptr;
153: return EvalFunctionGradient(snes,x,NULL,g,GradientEval,user);
154: }
155: /* -------------------------------------------------------------------- */
156: /*
157: FormHessian - Forms Hessian matrix by computing a column at a time.
158: */
159: int FormHessian(SNES snes,Vec X,Mat *H,Mat *PrecH,MatStructure *flag,
160: void *ptr)
161: {
162: AppCtx *user = (AppCtx*)ptr;
163: int i,j,ierr,ndim,xs,ys, xm,ym,rstart,rend,ldim,iglob;
164: Scalar *y,zero = 0.0,one = 1.0;
166: MatZeroEntries(*H);
167: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
169: ndim = user->ndim;
170: VecSet(&zero,user->s);
171: user->xvec = X; /* Set location of vector */
172: VecGetOwnershipRange(user->y,&rstart,&rend);
173: VecGetLocalSize(user->y,&ldim);
175: for (j=0; j<ndim; j++) { /* loop over columns */
177: VecSetValues(user->s,1,&j,&one,INSERT_VALUES);
178: VecAssemblyBegin(user->s);
179: VecAssemblyEnd(user->s);
181: HessianProduct(ptr,user->s,user->y);
183: VecSetValues(user->s,1,&j,&zero,INSERT_VALUES);
184: VecAssemblyBegin(user->s);
185: VecAssemblyEnd(user->s);
187: VecGetArray(user->y,&y);
188: for (i=0; i<ldim; i++) {
189: if (y[i] != zero) {
190: iglob = i+rstart;
191: MatSetValues(*H,1,&iglob,1,&j,&y[i],ADD_VALUES);
192: }
193: }
194: VecRestoreArray(user->y,&y);
195: }
196: MatAssemblyBegin(*H,MAT_FINAL_ASSEMBLY);
197: MatAssemblyEnd(*H,MAT_FINAL_ASSEMBLY);
199: return 0;
200: }
201: /* -------------------------------------------------------------------- */
202: /*
203: MatrixFreeHessian
204: */
205: int MatrixFreeHessian(SNES snes,Vec X,Mat *H,Mat *PrecH,MatStructure *flag,
206: void *ptr)
207: {
208: AppCtx *user = (AppCtx*)ptr;
210: /* Sets location of vector for use in computing matrix-vector products
211: of the form H(X)*y */
212: user->xvec = X;
213: return 0;
214: }
216: /* ------------------------------------------------------------------ */
217: /* Elastic-Plastic Torsion Test Problem */
218: /* ------------------------------------------------------------------ */
220: /* -------------------- Form initial approximation ----------------- */
222: int FormInitialGuess(AppCtx *user,Vec X)
223: {
224: int ierr,i,j,k,nx = user->mx,ny = user->my;
225: Scalar hx = user->hx,hy = user->hy,temp,*x;
226: int xs,ys,xm,ym,Xm,Ym,Xs,Ys,xe,ye;
228: /* Get local vector (including ghost points) */
229: VecGetArray(user->localX,&x);
230: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
231: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
233: xe = xs+xm;
234: ye = ys+ym;
235: for (j=ys; j<ye; j++) { /* for (j=0; j<ny; j++) */
236: temp = (double)PetscMin(j+1,ny-j)*hy;
237: for (i=xs; i<xe; i++) { /* for (i=0; i<nx; i++) */
238: k = (j-Ys)*Xm + i-Xs;
239: #if !defined(PETSC_USE_COMPLEX)
240: x[k] = PetscMin((PetscMin(i+1,nx-i))*hx,temp);
241: #else
242: x[k] = PetscMin(PetscRealPart((double)(PetscMin(i+1,nx-i))*hx),PetscRealPart(temp));
243: #endif
244: }
245: }
246: VecRestoreArray(user->localX,&x);
248: /* Insert values into global vector */
249: DALocalToGlobal(user->da,user->localX,INSERT_VALUES,X);
250: return 0;
251: }
252: /* ---------- Evaluate function f(x) and/or gradient g(x) ----------- */
254: int EvalFunctionGradient(SNES snes,Vec X,double *f,Vec gvec,FctGradFlag fg,AppCtx *user)
255: {
256: Scalar hx = user->hx,hy = user->hy,area,three = 3.0,p5 = 0.5,cdiv3;
257: Scalar zero = 0.0,v,vb,vl,vr,vt,dvdx,dvdy,flin = 0.0,fquad = 0.0;
258: Scalar val,*x,szero = 0.0,floc;
259: Vec localX = user->localX;
260: int xs,ys,xm,ym,Xm,Ym,Xs,Ys,xe,ye,xsm,ysm,xep,yep;
261: int ierr,nx = user->mx,ny = user->my,ind,i,j,k,*ltog,nloc;
263: cdiv3 = user->param/three;
265: /* Get ghost points */
266: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
267: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
268: VecGetArray(localX,&x);
269: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
270: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
271: DAGetGlobalIndices(user->da,&nloc,<og);
272: xe = xs+xm;
273: ye = ys+ym;
274: if (xs == 0) xsm = xs-1;
275: else xsm = xs;
276: if (ys == 0) ysm = ys-1;
277: else ysm = ys;
278: if (xe == nx) xep = xe+1;
279: else xep = xe;
280: if (ye == ny) yep = ye+1;
281: else yep = ye;
283: if (fg & GradientEval) {
284: VecSet(&szero,gvec);
285: }
287: /* Compute function and gradient over the lower triangular elements */
288: for (j=ysm; j<ye; j++) { /* for (j=-1; j<ny; j++) */
289: for (i=xsm; i<xe; i++) { /* for (i=-1; i<nx; i++) */
290: k = (j-Ys)*Xm + i-Xs;
291: v = zero;
292: vr = zero;
293: vt = zero;
294: if (i >= 0 && j >= 0) v = x[k];
295: if (i < nx-1 && j > -1) vr = x[k+1];
296: if (i > -1 && j < ny-1) vt = x[k+Xm];
297: dvdx = (vr-v)/hx;
298: dvdy = (vt-v)/hy;
299: if (fg & FunctionEval) {
300: fquad += dvdx*dvdx + dvdy*dvdy;
301: flin -= cdiv3*(v+vr+vt);
302: }
303: if (fg & GradientEval) {
304: if (i != -1 && j != -1) {
305: ind = ltog[k]; val = - dvdx/hx - dvdy/hy - cdiv3;
306: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
307: }
308: if (i != nx-1 && j != -1) {
309: ind = ltog[k+1]; val = dvdx/hx - cdiv3;
310: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
311: }
312: if (i != -1 && j != ny-1) {
313: ind = ltog[k+Xm]; val = dvdy/hy - cdiv3;
314: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
315: }
316: }
317: }
318: }
320: /* Compute function and gradient over the upper triangular elements */
321: for (j=ys; j<yep; j++) { /* for (j=0; j<=ny; j++) */
322: for (i=xs; i<xep; i++) { /* for (i=0; i<=nx; i++) */
323: k = (j-Ys)*Xm + i-Xs;
324: vb = zero;
325: vl = zero;
326: v = zero;
327: if (i < nx && j > 0) vb = x[k-Xm];
328: if (i > 0 && j < ny) vl = x[k-1];
329: if (i < nx && j < ny) v = x[k];
330: dvdx = (v-vl)/hx;
331: dvdy = (v-vb)/hy;
332: if (fg & FunctionEval) {
333: fquad = fquad + dvdx*dvdx + dvdy*dvdy;
334: flin = flin - cdiv3*(vb+vl+v);
335: } if (fg & GradientEval) {
336: if (i != nx && j != 0) {
337: ind = ltog[k-Xm]; val = - dvdy/hy - cdiv3;
338: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
339: }
340: if (i != 0 && j != ny) {
341: ind = ltog[k-1]; val = - dvdx/hx - cdiv3;
342: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
343: }
344: if (i != nx && j != ny) {
345: ind = ltog[k]; val = dvdx/hx + dvdy/hy - cdiv3;
346: VecSetValues(gvec,1,&ind,&val,ADD_VALUES);
347: }
348: }
349: }
350: }
351: VecRestoreArray(localX,&x);
352: area = p5*hx*hy;
353: if (fg & FunctionEval) { /* Scale the function */
354: #if !defined(PETSC_USE_COMPLEX)
355: floc = area*(p5*fquad+flin);
356: #else
357: floc = PetscRealPart(area*(p5*fquad+flin));
358: #endif
359: MPI_Allreduce((void*)&floc,(void*)f,1,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD);
360: } if (fg & GradientEval) { /* Scale the gradient */
361: VecAssemblyBegin(gvec);
362: VecAssemblyEnd(gvec);
363: VecScale((Scalar*)&area,gvec);
364: }
365: return 0;
366: }
367: /* --------------------------------------------------------------------- */
368: int HessianProductMat(Mat mat,Vec svec,Vec y)
369: {
370: void *ptr;
371: MatShellGetContext(mat,&ptr);
372: HessianProduct(ptr,svec,y);
373: return 0;
374: }
375: /*
376: HessianProduct - Computes the matrix-vector product: y = f''(x)*s
377: */
378: int HessianProduct(void *ptr,Vec svec,Vec y)
379: {
380: AppCtx *user = (AppCtx*)ptr;
381: Scalar p5 = 0.5,one = 1.0,zero = 0.0,hx,hy;
382: Scalar val,area,*x,*s,szero = 0.0,v,vb,vl,vr,vt,hxhx,hyhy;
383: Vec localX,localS;
384: int xs,ys,xm,ym,Xm,Ym,Xs,Ys,xe,ye,xsm,ysm,xep,yep;
385: int nx,ny,i,j,k,ierr,ind,nloc,*ltog;
387: hx = user->hx;
388: hy = user->hy;
389: localX = user->localX;
390: localS = user->localS;
391: nx = user->mx;
392: ny = user->my;
394: hxhx = one/(hx*hx);
395: hyhy = one/(hy*hy);
397: /* Get ghost points */
398: DAGlobalToLocalBegin(user->da,user->xvec,INSERT_VALUES,localX);
399: DAGlobalToLocalEnd(user->da,user->xvec,INSERT_VALUES,localX);
400: DAGlobalToLocalBegin(user->da,svec,INSERT_VALUES,localS);
401: DAGlobalToLocalEnd(user->da,svec,INSERT_VALUES,localS);
402: VecGetArray(localS,&s);
403: VecGetArray(localX,&x);
404: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
405: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
406: DAGetGlobalIndices(user->da,&nloc,<og);
407: xe = xs+xm;
408: ye = ys+ym;
409: if (xs == 0) xsm = xs-1;
410: else xsm = xs;
411: if (ys == 0) ysm = ys-1;
412: else ysm = ys;
413: if (xe == nx) xep = xe+1;
414: else xep = xe;
415: if (ye == ny) yep = ye+1;
416: else yep = ye;
418: VecSet(&szero,y);
420: /* Compute f''(x)*s over the lower triangular elements */
421: for (j=ysm; j<ye; j++) { /* for (j=-1; j<ny; j++) */
422: for (i=xsm; i<xe; i++) { /* for (i=-1; i<nx; i++) */
423: k = (j-Ys)*Xm + i-Xs;
424: v = zero;
425: vr = zero;
426: vt = zero;
427: if (i != -1 && j != -1) v = s[k];
428: if (i != nx-1 && j != -1) {
429: vr = s[k+1];
430: ind = ltog[k+1]; val = hxhx*(vr-v);
431: VecSetValues(y,1,&ind,&val,ADD_VALUES);
432: }
433: if (i != -1 && j != ny-1) {
434: vt = s[k+Xm];
435: ind = ltog[k+Xm]; val = hyhy*(vt-v);
436: VecSetValues(y,1,&ind,&val,ADD_VALUES);
437: }
438: if (i != -1 && j != -1) {
439: ind = ltog[k]; val = hxhx*(v-vr) + hyhy*(v-vt);
440: VecSetValues(y,1,&ind,&val,ADD_VALUES);
441: }
442: }
443: }
445: /* Compute f''(x)*s over the upper triangular elements */
446: for (j=ys; j<yep; j++) { /* for (j=0; j<=ny; j++) */
447: for (i=xs; i<xep; i++) { /* for (i=0; i<=nx; i++) */
448: k = (j-Ys)*Xm + i-Xs;
449: v = zero;
450: vl = zero;
451: vb = zero;
452: if (i != nx && j != ny) v = s[k];
453: if (i != nx && j != 0) {
454: vb = s[k-Xm];
455: ind = ltog[k-Xm]; val = hyhy*(vb-v);
456: VecSetValues(y,1,&ind,&val,ADD_VALUES);
457: }
458: if (i != 0 && j != ny) {
459: vl = s[k-1];
460: ind = ltog[k-1]; val = hxhx*(vl-v);
461: VecSetValues(y,1,&ind,&val,ADD_VALUES);
462: }
463: if (i != nx && j != ny) {
464: ind = ltog[k]; val = hxhx*(v-vl) + hyhy*(v-vb);
465: VecSetValues(y,1,&ind,&val,ADD_VALUES);
466: }
467: }
468: }
469: VecAssemblyBegin(y);
470: VecRestoreArray(localX,&x);
471: VecRestoreArray(localS,&x);
472: VecAssemblyEnd(y);
474: /* Scale result by area */
475: area = p5*hx*hy;
476: VecScale(&area,y);
477: return 0;
478: }