xref: /petsc/src/tao/unconstrained/tutorials/minsurf2.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown /* Program usage: mpiexec -n <proc> minsurf2 [-help] [all TAO options] */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   Include "petsctao.h" so we can use TAO solvers.
5c4762a1bSJed Brown   petscdmda.h for distributed array
6c4762a1bSJed Brown */
7c4762a1bSJed Brown #include <petsctao.h>
8c4762a1bSJed Brown #include <petscdmda.h>
9c4762a1bSJed Brown 
109371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to \n\
11c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
12c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and \n\
13c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
14c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
15c4762a1bSJed Brown The command line options are:\n\
16c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
17c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
18c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
19c4762a1bSJed Brown                for an average of the boundary conditions\n\n";
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown    User-defined application context - contains data needed by the
23c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient()
24c4762a1bSJed Brown    and FormHessian().
25c4762a1bSJed Brown */
26c4762a1bSJed Brown typedef struct {
27c4762a1bSJed Brown   PetscInt   mx, my;                      /* discretization in x, y directions */
28c4762a1bSJed Brown   PetscReal *bottom, *top, *left, *right; /* boundary values */
29c4762a1bSJed Brown   DM         dm;                          /* distributed array data structure */
30c4762a1bSJed Brown   Mat        H;                           /* Hessian */
31c4762a1bSJed Brown } AppCtx;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown /* -------- User-defined Routines --------- */
34c4762a1bSJed Brown 
35c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
36c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
37c4762a1bSJed Brown PetscErrorCode        QuadraticH(AppCtx *, Vec, Mat);
38c4762a1bSJed Brown PetscErrorCode        FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
39c4762a1bSJed Brown PetscErrorCode        FormGradient(Tao, Vec, Vec, void *);
40c4762a1bSJed Brown PetscErrorCode        FormHessian(Tao, Vec, Mat, Mat, void *);
41c4762a1bSJed Brown PetscErrorCode        My_Monitor(Tao, void *);
42c4762a1bSJed Brown 
439371c9d4SSatish Balay int main(int argc, char **argv) {
44c4762a1bSJed Brown   PetscInt      Nx, Ny;                /* number of processors in x- and y- directions */
45c4762a1bSJed Brown   Vec           x;                     /* solution, gradient vectors */
46c4762a1bSJed Brown   PetscBool     flg, viewmat;          /* flags */
47c4762a1bSJed Brown   PetscBool     fddefault, fdcoloring; /* flags */
48c4762a1bSJed Brown   Tao           tao;                   /* TAO solver context */
49c4762a1bSJed Brown   AppCtx        user;                  /* user-defined work context */
50c4762a1bSJed Brown   ISColoring    iscoloring;
51c4762a1bSJed Brown   MatFDColoring matfdcoloring;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Initialize TAO */
54327415f7SBarry Smith   PetscFunctionBeginUser;
559566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* Specify dimension of the problem */
589371c9d4SSatish Balay   user.mx = 10;
599371c9d4SSatish Balay   user.my = 10;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(MPI_COMM_WORLD, "\n---- Minimum Surface Area Problem -----\n"));
6663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(MPI_COMM_WORLD, "mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n", user.mx, user.my));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Let PETSc determine the vector distribution */
699371c9d4SSatish Balay   Nx = PETSC_DECIDE;
709371c9d4SSatish Balay   Ny = PETSC_DECIDE;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Create distributed array (DM) to manage parallel grid and vectors  */
739566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, user.mx, user.my, Nx, Ny, 1, 1, NULL, NULL, &user.dm));
749566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
759566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Create TAO solver and set desired solution method.*/
789566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
799566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOCG));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /*
82c4762a1bSJed Brown      Extract global vector from DA for the vector of variables --  PETSC routine
83c4762a1bSJed Brown      Compute the initial solution                              --  application specific, see below
84c4762a1bSJed Brown      Set this vector for use by TAO                            --  TAO routine
85c4762a1bSJed Brown   */
869566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm, &x));
879566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user));
889566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user, x));
899566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /*
92c4762a1bSJed Brown      Initialize the Application context for use in function evaluations  --  application specific, see below.
93c4762a1bSJed Brown      Set routines for function and gradient evaluation
94c4762a1bSJed Brown   */
959566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /*
98c4762a1bSJed Brown      Given the command line arguments, calculate the hessian with either the user-
99c4762a1bSJed Brown      provided function FormHessian, or the default finite-difference driven Hessian
100c4762a1bSJed Brown      functions
101c4762a1bSJed Brown   */
1029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fddefault", &fddefault));
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fdcoloring", &fdcoloring));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /*
106c4762a1bSJed Brown      Create a matrix data structure to store the Hessian and set
107c4762a1bSJed Brown      the Hessian evalution routine.
108c4762a1bSJed Brown      Set the matrix structure to be used for Hessian evalutions
109c4762a1bSJed Brown   */
1109566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.dm, &user.H));
1119566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   if (fdcoloring) {
1149566063dSJacob Faibussowitsch     PetscCall(DMCreateColoring(user.dm, IS_COLORING_GLOBAL, &iscoloring));
1159566063dSJacob Faibussowitsch     PetscCall(MatFDColoringCreate(user.H, iscoloring, &matfdcoloring));
1169566063dSJacob Faibussowitsch     PetscCall(ISColoringDestroy(&iscoloring));
1179566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))FormGradient, (void *)&user));
1189566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1199566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessianColor, (void *)matfdcoloring));
120c4762a1bSJed Brown   } else if (fddefault) {
1219566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessian, (void *)NULL));
122c4762a1bSJed Brown   } else {
1239566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /*
127c4762a1bSJed Brown      If my_monitor option is in command line, then use the user-provided
128c4762a1bSJed Brown      monitoring function
129c4762a1bSJed Brown   */
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-my_monitor", &viewmat));
1311baa6e33SBarry Smith   if (viewmat) PetscCall(TaoSetMonitor(tao, My_Monitor, NULL, NULL));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* Check for any tao command line options */
1349566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1379566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
138c4762a1bSJed Brown 
1399566063dSJacob Faibussowitsch   PetscCall(TaoView(tao, PETSC_VIEWER_STDOUT_WORLD));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Free TAO data structures */
1429566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Free PETSc data structures */
1459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
147*48a46eb9SPierre Jolivet   if (fdcoloring) PetscCall(MatFDColoringDestroy(&matfdcoloring));
1489566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.bottom));
1499566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.top));
1509566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.left));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.right));
1529566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
1539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
154b122ec5aSJacob Faibussowitsch   return 0;
155c4762a1bSJed Brown }
156c4762a1bSJed Brown 
1579371c9d4SSatish Balay PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *userCtx) {
158c4762a1bSJed Brown   PetscReal fcn;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBegin;
1619566063dSJacob Faibussowitsch   PetscCall(FormFunctionGradient(tao, X, &fcn, G, userCtx));
162c4762a1bSJed Brown   PetscFunctionReturn(0);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown /* -------------------------------------------------------------------- */
166c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates the function and corresponding gradient.
167c4762a1bSJed Brown 
168c4762a1bSJed Brown     Input Parameters:
169c4762a1bSJed Brown .   tao     - the Tao context
170c4762a1bSJed Brown .   XX      - input vector
171a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
172c4762a1bSJed Brown 
173c4762a1bSJed Brown     Output Parameters:
174c4762a1bSJed Brown .   fcn     - the newly evaluated function
175c4762a1bSJed Brown .   GG       - vector containing the newly evaluated gradient
176c4762a1bSJed Brown */
1779371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx) {
178c4762a1bSJed Brown   AppCtx     *user = (AppCtx *)userCtx;
179c4762a1bSJed Brown   PetscInt    i, j;
180c4762a1bSJed Brown   PetscInt    mx = user->mx, my = user->my;
181c4762a1bSJed Brown   PetscInt    xs, xm, gxs, gxm, ys, ym, gys, gym;
182c4762a1bSJed Brown   PetscReal   ft = 0.0;
183c4762a1bSJed Brown   PetscReal   hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy;
184c4762a1bSJed Brown   PetscReal   rhx = mx + 1, rhy = my + 1;
185c4762a1bSJed Brown   PetscReal   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
186c4762a1bSJed Brown   PetscReal   df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
187c4762a1bSJed Brown   PetscReal **g, **x;
188c4762a1bSJed Brown   Vec         localX;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBegin;
191c4762a1bSJed Brown   /* Get local mesh boundaries */
1929566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->dm, &localX));
1939566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
1949566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /* Scatter ghost points to local vector */
1979566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
1989566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* Get pointers to vector data */
2019566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x));
2029566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->dm, G, (void **)&g));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* Compute function and gradient over the locally owned part of the mesh */
205c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
206c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
207c4762a1bSJed Brown       xc  = x[j][i];
208c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
209c4762a1bSJed Brown 
210c4762a1bSJed Brown       if (i == 0) { /* left side */
211c4762a1bSJed Brown         xl  = user->left[j - ys + 1];
212c4762a1bSJed Brown         xlt = user->left[j - ys + 2];
213c4762a1bSJed Brown       } else {
214c4762a1bSJed Brown         xl = x[j][i - 1];
215c4762a1bSJed Brown       }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown       if (j == 0) { /* bottom side */
218c4762a1bSJed Brown         xb  = user->bottom[i - xs + 1];
219c4762a1bSJed Brown         xrb = user->bottom[i - xs + 2];
220c4762a1bSJed Brown       } else {
221c4762a1bSJed Brown         xb = x[j - 1][i];
222c4762a1bSJed Brown       }
223c4762a1bSJed Brown 
224c4762a1bSJed Brown       if (i + 1 == gxs + gxm) { /* right side */
225c4762a1bSJed Brown         xr  = user->right[j - ys + 1];
226c4762a1bSJed Brown         xrb = user->right[j - ys];
227c4762a1bSJed Brown       } else {
228c4762a1bSJed Brown         xr = x[j][i + 1];
229c4762a1bSJed Brown       }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown       if (j + 1 == gys + gym) { /* top side */
232c4762a1bSJed Brown         xt  = user->top[i - xs + 1];
233c4762a1bSJed Brown         xlt = user->top[i - xs];
234c4762a1bSJed Brown       } else {
235c4762a1bSJed Brown         xt = x[j + 1][i];
236c4762a1bSJed Brown       }
237c4762a1bSJed Brown 
2389371c9d4SSatish Balay       if (i > gxs && j + 1 < gys + gym) { xlt = x[j + 1][i - 1]; }
2399371c9d4SSatish Balay       if (j > gys && i + 1 < gxs + gxm) { xrb = x[j - 1][i + 1]; }
240c4762a1bSJed Brown 
241c4762a1bSJed Brown       d1 = (xc - xl);
242c4762a1bSJed Brown       d2 = (xc - xr);
243c4762a1bSJed Brown       d3 = (xc - xt);
244c4762a1bSJed Brown       d4 = (xc - xb);
245c4762a1bSJed Brown       d5 = (xr - xrb);
246c4762a1bSJed Brown       d6 = (xrb - xb);
247c4762a1bSJed Brown       d7 = (xlt - xl);
248c4762a1bSJed Brown       d8 = (xt - xlt);
249c4762a1bSJed Brown 
250c4762a1bSJed Brown       df1dxc = d1 * hydhx;
251c4762a1bSJed Brown       df2dxc = (d1 * hydhx + d4 * hxdhy);
252c4762a1bSJed Brown       df3dxc = d3 * hxdhy;
253c4762a1bSJed Brown       df4dxc = (d2 * hydhx + d3 * hxdhy);
254c4762a1bSJed Brown       df5dxc = d2 * hydhx;
255c4762a1bSJed Brown       df6dxc = d4 * hxdhy;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown       d1 *= rhx;
258c4762a1bSJed Brown       d2 *= rhx;
259c4762a1bSJed Brown       d3 *= rhy;
260c4762a1bSJed Brown       d4 *= rhy;
261c4762a1bSJed Brown       d5 *= rhy;
262c4762a1bSJed Brown       d6 *= rhx;
263c4762a1bSJed Brown       d7 *= rhy;
264c4762a1bSJed Brown       d8 *= rhx;
265c4762a1bSJed Brown 
266c4762a1bSJed Brown       f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7);
267c4762a1bSJed Brown       f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
268c4762a1bSJed Brown       f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8);
269c4762a1bSJed Brown       f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
270c4762a1bSJed Brown       f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5);
271c4762a1bSJed Brown       f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6);
272c4762a1bSJed Brown 
273c4762a1bSJed Brown       ft = ft + (f2 + f4);
274c4762a1bSJed Brown 
275c4762a1bSJed Brown       df1dxc /= f1;
276c4762a1bSJed Brown       df2dxc /= f2;
277c4762a1bSJed Brown       df3dxc /= f3;
278c4762a1bSJed Brown       df4dxc /= f4;
279c4762a1bSJed Brown       df5dxc /= f5;
280c4762a1bSJed Brown       df6dxc /= f6;
281c4762a1bSJed Brown 
282c4762a1bSJed Brown       g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5;
283c4762a1bSJed Brown     }
284c4762a1bSJed Brown   }
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
287c4762a1bSJed Brown   if (xs == 0) { /* left side */
288c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
289c4762a1bSJed Brown       d3 = (user->left[j - ys + 1] - user->left[j - ys + 2]) * rhy;
290c4762a1bSJed Brown       d2 = (user->left[j - ys + 1] - x[j][0]) * rhx;
291c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
292c4762a1bSJed Brown     }
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown   if (ys == 0) { /* bottom side */
295c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
296c4762a1bSJed Brown       d2 = (user->bottom[i + 1 - xs] - user->bottom[i - xs + 2]) * rhx;
297c4762a1bSJed Brown       d3 = (user->bottom[i - xs + 1] - x[0][i]) * rhy;
298c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
299c4762a1bSJed Brown     }
300c4762a1bSJed Brown   }
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   if (xs + xm == mx) { /* right side */
303c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
304c4762a1bSJed Brown       d1 = (x[j][mx - 1] - user->right[j - ys + 1]) * rhx;
305c4762a1bSJed Brown       d4 = (user->right[j - ys] - user->right[j - ys + 1]) * rhy;
306c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
307c4762a1bSJed Brown     }
308c4762a1bSJed Brown   }
309c4762a1bSJed Brown   if (ys + ym == my) { /* top side */
310c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
311c4762a1bSJed Brown       d1 = (x[my - 1][i] - user->top[i - xs + 1]) * rhy;
312c4762a1bSJed Brown       d4 = (user->top[i - xs + 1] - user->top[i - xs]) * rhx;
313c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
314c4762a1bSJed Brown     }
315c4762a1bSJed Brown   }
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   if (ys == 0 && xs == 0) {
318c4762a1bSJed Brown     d1 = (user->left[0] - user->left[1]) / hy;
319c4762a1bSJed Brown     d2 = (user->bottom[0] - user->bottom[1]) * rhx;
320c4762a1bSJed Brown     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
321c4762a1bSJed Brown   }
322c4762a1bSJed Brown   if (ys + ym == my && xs + xm == mx) {
323c4762a1bSJed Brown     d1 = (user->right[ym + 1] - user->right[ym]) * rhy;
324c4762a1bSJed Brown     d2 = (user->top[xm + 1] - user->top[xm]) * rhx;
325c4762a1bSJed Brown     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
326c4762a1bSJed Brown   }
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   ft = ft * area;
3299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   /* Restore vectors */
3329566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x));
3339566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->dm, G, (void **)&g));
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   /* Scatter values to global vector */
3369566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->dm, &localX));
3379566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(67.0 * xm * ym));
338c4762a1bSJed Brown   PetscFunctionReturn(0);
339c4762a1bSJed Brown }
340c4762a1bSJed Brown 
341c4762a1bSJed Brown /* ------------------------------------------------------------------- */
342c4762a1bSJed Brown /*
343c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
344c4762a1bSJed Brown 
345c4762a1bSJed Brown    Input Parameters:
346c4762a1bSJed Brown .  tao  - the Tao context
347c4762a1bSJed Brown .  x    - input vector
348a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
349c4762a1bSJed Brown 
350c4762a1bSJed Brown    Output Parameters:
351c4762a1bSJed Brown .  H    - Hessian matrix
352c4762a1bSJed Brown .  Hpre - optionally different preconditioning matrix
353c4762a1bSJed Brown .  flg  - flag indicating matrix structure
354c4762a1bSJed Brown 
355c4762a1bSJed Brown */
3569371c9d4SSatish Balay PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) {
357c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   PetscFunctionBegin;
360c4762a1bSJed Brown   /* Evaluate the Hessian entries*/
3619566063dSJacob Faibussowitsch   PetscCall(QuadraticH(user, X, H));
362c4762a1bSJed Brown   PetscFunctionReturn(0);
363c4762a1bSJed Brown }
364c4762a1bSJed Brown 
365c4762a1bSJed Brown /* ------------------------------------------------------------------- */
366c4762a1bSJed Brown /*
367c4762a1bSJed Brown    QuadraticH - Evaluates Hessian matrix.
368c4762a1bSJed Brown 
369c4762a1bSJed Brown    Input Parameters:
370c4762a1bSJed Brown .  user - user-defined context, as set by TaoSetHessian()
371c4762a1bSJed Brown .  X    - input vector
372c4762a1bSJed Brown 
373c4762a1bSJed Brown    Output Parameter:
374c4762a1bSJed Brown .  H    - Hessian matrix
375c4762a1bSJed Brown */
3769371c9d4SSatish Balay PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) {
377c4762a1bSJed Brown   PetscInt    i, j, k;
378c4762a1bSJed Brown   PetscInt    mx = user->mx, my = user->my;
379c4762a1bSJed Brown   PetscInt    xs, xm, gxs, gxm, ys, ym, gys, gym;
380c4762a1bSJed Brown   PetscReal   hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
381c4762a1bSJed Brown   PetscReal   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
382c4762a1bSJed Brown   PetscReal   hl, hr, ht, hb, hc, htl, hbr;
383c4762a1bSJed Brown   PetscReal **x, v[7];
384c4762a1bSJed Brown   MatStencil  col[7], row;
385c4762a1bSJed Brown   Vec         localX;
386c4762a1bSJed Brown   PetscBool   assembled;
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   PetscFunctionBegin;
389c4762a1bSJed Brown   /* Get local mesh boundaries */
3909566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->dm, &localX));
391c4762a1bSJed Brown 
3929566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
3939566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   /* Scatter ghost points to local vector */
3969566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
3979566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
398c4762a1bSJed Brown 
399c4762a1bSJed Brown   /* Get pointers to vector data */
4009566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x));
401c4762a1bSJed Brown 
402c4762a1bSJed Brown   /* Initialize matrix entries to zero */
4039566063dSJacob Faibussowitsch   PetscCall(MatAssembled(Hessian, &assembled));
4049566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(Hessian));
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   /* Set various matrix options */
4079566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
412c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
413c4762a1bSJed Brown       xc  = x[j][i];
414c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
415c4762a1bSJed Brown 
416c4762a1bSJed Brown       /* Left side */
417c4762a1bSJed Brown       if (i == 0) {
418c4762a1bSJed Brown         xl  = user->left[j - ys + 1];
419c4762a1bSJed Brown         xlt = user->left[j - ys + 2];
420c4762a1bSJed Brown       } else {
421c4762a1bSJed Brown         xl = x[j][i - 1];
422c4762a1bSJed Brown       }
423c4762a1bSJed Brown 
424c4762a1bSJed Brown       if (j == 0) {
425c4762a1bSJed Brown         xb  = user->bottom[i - xs + 1];
426c4762a1bSJed Brown         xrb = user->bottom[i - xs + 2];
427c4762a1bSJed Brown       } else {
428c4762a1bSJed Brown         xb = x[j - 1][i];
429c4762a1bSJed Brown       }
430c4762a1bSJed Brown 
431c4762a1bSJed Brown       if (i + 1 == mx) {
432c4762a1bSJed Brown         xr  = user->right[j - ys + 1];
433c4762a1bSJed Brown         xrb = user->right[j - ys];
434c4762a1bSJed Brown       } else {
435c4762a1bSJed Brown         xr = x[j][i + 1];
436c4762a1bSJed Brown       }
437c4762a1bSJed Brown 
438c4762a1bSJed Brown       if (j + 1 == my) {
439c4762a1bSJed Brown         xt  = user->top[i - xs + 1];
440c4762a1bSJed Brown         xlt = user->top[i - xs];
441c4762a1bSJed Brown       } else {
442c4762a1bSJed Brown         xt = x[j + 1][i];
443c4762a1bSJed Brown       }
444c4762a1bSJed Brown 
4459371c9d4SSatish Balay       if (i > 0 && j + 1 < my) { xlt = x[j + 1][i - 1]; }
4469371c9d4SSatish Balay       if (j > 0 && i + 1 < mx) { xrb = x[j - 1][i + 1]; }
447c4762a1bSJed Brown 
448c4762a1bSJed Brown       d1 = (xc - xl) / hx;
449c4762a1bSJed Brown       d2 = (xc - xr) / hx;
450c4762a1bSJed Brown       d3 = (xc - xt) / hy;
451c4762a1bSJed Brown       d4 = (xc - xb) / hy;
452c4762a1bSJed Brown       d5 = (xrb - xr) / hy;
453c4762a1bSJed Brown       d6 = (xrb - xb) / hx;
454c4762a1bSJed Brown       d7 = (xlt - xl) / hy;
455c4762a1bSJed Brown       d8 = (xlt - xt) / hx;
456c4762a1bSJed Brown 
457c4762a1bSJed Brown       f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7);
458c4762a1bSJed Brown       f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
459c4762a1bSJed Brown       f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8);
460c4762a1bSJed Brown       f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
461c4762a1bSJed Brown       f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5);
462c4762a1bSJed Brown       f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6);
463c4762a1bSJed Brown 
4649371c9d4SSatish Balay       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
4659371c9d4SSatish Balay       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
4669371c9d4SSatish Balay       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
4679371c9d4SSatish Balay       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
468c4762a1bSJed Brown 
469c4762a1bSJed Brown       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
470c4762a1bSJed Brown       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
471c4762a1bSJed Brown 
4729371c9d4SSatish 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);
473c4762a1bSJed Brown 
4749371c9d4SSatish Balay       hl /= 2.0;
4759371c9d4SSatish Balay       hr /= 2.0;
4769371c9d4SSatish Balay       ht /= 2.0;
4779371c9d4SSatish Balay       hb /= 2.0;
4789371c9d4SSatish Balay       hbr /= 2.0;
4799371c9d4SSatish Balay       htl /= 2.0;
4809371c9d4SSatish Balay       hc /= 2.0;
481c4762a1bSJed Brown 
4829371c9d4SSatish Balay       row.j = j;
4839371c9d4SSatish Balay       row.i = i;
484c4762a1bSJed Brown       k     = 0;
485c4762a1bSJed Brown       if (j > 0) {
486c4762a1bSJed Brown         v[k]     = hb;
4879371c9d4SSatish Balay         col[k].j = j - 1;
4889371c9d4SSatish Balay         col[k].i = i;
489c4762a1bSJed Brown         k++;
490c4762a1bSJed Brown       }
491c4762a1bSJed Brown 
492c4762a1bSJed Brown       if (j > 0 && i < mx - 1) {
493c4762a1bSJed Brown         v[k]     = hbr;
4949371c9d4SSatish Balay         col[k].j = j - 1;
4959371c9d4SSatish Balay         col[k].i = i + 1;
496c4762a1bSJed Brown         k++;
497c4762a1bSJed Brown       }
498c4762a1bSJed Brown 
499c4762a1bSJed Brown       if (i > 0) {
500c4762a1bSJed Brown         v[k]     = hl;
5019371c9d4SSatish Balay         col[k].j = j;
5029371c9d4SSatish Balay         col[k].i = i - 1;
503c4762a1bSJed Brown         k++;
504c4762a1bSJed Brown       }
505c4762a1bSJed Brown 
506c4762a1bSJed Brown       v[k]     = hc;
5079371c9d4SSatish Balay       col[k].j = j;
5089371c9d4SSatish Balay       col[k].i = i;
509c4762a1bSJed Brown       k++;
510c4762a1bSJed Brown 
511c4762a1bSJed Brown       if (i < mx - 1) {
512c4762a1bSJed Brown         v[k]     = hr;
5139371c9d4SSatish Balay         col[k].j = j;
5149371c9d4SSatish Balay         col[k].i = i + 1;
515c4762a1bSJed Brown         k++;
516c4762a1bSJed Brown       }
517c4762a1bSJed Brown 
518c4762a1bSJed Brown       if (i > 0 && j < my - 1) {
519c4762a1bSJed Brown         v[k]     = htl;
5209371c9d4SSatish Balay         col[k].j = j + 1;
5219371c9d4SSatish Balay         col[k].i = i - 1;
522c4762a1bSJed Brown         k++;
523c4762a1bSJed Brown       }
524c4762a1bSJed Brown 
525c4762a1bSJed Brown       if (j < my - 1) {
526c4762a1bSJed Brown         v[k]     = ht;
5279371c9d4SSatish Balay         col[k].j = j + 1;
5289371c9d4SSatish Balay         col[k].i = i;
529c4762a1bSJed Brown         k++;
530c4762a1bSJed Brown       }
531c4762a1bSJed Brown 
532c4762a1bSJed Brown       /*
533c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
534c4762a1bSJed Brown          earlier, in the main routine.
535c4762a1bSJed Brown       */
5369566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(Hessian, 1, &row, k, col, v, INSERT_VALUES));
537c4762a1bSJed Brown     }
538c4762a1bSJed Brown   }
539c4762a1bSJed Brown 
5409566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x));
5419566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->dm, &localX));
542c4762a1bSJed Brown 
5439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY));
5449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY));
545c4762a1bSJed Brown 
5469566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0 * xm * ym));
547c4762a1bSJed Brown   PetscFunctionReturn(0);
548c4762a1bSJed Brown }
549c4762a1bSJed Brown 
550c4762a1bSJed Brown /* ------------------------------------------------------------------- */
551c4762a1bSJed Brown /*
552c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
553c4762a1bSJed Brown    the region.
554c4762a1bSJed Brown 
555c4762a1bSJed Brown    Input Parameter:
556c4762a1bSJed Brown .  user - user-defined application context
557c4762a1bSJed Brown 
558c4762a1bSJed Brown    Output Parameter:
559c4762a1bSJed Brown .  user - user-defined application context
560c4762a1bSJed Brown */
5619371c9d4SSatish Balay static PetscErrorCode MSA_BoundaryConditions(AppCtx *user) {
562c4762a1bSJed Brown   PetscInt   i, j, k, limit = 0, maxits = 5;
563c4762a1bSJed Brown   PetscInt   xs, ys, xm, ym, gxs, gys, gxm, gym;
564c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
565c4762a1bSJed Brown   PetscInt   bsize = 0, lsize = 0, tsize = 0, rsize = 0;
566c4762a1bSJed Brown   PetscReal  one = 1.0, two = 2.0, three = 3.0, tol = 1e-10;
567c4762a1bSJed Brown   PetscReal  fnorm, det, hx, hy, xt = 0, yt = 0;
568c4762a1bSJed Brown   PetscReal  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
569c4762a1bSJed Brown   PetscReal  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
570c4762a1bSJed Brown   PetscReal *boundary;
571c4762a1bSJed Brown   PetscBool  flg;
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   PetscFunctionBegin;
574c4762a1bSJed Brown   /* Get local mesh boundaries */
5759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
5769566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
577c4762a1bSJed Brown 
578c4762a1bSJed Brown   bsize = xm + 2;
579c4762a1bSJed Brown   lsize = ym + 2;
580c4762a1bSJed Brown   rsize = ym + 2;
581c4762a1bSJed Brown   tsize = xm + 2;
582c4762a1bSJed Brown 
5839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bsize, &user->bottom));
5849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsize, &user->top));
5859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lsize, &user->left));
5869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(rsize, &user->right));
587c4762a1bSJed Brown 
5889371c9d4SSatish Balay   hx = (r - l) / (mx + 1);
5899371c9d4SSatish Balay   hy = (t - b) / (my + 1);
590c4762a1bSJed Brown 
591c4762a1bSJed Brown   for (j = 0; j < 4; j++) {
592c4762a1bSJed Brown     if (j == 0) {
593c4762a1bSJed Brown       yt       = b;
594c4762a1bSJed Brown       xt       = l + hx * xs;
595c4762a1bSJed Brown       limit    = bsize;
596c4762a1bSJed Brown       boundary = user->bottom;
597c4762a1bSJed Brown     } else if (j == 1) {
598c4762a1bSJed Brown       yt       = t;
599c4762a1bSJed Brown       xt       = l + hx * xs;
600c4762a1bSJed Brown       limit    = tsize;
601c4762a1bSJed Brown       boundary = user->top;
602c4762a1bSJed Brown     } else if (j == 2) {
603c4762a1bSJed Brown       yt       = b + hy * ys;
604c4762a1bSJed Brown       xt       = l;
605c4762a1bSJed Brown       limit    = lsize;
606c4762a1bSJed Brown       boundary = user->left;
607c4762a1bSJed Brown     } else { /* if (j==3) */
608c4762a1bSJed Brown       yt       = b + hy * ys;
609c4762a1bSJed Brown       xt       = r;
610c4762a1bSJed Brown       limit    = rsize;
611c4762a1bSJed Brown       boundary = user->right;
612c4762a1bSJed Brown     }
613c4762a1bSJed Brown 
614c4762a1bSJed Brown     for (i = 0; i < limit; i++) {
615c4762a1bSJed Brown       u1 = xt;
616c4762a1bSJed Brown       u2 = -yt;
617c4762a1bSJed Brown       for (k = 0; k < maxits; k++) {
618c4762a1bSJed Brown         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
619c4762a1bSJed Brown         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
620c4762a1bSJed Brown         fnorm = PetscSqrtReal(nf1 * nf1 + nf2 * nf2);
621c4762a1bSJed Brown         if (fnorm <= tol) break;
622c4762a1bSJed Brown         njac11 = one + u2 * u2 - u1 * u1;
623c4762a1bSJed Brown         njac12 = two * u1 * u2;
624c4762a1bSJed Brown         njac21 = -two * u1 * u2;
625c4762a1bSJed Brown         njac22 = -one - u1 * u1 + u2 * u2;
626c4762a1bSJed Brown         det    = njac11 * njac22 - njac21 * njac12;
627c4762a1bSJed Brown         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
628c4762a1bSJed Brown         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
629c4762a1bSJed Brown       }
630c4762a1bSJed Brown 
631c4762a1bSJed Brown       boundary[i] = u1 * u1 - u2 * u2;
632c4762a1bSJed Brown       if (j == 0 || j == 1) {
633c4762a1bSJed Brown         xt = xt + hx;
634c4762a1bSJed Brown       } else { /*  if (j==2 || j==3) */
635c4762a1bSJed Brown         yt = yt + hy;
636c4762a1bSJed Brown       }
637c4762a1bSJed Brown     }
638c4762a1bSJed Brown   }
639c4762a1bSJed Brown 
640c4762a1bSJed Brown   /* Scale the boundary if desired */
641c4762a1bSJed Brown   if (1 == 1) {
642c4762a1bSJed Brown     PetscReal scl = 1.0;
643c4762a1bSJed Brown 
6449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg));
645c4762a1bSJed Brown     if (flg) {
646c4762a1bSJed Brown       for (i = 0; i < bsize; i++) user->bottom[i] *= scl;
647c4762a1bSJed Brown     }
648c4762a1bSJed Brown 
6499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg));
650c4762a1bSJed Brown     if (flg) {
651c4762a1bSJed Brown       for (i = 0; i < tsize; i++) user->top[i] *= scl;
652c4762a1bSJed Brown     }
653c4762a1bSJed Brown 
6549566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg));
655c4762a1bSJed Brown     if (flg) {
656c4762a1bSJed Brown       for (i = 0; i < rsize; i++) user->right[i] *= scl;
657c4762a1bSJed Brown     }
658c4762a1bSJed Brown 
6599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg));
660c4762a1bSJed Brown     if (flg) {
661c4762a1bSJed Brown       for (i = 0; i < lsize; i++) user->left[i] *= scl;
662c4762a1bSJed Brown     }
663c4762a1bSJed Brown   }
664c4762a1bSJed Brown   PetscFunctionReturn(0);
665c4762a1bSJed Brown }
666c4762a1bSJed Brown 
667c4762a1bSJed Brown /* ------------------------------------------------------------------- */
668c4762a1bSJed Brown /*
669c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
670c4762a1bSJed Brown 
671c4762a1bSJed Brown    Input Parameters:
672c4762a1bSJed Brown .  user - user-defined application context
673c4762a1bSJed Brown .  X - vector for initial guess
674c4762a1bSJed Brown 
675c4762a1bSJed Brown    Output Parameters:
676c4762a1bSJed Brown .  X - newly computed initial guess
677c4762a1bSJed Brown */
6789371c9d4SSatish Balay static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X) {
679c4762a1bSJed Brown   PetscInt  start2 = -1, i, j;
680c4762a1bSJed Brown   PetscReal start1 = 0;
681c4762a1bSJed Brown   PetscBool flg1, flg2;
682c4762a1bSJed Brown 
683c4762a1bSJed Brown   PetscFunctionBegin;
6849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-start", &start1, &flg1));
6859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-random", &start2, &flg2));
686c4762a1bSJed Brown 
687c4762a1bSJed Brown   if (flg1) { /* The zero vector is reasonable */
688c4762a1bSJed Brown 
6899566063dSJacob Faibussowitsch     PetscCall(VecSet(X, start1));
690c4762a1bSJed Brown 
691c4762a1bSJed Brown   } else if (flg2 && start2 > 0) { /* Try a random start between -0.5 and 0.5 */
6929371c9d4SSatish Balay     PetscRandom rctx;
6939371c9d4SSatish Balay     PetscReal   np5 = -0.5;
694c4762a1bSJed Brown 
6959566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
6969566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(X, rctx));
6979566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rctx));
6989566063dSJacob Faibussowitsch     PetscCall(VecShift(X, np5));
699c4762a1bSJed Brown 
700c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
701c4762a1bSJed Brown     PetscInt    xs, xm, ys, ym;
702c4762a1bSJed Brown     PetscInt    mx = user->mx, my = user->my;
703c4762a1bSJed Brown     PetscReal **x;
704c4762a1bSJed Brown 
705c4762a1bSJed Brown     /* Get local mesh boundaries */
7069566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
707c4762a1bSJed Brown 
708c4762a1bSJed Brown     /* Get pointers to vector data */
7099566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(user->dm, X, (void **)&x));
710c4762a1bSJed Brown 
711c4762a1bSJed Brown     /* Perform local computations */
712c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
7139371c9d4SSatish Balay       for (i = xs; i < xs + xm; i++) { x[j][i] = (((j + 1) * user->bottom[i - xs + 1] + (my - j + 1) * user->top[i - xs + 1]) / (my + 2) + ((i + 1) * user->left[j - ys + 1] + (mx - i + 1) * user->right[j - ys + 1]) / (mx + 2)) / 2.0; }
714c4762a1bSJed Brown     }
7159566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(user->dm, X, (void **)&x));
7169566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(9.0 * xm * ym));
717c4762a1bSJed Brown   }
718c4762a1bSJed Brown   PetscFunctionReturn(0);
719c4762a1bSJed Brown }
720c4762a1bSJed Brown 
721c4762a1bSJed Brown /*-----------------------------------------------------------------------*/
7229371c9d4SSatish Balay PetscErrorCode My_Monitor(Tao tao, void *ctx) {
723c4762a1bSJed Brown   Vec X;
724c4762a1bSJed Brown 
725c4762a1bSJed Brown   PetscFunctionBegin;
7269566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao, &X));
7279566063dSJacob Faibussowitsch   PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
728c4762a1bSJed Brown   PetscFunctionReturn(0);
729c4762a1bSJed Brown }
730c4762a1bSJed Brown 
731c4762a1bSJed Brown /*TEST
732c4762a1bSJed Brown 
733c4762a1bSJed Brown    build:
734c4762a1bSJed Brown       requires: !complex
735c4762a1bSJed Brown 
736c4762a1bSJed Brown    test:
737c4762a1bSJed Brown       args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3
738c4762a1bSJed Brown       requires: !single
739c4762a1bSJed Brown 
740c4762a1bSJed Brown    test:
741c4762a1bSJed Brown       suffix: 2
742c4762a1bSJed Brown       nsize: 2
743c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4
744c4762a1bSJed Brown       filter: grep -v "nls ksp"
745c4762a1bSJed Brown       requires: !single
746c4762a1bSJed Brown 
747c4762a1bSJed Brown    test:
748c4762a1bSJed Brown       suffix: 3
749c4762a1bSJed Brown       nsize: 3
750c4762a1bSJed Brown       args: -tao_smonitor -tao_type cg -tao_cg_type fr -mx 10 -my 10 -tao_gatol 1.e-3
751c4762a1bSJed Brown       requires: !single
752c4762a1bSJed Brown 
753c4762a1bSJed Brown    test:
754c4762a1bSJed Brown       suffix: 5
755c4762a1bSJed Brown       nsize: 2
756c4762a1bSJed Brown       args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3
757c4762a1bSJed Brown       requires: !single
758c4762a1bSJed Brown 
759c4762a1bSJed Brown TEST*/
760