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
4501c1178eSBarry Smith /* For testing matrix-free submatrices */
46c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *);
47c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat, Vec, Vec);
48c4762a1bSJed Brown
main(int argc,char ** argv)49d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
50d71ae5a4SJacob 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;
62c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, 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
85f0b74427SPierre Jolivet /* 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));
12977433607SBarry Smith PetscCall(DMCreateMatrix(user.dm, &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));
13757d50842SBarry Smith PetscCall(MatShellSetOperation(H_shell, MATOP_MULT, (PetscErrorCodeFn *)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 */
FormFunctionGradient(Tao tao,Vec X,PetscReal * fcn,Vec G,void * userCtx)197d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx)
198d71ae5a4SJacob 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
2123ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
213c4762a1bSJed Brown /* Get local mesh boundaries */
2149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
2159566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
216c4762a1bSJed Brown
217c4762a1bSJed Brown /* Scatter ghost points to local vector */
2189566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
2199566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
220c4762a1bSJed Brown
221c4762a1bSJed Brown /* Initialize vector to zero */
2229566063dSJacob Faibussowitsch PetscCall(VecSet(localG, zero));
223c4762a1bSJed Brown
224c4762a1bSJed Brown /* Get pointers to vector data */
2259566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x));
2269566063dSJacob Faibussowitsch PetscCall(VecGetArray(localG, &g));
2279566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Top, &top));
2289566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Bottom, &bottom));
2299566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Left, &left));
2309566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Right, &right));
231c4762a1bSJed Brown
232c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */
233c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
234c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
235c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs);
236c4762a1bSJed Brown
237c4762a1bSJed Brown xc = x[row];
238c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc;
239c4762a1bSJed Brown
240c4762a1bSJed Brown if (i == 0) { /* left side */
241c4762a1bSJed Brown xl = left[j - ys + 1];
242c4762a1bSJed Brown xlt = left[j - ys + 2];
243c4762a1bSJed Brown } else {
244c4762a1bSJed Brown xl = x[row - 1];
245c4762a1bSJed Brown }
246c4762a1bSJed Brown
247c4762a1bSJed Brown if (j == 0) { /* bottom side */
248c4762a1bSJed Brown xb = bottom[i - xs + 1];
249c4762a1bSJed Brown xrb = bottom[i - xs + 2];
250c4762a1bSJed Brown } else {
251c4762a1bSJed Brown xb = x[row - gxm];
252c4762a1bSJed Brown }
253c4762a1bSJed Brown
254c4762a1bSJed Brown if (i + 1 == gxs + gxm) { /* right side */
255c4762a1bSJed Brown xr = right[j - ys + 1];
256c4762a1bSJed Brown xrb = right[j - ys];
257c4762a1bSJed Brown } else {
258c4762a1bSJed Brown xr = x[row + 1];
259c4762a1bSJed Brown }
260c4762a1bSJed Brown
261c4762a1bSJed Brown if (j + 1 == gys + gym) { /* top side */
262c4762a1bSJed Brown xt = top[i - xs + 1];
263c4762a1bSJed Brown xlt = top[i - xs];
264c4762a1bSJed Brown } else {
265c4762a1bSJed Brown xt = x[row + gxm];
266c4762a1bSJed Brown }
267c4762a1bSJed Brown
268ad540459SPierre Jolivet if (i > gxs && j + 1 < gys + gym) xlt = x[row - 1 + gxm];
269ad540459SPierre Jolivet if (j > gys && i + 1 < gxs + gxm) xrb = x[row + 1 - gxm];
270c4762a1bSJed Brown
271c4762a1bSJed Brown d1 = (xc - xl);
272c4762a1bSJed Brown d2 = (xc - xr);
273c4762a1bSJed Brown d3 = (xc - xt);
274c4762a1bSJed Brown d4 = (xc - xb);
275c4762a1bSJed Brown d5 = (xr - xrb);
276c4762a1bSJed Brown d6 = (xrb - xb);
277c4762a1bSJed Brown d7 = (xlt - xl);
278c4762a1bSJed Brown d8 = (xt - xlt);
279c4762a1bSJed Brown
280c4762a1bSJed Brown df1dxc = d1 * hydhx;
281c4762a1bSJed Brown df2dxc = (d1 * hydhx + d4 * hxdhy);
282c4762a1bSJed Brown df3dxc = d3 * hxdhy;
283c4762a1bSJed Brown df4dxc = (d2 * hydhx + d3 * hxdhy);
284c4762a1bSJed Brown df5dxc = d2 * hydhx;
285c4762a1bSJed Brown df6dxc = d4 * hxdhy;
286c4762a1bSJed Brown
287c4762a1bSJed Brown d1 *= rhx;
288c4762a1bSJed Brown d2 *= rhx;
289c4762a1bSJed Brown d3 *= rhy;
290c4762a1bSJed Brown d4 *= rhy;
291c4762a1bSJed Brown d5 *= rhy;
292c4762a1bSJed Brown d6 *= rhx;
293c4762a1bSJed Brown d7 *= rhy;
294c4762a1bSJed Brown d8 *= rhx;
295c4762a1bSJed Brown
296c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
297c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
298c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
299c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
300c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
301c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
302c4762a1bSJed Brown
303c4762a1bSJed Brown ft = ft + (f2 + f4);
304c4762a1bSJed Brown
305c4762a1bSJed Brown df1dxc /= f1;
306c4762a1bSJed Brown df2dxc /= f2;
307c4762a1bSJed Brown df3dxc /= f3;
308c4762a1bSJed Brown df4dxc /= f4;
309c4762a1bSJed Brown df5dxc /= f5;
310c4762a1bSJed Brown df6dxc /= f6;
311c4762a1bSJed Brown
312c4762a1bSJed Brown g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5;
313c4762a1bSJed Brown }
314c4762a1bSJed Brown }
315c4762a1bSJed Brown
316c4762a1bSJed Brown /* Compute triangular areas along the border of the domain. */
317c4762a1bSJed Brown if (xs == 0) { /* left side */
318c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
319c4762a1bSJed Brown d3 = (left[j - ys + 1] - left[j - ys + 2]) * rhy;
320c4762a1bSJed Brown d2 = (left[j - ys + 1] - x[(j - gys) * gxm]) * rhx;
321c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
322c4762a1bSJed Brown }
323c4762a1bSJed Brown }
324c4762a1bSJed Brown if (ys == 0) { /* bottom side */
325c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
326c4762a1bSJed Brown d2 = (bottom[i + 1 - xs] - bottom[i - xs + 2]) * rhx;
327c4762a1bSJed Brown d3 = (bottom[i - xs + 1] - x[i - gxs]) * rhy;
328c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
329c4762a1bSJed Brown }
330c4762a1bSJed Brown }
331c4762a1bSJed Brown
332c4762a1bSJed Brown if (xs + xm == mx) { /* right side */
333c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
334c4762a1bSJed Brown d1 = (x[(j + 1 - gys) * gxm - 1] - right[j - ys + 1]) * rhx;
335c4762a1bSJed Brown d4 = (right[j - ys] - right[j - ys + 1]) * rhy;
336c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown }
339c4762a1bSJed Brown if (ys + ym == my) { /* top side */
340c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
341c4762a1bSJed Brown d1 = (x[(gym - 1) * gxm + i - gxs] - top[i - xs + 1]) * rhy;
342c4762a1bSJed Brown d4 = (top[i - xs + 1] - top[i - xs]) * rhx;
343c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
344c4762a1bSJed Brown }
345c4762a1bSJed Brown }
346c4762a1bSJed Brown
347c4762a1bSJed Brown if (ys == 0 && xs == 0) {
348c4762a1bSJed Brown d1 = (left[0] - left[1]) * rhy;
349c4762a1bSJed Brown d2 = (bottom[0] - bottom[1]) * rhx;
350c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown if (ys + ym == my && xs + xm == mx) {
353c4762a1bSJed Brown d1 = (right[ym + 1] - right[ym]) * rhy;
354c4762a1bSJed Brown d2 = (top[xm + 1] - top[xm]) * rhx;
355c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
356c4762a1bSJed Brown }
357c4762a1bSJed Brown
358c4762a1bSJed Brown ft = ft * area;
359462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
360c4762a1bSJed Brown
361c4762a1bSJed Brown /* Restore vectors */
3629566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x));
3639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localG, &g));
3649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Left, &left));
3659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Top, &top));
3669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Bottom, &bottom));
3679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Right, &right));
368c4762a1bSJed Brown
369c4762a1bSJed Brown /* Scatter values to global vector */
3709566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(user->dm, localG, INSERT_VALUES, G));
3719566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(user->dm, localG, INSERT_VALUES, G));
372c4762a1bSJed Brown
3739566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(70.0 * xm * ym));
3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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
3887addb90fSBarry Smith . B - optionally different matrix used to construct the preconditioner
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 */
FormHessian(Tao tao,Vec X,Mat Hptr,Mat Hessian,void * ptr)415d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat Hptr, Mat Hessian, void *ptr)
416d71ae5a4SJacob 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
4303ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
431c4762a1bSJed Brown /* Set various matrix options */
4329566063dSJacob Faibussowitsch PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
433c4762a1bSJed Brown
434c4762a1bSJed Brown /* Initialize matrix entries to zero */
4359566063dSJacob Faibussowitsch PetscCall(MatAssembled(Hessian, &assembled));
4369566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(Hessian));
437c4762a1bSJed Brown
438c4762a1bSJed Brown /* Get local mesh boundaries */
4399566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
4409566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
441c4762a1bSJed Brown
442c4762a1bSJed Brown /* Scatter ghost points to local vector */
4439566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
4449566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
445c4762a1bSJed Brown
446c4762a1bSJed Brown /* Get pointers to vector data */
4479566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x));
4489566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Top, &top));
4499566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Bottom, &bottom));
4509566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Left, &left));
4519566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Right, &right));
452c4762a1bSJed Brown
453c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */
454c4762a1bSJed Brown
455c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
456c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
457c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs);
458c4762a1bSJed Brown
459c4762a1bSJed Brown xc = x[row];
460c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc;
461c4762a1bSJed Brown
462c4762a1bSJed Brown /* Left side */
463c4762a1bSJed Brown if (i == gxs) {
464c4762a1bSJed Brown xl = left[j - ys + 1];
465c4762a1bSJed Brown xlt = left[j - ys + 2];
466c4762a1bSJed Brown } else {
467c4762a1bSJed Brown xl = x[row - 1];
468c4762a1bSJed Brown }
469c4762a1bSJed Brown
470c4762a1bSJed Brown if (j == gys) {
471c4762a1bSJed Brown xb = bottom[i - xs + 1];
472c4762a1bSJed Brown xrb = bottom[i - xs + 2];
473c4762a1bSJed Brown } else {
474c4762a1bSJed Brown xb = x[row - gxm];
475c4762a1bSJed Brown }
476c4762a1bSJed Brown
477c4762a1bSJed Brown if (i + 1 == gxs + gxm) {
478c4762a1bSJed Brown xr = right[j - ys + 1];
479c4762a1bSJed Brown xrb = right[j - ys];
480c4762a1bSJed Brown } else {
481c4762a1bSJed Brown xr = x[row + 1];
482c4762a1bSJed Brown }
483c4762a1bSJed Brown
484c4762a1bSJed Brown if (j + 1 == gys + gym) {
485c4762a1bSJed Brown xt = top[i - xs + 1];
486c4762a1bSJed Brown xlt = top[i - xs];
487c4762a1bSJed Brown } else {
488c4762a1bSJed Brown xt = x[row + gxm];
489c4762a1bSJed Brown }
490c4762a1bSJed Brown
491ad540459SPierre Jolivet if (i > gxs && j + 1 < gys + gym) xlt = x[row - 1 + gxm];
492ad540459SPierre Jolivet if (j > gys && i + 1 < gxs + gxm) xrb = x[row + 1 - gxm];
493c4762a1bSJed Brown
494c4762a1bSJed Brown d1 = (xc - xl) * rhx;
495c4762a1bSJed Brown d2 = (xc - xr) * rhx;
496c4762a1bSJed Brown d3 = (xc - xt) * rhy;
497c4762a1bSJed Brown d4 = (xc - xb) * rhy;
498c4762a1bSJed Brown d5 = (xrb - xr) * rhy;
499c4762a1bSJed Brown d6 = (xrb - xb) * rhx;
500c4762a1bSJed Brown d7 = (xlt - xl) * rhy;
501c4762a1bSJed Brown d8 = (xlt - xt) * rhx;
502c4762a1bSJed Brown
503c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
504c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
505c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
506c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
507c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
508c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
509c4762a1bSJed Brown
5109371c9d4SSatish Balay hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
5119371c9d4SSatish Balay hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
5129371c9d4SSatish Balay ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
5139371c9d4SSatish Balay hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
514c4762a1bSJed Brown
515c4762a1bSJed Brown hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
516c4762a1bSJed Brown htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
517c4762a1bSJed Brown
5189371c9d4SSatish 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);
519c4762a1bSJed Brown
5209371c9d4SSatish Balay hl *= 0.5;
5219371c9d4SSatish Balay hr *= 0.5;
5229371c9d4SSatish Balay ht *= 0.5;
5239371c9d4SSatish Balay hb *= 0.5;
5249371c9d4SSatish Balay hbr *= 0.5;
5259371c9d4SSatish Balay htl *= 0.5;
5269371c9d4SSatish Balay hc *= 0.5;
527c4762a1bSJed Brown
528c4762a1bSJed Brown k = 0;
529c4762a1bSJed Brown if (j > 0) {
5309371c9d4SSatish Balay v[k] = hb;
5319371c9d4SSatish Balay col[k] = row - gxm;
5329371c9d4SSatish Balay k++;
533c4762a1bSJed Brown }
534c4762a1bSJed Brown
535c4762a1bSJed Brown if (j > 0 && i < mx - 1) {
5369371c9d4SSatish Balay v[k] = hbr;
5379371c9d4SSatish Balay col[k] = row - gxm + 1;
5389371c9d4SSatish Balay k++;
539c4762a1bSJed Brown }
540c4762a1bSJed Brown
541c4762a1bSJed Brown if (i > 0) {
5429371c9d4SSatish Balay v[k] = hl;
5439371c9d4SSatish Balay col[k] = row - 1;
5449371c9d4SSatish Balay k++;
545c4762a1bSJed Brown }
546c4762a1bSJed Brown
5479371c9d4SSatish Balay v[k] = hc;
5489371c9d4SSatish Balay col[k] = row;
5499371c9d4SSatish Balay k++;
550c4762a1bSJed Brown
551c4762a1bSJed Brown if (i < mx - 1) {
5529371c9d4SSatish Balay v[k] = hr;
5539371c9d4SSatish Balay col[k] = row + 1;
5549371c9d4SSatish Balay k++;
555c4762a1bSJed Brown }
556c4762a1bSJed Brown
557c4762a1bSJed Brown if (i > 0 && j < my - 1) {
5589371c9d4SSatish Balay v[k] = htl;
5599371c9d4SSatish Balay col[k] = row + gxm - 1;
5609371c9d4SSatish Balay k++;
561c4762a1bSJed Brown }
562c4762a1bSJed Brown
563c4762a1bSJed Brown if (j < my - 1) {
5649371c9d4SSatish Balay v[k] = ht;
5659371c9d4SSatish Balay col[k] = row + gxm;
5669371c9d4SSatish Balay k++;
567c4762a1bSJed Brown }
568c4762a1bSJed Brown
569c4762a1bSJed Brown /*
570c4762a1bSJed Brown Set matrix values using local numbering, which was defined
571c4762a1bSJed Brown earlier, in the main routine.
572c4762a1bSJed Brown */
5739566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Hessian, 1, &row, k, col, v, INSERT_VALUES));
574c4762a1bSJed Brown }
575c4762a1bSJed Brown }
576c4762a1bSJed Brown
577c4762a1bSJed Brown /* Restore vectors */
5789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x));
5799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Left, &left));
5809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Top, &top));
5819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Bottom, &bottom));
5829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Right, &right));
583c4762a1bSJed Brown
584c4762a1bSJed Brown /* Assemble the matrix */
5859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY));
5869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY));
587c4762a1bSJed Brown
5889566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199.0 * xm * ym));
5893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
590c4762a1bSJed Brown }
591c4762a1bSJed Brown
592c4762a1bSJed Brown /* ------------------------------------------------------------------- */
593c4762a1bSJed Brown /*
594c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for
595c4762a1bSJed Brown the region.
596c4762a1bSJed Brown
597c4762a1bSJed Brown Input Parameter:
598c4762a1bSJed Brown . user - user-defined application context
599c4762a1bSJed Brown
600c4762a1bSJed Brown Output Parameter:
601c4762a1bSJed Brown . user - user-defined application context
602c4762a1bSJed Brown */
MSA_BoundaryConditions(AppCtx * user)603d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user)
604d71ae5a4SJacob Faibussowitsch {
605c4762a1bSJed Brown PetscInt i, j, k, maxits = 5, limit = 0;
606c4762a1bSJed Brown PetscInt xs, ys, xm, ym, gxs, gys, gxm, gym;
607c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
608c4762a1bSJed Brown PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0;
609c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, three = 3.0, scl = 1.0, tol = 1e-10;
610c4762a1bSJed Brown PetscReal fnorm, det, hx, hy, xt = 0, yt = 0;
611c4762a1bSJed Brown PetscReal u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
612c4762a1bSJed Brown PetscReal b = -0.5, t = 0.5, l = -0.5, r = 0.5;
613c4762a1bSJed Brown PetscReal *boundary;
614c4762a1bSJed Brown PetscBool flg;
615c4762a1bSJed Brown Vec Bottom, Top, Right, Left;
616c4762a1bSJed Brown
6173ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
618c4762a1bSJed Brown /* Get local mesh boundaries */
6199566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
6209566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
621c4762a1bSJed Brown
622c4762a1bSJed Brown bsize = xm + 2;
623c4762a1bSJed Brown lsize = ym + 2;
624c4762a1bSJed Brown rsize = ym + 2;
625c4762a1bSJed Brown tsize = xm + 2;
626c4762a1bSJed Brown
6279566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(MPI_COMM_WORLD, bsize, PETSC_DECIDE, &Bottom));
6289566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(MPI_COMM_WORLD, tsize, PETSC_DECIDE, &Top));
6299566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(MPI_COMM_WORLD, lsize, PETSC_DECIDE, &Left));
6309566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(MPI_COMM_WORLD, rsize, PETSC_DECIDE, &Right));
631c4762a1bSJed Brown
632c4762a1bSJed Brown user->Top = Top;
633c4762a1bSJed Brown user->Left = Left;
634c4762a1bSJed Brown user->Bottom = Bottom;
635c4762a1bSJed Brown user->Right = Right;
636c4762a1bSJed Brown
6379371c9d4SSatish Balay hx = (r - l) / (mx + 1);
6389371c9d4SSatish Balay hy = (t - b) / (my + 1);
639c4762a1bSJed Brown
640c4762a1bSJed Brown for (j = 0; j < 4; j++) {
641c4762a1bSJed Brown if (j == 0) {
642c4762a1bSJed Brown yt = b;
643c4762a1bSJed Brown xt = l + hx * xs;
644c4762a1bSJed Brown limit = bsize;
6453ba16761SJacob Faibussowitsch PetscCall(VecGetArray(Bottom, &boundary));
646c4762a1bSJed Brown } else if (j == 1) {
647c4762a1bSJed Brown yt = t;
648c4762a1bSJed Brown xt = l + hx * xs;
649c4762a1bSJed Brown limit = tsize;
6503ba16761SJacob Faibussowitsch PetscCall(VecGetArray(Top, &boundary));
651c4762a1bSJed Brown } else if (j == 2) {
652c4762a1bSJed Brown yt = b + hy * ys;
653c4762a1bSJed Brown xt = l;
654c4762a1bSJed Brown limit = lsize;
6553ba16761SJacob Faibussowitsch PetscCall(VecGetArray(Left, &boundary));
656c4762a1bSJed Brown } else if (j == 3) {
657c4762a1bSJed Brown yt = b + hy * ys;
658c4762a1bSJed Brown xt = r;
659c4762a1bSJed Brown limit = rsize;
6603ba16761SJacob Faibussowitsch PetscCall(VecGetArray(Right, &boundary));
661c4762a1bSJed Brown }
662c4762a1bSJed Brown
663c4762a1bSJed Brown for (i = 0; i < limit; i++) {
664c4762a1bSJed Brown u1 = xt;
665c4762a1bSJed Brown u2 = -yt;
666c4762a1bSJed Brown for (k = 0; k < maxits; k++) {
667c4762a1bSJed Brown nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
668c4762a1bSJed Brown nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
669c4762a1bSJed Brown fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2);
670c4762a1bSJed Brown if (fnorm <= tol) break;
671c4762a1bSJed Brown njac11 = one + u2 * u2 - u1 * u1;
672c4762a1bSJed Brown njac12 = two * u1 * u2;
673c4762a1bSJed Brown njac21 = -two * u1 * u2;
674c4762a1bSJed Brown njac22 = -one - u1 * u1 + u2 * u2;
675c4762a1bSJed Brown det = njac11 * njac22 - njac21 * njac12;
676c4762a1bSJed Brown u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det;
677c4762a1bSJed Brown u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det;
678c4762a1bSJed Brown }
679c4762a1bSJed Brown
680c4762a1bSJed Brown boundary[i] = u1 * u1 - u2 * u2;
681c4762a1bSJed Brown if (j == 0 || j == 1) {
682c4762a1bSJed Brown xt = xt + hx;
683c4762a1bSJed Brown } else if (j == 2 || j == 3) {
684c4762a1bSJed Brown yt = yt + hy;
685c4762a1bSJed Brown }
686c4762a1bSJed Brown }
687c4762a1bSJed Brown if (j == 0) {
6889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Bottom, &boundary));
689c4762a1bSJed Brown } else if (j == 1) {
6909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Top, &boundary));
691c4762a1bSJed Brown } else if (j == 2) {
6929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Left, &boundary));
693c4762a1bSJed Brown } else if (j == 3) {
6949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Right, &boundary));
695c4762a1bSJed Brown }
696c4762a1bSJed Brown }
697c4762a1bSJed Brown
698c4762a1bSJed Brown /* Scale the boundary if desired */
699c4762a1bSJed Brown
7009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg));
7011baa6e33SBarry Smith if (flg) PetscCall(VecScale(Bottom, scl));
7029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg));
7031baa6e33SBarry Smith if (flg) PetscCall(VecScale(Top, scl));
7049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg));
7051baa6e33SBarry Smith if (flg) PetscCall(VecScale(Right, scl));
706c4762a1bSJed Brown
7079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg));
7081baa6e33SBarry Smith if (flg) PetscCall(VecScale(Left, scl));
7093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
710c4762a1bSJed Brown }
711c4762a1bSJed Brown
712c4762a1bSJed Brown /* ------------------------------------------------------------------- */
713c4762a1bSJed Brown /*
714c4762a1bSJed Brown MSA_Plate - Calculates an obstacle for surface to stretch over.
715c4762a1bSJed Brown
716c4762a1bSJed Brown Input Parameter:
717c4762a1bSJed Brown . user - user-defined application context
718c4762a1bSJed Brown
719c4762a1bSJed Brown Output Parameter:
720c4762a1bSJed Brown . user - user-defined application context
721c4762a1bSJed Brown */
MSA_Plate(Vec XL,Vec XU,PetscCtx ctx)7222a8381b2SBarry Smith static PetscErrorCode MSA_Plate(Vec XL, Vec XU, PetscCtx ctx)
723d71ae5a4SJacob Faibussowitsch {
724c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
725c4762a1bSJed Brown PetscInt i, j, row;
726c4762a1bSJed Brown PetscInt xs, ys, xm, ym;
727c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my, bmy, bmx;
728c4762a1bSJed Brown PetscReal t1, t2, t3;
729c4762a1bSJed Brown PetscReal *xl, lb = PETSC_NINFINITY, ub = PETSC_INFINITY;
730c4762a1bSJed Brown PetscBool cylinder;
731c4762a1bSJed Brown
7323ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
7339371c9d4SSatish Balay user->bmy = PetscMax(0, user->bmy);
7349371c9d4SSatish Balay user->bmy = PetscMin(my, user->bmy);
7359371c9d4SSatish Balay user->bmx = PetscMax(0, user->bmx);
7369371c9d4SSatish Balay user->bmx = PetscMin(mx, user->bmx);
7379371c9d4SSatish Balay bmy = user->bmy;
7389371c9d4SSatish Balay bmx = user->bmx;
739c4762a1bSJed Brown
7409566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
741c4762a1bSJed Brown
7429566063dSJacob Faibussowitsch PetscCall(VecSet(XL, lb));
7439566063dSJacob Faibussowitsch PetscCall(VecSet(XU, ub));
744c4762a1bSJed Brown
7459566063dSJacob Faibussowitsch PetscCall(VecGetArray(XL, &xl));
746c4762a1bSJed Brown
7479566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-cylinder", &cylinder));
748c4762a1bSJed Brown /* Compute the optional lower box */
749c4762a1bSJed Brown if (cylinder) {
750c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
751c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
752c4762a1bSJed Brown row = (j - ys) * xm + (i - xs);
753c4762a1bSJed Brown t1 = (2.0 * i - mx) * bmy;
754c4762a1bSJed Brown t2 = (2.0 * j - my) * bmx;
755c4762a1bSJed Brown t3 = bmx * bmx * bmy * bmy;
756ad540459SPierre Jolivet if (t1 * t1 + t2 * t2 <= t3) xl[row] = user->bheight;
757c4762a1bSJed Brown }
758c4762a1bSJed Brown }
759c4762a1bSJed Brown } else {
760c4762a1bSJed Brown /* Compute the optional lower box */
761c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
762c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
763c4762a1bSJed Brown row = (j - ys) * xm + (i - xs);
764ad540459SPierre Jolivet if (i >= (mx - bmx) / 2 && i < mx - (mx - bmx) / 2 && j >= (my - bmy) / 2 && j < my - (my - bmy) / 2) xl[row] = user->bheight;
765c4762a1bSJed Brown }
766c4762a1bSJed Brown }
767c4762a1bSJed Brown }
7689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(XL, &xl));
7693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
770c4762a1bSJed Brown }
771c4762a1bSJed Brown
772c4762a1bSJed Brown /* ------------------------------------------------------------------- */
773c4762a1bSJed Brown /*
774c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways.
775c4762a1bSJed Brown
776c4762a1bSJed Brown Input Parameters:
777c4762a1bSJed Brown . user - user-defined application context
778c4762a1bSJed Brown . X - vector for initial guess
779c4762a1bSJed Brown
780c4762a1bSJed Brown Output Parameters:
781c4762a1bSJed Brown . X - newly computed initial guess
782c4762a1bSJed Brown */
MSA_InitialPoint(AppCtx * user,Vec X)783d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X)
784d71ae5a4SJacob Faibussowitsch {
785c4762a1bSJed Brown PetscInt start = -1, i, j;
786c4762a1bSJed Brown PetscReal zero = 0.0;
787c4762a1bSJed Brown PetscBool flg;
788c4762a1bSJed Brown
7893ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
7909566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg));
791c4762a1bSJed Brown if (flg && start == 0) { /* The zero vector is reasonable */
7929566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero));
793c4762a1bSJed Brown } else if (flg && start > 0) { /* Try a random start between -0.5 and 0.5 */
7949371c9d4SSatish Balay PetscRandom rctx;
7959371c9d4SSatish Balay PetscReal np5 = -0.5;
796c4762a1bSJed Brown
7979566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(MPI_COMM_WORLD, &rctx));
79848a46eb9SPierre Jolivet for (i = 0; i < start; i++) PetscCall(VecSetRandom(X, rctx));
7999566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx));
8009566063dSJacob Faibussowitsch PetscCall(VecShift(X, np5));
801c4762a1bSJed Brown
802c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */
803c4762a1bSJed Brown
804c4762a1bSJed Brown PetscInt row, xs, xm, gxs, gxm, ys, ym, gys, gym;
805c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
806c4762a1bSJed Brown PetscReal *x, *left, *right, *bottom, *top;
807c4762a1bSJed Brown Vec localX = user->localX;
808c4762a1bSJed Brown
809c4762a1bSJed Brown /* Get local mesh boundaries */
8109566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
8119566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
812c4762a1bSJed Brown
813c4762a1bSJed Brown /* Get pointers to vector data */
8149566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Top, &top));
8159566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Bottom, &bottom));
8169566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Left, &left));
8179566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->Right, &right));
818c4762a1bSJed Brown
8199566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x));
820c4762a1bSJed Brown /* Perform local computations */
821c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
822c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
823c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs);
8242da392ccSBarry 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;
825c4762a1bSJed Brown }
826c4762a1bSJed Brown }
827c4762a1bSJed Brown
828c4762a1bSJed Brown /* Restore vectors */
8299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x));
830c4762a1bSJed Brown
8319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Left, &left));
8329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Top, &top));
8339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Bottom, &bottom));
8349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->Right, &right));
835c4762a1bSJed Brown
836c4762a1bSJed Brown /* Scatter values into global vector */
8379566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(user->dm, localX, INSERT_VALUES, X));
8389566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(user->dm, localX, INSERT_VALUES, X));
839c4762a1bSJed Brown }
8403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
841c4762a1bSJed Brown }
842c4762a1bSJed Brown
84301c1178eSBarry Smith /* For testing matrix-free submatrices */
MatrixFreeHessian(Tao tao,Vec x,Mat H,Mat Hpre,void * ptr)844d71ae5a4SJacob Faibussowitsch PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
845d71ae5a4SJacob Faibussowitsch {
846c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
8473ba16761SJacob Faibussowitsch
848c4762a1bSJed Brown PetscFunctionBegin;
8499566063dSJacob Faibussowitsch PetscCall(FormHessian(tao, x, user->H, user->H, ptr));
8503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
851c4762a1bSJed Brown }
8523ba16761SJacob Faibussowitsch
MyMatMult(Mat H_shell,Vec X,Vec Y)853d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
854d71ae5a4SJacob Faibussowitsch {
855c4762a1bSJed Brown void *ptr;
856c4762a1bSJed Brown AppCtx *user;
8573ba16761SJacob Faibussowitsch
858c4762a1bSJed Brown PetscFunctionBegin;
8599566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(H_shell, &ptr));
860c4762a1bSJed Brown user = (AppCtx *)ptr;
8619566063dSJacob Faibussowitsch PetscCall(MatMult(user->H, X, Y));
8623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
863c4762a1bSJed Brown }
864c4762a1bSJed Brown
865c4762a1bSJed Brown /*TEST
866c4762a1bSJed Brown
867c4762a1bSJed Brown build:
868c4762a1bSJed Brown requires: !complex
869c4762a1bSJed Brown
870c4762a1bSJed Brown test:
87110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
872c4762a1bSJed Brown requires: !single
873c4762a1bSJed Brown
874c4762a1bSJed Brown test:
875c4762a1bSJed Brown suffix: 2
876c4762a1bSJed Brown nsize: 2
87710978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
878c4762a1bSJed Brown requires: !single
879c4762a1bSJed Brown
880c4762a1bSJed Brown test:
881c4762a1bSJed Brown suffix: 3
882c4762a1bSJed Brown nsize: 3
88310978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
884c4762a1bSJed Brown requires: !single
885c4762a1bSJed Brown
886c4762a1bSJed Brown test:
887c4762a1bSJed Brown suffix: 4
888c4762a1bSJed Brown nsize: 3
88910978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type mask -tao_type tron -tao_gatol 1.e-5
890c4762a1bSJed Brown requires: !single
891c4762a1bSJed Brown
892c4762a1bSJed Brown test:
893c4762a1bSJed Brown suffix: 5
894c4762a1bSJed Brown nsize: 3
89510978b7dSBarry Smith args: -tao_monitor_short -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
896c4762a1bSJed Brown requires: !single
897c4762a1bSJed Brown
898c4762a1bSJed Brown test:
899c4762a1bSJed Brown suffix: 6
900c4762a1bSJed Brown nsize: 3
90110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -tao_type blmvm -tao_gatol 1.e-4
902c4762a1bSJed Brown requires: !single
903c4762a1bSJed Brown
904c4762a1bSJed Brown test:
905c4762a1bSJed Brown suffix: 7
906c4762a1bSJed Brown nsize: 3
90710978b7dSBarry Smith args: -tao_monitor_short -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
908c4762a1bSJed Brown requires: !single
909c4762a1bSJed Brown
910c4762a1bSJed Brown test:
911c4762a1bSJed Brown suffix: 8
91210978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
913c4762a1bSJed Brown requires: !single
914c4762a1bSJed Brown
915c4762a1bSJed Brown test:
916c4762a1bSJed Brown suffix: 9
91710978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
918c4762a1bSJed Brown requires: !single
919c4762a1bSJed Brown
920c4762a1bSJed Brown test:
921c4762a1bSJed Brown suffix: 10
92210978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
923c4762a1bSJed Brown requires: !single
924c4762a1bSJed Brown
925c4762a1bSJed Brown test:
926c4762a1bSJed Brown suffix: 11
92710978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
928c4762a1bSJed Brown requires: !single
929c4762a1bSJed Brown
930c4762a1bSJed Brown test:
931c4762a1bSJed Brown suffix: 12
93210978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
933c4762a1bSJed Brown requires: !single
934c4762a1bSJed Brown
935c4762a1bSJed Brown test:
936c4762a1bSJed Brown suffix: 13
937*a336c150SZach Atkins args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 -tao_bnk_cg_tao_monitor_short
938c4762a1bSJed Brown requires: !single
939c4762a1bSJed Brown
940c4762a1bSJed Brown test:
941c4762a1bSJed Brown suffix: 14
942*a336c150SZach Atkins args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 -tao_bnk_cg_tao_monitor_short
943c4762a1bSJed Brown requires: !single
944c4762a1bSJed Brown
945c4762a1bSJed Brown test:
946c4762a1bSJed Brown suffix: 15
947*a336c150SZach Atkins args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 -tao_bnk_cg_tao_monitor_short
948c4762a1bSJed Brown requires: !single
949c4762a1bSJed Brown
950c4762a1bSJed Brown test:
951c4762a1bSJed Brown suffix: 16
95210978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
953c4762a1bSJed Brown requires: !single
954c4762a1bSJed Brown
955c4762a1bSJed Brown test:
956c4762a1bSJed Brown suffix: 17
95710978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
958c4762a1bSJed Brown requires: !single
959c4762a1bSJed Brown
96034ad9904SAlp Dener test:
96134ad9904SAlp Dener suffix: 18
96210978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
96334ad9904SAlp Dener requires: !single
96434ad9904SAlp Dener
96534ad9904SAlp Dener test:
96634ad9904SAlp Dener suffix: 19
96710978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
96834ad9904SAlp Dener requires: !single
96934ad9904SAlp Dener
97034ad9904SAlp Dener test:
97134ad9904SAlp Dener suffix: 20
97210978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
973910b42cbSStefano Zampini requires: !single
97434ad9904SAlp Dener
975c4762a1bSJed Brown TEST*/
976