xref: /petsc/src/tao/bound/tutorials/plate2.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
499371c9d4SSatish Balay int main(int argc, char **argv) {
50c4762a1bSJed Brown   PetscInt               Nx, Ny;    /* number of processors in x- and y- directions */
51c4762a1bSJed Brown   PetscInt               m, N;      /* number of local and global elements in vectors */
52c4762a1bSJed Brown   Vec                    x, xl, xu; /* solution vector  and bounds*/
53c4762a1bSJed Brown   PetscBool              flg;       /* A return variable when checking for user options */
54c4762a1bSJed Brown   Tao                    tao;       /* Tao solver context */
55c4762a1bSJed Brown   ISLocalToGlobalMapping isltog;    /* local-to-global mapping object */
56c4762a1bSJed Brown   Mat                    H_shell;   /* to test matrix-free submatrices */
57c4762a1bSJed Brown   AppCtx                 user;      /* user-defined work context */
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Initialize PETSc, TAO */
60327415f7SBarry Smith   PetscFunctionBeginUser;
619566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* Specify default dimension of the problem */
649371c9d4SSatish Balay   user.mx      = 10;
659371c9d4SSatish Balay   user.my      = 10;
669371c9d4SSatish Balay   user.bheight = 0.1;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-bheight", &user.bheight, &flg));
72c4762a1bSJed Brown 
739371c9d4SSatish Balay   user.bmx = user.mx / 2;
749371c9d4SSatish Balay   user.bmy = user.my / 2;
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bmx", &user.bmx, &flg));
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bmy", &user.bmy, &flg));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n---- Minimum Surface Area With Plate Problem -----\n"));
7963a3b9bcSJacob 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));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Calculate any derived values from parameters */
82c4762a1bSJed Brown   N = user.mx * user.my;
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* Let Petsc determine the dimensions of the local vectors */
859371c9d4SSatish Balay   Nx = PETSC_DECIDE;
869371c9d4SSatish Balay   Ny = PETSC_DECIDE;
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /*
89c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
90c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
91c4762a1bSJed Brown      the distributed array, Create the vectors.
92c4762a1bSJed Brown   */
939566063dSJacob 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));
949566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
959566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
96c4762a1bSJed Brown   /*
97c4762a1bSJed Brown      Extract global and local vectors from DM; The local vectors are
98c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
99c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
100c4762a1bSJed Brown      the same types.
101c4762a1bSJed Brown   */
1029566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm, &x)); /* Solution */
1039566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(user.dm, &user.localX));
1049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.localX, &user.localV));
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &xl));
1079566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &xu));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* The TAO code begins here */
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /*
112c4762a1bSJed Brown      Create TAO solver and set desired solution method
113c4762a1bSJed Brown      The method must either be TAOTRON or TAOBLMVM
114c4762a1bSJed Brown      If TAOBLMVM is used, then hessian function is not called.
115c4762a1bSJed Brown   */
1169566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1179566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBLMVM));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Set initial solution guess; */
1209566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user));
1219566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user, x));
1229566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* Set routines for function, gradient and hessian evaluation */
1259566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &m));
1289566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(MPI_COMM_WORLD, m, m, N, N, 7, NULL, 3, NULL, &(user.H)));
1299566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
130c4762a1bSJed Brown 
1319566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(user.dm, &isltog));
1329566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(user.H, isltog, isltog));
1339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-matrixfree", &flg));
134c4762a1bSJed Brown   if (flg) {
1359566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PETSC_COMM_WORLD, m, m, N, N, (void *)&user, &H_shell));
1369566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(H_shell, MATOP_MULT, (void (*)(void))MyMatMult));
1379566063dSJacob Faibussowitsch     PetscCall(MatSetOption(H_shell, MAT_SYMMETRIC, PETSC_TRUE));
1389566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, H_shell, H_shell, MatrixFreeHessian, (void *)&user));
139c4762a1bSJed Brown   } else {
1409566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
141c4762a1bSJed Brown   }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* Set Variable bounds */
1449566063dSJacob Faibussowitsch   PetscCall(MSA_Plate(xl, xu, (void *)&user));
1459566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, xl, xu));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* Check for any tao command line options */
1489566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1519566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(TaoView(tao, PETSC_VIEWER_STDOUT_WORLD));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Free TAO data structures */
1569566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* Free PETSc data structures */
1599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
1619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
1629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
1639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localX));
1649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localV));
1659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Bottom));
1669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Top));
1679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Left));
1689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Right));
1699566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
170*48a46eb9SPierre Jolivet   if (flg) PetscCall(MatDestroy(&H_shell));
1719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
172b122ec5aSJacob Faibussowitsch   return 0;
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).
176c4762a1bSJed Brown 
177c4762a1bSJed Brown     Input Parameters:
178c4762a1bSJed Brown .   tao     - the Tao context
179c4762a1bSJed Brown .   X      - input vector
180a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
181c4762a1bSJed Brown 
182c4762a1bSJed Brown     Output Parameters:
183c4762a1bSJed Brown .   fcn     - the function value
184c4762a1bSJed Brown .   G      - vector containing the newly evaluated gradient
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    Notes:
187c4762a1bSJed Brown    In this case, we discretize the domain and Create triangles. The
188c4762a1bSJed Brown    surface of each triangle is planar, whose surface area can be easily
189c4762a1bSJed Brown    computed.  The total surface area is found by sweeping through the grid
190c4762a1bSJed Brown    and computing the surface area of the two triangles that have their
191c4762a1bSJed Brown    right angle at the grid point.  The diagonal line segments on the
192c4762a1bSJed Brown    grid that define the triangles run from top left to lower right.
193c4762a1bSJed Brown    The numbering of points starts at the lower left and runs left to
194c4762a1bSJed Brown    right, then bottom to top.
195c4762a1bSJed Brown */
1969371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx) {
197c4762a1bSJed Brown   AppCtx    *user = (AppCtx *)userCtx;
198c4762a1bSJed Brown   PetscInt   i, j, row;
199c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
200c4762a1bSJed Brown   PetscInt   xs, xm, gxs, gxm, ys, ym, gys, gym;
201c4762a1bSJed Brown   PetscReal  ft   = 0;
202c4762a1bSJed Brown   PetscReal  zero = 0.0;
203c4762a1bSJed Brown   PetscReal  hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy;
204c4762a1bSJed Brown   PetscReal  rhx = mx + 1, rhy = my + 1;
205c4762a1bSJed Brown   PetscReal  f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
206c4762a1bSJed Brown   PetscReal  df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
207c4762a1bSJed Brown   PetscReal *g, *x, *left, *right, *bottom, *top;
208c4762a1bSJed Brown   Vec        localX = user->localX, localG = user->localV;
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* Get local mesh boundaries */
2119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
2129566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* Scatter ghost points to local vector */
2159566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
2169566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* Initialize vector to zero */
2199566063dSJacob Faibussowitsch   PetscCall(VecSet(localG, zero));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Get pointers to vector data */
2229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &x));
2239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localG, &g));
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top, &top));
2259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom, &bottom));
2269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left, &left));
2279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right, &right));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
230c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
231c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
232c4762a1bSJed Brown       row = (j - gys) * gxm + (i - gxs);
233c4762a1bSJed Brown 
234c4762a1bSJed Brown       xc  = x[row];
235c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
236c4762a1bSJed Brown 
237c4762a1bSJed Brown       if (i == 0) { /* left side */
238c4762a1bSJed Brown         xl  = left[j - ys + 1];
239c4762a1bSJed Brown         xlt = left[j - ys + 2];
240c4762a1bSJed Brown       } else {
241c4762a1bSJed Brown         xl = x[row - 1];
242c4762a1bSJed Brown       }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown       if (j == 0) { /* bottom side */
245c4762a1bSJed Brown         xb  = bottom[i - xs + 1];
246c4762a1bSJed Brown         xrb = bottom[i - xs + 2];
247c4762a1bSJed Brown       } else {
248c4762a1bSJed Brown         xb = x[row - gxm];
249c4762a1bSJed Brown       }
250c4762a1bSJed Brown 
251c4762a1bSJed Brown       if (i + 1 == gxs + gxm) { /* right side */
252c4762a1bSJed Brown         xr  = right[j - ys + 1];
253c4762a1bSJed Brown         xrb = right[j - ys];
254c4762a1bSJed Brown       } else {
255c4762a1bSJed Brown         xr = x[row + 1];
256c4762a1bSJed Brown       }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown       if (j + 1 == gys + gym) { /* top side */
259c4762a1bSJed Brown         xt  = top[i - xs + 1];
260c4762a1bSJed Brown         xlt = top[i - xs];
261c4762a1bSJed Brown       } else {
262c4762a1bSJed Brown         xt = x[row + gxm];
263c4762a1bSJed Brown       }
264c4762a1bSJed Brown 
2659371c9d4SSatish Balay       if (i > gxs && j + 1 < gys + gym) { xlt = x[row - 1 + gxm]; }
2669371c9d4SSatish Balay       if (j > gys && i + 1 < gxs + gxm) { xrb = x[row + 1 - gxm]; }
267c4762a1bSJed Brown 
268c4762a1bSJed Brown       d1 = (xc - xl);
269c4762a1bSJed Brown       d2 = (xc - xr);
270c4762a1bSJed Brown       d3 = (xc - xt);
271c4762a1bSJed Brown       d4 = (xc - xb);
272c4762a1bSJed Brown       d5 = (xr - xrb);
273c4762a1bSJed Brown       d6 = (xrb - xb);
274c4762a1bSJed Brown       d7 = (xlt - xl);
275c4762a1bSJed Brown       d8 = (xt - xlt);
276c4762a1bSJed Brown 
277c4762a1bSJed Brown       df1dxc = d1 * hydhx;
278c4762a1bSJed Brown       df2dxc = (d1 * hydhx + d4 * hxdhy);
279c4762a1bSJed Brown       df3dxc = d3 * hxdhy;
280c4762a1bSJed Brown       df4dxc = (d2 * hydhx + d3 * hxdhy);
281c4762a1bSJed Brown       df5dxc = d2 * hydhx;
282c4762a1bSJed Brown       df6dxc = d4 * hxdhy;
283c4762a1bSJed Brown 
284c4762a1bSJed Brown       d1 *= rhx;
285c4762a1bSJed Brown       d2 *= rhx;
286c4762a1bSJed Brown       d3 *= rhy;
287c4762a1bSJed Brown       d4 *= rhy;
288c4762a1bSJed Brown       d5 *= rhy;
289c4762a1bSJed Brown       d6 *= rhx;
290c4762a1bSJed Brown       d7 *= rhy;
291c4762a1bSJed Brown       d8 *= rhx;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
294c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
295c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
296c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
297c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
298c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
299c4762a1bSJed Brown 
300c4762a1bSJed Brown       ft = ft + (f2 + f4);
301c4762a1bSJed Brown 
302c4762a1bSJed Brown       df1dxc /= f1;
303c4762a1bSJed Brown       df2dxc /= f2;
304c4762a1bSJed Brown       df3dxc /= f3;
305c4762a1bSJed Brown       df4dxc /= f4;
306c4762a1bSJed Brown       df5dxc /= f5;
307c4762a1bSJed Brown       df6dxc /= f6;
308c4762a1bSJed Brown 
309c4762a1bSJed Brown       g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5;
310c4762a1bSJed Brown     }
311c4762a1bSJed Brown   }
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
314c4762a1bSJed Brown   if (xs == 0) { /* left side */
315c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
316c4762a1bSJed Brown       d3 = (left[j - ys + 1] - left[j - ys + 2]) * rhy;
317c4762a1bSJed Brown       d2 = (left[j - ys + 1] - x[(j - gys) * gxm]) * rhx;
318c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
319c4762a1bSJed Brown     }
320c4762a1bSJed Brown   }
321c4762a1bSJed Brown   if (ys == 0) { /* bottom side */
322c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
323c4762a1bSJed Brown       d2 = (bottom[i + 1 - xs] - bottom[i - xs + 2]) * rhx;
324c4762a1bSJed Brown       d3 = (bottom[i - xs + 1] - x[i - gxs]) * rhy;
325c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
326c4762a1bSJed Brown     }
327c4762a1bSJed Brown   }
328c4762a1bSJed Brown 
329c4762a1bSJed Brown   if (xs + xm == mx) { /* right side */
330c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
331c4762a1bSJed Brown       d1 = (x[(j + 1 - gys) * gxm - 1] - right[j - ys + 1]) * rhx;
332c4762a1bSJed Brown       d4 = (right[j - ys] - right[j - ys + 1]) * rhy;
333c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
334c4762a1bSJed Brown     }
335c4762a1bSJed Brown   }
336c4762a1bSJed Brown   if (ys + ym == my) { /* top side */
337c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
338c4762a1bSJed Brown       d1 = (x[(gym - 1) * gxm + i - gxs] - top[i - xs + 1]) * rhy;
339c4762a1bSJed Brown       d4 = (top[i - xs + 1] - top[i - xs]) * rhx;
340c4762a1bSJed Brown       ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
341c4762a1bSJed Brown     }
342c4762a1bSJed Brown   }
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   if (ys == 0 && xs == 0) {
345c4762a1bSJed Brown     d1 = (left[0] - left[1]) * rhy;
346c4762a1bSJed Brown     d2 = (bottom[0] - bottom[1]) * rhx;
347c4762a1bSJed Brown     ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
348c4762a1bSJed Brown   }
349c4762a1bSJed Brown   if (ys + ym == my && xs + xm == mx) {
350c4762a1bSJed Brown     d1 = (right[ym + 1] - right[ym]) * rhy;
351c4762a1bSJed Brown     d2 = (top[xm + 1] - top[xm]) * rhx;
352c4762a1bSJed Brown     ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
353c4762a1bSJed Brown   }
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   ft = ft * area;
3569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   /* Restore vectors */
3599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &x));
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localG, &g));
3619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left, &left));
3629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top, &top));
3639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom, &bottom));
3649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right, &right));
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   /* Scatter values to global vector */
3679566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(user->dm, localG, INSERT_VALUES, G));
3689566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(user->dm, localG, INSERT_VALUES, G));
369c4762a1bSJed Brown 
3709566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(70.0 * xm * ym));
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   return 0;
373c4762a1bSJed Brown }
374c4762a1bSJed Brown 
375c4762a1bSJed Brown /* ------------------------------------------------------------------- */
376c4762a1bSJed Brown /*
377c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
378c4762a1bSJed Brown 
379c4762a1bSJed Brown    Input Parameters:
380c4762a1bSJed Brown .  tao  - the Tao context
381c4762a1bSJed Brown .  x    - input vector
382a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
383c4762a1bSJed Brown 
384c4762a1bSJed Brown    Output Parameters:
385c4762a1bSJed Brown .  A    - Hessian matrix
386c4762a1bSJed Brown .  B    - optionally different preconditioning matrix
387c4762a1bSJed Brown 
388c4762a1bSJed Brown    Notes:
389c4762a1bSJed Brown    Due to mesh point reordering with DMs, we must always work
390c4762a1bSJed Brown    with the local mesh points, and then transform them to the new
391c4762a1bSJed Brown    global numbering with the local-to-global mapping.  We cannot work
392c4762a1bSJed Brown    directly with the global numbers for the original uniprocessor mesh!
393c4762a1bSJed Brown 
394c4762a1bSJed Brown    Two methods are available for imposing this transformation
395c4762a1bSJed Brown    when setting matrix entries:
396c4762a1bSJed Brown      (A) MatSetValuesLocal(), using the local ordering (including
397c4762a1bSJed Brown          ghost points!)
398c4762a1bSJed Brown          - Do the following two steps once, before calling TaoSolve()
399c4762a1bSJed Brown            - Use DMGetISLocalToGlobalMapping() to extract the
400c4762a1bSJed Brown              local-to-global map from the DM
401c4762a1bSJed Brown            - Associate this map with the matrix by calling
402c4762a1bSJed Brown              MatSetLocalToGlobalMapping()
403c4762a1bSJed Brown          - Then set matrix entries using the local ordering
404c4762a1bSJed Brown            by calling MatSetValuesLocal()
405c4762a1bSJed Brown      (B) MatSetValues(), using the global ordering
406c4762a1bSJed Brown          - Use DMGetGlobalIndices() to extract the local-to-global map
407c4762a1bSJed Brown          - Then apply this map explicitly yourself
408c4762a1bSJed Brown          - Set matrix entries using the global ordering by calling
409c4762a1bSJed Brown            MatSetValues()
410c4762a1bSJed Brown    Option (A) seems cleaner/easier in many cases, and is the procedure
411c4762a1bSJed Brown    used in this example.
412c4762a1bSJed Brown */
4139371c9d4SSatish Balay PetscErrorCode FormHessian(Tao tao, Vec X, Mat Hptr, Mat Hessian, void *ptr) {
414c4762a1bSJed Brown   AppCtx    *user = (AppCtx *)ptr;
415c4762a1bSJed Brown   PetscInt   i, j, k, row;
416c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
417c4762a1bSJed Brown   PetscInt   xs, xm, gxs, gxm, ys, ym, gys, gym, col[7];
418c4762a1bSJed Brown   PetscReal  hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
419c4762a1bSJed Brown   PetscReal  rhx = mx + 1, rhy = my + 1;
420c4762a1bSJed Brown   PetscReal  f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
421c4762a1bSJed Brown   PetscReal  hl, hr, ht, hb, hc, htl, hbr;
422c4762a1bSJed Brown   PetscReal *x, *left, *right, *bottom, *top;
423c4762a1bSJed Brown   PetscReal  v[7];
424c4762a1bSJed Brown   Vec        localX = user->localX;
425c4762a1bSJed Brown   PetscBool  assembled;
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   /* Set various matrix options */
4289566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
429c4762a1bSJed Brown 
430c4762a1bSJed Brown   /* Initialize matrix entries to zero */
4319566063dSJacob Faibussowitsch   PetscCall(MatAssembled(Hessian, &assembled));
4329566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(Hessian));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /* Get local mesh boundaries */
4359566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
4369566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   /* Scatter ghost points to local vector */
4399566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
4409566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   /* Get pointers to vector data */
4439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &x));
4449566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top, &top));
4459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom, &bottom));
4469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left, &left));
4479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right, &right));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
452c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
453c4762a1bSJed Brown       row = (j - gys) * gxm + (i - gxs);
454c4762a1bSJed Brown 
455c4762a1bSJed Brown       xc  = x[row];
456c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
457c4762a1bSJed Brown 
458c4762a1bSJed Brown       /* Left side */
459c4762a1bSJed Brown       if (i == gxs) {
460c4762a1bSJed Brown         xl  = left[j - ys + 1];
461c4762a1bSJed Brown         xlt = left[j - ys + 2];
462c4762a1bSJed Brown       } else {
463c4762a1bSJed Brown         xl = x[row - 1];
464c4762a1bSJed Brown       }
465c4762a1bSJed Brown 
466c4762a1bSJed Brown       if (j == gys) {
467c4762a1bSJed Brown         xb  = bottom[i - xs + 1];
468c4762a1bSJed Brown         xrb = bottom[i - xs + 2];
469c4762a1bSJed Brown       } else {
470c4762a1bSJed Brown         xb = x[row - gxm];
471c4762a1bSJed Brown       }
472c4762a1bSJed Brown 
473c4762a1bSJed Brown       if (i + 1 == gxs + gxm) {
474c4762a1bSJed Brown         xr  = right[j - ys + 1];
475c4762a1bSJed Brown         xrb = right[j - ys];
476c4762a1bSJed Brown       } else {
477c4762a1bSJed Brown         xr = x[row + 1];
478c4762a1bSJed Brown       }
479c4762a1bSJed Brown 
480c4762a1bSJed Brown       if (j + 1 == gys + gym) {
481c4762a1bSJed Brown         xt  = top[i - xs + 1];
482c4762a1bSJed Brown         xlt = top[i - xs];
483c4762a1bSJed Brown       } else {
484c4762a1bSJed Brown         xt = x[row + gxm];
485c4762a1bSJed Brown       }
486c4762a1bSJed Brown 
4879371c9d4SSatish Balay       if (i > gxs && j + 1 < gys + gym) { xlt = x[row - 1 + gxm]; }
4889371c9d4SSatish Balay       if (j > gys && i + 1 < gxs + gxm) { xrb = x[row + 1 - gxm]; }
489c4762a1bSJed Brown 
490c4762a1bSJed Brown       d1 = (xc - xl) * rhx;
491c4762a1bSJed Brown       d2 = (xc - xr) * rhx;
492c4762a1bSJed Brown       d3 = (xc - xt) * rhy;
493c4762a1bSJed Brown       d4 = (xc - xb) * rhy;
494c4762a1bSJed Brown       d5 = (xrb - xr) * rhy;
495c4762a1bSJed Brown       d6 = (xrb - xb) * rhx;
496c4762a1bSJed Brown       d7 = (xlt - xl) * rhy;
497c4762a1bSJed Brown       d8 = (xlt - xt) * rhx;
498c4762a1bSJed Brown 
499c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
500c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
501c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
502c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
503c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
504c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
505c4762a1bSJed Brown 
5069371c9d4SSatish Balay       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
5079371c9d4SSatish Balay       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
5089371c9d4SSatish Balay       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
5099371c9d4SSatish Balay       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
510c4762a1bSJed Brown 
511c4762a1bSJed Brown       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
512c4762a1bSJed Brown       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
513c4762a1bSJed Brown 
5149371c9d4SSatish 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);
515c4762a1bSJed Brown 
5169371c9d4SSatish Balay       hl *= 0.5;
5179371c9d4SSatish Balay       hr *= 0.5;
5189371c9d4SSatish Balay       ht *= 0.5;
5199371c9d4SSatish Balay       hb *= 0.5;
5209371c9d4SSatish Balay       hbr *= 0.5;
5219371c9d4SSatish Balay       htl *= 0.5;
5229371c9d4SSatish Balay       hc *= 0.5;
523c4762a1bSJed Brown 
524c4762a1bSJed Brown       k = 0;
525c4762a1bSJed Brown       if (j > 0) {
5269371c9d4SSatish Balay         v[k]   = hb;
5279371c9d4SSatish Balay         col[k] = row - gxm;
5289371c9d4SSatish Balay         k++;
529c4762a1bSJed Brown       }
530c4762a1bSJed Brown 
531c4762a1bSJed Brown       if (j > 0 && i < mx - 1) {
5329371c9d4SSatish Balay         v[k]   = hbr;
5339371c9d4SSatish Balay         col[k] = row - gxm + 1;
5349371c9d4SSatish Balay         k++;
535c4762a1bSJed Brown       }
536c4762a1bSJed Brown 
537c4762a1bSJed Brown       if (i > 0) {
5389371c9d4SSatish Balay         v[k]   = hl;
5399371c9d4SSatish Balay         col[k] = row - 1;
5409371c9d4SSatish Balay         k++;
541c4762a1bSJed Brown       }
542c4762a1bSJed Brown 
5439371c9d4SSatish Balay       v[k]   = hc;
5449371c9d4SSatish Balay       col[k] = row;
5459371c9d4SSatish Balay       k++;
546c4762a1bSJed Brown 
547c4762a1bSJed Brown       if (i < mx - 1) {
5489371c9d4SSatish Balay         v[k]   = hr;
5499371c9d4SSatish Balay         col[k] = row + 1;
5509371c9d4SSatish Balay         k++;
551c4762a1bSJed Brown       }
552c4762a1bSJed Brown 
553c4762a1bSJed Brown       if (i > 0 && j < my - 1) {
5549371c9d4SSatish Balay         v[k]   = htl;
5559371c9d4SSatish Balay         col[k] = row + gxm - 1;
5569371c9d4SSatish Balay         k++;
557c4762a1bSJed Brown       }
558c4762a1bSJed Brown 
559c4762a1bSJed Brown       if (j < my - 1) {
5609371c9d4SSatish Balay         v[k]   = ht;
5619371c9d4SSatish Balay         col[k] = row + gxm;
5629371c9d4SSatish Balay         k++;
563c4762a1bSJed Brown       }
564c4762a1bSJed Brown 
565c4762a1bSJed Brown       /*
566c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
567c4762a1bSJed Brown          earlier, in the main routine.
568c4762a1bSJed Brown       */
5699566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Hessian, 1, &row, k, col, v, INSERT_VALUES));
570c4762a1bSJed Brown     }
571c4762a1bSJed Brown   }
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   /* Restore vectors */
5749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &x));
5759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left, &left));
5769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top, &top));
5779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom, &bottom));
5789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right, &right));
579c4762a1bSJed Brown 
580c4762a1bSJed Brown   /* Assemble the matrix */
5819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY));
5829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY));
583c4762a1bSJed Brown 
5849566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0 * xm * ym));
585c4762a1bSJed Brown   return 0;
586c4762a1bSJed Brown }
587c4762a1bSJed Brown 
588c4762a1bSJed Brown /* ------------------------------------------------------------------- */
589c4762a1bSJed Brown /*
590c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
591c4762a1bSJed Brown    the region.
592c4762a1bSJed Brown 
593c4762a1bSJed Brown    Input Parameter:
594c4762a1bSJed Brown .  user - user-defined application context
595c4762a1bSJed Brown 
596c4762a1bSJed Brown    Output Parameter:
597c4762a1bSJed Brown .  user - user-defined application context
598c4762a1bSJed Brown */
5999371c9d4SSatish Balay static PetscErrorCode MSA_BoundaryConditions(AppCtx *user) {
600c4762a1bSJed Brown   PetscInt   i, j, k, maxits = 5, limit = 0;
601c4762a1bSJed Brown   PetscInt   xs, ys, xm, ym, gxs, gys, gxm, gym;
602c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
603c4762a1bSJed Brown   PetscInt   bsize = 0, lsize = 0, tsize = 0, rsize = 0;
604c4762a1bSJed Brown   PetscReal  one = 1.0, two = 2.0, three = 3.0, scl = 1.0, tol = 1e-10;
605c4762a1bSJed Brown   PetscReal  fnorm, det, hx, hy, xt = 0, yt = 0;
606c4762a1bSJed Brown   PetscReal  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
607c4762a1bSJed Brown   PetscReal  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
608c4762a1bSJed Brown   PetscReal *boundary;
609c4762a1bSJed Brown   PetscBool  flg;
610c4762a1bSJed Brown   Vec        Bottom, Top, Right, Left;
611c4762a1bSJed Brown 
612c4762a1bSJed Brown   /* Get local mesh boundaries */
6139566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
6149566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
615c4762a1bSJed Brown 
616c4762a1bSJed Brown   bsize = xm + 2;
617c4762a1bSJed Brown   lsize = ym + 2;
618c4762a1bSJed Brown   rsize = ym + 2;
619c4762a1bSJed Brown   tsize = xm + 2;
620c4762a1bSJed Brown 
6219566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, bsize, PETSC_DECIDE, &Bottom));
6229566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, tsize, PETSC_DECIDE, &Top));
6239566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, lsize, PETSC_DECIDE, &Left));
6249566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD, rsize, PETSC_DECIDE, &Right));
625c4762a1bSJed Brown 
626c4762a1bSJed Brown   user->Top    = Top;
627c4762a1bSJed Brown   user->Left   = Left;
628c4762a1bSJed Brown   user->Bottom = Bottom;
629c4762a1bSJed Brown   user->Right  = Right;
630c4762a1bSJed Brown 
6319371c9d4SSatish Balay   hx = (r - l) / (mx + 1);
6329371c9d4SSatish Balay   hy = (t - b) / (my + 1);
633c4762a1bSJed Brown 
634c4762a1bSJed Brown   for (j = 0; j < 4; j++) {
635c4762a1bSJed Brown     if (j == 0) {
636c4762a1bSJed Brown       yt    = b;
637c4762a1bSJed Brown       xt    = l + hx * xs;
638c4762a1bSJed Brown       limit = bsize;
639c4762a1bSJed Brown       VecGetArray(Bottom, &boundary);
640c4762a1bSJed Brown     } else if (j == 1) {
641c4762a1bSJed Brown       yt    = t;
642c4762a1bSJed Brown       xt    = l + hx * xs;
643c4762a1bSJed Brown       limit = tsize;
644c4762a1bSJed Brown       VecGetArray(Top, &boundary);
645c4762a1bSJed Brown     } else if (j == 2) {
646c4762a1bSJed Brown       yt    = b + hy * ys;
647c4762a1bSJed Brown       xt    = l;
648c4762a1bSJed Brown       limit = lsize;
649c4762a1bSJed Brown       VecGetArray(Left, &boundary);
650c4762a1bSJed Brown     } else if (j == 3) {
651c4762a1bSJed Brown       yt    = b + hy * ys;
652c4762a1bSJed Brown       xt    = r;
653c4762a1bSJed Brown       limit = rsize;
654c4762a1bSJed Brown       VecGetArray(Right, &boundary);
655c4762a1bSJed Brown     }
656c4762a1bSJed Brown 
657c4762a1bSJed Brown     for (i = 0; i < limit; i++) {
658c4762a1bSJed Brown       u1 = xt;
659c4762a1bSJed Brown       u2 = -yt;
660c4762a1bSJed Brown       for (k = 0; k < maxits; k++) {
661c4762a1bSJed Brown         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
662c4762a1bSJed Brown         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
663c4762a1bSJed Brown         fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2);
664c4762a1bSJed Brown         if (fnorm <= tol) break;
665c4762a1bSJed Brown         njac11 = one + u2 * u2 - u1 * u1;
666c4762a1bSJed Brown         njac12 = two * u1 * u2;
667c4762a1bSJed Brown         njac21 = -two * u1 * u2;
668c4762a1bSJed Brown         njac22 = -one - u1 * u1 + u2 * u2;
669c4762a1bSJed Brown         det    = njac11 * njac22 - njac21 * njac12;
670c4762a1bSJed Brown         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
671c4762a1bSJed Brown         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
672c4762a1bSJed Brown       }
673c4762a1bSJed Brown 
674c4762a1bSJed Brown       boundary[i] = u1 * u1 - u2 * u2;
675c4762a1bSJed Brown       if (j == 0 || j == 1) {
676c4762a1bSJed Brown         xt = xt + hx;
677c4762a1bSJed Brown       } else if (j == 2 || j == 3) {
678c4762a1bSJed Brown         yt = yt + hy;
679c4762a1bSJed Brown       }
680c4762a1bSJed Brown     }
681c4762a1bSJed Brown     if (j == 0) {
6829566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Bottom, &boundary));
683c4762a1bSJed Brown     } else if (j == 1) {
6849566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Top, &boundary));
685c4762a1bSJed Brown     } else if (j == 2) {
6869566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Left, &boundary));
687c4762a1bSJed Brown     } else if (j == 3) {
6889566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Right, &boundary));
689c4762a1bSJed Brown     }
690c4762a1bSJed Brown   }
691c4762a1bSJed Brown 
692c4762a1bSJed Brown   /* Scale the boundary if desired */
693c4762a1bSJed Brown 
6949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg));
6951baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Bottom, scl));
6969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg));
6971baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Top, scl));
6989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg));
6991baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Right, scl));
700c4762a1bSJed Brown 
7019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg));
7021baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Left, scl));
703c4762a1bSJed Brown   return 0;
704c4762a1bSJed Brown }
705c4762a1bSJed Brown 
706c4762a1bSJed Brown /* ------------------------------------------------------------------- */
707c4762a1bSJed Brown /*
708c4762a1bSJed Brown    MSA_Plate -  Calculates an obstacle for surface to stretch over.
709c4762a1bSJed Brown 
710c4762a1bSJed Brown    Input Parameter:
711c4762a1bSJed Brown .  user - user-defined application context
712c4762a1bSJed Brown 
713c4762a1bSJed Brown    Output Parameter:
714c4762a1bSJed Brown .  user - user-defined application context
715c4762a1bSJed Brown */
7169371c9d4SSatish Balay static PetscErrorCode MSA_Plate(Vec XL, Vec XU, void *ctx) {
717c4762a1bSJed Brown   AppCtx    *user = (AppCtx *)ctx;
718c4762a1bSJed Brown   PetscInt   i, j, row;
719c4762a1bSJed Brown   PetscInt   xs, ys, xm, ym;
720c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my, bmy, bmx;
721c4762a1bSJed Brown   PetscReal  t1, t2, t3;
722c4762a1bSJed Brown   PetscReal *xl, lb = PETSC_NINFINITY, ub = PETSC_INFINITY;
723c4762a1bSJed Brown   PetscBool  cylinder;
724c4762a1bSJed Brown 
7259371c9d4SSatish Balay   user->bmy = PetscMax(0, user->bmy);
7269371c9d4SSatish Balay   user->bmy = PetscMin(my, user->bmy);
7279371c9d4SSatish Balay   user->bmx = PetscMax(0, user->bmx);
7289371c9d4SSatish Balay   user->bmx = PetscMin(mx, user->bmx);
7299371c9d4SSatish Balay   bmy       = user->bmy;
7309371c9d4SSatish Balay   bmx       = user->bmx;
731c4762a1bSJed Brown 
7329566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
733c4762a1bSJed Brown 
7349566063dSJacob Faibussowitsch   PetscCall(VecSet(XL, lb));
7359566063dSJacob Faibussowitsch   PetscCall(VecSet(XU, ub));
736c4762a1bSJed Brown 
7379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(XL, &xl));
738c4762a1bSJed Brown 
7399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-cylinder", &cylinder));
740c4762a1bSJed Brown   /* Compute the optional lower box */
741c4762a1bSJed Brown   if (cylinder) {
742c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
743c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
744c4762a1bSJed Brown         row = (j - ys) * xm + (i - xs);
745c4762a1bSJed Brown         t1  = (2.0 * i - mx) * bmy;
746c4762a1bSJed Brown         t2  = (2.0 * j - my) * bmx;
747c4762a1bSJed Brown         t3  = bmx * bmx * bmy * bmy;
7489371c9d4SSatish Balay         if (t1 * t1 + t2 * t2 <= t3) { xl[row] = user->bheight; }
749c4762a1bSJed Brown       }
750c4762a1bSJed Brown     }
751c4762a1bSJed Brown   } else {
752c4762a1bSJed Brown     /* Compute the optional lower box */
753c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
754c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
755c4762a1bSJed Brown         row = (j - ys) * xm + (i - xs);
7569371c9d4SSatish Balay         if (i >= (mx - bmx) / 2 && i < mx - (mx - bmx) / 2 && j >= (my - bmy) / 2 && j < my - (my - bmy) / 2) { xl[row] = user->bheight; }
757c4762a1bSJed Brown       }
758c4762a1bSJed Brown     }
759c4762a1bSJed Brown   }
7609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(XL, &xl));
761c4762a1bSJed Brown 
762c4762a1bSJed Brown   return 0;
763c4762a1bSJed Brown }
764c4762a1bSJed Brown 
765c4762a1bSJed Brown /* ------------------------------------------------------------------- */
766c4762a1bSJed Brown /*
767c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
768c4762a1bSJed Brown 
769c4762a1bSJed Brown    Input Parameters:
770c4762a1bSJed Brown .  user - user-defined application context
771c4762a1bSJed Brown .  X - vector for initial guess
772c4762a1bSJed Brown 
773c4762a1bSJed Brown    Output Parameters:
774c4762a1bSJed Brown .  X - newly computed initial guess
775c4762a1bSJed Brown */
7769371c9d4SSatish Balay static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X) {
777c4762a1bSJed Brown   PetscInt  start = -1, i, j;
778c4762a1bSJed Brown   PetscReal zero  = 0.0;
779c4762a1bSJed Brown   PetscBool flg;
780c4762a1bSJed Brown 
7819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg));
782c4762a1bSJed Brown   if (flg && start == 0) { /* The zero vector is reasonable */
7839566063dSJacob Faibussowitsch     PetscCall(VecSet(X, zero));
784c4762a1bSJed Brown   } else if (flg && start > 0) { /* Try a random start between -0.5 and 0.5 */
7859371c9d4SSatish Balay     PetscRandom rctx;
7869371c9d4SSatish Balay     PetscReal   np5 = -0.5;
787c4762a1bSJed Brown 
7889566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(MPI_COMM_WORLD, &rctx));
789*48a46eb9SPierre Jolivet     for (i = 0; i < start; i++) PetscCall(VecSetRandom(X, rctx));
7909566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rctx));
7919566063dSJacob Faibussowitsch     PetscCall(VecShift(X, np5));
792c4762a1bSJed Brown 
793c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
794c4762a1bSJed Brown 
795c4762a1bSJed Brown     PetscInt   row, xs, xm, gxs, gxm, ys, ym, gys, gym;
796c4762a1bSJed Brown     PetscInt   mx = user->mx, my = user->my;
797c4762a1bSJed Brown     PetscReal *x, *left, *right, *bottom, *top;
798c4762a1bSJed Brown     Vec        localX = user->localX;
799c4762a1bSJed Brown 
800c4762a1bSJed Brown     /* Get local mesh boundaries */
8019566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
8029566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
803c4762a1bSJed Brown 
804c4762a1bSJed Brown     /* Get pointers to vector data */
8059566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Top, &top));
8069566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Bottom, &bottom));
8079566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Left, &left));
8089566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Right, &right));
809c4762a1bSJed Brown 
8109566063dSJacob Faibussowitsch     PetscCall(VecGetArray(localX, &x));
811c4762a1bSJed Brown     /* Perform local computations */
812c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
813c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
814c4762a1bSJed Brown         row    = (j - gys) * gxm + (i - gxs);
8152da392ccSBarry 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;
816c4762a1bSJed Brown       }
817c4762a1bSJed Brown     }
818c4762a1bSJed Brown 
819c4762a1bSJed Brown     /* Restore vectors */
8209566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(localX, &x));
821c4762a1bSJed Brown 
8229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Left, &left));
8239566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Top, &top));
8249566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Bottom, &bottom));
8259566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Right, &right));
826c4762a1bSJed Brown 
827c4762a1bSJed Brown     /* Scatter values into global vector */
8289566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(user->dm, localX, INSERT_VALUES, X));
8299566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(user->dm, localX, INSERT_VALUES, X));
830c4762a1bSJed Brown   }
831c4762a1bSJed Brown   return 0;
832c4762a1bSJed Brown }
833c4762a1bSJed Brown 
834c4762a1bSJed Brown /* For testing matrix free submatrices */
8359371c9d4SSatish Balay PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) {
836c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
837c4762a1bSJed Brown   PetscFunctionBegin;
8389566063dSJacob Faibussowitsch   PetscCall(FormHessian(tao, x, user->H, user->H, ptr));
839c4762a1bSJed Brown   PetscFunctionReturn(0);
840c4762a1bSJed Brown }
8419371c9d4SSatish Balay PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y) {
842c4762a1bSJed Brown   void   *ptr;
843c4762a1bSJed Brown   AppCtx *user;
844c4762a1bSJed Brown   PetscFunctionBegin;
8459566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(H_shell, &ptr));
846c4762a1bSJed Brown   user = (AppCtx *)ptr;
8479566063dSJacob Faibussowitsch   PetscCall(MatMult(user->H, X, Y));
848c4762a1bSJed Brown   PetscFunctionReturn(0);
849c4762a1bSJed Brown }
850c4762a1bSJed Brown 
851c4762a1bSJed Brown /*TEST
852c4762a1bSJed Brown 
853c4762a1bSJed Brown    build:
854c4762a1bSJed Brown       requires: !complex
855c4762a1bSJed Brown 
856c4762a1bSJed Brown    test:
857c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
858c4762a1bSJed Brown       requires: !single
859c4762a1bSJed Brown 
860c4762a1bSJed Brown    test:
861c4762a1bSJed Brown       suffix: 2
862c4762a1bSJed Brown       nsize: 2
863c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
864c4762a1bSJed Brown       requires: !single
865c4762a1bSJed Brown 
866c4762a1bSJed Brown    test:
867c4762a1bSJed Brown       suffix: 3
868c4762a1bSJed Brown       nsize: 3
869c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
870c4762a1bSJed Brown       requires: !single
871c4762a1bSJed Brown 
872c4762a1bSJed Brown    test:
873c4762a1bSJed Brown       suffix: 4
874c4762a1bSJed Brown       nsize: 3
875c4762a1bSJed 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
876c4762a1bSJed Brown       requires: !single
877c4762a1bSJed Brown 
878c4762a1bSJed Brown    test:
879c4762a1bSJed Brown       suffix: 5
880c4762a1bSJed Brown       nsize: 3
881c4762a1bSJed 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
882c4762a1bSJed Brown       requires: !single
883c4762a1bSJed Brown 
884c4762a1bSJed Brown    test:
885c4762a1bSJed Brown       suffix: 6
886c4762a1bSJed Brown       nsize: 3
887c4762a1bSJed 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
888c4762a1bSJed Brown       requires: !single
889c4762a1bSJed Brown 
890c4762a1bSJed Brown    test:
891c4762a1bSJed Brown       suffix: 7
892c4762a1bSJed Brown       nsize: 3
893c4762a1bSJed 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
894c4762a1bSJed Brown       requires: !single
895c4762a1bSJed Brown 
896c4762a1bSJed Brown    test:
897c4762a1bSJed Brown       suffix: 8
898c4762a1bSJed 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
899c4762a1bSJed Brown       requires: !single
900c4762a1bSJed Brown 
901c4762a1bSJed Brown    test:
902c4762a1bSJed Brown       suffix: 9
903c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
904c4762a1bSJed Brown       requires: !single
905c4762a1bSJed Brown 
906c4762a1bSJed Brown    test:
907c4762a1bSJed Brown       suffix: 10
908c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
909c4762a1bSJed Brown       requires: !single
910c4762a1bSJed Brown 
911c4762a1bSJed Brown    test:
912c4762a1bSJed Brown       suffix: 11
913c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
914c4762a1bSJed Brown       requires: !single
915c4762a1bSJed Brown 
916c4762a1bSJed Brown    test:
917c4762a1bSJed Brown       suffix: 12
918c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
919c4762a1bSJed Brown       requires: !single
920c4762a1bSJed Brown 
921c4762a1bSJed Brown    test:
922c4762a1bSJed Brown       suffix: 13
923c4762a1bSJed 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
924c4762a1bSJed Brown       requires: !single
925c4762a1bSJed Brown 
926c4762a1bSJed Brown    test:
927c4762a1bSJed Brown       suffix: 14
928c4762a1bSJed 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
929c4762a1bSJed Brown       requires: !single
930c4762a1bSJed Brown 
931c4762a1bSJed Brown    test:
932c4762a1bSJed Brown       suffix: 15
933c4762a1bSJed 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
934c4762a1bSJed Brown       requires: !single
935c4762a1bSJed Brown 
936c4762a1bSJed Brown    test:
937c4762a1bSJed Brown      suffix: 16
938c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
939c4762a1bSJed Brown      requires: !single
940c4762a1bSJed Brown 
941c4762a1bSJed Brown    test:
942c4762a1bSJed Brown      suffix: 17
943c4762a1bSJed 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
944c4762a1bSJed Brown      requires: !single
945c4762a1bSJed Brown 
94634ad9904SAlp Dener    test:
94734ad9904SAlp Dener      suffix: 18
94834ad9904SAlp 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
94934ad9904SAlp Dener      requires: !single
95034ad9904SAlp Dener 
95134ad9904SAlp Dener    test:
95234ad9904SAlp Dener      suffix: 19
95334ad9904SAlp 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
95434ad9904SAlp Dener      requires: !single
95534ad9904SAlp Dener 
95634ad9904SAlp Dener    test:
95734ad9904SAlp Dener      suffix: 20
95834ad9904SAlp 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
95934ad9904SAlp Dener 
960c4762a1bSJed Brown TEST*/
961