xref: /petsc/src/tao/unconstrained/tutorials/minsurf1.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 minsurf1 [-help] [all TAO options] */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*  Include "petsctao.h" so we can use TAO solvers.  */
4c4762a1bSJed Brown #include <petsctao.h>
5c4762a1bSJed Brown 
69371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to \n\
7c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
8c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and \n\
9c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
10c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
11c4762a1bSJed Brown The command line options are:\n\
12c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
13c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
14c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";
15c4762a1bSJed Brown 
16c4762a1bSJed Brown /*
17c4762a1bSJed Brown    User-defined application context - contains data needed by the
18c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient()
19c4762a1bSJed Brown    and FormHessian().
20c4762a1bSJed Brown */
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscInt   mx, my;                      /* discretization in x, y directions */
23c4762a1bSJed Brown   PetscReal *bottom, *top, *left, *right; /* boundary values */
24c4762a1bSJed Brown   Mat        H;
25c4762a1bSJed Brown } AppCtx;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /* -------- User-defined Routines --------- */
28c4762a1bSJed Brown 
29c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
30c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
31c4762a1bSJed Brown static PetscErrorCode QuadraticH(AppCtx *, Vec, Mat);
32c4762a1bSJed Brown PetscErrorCode        FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
33c4762a1bSJed Brown PetscErrorCode        FormHessian(Tao, Vec, Mat, Mat, void *);
34c4762a1bSJed Brown 
35d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
36d71ae5a4SJacob Faibussowitsch {
37c4762a1bSJed Brown   PetscInt    N;    /* Size of vector */
38c4762a1bSJed Brown   PetscMPIInt size; /* Number of processors */
39c4762a1bSJed Brown   Vec         x;    /* solution, gradient vectors */
40c4762a1bSJed Brown   KSP         ksp;  /*  PETSc Krylov subspace method */
41c4762a1bSJed Brown   PetscBool   flg;  /* A return value when checking for user options */
42c4762a1bSJed Brown   Tao         tao;  /* Tao solver context */
43c4762a1bSJed Brown   AppCtx      user; /* user-defined work context */
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Initialize TAO,PETSc */
46327415f7SBarry Smith   PetscFunctionBeginUser;
479566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
503c859ba3SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors");
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Specify default dimension of the problem */
539371c9d4SSatish Balay   user.mx = 4;
549371c9d4SSatish Balay   user.my = 4;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg));
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Minimum Surface Area Problem -----\n"));
6163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n", user.mx, user.my));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* Calculate any derived values from parameters */
64c4762a1bSJed Brown   N = user.mx * user.my;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* Create TAO solver and set desired solution method  */
679566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
689566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOLMVM));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* Initialize minsurf application data structure for use in the function evaluations  */
719566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user)); /* Application specific routine */
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /*
74c4762a1bSJed Brown      Create a vector to hold the variables.  Compute an initial solution.
75c4762a1bSJed Brown      Set this vector, which will be used by TAO.
76c4762a1bSJed Brown   */
779566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));
789566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user, x)); /* Application specific routine */
799566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));     /* A TAO routine                */
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Provide TAO routines for function, gradient, and Hessian evaluation */
829566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* Create a matrix data structure to store the Hessian.  This structure will be used by TAO */
85*f4f49eeaSPierre Jolivet   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 7, NULL, &user.H));
869566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
879566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* Check for any TAO command line options */
909566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Limit the number of iterations in the KSP linear solver */
939566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
941baa6e33SBarry Smith   if (ksp) PetscCall(KSPSetTolerances(ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, user.mx * user.my));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
979566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
1029566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.bottom));
1039566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.top));
1049566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.left));
1059566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.right));
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
108b122ec5aSJacob Faibussowitsch   return 0;
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111c4762a1bSJed Brown /* -------------------------------------------------------------------- */
112c4762a1bSJed Brown 
113c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates function and corresponding gradient.
114c4762a1bSJed Brown 
115c4762a1bSJed Brown     Input Parameters:
116c4762a1bSJed Brown .   tao     - the Tao context
117c4762a1bSJed Brown .   X       - input vector
118c4762a1bSJed Brown .   userCtx - optional user-defined context, as set by TaoSetFunctionGradient()
119c4762a1bSJed Brown 
120c4762a1bSJed Brown     Output Parameters:
121c4762a1bSJed Brown .   fcn     - the newly evaluated function
122c4762a1bSJed Brown .   G       - vector containing the newly evaluated gradient
123c4762a1bSJed Brown */
124d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx)
125d71ae5a4SJacob Faibussowitsch {
126c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)userCtx;
127c4762a1bSJed Brown   PetscInt           i, j, row;
128c4762a1bSJed Brown   PetscInt           mx = user->mx, my = user->my;
129c4762a1bSJed Brown   PetscReal          rhx = mx + 1, rhy = my + 1;
130c4762a1bSJed Brown   PetscReal          hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy, ft = 0;
131c4762a1bSJed Brown   PetscReal          f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
132c4762a1bSJed Brown   PetscReal          df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
133c4762a1bSJed Brown   PetscReal          zero = 0.0;
134c4762a1bSJed Brown   PetscScalar       *g;
135c4762a1bSJed Brown   const PetscScalar *x;
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   PetscFunctionBeginUser;
1389566063dSJacob Faibussowitsch   PetscCall(VecSet(G, zero));
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
144c4762a1bSJed Brown   for (j = 0; j < my; j++) {
145c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
146c4762a1bSJed Brown       row = (j)*mx + (i);
147c4762a1bSJed Brown       xc  = x[row];
148c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
149c4762a1bSJed Brown       if (i == 0) { /* left side */
150c4762a1bSJed Brown         xl  = user->left[j + 1];
151c4762a1bSJed Brown         xlt = user->left[j + 2];
152c4762a1bSJed Brown       } else {
153c4762a1bSJed Brown         xl = x[row - 1];
154c4762a1bSJed Brown       }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown       if (j == 0) { /* bottom side */
157c4762a1bSJed Brown         xb  = user->bottom[i + 1];
158c4762a1bSJed Brown         xrb = user->bottom[i + 2];
159c4762a1bSJed Brown       } else {
160c4762a1bSJed Brown         xb = x[row - mx];
161c4762a1bSJed Brown       }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown       if (i + 1 == mx) { /* right side */
164c4762a1bSJed Brown         xr  = user->right[j + 1];
165c4762a1bSJed Brown         xrb = user->right[j];
166c4762a1bSJed Brown       } else {
167c4762a1bSJed Brown         xr = x[row + 1];
168c4762a1bSJed Brown       }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown       if (j + 1 == 0 + my) { /* top side */
171c4762a1bSJed Brown         xt  = user->top[i + 1];
172c4762a1bSJed Brown         xlt = user->top[i];
173c4762a1bSJed Brown       } else {
174c4762a1bSJed Brown         xt = x[row + mx];
175c4762a1bSJed Brown       }
176c4762a1bSJed Brown 
177ad540459SPierre Jolivet       if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx];
178ad540459SPierre Jolivet       if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx];
179c4762a1bSJed Brown 
180c4762a1bSJed Brown       d1 = (xc - xl);
181c4762a1bSJed Brown       d2 = (xc - xr);
182c4762a1bSJed Brown       d3 = (xc - xt);
183c4762a1bSJed Brown       d4 = (xc - xb);
184c4762a1bSJed Brown       d5 = (xr - xrb);
185c4762a1bSJed Brown       d6 = (xrb - xb);
186c4762a1bSJed Brown       d7 = (xlt - xl);
187c4762a1bSJed Brown       d8 = (xt - xlt);
188c4762a1bSJed Brown 
189c4762a1bSJed Brown       df1dxc = d1 * hydhx;
190c4762a1bSJed Brown       df2dxc = (d1 * hydhx + d4 * hxdhy);
191c4762a1bSJed Brown       df3dxc = d3 * hxdhy;
192c4762a1bSJed Brown       df4dxc = (d2 * hydhx + d3 * hxdhy);
193c4762a1bSJed Brown       df5dxc = d2 * hydhx;
194c4762a1bSJed Brown       df6dxc = d4 * hxdhy;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown       d1 *= rhx;
197c4762a1bSJed Brown       d2 *= rhx;
198c4762a1bSJed Brown       d3 *= rhy;
199c4762a1bSJed Brown       d4 *= rhy;
200c4762a1bSJed Brown       d5 *= rhy;
201c4762a1bSJed Brown       d6 *= rhx;
202c4762a1bSJed Brown       d7 *= rhy;
203c4762a1bSJed Brown       d8 *= rhx;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
206c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
207c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
208c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
209c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
210c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
211c4762a1bSJed Brown 
212c4762a1bSJed Brown       ft = ft + (f2 + f4);
213c4762a1bSJed Brown 
214c4762a1bSJed Brown       df1dxc /= f1;
215c4762a1bSJed Brown       df2dxc /= f2;
216c4762a1bSJed Brown       df3dxc /= f3;
217c4762a1bSJed Brown       df4dxc /= f4;
218c4762a1bSJed Brown       df5dxc /= f5;
219c4762a1bSJed Brown       df6dxc /= f6;
220c4762a1bSJed Brown 
221c4762a1bSJed Brown       g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0;
222c4762a1bSJed Brown     }
223c4762a1bSJed Brown   }
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   for (j = 0; j < my; j++) { /* left side */
226c4762a1bSJed Brown     d3 = (user->left[j + 1] - user->left[j + 2]) * rhy;
227c4762a1bSJed Brown     d2 = (user->left[j + 1] - x[j * mx]) * rhx;
228c4762a1bSJed Brown     ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
229c4762a1bSJed Brown   }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   for (i = 0; i < mx; i++) { /* bottom */
232c4762a1bSJed Brown     d2 = (user->bottom[i + 1] - user->bottom[i + 2]) * rhx;
233c4762a1bSJed Brown     d3 = (user->bottom[i + 1] - x[i]) * rhy;
234c4762a1bSJed Brown     ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   for (j = 0; j < my; j++) { /* right side */
238c4762a1bSJed Brown     d1 = (x[(j + 1) * mx - 1] - user->right[j + 1]) * rhx;
239c4762a1bSJed Brown     d4 = (user->right[j] - user->right[j + 1]) * rhy;
240c4762a1bSJed Brown     ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
241c4762a1bSJed Brown   }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   for (i = 0; i < mx; i++) { /* top side */
244c4762a1bSJed Brown     d1 = (x[(my - 1) * mx + i] - user->top[i + 1]) * rhy;
245c4762a1bSJed Brown     d4 = (user->top[i + 1] - user->top[i]) * rhx;
246c4762a1bSJed Brown     ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Bottom left corner */
250c4762a1bSJed Brown   d1 = (user->left[0] - user->left[1]) * rhy;
251c4762a1bSJed Brown   d2 = (user->bottom[0] - user->bottom[1]) * rhx;
252c4762a1bSJed Brown   ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /* Top right corner */
255c4762a1bSJed Brown   d1 = (user->right[my + 1] - user->right[my]) * rhy;
256c4762a1bSJed Brown   d2 = (user->top[mx + 1] - user->top[mx]) * rhx;
257c4762a1bSJed Brown   ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   (*fcn) = ft * area;
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /* Restore vectors */
2629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
2649566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(67.0 * mx * my));
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266c4762a1bSJed Brown }
267c4762a1bSJed Brown 
268c4762a1bSJed Brown /* ------------------------------------------------------------------- */
269c4762a1bSJed Brown /*
270c4762a1bSJed Brown    FormHessian - Evaluates the Hessian matrix.
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    Input Parameters:
273c4762a1bSJed Brown .  tao  - the Tao context
274c4762a1bSJed Brown .  x    - input vector
275c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetHessian()
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    Output Parameters:
278c4762a1bSJed Brown .  H    - Hessian matrix
279c4762a1bSJed Brown .  Hpre - optionally different preconditioning matrix
280c4762a1bSJed Brown .  flg  - flag indicating matrix structure
281c4762a1bSJed Brown 
282c4762a1bSJed Brown */
283d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
284d71ae5a4SJacob Faibussowitsch {
285c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   PetscFunctionBeginUser;
288c4762a1bSJed Brown   /* Evaluate the Hessian entries*/
2899566063dSJacob Faibussowitsch   PetscCall(QuadraticH(user, X, H));
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
291c4762a1bSJed Brown }
292c4762a1bSJed Brown 
293c4762a1bSJed Brown /* ------------------------------------------------------------------- */
294c4762a1bSJed Brown /*
295c4762a1bSJed Brown    QuadraticH - Evaluates the Hessian matrix.
296c4762a1bSJed Brown 
297c4762a1bSJed Brown    Input Parameters:
298c4762a1bSJed Brown .  user - user-defined context, as set by TaoSetHessian()
299c4762a1bSJed Brown .  X    - input vector
300c4762a1bSJed Brown 
301c4762a1bSJed Brown    Output Parameter:
302c4762a1bSJed Brown .  H    - Hessian matrix
303c4762a1bSJed Brown */
304d71ae5a4SJacob Faibussowitsch PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian)
305d71ae5a4SJacob Faibussowitsch {
306c4762a1bSJed Brown   PetscInt           i, j, k, row;
307c4762a1bSJed Brown   PetscInt           mx = user->mx, my = user->my;
308c4762a1bSJed Brown   PetscInt           col[7];
309c4762a1bSJed Brown   PetscReal          hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
310c4762a1bSJed Brown   PetscReal          rhx = mx + 1, rhy = my + 1;
311c4762a1bSJed Brown   PetscReal          f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
312c4762a1bSJed Brown   PetscReal          hl, hr, ht, hb, hc, htl, hbr;
313c4762a1bSJed Brown   const PetscScalar *x;
314c4762a1bSJed Brown   PetscReal          v[7];
315c4762a1bSJed Brown 
316362febeeSStefano Zampini   PetscFunctionBeginUser;
317c4762a1bSJed Brown   /* Get pointers to vector data */
3189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   /* Initialize matrix entries to zero */
3219566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Hessian));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /* Set various matrix options */
3249566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
327c4762a1bSJed Brown   for (i = 0; i < mx; i++) {
328c4762a1bSJed Brown     for (j = 0; j < my; j++) {
329c4762a1bSJed Brown       row = (j)*mx + (i);
330c4762a1bSJed Brown 
331c4762a1bSJed Brown       xc  = x[row];
332c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
333c4762a1bSJed Brown 
334c4762a1bSJed Brown       /* Left side */
335c4762a1bSJed Brown       if (i == 0) {
336c4762a1bSJed Brown         xl  = user->left[j + 1];
337c4762a1bSJed Brown         xlt = user->left[j + 2];
338c4762a1bSJed Brown       } else {
339c4762a1bSJed Brown         xl = x[row - 1];
340c4762a1bSJed Brown       }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown       if (j == 0) {
343c4762a1bSJed Brown         xb  = user->bottom[i + 1];
344c4762a1bSJed Brown         xrb = user->bottom[i + 2];
345c4762a1bSJed Brown       } else {
346c4762a1bSJed Brown         xb = x[row - mx];
347c4762a1bSJed Brown       }
348c4762a1bSJed Brown 
349c4762a1bSJed Brown       if (i + 1 == mx) {
350c4762a1bSJed Brown         xr  = user->right[j + 1];
351c4762a1bSJed Brown         xrb = user->right[j];
352c4762a1bSJed Brown       } else {
353c4762a1bSJed Brown         xr = x[row + 1];
354c4762a1bSJed Brown       }
355c4762a1bSJed Brown 
356c4762a1bSJed Brown       if (j + 1 == my) {
357c4762a1bSJed Brown         xt  = user->top[i + 1];
358c4762a1bSJed Brown         xlt = user->top[i];
359c4762a1bSJed Brown       } else {
360c4762a1bSJed Brown         xt = x[row + mx];
361c4762a1bSJed Brown       }
362c4762a1bSJed Brown 
363ad540459SPierre Jolivet       if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx];
364ad540459SPierre Jolivet       if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx];
365c4762a1bSJed Brown 
366c4762a1bSJed Brown       d1 = (xc - xl) * rhx;
367c4762a1bSJed Brown       d2 = (xc - xr) * rhx;
368c4762a1bSJed Brown       d3 = (xc - xt) * rhy;
369c4762a1bSJed Brown       d4 = (xc - xb) * rhy;
370c4762a1bSJed Brown       d5 = (xrb - xr) * rhy;
371c4762a1bSJed Brown       d6 = (xrb - xb) * rhx;
372c4762a1bSJed Brown       d7 = (xlt - xl) * rhy;
373c4762a1bSJed Brown       d8 = (xlt - xt) * rhx;
374c4762a1bSJed Brown 
375c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
376c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
377c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
378c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
379c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
380c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
381c4762a1bSJed Brown 
382c4762a1bSJed Brown       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
383c4762a1bSJed Brown       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
384c4762a1bSJed Brown       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
385c4762a1bSJed Brown       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
386c4762a1bSJed Brown 
387c4762a1bSJed Brown       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
388c4762a1bSJed Brown       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
389c4762a1bSJed Brown 
3909371c9d4SSatish 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);
391c4762a1bSJed Brown 
3929371c9d4SSatish Balay       hl *= 0.5;
3939371c9d4SSatish Balay       hr *= 0.5;
3949371c9d4SSatish Balay       ht *= 0.5;
3959371c9d4SSatish Balay       hb *= 0.5;
3969371c9d4SSatish Balay       hbr *= 0.5;
3979371c9d4SSatish Balay       htl *= 0.5;
3989371c9d4SSatish Balay       hc *= 0.5;
399c4762a1bSJed Brown 
400c4762a1bSJed Brown       k = 0;
401c4762a1bSJed Brown       if (j > 0) {
4029371c9d4SSatish Balay         v[k]   = hb;
4039371c9d4SSatish Balay         col[k] = row - mx;
4049371c9d4SSatish Balay         k++;
405c4762a1bSJed Brown       }
406c4762a1bSJed Brown 
407c4762a1bSJed Brown       if (j > 0 && i < mx - 1) {
4089371c9d4SSatish Balay         v[k]   = hbr;
4099371c9d4SSatish Balay         col[k] = row - mx + 1;
4109371c9d4SSatish Balay         k++;
411c4762a1bSJed Brown       }
412c4762a1bSJed Brown 
413c4762a1bSJed Brown       if (i > 0) {
4149371c9d4SSatish Balay         v[k]   = hl;
4159371c9d4SSatish Balay         col[k] = row - 1;
4169371c9d4SSatish Balay         k++;
417c4762a1bSJed Brown       }
418c4762a1bSJed Brown 
4199371c9d4SSatish Balay       v[k]   = hc;
4209371c9d4SSatish Balay       col[k] = row;
4219371c9d4SSatish Balay       k++;
422c4762a1bSJed Brown 
423c4762a1bSJed Brown       if (i < mx - 1) {
4249371c9d4SSatish Balay         v[k]   = hr;
4259371c9d4SSatish Balay         col[k] = row + 1;
4269371c9d4SSatish Balay         k++;
427c4762a1bSJed Brown       }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown       if (i > 0 && j < my - 1) {
4309371c9d4SSatish Balay         v[k]   = htl;
4319371c9d4SSatish Balay         col[k] = row + mx - 1;
4329371c9d4SSatish Balay         k++;
433c4762a1bSJed Brown       }
434c4762a1bSJed Brown 
435c4762a1bSJed Brown       if (j < my - 1) {
4369371c9d4SSatish Balay         v[k]   = ht;
4379371c9d4SSatish Balay         col[k] = row + mx;
4389371c9d4SSatish Balay         k++;
439c4762a1bSJed Brown       }
440c4762a1bSJed Brown 
441c4762a1bSJed Brown       /*
442c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
443c4762a1bSJed Brown          earlier, in the main routine.
444c4762a1bSJed Brown       */
4459566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Hessian, 1, &row, k, col, v, INSERT_VALUES));
446c4762a1bSJed Brown     }
447c4762a1bSJed Brown   }
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /* Restore vectors */
4509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   /* Assemble the matrix */
4539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY));
4549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY));
455c4762a1bSJed Brown 
4569566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0 * mx * my));
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
458c4762a1bSJed Brown }
459c4762a1bSJed Brown 
460c4762a1bSJed Brown /* ------------------------------------------------------------------- */
461c4762a1bSJed Brown /*
462c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
463c4762a1bSJed Brown    the region.
464c4762a1bSJed Brown 
465c4762a1bSJed Brown    Input Parameter:
466c4762a1bSJed Brown .  user - user-defined application context
467c4762a1bSJed Brown 
468c4762a1bSJed Brown    Output Parameter:
469c4762a1bSJed Brown .  user - user-defined application context
470c4762a1bSJed Brown */
471d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user)
472d71ae5a4SJacob Faibussowitsch {
473c4762a1bSJed Brown   PetscInt   i, j, k, limit = 0;
474c4762a1bSJed Brown   PetscInt   maxits = 5;
475c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
476c4762a1bSJed Brown   PetscInt   bsize = 0, lsize = 0, tsize = 0, rsize = 0;
477c4762a1bSJed Brown   PetscReal  one = 1.0, two = 2.0, three = 3.0, tol = 1e-10;
478c4762a1bSJed Brown   PetscReal  fnorm, det, hx, hy, xt = 0, yt = 0;
479c4762a1bSJed Brown   PetscReal  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
480c4762a1bSJed Brown   PetscReal  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
481c4762a1bSJed Brown   PetscReal *boundary;
482c4762a1bSJed Brown 
483c4762a1bSJed Brown   PetscFunctionBeginUser;
4849371c9d4SSatish Balay   bsize = mx + 2;
4859371c9d4SSatish Balay   lsize = my + 2;
4869371c9d4SSatish Balay   rsize = my + 2;
4879371c9d4SSatish Balay   tsize = mx + 2;
488c4762a1bSJed Brown 
4899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bsize, &user->bottom));
4909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsize, &user->top));
4919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lsize, &user->left));
4929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(rsize, &user->right));
493c4762a1bSJed Brown 
4949371c9d4SSatish Balay   hx = (r - l) / (mx + 1);
4959371c9d4SSatish Balay   hy = (t - b) / (my + 1);
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   for (j = 0; j < 4; j++) {
498c4762a1bSJed Brown     if (j == 0) {
499c4762a1bSJed Brown       yt       = b;
500c4762a1bSJed Brown       xt       = l;
501c4762a1bSJed Brown       limit    = bsize;
502c4762a1bSJed Brown       boundary = user->bottom;
503c4762a1bSJed Brown     } else if (j == 1) {
504c4762a1bSJed Brown       yt       = t;
505c4762a1bSJed Brown       xt       = l;
506c4762a1bSJed Brown       limit    = tsize;
507c4762a1bSJed Brown       boundary = user->top;
508c4762a1bSJed Brown     } else if (j == 2) {
509c4762a1bSJed Brown       yt       = b;
510c4762a1bSJed Brown       xt       = l;
511c4762a1bSJed Brown       limit    = lsize;
512c4762a1bSJed Brown       boundary = user->left;
513c4762a1bSJed Brown     } else { /*  if (j==3) */
514c4762a1bSJed Brown       yt       = b;
515c4762a1bSJed Brown       xt       = r;
516c4762a1bSJed Brown       limit    = rsize;
517c4762a1bSJed Brown       boundary = user->right;
518c4762a1bSJed Brown     }
519c4762a1bSJed Brown 
520c4762a1bSJed Brown     for (i = 0; i < limit; i++) {
521c4762a1bSJed Brown       u1 = xt;
522c4762a1bSJed Brown       u2 = -yt;
523c4762a1bSJed Brown       for (k = 0; k < maxits; k++) {
524c4762a1bSJed Brown         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
525c4762a1bSJed Brown         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
526c4762a1bSJed Brown         fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2);
527c4762a1bSJed Brown         if (fnorm <= tol) break;
528c4762a1bSJed Brown         njac11 = one + u2 * u2 - u1 * u1;
529c4762a1bSJed Brown         njac12 = two * u1 * u2;
530c4762a1bSJed Brown         njac21 = -two * u1 * u2;
531c4762a1bSJed Brown         njac22 = -one - u1 * u1 + u2 * u2;
532c4762a1bSJed Brown         det    = njac11 * njac22 - njac21 * njac12;
533c4762a1bSJed Brown         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
534c4762a1bSJed Brown         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
535c4762a1bSJed Brown       }
536c4762a1bSJed Brown 
537c4762a1bSJed Brown       boundary[i] = u1 * u1 - u2 * u2;
538c4762a1bSJed Brown       if (j == 0 || j == 1) {
539c4762a1bSJed Brown         xt = xt + hx;
540c4762a1bSJed Brown       } else { /*  if (j==2 || j==3) */
541c4762a1bSJed Brown         yt = yt + hy;
542c4762a1bSJed Brown       }
543c4762a1bSJed Brown     }
544c4762a1bSJed Brown   }
5453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
546c4762a1bSJed Brown }
547c4762a1bSJed Brown 
548c4762a1bSJed Brown /* ------------------------------------------------------------------- */
549c4762a1bSJed Brown /*
550c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
551c4762a1bSJed Brown 
552c4762a1bSJed Brown    Input Parameters:
553c4762a1bSJed Brown .  user - user-defined application context
554c4762a1bSJed Brown .  X - vector for initial guess
555c4762a1bSJed Brown 
556c4762a1bSJed Brown    Output Parameters:
557c4762a1bSJed Brown .  X - newly computed initial guess
558c4762a1bSJed Brown */
559d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X)
560d71ae5a4SJacob Faibussowitsch {
561c4762a1bSJed Brown   PetscInt  start = -1, i, j;
562c4762a1bSJed Brown   PetscReal zero  = 0.0;
563c4762a1bSJed Brown   PetscBool flg;
564c4762a1bSJed Brown 
565c4762a1bSJed Brown   PetscFunctionBeginUser;
5669566063dSJacob Faibussowitsch   PetscCall(VecSet(X, zero));
5679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg));
568c4762a1bSJed Brown 
569c4762a1bSJed Brown   if (flg && start == 0) { /* The zero vector is reasonable */
5709566063dSJacob Faibussowitsch     PetscCall(VecSet(X, zero));
571c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
572c4762a1bSJed Brown     PetscInt   row;
573c4762a1bSJed Brown     PetscInt   mx = user->mx, my = user->my;
574c4762a1bSJed Brown     PetscReal *x;
575c4762a1bSJed Brown 
576c4762a1bSJed Brown     /* Get pointers to vector data */
5779566063dSJacob Faibussowitsch     PetscCall(VecGetArray(X, &x));
578c4762a1bSJed Brown     /* Perform local computations */
579c4762a1bSJed Brown     for (j = 0; j < my; j++) {
580c4762a1bSJed Brown       for (i = 0; i < mx; i++) {
581c4762a1bSJed Brown         row    = (j)*mx + (i);
582c4762a1bSJed Brown         x[row] = (((j + 1) * user->bottom[i + 1] + (my - j + 1) * user->top[i + 1]) / (my + 2) + ((i + 1) * user->left[j + 1] + (mx - i + 1) * user->right[j + 1]) / (mx + 2)) / 2.0;
583c4762a1bSJed Brown       }
584c4762a1bSJed Brown     }
585c4762a1bSJed Brown     /* Restore vectors */
5869566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(X, &x));
587c4762a1bSJed Brown   }
5883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
589c4762a1bSJed Brown }
590c4762a1bSJed Brown 
591c4762a1bSJed Brown /*TEST
592c4762a1bSJed Brown 
593c4762a1bSJed Brown    build:
594c4762a1bSJed Brown       requires: !complex
595c4762a1bSJed Brown 
596c4762a1bSJed Brown    test:
59710978b7dSBarry Smith       args: -tao_monitor_short -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4
598c4762a1bSJed Brown       requires: !single
599c4762a1bSJed Brown 
600c4762a1bSJed Brown    test:
601c4762a1bSJed Brown       suffix: 2
60210978b7dSBarry Smith       args: -tao_monitor_short -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3
603c4762a1bSJed Brown       requires: !single
604c4762a1bSJed Brown 
605c4762a1bSJed Brown    test:
606c4762a1bSJed Brown       suffix: 3
60710978b7dSBarry Smith       args: -tao_monitor_short -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3
608c4762a1bSJed Brown       requires: !single
609c4762a1bSJed Brown 
610c4762a1bSJed Brown    test:
611c4762a1bSJed Brown       suffix: 4
61210978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4
613c4762a1bSJed Brown       requires: !single
614c4762a1bSJed Brown 
615c4762a1bSJed Brown    test:
6160f0abf79SStefano Zampini       suffix: 4_ew
61710978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntr -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4
6180f0abf79SStefano Zampini       requires: !single
6190f0abf79SStefano Zampini 
6200f0abf79SStefano Zampini    test:
621c4762a1bSJed Brown       suffix: 5
62210978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4
623c4762a1bSJed Brown       requires: !single
624c4762a1bSJed Brown 
625c4762a1bSJed Brown    test:
6260f0abf79SStefano Zampini       suffix: 5_ew
62710978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntl -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4
6280f0abf79SStefano Zampini       requires: !single
6290f0abf79SStefano Zampini 
6300f0abf79SStefano Zampini    test:
631c4762a1bSJed Brown       suffix: 6
63210978b7dSBarry Smith       args: -tao_monitor_short -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4
633c4762a1bSJed Brown       requires: !single
634c4762a1bSJed Brown 
635c4762a1bSJed Brown    test:
6360f0abf79SStefano Zampini       suffix: 6_ew
63710978b7dSBarry Smith       args: -tao_monitor_short -tao_type bnls -tao_ksp_ew -tao_bnk_ksp_ew_version 3 -mx 10 -my 8 -tao_gatol 1.e-4
6380f0abf79SStefano Zampini       requires: !single
6390f0abf79SStefano Zampini 
6400f0abf79SStefano Zampini    test:
641c4762a1bSJed Brown       suffix: 7
64210978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
643c4762a1bSJed Brown       requires: !single
644c4762a1bSJed Brown 
645c4762a1bSJed Brown    test:
646c4762a1bSJed Brown       suffix: 8
64710978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
648c4762a1bSJed Brown       requires: !single
649c4762a1bSJed Brown 
650c4762a1bSJed Brown    test:
651c4762a1bSJed Brown       suffix: 9
65210978b7dSBarry Smith       args: -tao_monitor_short -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
653c4762a1bSJed Brown       requires: !single
654c4762a1bSJed Brown 
65534ad9904SAlp Dener    test:
65634ad9904SAlp Dener       suffix: 10
65710978b7dSBarry Smith       args: -tao_monitor_short -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_mf_hessian
65834ad9904SAlp Dener 
65934ad9904SAlp Dener    test:
66034ad9904SAlp Dener       suffix: 11
66110978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
66234ad9904SAlp Dener       requires: !single
66334ad9904SAlp Dener 
66434ad9904SAlp Dener    test:
66534ad9904SAlp Dener       suffix: 12
66610978b7dSBarry Smith       args: -tao_monitor_short -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
66734ad9904SAlp Dener       requires: !single
66834ad9904SAlp Dener 
669c4762a1bSJed Brown TEST*/
670