xref: /petsc/src/tao/unconstrained/tutorials/eptorsion1.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 eptorsion1 [-help] [all TAO options] */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /* ----------------------------------------------------------------------
4c4762a1bSJed Brown 
5c4762a1bSJed Brown   Elastic-plastic torsion problem.
6c4762a1bSJed Brown 
7c4762a1bSJed Brown   The elastic plastic torsion problem arises from the determination
8c4762a1bSJed Brown   of the stress field on an infinitely long cylindrical bar, which is
9c4762a1bSJed Brown   equivalent to the solution of the following problem:
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)}
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   where C is the torsion angle per unit length.
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   The multiprocessor version of this code is eptorsion2.c.
16c4762a1bSJed Brown 
17c4762a1bSJed Brown ---------------------------------------------------------------------- */
18c4762a1bSJed Brown 
19c4762a1bSJed Brown /*
20c4762a1bSJed Brown   Include "petsctao.h" so that we can use TAO solvers.  Note that this
21c4762a1bSJed Brown   file automatically includes files for lower-level support, such as those
22c4762a1bSJed Brown   provided by the PETSc library:
23c4762a1bSJed Brown      petsc.h       - base PETSc routines   petscvec.h - vectors
24a5b23f4aSJose E. Roman      petscsys.h    - system routines        petscmat.h - matrices
25c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
26c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
27c4762a1bSJed Brown */
28c4762a1bSJed Brown 
29c4762a1bSJed Brown #include <petsctao.h>
30c4762a1bSJed Brown 
319371c9d4SSatish Balay static char help[] = "Demonstrates use of the TAO package to solve \n\
32c4762a1bSJed Brown unconstrained minimization problems on a single processor.  This example \n\
33c4762a1bSJed Brown is based on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 \n\
34c4762a1bSJed Brown test suite.\n\
35c4762a1bSJed Brown The command line options are:\n\
36c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
37c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
38c4762a1bSJed Brown   -par <param>, where <param> = angle of twist per unit length\n\n";
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /*
41c4762a1bSJed Brown    User-defined application context - contains data needed by the
42c4762a1bSJed Brown    application-provided call-back routines, FormFunction(),
43c4762a1bSJed Brown    FormGradient(), and FormHessian().
44c4762a1bSJed Brown */
45c4762a1bSJed Brown 
46c4762a1bSJed Brown typedef struct {
47c4762a1bSJed Brown   PetscReal param;      /* nonlinearity parameter */
48c4762a1bSJed Brown   PetscInt  mx, my;     /* discretization in x- and y-directions */
49c4762a1bSJed Brown   PetscInt  ndim;       /* problem dimension */
50c4762a1bSJed Brown   Vec       s, y, xvec; /* work space for computing Hessian */
51c4762a1bSJed Brown   PetscReal hx, hy;     /* mesh spacing in x- and y-directions */
52c4762a1bSJed Brown } AppCtx;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /* -------- User-defined Routines --------- */
55c4762a1bSJed Brown 
56c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *, Vec);
57c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
58c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
59c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
60c4762a1bSJed Brown PetscErrorCode HessianProductMat(Mat, Vec, Vec);
61c4762a1bSJed Brown PetscErrorCode HessianProduct(void *, Vec, Vec);
62c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *);
63c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
64c4762a1bSJed Brown 
659371c9d4SSatish Balay PetscErrorCode main(int argc, char **argv) {
66c4762a1bSJed Brown   PetscInt    mx = 10; /* discretization in x-direction */
67c4762a1bSJed Brown   PetscInt    my = 10; /* discretization in y-direction */
68c4762a1bSJed Brown   Vec         x;       /* solution, gradient vectors */
69c4762a1bSJed Brown   PetscBool   flg;     /* A return value when checking for use options */
70c4762a1bSJed Brown   Tao         tao;     /* Tao solver context */
71c4762a1bSJed Brown   Mat         H;       /* Hessian matrix */
72c4762a1bSJed Brown   AppCtx      user;    /* application context */
73c4762a1bSJed Brown   PetscMPIInt size;    /* number of processes */
74c4762a1bSJed Brown   PetscReal   one = 1.0;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   PetscBool test_lmvm = PETSC_FALSE;
77c4762a1bSJed Brown   KSP       ksp;
78c4762a1bSJed Brown   PC        pc;
79c4762a1bSJed Brown   Mat       M;
80c4762a1bSJed Brown   Vec       in, out, out2;
81c4762a1bSJed Brown   PetscReal mult_solve_dist;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Initialize TAO,PETSc */
84327415f7SBarry Smith   PetscFunctionBeginUser;
859566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
873c859ba3SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors");
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* Specify default parameters for the problem, check for command-line overrides */
90c4762a1bSJed Brown   user.param = 5.0;
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, &flg));
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, &flg));
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, &flg));
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg));
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Elastic-Plastic Torsion Problem -----\n"));
9763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n", mx, my));
989371c9d4SSatish Balay   user.ndim = mx * my;
999371c9d4SSatish Balay   user.mx   = mx;
1009371c9d4SSatish Balay   user.my   = my;
1019371c9d4SSatish Balay   user.hx   = one / (mx + 1);
1029371c9d4SSatish Balay   user.hy   = one / (my + 1);
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Allocate vectors */
1059566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.ndim, &user.y));
1069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.y, &user.s));
1079566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.y, &x));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1109566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
1119566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOLMVM));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* Set solution vector with an initial guess */
1149566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user, x));
1159566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
1189566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* From command line options, determine if using matrix-free hessian */
1219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-my_tao_mf", &flg));
122c4762a1bSJed Brown   if (flg) {
1239566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PETSC_COMM_SELF, user.ndim, user.ndim, user.ndim, user.ndim, (void *)&user, &H));
1249566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(H, MATOP_MULT, (void (*)(void))HessianProductMat));
1259566063dSJacob Faibussowitsch     PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE));
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, H, H, MatrixFreeHessian, (void *)&user));
128c4762a1bSJed Brown   } else {
1299566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, user.ndim, user.ndim, 5, NULL, &H));
1309566063dSJacob Faibussowitsch     PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE));
1319566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, H, H, FormHessian, (void *)&user));
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Test the LMVM matrix */
135c4762a1bSJed Brown   if (test_lmvm) {
1369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bntr"));
1379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsSetValue(NULL, "-tao_bnk_pc_type", "lmvm"));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* Check for any TAO command line options */
1419566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1449566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Test the LMVM matrix */
147c4762a1bSJed Brown   if (test_lmvm) {
1489566063dSJacob Faibussowitsch     PetscCall(TaoGetKSP(tao, &ksp));
1499566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
1509566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(pc, &M));
1519566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &in));
1529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &out));
1539566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &out2));
1549566063dSJacob Faibussowitsch     PetscCall(VecSet(in, 5.0));
1559566063dSJacob Faibussowitsch     PetscCall(MatMult(M, in, out));
1569566063dSJacob Faibussowitsch     PetscCall(MatSolve(M, out, out2));
1579566063dSJacob Faibussowitsch     PetscCall(VecAXPY(out2, -1.0, in));
1589566063dSJacob Faibussowitsch     PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist));
15963a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between MatMult and MatSolve: %e\n", (double)mult_solve_dist));
1609566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&in));
1619566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&out));
1629566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&out2));
163c4762a1bSJed Brown   }
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.s));
1679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.y));
1689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&H));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
172b122ec5aSJacob Faibussowitsch   return 0;
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown /* ------------------------------------------------------------------- */
176c4762a1bSJed Brown /*
177c4762a1bSJed Brown     FormInitialGuess - Computes an initial approximation to the solution.
178c4762a1bSJed Brown 
179c4762a1bSJed Brown     Input Parameters:
180c4762a1bSJed Brown .   user - user-defined application context
181c4762a1bSJed Brown .   X    - vector
182c4762a1bSJed Brown 
183c4762a1bSJed Brown     Output Parameters:
184c4762a1bSJed Brown .   X    - vector
185c4762a1bSJed Brown */
1869371c9d4SSatish Balay PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) {
187c4762a1bSJed Brown   PetscReal hx = user->hx, hy = user->hy, temp;
188c4762a1bSJed Brown   PetscReal val;
189c4762a1bSJed Brown   PetscInt  i, j, k, nx = user->mx, ny = user->my;
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* Compute initial guess */
192c4762a1bSJed Brown   PetscFunctionBeginUser;
193c4762a1bSJed Brown   for (j = 0; j < ny; j++) {
194c4762a1bSJed Brown     temp = PetscMin(j + 1, ny - j) * hy;
195c4762a1bSJed Brown     for (i = 0; i < nx; i++) {
196c4762a1bSJed Brown       k   = nx * j + i;
197c4762a1bSJed Brown       val = PetscMin((PetscMin(i + 1, nx - i)) * hx, temp);
1989566063dSJacob Faibussowitsch       PetscCall(VecSetValues(X, 1, &k, &val, ADD_VALUES));
199c4762a1bSJed Brown     }
200c4762a1bSJed Brown   }
2019566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
2029566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
203c4762a1bSJed Brown   PetscFunctionReturn(0);
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /* ------------------------------------------------------------------- */
207c4762a1bSJed Brown /*
208c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
209c4762a1bSJed Brown 
210c4762a1bSJed Brown    Input Parameters:
211c4762a1bSJed Brown    tao - the Tao context
212c4762a1bSJed Brown    X   - the input vector
213c4762a1bSJed Brown    ptr - optional user-defined context, as set by TaoSetFunction()
214c4762a1bSJed Brown 
215c4762a1bSJed Brown    Output Parameters:
216c4762a1bSJed Brown    f   - the newly evaluated function
217c4762a1bSJed Brown    G   - the newly evaluated gradient
218c4762a1bSJed Brown */
2199371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) {
220c4762a1bSJed Brown   PetscFunctionBeginUser;
2219566063dSJacob Faibussowitsch   PetscCall(FormFunction(tao, X, f, ptr));
2229566063dSJacob Faibussowitsch   PetscCall(FormGradient(tao, X, G, ptr));
223c4762a1bSJed Brown   PetscFunctionReturn(0);
224c4762a1bSJed Brown }
225c4762a1bSJed Brown 
226c4762a1bSJed Brown /* ------------------------------------------------------------------- */
227c4762a1bSJed Brown /*
228c4762a1bSJed Brown    FormFunction - Evaluates the function, f(X).
229c4762a1bSJed Brown 
230c4762a1bSJed Brown    Input Parameters:
231c4762a1bSJed Brown .  tao - the Tao context
232c4762a1bSJed Brown .  X   - the input vector
233c4762a1bSJed Brown .  ptr - optional user-defined context, as set by TaoSetFunction()
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    Output Parameters:
236c4762a1bSJed Brown .  f    - the newly evaluated function
237c4762a1bSJed Brown */
2389371c9d4SSatish Balay PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) {
239c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
240c4762a1bSJed Brown   PetscReal          hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5;
241c4762a1bSJed Brown   PetscReal          zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0;
242c4762a1bSJed Brown   PetscReal          v, cdiv3 = user->param / three;
243c4762a1bSJed Brown   const PetscScalar *x;
244c4762a1bSJed Brown   PetscInt           nx = user->mx, ny = user->my, i, j, k;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   PetscFunctionBeginUser;
247c4762a1bSJed Brown   /* Get pointer to vector data */
2489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   /* Compute function contributions over the lower triangular elements */
251c4762a1bSJed Brown   for (j = -1; j < ny; j++) {
252c4762a1bSJed Brown     for (i = -1; i < nx; i++) {
253c4762a1bSJed Brown       k  = nx * j + i;
254c4762a1bSJed Brown       v  = zero;
255c4762a1bSJed Brown       vr = zero;
256c4762a1bSJed Brown       vt = zero;
257c4762a1bSJed Brown       if (i >= 0 && j >= 0) v = x[k];
258c4762a1bSJed Brown       if (i < nx - 1 && j > -1) vr = x[k + 1];
259c4762a1bSJed Brown       if (i > -1 && j < ny - 1) vt = x[k + nx];
260c4762a1bSJed Brown       dvdx = (vr - v) / hx;
261c4762a1bSJed Brown       dvdy = (vt - v) / hy;
262c4762a1bSJed Brown       fquad += dvdx * dvdx + dvdy * dvdy;
263c4762a1bSJed Brown       flin -= cdiv3 * (v + vr + vt);
264c4762a1bSJed Brown     }
265c4762a1bSJed Brown   }
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Compute function contributions over the upper triangular elements */
268c4762a1bSJed Brown   for (j = 0; j <= ny; j++) {
269c4762a1bSJed Brown     for (i = 0; i <= nx; i++) {
270c4762a1bSJed Brown       k  = nx * j + i;
271c4762a1bSJed Brown       vb = zero;
272c4762a1bSJed Brown       vl = zero;
273c4762a1bSJed Brown       v  = zero;
274c4762a1bSJed Brown       if (i < nx && j > 0) vb = x[k - nx];
275c4762a1bSJed Brown       if (i > 0 && j < ny) vl = x[k - 1];
276c4762a1bSJed Brown       if (i < nx && j < ny) v = x[k];
277c4762a1bSJed Brown       dvdx  = (v - vl) / hx;
278c4762a1bSJed Brown       dvdy  = (v - vb) / hy;
279c4762a1bSJed Brown       fquad = fquad + dvdx * dvdx + dvdy * dvdy;
280c4762a1bSJed Brown       flin  = flin - cdiv3 * (vb + vl + v);
281c4762a1bSJed Brown     }
282c4762a1bSJed Brown   }
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /* Restore vector */
2859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /* Scale the function */
288c4762a1bSJed Brown   area = p5 * hx * hy;
289c4762a1bSJed Brown   *f   = area * (p5 * fquad + flin);
290c4762a1bSJed Brown 
2919566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(24.0 * nx * ny));
292c4762a1bSJed Brown   PetscFunctionReturn(0);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown /* ------------------------------------------------------------------- */
296c4762a1bSJed Brown /*
297c4762a1bSJed Brown     FormGradient - Evaluates the gradient, G(X).
298c4762a1bSJed Brown 
299c4762a1bSJed Brown     Input Parameters:
300c4762a1bSJed Brown .   tao  - the Tao context
301c4762a1bSJed Brown .   X    - input vector
302c4762a1bSJed Brown .   ptr  - optional user-defined context
303c4762a1bSJed Brown 
304c4762a1bSJed Brown     Output Parameters:
305c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
306c4762a1bSJed Brown */
3079371c9d4SSatish Balay PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) {
308c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
309c4762a1bSJed Brown   PetscReal          zero = 0.0, p5 = 0.5, three = 3.0, area, val;
310c4762a1bSJed Brown   PetscInt           nx = user->mx, ny = user->my, ind, i, j, k;
311c4762a1bSJed Brown   PetscReal          hx = user->hx, hy = user->hy;
312c4762a1bSJed Brown   PetscReal          vb, vl, vr, vt, dvdx, dvdy;
313c4762a1bSJed Brown   PetscReal          v, cdiv3 = user->param / three;
314c4762a1bSJed Brown   const PetscScalar *x;
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   PetscFunctionBeginUser;
317c4762a1bSJed Brown   /* Initialize gradient to zero */
3189566063dSJacob Faibussowitsch   PetscCall(VecSet(G, zero));
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   /* Get pointer to vector data */
3219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /* Compute gradient contributions over the lower triangular elements */
324c4762a1bSJed Brown   for (j = -1; j < ny; j++) {
325c4762a1bSJed Brown     for (i = -1; i < nx; i++) {
326c4762a1bSJed Brown       k  = nx * j + i;
327c4762a1bSJed Brown       v  = zero;
328c4762a1bSJed Brown       vr = zero;
329c4762a1bSJed Brown       vt = zero;
330c4762a1bSJed Brown       if (i >= 0 && j >= 0) v = x[k];
331c4762a1bSJed Brown       if (i < nx - 1 && j > -1) vr = x[k + 1];
332c4762a1bSJed Brown       if (i > -1 && j < ny - 1) vt = x[k + nx];
333c4762a1bSJed Brown       dvdx = (vr - v) / hx;
334c4762a1bSJed Brown       dvdy = (vt - v) / hy;
335c4762a1bSJed Brown       if (i != -1 && j != -1) {
3369371c9d4SSatish Balay         ind = k;
3379371c9d4SSatish Balay         val = -dvdx / hx - dvdy / hy - cdiv3;
3389566063dSJacob Faibussowitsch         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
339c4762a1bSJed Brown       }
340c4762a1bSJed Brown       if (i != nx - 1 && j != -1) {
3419371c9d4SSatish Balay         ind = k + 1;
3429371c9d4SSatish Balay         val = dvdx / hx - cdiv3;
3439566063dSJacob Faibussowitsch         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
344c4762a1bSJed Brown       }
345c4762a1bSJed Brown       if (i != -1 && j != ny - 1) {
3469371c9d4SSatish Balay         ind = k + nx;
3479371c9d4SSatish Balay         val = dvdy / hy - cdiv3;
3489566063dSJacob Faibussowitsch         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
349c4762a1bSJed Brown       }
350c4762a1bSJed Brown     }
351c4762a1bSJed Brown   }
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   /* Compute gradient contributions over the upper triangular elements */
354c4762a1bSJed Brown   for (j = 0; j <= ny; j++) {
355c4762a1bSJed Brown     for (i = 0; i <= nx; i++) {
356c4762a1bSJed Brown       k  = nx * j + i;
357c4762a1bSJed Brown       vb = zero;
358c4762a1bSJed Brown       vl = zero;
359c4762a1bSJed Brown       v  = zero;
360c4762a1bSJed Brown       if (i < nx && j > 0) vb = x[k - nx];
361c4762a1bSJed Brown       if (i > 0 && j < ny) vl = x[k - 1];
362c4762a1bSJed Brown       if (i < nx && j < ny) v = x[k];
363c4762a1bSJed Brown       dvdx = (v - vl) / hx;
364c4762a1bSJed Brown       dvdy = (v - vb) / hy;
365c4762a1bSJed Brown       if (i != nx && j != 0) {
3669371c9d4SSatish Balay         ind = k - nx;
3679371c9d4SSatish Balay         val = -dvdy / hy - cdiv3;
3689566063dSJacob Faibussowitsch         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
369c4762a1bSJed Brown       }
370c4762a1bSJed Brown       if (i != 0 && j != ny) {
3719371c9d4SSatish Balay         ind = k - 1;
3729371c9d4SSatish Balay         val = -dvdx / hx - cdiv3;
3739566063dSJacob Faibussowitsch         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
374c4762a1bSJed Brown       }
375c4762a1bSJed Brown       if (i != nx && j != ny) {
3769371c9d4SSatish Balay         ind = k;
3779371c9d4SSatish Balay         val = dvdx / hx + dvdy / hy - cdiv3;
3789566063dSJacob Faibussowitsch         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
379c4762a1bSJed Brown       }
380c4762a1bSJed Brown     }
381c4762a1bSJed Brown   }
3829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
383c4762a1bSJed Brown 
384c4762a1bSJed Brown   /* Assemble gradient vector */
3859566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(G));
3869566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(G));
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   /* Scale the gradient */
389c4762a1bSJed Brown   area = p5 * hx * hy;
3909566063dSJacob Faibussowitsch   PetscCall(VecScale(G, area));
3919566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(24.0 * nx * ny));
392c4762a1bSJed Brown   PetscFunctionReturn(0);
393c4762a1bSJed Brown }
394c4762a1bSJed Brown 
395c4762a1bSJed Brown /* ------------------------------------------------------------------- */
396c4762a1bSJed Brown /*
397c4762a1bSJed Brown    FormHessian - Forms the Hessian matrix.
398c4762a1bSJed Brown 
399c4762a1bSJed Brown    Input Parameters:
400c4762a1bSJed Brown .  tao - the Tao context
401c4762a1bSJed Brown .  X    - the input vector
402c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetHessian()
403c4762a1bSJed Brown 
404c4762a1bSJed Brown    Output Parameters:
405c4762a1bSJed Brown .  H     - Hessian matrix
406c4762a1bSJed Brown .  PrecH - optionally different preconditioning Hessian
407c4762a1bSJed Brown .  flag  - flag indicating matrix structure
408c4762a1bSJed Brown 
409c4762a1bSJed Brown    Notes:
410c4762a1bSJed Brown    This routine is intended simply as an example of the interface
411c4762a1bSJed Brown    to a Hessian evaluation routine.  Since this example compute the
412c4762a1bSJed Brown    Hessian a column at a time, it is not particularly efficient and
413c4762a1bSJed Brown    is not recommended.
414c4762a1bSJed Brown */
4159371c9d4SSatish Balay PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) {
416c4762a1bSJed Brown   AppCtx    *user = (AppCtx *)ptr;
417c4762a1bSJed Brown   PetscInt   i, j, ndim = user->ndim;
418c4762a1bSJed Brown   PetscReal *y, zero = 0.0, one = 1.0;
419c4762a1bSJed Brown   PetscBool  assembled;
420c4762a1bSJed Brown 
421c4762a1bSJed Brown   PetscFunctionBeginUser;
422c4762a1bSJed Brown   user->xvec = X;
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /* Initialize Hessian entries and work vector to zero */
4259566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H, &assembled));
4269566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(H));
427c4762a1bSJed Brown 
4289566063dSJacob Faibussowitsch   PetscCall(VecSet(user->s, zero));
429c4762a1bSJed Brown 
430c4762a1bSJed Brown   /* Loop over matrix columns to compute entries of the Hessian */
431c4762a1bSJed Brown   for (j = 0; j < ndim; j++) {
4329566063dSJacob Faibussowitsch     PetscCall(VecSetValues(user->s, 1, &j, &one, INSERT_VALUES));
4339566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(user->s));
4349566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(user->s));
435c4762a1bSJed Brown 
4369566063dSJacob Faibussowitsch     PetscCall(HessianProduct(ptr, user->s, user->y));
437c4762a1bSJed Brown 
4389566063dSJacob Faibussowitsch     PetscCall(VecSetValues(user->s, 1, &j, &zero, INSERT_VALUES));
4399566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(user->s));
4409566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(user->s));
441c4762a1bSJed Brown 
4429566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->y, &y));
443c4762a1bSJed Brown     for (i = 0; i < ndim; i++) {
444*48a46eb9SPierre Jolivet       if (y[i] != zero) PetscCall(MatSetValues(H, 1, &i, 1, &j, &y[i], ADD_VALUES));
445c4762a1bSJed Brown     }
4469566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->y, &y));
447c4762a1bSJed Brown   }
4489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
4499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
450c4762a1bSJed Brown   PetscFunctionReturn(0);
451c4762a1bSJed Brown }
452c4762a1bSJed Brown 
453c4762a1bSJed Brown /* ------------------------------------------------------------------- */
454c4762a1bSJed Brown /*
455c4762a1bSJed Brown    MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector
456c4762a1bSJed Brown    products.
457c4762a1bSJed Brown 
458c4762a1bSJed Brown    Input Parameters:
459c4762a1bSJed Brown .  tao - the Tao context
460c4762a1bSJed Brown .  X    - the input vector
461c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetHessian()
462c4762a1bSJed Brown 
463c4762a1bSJed Brown    Output Parameters:
464c4762a1bSJed Brown .  H     - Hessian matrix
465c4762a1bSJed Brown .  PrecH - optionally different preconditioning Hessian
466c4762a1bSJed Brown .  flag  - flag indicating matrix structure
467c4762a1bSJed Brown */
4689371c9d4SSatish Balay PetscErrorCode MatrixFreeHessian(Tao tao, Vec X, Mat H, Mat PrecH, void *ptr) {
469c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
470c4762a1bSJed Brown 
471c4762a1bSJed Brown   /* Sets location of vector for use in computing matrix-vector products  of the form H(X)*y  */
472362febeeSStefano Zampini   PetscFunctionBeginUser;
473c4762a1bSJed Brown   user->xvec = X;
474c4762a1bSJed Brown   PetscFunctionReturn(0);
475c4762a1bSJed Brown }
476c4762a1bSJed Brown 
477c4762a1bSJed Brown /* ------------------------------------------------------------------- */
478c4762a1bSJed Brown /*
479c4762a1bSJed Brown    HessianProductMat - Computes the matrix-vector product
480c4762a1bSJed Brown    y = mat*svec.
481c4762a1bSJed Brown 
482c4762a1bSJed Brown    Input Parameters:
483c4762a1bSJed Brown .  mat  - input matrix
484c4762a1bSJed Brown .  svec - input vector
485c4762a1bSJed Brown 
486c4762a1bSJed Brown    Output Parameters:
487c4762a1bSJed Brown .  y    - solution vector
488c4762a1bSJed Brown */
4899371c9d4SSatish Balay PetscErrorCode HessianProductMat(Mat mat, Vec svec, Vec y) {
490c4762a1bSJed Brown   void *ptr;
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   PetscFunctionBeginUser;
4939566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ptr));
4949566063dSJacob Faibussowitsch   PetscCall(HessianProduct(ptr, svec, y));
495c4762a1bSJed Brown   PetscFunctionReturn(0);
496c4762a1bSJed Brown }
497c4762a1bSJed Brown 
498c4762a1bSJed Brown /* ------------------------------------------------------------------- */
499c4762a1bSJed Brown /*
500c4762a1bSJed Brown    Hessian Product - Computes the matrix-vector product:
501c4762a1bSJed Brown    y = f''(x)*svec.
502c4762a1bSJed Brown 
5037a7aea1fSJed Brown    Input Parameters:
504c4762a1bSJed Brown .  ptr  - pointer to the user-defined context
505c4762a1bSJed Brown .  svec - input vector
506c4762a1bSJed Brown 
507c4762a1bSJed Brown    Output Parameters:
508c4762a1bSJed Brown .  y    - product vector
509c4762a1bSJed Brown */
5109371c9d4SSatish Balay PetscErrorCode HessianProduct(void *ptr, Vec svec, Vec y) {
511c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
512c4762a1bSJed Brown   PetscReal          p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area;
513c4762a1bSJed Brown   const PetscScalar *x, *s;
514c4762a1bSJed Brown   PetscReal          v, vb, vl, vr, vt, hxhx, hyhy;
515c4762a1bSJed Brown   PetscInt           nx, ny, i, j, k, ind;
516c4762a1bSJed Brown 
517c4762a1bSJed Brown   PetscFunctionBeginUser;
518c4762a1bSJed Brown   nx   = user->mx;
519c4762a1bSJed Brown   ny   = user->my;
520c4762a1bSJed Brown   hx   = user->hx;
521c4762a1bSJed Brown   hy   = user->hy;
522c4762a1bSJed Brown   hxhx = one / (hx * hx);
523c4762a1bSJed Brown   hyhy = one / (hy * hy);
524c4762a1bSJed Brown 
525c4762a1bSJed Brown   /* Get pointers to vector data */
5269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user->xvec, &x));
5279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(svec, &s));
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   /* Initialize product vector to zero */
5309566063dSJacob Faibussowitsch   PetscCall(VecSet(y, zero));
531c4762a1bSJed Brown 
532c4762a1bSJed Brown   /* Compute f''(x)*s over the lower triangular elements */
533c4762a1bSJed Brown   for (j = -1; j < ny; j++) {
534c4762a1bSJed Brown     for (i = -1; i < nx; i++) {
535c4762a1bSJed Brown       k  = nx * j + i;
536c4762a1bSJed Brown       v  = zero;
537c4762a1bSJed Brown       vr = zero;
538c4762a1bSJed Brown       vt = zero;
539c4762a1bSJed Brown       if (i != -1 && j != -1) v = s[k];
540c4762a1bSJed Brown       if (i != nx - 1 && j != -1) {
541c4762a1bSJed Brown         vr  = s[k + 1];
5429371c9d4SSatish Balay         ind = k + 1;
5439371c9d4SSatish Balay         val = hxhx * (vr - v);
5449566063dSJacob Faibussowitsch         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
545c4762a1bSJed Brown       }
546c4762a1bSJed Brown       if (i != -1 && j != ny - 1) {
547c4762a1bSJed Brown         vt  = s[k + nx];
5489371c9d4SSatish Balay         ind = k + nx;
5499371c9d4SSatish Balay         val = hyhy * (vt - v);
5509566063dSJacob Faibussowitsch         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
551c4762a1bSJed Brown       }
552c4762a1bSJed Brown       if (i != -1 && j != -1) {
5539371c9d4SSatish Balay         ind = k;
5549371c9d4SSatish Balay         val = hxhx * (v - vr) + hyhy * (v - vt);
5559566063dSJacob Faibussowitsch         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
556c4762a1bSJed Brown       }
557c4762a1bSJed Brown     }
558c4762a1bSJed Brown   }
559c4762a1bSJed Brown 
560c4762a1bSJed Brown   /* Compute f''(x)*s over the upper triangular elements */
561c4762a1bSJed Brown   for (j = 0; j <= ny; j++) {
562c4762a1bSJed Brown     for (i = 0; i <= nx; i++) {
563c4762a1bSJed Brown       k  = nx * j + i;
564c4762a1bSJed Brown       v  = zero;
565c4762a1bSJed Brown       vl = zero;
566c4762a1bSJed Brown       vb = zero;
567c4762a1bSJed Brown       if (i != nx && j != ny) v = s[k];
568c4762a1bSJed Brown       if (i != nx && j != 0) {
569c4762a1bSJed Brown         vb  = s[k - nx];
5709371c9d4SSatish Balay         ind = k - nx;
5719371c9d4SSatish Balay         val = hyhy * (vb - v);
5729566063dSJacob Faibussowitsch         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
573c4762a1bSJed Brown       }
574c4762a1bSJed Brown       if (i != 0 && j != ny) {
575c4762a1bSJed Brown         vl  = s[k - 1];
5769371c9d4SSatish Balay         ind = k - 1;
5779371c9d4SSatish Balay         val = hxhx * (vl - v);
5789566063dSJacob Faibussowitsch         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
579c4762a1bSJed Brown       }
580c4762a1bSJed Brown       if (i != nx && j != ny) {
5819371c9d4SSatish Balay         ind = k;
5829371c9d4SSatish Balay         val = hxhx * (v - vl) + hyhy * (v - vb);
5839566063dSJacob Faibussowitsch         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
584c4762a1bSJed Brown       }
585c4762a1bSJed Brown     }
586c4762a1bSJed Brown   }
587c4762a1bSJed Brown   /* Restore vector data */
5889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(svec, &s));
5899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user->xvec, &x));
590c4762a1bSJed Brown 
591c4762a1bSJed Brown   /* Assemble vector */
5929566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(y));
5939566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(y));
594c4762a1bSJed Brown 
595c4762a1bSJed Brown   /* Scale resulting vector by area */
596c4762a1bSJed Brown   area = p5 * hx * hy;
5979566063dSJacob Faibussowitsch   PetscCall(VecScale(y, area));
5989566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(18.0 * nx * ny));
599c4762a1bSJed Brown   PetscFunctionReturn(0);
600c4762a1bSJed Brown }
601c4762a1bSJed Brown 
602c4762a1bSJed Brown /*TEST
603c4762a1bSJed Brown 
604c4762a1bSJed Brown    build:
605c4762a1bSJed Brown       requires: !complex
606c4762a1bSJed Brown 
607c4762a1bSJed Brown    test:
608c4762a1bSJed Brown       suffix: 1
609c4762a1bSJed Brown       args: -tao_smonitor -tao_type ntl -tao_gatol 1.e-4
610c4762a1bSJed Brown 
611c4762a1bSJed Brown    test:
612c4762a1bSJed Brown       suffix: 2
613c4762a1bSJed Brown       args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4
614c4762a1bSJed Brown 
615c4762a1bSJed Brown    test:
616c4762a1bSJed Brown       suffix: 3
617c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4 -my_tao_mf -tao_test_hessian
618c4762a1bSJed Brown 
619c4762a1bSJed Brown    test:
620c4762a1bSJed Brown      suffix: 4
621c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnls
622c4762a1bSJed Brown 
623c4762a1bSJed Brown    test:
624c4762a1bSJed Brown      suffix: 5
625c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-3 -tao_type blmvm
626c4762a1bSJed Brown 
627c4762a1bSJed Brown    test:
628c4762a1bSJed Brown      suffix: 6
629c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
630c4762a1bSJed Brown 
631c4762a1bSJed Brown TEST*/
632