xref: /petsc/src/tao/bound/tutorials/plate2.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown #include <petscdmda.h>
2c4762a1bSJed Brown #include <petsctao.h>
3c4762a1bSJed Brown 
49371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to \n\
5c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
6c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain, \n\
7c4762a1bSJed Brown boundary values along the edges of the domain, and a plate represented by \n\
8c4762a1bSJed Brown lower boundary conditions, the objective is to find the\n\
9c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
10c4762a1bSJed Brown The command line options are:\n\
11c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
12c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
13c4762a1bSJed Brown   -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\
14c4762a1bSJed Brown   -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\
15c4762a1bSJed Brown   -bheight <ht>, where <ht> = height of the plate\n\
16c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
17c4762a1bSJed Brown                for an average of the boundary conditions\n\n";
18c4762a1bSJed Brown 
19c4762a1bSJed Brown /*
20c4762a1bSJed Brown    User-defined application context - contains data needed by the
21c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
22c4762a1bSJed Brown    FormHessian().
23c4762a1bSJed Brown */
24c4762a1bSJed Brown typedef struct {
25c4762a1bSJed Brown   /* problem parameters */
26c4762a1bSJed Brown   PetscReal bheight;                  /* Height of plate under the surface */
27c4762a1bSJed Brown   PetscInt  mx, my;                   /* discretization in x, y directions */
28c4762a1bSJed Brown   PetscInt  bmx, bmy;                 /* Size of plate under the surface */
29c4762a1bSJed Brown   Vec       Bottom, Top, Left, Right; /* boundary values */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* Working space */
32c4762a1bSJed Brown   Vec localX, localV; /* ghosted local vector */
33c4762a1bSJed Brown   DM  dm;             /* distributed array data structure */
34c4762a1bSJed Brown   Mat H;
35c4762a1bSJed Brown } AppCtx;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown /* -------- User-defined Routines --------- */
38c4762a1bSJed Brown 
39c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
40c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
41c4762a1bSJed Brown static PetscErrorCode MSA_Plate(Vec, Vec, void *);
42c4762a1bSJed Brown PetscErrorCode        FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
43c4762a1bSJed Brown PetscErrorCode        FormHessian(Tao, Vec, Mat, Mat, void *);
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /* For testing matrix free submatrices */
46c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *);
47c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat, Vec, Vec);
48c4762a1bSJed Brown 
49*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
50*d71ae5a4SJacob Faibussowitsch {
51c4762a1bSJed Brown   PetscInt               Nx, Ny;    /* number of processors in x- and y- directions */
52c4762a1bSJed Brown   PetscInt               m, N;      /* number of local and global elements in vectors */
53c4762a1bSJed Brown   Vec                    x, xl, xu; /* solution vector  and bounds*/
54c4762a1bSJed Brown   PetscBool              flg;       /* A return variable when checking for user options */
55c4762a1bSJed Brown   Tao                    tao;       /* Tao solver context */
56c4762a1bSJed Brown   ISLocalToGlobalMapping isltog;    /* local-to-global mapping object */
57c4762a1bSJed Brown   Mat                    H_shell;   /* to test matrix-free submatrices */
58c4762a1bSJed Brown   AppCtx                 user;      /* user-defined work context */
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Initialize PETSc, TAO */
61327415f7SBarry Smith   PetscFunctionBeginUser;
629566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Specify default dimension of the problem */
659371c9d4SSatish Balay   user.mx      = 10;
669371c9d4SSatish Balay   user.my      = 10;
679371c9d4SSatish Balay   user.bheight = 0.1;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-bheight", &user.bheight, &flg));
73c4762a1bSJed Brown 
749371c9d4SSatish Balay   user.bmx = user.mx / 2;
759371c9d4SSatish Balay   user.bmy = user.my / 2;
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bmx", &user.bmx, &flg));
779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bmy", &user.bmy, &flg));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n---- Minimum Surface Area With Plate Problem -----\n"));
8063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mx:%" PetscInt_FMT ", my:%" PetscInt_FMT ", bmx:%" PetscInt_FMT ", bmy:%" PetscInt_FMT ", height:%g\n", user.mx, user.my, user.bmx, user.bmy, (double)user.bheight));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Calculate any derived values from parameters */
83c4762a1bSJed Brown   N = user.mx * user.my;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* Let Petsc determine the dimensions of the local vectors */
869371c9d4SSatish Balay   Nx = PETSC_DECIDE;
879371c9d4SSatish Balay   Ny = PETSC_DECIDE;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /*
90c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
91c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
92c4762a1bSJed Brown      the distributed array, Create the vectors.
93c4762a1bSJed Brown   */
949566063dSJacob Faibussowitsch   PetscCall(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));
959566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
969566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
97c4762a1bSJed Brown   /*
98c4762a1bSJed Brown      Extract global and local vectors from DM; The local vectors are
99c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
100c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
101c4762a1bSJed Brown      the same types.
102c4762a1bSJed Brown   */
1039566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm, &x)); /* Solution */
1049566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(user.dm, &user.localX));
1059566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.localX, &user.localV));
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &xl));
1089566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &xu));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* The TAO code begins here */
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /*
113c4762a1bSJed Brown      Create TAO solver and set desired solution method
114c4762a1bSJed Brown      The method must either be TAOTRON or TAOBLMVM
115c4762a1bSJed Brown      If TAOBLMVM is used, then hessian function is not called.
116c4762a1bSJed Brown   */
1179566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1189566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBLMVM));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* Set initial solution guess; */
1219566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user));
1229566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user, x));
1239566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Set routines for function, gradient and hessian evaluation */
1269566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
127c4762a1bSJed Brown 
1289566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &m));
1299566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(MPI_COMM_WORLD, m, m, N, N, 7, NULL, 3, NULL, &(user.H)));
1309566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
131c4762a1bSJed Brown 
1329566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(user.dm, &isltog));
1339566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(user.H, isltog, isltog));
1349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-matrixfree", &flg));
135c4762a1bSJed Brown   if (flg) {
1369566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PETSC_COMM_WORLD, m, m, N, N, (void *)&user, &H_shell));
1379566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(H_shell, MATOP_MULT, (void (*)(void))MyMatMult));
1389566063dSJacob Faibussowitsch     PetscCall(MatSetOption(H_shell, MAT_SYMMETRIC, PETSC_TRUE));
1399566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, H_shell, H_shell, MatrixFreeHessian, (void *)&user));
140c4762a1bSJed Brown   } else {
1419566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
142c4762a1bSJed Brown   }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Set Variable bounds */
1459566063dSJacob Faibussowitsch   PetscCall(MSA_Plate(xl, xu, (void *)&user));
1469566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, xl, xu));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Check for any tao command line options */
1499566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1529566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(TaoView(tao, PETSC_VIEWER_STDOUT_WORLD));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Free TAO data structures */
1579566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Free PETSc data structures */
1609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
1629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
1639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
1649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localX));
1659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localV));
1669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Bottom));
1679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Top));
1689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Left));
1699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Right));
1709566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
17148a46eb9SPierre Jolivet   if (flg) PetscCall(MatDestroy(&H_shell));
1729566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
173b122ec5aSJacob Faibussowitsch   return 0;
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).
177c4762a1bSJed Brown 
178c4762a1bSJed Brown     Input Parameters:
179c4762a1bSJed Brown .   tao     - the Tao context
180c4762a1bSJed Brown .   X      - input vector
181a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
182c4762a1bSJed Brown 
183c4762a1bSJed Brown     Output Parameters:
184c4762a1bSJed Brown .   fcn     - the function value
185c4762a1bSJed Brown .   G      - vector containing the newly evaluated gradient
186c4762a1bSJed Brown 
187c4762a1bSJed Brown    Notes:
188c4762a1bSJed Brown    In this case, we discretize the domain and Create triangles. The
189c4762a1bSJed Brown    surface of each triangle is planar, whose surface area can be easily
190c4762a1bSJed Brown    computed.  The total surface area is found by sweeping through the grid
191c4762a1bSJed Brown    and computing the surface area of the two triangles that have their
192c4762a1bSJed Brown    right angle at the grid point.  The diagonal line segments on the
193c4762a1bSJed Brown    grid that define the triangles run from top left to lower right.
194c4762a1bSJed Brown    The numbering of points starts at the lower left and runs left to
195c4762a1bSJed Brown    right, then bottom to top.
196c4762a1bSJed Brown */
197*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx)
198*d71ae5a4SJacob Faibussowitsch {
199c4762a1bSJed Brown   AppCtx    *user = (AppCtx *)userCtx;
200c4762a1bSJed Brown   PetscInt   i, j, row;
201c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
202c4762a1bSJed Brown   PetscInt   xs, xm, gxs, gxm, ys, ym, gys, gym;
203c4762a1bSJed Brown   PetscReal  ft   = 0;
204c4762a1bSJed Brown   PetscReal  zero = 0.0;
205c4762a1bSJed Brown   PetscReal  hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy;
206c4762a1bSJed Brown   PetscReal  rhx = mx + 1, rhy = my + 1;
207c4762a1bSJed Brown   PetscReal  f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
208c4762a1bSJed Brown   PetscReal  df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
209c4762a1bSJed Brown   PetscReal *g, *x, *left, *right, *bottom, *top;
210c4762a1bSJed Brown   Vec        localX = user->localX, localG = user->localV;
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /* Get local mesh boundaries */
2139566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
2149566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   /* Scatter ghost points to local vector */
2179566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
2189566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* Initialize vector to zero */
2219566063dSJacob Faibussowitsch   PetscCall(VecSet(localG, zero));
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* Get pointers to vector data */
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &x));
2259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localG, &g));
2269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top, &top));
2279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom, &bottom));
2289566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left, &left));
2299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right, &right));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
232c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
233c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
234c4762a1bSJed Brown       row = (j - gys) * gxm + (i - gxs);
235c4762a1bSJed Brown 
236c4762a1bSJed Brown       xc  = x[row];
237c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
238c4762a1bSJed Brown 
239c4762a1bSJed Brown       if (i == 0) { /* left side */
240c4762a1bSJed Brown         xl  = left[j - ys + 1];
241c4762a1bSJed Brown         xlt = left[j - ys + 2];
242c4762a1bSJed Brown       } else {
243c4762a1bSJed Brown         xl = x[row - 1];
244c4762a1bSJed Brown       }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown       if (j == 0) { /* bottom side */
247c4762a1bSJed Brown         xb  = bottom[i - xs + 1];
248c4762a1bSJed Brown         xrb = bottom[i - xs + 2];
249c4762a1bSJed Brown       } else {
250c4762a1bSJed Brown         xb = x[row - gxm];
251c4762a1bSJed Brown       }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown       if (i + 1 == gxs + gxm) { /* right side */
254c4762a1bSJed Brown         xr  = right[j - ys + 1];
255c4762a1bSJed Brown         xrb = right[j - ys];
256c4762a1bSJed Brown       } else {
257c4762a1bSJed Brown         xr = x[row + 1];
258c4762a1bSJed Brown       }
259c4762a1bSJed Brown 
260c4762a1bSJed Brown       if (j + 1 == gys + gym) { /* top side */
261c4762a1bSJed Brown         xt  = top[i - xs + 1];
262c4762a1bSJed Brown         xlt = top[i - xs];
263c4762a1bSJed Brown       } else {
264c4762a1bSJed Brown         xt = x[row + gxm];
265c4762a1bSJed Brown       }
266c4762a1bSJed Brown 
267ad540459SPierre Jolivet       if (i > gxs && j + 1 < gys + gym) xlt = x[row - 1 + gxm];
268ad540459SPierre Jolivet       if (j > gys && i + 1 < gxs + gxm) xrb = x[row + 1 - gxm];
269c4762a1bSJed Brown 
270c4762a1bSJed Brown       d1 = (xc - xl);
271c4762a1bSJed Brown       d2 = (xc - xr);
272c4762a1bSJed Brown       d3 = (xc - xt);
273c4762a1bSJed Brown       d4 = (xc - xb);
274c4762a1bSJed Brown       d5 = (xr - xrb);
275c4762a1bSJed Brown       d6 = (xrb - xb);
276c4762a1bSJed Brown       d7 = (xlt - xl);
277c4762a1bSJed Brown       d8 = (xt - xlt);
278c4762a1bSJed Brown 
279c4762a1bSJed Brown       df1dxc = d1 * hydhx;
280c4762a1bSJed Brown       df2dxc = (d1 * hydhx + d4 * hxdhy);
281c4762a1bSJed Brown       df3dxc = d3 * hxdhy;
282c4762a1bSJed Brown       df4dxc = (d2 * hydhx + d3 * hxdhy);
283c4762a1bSJed Brown       df5dxc = d2 * hydhx;
284c4762a1bSJed Brown       df6dxc = d4 * hxdhy;
285c4762a1bSJed Brown 
286c4762a1bSJed Brown       d1 *= rhx;
287c4762a1bSJed Brown       d2 *= rhx;
288c4762a1bSJed Brown       d3 *= rhy;
289c4762a1bSJed Brown       d4 *= rhy;
290c4762a1bSJed Brown       d5 *= rhy;
291c4762a1bSJed Brown       d6 *= rhx;
292c4762a1bSJed Brown       d7 *= rhy;
293c4762a1bSJed Brown       d8 *= rhx;
294c4762a1bSJed Brown 
295c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
296c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
297c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
298c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
299c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
300c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
301c4762a1bSJed Brown 
302c4762a1bSJed Brown       ft = ft + (f2 + f4);
303c4762a1bSJed Brown 
304c4762a1bSJed Brown       df1dxc /= f1;
305c4762a1bSJed Brown       df2dxc /= f2;
306c4762a1bSJed Brown       df3dxc /= f3;
307c4762a1bSJed Brown       df4dxc /= f4;
308c4762a1bSJed Brown       df5dxc /= f5;
309c4762a1bSJed Brown       df6dxc /= f6;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown       g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5;
312c4762a1bSJed Brown     }
313c4762a1bSJed Brown   }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
316c4762a1bSJed Brown   if (xs == 0) { /* left side */
317c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
318c4762a1bSJed Brown       d3 = (left[j - ys + 1] - left[j - ys + 2]) * rhy;
319c4762a1bSJed Brown       d2 = (left[j - ys + 1] - x[(j - gys) * gxm]) * rhx;
320c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown   }
323c4762a1bSJed Brown   if (ys == 0) { /* bottom side */
324c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
325c4762a1bSJed Brown       d2 = (bottom[i + 1 - xs] - bottom[i - xs + 2]) * rhx;
326c4762a1bSJed Brown       d3 = (bottom[i - xs + 1] - x[i - gxs]) * rhy;
327c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
328c4762a1bSJed Brown     }
329c4762a1bSJed Brown   }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   if (xs + xm == mx) { /* right side */
332c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
333c4762a1bSJed Brown       d1 = (x[(j + 1 - gys) * gxm - 1] - right[j - ys + 1]) * rhx;
334c4762a1bSJed Brown       d4 = (right[j - ys] - right[j - ys + 1]) * rhy;
335c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
336c4762a1bSJed Brown     }
337c4762a1bSJed Brown   }
338c4762a1bSJed Brown   if (ys + ym == my) { /* top side */
339c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
340c4762a1bSJed Brown       d1 = (x[(gym - 1) * gxm + i - gxs] - top[i - xs + 1]) * rhy;
341c4762a1bSJed Brown       d4 = (top[i - xs + 1] - top[i - xs]) * rhx;
342c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
343c4762a1bSJed Brown     }
344c4762a1bSJed Brown   }
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   if (ys == 0 && xs == 0) {
347c4762a1bSJed Brown     d1 = (left[0] - left[1]) * rhy;
348c4762a1bSJed Brown     d2 = (bottom[0] - bottom[1]) * rhx;
349c4762a1bSJed Brown     ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
350c4762a1bSJed Brown   }
351c4762a1bSJed Brown   if (ys + ym == my && xs + xm == mx) {
352c4762a1bSJed Brown     d1 = (right[ym + 1] - right[ym]) * rhy;
353c4762a1bSJed Brown     d2 = (top[xm + 1] - top[xm]) * rhx;
354c4762a1bSJed Brown     ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
355c4762a1bSJed Brown   }
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   ft = ft * area;
3589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
359c4762a1bSJed Brown 
360c4762a1bSJed Brown   /* Restore vectors */
3619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &x));
3629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localG, &g));
3639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left, &left));
3649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top, &top));
3659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom, &bottom));
3669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right, &right));
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   /* Scatter values to global vector */
3699566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(user->dm, localG, INSERT_VALUES, G));
3709566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(user->dm, localG, INSERT_VALUES, G));
371c4762a1bSJed Brown 
3729566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(70.0 * xm * ym));
373c4762a1bSJed Brown 
374c4762a1bSJed Brown   return 0;
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown /* ------------------------------------------------------------------- */
378c4762a1bSJed Brown /*
379c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
380c4762a1bSJed Brown 
381c4762a1bSJed Brown    Input Parameters:
382c4762a1bSJed Brown .  tao  - the Tao context
383c4762a1bSJed Brown .  x    - input vector
384a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
385c4762a1bSJed Brown 
386c4762a1bSJed Brown    Output Parameters:
387c4762a1bSJed Brown .  A    - Hessian matrix
388c4762a1bSJed Brown .  B    - optionally different preconditioning matrix
389c4762a1bSJed Brown 
390c4762a1bSJed Brown    Notes:
391c4762a1bSJed Brown    Due to mesh point reordering with DMs, we must always work
392c4762a1bSJed Brown    with the local mesh points, and then transform them to the new
393c4762a1bSJed Brown    global numbering with the local-to-global mapping.  We cannot work
394c4762a1bSJed Brown    directly with the global numbers for the original uniprocessor mesh!
395c4762a1bSJed Brown 
396c4762a1bSJed Brown    Two methods are available for imposing this transformation
397c4762a1bSJed Brown    when setting matrix entries:
398c4762a1bSJed Brown      (A) MatSetValuesLocal(), using the local ordering (including
399c4762a1bSJed Brown          ghost points!)
400c4762a1bSJed Brown          - Do the following two steps once, before calling TaoSolve()
401c4762a1bSJed Brown            - Use DMGetISLocalToGlobalMapping() to extract the
402c4762a1bSJed Brown              local-to-global map from the DM
403c4762a1bSJed Brown            - Associate this map with the matrix by calling
404c4762a1bSJed Brown              MatSetLocalToGlobalMapping()
405c4762a1bSJed Brown          - Then set matrix entries using the local ordering
406c4762a1bSJed Brown            by calling MatSetValuesLocal()
407c4762a1bSJed Brown      (B) MatSetValues(), using the global ordering
408c4762a1bSJed Brown          - Use DMGetGlobalIndices() to extract the local-to-global map
409c4762a1bSJed Brown          - Then apply this map explicitly yourself
410c4762a1bSJed Brown          - Set matrix entries using the global ordering by calling
411c4762a1bSJed Brown            MatSetValues()
412c4762a1bSJed Brown    Option (A) seems cleaner/easier in many cases, and is the procedure
413c4762a1bSJed Brown    used in this example.
414c4762a1bSJed Brown */
415*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat Hptr, Mat Hessian, void *ptr)
416*d71ae5a4SJacob Faibussowitsch {
417c4762a1bSJed Brown   AppCtx    *user = (AppCtx *)ptr;
418c4762a1bSJed Brown   PetscInt   i, j, k, row;
419c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
420c4762a1bSJed Brown   PetscInt   xs, xm, gxs, gxm, ys, ym, gys, gym, col[7];
421c4762a1bSJed Brown   PetscReal  hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
422c4762a1bSJed Brown   PetscReal  rhx = mx + 1, rhy = my + 1;
423c4762a1bSJed Brown   PetscReal  f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
424c4762a1bSJed Brown   PetscReal  hl, hr, ht, hb, hc, htl, hbr;
425c4762a1bSJed Brown   PetscReal *x, *left, *right, *bottom, *top;
426c4762a1bSJed Brown   PetscReal  v[7];
427c4762a1bSJed Brown   Vec        localX = user->localX;
428c4762a1bSJed Brown   PetscBool  assembled;
429c4762a1bSJed Brown 
430c4762a1bSJed Brown   /* Set various matrix options */
4319566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   /* Initialize matrix entries to zero */
4349566063dSJacob Faibussowitsch   PetscCall(MatAssembled(Hessian, &assembled));
4359566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(Hessian));
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   /* Get local mesh boundaries */
4389566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
4399566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /* Scatter ghost points to local vector */
4429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
4439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /* Get pointers to vector data */
4469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &x));
4479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top, &top));
4489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom, &bottom));
4499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left, &left));
4509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right, &right));
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
455c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
456c4762a1bSJed Brown       row = (j - gys) * gxm + (i - gxs);
457c4762a1bSJed Brown 
458c4762a1bSJed Brown       xc  = x[row];
459c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
460c4762a1bSJed Brown 
461c4762a1bSJed Brown       /* Left side */
462c4762a1bSJed Brown       if (i == gxs) {
463c4762a1bSJed Brown         xl  = left[j - ys + 1];
464c4762a1bSJed Brown         xlt = left[j - ys + 2];
465c4762a1bSJed Brown       } else {
466c4762a1bSJed Brown         xl = x[row - 1];
467c4762a1bSJed Brown       }
468c4762a1bSJed Brown 
469c4762a1bSJed Brown       if (j == gys) {
470c4762a1bSJed Brown         xb  = bottom[i - xs + 1];
471c4762a1bSJed Brown         xrb = bottom[i - xs + 2];
472c4762a1bSJed Brown       } else {
473c4762a1bSJed Brown         xb = x[row - gxm];
474c4762a1bSJed Brown       }
475c4762a1bSJed Brown 
476c4762a1bSJed Brown       if (i + 1 == gxs + gxm) {
477c4762a1bSJed Brown         xr  = right[j - ys + 1];
478c4762a1bSJed Brown         xrb = right[j - ys];
479c4762a1bSJed Brown       } else {
480c4762a1bSJed Brown         xr = x[row + 1];
481c4762a1bSJed Brown       }
482c4762a1bSJed Brown 
483c4762a1bSJed Brown       if (j + 1 == gys + gym) {
484c4762a1bSJed Brown         xt  = top[i - xs + 1];
485c4762a1bSJed Brown         xlt = top[i - xs];
486c4762a1bSJed Brown       } else {
487c4762a1bSJed Brown         xt = x[row + gxm];
488c4762a1bSJed Brown       }
489c4762a1bSJed Brown 
490ad540459SPierre Jolivet       if (i > gxs && j + 1 < gys + gym) xlt = x[row - 1 + gxm];
491ad540459SPierre Jolivet       if (j > gys && i + 1 < gxs + gxm) xrb = x[row + 1 - gxm];
492c4762a1bSJed Brown 
493c4762a1bSJed Brown       d1 = (xc - xl) * rhx;
494c4762a1bSJed Brown       d2 = (xc - xr) * rhx;
495c4762a1bSJed Brown       d3 = (xc - xt) * rhy;
496c4762a1bSJed Brown       d4 = (xc - xb) * rhy;
497c4762a1bSJed Brown       d5 = (xrb - xr) * rhy;
498c4762a1bSJed Brown       d6 = (xrb - xb) * rhx;
499c4762a1bSJed Brown       d7 = (xlt - xl) * rhy;
500c4762a1bSJed Brown       d8 = (xlt - xt) * rhx;
501c4762a1bSJed Brown 
502c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
503c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
504c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
505c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
506c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
507c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
508c4762a1bSJed Brown 
5099371c9d4SSatish Balay       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
5109371c9d4SSatish Balay       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
5119371c9d4SSatish Balay       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
5129371c9d4SSatish Balay       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
513c4762a1bSJed Brown 
514c4762a1bSJed Brown       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
515c4762a1bSJed Brown       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
516c4762a1bSJed Brown 
5179371c9d4SSatish Balay       hc = hydhx * (1.0 + d7 * d7) / (f1 * f1 * f1) + hxdhy * (1.0 + d8 * d8) / (f3 * f3 * f3) + hydhx * (1.0 + d5 * d5) / (f5 * f5 * f5) + hxdhy * (1.0 + d6 * d6) / (f6 * f6 * f6) + (hxdhy * (1.0 + d1 * d1) + hydhx * (1.0 + d4 * d4) - 2 * d1 * d4) / (f2 * f2 * f2) + (hxdhy * (1.0 + d2 * d2) + hydhx * (1.0 + d3 * d3) - 2 * d2 * d3) / (f4 * f4 * f4);
518c4762a1bSJed Brown 
5199371c9d4SSatish Balay       hl *= 0.5;
5209371c9d4SSatish Balay       hr *= 0.5;
5219371c9d4SSatish Balay       ht *= 0.5;
5229371c9d4SSatish Balay       hb *= 0.5;
5239371c9d4SSatish Balay       hbr *= 0.5;
5249371c9d4SSatish Balay       htl *= 0.5;
5259371c9d4SSatish Balay       hc *= 0.5;
526c4762a1bSJed Brown 
527c4762a1bSJed Brown       k = 0;
528c4762a1bSJed Brown       if (j > 0) {
5299371c9d4SSatish Balay         v[k]   = hb;
5309371c9d4SSatish Balay         col[k] = row - gxm;
5319371c9d4SSatish Balay         k++;
532c4762a1bSJed Brown       }
533c4762a1bSJed Brown 
534c4762a1bSJed Brown       if (j > 0 && i < mx - 1) {
5359371c9d4SSatish Balay         v[k]   = hbr;
5369371c9d4SSatish Balay         col[k] = row - gxm + 1;
5379371c9d4SSatish Balay         k++;
538c4762a1bSJed Brown       }
539c4762a1bSJed Brown 
540c4762a1bSJed Brown       if (i > 0) {
5419371c9d4SSatish Balay         v[k]   = hl;
5429371c9d4SSatish Balay         col[k] = row - 1;
5439371c9d4SSatish Balay         k++;
544c4762a1bSJed Brown       }
545c4762a1bSJed Brown 
5469371c9d4SSatish Balay       v[k]   = hc;
5479371c9d4SSatish Balay       col[k] = row;
5489371c9d4SSatish Balay       k++;
549c4762a1bSJed Brown 
550c4762a1bSJed Brown       if (i < mx - 1) {
5519371c9d4SSatish Balay         v[k]   = hr;
5529371c9d4SSatish Balay         col[k] = row + 1;
5539371c9d4SSatish Balay         k++;
554c4762a1bSJed Brown       }
555c4762a1bSJed Brown 
556c4762a1bSJed Brown       if (i > 0 && j < my - 1) {
5579371c9d4SSatish Balay         v[k]   = htl;
5589371c9d4SSatish Balay         col[k] = row + gxm - 1;
5599371c9d4SSatish Balay         k++;
560c4762a1bSJed Brown       }
561c4762a1bSJed Brown 
562c4762a1bSJed Brown       if (j < my - 1) {
5639371c9d4SSatish Balay         v[k]   = ht;
5649371c9d4SSatish Balay         col[k] = row + gxm;
5659371c9d4SSatish Balay         k++;
566c4762a1bSJed Brown       }
567c4762a1bSJed Brown 
568c4762a1bSJed Brown       /*
569c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
570c4762a1bSJed Brown          earlier, in the main routine.
571c4762a1bSJed Brown       */
5729566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Hessian, 1, &row, k, col, v, INSERT_VALUES));
573c4762a1bSJed Brown     }
574c4762a1bSJed Brown   }
575c4762a1bSJed Brown 
576c4762a1bSJed Brown   /* Restore vectors */
5779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &x));
5789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left, &left));
5799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top, &top));
5809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom, &bottom));
5819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right, &right));
582c4762a1bSJed Brown 
583c4762a1bSJed Brown   /* Assemble the matrix */
5849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY));
5859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY));
586c4762a1bSJed Brown 
5879566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0 * xm * ym));
588c4762a1bSJed Brown   return 0;
589c4762a1bSJed Brown }
590c4762a1bSJed Brown 
591c4762a1bSJed Brown /* ------------------------------------------------------------------- */
592c4762a1bSJed Brown /*
593c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
594c4762a1bSJed Brown    the region.
595c4762a1bSJed Brown 
596c4762a1bSJed Brown    Input Parameter:
597c4762a1bSJed Brown .  user - user-defined application context
598c4762a1bSJed Brown 
599c4762a1bSJed Brown    Output Parameter:
600c4762a1bSJed Brown .  user - user-defined application context
601c4762a1bSJed Brown */
602*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user)
603*d71ae5a4SJacob Faibussowitsch {
604c4762a1bSJed Brown   PetscInt   i, j, k, maxits = 5, limit = 0;
605c4762a1bSJed Brown   PetscInt   xs, ys, xm, ym, gxs, gys, gxm, gym;
606c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
607c4762a1bSJed Brown   PetscInt   bsize = 0, lsize = 0, tsize = 0, rsize = 0;
608c4762a1bSJed Brown   PetscReal  one = 1.0, two = 2.0, three = 3.0, scl = 1.0, tol = 1e-10;
609c4762a1bSJed Brown   PetscReal  fnorm, det, hx, hy, xt = 0, yt = 0;
610c4762a1bSJed Brown   PetscReal  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
611c4762a1bSJed Brown   PetscReal  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
612c4762a1bSJed Brown   PetscReal *boundary;
613c4762a1bSJed Brown   PetscBool  flg;
614c4762a1bSJed Brown   Vec        Bottom, Top, Right, Left;
615c4762a1bSJed Brown 
616c4762a1bSJed Brown   /* Get local mesh boundaries */
6179566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
6189566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
619c4762a1bSJed Brown 
620c4762a1bSJed Brown   bsize = xm + 2;
621c4762a1bSJed Brown   lsize = ym + 2;
622c4762a1bSJed Brown   rsize = ym + 2;
623c4762a1bSJed Brown   tsize = xm + 2;
624c4762a1bSJed Brown 
6259566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, bsize, PETSC_DECIDE, &Bottom));
6269566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, tsize, PETSC_DECIDE, &Top));
6279566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, lsize, PETSC_DECIDE, &Left));
6289566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, rsize, PETSC_DECIDE, &Right));
629c4762a1bSJed Brown 
630c4762a1bSJed Brown   user->Top    = Top;
631c4762a1bSJed Brown   user->Left   = Left;
632c4762a1bSJed Brown   user->Bottom = Bottom;
633c4762a1bSJed Brown   user->Right  = Right;
634c4762a1bSJed Brown 
6359371c9d4SSatish Balay   hx = (r - l) / (mx + 1);
6369371c9d4SSatish Balay   hy = (t - b) / (my + 1);
637c4762a1bSJed Brown 
638c4762a1bSJed Brown   for (j = 0; j < 4; j++) {
639c4762a1bSJed Brown     if (j == 0) {
640c4762a1bSJed Brown       yt    = b;
641c4762a1bSJed Brown       xt    = l + hx * xs;
642c4762a1bSJed Brown       limit = bsize;
643c4762a1bSJed Brown       VecGetArray(Bottom, &boundary);
644c4762a1bSJed Brown     } else if (j == 1) {
645c4762a1bSJed Brown       yt    = t;
646c4762a1bSJed Brown       xt    = l + hx * xs;
647c4762a1bSJed Brown       limit = tsize;
648c4762a1bSJed Brown       VecGetArray(Top, &boundary);
649c4762a1bSJed Brown     } else if (j == 2) {
650c4762a1bSJed Brown       yt    = b + hy * ys;
651c4762a1bSJed Brown       xt    = l;
652c4762a1bSJed Brown       limit = lsize;
653c4762a1bSJed Brown       VecGetArray(Left, &boundary);
654c4762a1bSJed Brown     } else if (j == 3) {
655c4762a1bSJed Brown       yt    = b + hy * ys;
656c4762a1bSJed Brown       xt    = r;
657c4762a1bSJed Brown       limit = rsize;
658c4762a1bSJed Brown       VecGetArray(Right, &boundary);
659c4762a1bSJed Brown     }
660c4762a1bSJed Brown 
661c4762a1bSJed Brown     for (i = 0; i < limit; i++) {
662c4762a1bSJed Brown       u1 = xt;
663c4762a1bSJed Brown       u2 = -yt;
664c4762a1bSJed Brown       for (k = 0; k < maxits; k++) {
665c4762a1bSJed Brown         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
666c4762a1bSJed Brown         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
667c4762a1bSJed Brown         fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2);
668c4762a1bSJed Brown         if (fnorm <= tol) break;
669c4762a1bSJed Brown         njac11 = one + u2 * u2 - u1 * u1;
670c4762a1bSJed Brown         njac12 = two * u1 * u2;
671c4762a1bSJed Brown         njac21 = -two * u1 * u2;
672c4762a1bSJed Brown         njac22 = -one - u1 * u1 + u2 * u2;
673c4762a1bSJed Brown         det    = njac11 * njac22 - njac21 * njac12;
674c4762a1bSJed Brown         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
675c4762a1bSJed Brown         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
676c4762a1bSJed Brown       }
677c4762a1bSJed Brown 
678c4762a1bSJed Brown       boundary[i] = u1 * u1 - u2 * u2;
679c4762a1bSJed Brown       if (j == 0 || j == 1) {
680c4762a1bSJed Brown         xt = xt + hx;
681c4762a1bSJed Brown       } else if (j == 2 || j == 3) {
682c4762a1bSJed Brown         yt = yt + hy;
683c4762a1bSJed Brown       }
684c4762a1bSJed Brown     }
685c4762a1bSJed Brown     if (j == 0) {
6869566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Bottom, &boundary));
687c4762a1bSJed Brown     } else if (j == 1) {
6889566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Top, &boundary));
689c4762a1bSJed Brown     } else if (j == 2) {
6909566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Left, &boundary));
691c4762a1bSJed Brown     } else if (j == 3) {
6929566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Right, &boundary));
693c4762a1bSJed Brown     }
694c4762a1bSJed Brown   }
695c4762a1bSJed Brown 
696c4762a1bSJed Brown   /* Scale the boundary if desired */
697c4762a1bSJed Brown 
6989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg));
6991baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Bottom, scl));
7009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg));
7011baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Top, scl));
7029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg));
7031baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Right, scl));
704c4762a1bSJed Brown 
7059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg));
7061baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Left, scl));
707c4762a1bSJed Brown   return 0;
708c4762a1bSJed Brown }
709c4762a1bSJed Brown 
710c4762a1bSJed Brown /* ------------------------------------------------------------------- */
711c4762a1bSJed Brown /*
712c4762a1bSJed Brown    MSA_Plate -  Calculates an obstacle for surface to stretch over.
713c4762a1bSJed Brown 
714c4762a1bSJed Brown    Input Parameter:
715c4762a1bSJed Brown .  user - user-defined application context
716c4762a1bSJed Brown 
717c4762a1bSJed Brown    Output Parameter:
718c4762a1bSJed Brown .  user - user-defined application context
719c4762a1bSJed Brown */
720*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_Plate(Vec XL, Vec XU, void *ctx)
721*d71ae5a4SJacob Faibussowitsch {
722c4762a1bSJed Brown   AppCtx    *user = (AppCtx *)ctx;
723c4762a1bSJed Brown   PetscInt   i, j, row;
724c4762a1bSJed Brown   PetscInt   xs, ys, xm, ym;
725c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my, bmy, bmx;
726c4762a1bSJed Brown   PetscReal  t1, t2, t3;
727c4762a1bSJed Brown   PetscReal *xl, lb = PETSC_NINFINITY, ub = PETSC_INFINITY;
728c4762a1bSJed Brown   PetscBool  cylinder;
729c4762a1bSJed Brown 
7309371c9d4SSatish Balay   user->bmy = PetscMax(0, user->bmy);
7319371c9d4SSatish Balay   user->bmy = PetscMin(my, user->bmy);
7329371c9d4SSatish Balay   user->bmx = PetscMax(0, user->bmx);
7339371c9d4SSatish Balay   user->bmx = PetscMin(mx, user->bmx);
7349371c9d4SSatish Balay   bmy       = user->bmy;
7359371c9d4SSatish Balay   bmx       = user->bmx;
736c4762a1bSJed Brown 
7379566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
738c4762a1bSJed Brown 
7399566063dSJacob Faibussowitsch   PetscCall(VecSet(XL, lb));
7409566063dSJacob Faibussowitsch   PetscCall(VecSet(XU, ub));
741c4762a1bSJed Brown 
7429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(XL, &xl));
743c4762a1bSJed Brown 
7449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-cylinder", &cylinder));
745c4762a1bSJed Brown   /* Compute the optional lower box */
746c4762a1bSJed Brown   if (cylinder) {
747c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
748c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
749c4762a1bSJed Brown         row = (j - ys) * xm + (i - xs);
750c4762a1bSJed Brown         t1  = (2.0 * i - mx) * bmy;
751c4762a1bSJed Brown         t2  = (2.0 * j - my) * bmx;
752c4762a1bSJed Brown         t3  = bmx * bmx * bmy * bmy;
753ad540459SPierre Jolivet         if (t1 * t1 + t2 * t2 <= t3) xl[row] = user->bheight;
754c4762a1bSJed Brown       }
755c4762a1bSJed Brown     }
756c4762a1bSJed Brown   } else {
757c4762a1bSJed Brown     /* Compute the optional lower box */
758c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
759c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
760c4762a1bSJed Brown         row = (j - ys) * xm + (i - xs);
761ad540459SPierre Jolivet         if (i >= (mx - bmx) / 2 && i < mx - (mx - bmx) / 2 && j >= (my - bmy) / 2 && j < my - (my - bmy) / 2) xl[row] = user->bheight;
762c4762a1bSJed Brown       }
763c4762a1bSJed Brown     }
764c4762a1bSJed Brown   }
7659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(XL, &xl));
766c4762a1bSJed Brown 
767c4762a1bSJed Brown   return 0;
768c4762a1bSJed Brown }
769c4762a1bSJed Brown 
770c4762a1bSJed Brown /* ------------------------------------------------------------------- */
771c4762a1bSJed Brown /*
772c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
773c4762a1bSJed Brown 
774c4762a1bSJed Brown    Input Parameters:
775c4762a1bSJed Brown .  user - user-defined application context
776c4762a1bSJed Brown .  X - vector for initial guess
777c4762a1bSJed Brown 
778c4762a1bSJed Brown    Output Parameters:
779c4762a1bSJed Brown .  X - newly computed initial guess
780c4762a1bSJed Brown */
781*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X)
782*d71ae5a4SJacob Faibussowitsch {
783c4762a1bSJed Brown   PetscInt  start = -1, i, j;
784c4762a1bSJed Brown   PetscReal zero  = 0.0;
785c4762a1bSJed Brown   PetscBool flg;
786c4762a1bSJed Brown 
7879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg));
788c4762a1bSJed Brown   if (flg && start == 0) { /* The zero vector is reasonable */
7899566063dSJacob Faibussowitsch     PetscCall(VecSet(X, zero));
790c4762a1bSJed Brown   } else if (flg && start > 0) { /* Try a random start between -0.5 and 0.5 */
7919371c9d4SSatish Balay     PetscRandom rctx;
7929371c9d4SSatish Balay     PetscReal   np5 = -0.5;
793c4762a1bSJed Brown 
7949566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(MPI_COMM_WORLD, &rctx));
79548a46eb9SPierre Jolivet     for (i = 0; i < start; i++) PetscCall(VecSetRandom(X, rctx));
7969566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rctx));
7979566063dSJacob Faibussowitsch     PetscCall(VecShift(X, np5));
798c4762a1bSJed Brown 
799c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
800c4762a1bSJed Brown 
801c4762a1bSJed Brown     PetscInt   row, xs, xm, gxs, gxm, ys, ym, gys, gym;
802c4762a1bSJed Brown     PetscInt   mx = user->mx, my = user->my;
803c4762a1bSJed Brown     PetscReal *x, *left, *right, *bottom, *top;
804c4762a1bSJed Brown     Vec        localX = user->localX;
805c4762a1bSJed Brown 
806c4762a1bSJed Brown     /* Get local mesh boundaries */
8079566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
8089566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
809c4762a1bSJed Brown 
810c4762a1bSJed Brown     /* Get pointers to vector data */
8119566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Top, &top));
8129566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Bottom, &bottom));
8139566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Left, &left));
8149566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Right, &right));
815c4762a1bSJed Brown 
8169566063dSJacob Faibussowitsch     PetscCall(VecGetArray(localX, &x));
817c4762a1bSJed Brown     /* Perform local computations */
818c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
819c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
820c4762a1bSJed Brown         row    = (j - gys) * gxm + (i - gxs);
8212da392ccSBarry 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;
822c4762a1bSJed Brown       }
823c4762a1bSJed Brown     }
824c4762a1bSJed Brown 
825c4762a1bSJed Brown     /* Restore vectors */
8269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(localX, &x));
827c4762a1bSJed Brown 
8289566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Left, &left));
8299566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Top, &top));
8309566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Bottom, &bottom));
8319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Right, &right));
832c4762a1bSJed Brown 
833c4762a1bSJed Brown     /* Scatter values into global vector */
8349566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(user->dm, localX, INSERT_VALUES, X));
8359566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(user->dm, localX, INSERT_VALUES, X));
836c4762a1bSJed Brown   }
837c4762a1bSJed Brown   return 0;
838c4762a1bSJed Brown }
839c4762a1bSJed Brown 
840c4762a1bSJed Brown /* For testing matrix free submatrices */
841*d71ae5a4SJacob Faibussowitsch PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
842*d71ae5a4SJacob Faibussowitsch {
843c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
844c4762a1bSJed Brown   PetscFunctionBegin;
8459566063dSJacob Faibussowitsch   PetscCall(FormHessian(tao, x, user->H, user->H, ptr));
846c4762a1bSJed Brown   PetscFunctionReturn(0);
847c4762a1bSJed Brown }
848*d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
849*d71ae5a4SJacob Faibussowitsch {
850c4762a1bSJed Brown   void   *ptr;
851c4762a1bSJed Brown   AppCtx *user;
852c4762a1bSJed Brown   PetscFunctionBegin;
8539566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(H_shell, &ptr));
854c4762a1bSJed Brown   user = (AppCtx *)ptr;
8559566063dSJacob Faibussowitsch   PetscCall(MatMult(user->H, X, Y));
856c4762a1bSJed Brown   PetscFunctionReturn(0);
857c4762a1bSJed Brown }
858c4762a1bSJed Brown 
859c4762a1bSJed Brown /*TEST
860c4762a1bSJed Brown 
861c4762a1bSJed Brown    build:
862c4762a1bSJed Brown       requires: !complex
863c4762a1bSJed Brown 
864c4762a1bSJed Brown    test:
865c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
866c4762a1bSJed Brown       requires: !single
867c4762a1bSJed Brown 
868c4762a1bSJed Brown    test:
869c4762a1bSJed Brown       suffix: 2
870c4762a1bSJed Brown       nsize: 2
871c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
872c4762a1bSJed Brown       requires: !single
873c4762a1bSJed Brown 
874c4762a1bSJed Brown    test:
875c4762a1bSJed Brown       suffix: 3
876c4762a1bSJed Brown       nsize: 3
877c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
878c4762a1bSJed Brown       requires: !single
879c4762a1bSJed Brown 
880c4762a1bSJed Brown    test:
881c4762a1bSJed Brown       suffix: 4
882c4762a1bSJed Brown       nsize: 3
883c4762a1bSJed 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
884c4762a1bSJed Brown       requires: !single
885c4762a1bSJed Brown 
886c4762a1bSJed Brown    test:
887c4762a1bSJed Brown       suffix: 5
888c4762a1bSJed Brown       nsize: 3
889c4762a1bSJed 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
890c4762a1bSJed Brown       requires: !single
891c4762a1bSJed Brown 
892c4762a1bSJed Brown    test:
893c4762a1bSJed Brown       suffix: 6
894c4762a1bSJed Brown       nsize: 3
895c4762a1bSJed 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
896c4762a1bSJed Brown       requires: !single
897c4762a1bSJed Brown 
898c4762a1bSJed Brown    test:
899c4762a1bSJed Brown       suffix: 7
900c4762a1bSJed Brown       nsize: 3
901c4762a1bSJed 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
902c4762a1bSJed Brown       requires: !single
903c4762a1bSJed Brown 
904c4762a1bSJed Brown    test:
905c4762a1bSJed Brown       suffix: 8
906c4762a1bSJed 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
907c4762a1bSJed Brown       requires: !single
908c4762a1bSJed Brown 
909c4762a1bSJed Brown    test:
910c4762a1bSJed Brown       suffix: 9
911c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
912c4762a1bSJed Brown       requires: !single
913c4762a1bSJed Brown 
914c4762a1bSJed Brown    test:
915c4762a1bSJed Brown       suffix: 10
916c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
917c4762a1bSJed Brown       requires: !single
918c4762a1bSJed Brown 
919c4762a1bSJed Brown    test:
920c4762a1bSJed Brown       suffix: 11
921c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
922c4762a1bSJed Brown       requires: !single
923c4762a1bSJed Brown 
924c4762a1bSJed Brown    test:
925c4762a1bSJed Brown       suffix: 12
926c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
927c4762a1bSJed Brown       requires: !single
928c4762a1bSJed Brown 
929c4762a1bSJed Brown    test:
930c4762a1bSJed Brown       suffix: 13
931c4762a1bSJed 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
932c4762a1bSJed Brown       requires: !single
933c4762a1bSJed Brown 
934c4762a1bSJed Brown    test:
935c4762a1bSJed Brown       suffix: 14
936c4762a1bSJed 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
937c4762a1bSJed Brown       requires: !single
938c4762a1bSJed Brown 
939c4762a1bSJed Brown    test:
940c4762a1bSJed Brown       suffix: 15
941c4762a1bSJed 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
942c4762a1bSJed Brown       requires: !single
943c4762a1bSJed Brown 
944c4762a1bSJed Brown    test:
945c4762a1bSJed Brown      suffix: 16
946c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
947c4762a1bSJed Brown      requires: !single
948c4762a1bSJed Brown 
949c4762a1bSJed Brown    test:
950c4762a1bSJed Brown      suffix: 17
951c4762a1bSJed 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
952c4762a1bSJed Brown      requires: !single
953c4762a1bSJed Brown 
95434ad9904SAlp Dener    test:
95534ad9904SAlp Dener      suffix: 18
95634ad9904SAlp 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
95734ad9904SAlp Dener      requires: !single
95834ad9904SAlp Dener 
95934ad9904SAlp Dener    test:
96034ad9904SAlp Dener      suffix: 19
96134ad9904SAlp 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
96234ad9904SAlp Dener      requires: !single
96334ad9904SAlp Dener 
96434ad9904SAlp Dener    test:
96534ad9904SAlp Dener      suffix: 20
96634ad9904SAlp 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
96734ad9904SAlp Dener 
968c4762a1bSJed Brown TEST*/
969