xref: /petsc/src/tao/bound/tutorials/plate2.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown #include <petscdmda.h>
2c4762a1bSJed Brown #include <petsctao.h>
3c4762a1bSJed Brown 
4c4762a1bSJed Brown static  char help[] =
5c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\
6c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
7c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain, \n\
8c4762a1bSJed Brown boundary values along the edges of the domain, and a plate represented by \n\
9c4762a1bSJed Brown lower boundary conditions, the objective is to find the\n\
10c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
11c4762a1bSJed Brown The command line options are:\n\
12c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
13c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
14c4762a1bSJed Brown   -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\
15c4762a1bSJed Brown   -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\
16c4762a1bSJed Brown   -bheight <ht>, where <ht> = height of the plate\n\
17c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
18c4762a1bSJed Brown                for an average of the boundary conditions\n\n";
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /*T
21c4762a1bSJed Brown    Concepts: TAO^Solving a bound constrained minimization problem
22c4762a1bSJed Brown    Routines: TaoCreate();
23a82e8c82SStefano Zampini    Routines: TaoSetType(); TaoSetObjectiveAndGradient();
24a82e8c82SStefano Zampini    Routines: TaoSetHessian();
25a82e8c82SStefano Zampini    Routines: TaoSetSolution();
26c4762a1bSJed Brown    Routines: TaoSetVariableBounds();
27c4762a1bSJed Brown    Routines: TaoSetFromOptions();
28c4762a1bSJed Brown    Routines: TaoSolve(); TaoView();
29c4762a1bSJed Brown    Routines: TaoDestroy();
30c4762a1bSJed Brown    Processors: n
31c4762a1bSJed Brown T*/
32c4762a1bSJed Brown 
33c4762a1bSJed Brown /*
34c4762a1bSJed Brown    User-defined application context - contains data needed by the
35c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
36c4762a1bSJed Brown    FormHessian().
37c4762a1bSJed Brown */
38c4762a1bSJed Brown typedef struct {
39c4762a1bSJed Brown   /* problem parameters */
40c4762a1bSJed Brown   PetscReal      bheight;                  /* Height of plate under the surface */
41c4762a1bSJed Brown   PetscInt       mx, my;                   /* discretization in x, y directions */
42c4762a1bSJed Brown   PetscInt       bmx,bmy;                  /* Size of plate under the surface */
43c4762a1bSJed Brown   Vec            Bottom, Top, Left, Right; /* boundary values */
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Working space */
46c4762a1bSJed Brown   Vec         localX, localV;           /* ghosted local vector */
47c4762a1bSJed Brown   DM          dm;                       /* distributed array data structure */
48c4762a1bSJed Brown   Mat         H;
49c4762a1bSJed Brown } AppCtx;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown /* -------- User-defined Routines --------- */
52c4762a1bSJed Brown 
53c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
54c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
55c4762a1bSJed Brown static PetscErrorCode MSA_Plate(Vec,Vec,void*);
56c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
57c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /* For testing matrix free submatrices */
60c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao,Vec,Mat, Mat,void*);
61c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat,Vec,Vec);
62c4762a1bSJed Brown 
63c4762a1bSJed Brown int main(int argc, char **argv)
64c4762a1bSJed Brown {
65c4762a1bSJed Brown   PetscInt               Nx, Ny;               /* number of processors in x- and y- directions */
66c4762a1bSJed Brown   PetscInt               m, N;                 /* number of local and global elements in vectors */
67c4762a1bSJed Brown   Vec                    x,xl,xu;               /* solution vector  and bounds*/
68c4762a1bSJed Brown   PetscBool              flg;                /* A return variable when checking for user options */
69c4762a1bSJed Brown   Tao                    tao;                  /* Tao solver context */
70c4762a1bSJed Brown   ISLocalToGlobalMapping isltog;   /* local-to-global mapping object */
71c4762a1bSJed Brown   Mat                    H_shell;                  /* to test matrix-free submatrices */
72c4762a1bSJed Brown   AppCtx                 user;                 /* user-defined work context */
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Initialize PETSc, TAO */
75*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv,(char *)0,help));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Specify default dimension of the problem */
78c4762a1bSJed Brown   user.mx = 10; user.my = 10; user.bheight=0.1;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-bheight",&user.bheight,&flg));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   user.bmx = user.mx/2; user.bmy = user.my/2;
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bmx",&user.bmx,&flg));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bmy",&user.bmy,&flg));
88c4762a1bSJed Brown 
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n"));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mx:%D, my:%D, bmx:%D, bmy:%D, height:%g\n",user.mx,user.my,user.bmx,user.bmy,(double)user.bheight));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Calculate any derived values from parameters */
93c4762a1bSJed Brown   N    = user.mx*user.my;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Let Petsc determine the dimensions of the local vectors */
96c4762a1bSJed Brown   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /*
99c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
100c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
101c4762a1bSJed Brown      the distributed array, Create the vectors.
102c4762a1bSJed Brown   */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(MPI_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(user.dm));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(user.dm));
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Extract global and local vectors from DM; The local vectors are
108c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
109c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
110c4762a1bSJed Brown      the same types.
111c4762a1bSJed Brown   */
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(user.dm,&x)); /* Solution */
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(user.dm,&user.localX));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.localX,&user.localV));
115c4762a1bSJed Brown 
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&xl));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&xu));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* The TAO code begins here */
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /*
122c4762a1bSJed Brown      Create TAO solver and set desired solution method
123c4762a1bSJed Brown      The method must either be TAOTRON or TAOBLMVM
124c4762a1bSJed Brown      If TAOBLMVM is used, then hessian function is not called.
125c4762a1bSJed Brown   */
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBLMVM));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Set initial solution guess; */
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(MSA_BoundaryConditions(&user));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(MSA_InitialPoint(&user,x));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Set routines for function, gradient and hessian evaluation */
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user));
136c4762a1bSJed Brown 
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(x,&m));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H)));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
140c4762a1bSJed Brown 
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalToGlobalMapping(user.dm,&isltog));
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetLocalToGlobalMapping(user.H,isltog,isltog));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-matrixfree",&flg));
144c4762a1bSJed Brown   if (flg) {
1455f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,m,m,N,N,(void*)&user,&H_shell));
1465f80ce2aSJacob Faibussowitsch       CHKERRQ(MatShellSetOperation(H_shell,MATOP_MULT,(void(*)(void))MyMatMult));
1475f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(H_shell,MAT_SYMMETRIC,PETSC_TRUE));
1485f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetHessian(tao,H_shell,H_shell,MatrixFreeHessian,(void*)&user));
149c4762a1bSJed Brown   } else {
1505f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
151c4762a1bSJed Brown   }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* Set Variable bounds */
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(MSA_Plate(xl,xu,(void*)&user));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetVariableBounds(tao,xl,xu));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* Check for any tao command line options */
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
162c4762a1bSJed Brown 
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* Free TAO data structures */
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* Free PETSc data structures */
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&xl));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&xu));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.H));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.localX));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.localV));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Bottom));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Top));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Left));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Right));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.dm));
180c4762a1bSJed Brown   if (flg) {
1815f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&H_shell));
182c4762a1bSJed Brown   }
183*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
184*b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).
188c4762a1bSJed Brown 
189c4762a1bSJed Brown     Input Parameters:
190c4762a1bSJed Brown .   tao     - the Tao context
191c4762a1bSJed Brown .   X      - input vector
192a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
193c4762a1bSJed Brown 
194c4762a1bSJed Brown     Output Parameters:
195c4762a1bSJed Brown .   fcn     - the function value
196c4762a1bSJed Brown .   G      - vector containing the newly evaluated gradient
197c4762a1bSJed Brown 
198c4762a1bSJed Brown    Notes:
199c4762a1bSJed Brown    In this case, we discretize the domain and Create triangles. The
200c4762a1bSJed Brown    surface of each triangle is planar, whose surface area can be easily
201c4762a1bSJed Brown    computed.  The total surface area is found by sweeping through the grid
202c4762a1bSJed Brown    and computing the surface area of the two triangles that have their
203c4762a1bSJed Brown    right angle at the grid point.  The diagonal line segments on the
204c4762a1bSJed Brown    grid that define the triangles run from top left to lower right.
205c4762a1bSJed Brown    The numbering of points starts at the lower left and runs left to
206c4762a1bSJed Brown    right, then bottom to top.
207c4762a1bSJed Brown */
208c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G,void *userCtx)
209c4762a1bSJed Brown {
210c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) userCtx;
211c4762a1bSJed Brown   PetscInt       i,j,row;
212c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
213c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
214c4762a1bSJed Brown   PetscReal      ft=0;
215c4762a1bSJed Brown   PetscReal      zero=0.0;
216c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
217c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
218c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
219c4762a1bSJed Brown   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
220c4762a1bSJed Brown   PetscReal      *g, *x,*left,*right,*bottom,*top;
221c4762a1bSJed Brown   Vec            localX = user->localX, localG = user->localV;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* Get local mesh boundaries */
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* Scatter ghost points to local vector */
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* Initialize vector to zero */
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(localG, zero));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   /* Get pointers to vector data */
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localX,&x));
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localG,&g));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Top,&top));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Bottom,&bottom));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Left,&left));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Right,&right));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
243c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
244c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
245c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
246c4762a1bSJed Brown 
247c4762a1bSJed Brown       xc = x[row];
248c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
249c4762a1bSJed Brown 
250c4762a1bSJed Brown       if (i==0) { /* left side */
251c4762a1bSJed Brown         xl= left[j-ys+1];
252c4762a1bSJed Brown         xlt = left[j-ys+2];
253c4762a1bSJed Brown       } else {
254c4762a1bSJed Brown         xl = x[row-1];
255c4762a1bSJed Brown       }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown       if (j==0) { /* bottom side */
258c4762a1bSJed Brown         xb=bottom[i-xs+1];
259c4762a1bSJed Brown         xrb = bottom[i-xs+2];
260c4762a1bSJed Brown       } else {
261c4762a1bSJed Brown         xb = x[row-gxm];
262c4762a1bSJed Brown       }
263c4762a1bSJed Brown 
264c4762a1bSJed Brown       if (i+1 == gxs+gxm) { /* right side */
265c4762a1bSJed Brown         xr=right[j-ys+1];
266c4762a1bSJed Brown         xrb = right[j-ys];
267c4762a1bSJed Brown       } else {
268c4762a1bSJed Brown         xr = x[row+1];
269c4762a1bSJed Brown       }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown       if (j+1==gys+gym) { /* top side */
272c4762a1bSJed Brown         xt=top[i-xs+1];
273c4762a1bSJed Brown         xlt = top[i-xs];
274c4762a1bSJed Brown       }else {
275c4762a1bSJed Brown         xt = x[row+gxm];
276c4762a1bSJed Brown       }
277c4762a1bSJed Brown 
278c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym) {
279c4762a1bSJed Brown         xlt = x[row-1+gxm];
280c4762a1bSJed Brown       }
281c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm) {
282c4762a1bSJed Brown         xrb = x[row+1-gxm];
283c4762a1bSJed Brown       }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown       d1 = (xc-xl);
286c4762a1bSJed Brown       d2 = (xc-xr);
287c4762a1bSJed Brown       d3 = (xc-xt);
288c4762a1bSJed Brown       d4 = (xc-xb);
289c4762a1bSJed Brown       d5 = (xr-xrb);
290c4762a1bSJed Brown       d6 = (xrb-xb);
291c4762a1bSJed Brown       d7 = (xlt-xl);
292c4762a1bSJed Brown       d8 = (xt-xlt);
293c4762a1bSJed Brown 
294c4762a1bSJed Brown       df1dxc = d1*hydhx;
295c4762a1bSJed Brown       df2dxc = (d1*hydhx + d4*hxdhy);
296c4762a1bSJed Brown       df3dxc = d3*hxdhy;
297c4762a1bSJed Brown       df4dxc = (d2*hydhx + d3*hxdhy);
298c4762a1bSJed Brown       df5dxc = d2*hydhx;
299c4762a1bSJed Brown       df6dxc = d4*hxdhy;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown       d1 *= rhx;
302c4762a1bSJed Brown       d2 *= rhx;
303c4762a1bSJed Brown       d3 *= rhy;
304c4762a1bSJed Brown       d4 *= rhy;
305c4762a1bSJed Brown       d5 *= rhy;
306c4762a1bSJed Brown       d6 *= rhx;
307c4762a1bSJed Brown       d7 *= rhy;
308c4762a1bSJed Brown       d8 *= rhx;
309c4762a1bSJed Brown 
310c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
311c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
312c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
313c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
314c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
315c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
316c4762a1bSJed Brown 
317c4762a1bSJed Brown       ft = ft + (f2 + f4);
318c4762a1bSJed Brown 
319c4762a1bSJed Brown       df1dxc /= f1;
320c4762a1bSJed Brown       df2dxc /= f2;
321c4762a1bSJed Brown       df3dxc /= f3;
322c4762a1bSJed Brown       df4dxc /= f4;
323c4762a1bSJed Brown       df5dxc /= f5;
324c4762a1bSJed Brown       df6dxc /= f6;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5;
327c4762a1bSJed Brown 
328c4762a1bSJed Brown     }
329c4762a1bSJed Brown   }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
332c4762a1bSJed Brown   if (xs==0) { /* left side */
333c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
334c4762a1bSJed Brown       d3=(left[j-ys+1] - left[j-ys+2])*rhy;
335c4762a1bSJed Brown       d2=(left[j-ys+1] - x[(j-gys)*gxm])*rhx;
336c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
337c4762a1bSJed Brown     }
338c4762a1bSJed Brown   }
339c4762a1bSJed Brown   if (ys==0) { /* bottom side */
340c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
341c4762a1bSJed Brown       d2=(bottom[i+1-xs]-bottom[i-xs+2])*rhx;
342c4762a1bSJed Brown       d3=(bottom[i-xs+1]-x[i-gxs])*rhy;
343c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
344c4762a1bSJed Brown     }
345c4762a1bSJed Brown   }
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   if (xs+xm==mx) { /* right side */
348c4762a1bSJed Brown     for (j=ys; j< ys+ym; j++) {
349c4762a1bSJed Brown       d1=(x[(j+1-gys)*gxm-1]-right[j-ys+1])*rhx;
350c4762a1bSJed Brown       d4=(right[j-ys]-right[j-ys+1])*rhy;
351c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
352c4762a1bSJed Brown     }
353c4762a1bSJed Brown   }
354c4762a1bSJed Brown   if (ys+ym==my) { /* top side */
355c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
356c4762a1bSJed Brown       d1=(x[(gym-1)*gxm + i-gxs] - top[i-xs+1])*rhy;
357c4762a1bSJed Brown       d4=(top[i-xs+1] - top[i-xs])*rhx;
358c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
359c4762a1bSJed Brown     }
360c4762a1bSJed Brown   }
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   if (ys==0 && xs==0) {
363c4762a1bSJed Brown     d1=(left[0]-left[1])*rhy;
364c4762a1bSJed Brown     d2=(bottom[0]-bottom[1])*rhx;
365c4762a1bSJed Brown     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
366c4762a1bSJed Brown   }
367c4762a1bSJed Brown   if (ys+ym == my && xs+xm == mx) {
368c4762a1bSJed Brown     d1=(right[ym+1] - right[ym])*rhy;
369c4762a1bSJed Brown     d2=(top[xm+1] - top[xm])*rhx;
370c4762a1bSJed Brown     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
371c4762a1bSJed Brown   }
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   ft=ft*area;
3745f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD));
375c4762a1bSJed Brown 
376c4762a1bSJed Brown   /* Restore vectors */
3775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localX,&x));
3785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localG,&g));
3795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Left,&left));
3805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Top,&top));
3815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Bottom,&bottom));
3825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Right,&right));
383c4762a1bSJed Brown 
384c4762a1bSJed Brown   /* Scatter values to global vector */
3855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G));
3865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G));
387c4762a1bSJed Brown 
3885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(70.0*xm*ym));
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   return 0;
391c4762a1bSJed Brown }
392c4762a1bSJed Brown 
393c4762a1bSJed Brown /* ------------------------------------------------------------------- */
394c4762a1bSJed Brown /*
395c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
396c4762a1bSJed Brown 
397c4762a1bSJed Brown    Input Parameters:
398c4762a1bSJed Brown .  tao  - the Tao context
399c4762a1bSJed Brown .  x    - input vector
400a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
401c4762a1bSJed Brown 
402c4762a1bSJed Brown    Output Parameters:
403c4762a1bSJed Brown .  A    - Hessian matrix
404c4762a1bSJed Brown .  B    - optionally different preconditioning matrix
405c4762a1bSJed Brown 
406c4762a1bSJed Brown    Notes:
407c4762a1bSJed Brown    Due to mesh point reordering with DMs, we must always work
408c4762a1bSJed Brown    with the local mesh points, and then transform them to the new
409c4762a1bSJed Brown    global numbering with the local-to-global mapping.  We cannot work
410c4762a1bSJed Brown    directly with the global numbers for the original uniprocessor mesh!
411c4762a1bSJed Brown 
412c4762a1bSJed Brown    Two methods are available for imposing this transformation
413c4762a1bSJed Brown    when setting matrix entries:
414c4762a1bSJed Brown      (A) MatSetValuesLocal(), using the local ordering (including
415c4762a1bSJed Brown          ghost points!)
416c4762a1bSJed Brown          - Do the following two steps once, before calling TaoSolve()
417c4762a1bSJed Brown            - Use DMGetISLocalToGlobalMapping() to extract the
418c4762a1bSJed Brown              local-to-global map from the DM
419c4762a1bSJed Brown            - Associate this map with the matrix by calling
420c4762a1bSJed Brown              MatSetLocalToGlobalMapping()
421c4762a1bSJed Brown          - Then set matrix entries using the local ordering
422c4762a1bSJed Brown            by calling MatSetValuesLocal()
423c4762a1bSJed Brown      (B) MatSetValues(), using the global ordering
424c4762a1bSJed Brown          - Use DMGetGlobalIndices() to extract the local-to-global map
425c4762a1bSJed Brown          - Then apply this map explicitly yourself
426c4762a1bSJed Brown          - Set matrix entries using the global ordering by calling
427c4762a1bSJed Brown            MatSetValues()
428c4762a1bSJed Brown    Option (A) seems cleaner/easier in many cases, and is the procedure
429c4762a1bSJed Brown    used in this example.
430c4762a1bSJed Brown */
431c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat Hptr, Mat Hessian, void *ptr)
432c4762a1bSJed Brown {
433c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
434c4762a1bSJed Brown   PetscInt       i,j,k,row;
435c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
436c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym,col[7];
437c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
438c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
439c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
440c4762a1bSJed Brown   PetscReal      hl,hr,ht,hb,hc,htl,hbr;
441c4762a1bSJed Brown   PetscReal      *x,*left,*right,*bottom,*top;
442c4762a1bSJed Brown   PetscReal      v[7];
443c4762a1bSJed Brown   Vec            localX = user->localX;
444c4762a1bSJed Brown   PetscBool      assembled;
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   /* Set various matrix options */
4475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /* Initialize matrix entries to zero */
4505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssembled(Hessian,&assembled));
4515f80ce2aSJacob Faibussowitsch   if (assembled) CHKERRQ(MatZeroEntries(Hessian));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   /* Get local mesh boundaries */
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
4555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   /* Scatter ghost points to local vector */
4585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
4595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
460c4762a1bSJed Brown 
461c4762a1bSJed Brown   /* Get pointers to vector data */
4625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localX,&x));
4635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Top,&top));
4645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Bottom,&bottom));
4655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Left,&left));
4665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->Right,&right));
467c4762a1bSJed Brown 
468c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
469c4762a1bSJed Brown 
470c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++) {
471c4762a1bSJed Brown 
472c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
473c4762a1bSJed Brown 
474c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
475c4762a1bSJed Brown 
476c4762a1bSJed Brown       xc = x[row];
477c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
478c4762a1bSJed Brown 
479c4762a1bSJed Brown       /* Left side */
480c4762a1bSJed Brown       if (i==gxs) {
481c4762a1bSJed Brown         xl= left[j-ys+1];
482c4762a1bSJed Brown         xlt = left[j-ys+2];
483c4762a1bSJed Brown       } else {
484c4762a1bSJed Brown         xl = x[row-1];
485c4762a1bSJed Brown       }
486c4762a1bSJed Brown 
487c4762a1bSJed Brown       if (j==gys) {
488c4762a1bSJed Brown         xb=bottom[i-xs+1];
489c4762a1bSJed Brown         xrb = bottom[i-xs+2];
490c4762a1bSJed Brown       } else {
491c4762a1bSJed Brown         xb = x[row-gxm];
492c4762a1bSJed Brown       }
493c4762a1bSJed Brown 
494c4762a1bSJed Brown       if (i+1 == gxs+gxm) {
495c4762a1bSJed Brown         xr=right[j-ys+1];
496c4762a1bSJed Brown         xrb = right[j-ys];
497c4762a1bSJed Brown       } else {
498c4762a1bSJed Brown         xr = x[row+1];
499c4762a1bSJed Brown       }
500c4762a1bSJed Brown 
501c4762a1bSJed Brown       if (j+1==gys+gym) {
502c4762a1bSJed Brown         xt=top[i-xs+1];
503c4762a1bSJed Brown         xlt = top[i-xs];
504c4762a1bSJed Brown       }else {
505c4762a1bSJed Brown         xt = x[row+gxm];
506c4762a1bSJed Brown       }
507c4762a1bSJed Brown 
508c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym) {
509c4762a1bSJed Brown         xlt = x[row-1+gxm];
510c4762a1bSJed Brown       }
511c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm) {
512c4762a1bSJed Brown         xrb = x[row+1-gxm];
513c4762a1bSJed Brown       }
514c4762a1bSJed Brown 
515c4762a1bSJed Brown       d1 = (xc-xl)*rhx;
516c4762a1bSJed Brown       d2 = (xc-xr)*rhx;
517c4762a1bSJed Brown       d3 = (xc-xt)*rhy;
518c4762a1bSJed Brown       d4 = (xc-xb)*rhy;
519c4762a1bSJed Brown       d5 = (xrb-xr)*rhy;
520c4762a1bSJed Brown       d6 = (xrb-xb)*rhx;
521c4762a1bSJed Brown       d7 = (xlt-xl)*rhy;
522c4762a1bSJed Brown       d8 = (xlt-xt)*rhx;
523c4762a1bSJed Brown 
524c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
525c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
526c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
527c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
528c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
529c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
530c4762a1bSJed Brown 
531c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
532c4762a1bSJed Brown         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
533c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
534c4762a1bSJed Brown         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
535c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
536c4762a1bSJed Brown         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
537c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
538c4762a1bSJed Brown         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
539c4762a1bSJed Brown 
540c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
541c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
542c4762a1bSJed Brown 
543c4762a1bSJed Brown       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
544c4762a1bSJed Brown         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
545c4762a1bSJed Brown         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
546c4762a1bSJed Brown         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
547c4762a1bSJed Brown 
548c4762a1bSJed Brown       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;
549c4762a1bSJed Brown 
550c4762a1bSJed Brown       k=0;
551c4762a1bSJed Brown       if (j>0) {
552c4762a1bSJed Brown         v[k]=hb; col[k]=row - gxm; k++;
553c4762a1bSJed Brown       }
554c4762a1bSJed Brown 
555c4762a1bSJed Brown       if (j>0 && i < mx -1) {
556c4762a1bSJed Brown         v[k]=hbr; col[k]=row - gxm+1; k++;
557c4762a1bSJed Brown       }
558c4762a1bSJed Brown 
559c4762a1bSJed Brown       if (i>0) {
560c4762a1bSJed Brown         v[k]= hl; col[k]=row - 1; k++;
561c4762a1bSJed Brown       }
562c4762a1bSJed Brown 
563c4762a1bSJed Brown       v[k]= hc; col[k]=row; k++;
564c4762a1bSJed Brown 
565c4762a1bSJed Brown       if (i < mx-1) {
566c4762a1bSJed Brown         v[k]= hr; col[k]=row+1; k++;
567c4762a1bSJed Brown       }
568c4762a1bSJed Brown 
569c4762a1bSJed Brown       if (i>0 && j < my-1) {
570c4762a1bSJed Brown         v[k]= htl; col[k] = row+gxm-1; k++;
571c4762a1bSJed Brown       }
572c4762a1bSJed Brown 
573c4762a1bSJed Brown       if (j < my-1) {
574c4762a1bSJed Brown         v[k]= ht; col[k] = row+gxm; k++;
575c4762a1bSJed Brown       }
576c4762a1bSJed Brown 
577c4762a1bSJed Brown       /*
578c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
579c4762a1bSJed Brown          earlier, in the main routine.
580c4762a1bSJed Brown       */
5815f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesLocal(Hessian,1,&row,k,col,v,INSERT_VALUES));
582c4762a1bSJed Brown 
583c4762a1bSJed Brown     }
584c4762a1bSJed Brown   }
585c4762a1bSJed Brown 
586c4762a1bSJed Brown   /* Restore vectors */
5875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localX,&x));
5885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Left,&left));
5895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Top,&top));
5905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Bottom,&bottom));
5915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->Right,&right));
592c4762a1bSJed Brown 
593c4762a1bSJed Brown   /* Assemble the matrix */
5945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
5955f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
596c4762a1bSJed Brown 
5975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(199.0*xm*ym));
598c4762a1bSJed Brown   return 0;
599c4762a1bSJed Brown }
600c4762a1bSJed Brown 
601c4762a1bSJed Brown /* ------------------------------------------------------------------- */
602c4762a1bSJed Brown /*
603c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
604c4762a1bSJed Brown    the region.
605c4762a1bSJed Brown 
606c4762a1bSJed Brown    Input Parameter:
607c4762a1bSJed Brown .  user - user-defined application context
608c4762a1bSJed Brown 
609c4762a1bSJed Brown    Output Parameter:
610c4762a1bSJed Brown .  user - user-defined application context
611c4762a1bSJed Brown */
612c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
613c4762a1bSJed Brown {
614c4762a1bSJed Brown   PetscInt   i,j,k,maxits=5,limit=0;
615c4762a1bSJed Brown   PetscInt   xs,ys,xm,ym,gxs,gys,gxm,gym;
616c4762a1bSJed Brown   PetscInt   mx=user->mx,my=user->my;
617c4762a1bSJed Brown   PetscInt   bsize=0, lsize=0, tsize=0, rsize=0;
618c4762a1bSJed Brown   PetscReal  one=1.0, two=2.0, three=3.0, scl=1.0, tol=1e-10;
619c4762a1bSJed Brown   PetscReal  fnorm,det,hx,hy,xt=0,yt=0;
620c4762a1bSJed Brown   PetscReal  u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
621c4762a1bSJed Brown   PetscReal  b=-0.5, t=0.5, l=-0.5, r=0.5;
622c4762a1bSJed Brown   PetscReal  *boundary;
623c4762a1bSJed Brown   PetscBool  flg;
624c4762a1bSJed Brown   Vec        Bottom,Top,Right,Left;
625c4762a1bSJed Brown 
626c4762a1bSJed Brown   /* Get local mesh boundaries */
6275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
6285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
629c4762a1bSJed Brown 
630c4762a1bSJed Brown   bsize=xm+2;
631c4762a1bSJed Brown   lsize=ym+2;
632c4762a1bSJed Brown   rsize=ym+2;
633c4762a1bSJed Brown   tsize=xm+2;
634c4762a1bSJed Brown 
6355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom));
6365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top));
6375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left));
6385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right));
639c4762a1bSJed Brown 
640c4762a1bSJed Brown   user->Top=Top;
641c4762a1bSJed Brown   user->Left=Left;
642c4762a1bSJed Brown   user->Bottom=Bottom;
643c4762a1bSJed Brown   user->Right=Right;
644c4762a1bSJed Brown 
645c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
646c4762a1bSJed Brown 
647c4762a1bSJed Brown   for (j=0; j<4; j++) {
648c4762a1bSJed Brown     if (j==0) {
649c4762a1bSJed Brown       yt=b;
650c4762a1bSJed Brown       xt=l+hx*xs;
651c4762a1bSJed Brown       limit=bsize;
652c4762a1bSJed Brown       VecGetArray(Bottom,&boundary);
653c4762a1bSJed Brown     } else if (j==1) {
654c4762a1bSJed Brown       yt=t;
655c4762a1bSJed Brown       xt=l+hx*xs;
656c4762a1bSJed Brown       limit=tsize;
657c4762a1bSJed Brown       VecGetArray(Top,&boundary);
658c4762a1bSJed Brown     } else if (j==2) {
659c4762a1bSJed Brown       yt=b+hy*ys;
660c4762a1bSJed Brown       xt=l;
661c4762a1bSJed Brown       limit=lsize;
662c4762a1bSJed Brown       VecGetArray(Left,&boundary);
663c4762a1bSJed Brown     } else if (j==3) {
664c4762a1bSJed Brown       yt=b+hy*ys;
665c4762a1bSJed Brown       xt=r;
666c4762a1bSJed Brown       limit=rsize;
667c4762a1bSJed Brown       VecGetArray(Right,&boundary);
668c4762a1bSJed Brown     }
669c4762a1bSJed Brown 
670c4762a1bSJed Brown     for (i=0; i<limit; i++) {
671c4762a1bSJed Brown       u1=xt;
672c4762a1bSJed Brown       u2=-yt;
673c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
674c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
675c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
676c4762a1bSJed Brown         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
677c4762a1bSJed Brown         if (fnorm <= tol) break;
678c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
679c4762a1bSJed Brown         njac12=two*u1*u2;
680c4762a1bSJed Brown         njac21=-two*u1*u2;
681c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
682c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
683c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
684c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
685c4762a1bSJed Brown       }
686c4762a1bSJed Brown 
687c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
688c4762a1bSJed Brown       if (j==0 || j==1) {
689c4762a1bSJed Brown         xt=xt+hx;
690c4762a1bSJed Brown       } else if (j==2 || j==3) {
691c4762a1bSJed Brown         yt=yt+hy;
692c4762a1bSJed Brown       }
693c4762a1bSJed Brown     }
694c4762a1bSJed Brown     if (j==0) {
6955f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(Bottom,&boundary));
696c4762a1bSJed Brown     } else if (j==1) {
6975f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(Top,&boundary));
698c4762a1bSJed Brown     } else if (j==2) {
6995f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(Left,&boundary));
700c4762a1bSJed Brown     } else if (j==3) {
7015f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(Right,&boundary));
702c4762a1bSJed Brown     }
703c4762a1bSJed Brown   }
704c4762a1bSJed Brown 
705c4762a1bSJed Brown   /* Scale the boundary if desired */
706c4762a1bSJed Brown 
7075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg));
708c4762a1bSJed Brown   if (flg) {
7095f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(Bottom, scl));
710c4762a1bSJed Brown   }
7115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg));
712c4762a1bSJed Brown   if (flg) {
7135f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(Top, scl));
714c4762a1bSJed Brown   }
7155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg));
716c4762a1bSJed Brown   if (flg) {
7175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(Right, scl));
718c4762a1bSJed Brown   }
719c4762a1bSJed Brown 
7205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg));
721c4762a1bSJed Brown   if (flg) {
7225f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(Left, scl));
723c4762a1bSJed Brown   }
724c4762a1bSJed Brown   return 0;
725c4762a1bSJed Brown }
726c4762a1bSJed Brown 
727c4762a1bSJed Brown /* ------------------------------------------------------------------- */
728c4762a1bSJed Brown /*
729c4762a1bSJed Brown    MSA_Plate -  Calculates an obstacle for surface to stretch over.
730c4762a1bSJed Brown 
731c4762a1bSJed Brown    Input Parameter:
732c4762a1bSJed Brown .  user - user-defined application context
733c4762a1bSJed Brown 
734c4762a1bSJed Brown    Output Parameter:
735c4762a1bSJed Brown .  user - user-defined application context
736c4762a1bSJed Brown */
73770a7d78aSStefano Zampini static PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx)
73870a7d78aSStefano Zampini {
739c4762a1bSJed Brown   AppCtx         *user=(AppCtx *)ctx;
740c4762a1bSJed Brown   PetscInt       i,j,row;
741c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym;
742c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my, bmy, bmx;
743c4762a1bSJed Brown   PetscReal      t1,t2,t3;
744c4762a1bSJed Brown   PetscReal      *xl, lb=PETSC_NINFINITY, ub=PETSC_INFINITY;
745c4762a1bSJed Brown   PetscBool      cylinder;
746c4762a1bSJed Brown 
747c4762a1bSJed Brown   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
748c4762a1bSJed Brown   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
749c4762a1bSJed Brown   bmy=user->bmy; bmx=user->bmx;
750c4762a1bSJed Brown 
7515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
752c4762a1bSJed Brown 
7535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(XL, lb));
7545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(XU, ub));
755c4762a1bSJed Brown 
7565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(XL,&xl));
757c4762a1bSJed Brown 
7585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-cylinder",&cylinder));
759c4762a1bSJed Brown   /* Compute the optional lower box */
760c4762a1bSJed Brown   if (cylinder) {
761c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
762c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
763c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
764c4762a1bSJed Brown         t1=(2.0*i-mx)*bmy;
765c4762a1bSJed Brown         t2=(2.0*j-my)*bmx;
766c4762a1bSJed Brown         t3=bmx*bmx*bmy*bmy;
767c4762a1bSJed Brown         if (t1*t1 + t2*t2 <= t3) {
768c4762a1bSJed Brown           xl[row] = user->bheight;
769c4762a1bSJed Brown         }
770c4762a1bSJed Brown       }
771c4762a1bSJed Brown     }
772c4762a1bSJed Brown   } else {
773c4762a1bSJed Brown     /* Compute the optional lower box */
774c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
775c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
776c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
777c4762a1bSJed Brown         if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
778c4762a1bSJed Brown             j>=(my-bmy)/2 && j<my-(my-bmy)/2) {
779c4762a1bSJed Brown           xl[row] = user->bheight;
780c4762a1bSJed Brown         }
781c4762a1bSJed Brown       }
782c4762a1bSJed Brown     }
783c4762a1bSJed Brown   }
7845f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(XL,&xl));
785c4762a1bSJed Brown 
786c4762a1bSJed Brown   return 0;
787c4762a1bSJed Brown }
788c4762a1bSJed Brown 
789c4762a1bSJed Brown /* ------------------------------------------------------------------- */
790c4762a1bSJed Brown /*
791c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
792c4762a1bSJed Brown 
793c4762a1bSJed Brown    Input Parameters:
794c4762a1bSJed Brown .  user - user-defined application context
795c4762a1bSJed Brown .  X - vector for initial guess
796c4762a1bSJed Brown 
797c4762a1bSJed Brown    Output Parameters:
798c4762a1bSJed Brown .  X - newly computed initial guess
799c4762a1bSJed Brown */
800c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
801c4762a1bSJed Brown {
802c4762a1bSJed Brown   PetscInt       start=-1,i,j;
803c4762a1bSJed Brown   PetscReal      zero=0.0;
804c4762a1bSJed Brown   PetscBool      flg;
805c4762a1bSJed Brown 
8065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
807c4762a1bSJed Brown   if (flg && start==0) { /* The zero vector is reasonable */
8085f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(X, zero));
809c4762a1bSJed Brown   } else if (flg && start>0) { /* Try a random start between -0.5 and 0.5 */
810c4762a1bSJed Brown     PetscRandom rctx;  PetscReal np5=-0.5;
811c4762a1bSJed Brown 
8125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomCreate(MPI_COMM_WORLD,&rctx));
813c4762a1bSJed Brown     for (i=0; i<start; i++) {
8145f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(X, rctx));
815c4762a1bSJed Brown     }
8165f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomDestroy(&rctx));
8175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecShift(X, np5));
818c4762a1bSJed Brown 
819c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
820c4762a1bSJed Brown 
821c4762a1bSJed Brown     PetscInt row,xs,xm,gxs,gxm,ys,ym,gys,gym;
822c4762a1bSJed Brown     PetscInt mx=user->mx,my=user->my;
823c4762a1bSJed Brown     PetscReal *x,*left,*right,*bottom,*top;
824c4762a1bSJed Brown     Vec    localX = user->localX;
825c4762a1bSJed Brown 
826c4762a1bSJed Brown     /* Get local mesh boundaries */
8275f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
8285f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
829c4762a1bSJed Brown 
830c4762a1bSJed Brown     /* Get pointers to vector data */
8315f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(user->Top,&top));
8325f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(user->Bottom,&bottom));
8335f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(user->Left,&left));
8345f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(user->Right,&right));
835c4762a1bSJed Brown 
8365f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(localX,&x));
837c4762a1bSJed Brown     /* Perform local computations */
838c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
839c4762a1bSJed Brown       for (i=xs; i< xs+xm; i++) {
840c4762a1bSJed Brown         row=(j-gys)*gxm + (i-gxs);
8412da392ccSBarry Smith         x[row] = ((j+1)*bottom[i-xs+1]/my + (my-j+1)*top[i-xs+1]/(my+2)+(i+1)*left[j-ys+1]/mx + (mx-i+1)*right[j-ys+1]/(mx+2))/2.0;
842c4762a1bSJed Brown       }
843c4762a1bSJed Brown     }
844c4762a1bSJed Brown 
845c4762a1bSJed Brown     /* Restore vectors */
8465f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(localX,&x));
847c4762a1bSJed Brown 
8485f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(user->Left,&left));
8495f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(user->Top,&top));
8505f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(user->Bottom,&bottom));
8515f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(user->Right,&right));
852c4762a1bSJed Brown 
853c4762a1bSJed Brown     /* Scatter values into global vector */
8545f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X));
8555f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X));
856c4762a1bSJed Brown 
857c4762a1bSJed Brown   }
858c4762a1bSJed Brown   return 0;
859c4762a1bSJed Brown }
860c4762a1bSJed Brown 
861c4762a1bSJed Brown /* For testing matrix free submatrices */
862c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
863c4762a1bSJed Brown {
864c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
865c4762a1bSJed Brown   PetscFunctionBegin;
8665f80ce2aSJacob Faibussowitsch   CHKERRQ(FormHessian(tao,x,user->H,user->H,ptr));
867c4762a1bSJed Brown   PetscFunctionReturn(0);
868c4762a1bSJed Brown }
869c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
870c4762a1bSJed Brown {
871c4762a1bSJed Brown   void           *ptr;
872c4762a1bSJed Brown   AppCtx         *user;
873c4762a1bSJed Brown   PetscFunctionBegin;
8745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(H_shell,&ptr));
875c4762a1bSJed Brown   user = (AppCtx*)ptr;
8765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->H,X,Y));
877c4762a1bSJed Brown   PetscFunctionReturn(0);
878c4762a1bSJed Brown }
879c4762a1bSJed Brown 
880c4762a1bSJed Brown /*TEST
881c4762a1bSJed Brown 
882c4762a1bSJed Brown    build:
883c4762a1bSJed Brown       requires: !complex
884c4762a1bSJed Brown 
885c4762a1bSJed Brown    test:
886c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
887c4762a1bSJed Brown       requires: !single
888c4762a1bSJed Brown 
889c4762a1bSJed Brown    test:
890c4762a1bSJed Brown       suffix: 2
891c4762a1bSJed Brown       nsize: 2
892c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
893c4762a1bSJed Brown       requires: !single
894c4762a1bSJed Brown 
895c4762a1bSJed Brown    test:
896c4762a1bSJed Brown       suffix: 3
897c4762a1bSJed Brown       nsize: 3
898c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
899c4762a1bSJed Brown       requires: !single
900c4762a1bSJed Brown 
901c4762a1bSJed Brown    test:
902c4762a1bSJed Brown       suffix: 4
903c4762a1bSJed Brown       nsize: 3
904c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type mask -tao_type tron -tao_gatol 1.e-5
905c4762a1bSJed Brown       requires: !single
906c4762a1bSJed Brown 
907c4762a1bSJed Brown    test:
908c4762a1bSJed Brown       suffix: 5
909c4762a1bSJed Brown       nsize: 3
910c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -pc_type none -tao_type tron -tao_gatol 1.e-5
911c4762a1bSJed Brown       requires: !single
912c4762a1bSJed Brown 
913c4762a1bSJed Brown    test:
914c4762a1bSJed Brown       suffix: 6
915c4762a1bSJed Brown       nsize: 3
916c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -tao_type blmvm -tao_gatol 1.e-4
917c4762a1bSJed Brown       requires: !single
918c4762a1bSJed Brown 
919c4762a1bSJed Brown    test:
920c4762a1bSJed Brown       suffix: 7
921c4762a1bSJed Brown       nsize: 3
922c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -pc_type none -tao_type gpcg -tao_gatol 1.e-5
923c4762a1bSJed Brown       requires: !single
924c4762a1bSJed Brown 
925c4762a1bSJed Brown    test:
926c4762a1bSJed Brown       suffix: 8
927c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
928c4762a1bSJed Brown       requires: !single
929c4762a1bSJed Brown 
930c4762a1bSJed Brown    test:
931c4762a1bSJed Brown       suffix: 9
932c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
933c4762a1bSJed Brown       requires: !single
934c4762a1bSJed Brown 
935c4762a1bSJed Brown    test:
936c4762a1bSJed Brown       suffix: 10
937c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
938c4762a1bSJed Brown       requires: !single
939c4762a1bSJed Brown 
940c4762a1bSJed Brown    test:
941c4762a1bSJed Brown       suffix: 11
942c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
943c4762a1bSJed Brown       requires: !single
944c4762a1bSJed Brown 
945c4762a1bSJed Brown    test:
946c4762a1bSJed Brown       suffix: 12
947c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
948c4762a1bSJed Brown       requires: !single
949c4762a1bSJed Brown 
950c4762a1bSJed Brown    test:
951c4762a1bSJed Brown       suffix: 13
952c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
953c4762a1bSJed Brown       requires: !single
954c4762a1bSJed Brown 
955c4762a1bSJed Brown    test:
956c4762a1bSJed Brown       suffix: 14
957c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
958c4762a1bSJed Brown       requires: !single
959c4762a1bSJed Brown 
960c4762a1bSJed Brown    test:
961c4762a1bSJed Brown       suffix: 15
962c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
963c4762a1bSJed Brown       requires: !single
964c4762a1bSJed Brown 
965c4762a1bSJed Brown    test:
966c4762a1bSJed Brown      suffix: 16
967c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
968c4762a1bSJed Brown      requires: !single
969c4762a1bSJed Brown 
970c4762a1bSJed Brown    test:
971c4762a1bSJed Brown      suffix: 17
972c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
973c4762a1bSJed Brown      requires: !single
974c4762a1bSJed Brown 
97534ad9904SAlp Dener    test:
97634ad9904SAlp Dener      suffix: 18
97734ad9904SAlp Dener      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
97834ad9904SAlp Dener      requires: !single
97934ad9904SAlp Dener 
98034ad9904SAlp Dener    test:
98134ad9904SAlp Dener      suffix: 19
98234ad9904SAlp Dener      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
98334ad9904SAlp Dener      requires: !single
98434ad9904SAlp Dener 
98534ad9904SAlp Dener    test:
98634ad9904SAlp Dener      suffix: 20
98734ad9904SAlp Dener      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
98834ad9904SAlp Dener 
989c4762a1bSJed Brown TEST*/
990