xref: /petsc/src/tao/unconstrained/tutorials/minsurf2.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown /* Program usage: mpiexec -n <proc> minsurf2 [-help] [all TAO options] */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   Include "petsctao.h" so we can use TAO solvers.
5c4762a1bSJed Brown   petscdmda.h for distributed array
6c4762a1bSJed Brown */
7c4762a1bSJed Brown #include <petsctao.h>
8c4762a1bSJed Brown #include <petscdmda.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown static  char help[] =
11c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\
12c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
13c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and \n\
14c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
15c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
16c4762a1bSJed Brown The command line options are:\n\
17c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
18c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
19c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
20c4762a1bSJed Brown                for an average of the boundary conditions\n\n";
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /*T
23c4762a1bSJed Brown    Concepts: TAO^Solving an unconstrained minimization problem
24c4762a1bSJed Brown    Routines: TaoCreate(); TaoSetType();
25a82e8c82SStefano Zampini    Routines: TaoSetSolution();
26a82e8c82SStefano Zampini    Routines: TaoSetObjectiveAndGradient();
27a82e8c82SStefano Zampini    Routines: TaoSetHessian(); TaoSetFromOptions();
28c4762a1bSJed Brown    Routines: TaoSetMonitor();
29c4762a1bSJed Brown    Routines: TaoSolve(); TaoView();
30c4762a1bSJed Brown    Routines: TaoDestroy();
31c4762a1bSJed Brown    Processors: n
32c4762a1bSJed Brown T*/
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /*
35c4762a1bSJed Brown    User-defined application context - contains data needed by the
36c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient()
37c4762a1bSJed Brown    and FormHessian().
38c4762a1bSJed Brown */
39c4762a1bSJed Brown typedef struct {
40c4762a1bSJed Brown   PetscInt    mx, my;                 /* discretization in x, y directions */
41c4762a1bSJed Brown   PetscReal   *bottom, *top, *left, *right;             /* boundary values */
42c4762a1bSJed Brown   DM          dm;                      /* distributed array data structure */
43c4762a1bSJed Brown   Mat         H;                       /* Hessian */
44c4762a1bSJed Brown } AppCtx;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /* -------- User-defined Routines --------- */
47c4762a1bSJed Brown 
48c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
49c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
50c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx*,Vec,Mat);
51c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void*);
52c4762a1bSJed Brown PetscErrorCode FormGradient(Tao,Vec,Vec,void*);
53c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
54c4762a1bSJed Brown PetscErrorCode My_Monitor(Tao, void *);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown int main(int argc, char **argv)
57c4762a1bSJed Brown {
58c4762a1bSJed Brown   PetscInt           Nx, Ny;              /* number of processors in x- and y- directions */
59c4762a1bSJed Brown   Vec                x;                   /* solution, gradient vectors */
60c4762a1bSJed Brown   PetscBool          flg, viewmat;        /* flags */
61c4762a1bSJed Brown   PetscBool          fddefault, fdcoloring;   /* flags */
62c4762a1bSJed Brown   Tao                tao;                 /* TAO solver context */
63c4762a1bSJed Brown   AppCtx             user;                /* user-defined work context */
64c4762a1bSJed Brown   ISColoring         iscoloring;
65c4762a1bSJed Brown   MatFDColoring      matfdcoloring;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Initialize TAO */
68*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv,(char *)0,help));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* Specify dimension of the problem */
71c4762a1bSJed Brown   user.mx = 10; user.my = 10;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
76c4762a1bSJed Brown 
775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(MPI_COMM_WORLD,"\n---- Minimum Surface Area Problem -----\n"));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(MPI_COMM_WORLD,"mx: %D     my: %D   \n\n",user.mx,user.my));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Let PETSc determine the vector distribution */
81c4762a1bSJed Brown   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Create distributed array (DM) to manage parallel grid and vectors  */
845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,user.mx, user.my,Nx,Ny,1,1,NULL,NULL,&user.dm));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(user.dm));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(user.dm));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Create TAO solver and set desired solution method.*/
895f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOCG));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /*
93c4762a1bSJed Brown      Extract global vector from DA for the vector of variables --  PETSC routine
94c4762a1bSJed Brown      Compute the initial solution                              --  application specific, see below
95c4762a1bSJed Brown      Set this vector for use by TAO                            --  TAO routine
96c4762a1bSJed Brown   */
975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(user.dm,&x));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MSA_BoundaryConditions(&user));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(MSA_InitialPoint(&user,x));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /*
103c4762a1bSJed Brown      Initialize the Application context for use in function evaluations  --  application specific, see below.
104c4762a1bSJed Brown      Set routines for function and gradient evaluation
105c4762a1bSJed Brown   */
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /*
109c4762a1bSJed Brown      Given the command line arguments, calculate the hessian with either the user-
110c4762a1bSJed Brown      provided function FormHessian, or the default finite-difference driven Hessian
111c4762a1bSJed Brown      functions
112c4762a1bSJed Brown   */
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-tao_fddefault",&fddefault));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-tao_fdcoloring",&fdcoloring));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /*
117c4762a1bSJed Brown      Create a matrix data structure to store the Hessian and set
118c4762a1bSJed Brown      the Hessian evalution routine.
119c4762a1bSJed Brown      Set the matrix structure to be used for Hessian evalutions
120c4762a1bSJed Brown   */
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(user.dm,&user.H));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   if (fdcoloring) {
1255f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateColoring(user.dm,IS_COLORING_GLOBAL,&iscoloring));
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFDColoringCreate(user.H,iscoloring,&matfdcoloring));
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(ISColoringDestroy(&iscoloring));
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormGradient,(void*)&user));
1295f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFDColoringSetFromOptions(matfdcoloring));
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessianColor,(void *)matfdcoloring));
131c4762a1bSJed Brown   } else if (fddefault) {
1325f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessian,(void *)NULL));
133c4762a1bSJed Brown   } else {
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user));
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown      If my_monitor option is in command line, then use the user-provided
139c4762a1bSJed Brown      monitoring function
140c4762a1bSJed Brown   */
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-my_monitor",&viewmat));
142c4762a1bSJed Brown   if (viewmat) {
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetMonitor(tao,My_Monitor,NULL,NULL));
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Check for any tao command line options */
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
151c4762a1bSJed Brown 
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /* Free TAO data structures */
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* Free PETSc data structures */
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.H));
160c4762a1bSJed Brown   if (fdcoloring) {
1615f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFDColoringDestroy(&matfdcoloring));
162c4762a1bSJed Brown   }
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.bottom));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.top));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.left));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.right));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.dm));
168*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
169*b122ec5aSJacob Faibussowitsch   return 0;
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao, Vec X, Vec G,void *userCtx)
173c4762a1bSJed Brown {
174c4762a1bSJed Brown   PetscReal      fcn;
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   PetscFunctionBegin;
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(FormFunctionGradient(tao,X,&fcn,G,userCtx));
178c4762a1bSJed Brown   PetscFunctionReturn(0);
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown /* -------------------------------------------------------------------- */
182c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates the function and corresponding gradient.
183c4762a1bSJed Brown 
184c4762a1bSJed Brown     Input Parameters:
185c4762a1bSJed Brown .   tao     - the Tao context
186c4762a1bSJed Brown .   XX      - input vector
187a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
188c4762a1bSJed Brown 
189c4762a1bSJed Brown     Output Parameters:
190c4762a1bSJed Brown .   fcn     - the newly evaluated function
191c4762a1bSJed Brown .   GG       - vector containing the newly evaluated gradient
192c4762a1bSJed Brown */
19370a7d78aSStefano Zampini PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *userCtx)
19470a7d78aSStefano Zampini {
195c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) userCtx;
196c4762a1bSJed Brown   PetscInt       i,j;
197c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
198c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
199c4762a1bSJed Brown   PetscReal      ft=0.0;
200c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
201c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
202c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
203c4762a1bSJed Brown   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
204c4762a1bSJed Brown   PetscReal      **g, **x;
205c4762a1bSJed Brown   Vec            localX;
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   PetscFunctionBegin;
208c4762a1bSJed Brown   /* Get local mesh boundaries */
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(user->dm,&localX));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* Scatter ghost points to local vector */
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* Get pointers to vector data */
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->dm,localX,(void**)&x));
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->dm,G,(void**)&g));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Compute function and gradient over the locally owned part of the mesh */
222c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
223c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
224c4762a1bSJed Brown 
225c4762a1bSJed Brown       xc = x[j][i];
226c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
227c4762a1bSJed Brown 
228c4762a1bSJed Brown       if (i==0) { /* left side */
229c4762a1bSJed Brown         xl= user->left[j-ys+1];
230c4762a1bSJed Brown         xlt = user->left[j-ys+2];
231c4762a1bSJed Brown       } else {
232c4762a1bSJed Brown         xl = x[j][i-1];
233c4762a1bSJed Brown       }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown       if (j==0) { /* bottom side */
236c4762a1bSJed Brown         xb=user->bottom[i-xs+1];
237c4762a1bSJed Brown         xrb =user->bottom[i-xs+2];
238c4762a1bSJed Brown       } else {
239c4762a1bSJed Brown         xb = x[j-1][i];
240c4762a1bSJed Brown       }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown       if (i+1 == gxs+gxm) { /* right side */
243c4762a1bSJed Brown         xr=user->right[j-ys+1];
244c4762a1bSJed Brown         xrb = user->right[j-ys];
245c4762a1bSJed Brown       } else {
246c4762a1bSJed Brown         xr = x[j][i+1];
247c4762a1bSJed Brown       }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown       if (j+1==gys+gym) { /* top side */
250c4762a1bSJed Brown         xt=user->top[i-xs+1];
251c4762a1bSJed Brown         xlt = user->top[i-xs];
252c4762a1bSJed Brown       }else {
253c4762a1bSJed Brown         xt = x[j+1][i];
254c4762a1bSJed Brown       }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym) {
257c4762a1bSJed Brown         xlt = x[j+1][i-1];
258c4762a1bSJed Brown       }
259c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm) {
260c4762a1bSJed Brown         xrb = x[j-1][i+1];
261c4762a1bSJed Brown       }
262c4762a1bSJed Brown 
263c4762a1bSJed Brown       d1 = (xc-xl);
264c4762a1bSJed Brown       d2 = (xc-xr);
265c4762a1bSJed Brown       d3 = (xc-xt);
266c4762a1bSJed Brown       d4 = (xc-xb);
267c4762a1bSJed Brown       d5 = (xr-xrb);
268c4762a1bSJed Brown       d6 = (xrb-xb);
269c4762a1bSJed Brown       d7 = (xlt-xl);
270c4762a1bSJed Brown       d8 = (xt-xlt);
271c4762a1bSJed Brown 
272c4762a1bSJed Brown       df1dxc = d1*hydhx;
273c4762a1bSJed Brown       df2dxc = (d1*hydhx + d4*hxdhy);
274c4762a1bSJed Brown       df3dxc = d3*hxdhy;
275c4762a1bSJed Brown       df4dxc = (d2*hydhx + d3*hxdhy);
276c4762a1bSJed Brown       df5dxc = d2*hydhx;
277c4762a1bSJed Brown       df6dxc = d4*hxdhy;
278c4762a1bSJed Brown 
279c4762a1bSJed Brown       d1 *= rhx;
280c4762a1bSJed Brown       d2 *= rhx;
281c4762a1bSJed Brown       d3 *= rhy;
282c4762a1bSJed Brown       d4 *= rhy;
283c4762a1bSJed Brown       d5 *= rhy;
284c4762a1bSJed Brown       d6 *= rhx;
285c4762a1bSJed Brown       d7 *= rhy;
286c4762a1bSJed Brown       d8 *= rhx;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
289c4762a1bSJed Brown       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
290c4762a1bSJed Brown       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
291c4762a1bSJed Brown       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
292c4762a1bSJed Brown       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
293c4762a1bSJed Brown       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
294c4762a1bSJed Brown 
295c4762a1bSJed Brown       ft = ft + (f2 + f4);
296c4762a1bSJed Brown 
297c4762a1bSJed Brown       df1dxc /= f1;
298c4762a1bSJed Brown       df2dxc /= f2;
299c4762a1bSJed Brown       df3dxc /= f3;
300c4762a1bSJed Brown       df4dxc /= f4;
301c4762a1bSJed Brown       df5dxc /= f5;
302c4762a1bSJed Brown       df6dxc /= f6;
303c4762a1bSJed Brown 
304c4762a1bSJed Brown       g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5;
305c4762a1bSJed Brown 
306c4762a1bSJed Brown     }
307c4762a1bSJed Brown   }
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
310c4762a1bSJed Brown   if (xs==0) { /* left side */
311c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
312c4762a1bSJed Brown       d3=(user->left[j-ys+1] - user->left[j-ys+2])*rhy;
313c4762a1bSJed Brown       d2=(user->left[j-ys+1] - x[j][0]) *rhx;
314c4762a1bSJed Brown       ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2);
315c4762a1bSJed Brown     }
316c4762a1bSJed Brown   }
317c4762a1bSJed Brown   if (ys==0) { /* bottom side */
318c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
319c4762a1bSJed Brown       d2=(user->bottom[i+1-xs]-user->bottom[i-xs+2])*rhx;
320c4762a1bSJed Brown       d3=(user->bottom[i-xs+1]-x[0][i])*rhy;
321c4762a1bSJed Brown       ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2);
322c4762a1bSJed Brown     }
323c4762a1bSJed Brown   }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   if (xs+xm==mx) { /* right side */
326c4762a1bSJed Brown     for (j=ys; j< ys+ym; j++) {
327c4762a1bSJed Brown       d1=(x[j][mx-1] - user->right[j-ys+1])*rhx;
328c4762a1bSJed Brown       d4=(user->right[j-ys]-user->right[j-ys+1])*rhy;
329c4762a1bSJed Brown       ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4);
330c4762a1bSJed Brown     }
331c4762a1bSJed Brown   }
332c4762a1bSJed Brown   if (ys+ym==my) { /* top side */
333c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
334c4762a1bSJed Brown       d1=(x[my-1][i] - user->top[i-xs+1])*rhy;
335c4762a1bSJed Brown       d4=(user->top[i-xs+1] - user->top[i-xs])*rhx;
336c4762a1bSJed Brown       ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4);
337c4762a1bSJed Brown     }
338c4762a1bSJed Brown   }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   if (ys==0 && xs==0) {
341c4762a1bSJed Brown     d1=(user->left[0]-user->left[1])/hy;
342c4762a1bSJed Brown     d2=(user->bottom[0]-user->bottom[1])*rhx;
343c4762a1bSJed Brown     ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2);
344c4762a1bSJed Brown   }
345c4762a1bSJed Brown   if (ys+ym == my && xs+xm == mx) {
346c4762a1bSJed Brown     d1=(user->right[ym+1] - user->right[ym])*rhy;
347c4762a1bSJed Brown     d2=(user->top[xm+1] - user->top[xm])*rhx;
348c4762a1bSJed Brown     ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2);
349c4762a1bSJed Brown   }
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   ft=ft*area;
3525f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD));
353c4762a1bSJed Brown 
354c4762a1bSJed Brown   /* Restore vectors */
3555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->dm,localX,(void**)&x));
3565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->dm,G,(void**)&g));
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   /* Scatter values to global vector */
3595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(user->dm,&localX));
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(67.0*xm*ym));
361c4762a1bSJed Brown   PetscFunctionReturn(0);
362c4762a1bSJed Brown }
363c4762a1bSJed Brown 
364c4762a1bSJed Brown /* ------------------------------------------------------------------- */
365c4762a1bSJed Brown /*
366c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
367c4762a1bSJed Brown 
368c4762a1bSJed Brown    Input Parameters:
369c4762a1bSJed Brown .  tao  - the Tao context
370c4762a1bSJed Brown .  x    - input vector
371a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
372c4762a1bSJed Brown 
373c4762a1bSJed Brown    Output Parameters:
374c4762a1bSJed Brown .  H    - Hessian matrix
375c4762a1bSJed Brown .  Hpre - optionally different preconditioning matrix
376c4762a1bSJed Brown .  flg  - flag indicating matrix structure
377c4762a1bSJed Brown 
378c4762a1bSJed Brown */
379c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
380c4762a1bSJed Brown {
381c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   PetscFunctionBegin;
384c4762a1bSJed Brown   /* Evaluate the Hessian entries*/
3855f80ce2aSJacob Faibussowitsch   CHKERRQ(QuadraticH(user,X,H));
386c4762a1bSJed Brown   PetscFunctionReturn(0);
387c4762a1bSJed Brown }
388c4762a1bSJed Brown 
389c4762a1bSJed Brown /* ------------------------------------------------------------------- */
390c4762a1bSJed Brown /*
391c4762a1bSJed Brown    QuadraticH - Evaluates Hessian matrix.
392c4762a1bSJed Brown 
393c4762a1bSJed Brown    Input Parameters:
394c4762a1bSJed Brown .  user - user-defined context, as set by TaoSetHessian()
395c4762a1bSJed Brown .  X    - input vector
396c4762a1bSJed Brown 
397c4762a1bSJed Brown    Output Parameter:
398c4762a1bSJed Brown .  H    - Hessian matrix
399c4762a1bSJed Brown */
400c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian)
401c4762a1bSJed Brown {
402c4762a1bSJed Brown   PetscInt i,j,k;
403c4762a1bSJed Brown   PetscInt mx=user->mx, my=user->my;
404c4762a1bSJed Brown   PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym;
405c4762a1bSJed Brown   PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
406c4762a1bSJed Brown   PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
407c4762a1bSJed Brown   PetscReal hl,hr,ht,hb,hc,htl,hbr;
408c4762a1bSJed Brown   PetscReal **x, v[7];
409c4762a1bSJed Brown   MatStencil col[7],row;
410c4762a1bSJed Brown   Vec    localX;
411c4762a1bSJed Brown   PetscBool assembled;
412c4762a1bSJed Brown 
413c4762a1bSJed Brown   PetscFunctionBegin;
414c4762a1bSJed Brown   /* Get local mesh boundaries */
4155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(user->dm,&localX));
416c4762a1bSJed Brown 
4175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
4185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   /* Scatter ghost points to local vector */
4215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
4225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /* Get pointers to vector data */
4255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->dm,localX,(void**)&x));
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   /* Initialize matrix entries to zero */
4285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssembled(Hessian,&assembled));
4295f80ce2aSJacob Faibussowitsch   if (assembled) CHKERRQ(MatZeroEntries(Hessian));
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   /* Set various matrix options */
4325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
437c4762a1bSJed Brown 
438c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
439c4762a1bSJed Brown 
440c4762a1bSJed Brown       xc = x[j][i];
441c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
442c4762a1bSJed Brown 
443c4762a1bSJed Brown       /* Left side */
444c4762a1bSJed Brown       if (i==0) {
445c4762a1bSJed Brown         xl  = user->left[j-ys+1];
446c4762a1bSJed Brown         xlt = user->left[j-ys+2];
447c4762a1bSJed Brown       } else {
448c4762a1bSJed Brown         xl  = x[j][i-1];
449c4762a1bSJed Brown       }
450c4762a1bSJed Brown 
451c4762a1bSJed Brown       if (j==0) {
452c4762a1bSJed Brown         xb  = user->bottom[i-xs+1];
453c4762a1bSJed Brown         xrb = user->bottom[i-xs+2];
454c4762a1bSJed Brown       } else {
455c4762a1bSJed Brown         xb  = x[j-1][i];
456c4762a1bSJed Brown       }
457c4762a1bSJed Brown 
458c4762a1bSJed Brown       if (i+1 == mx) {
459c4762a1bSJed Brown         xr  = user->right[j-ys+1];
460c4762a1bSJed Brown         xrb = user->right[j-ys];
461c4762a1bSJed Brown       } else {
462c4762a1bSJed Brown         xr  = x[j][i+1];
463c4762a1bSJed Brown       }
464c4762a1bSJed Brown 
465c4762a1bSJed Brown       if (j+1==my) {
466c4762a1bSJed Brown         xt  = user->top[i-xs+1];
467c4762a1bSJed Brown         xlt = user->top[i-xs];
468c4762a1bSJed Brown       }else {
469c4762a1bSJed Brown         xt  = x[j+1][i];
470c4762a1bSJed Brown       }
471c4762a1bSJed Brown 
472c4762a1bSJed Brown       if (i>0 && j+1<my) {
473c4762a1bSJed Brown         xlt = x[j+1][i-1];
474c4762a1bSJed Brown       }
475c4762a1bSJed Brown       if (j>0 && i+1<mx) {
476c4762a1bSJed Brown         xrb = x[j-1][i+1];
477c4762a1bSJed Brown       }
478c4762a1bSJed Brown 
479c4762a1bSJed Brown       d1 = (xc-xl)/hx;
480c4762a1bSJed Brown       d2 = (xc-xr)/hx;
481c4762a1bSJed Brown       d3 = (xc-xt)/hy;
482c4762a1bSJed Brown       d4 = (xc-xb)/hy;
483c4762a1bSJed Brown       d5 = (xrb-xr)/hy;
484c4762a1bSJed Brown       d6 = (xrb-xb)/hx;
485c4762a1bSJed Brown       d7 = (xlt-xl)/hy;
486c4762a1bSJed Brown       d8 = (xlt-xt)/hx;
487c4762a1bSJed Brown 
488c4762a1bSJed Brown       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
489c4762a1bSJed Brown       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
490c4762a1bSJed Brown       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
491c4762a1bSJed Brown       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
492c4762a1bSJed Brown       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
493c4762a1bSJed Brown       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
494c4762a1bSJed Brown 
495c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
496c4762a1bSJed Brown         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
497c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
498c4762a1bSJed Brown         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
499c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
500c4762a1bSJed Brown         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
501c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
502c4762a1bSJed Brown         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
503c4762a1bSJed Brown 
504c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
505c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
506c4762a1bSJed Brown 
507c4762a1bSJed Brown       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
508c4762a1bSJed Brown         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
509c4762a1bSJed Brown         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
510c4762a1bSJed Brown         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
511c4762a1bSJed Brown 
512c4762a1bSJed Brown       hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0;  hc/=2.0;
513c4762a1bSJed Brown 
514c4762a1bSJed Brown       row.j = j; row.i = i;
515c4762a1bSJed Brown       k=0;
516c4762a1bSJed Brown       if (j>0) {
517c4762a1bSJed Brown         v[k]=hb;
518c4762a1bSJed Brown         col[k].j = j - 1; col[k].i = i;
519c4762a1bSJed Brown         k++;
520c4762a1bSJed Brown       }
521c4762a1bSJed Brown 
522c4762a1bSJed Brown       if (j>0 && i < mx -1) {
523c4762a1bSJed Brown         v[k]=hbr;
524c4762a1bSJed Brown         col[k].j = j - 1; col[k].i = i+1;
525c4762a1bSJed Brown         k++;
526c4762a1bSJed Brown       }
527c4762a1bSJed Brown 
528c4762a1bSJed Brown       if (i>0) {
529c4762a1bSJed Brown         v[k]= hl;
530c4762a1bSJed Brown         col[k].j = j; col[k].i = i-1;
531c4762a1bSJed Brown         k++;
532c4762a1bSJed Brown       }
533c4762a1bSJed Brown 
534c4762a1bSJed Brown       v[k]= hc;
535c4762a1bSJed Brown       col[k].j = j; col[k].i = i;
536c4762a1bSJed Brown       k++;
537c4762a1bSJed Brown 
538c4762a1bSJed Brown       if (i < mx-1) {
539c4762a1bSJed Brown         v[k]= hr;
540c4762a1bSJed Brown         col[k].j = j; col[k].i = i+1;
541c4762a1bSJed Brown         k++;
542c4762a1bSJed Brown       }
543c4762a1bSJed Brown 
544c4762a1bSJed Brown       if (i>0 && j < my-1) {
545c4762a1bSJed Brown         v[k]= htl;
546c4762a1bSJed Brown         col[k].j = j+1; col[k].i = i-1;
547c4762a1bSJed Brown         k++;
548c4762a1bSJed Brown       }
549c4762a1bSJed Brown 
550c4762a1bSJed Brown       if (j < my-1) {
551c4762a1bSJed Brown         v[k]= ht;
552c4762a1bSJed Brown         col[k].j = j+1; col[k].i = i;
553c4762a1bSJed Brown         k++;
554c4762a1bSJed Brown       }
555c4762a1bSJed Brown 
556c4762a1bSJed Brown       /*
557c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
558c4762a1bSJed Brown          earlier, in the main routine.
559c4762a1bSJed Brown       */
5605f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(Hessian,1,&row,k,col,v,INSERT_VALUES));
561c4762a1bSJed Brown     }
562c4762a1bSJed Brown   }
563c4762a1bSJed Brown 
5645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->dm,localX,(void**)&x));
5655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(user->dm,&localX));
566c4762a1bSJed Brown 
5675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
5685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
569c4762a1bSJed Brown 
5705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(199.0*xm*ym));
571c4762a1bSJed Brown   PetscFunctionReturn(0);
572c4762a1bSJed Brown }
573c4762a1bSJed Brown 
574c4762a1bSJed Brown /* ------------------------------------------------------------------- */
575c4762a1bSJed Brown /*
576c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
577c4762a1bSJed Brown    the region.
578c4762a1bSJed Brown 
579c4762a1bSJed Brown    Input Parameter:
580c4762a1bSJed Brown .  user - user-defined application context
581c4762a1bSJed Brown 
582c4762a1bSJed Brown    Output Parameter:
583c4762a1bSJed Brown .  user - user-defined application context
584c4762a1bSJed Brown */
585c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
586c4762a1bSJed Brown {
587c4762a1bSJed Brown   PetscInt       i,j,k,limit=0,maxits=5;
588c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,gxs,gys,gxm,gym;
589c4762a1bSJed Brown   PetscInt       mx=user->mx,my=user->my;
590c4762a1bSJed Brown   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
591c4762a1bSJed Brown   PetscReal      one=1.0, two=2.0, three=3.0, tol=1e-10;
592c4762a1bSJed Brown   PetscReal      fnorm,det,hx,hy,xt=0,yt=0;
593c4762a1bSJed Brown   PetscReal      u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
594c4762a1bSJed Brown   PetscReal      b=-0.5, t=0.5, l=-0.5, r=0.5;
595c4762a1bSJed Brown   PetscReal      *boundary;
596c4762a1bSJed Brown   PetscBool      flg;
597c4762a1bSJed Brown 
598c4762a1bSJed Brown   PetscFunctionBegin;
599c4762a1bSJed Brown   /* Get local mesh boundaries */
6005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
6015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
602c4762a1bSJed Brown 
603c4762a1bSJed Brown   bsize=xm+2;
604c4762a1bSJed Brown   lsize=ym+2;
605c4762a1bSJed Brown   rsize=ym+2;
606c4762a1bSJed Brown   tsize=xm+2;
607c4762a1bSJed Brown 
6085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bsize,&user->bottom));
6095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(tsize,&user->top));
6105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(lsize,&user->left));
6115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(rsize,&user->right));
612c4762a1bSJed Brown 
613c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   for (j=0; j<4; j++) {
616c4762a1bSJed Brown     if (j==0) {
617c4762a1bSJed Brown       yt=b;
618c4762a1bSJed Brown       xt=l+hx*xs;
619c4762a1bSJed Brown       limit=bsize;
620c4762a1bSJed Brown       boundary=user->bottom;
621c4762a1bSJed Brown     } else if (j==1) {
622c4762a1bSJed Brown       yt=t;
623c4762a1bSJed Brown       xt=l+hx*xs;
624c4762a1bSJed Brown       limit=tsize;
625c4762a1bSJed Brown       boundary=user->top;
626c4762a1bSJed Brown     } else if (j==2) {
627c4762a1bSJed Brown       yt=b+hy*ys;
628c4762a1bSJed Brown       xt=l;
629c4762a1bSJed Brown       limit=lsize;
630c4762a1bSJed Brown       boundary=user->left;
631c4762a1bSJed Brown     } else { /* if (j==3) */
632c4762a1bSJed Brown       yt=b+hy*ys;
633c4762a1bSJed Brown       xt=r;
634c4762a1bSJed Brown       limit=rsize;
635c4762a1bSJed Brown       boundary=user->right;
636c4762a1bSJed Brown     }
637c4762a1bSJed Brown 
638c4762a1bSJed Brown     for (i=0; i<limit; i++) {
639c4762a1bSJed Brown       u1=xt;
640c4762a1bSJed Brown       u2=-yt;
641c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
642c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
643c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
644c4762a1bSJed Brown         fnorm=PetscSqrtReal(nf1*nf1+nf2*nf2);
645c4762a1bSJed Brown         if (fnorm <= tol) break;
646c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
647c4762a1bSJed Brown         njac12=two*u1*u2;
648c4762a1bSJed Brown         njac21=-two*u1*u2;
649c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
650c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
651c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
652c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
653c4762a1bSJed Brown       }
654c4762a1bSJed Brown 
655c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
656c4762a1bSJed Brown       if (j==0 || j==1) {
657c4762a1bSJed Brown         xt=xt+hx;
658c4762a1bSJed Brown       } else { /*  if (j==2 || j==3) */
659c4762a1bSJed Brown         yt=yt+hy;
660c4762a1bSJed Brown       }
661c4762a1bSJed Brown 
662c4762a1bSJed Brown     }
663c4762a1bSJed Brown 
664c4762a1bSJed Brown   }
665c4762a1bSJed Brown 
666c4762a1bSJed Brown   /* Scale the boundary if desired */
667c4762a1bSJed Brown   if (1==1) {
668c4762a1bSJed Brown     PetscReal scl = 1.0;
669c4762a1bSJed Brown 
6705f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg));
671c4762a1bSJed Brown     if (flg) {
672c4762a1bSJed Brown       for (i=0;i<bsize;i++) user->bottom[i]*=scl;
673c4762a1bSJed Brown     }
674c4762a1bSJed Brown 
6755f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg));
676c4762a1bSJed Brown     if (flg) {
677c4762a1bSJed Brown       for (i=0;i<tsize;i++) user->top[i]*=scl;
678c4762a1bSJed Brown     }
679c4762a1bSJed Brown 
6805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg));
681c4762a1bSJed Brown     if (flg) {
682c4762a1bSJed Brown       for (i=0;i<rsize;i++) user->right[i]*=scl;
683c4762a1bSJed Brown     }
684c4762a1bSJed Brown 
6855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg));
686c4762a1bSJed Brown     if (flg) {
687c4762a1bSJed Brown       for (i=0;i<lsize;i++) user->left[i]*=scl;
688c4762a1bSJed Brown     }
689c4762a1bSJed Brown   }
690c4762a1bSJed Brown   PetscFunctionReturn(0);
691c4762a1bSJed Brown }
692c4762a1bSJed Brown 
693c4762a1bSJed Brown /* ------------------------------------------------------------------- */
694c4762a1bSJed Brown /*
695c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
696c4762a1bSJed Brown 
697c4762a1bSJed Brown    Input Parameters:
698c4762a1bSJed Brown .  user - user-defined application context
699c4762a1bSJed Brown .  X - vector for initial guess
700c4762a1bSJed Brown 
701c4762a1bSJed Brown    Output Parameters:
702c4762a1bSJed Brown .  X - newly computed initial guess
703c4762a1bSJed Brown */
704c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
705c4762a1bSJed Brown {
706c4762a1bSJed Brown   PetscInt       start2=-1,i,j;
707c4762a1bSJed Brown   PetscReal      start1=0;
708c4762a1bSJed Brown   PetscBool      flg1,flg2;
709c4762a1bSJed Brown 
710c4762a1bSJed Brown   PetscFunctionBegin;
7115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-start",&start1,&flg1));
7125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-random",&start2,&flg2));
713c4762a1bSJed Brown 
714c4762a1bSJed Brown   if (flg1) { /* The zero vector is reasonable */
715c4762a1bSJed Brown 
7165f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(X, start1));
717c4762a1bSJed Brown 
718c4762a1bSJed Brown   } else if (flg2 && start2>0) { /* Try a random start between -0.5 and 0.5 */
719c4762a1bSJed Brown     PetscRandom rctx;  PetscReal np5=-0.5;
720c4762a1bSJed Brown 
7215f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
7225f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(X, rctx));
7235f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomDestroy(&rctx));
7245f80ce2aSJacob Faibussowitsch     CHKERRQ(VecShift(X, np5));
725c4762a1bSJed Brown 
726c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
727c4762a1bSJed Brown     PetscInt  xs,xm,ys,ym;
728c4762a1bSJed Brown     PetscInt  mx=user->mx,my=user->my;
729c4762a1bSJed Brown     PetscReal **x;
730c4762a1bSJed Brown 
731c4762a1bSJed Brown     /* Get local mesh boundaries */
7325f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
733c4762a1bSJed Brown 
734c4762a1bSJed Brown     /* Get pointers to vector data */
7355f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(user->dm,X,(void**)&x));
736c4762a1bSJed Brown 
737c4762a1bSJed Brown     /* Perform local computations */
738c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
739c4762a1bSJed Brown       for (i=xs; i< xs+xm; i++) {
740c4762a1bSJed Brown         x[j][i] = (((j+1)*user->bottom[i-xs+1]+(my-j+1)*user->top[i-xs+1])/(my+2)+((i+1)*user->left[j-ys+1]+(mx-i+1)*user->right[j-ys+1])/(mx+2))/2.0;
741c4762a1bSJed Brown       }
742c4762a1bSJed Brown     }
7435f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(user->dm,X,(void**)&x));
7445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogFlops(9.0*xm*ym));
745c4762a1bSJed Brown   }
746c4762a1bSJed Brown   PetscFunctionReturn(0);
747c4762a1bSJed Brown }
748c4762a1bSJed Brown 
749c4762a1bSJed Brown /*-----------------------------------------------------------------------*/
750c4762a1bSJed Brown PetscErrorCode My_Monitor(Tao tao, void *ctx)
751c4762a1bSJed Brown {
752c4762a1bSJed Brown   Vec            X;
753c4762a1bSJed Brown 
754c4762a1bSJed Brown   PetscFunctionBegin;
7555f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetSolution(tao,&X));
7565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD));
757c4762a1bSJed Brown   PetscFunctionReturn(0);
758c4762a1bSJed Brown }
759c4762a1bSJed Brown 
760c4762a1bSJed Brown /*TEST
761c4762a1bSJed Brown 
762c4762a1bSJed Brown    build:
763c4762a1bSJed Brown       requires: !complex
764c4762a1bSJed Brown 
765c4762a1bSJed Brown    test:
766c4762a1bSJed Brown       args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3
767c4762a1bSJed Brown       requires: !single
768c4762a1bSJed Brown 
769c4762a1bSJed Brown    test:
770c4762a1bSJed Brown       suffix: 2
771c4762a1bSJed Brown       nsize: 2
772c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4
773c4762a1bSJed Brown       filter: grep -v "nls ksp"
774c4762a1bSJed Brown       requires: !single
775c4762a1bSJed Brown 
776c4762a1bSJed Brown    test:
777c4762a1bSJed Brown       suffix: 3
778c4762a1bSJed Brown       nsize: 3
779c4762a1bSJed Brown       args: -tao_smonitor -tao_type cg -tao_cg_type fr -mx 10 -my 10 -tao_gatol 1.e-3
780c4762a1bSJed Brown       requires: !single
781c4762a1bSJed Brown 
782c4762a1bSJed Brown    test:
783c4762a1bSJed Brown       suffix: 5
784c4762a1bSJed Brown       nsize: 2
785c4762a1bSJed Brown       args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3
786c4762a1bSJed Brown       requires: !single
787c4762a1bSJed Brown 
788c4762a1bSJed Brown TEST*/
789