Actual source code: minsurf2.c
1: /*$Id$*/
3: /* Program usage: mpirun -np <proc> minsurf2 [-help] [all TAO options] */
5: /*
6: Include "tao.h" so we can use TAO solvers.
7: petscda.h for distributed array
8: */
9: #include "petscda.h"
10: #include "tao.h"
12: static char help[] =
13: "This example demonstrates use of the TAO package to \n\
14: solve an unconstrained minimization problem. This example is based on a \n\
15: problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\
16: boundary values along the edges of the domain, the objective is to find the\n\
17: surface with the minimal area that satisfies the boundary conditions.\n\
18: The command line options are:\n\
19: -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
20: -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
21: -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
22: for an average of the boundary conditions\n\n";
24: /*T
25: Concepts: TAO - Solving an unconstrained minimization problem
26: Routines: TaoInitialize(); TaoFinalize();
27: Routines: TaoCreate(); TaoDestroy();
28: Routines: TaoApplicationCreate(); TaoAppDestroy();
29: Routines: TaoAppSetInitialSolutionVec();
30: Routines: TaoAppSetObjectiveAndGradientRoutine();
31: Routines: TaoAppSetHessianMat(); TaoAppSetHessianRoutine();
32: Routines: TaoSetOptions();
33: Routines: TaoAppGetKSP(); TaoSolveApplication();
34: Routines: TaoAppSetMonitor(); TaoView();
35: Routines: TaoAppGetSolutionVec();
36: Processors: 1
37: T*/
39: /*
40: User-defined application context - contains data needed by the
41: application-provided call-back routines, FormFunctionGradient()
42: and FormHessian().
43: */
44: typedef struct {
45: int mx, my; /* discretization in x, y directions */
46: double *bottom, *top, *left, *right; /* boundary values */
47: DA da; /* distributed array data structure */
48: Mat H; /* Hessian */
49: ISColoring iscoloring;
50: } AppCtx;
53: /* -------- User-defined Routines --------- */
55: static int MSA_BoundaryConditions(AppCtx*);
56: static int MSA_InitialPoint(AppCtx*,Vec);
57: int QuadraticH(AppCtx*,Vec,Mat);
58: int FormFunctionGradient(TAO_APPLICATION,Vec,double *,Vec,void*);
59: int FormGradient(TAO_APPLICATION,Vec,Vec,void*);
60: int FormHessian(TAO_APPLICATION,Vec,Mat*,Mat*,MatStructure *,void*);
61: int My_Monitor(TAO_APPLICATION, void *);
65: int main( int argc, char **argv )
66: {
67: int info; /* used to check for functions returning nonzeros */
68: int Nx, Ny; /* number of processors in x- and y- directions */
69: int iter; /* iteration information */
70: double ff,gnorm;
71: Vec x; /* solution, gradient vectors */
72: PetscTruth flg, viewmat; /* flags */
73: PetscTruth fddefault, fdcoloring; /* flags */
74: KSP ksp; /* Krylov subspace method */
75: TaoMethod method = "tao_nls"; /* minimization method */
76: TaoTerminateReason reason;
77: TAO_SOLVER tao; /* TAO_SOLVER solver context */
78: TAO_APPLICATION minsurfapp; /* The PETSc application */
79: AppCtx user; /* user-defined work context */
81: /* Initialize TAO */
82: PetscInitialize( &argc, &argv,(char *)0,help );
83: TaoInitialize( &argc, &argv,(char *)0,help );
85: /* Specify dimension of the problem */
86: user.mx = 10; user.my = 10;
88: /* Check for any command line arguments that override defaults */
89: info = PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,&flg); CHKERRQ(info);
90: info = PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,&flg); CHKERRQ(info);
92: PetscPrintf(MPI_COMM_WORLD,"\n---- Minimum Surface Area Problem -----\n");
93: PetscPrintf(MPI_COMM_WORLD,"mx: %d my: %d \n\n",user.mx,user.my);
96: /* Let PETSc determine the vector distribution */
97: Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
99: /* Create distributed array (DA) to manage parallel grid and vectors */
100: info = DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,user.mx,
101: user.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.da); CHKERRQ(info);
102:
104: /* Create TAO solver and set desired solution method. Create an TAO application structure */
105: info = TaoCreate(PETSC_COMM_WORLD,method,&tao); CHKERRQ(info);
106: info = TaoApplicationCreate(PETSC_COMM_WORLD,&minsurfapp); CHKERRQ(info);
108: /*
109: Extract global vector from DA for the vector of variables -- PETSC routine
110: Compute the initial solution -- application specific, see below
111: Set this vector for use by TAO -- TAO routine
112: */
113: info = DACreateGlobalVector(user.da,&x); CHKERRQ(info);
114: info = MSA_BoundaryConditions(&user); CHKERRQ(info);
115: info = MSA_InitialPoint(&user,x); CHKERRQ(info);
116: info = TaoAppSetInitialSolutionVec(minsurfapp,x); CHKERRQ(info);
118: /*
119: Initialize the Application context for use in function evaluations -- application specific, see below.
120: Set routines for function and gradient evaluation
121: */
122: info = TaoAppSetObjectiveAndGradientRoutine(minsurfapp,FormFunctionGradient,(void *)&user); CHKERRQ(info);
124: /*
125: Given the command line arguments, calculate the hessian with either the user-
126: provided function FormHessian, or the default finite-difference driven Hessian
127: functions
128: */
129: info = PetscOptionsHasName(PETSC_NULL,"-tao_fddefault",&fddefault);CHKERRQ(info);
130: info = PetscOptionsHasName(PETSC_NULL,"-tao_fdcoloring",&fdcoloring);CHKERRQ(info);
132: if (fdcoloring) {
133: info = DAGetColoring(user.da,IS_COLORING_GLOBAL,&(user.iscoloring));
134: CHKERRQ(info);
135: info = TaoAppSetColoring(minsurfapp, user.iscoloring); CHKERRQ(info);
136: info = TaoAppSetHessianRoutine(minsurfapp,TaoAppDefaultComputeHessianColor,(void *)&user); CHKERRQ(info);
137: } else if (fddefault){
138: info = TaoAppSetHessianRoutine(minsurfapp,TaoAppDefaultComputeHessian,(void *)&user); CHKERRQ(info);
139: } else {
140: info = TaoAppSetHessianRoutine(minsurfapp,FormHessian,(void *)&user); CHKERRQ(info);
141: }
143: /*
144: Create a matrix data structure to store the Hessian and set
145: the Hessian evalution routine.
146: Set the matrix structure to be used for Hessian evalutions
147: */
148: info = DAGetMatrix(user.da,MATAIJ,&user.H);CHKERRQ(info);
149: info = MatSetOption(user.H,MAT_SYMMETRIC); CHKERRQ(info);
151: info = TaoAppSetHessianMat(minsurfapp,user.H,user.H); CHKERRQ(info);
153: /*
154: If my_monitor option is in command line, then use the user-provided
155: monitoring function
156: */
157: info = PetscOptionsHasName(PETSC_NULL,"-my_monitor",&viewmat); CHKERRQ(info);
158: if (viewmat){
159: info = TaoAppSetMonitor(minsurfapp,My_Monitor,TAO_NULL); CHKERRQ(info);
160: }
162: /* Check for any tao command line options */
163: info = TaoSetOptions(minsurfapp,tao); CHKERRQ(info);
165: /* Limit the number of iterations in the KSP linear solver */
166: info = TaoAppGetKSP(minsurfapp,&ksp); CHKERRQ(info);
167: if (ksp) { /* Modify the PETSc KSP structure */
168: info = KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,user.mx*user.my);
169: CHKERRQ(info);
170: }
172: /* SOLVE THE APPLICATION */
173: info = TaoSolveApplication(minsurfapp,tao); CHKERRQ(info);
175: /* Get information on termination */
176: info = TaoGetSolutionStatus(tao,&iter,&ff,&gnorm,0,0,&reason); CHKERRQ(info);
177: if (reason <= 0 ){
178: PetscPrintf(MPI_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");
179: PetscPrintf(MPI_COMM_WORLD," Iterations: %d, Function Value: %4.2e, Residual: %4.2e \n",iter,ff,gnorm);
180: }
182: /*
183: To View TAO solver information use
184: info = TaoView(tao); CHKERRQ(info);
185: */
187: /* Free TAO data structures */
188: info = TaoDestroy(tao); CHKERRQ(info);
189: info = TaoAppDestroy(minsurfapp); CHKERRQ(info);
191: /* Free PETSc data structures */
192: info = VecDestroy(x); CHKERRQ(info);
193: info = MatDestroy(user.H); CHKERRQ(info);
194: if (fdcoloring) {
195: info = ISColoringDestroy(user.iscoloring);CHKERRQ(info);
196: }
197: info = PetscFree(user.bottom); CHKERRQ(info);
198: info = PetscFree(user.top); CHKERRQ(info);
199: info = PetscFree(user.left); CHKERRQ(info);
200: info = PetscFree(user.right); CHKERRQ(info);
201: info = DADestroy(user.da); CHKERRQ(info);
203: /* Finalize TAO */
204: TaoFinalize();
205: PetscFinalize();
206:
207: return 0;
208: }
212: int FormGradient(TAO_APPLICATION taoapp, Vec X, Vec G,void *userCtx){
213: int info;
214: double fcn;
215: TaoFunctionBegin;
216: info = FormFunctionGradient(taoapp,X,&fcn,G,userCtx);CHKERRQ(info);
217: TaoFunctionReturn(0);
218: }
220: /* -------------------------------------------------------------------- */
223: /* FormFunctionGradient - Evaluates the function and corresponding gradient.
225: Input Parameters:
226: . taoapp - the TAO_APPLICATION context
227: . XX - input vector
228: . userCtx - optional user-defined context, as set by TaoSetFunctionGradient()
229:
230: Output Parameters:
231: . fcn - the newly evaluated function
232: . GG - vector containing the newly evaluated gradient
233: */
234: int FormFunctionGradient(TAO_APPLICATION taoapp, Vec X, double *fcn,Vec G,void *userCtx){
236: AppCtx * user = (AppCtx *) userCtx;
237: int info,i,j;
238: int mx=user->mx, my=user->my;
239: int xs,xm,gxs,gxm,ys,ym,gys,gym;
240: double ft=0;
241: double hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
242: double rhx=mx+1, rhy=my+1;
243: double f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
244: double df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
245: PetscScalar **g, **x;
246: Vec localX;
248: /* Get local mesh boundaries */
249: info = DAGetLocalVector(user->da,&localX);CHKERRQ(info);
251: info = DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL); CHKERRQ(info);
252: info = DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL); CHKERRQ(info);
254: /* Scatter ghost points to local vector */
255: info = DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX); CHKERRQ(info);
256: info = DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX); CHKERRQ(info);
258: /* Get pointers to vector data */
259: info = DAVecGetArray(user->da,localX,(void**)&x);
260: info = DAVecGetArray(user->da,G,(void**)&g);
262: /* Compute function and gradient over the locally owned part of the mesh */
263: for (j=ys; j<ys+ym; j++){
264: for (i=xs; i< xs+xm; i++){
265:
266: xc = x[j][i];
267: xlt=xrb=xl=xr=xb=xt=xc;
268:
269: if (i==0){ /* left side */
270: xl= user->left[j-ys+1];
271: xlt = user->left[j-ys+2];
272: } else {
273: xl = x[j][i-1];
274: }
276: if (j==0){ /* bottom side */
277: xb=user->bottom[i-xs+1];
278: xrb =user->bottom[i-xs+2];
279: } else {
280: xb = x[j-1][i];
281: }
282:
283: if (i+1 == gxs+gxm){ /* right side */
284: xr=user->right[j-ys+1];
285: xrb = user->right[j-ys];
286: } else {
287: xr = x[j][i+1];
288: }
290: if (j+1==gys+gym){ /* top side */
291: xt=user->top[i-xs+1];
292: xlt = user->top[i-xs];
293: }else {
294: xt = x[j+1][i];
295: }
297: if (i>gxs && j+1<gys+gym){
298: xlt = x[j+1][i-1];
299: }
300: if (j>gys && i+1<gxs+gxm){
301: xrb = x[j-1][i+1];
302: }
304: d1 = (xc-xl);
305: d2 = (xc-xr);
306: d3 = (xc-xt);
307: d4 = (xc-xb);
308: d5 = (xr-xrb);
309: d6 = (xrb-xb);
310: d7 = (xlt-xl);
311: d8 = (xt-xlt);
312:
313: df1dxc = d1*hydhx;
314: df2dxc = ( d1*hydhx + d4*hxdhy );
315: df3dxc = d3*hxdhy;
316: df4dxc = ( d2*hydhx + d3*hxdhy );
317: df5dxc = d2*hydhx;
318: df6dxc = d4*hxdhy;
320: d1 *= rhx;
321: d2 *= rhx;
322: d3 *= rhy;
323: d4 *= rhy;
324: d5 *= rhy;
325: d6 *= rhx;
326: d7 *= rhy;
327: d8 *= rhx;
329: f1 = sqrt( 1.0 + d1*d1 + d7*d7);
330: f2 = sqrt( 1.0 + d1*d1 + d4*d4);
331: f3 = sqrt( 1.0 + d3*d3 + d8*d8);
332: f4 = sqrt( 1.0 + d3*d3 + d2*d2);
333: f5 = sqrt( 1.0 + d2*d2 + d5*d5);
334: f6 = sqrt( 1.0 + d4*d4 + d6*d6);
335:
336: f2 = sqrt( 1.0 + d1*d1 + d4*d4);
337: f4 = sqrt( 1.0 + d3*d3 + d2*d2);
339: ft = ft + (f2 + f4);
341: df1dxc /= f1;
342: df2dxc /= f2;
343: df3dxc /= f3;
344: df4dxc /= f4;
345: df5dxc /= f5;
346: df6dxc /= f6;
348: g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc ) * 0.5;
349:
350: }
351: }
353: /* Compute triangular areas along the border of the domain. */
354: if (xs==0){ /* left side */
355: for (j=ys; j<ys+ym; j++){
356: d3=(user->left[j-ys+1] - user->left[j-ys+2])*rhy;
357: d2=(user->left[j-ys+1] - x[j][0]) *rhx;
358: ft = ft+sqrt( 1.0 + d3*d3 + d2*d2);
359: }
360: }
361: if (ys==0){ /* bottom side */
362: for (i=xs; i<xs+xm; i++){
363: d2=(user->bottom[i+1-xs]-user->bottom[i-xs+2])*rhx;
364: d3=(user->bottom[i-xs+1]-x[0][i])*rhy;
365: ft = ft+sqrt( 1.0 + d3*d3 + d2*d2);
366: }
367: }
369: if (xs+xm==mx){ /* right side */
370: for (j=ys; j< ys+ym; j++){
371: d1=(x[j][mx-1] - user->right[j-ys+1])*rhx;
372: d4=(user->right[j-ys]-user->right[j-ys+1])*rhy;
373: ft = ft+sqrt( 1.0 + d1*d1 + d4*d4);
374: }
375: }
376: if (ys+ym==my){ /* top side */
377: for (i=xs; i<xs+xm; i++){
378: d1=(x[my-1][i] - user->top[i-xs+1])*rhy;
379: d4=(user->top[i-xs+1] - user->top[i-xs])*rhx;
380: ft = ft+sqrt( 1.0 + d1*d1 + d4*d4);
381: }
382: }
384: if (ys==0 && xs==0){
385: d1=(user->left[0]-user->left[1])/hy;
386: d2=(user->bottom[0]-user->bottom[1])*rhx;
387: ft +=sqrt( 1.0 + d1*d1 + d2*d2);
388: }
389: if (ys+ym == my && xs+xm == mx){
390: d1=(user->right[ym+1] - user->right[ym])*rhy;
391: d2=(user->top[xm+1] - user->top[xm])*rhx;
392: ft +=sqrt( 1.0 + d1*d1 + d2*d2);
393: }
395: ft=ft*area;
396: info = MPI_Allreduce(&ft,fcn,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD);CHKERRQ(info);
398: /* Restore vectors */
399: info = DAVecRestoreArray(user->da,localX,(void**)&x);
400: info = DAVecRestoreArray(user->da,G,(void**)&g);
402: /* Scatter values to global vector */
403: info = DARestoreLocalVector(user->da,&localX); CHKERRQ(info);
405: info = PetscLogFlops(67*xm*ym); CHKERRQ(info);
407: return 0;
408: }
410: /* ------------------------------------------------------------------- */
413: /*
414: FormHessian - Evaluates Hessian matrix.
416: Input Parameters:
417: . taoapp - the TAO_APPLICATION context
418: . x - input vector
419: . ptr - optional user-defined context, as set by TaoSetHessian()
421: Output Parameters:
422: . H - Hessian matrix
423: . Hpre - optionally different preconditioning matrix
424: . flg - flag indicating matrix structure
426: */
427: int FormHessian(TAO_APPLICATION taoapp,Vec X,Mat *H, Mat *Hpre, MatStructure *flg, void *ptr)
428: {
429: int info;
430: AppCtx *user = (AppCtx *) ptr;
432: /* Evaluate the Hessian entries*/
433: info = QuadraticH(user,X,*H); CHKERRQ(info);
436: /* Indicate that this matrix has the same sparsity pattern during
437: successive iterations; setting this flag can save significant work
438: in computing the preconditioner for some methods. */
439: *flg=SAME_NONZERO_PATTERN;
441: return 0;
442: }
444: /* ------------------------------------------------------------------- */
447: /*
448: QuadraticH - Evaluates Hessian matrix.
450: Input Parameters:
451: . user - user-defined context, as set by TaoSetHessian()
452: . X - input vector
454: Output Parameter:
455: . H - Hessian matrix
456: */
457: int QuadraticH(AppCtx *user, Vec X, Mat Hessian)
458: {
459: int i,j,k,info;
460: int mx=user->mx, my=user->my;
461: int xs,xm,gxs,gxm,ys,ym,gys,gym;
462: double hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
463: double f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
464: double hl,hr,ht,hb,hc,htl,hbr;
465: PetscScalar **x, v[7];
466: MatStencil col[7],row;
467: Vec localX;
468: PetscTruth assembled;
470: /* Get local mesh boundaries */
471: info = DAGetLocalVector(user->da,&localX);CHKERRQ(info);
473: info = DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL); CHKERRQ(info);
474: info = DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL); CHKERRQ(info);
476: /* Scatter ghost points to local vector */
477: info = DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX); CHKERRQ(info);
478: info = DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX); CHKERRQ(info);
480: /* Get pointers to vector data */
481: info = DAVecGetArray(user->da,localX,(void**)&x);
483: /* Initialize matrix entries to zero */
484: info = MatAssembled(Hessian,&assembled); CHKERRQ(info);
485: if (assembled){info = MatZeroEntries(Hessian); CHKERRQ(info);}
488: /* Set various matrix options */
489: info = MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES); CHKERRQ(info);
490: info = MatSetOption(Hessian,MAT_COLUMNS_SORTED); CHKERRQ(info);
491: info = MatSetOption(Hessian,MAT_ROWS_SORTED); CHKERRQ(info);
493: /* Compute Hessian over the locally owned part of the mesh */
495: for (j=ys; j<ys+ym; j++){
496:
497: for (i=xs; i< xs+xm; i++){
499: xc = x[j][i];
500: xlt=xrb=xl=xr=xb=xt=xc;
502: /* Left side */
503: if (i==0){
504: xl = user->left[j-ys+1];
505: xlt = user->left[j-ys+2];
506: } else {
507: xl = x[j][i-1];
508: }
509:
510: if (j==0){
511: xb = user->bottom[i-xs+1];
512: xrb = user->bottom[i-xs+2];
513: } else {
514: xb = x[j-1][i];
515: }
516:
517: if (i+1 == mx){
518: xr = user->right[j-ys+1];
519: xrb = user->right[j-ys];
520: } else {
521: xr = x[j][i+1];
522: }
524: if (j+1==my){
525: xt = user->top[i-xs+1];
526: xlt = user->top[i-xs];
527: }else {
528: xt = x[j+1][i];
529: }
531: if (i>0 && j+1<my){
532: xlt = x[j+1][i-1];
533: }
534: if (j>0 && i+1<mx){
535: xrb = x[j-1][i+1];
536: }
539: d1 = (xc-xl)/hx;
540: d2 = (xc-xr)/hx;
541: d3 = (xc-xt)/hy;
542: d4 = (xc-xb)/hy;
543: d5 = (xrb-xr)/hy;
544: d6 = (xrb-xb)/hx;
545: d7 = (xlt-xl)/hy;
546: d8 = (xlt-xt)/hx;
547:
548: f1 = sqrt( 1.0 + d1*d1 + d7*d7);
549: f2 = sqrt( 1.0 + d1*d1 + d4*d4);
550: f3 = sqrt( 1.0 + d3*d3 + d8*d8);
551: f4 = sqrt( 1.0 + d3*d3 + d2*d2);
552: f5 = sqrt( 1.0 + d2*d2 + d5*d5);
553: f6 = sqrt( 1.0 + d4*d4 + d6*d6);
556: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
557: (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
558: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
559: (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
560: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
561: (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
562: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
563: (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
565: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
566: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
568: hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
569: hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
570: (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
571: (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
573: hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0;
575: row.j = j; row.i = i;
576: k=0;
577: if (j>0){
578: v[k]=hb;
579: col[k].j = j - 1; col[k].i = i;
580: k++;
581: }
582:
583: if (j>0 && i < mx -1){
584: v[k]=hbr;
585: col[k].j = j - 1; col[k].i = i+1;
586: k++;
587: }
588:
589: if (i>0){
590: v[k]= hl;
591: col[k].j = j; col[k].i = i-1;
592: k++;
593: }
594:
595: v[k]= hc;
596: col[k].j = j; col[k].i = i;
597: k++;
598:
599: if (i < mx-1 ){
600: v[k]= hr;
601: col[k].j = j; col[k].i = i+1;
602: k++;
603: }
604:
605: if (i>0 && j < my-1 ){
606: v[k]= htl;
607: col[k].j = j+1; col[k].i = i-1;
608: k++;
609: }
610:
611: if (j < my-1 ){
612: v[k]= ht;
613: col[k].j = j+1; col[k].i = i;
614: k++;
615: }
616:
617: /*
618: Set matrix values using local numbering, which was defined
619: earlier, in the main routine.
620: */
621: info = MatSetValuesStencil(Hessian,1,&row,k,col,v,INSERT_VALUES);
622: CHKERRQ(info);
623:
624: }
625: }
626:
627: /* Restore vectors */
628: info = DAVecRestoreArray(user->da,localX,(void**)&x);
630: info = DARestoreLocalVector(user->da,&localX); CHKERRQ(info);
632: /* Assemble the matrix */
633: info = MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
634: info = MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
636: info = PetscLogFlops(199*xm*ym); CHKERRQ(info);
637: return 0;
638: }
640: /* ------------------------------------------------------------------- */
643: /*
644: MSA_BoundaryConditions - Calculates the boundary conditions for
645: the region.
647: Input Parameter:
648: . user - user-defined application context
650: Output Parameter:
651: . user - user-defined application context
652: */
653: static int MSA_BoundaryConditions(AppCtx * user)
654: {
655: int i,j,k,limit=0,info,maxits=5;
656: int xs,ys,xm,ym,gxs,gys,gxm,gym;
657: int mx=user->mx,my=user->my;
658: int bsize=0, lsize=0, tsize=0, rsize=0;
659: double one=1.0, two=2.0, three=3.0, tol=1e-10;
660: double fnorm,det,hx,hy,xt=0,yt=0;
661: double u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
662: double b=-0.5, t=0.5, l=-0.5, r=0.5;
663: double *boundary;
664: PetscTruth flg;
666: /* Get local mesh boundaries */
667: info = DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL); CHKERRQ(info);
668: info = DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL); CHKERRQ(info);
670: bsize=xm+2;
671: lsize=ym+2;
672: rsize=ym+2;
673: tsize=xm+2;
675: info = PetscMalloc(bsize*sizeof(double),&user->bottom); CHKERRQ(info);
676: info = PetscMalloc(tsize*sizeof(double),&user->top); CHKERRQ(info);
677: info = PetscMalloc(lsize*sizeof(double),&user->left); CHKERRQ(info);
678: info = PetscMalloc(rsize*sizeof(double),&user->right); CHKERRQ(info);
680: hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
682: for (j=0; j<4; j++){
683: if (j==0){
684: yt=b;
685: xt=l+hx*xs;
686: limit=bsize;
687: boundary=user->bottom;
688: } else if (j==1){
689: yt=t;
690: xt=l+hx*xs;
691: limit=tsize;
692: boundary=user->top;
693: } else if (j==2){
694: yt=b+hy*ys;
695: xt=l;
696: limit=lsize;
697: boundary=user->left;
698: } else { //if (j==3)
699: yt=b+hy*ys;
700: xt=r;
701: limit=rsize;
702: boundary=user->right;
703: }
705: for (i=0; i<limit; i++){
706: u1=xt;
707: u2=-yt;
708: for (k=0; k<maxits; k++){
709: nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
710: nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
711: fnorm=sqrt(nf1*nf1+nf2*nf2);
712: if (fnorm <= tol) break;
713: njac11=one+u2*u2-u1*u1;
714: njac12=two*u1*u2;
715: njac21=-two*u1*u2;
716: njac22=-one - u1*u1 + u2*u2;
717: det = njac11*njac22-njac21*njac12;
718: u1 = u1-(njac22*nf1-njac12*nf2)/det;
719: u2 = u2-(njac11*nf2-njac21*nf1)/det;
720: }
722: boundary[i]=u1*u1-u2*u2;
723: if (j==0 || j==1) {
724: xt=xt+hx;
725: } else { // if (j==2 || j==3)
726: yt=yt+hy;
727: }
728:
729: }
731: }
733: /* Scale the boundary if desired */
734: if (1==1){
735: PetscReal scl = 1.0;
737: info = PetscOptionsGetReal(PETSC_NULL,"-bottom",&scl,&flg);
738: CHKERRQ(info);
739: if (flg){
740: for (i=0;i<bsize;i++) user->bottom[i]*=scl;
741: }
743: info = PetscOptionsGetReal(PETSC_NULL,"-top",&scl,&flg);
744: CHKERRQ(info);
745: if (flg){
746: for (i=0;i<tsize;i++) user->top[i]*=scl;
747: }
749: info = PetscOptionsGetReal(PETSC_NULL,"-right",&scl,&flg);
750: CHKERRQ(info);
751: if (flg){
752: for (i=0;i<rsize;i++) user->right[i]*=scl;
753: }
755: info = PetscOptionsGetReal(PETSC_NULL,"-left",&scl,&flg);
756: CHKERRQ(info);
757: if (flg){
758: for (i=0;i<lsize;i++) user->left[i]*=scl;
759: }
760: }
761:
762: return 0;
763: }
765: /* ------------------------------------------------------------------- */
768: /*
769: MSA_InitialPoint - Calculates the initial guess in one of three ways.
771: Input Parameters:
772: . user - user-defined application context
773: . X - vector for initial guess
775: Output Parameters:
776: . X - newly computed initial guess
777: */
778: static int MSA_InitialPoint(AppCtx * user, Vec X)
779: {
780: int start2=-1,i,j,info;
781: PetscReal start1=0;
782: PetscTruth flg1,flg2;
784: info = PetscOptionsGetReal(PETSC_NULL,"-start",&start1,&flg1); CHKERRQ(info);
785: info = PetscOptionsGetInt(PETSC_NULL,"-random",&start2,&flg2); CHKERRQ(info);
787: if (flg1){ /* The zero vector is reasonable */
788:
789: info = VecSet(X, start1); CHKERRQ(info);
791: } else if (flg2 && start2>0){ /* Try a random start between -0.5 and 0.5 */
793: PetscRandom rctx; PetscScalar np5=-0.5;
795: info = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
796: CHKERRQ(info);
797: for (i=0; i<start2; i++){
798: info = VecSetRandom(X, rctx); CHKERRQ(info);
799: }
800: info = PetscRandomDestroy(rctx); CHKERRQ(info);
801: info = VecShift(X, np5); CHKERRQ(info);
803: } else { /* Take an average of the boundary conditions */
805: int xs,xm,ys,ym;
806: int mx=user->mx,my=user->my;
807: PetscScalar **x;
808:
809: /* Get local mesh boundaries */
810: info = DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL); CHKERRQ(info);
811:
812: /* Get pointers to vector data */
813: info = DAVecGetArray(user->da,X,(void**)&x);
815: /* Perform local computations */
816: for (j=ys; j<ys+ym; j++){
817: for (i=xs; i< xs+xm; i++){
818: x[j][i] = ( ((j+1)*user->bottom[i-xs+1]+(my-j+1)*user->top[i-xs+1])/(my+2)+
819: ((i+1)*user->left[j-ys+1]+(mx-i+1)*user->right[j-ys+1])/(mx+2))/2.0;
820: }
821: }
822:
823: /* Restore vectors */
824: info = DAVecRestoreArray(user->da,X,(void**)&x); CHKERRQ(info);
826: info = PetscLogFlops(9*xm*ym); CHKERRQ(info);
827:
828: }
829: return 0;
830: }
832: /*-----------------------------------------------------------------------*/
835: int My_Monitor(TAO_APPLICATION minsurfapp, void *ctx){
836: int info;
837: Vec X;
839: info = TaoAppGetSolutionVec(minsurfapp,&X); CHKERRQ(info);
840: info = VecView(X,PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(info);
841: return 0;
842: }