xref: /petsc/src/tao/unconstrained/tutorials/eptorsion3.c (revision 2286efddd54511ab18e8e2adb1e023c4bf8f0b92)
1b3e6a353SBarry Smith /* Program usage: mpiexec -n 1 eptorsion1 [-help] [all TAO options] */
2b3e6a353SBarry Smith 
3b3e6a353SBarry Smith /* ----------------------------------------------------------------------
4b3e6a353SBarry Smith 
5b3e6a353SBarry Smith   Elastic-plastic torsion problem.
6b3e6a353SBarry Smith 
7b3e6a353SBarry Smith   The elastic plastic torsion problem arises from the determination
8b3e6a353SBarry Smith   of the stress field on an infinitely long cylindrical bar, which is
9b3e6a353SBarry Smith   equivalent to the solution of the following problem:
10b3e6a353SBarry Smith 
11b3e6a353SBarry Smith   min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)}
12b3e6a353SBarry Smith 
13b3e6a353SBarry Smith   where C is the torsion angle per unit length.
14b3e6a353SBarry Smith 
15b3e6a353SBarry Smith   The multiprocessor version of this code is eptorsion2.c.
16b3e6a353SBarry Smith 
17b3e6a353SBarry Smith ---------------------------------------------------------------------- */
18b3e6a353SBarry Smith 
19b3e6a353SBarry Smith /*
20b3e6a353SBarry Smith   Include "petsctao.h" so that we can use TAO solvers.  Note that this
21b3e6a353SBarry Smith   file automatically includes files for lower-level support, such as those
22b3e6a353SBarry Smith   provided by the PETSc library:
23b3e6a353SBarry Smith      petsc.h       - base PETSc routines   petscvec.h - vectors
24b3e6a353SBarry Smith      petscsys.h    - system routines        petscmat.h - matrices
25b3e6a353SBarry Smith      petscis.h     - index sets            petscksp.h - Krylov subspace methods
26b3e6a353SBarry Smith      petscviewer.h - viewers               petscpc.h  - preconditioners
27b3e6a353SBarry Smith */
28b3e6a353SBarry Smith 
29b3e6a353SBarry Smith #include <petsctao.h>
30b3e6a353SBarry Smith 
31b3e6a353SBarry Smith static char help[] = "Demonstrates use of the TAO package to solve \n\
32b3e6a353SBarry Smith unconstrained minimization problems on a single processor.  This example \n\
33b3e6a353SBarry Smith is based on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 \n\
34b3e6a353SBarry Smith test suite.\n\
35b3e6a353SBarry Smith The command line options are:\n\
36b3e6a353SBarry Smith   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
37b3e6a353SBarry Smith   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
38b3e6a353SBarry Smith   -par <param>, where <param> = angle of twist per unit length\n\n";
39b3e6a353SBarry Smith 
40b3e6a353SBarry Smith /*
41b3e6a353SBarry Smith    User-defined application context - contains data needed by the
42b3e6a353SBarry Smith    application-provided call-back routines, FormFunction(),
43b3e6a353SBarry Smith    FormGradient(), and FormHessian().
44b3e6a353SBarry Smith */
45b3e6a353SBarry Smith 
46b3e6a353SBarry Smith typedef struct {
47b3e6a353SBarry Smith   PetscReal param;      /* nonlinearity parameter */
48b3e6a353SBarry Smith   PetscInt  mx, my;     /* discretization in x- and y-directions */
49b3e6a353SBarry Smith   PetscInt  ndim;       /* problem dimension */
50b3e6a353SBarry Smith   Vec       s, y, xvec; /* work space for computing Hessian */
51b3e6a353SBarry Smith   PetscReal hx, hy;     /* mesh spacing in x- and y-directions */
52b3e6a353SBarry Smith } AppCtx;
53b3e6a353SBarry Smith 
54b3e6a353SBarry Smith /* -------- User-defined Routines --------- */
55b3e6a353SBarry Smith 
56b3e6a353SBarry Smith PetscErrorCode FormInitialGuess(AppCtx *, Vec);
57b3e6a353SBarry Smith PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
58b3e6a353SBarry Smith PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
59b3e6a353SBarry Smith PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
60b3e6a353SBarry Smith PetscErrorCode HessianProductMat(Mat, Vec, Vec);
61b3e6a353SBarry Smith PetscErrorCode HessianProduct(void *, Vec, Vec);
62b3e6a353SBarry Smith PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *);
63b3e6a353SBarry Smith PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
64b3e6a353SBarry Smith 
main(int argc,char ** argv)653ba16761SJacob Faibussowitsch int main(int argc, char **argv)
66b3e6a353SBarry Smith {
67b3e6a353SBarry Smith   PetscInt    mx = 10; /* discretization in x-direction */
68b3e6a353SBarry Smith   PetscInt    my = 10; /* discretization in y-direction */
69b3e6a353SBarry Smith   Vec         x;       /* solution, gradient vectors */
70b3e6a353SBarry Smith   PetscBool   flg;     /* A return value when checking for use options */
71b3e6a353SBarry Smith   Tao         tao;     /* Tao solver context */
72b3e6a353SBarry Smith   Mat         H;       /* Hessian matrix */
73b3e6a353SBarry Smith   AppCtx      user;    /* application context */
74b3e6a353SBarry Smith   PetscMPIInt size;    /* number of processes */
75b3e6a353SBarry Smith   PetscReal   one = 1.0;
76b3e6a353SBarry Smith 
77b3e6a353SBarry Smith   PetscBool test_lmvm = PETSC_FALSE;
78b3e6a353SBarry Smith   KSP       ksp;
79b3e6a353SBarry Smith   PC        pc;
80b3e6a353SBarry Smith   Mat       M;
81b3e6a353SBarry Smith   Vec       in, out, out2;
82b3e6a353SBarry Smith   PetscReal mult_solve_dist;
83b3e6a353SBarry Smith   Vec       ub, lb;
84b3e6a353SBarry Smith   PetscInt  i = 3;
85b3e6a353SBarry Smith 
86b3e6a353SBarry Smith   /* Initialize TAO,PETSc */
87b3e6a353SBarry Smith   PetscFunctionBeginUser;
88c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
89b3e6a353SBarry Smith   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
90b3e6a353SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors");
91b3e6a353SBarry Smith 
92b3e6a353SBarry Smith   /* Specify default parameters for the problem, check for command-line overrides */
93b3e6a353SBarry Smith   user.param = 5.0;
94b3e6a353SBarry Smith   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, &flg));
95b3e6a353SBarry Smith   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, &flg));
96b3e6a353SBarry Smith   PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, &flg));
97b3e6a353SBarry Smith   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg));
98b3e6a353SBarry Smith 
99b3e6a353SBarry Smith   PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Elastic-Plastic Torsion Problem -----\n"));
100b3e6a353SBarry Smith   PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n", mx, my));
101b3e6a353SBarry Smith   user.ndim = mx * my;
102b3e6a353SBarry Smith   user.mx   = mx;
103b3e6a353SBarry Smith   user.my   = my;
104b3e6a353SBarry Smith   user.hx   = one / (mx + 1);
105b3e6a353SBarry Smith   user.hy   = one / (my + 1);
106b3e6a353SBarry Smith 
107b3e6a353SBarry Smith   /* Allocate vectors */
108b3e6a353SBarry Smith   PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.ndim, &user.y));
109b3e6a353SBarry Smith   PetscCall(VecDuplicate(user.y, &user.s));
110b3e6a353SBarry Smith   PetscCall(VecDuplicate(user.y, &x));
111b3e6a353SBarry Smith   PetscCall(VecDuplicate(user.y, &lb));
112b3e6a353SBarry Smith   PetscCall(VecDuplicate(user.y, &ub));
113b3e6a353SBarry Smith 
114b3e6a353SBarry Smith   PetscCall(VecSet(ub, 0.1));
115b3e6a353SBarry Smith   PetscCall(VecSet(lb, -0.1));
116b3e6a353SBarry Smith   PetscCall(VecSetValue(ub, i, 0., INSERT_VALUES));
117b3e6a353SBarry Smith   PetscCall(VecSetValue(lb, i, 0., INSERT_VALUES));
118b3e6a353SBarry Smith   PetscCall(VecAssemblyBegin(ub));
119b3e6a353SBarry Smith   PetscCall(VecAssemblyBegin(lb));
120b3e6a353SBarry Smith   PetscCall(VecAssemblyEnd(ub));
121b3e6a353SBarry Smith   PetscCall(VecAssemblyEnd(lb));
122b3e6a353SBarry Smith 
123b3e6a353SBarry Smith   /* Create TAO solver and set desired solution method */
124b3e6a353SBarry Smith   PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
125b3e6a353SBarry Smith   PetscCall(TaoSetType(tao, TAOLMVM));
126b3e6a353SBarry Smith 
127b3e6a353SBarry Smith   /* Set solution vector with an initial guess */
128b3e6a353SBarry Smith   PetscCall(FormInitialGuess(&user, x));
129b3e6a353SBarry Smith   PetscCall(TaoSetSolution(tao, x));
130b3e6a353SBarry Smith   PetscCall(TaoSetVariableBounds(tao, lb, ub));
131b3e6a353SBarry Smith 
132b3e6a353SBarry Smith   /* Set routine for function and gradient evaluation */
133b3e6a353SBarry Smith   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
134b3e6a353SBarry Smith 
135b3e6a353SBarry Smith   /* From command line options, determine if using matrix-free hessian */
136b3e6a353SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-my_tao_mf", &flg));
137b3e6a353SBarry Smith   if (flg) {
138b3e6a353SBarry Smith     PetscCall(MatCreateShell(PETSC_COMM_SELF, user.ndim, user.ndim, user.ndim, user.ndim, (void *)&user, &H));
139*57d50842SBarry Smith     PetscCall(MatShellSetOperation(H, MATOP_MULT, (PetscErrorCodeFn *)HessianProductMat));
140b3e6a353SBarry Smith     PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE));
141b3e6a353SBarry Smith 
142b3e6a353SBarry Smith     PetscCall(TaoSetHessian(tao, H, H, MatrixFreeHessian, (void *)&user));
143b3e6a353SBarry Smith   } else {
144b3e6a353SBarry Smith     PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, user.ndim, user.ndim, 5, NULL, &H));
145b3e6a353SBarry Smith     PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE));
146b3e6a353SBarry Smith     PetscCall(TaoSetHessian(tao, H, H, FormHessian, (void *)&user));
147b3e6a353SBarry Smith   }
148b3e6a353SBarry Smith 
149b3e6a353SBarry Smith   /* Test the LMVM matrix */
150b3e6a353SBarry Smith   if (test_lmvm) {
151b3e6a353SBarry Smith     PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bntr"));
152b3e6a353SBarry Smith     PetscCall(PetscOptionsSetValue(NULL, "-tao_bnk_pc_type", "lmvm"));
153b3e6a353SBarry Smith   }
154b3e6a353SBarry Smith 
155b3e6a353SBarry Smith   /* Check for any TAO command line options */
156b3e6a353SBarry Smith   PetscCall(TaoSetFromOptions(tao));
157b3e6a353SBarry Smith 
158b3e6a353SBarry Smith   /* SOLVE THE APPLICATION */
159b3e6a353SBarry Smith   PetscCall(TaoSolve(tao));
160b3e6a353SBarry Smith 
161b3e6a353SBarry Smith   /* Test the LMVM matrix */
162b3e6a353SBarry Smith   if (test_lmvm) {
163b3e6a353SBarry Smith     PetscCall(TaoGetKSP(tao, &ksp));
164b3e6a353SBarry Smith     PetscCall(KSPGetPC(ksp, &pc));
165b3e6a353SBarry Smith     PetscCall(PCLMVMGetMatLMVM(pc, &M));
166b3e6a353SBarry Smith     PetscCall(VecDuplicate(x, &in));
167b3e6a353SBarry Smith     PetscCall(VecDuplicate(x, &out));
168b3e6a353SBarry Smith     PetscCall(VecDuplicate(x, &out2));
169b3e6a353SBarry Smith     PetscCall(VecSet(in, 5.0));
170b3e6a353SBarry Smith     PetscCall(MatMult(M, in, out));
171b3e6a353SBarry Smith     PetscCall(MatSolve(M, out, out2));
172b3e6a353SBarry Smith     PetscCall(VecAXPY(out2, -1.0, in));
173b3e6a353SBarry Smith     PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist));
174b3e6a353SBarry Smith     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between MatMult and MatSolve: %e\n", (double)mult_solve_dist));
175b3e6a353SBarry Smith     PetscCall(VecDestroy(&in));
176b3e6a353SBarry Smith     PetscCall(VecDestroy(&out));
177b3e6a353SBarry Smith     PetscCall(VecDestroy(&out2));
178b3e6a353SBarry Smith   }
179b3e6a353SBarry Smith 
180b3e6a353SBarry Smith   PetscCall(TaoDestroy(&tao));
181b3e6a353SBarry Smith   PetscCall(VecDestroy(&user.s));
182b3e6a353SBarry Smith   PetscCall(VecDestroy(&user.y));
183b3e6a353SBarry Smith   PetscCall(VecDestroy(&x));
184b3e6a353SBarry Smith   PetscCall(MatDestroy(&H));
185b3e6a353SBarry Smith   PetscCall(VecDestroy(&lb));
186b3e6a353SBarry Smith   PetscCall(VecDestroy(&ub));
187b3e6a353SBarry Smith 
188b3e6a353SBarry Smith   PetscCall(PetscFinalize());
189b3e6a353SBarry Smith   return 0;
190b3e6a353SBarry Smith }
191b3e6a353SBarry Smith 
192b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
193b3e6a353SBarry Smith /*
194b3e6a353SBarry Smith     FormInitialGuess - Computes an initial approximation to the solution.
195b3e6a353SBarry Smith 
196b3e6a353SBarry Smith     Input Parameters:
197b3e6a353SBarry Smith .   user - user-defined application context
198b3e6a353SBarry Smith .   X    - vector
199b3e6a353SBarry Smith 
200b3e6a353SBarry Smith     Output Parameters:
201b3e6a353SBarry Smith .   X    - vector
202b3e6a353SBarry Smith */
FormInitialGuess(AppCtx * user,Vec X)203b3e6a353SBarry Smith PetscErrorCode FormInitialGuess(AppCtx *user, Vec X)
204b3e6a353SBarry Smith {
205b3e6a353SBarry Smith   PetscReal hx = user->hx, hy = user->hy, temp;
206b3e6a353SBarry Smith   PetscReal val;
207b3e6a353SBarry Smith   PetscInt  i, j, k, nx = user->mx, ny = user->my;
208b3e6a353SBarry Smith 
209b3e6a353SBarry Smith   /* Compute initial guess */
210b3e6a353SBarry Smith   PetscFunctionBeginUser;
211b3e6a353SBarry Smith   for (j = 0; j < ny; j++) {
212b3e6a353SBarry Smith     temp = PetscMin(j + 1, ny - j) * hy;
213b3e6a353SBarry Smith     for (i = 0; i < nx; i++) {
214b3e6a353SBarry Smith       k   = nx * j + i;
215b3e6a353SBarry Smith       val = PetscMin((PetscMin(i + 1, nx - i)) * hx, temp);
216b3e6a353SBarry Smith       PetscCall(VecSetValues(X, 1, &k, &val, ADD_VALUES));
217b3e6a353SBarry Smith     }
218b3e6a353SBarry Smith   }
219b3e6a353SBarry Smith   PetscCall(VecAssemblyBegin(X));
220b3e6a353SBarry Smith   PetscCall(VecAssemblyEnd(X));
2213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
222b3e6a353SBarry Smith }
223b3e6a353SBarry Smith 
224b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
225b3e6a353SBarry Smith /*
226b3e6a353SBarry Smith    FormFunctionGradient - Evaluates the function and corresponding gradient.
227b3e6a353SBarry Smith 
228b3e6a353SBarry Smith    Input Parameters:
229b3e6a353SBarry Smith    tao - the Tao context
230b3e6a353SBarry Smith    X   - the input vector
231b3e6a353SBarry Smith    ptr - optional user-defined context, as set by TaoSetFunction()
232b3e6a353SBarry Smith 
233b3e6a353SBarry Smith    Output Parameters:
234b3e6a353SBarry Smith    f   - the newly evaluated function
235b3e6a353SBarry Smith    G   - the newly evaluated gradient
236b3e6a353SBarry Smith */
FormFunctionGradient(Tao tao,Vec X,PetscReal * f,Vec G,void * ptr)237b3e6a353SBarry Smith PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
238b3e6a353SBarry Smith {
239b3e6a353SBarry Smith   PetscFunctionBeginUser;
240b3e6a353SBarry Smith   PetscCall(FormFunction(tao, X, f, ptr));
241b3e6a353SBarry Smith   PetscCall(FormGradient(tao, X, G, ptr));
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243b3e6a353SBarry Smith }
244b3e6a353SBarry Smith 
245b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
246b3e6a353SBarry Smith /*
247b3e6a353SBarry Smith    FormFunction - Evaluates the function, f(X).
248b3e6a353SBarry Smith 
249b3e6a353SBarry Smith    Input Parameters:
250b3e6a353SBarry Smith .  tao - the Tao context
251b3e6a353SBarry Smith .  X   - the input vector
252b3e6a353SBarry Smith .  ptr - optional user-defined context, as set by TaoSetFunction()
253b3e6a353SBarry Smith 
254b3e6a353SBarry Smith    Output Parameters:
255b3e6a353SBarry Smith .  f    - the newly evaluated function
256b3e6a353SBarry Smith */
FormFunction(Tao tao,Vec X,PetscReal * f,void * ptr)257b3e6a353SBarry Smith PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr)
258b3e6a353SBarry Smith {
259b3e6a353SBarry Smith   AppCtx            *user = (AppCtx *)ptr;
260b3e6a353SBarry Smith   PetscReal          hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5;
261b3e6a353SBarry Smith   PetscReal          zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0;
262b3e6a353SBarry Smith   PetscReal          v, cdiv3 = user->param / three;
263b3e6a353SBarry Smith   const PetscScalar *x;
264b3e6a353SBarry Smith   PetscInt           nx = user->mx, ny = user->my, i, j, k;
265b3e6a353SBarry Smith 
266b3e6a353SBarry Smith   PetscFunctionBeginUser;
267b3e6a353SBarry Smith   /* Get pointer to vector data */
268b3e6a353SBarry Smith   PetscCall(VecGetArrayRead(X, &x));
269b3e6a353SBarry Smith 
270b3e6a353SBarry Smith   /* Compute function contributions over the lower triangular elements */
271b3e6a353SBarry Smith   for (j = -1; j < ny; j++) {
272b3e6a353SBarry Smith     for (i = -1; i < nx; i++) {
273b3e6a353SBarry Smith       k  = nx * j + i;
274b3e6a353SBarry Smith       v  = zero;
275b3e6a353SBarry Smith       vr = zero;
276b3e6a353SBarry Smith       vt = zero;
277b3e6a353SBarry Smith       if (i >= 0 && j >= 0) v = x[k];
278b3e6a353SBarry Smith       if (i < nx - 1 && j > -1) vr = x[k + 1];
279b3e6a353SBarry Smith       if (i > -1 && j < ny - 1) vt = x[k + nx];
280b3e6a353SBarry Smith       dvdx = (vr - v) / hx;
281b3e6a353SBarry Smith       dvdy = (vt - v) / hy;
282b3e6a353SBarry Smith       fquad += dvdx * dvdx + dvdy * dvdy;
283b3e6a353SBarry Smith       flin -= cdiv3 * (v + vr + vt);
284b3e6a353SBarry Smith     }
285b3e6a353SBarry Smith   }
286b3e6a353SBarry Smith 
287b3e6a353SBarry Smith   /* Compute function contributions over the upper triangular elements */
288b3e6a353SBarry Smith   for (j = 0; j <= ny; j++) {
289b3e6a353SBarry Smith     for (i = 0; i <= nx; i++) {
290b3e6a353SBarry Smith       k  = nx * j + i;
291b3e6a353SBarry Smith       vb = zero;
292b3e6a353SBarry Smith       vl = zero;
293b3e6a353SBarry Smith       v  = zero;
294b3e6a353SBarry Smith       if (i < nx && j > 0) vb = x[k - nx];
295b3e6a353SBarry Smith       if (i > 0 && j < ny) vl = x[k - 1];
296b3e6a353SBarry Smith       if (i < nx && j < ny) v = x[k];
297b3e6a353SBarry Smith       dvdx  = (v - vl) / hx;
298b3e6a353SBarry Smith       dvdy  = (v - vb) / hy;
299b3e6a353SBarry Smith       fquad = fquad + dvdx * dvdx + dvdy * dvdy;
300b3e6a353SBarry Smith       flin  = flin - cdiv3 * (vb + vl + v);
301b3e6a353SBarry Smith     }
302b3e6a353SBarry Smith   }
303b3e6a353SBarry Smith 
304b3e6a353SBarry Smith   /* Restore vector */
305b3e6a353SBarry Smith   PetscCall(VecRestoreArrayRead(X, &x));
306b3e6a353SBarry Smith 
307b3e6a353SBarry Smith   /* Scale the function */
308b3e6a353SBarry Smith   area = p5 * hx * hy;
309b3e6a353SBarry Smith   *f   = area * (p5 * fquad + flin);
310b3e6a353SBarry Smith 
311b3e6a353SBarry Smith   PetscCall(PetscLogFlops(24.0 * nx * ny));
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
313b3e6a353SBarry Smith }
314b3e6a353SBarry Smith 
315b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
316b3e6a353SBarry Smith /*
317b3e6a353SBarry Smith     FormGradient - Evaluates the gradient, G(X).
318b3e6a353SBarry Smith 
319b3e6a353SBarry Smith     Input Parameters:
320b3e6a353SBarry Smith .   tao  - the Tao context
321b3e6a353SBarry Smith .   X    - input vector
322b3e6a353SBarry Smith .   ptr  - optional user-defined context
323b3e6a353SBarry Smith 
324b3e6a353SBarry Smith     Output Parameters:
325b3e6a353SBarry Smith .   G - vector containing the newly evaluated gradient
326b3e6a353SBarry Smith */
FormGradient(Tao tao,Vec X,Vec G,void * ptr)327b3e6a353SBarry Smith PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr)
328b3e6a353SBarry Smith {
329b3e6a353SBarry Smith   AppCtx            *user = (AppCtx *)ptr;
330b3e6a353SBarry Smith   PetscReal          zero = 0.0, p5 = 0.5, three = 3.0, area, val;
331b3e6a353SBarry Smith   PetscInt           nx = user->mx, ny = user->my, ind, i, j, k;
332b3e6a353SBarry Smith   PetscReal          hx = user->hx, hy = user->hy;
333b3e6a353SBarry Smith   PetscReal          vb, vl, vr, vt, dvdx, dvdy;
334b3e6a353SBarry Smith   PetscReal          v, cdiv3 = user->param / three;
335b3e6a353SBarry Smith   const PetscScalar *x;
336b3e6a353SBarry Smith 
337b3e6a353SBarry Smith   PetscFunctionBeginUser;
338b3e6a353SBarry Smith   /* Initialize gradient to zero */
339b3e6a353SBarry Smith   PetscCall(VecSet(G, zero));
340b3e6a353SBarry Smith 
341b3e6a353SBarry Smith   /* Get pointer to vector data */
342b3e6a353SBarry Smith   PetscCall(VecGetArrayRead(X, &x));
343b3e6a353SBarry Smith 
344b3e6a353SBarry Smith   /* Compute gradient contributions over the lower triangular elements */
345b3e6a353SBarry Smith   for (j = -1; j < ny; j++) {
346b3e6a353SBarry Smith     for (i = -1; i < nx; i++) {
347b3e6a353SBarry Smith       k  = nx * j + i;
348b3e6a353SBarry Smith       v  = zero;
349b3e6a353SBarry Smith       vr = zero;
350b3e6a353SBarry Smith       vt = zero;
351b3e6a353SBarry Smith       if (i >= 0 && j >= 0) v = x[k];
352b3e6a353SBarry Smith       if (i < nx - 1 && j > -1) vr = x[k + 1];
353b3e6a353SBarry Smith       if (i > -1 && j < ny - 1) vt = x[k + nx];
354b3e6a353SBarry Smith       dvdx = (vr - v) / hx;
355b3e6a353SBarry Smith       dvdy = (vt - v) / hy;
356b3e6a353SBarry Smith       if (i != -1 && j != -1) {
357b3e6a353SBarry Smith         ind = k;
358b3e6a353SBarry Smith         val = -dvdx / hx - dvdy / hy - cdiv3;
359b3e6a353SBarry Smith         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
360b3e6a353SBarry Smith       }
361b3e6a353SBarry Smith       if (i != nx - 1 && j != -1) {
362b3e6a353SBarry Smith         ind = k + 1;
363b3e6a353SBarry Smith         val = dvdx / hx - cdiv3;
364b3e6a353SBarry Smith         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
365b3e6a353SBarry Smith       }
366b3e6a353SBarry Smith       if (i != -1 && j != ny - 1) {
367b3e6a353SBarry Smith         ind = k + nx;
368b3e6a353SBarry Smith         val = dvdy / hy - cdiv3;
369b3e6a353SBarry Smith         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
370b3e6a353SBarry Smith       }
371b3e6a353SBarry Smith     }
372b3e6a353SBarry Smith   }
373b3e6a353SBarry Smith 
374b3e6a353SBarry Smith   /* Compute gradient contributions over the upper triangular elements */
375b3e6a353SBarry Smith   for (j = 0; j <= ny; j++) {
376b3e6a353SBarry Smith     for (i = 0; i <= nx; i++) {
377b3e6a353SBarry Smith       k  = nx * j + i;
378b3e6a353SBarry Smith       vb = zero;
379b3e6a353SBarry Smith       vl = zero;
380b3e6a353SBarry Smith       v  = zero;
381b3e6a353SBarry Smith       if (i < nx && j > 0) vb = x[k - nx];
382b3e6a353SBarry Smith       if (i > 0 && j < ny) vl = x[k - 1];
383b3e6a353SBarry Smith       if (i < nx && j < ny) v = x[k];
384b3e6a353SBarry Smith       dvdx = (v - vl) / hx;
385b3e6a353SBarry Smith       dvdy = (v - vb) / hy;
386b3e6a353SBarry Smith       if (i != nx && j != 0) {
387b3e6a353SBarry Smith         ind = k - nx;
388b3e6a353SBarry Smith         val = -dvdy / hy - cdiv3;
389b3e6a353SBarry Smith         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
390b3e6a353SBarry Smith       }
391b3e6a353SBarry Smith       if (i != 0 && j != ny) {
392b3e6a353SBarry Smith         ind = k - 1;
393b3e6a353SBarry Smith         val = -dvdx / hx - cdiv3;
394b3e6a353SBarry Smith         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
395b3e6a353SBarry Smith       }
396b3e6a353SBarry Smith       if (i != nx && j != ny) {
397b3e6a353SBarry Smith         ind = k;
398b3e6a353SBarry Smith         val = dvdx / hx + dvdy / hy - cdiv3;
399b3e6a353SBarry Smith         PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES));
400b3e6a353SBarry Smith       }
401b3e6a353SBarry Smith     }
402b3e6a353SBarry Smith   }
403b3e6a353SBarry Smith   PetscCall(VecRestoreArrayRead(X, &x));
404b3e6a353SBarry Smith 
405b3e6a353SBarry Smith   /* Assemble gradient vector */
406b3e6a353SBarry Smith   PetscCall(VecAssemblyBegin(G));
407b3e6a353SBarry Smith   PetscCall(VecAssemblyEnd(G));
408b3e6a353SBarry Smith 
409b3e6a353SBarry Smith   /* Scale the gradient */
410b3e6a353SBarry Smith   area = p5 * hx * hy;
411b3e6a353SBarry Smith   PetscCall(VecScale(G, area));
412b3e6a353SBarry Smith   PetscCall(PetscLogFlops(24.0 * nx * ny));
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
414b3e6a353SBarry Smith }
415b3e6a353SBarry Smith 
416b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
417b3e6a353SBarry Smith /*
418b3e6a353SBarry Smith    FormHessian - Forms the Hessian matrix.
419b3e6a353SBarry Smith 
420b3e6a353SBarry Smith    Input Parameters:
421b3e6a353SBarry Smith .  tao - the Tao context
422b3e6a353SBarry Smith .  X    - the input vector
423b3e6a353SBarry Smith .  ptr  - optional user-defined context, as set by TaoSetHessian()
424b3e6a353SBarry Smith 
425b3e6a353SBarry Smith    Output Parameters:
426b3e6a353SBarry Smith .  H     - Hessian matrix
427b3e6a353SBarry Smith .  PrecH - optionally different preconditioning Hessian
428b3e6a353SBarry Smith 
429b3e6a353SBarry Smith    Notes:
430b3e6a353SBarry Smith    This routine is intended simply as an example of the interface
431b3e6a353SBarry Smith    to a Hessian evaluation routine.  Since this example compute the
432b3e6a353SBarry Smith    Hessian a column at a time, it is not particularly efficient and
433b3e6a353SBarry Smith    is not recommended.
434b3e6a353SBarry Smith */
FormHessian(Tao tao,Vec X,Mat H,Mat Hpre,void * ptr)435b3e6a353SBarry Smith PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
436b3e6a353SBarry Smith {
437b3e6a353SBarry Smith   AppCtx    *user = (AppCtx *)ptr;
438b3e6a353SBarry Smith   PetscInt   i, j, ndim = user->ndim;
439b3e6a353SBarry Smith   PetscReal *y, zero = 0.0, one = 1.0;
440b3e6a353SBarry Smith   PetscBool  assembled;
441b3e6a353SBarry Smith 
442b3e6a353SBarry Smith   PetscFunctionBeginUser;
443b3e6a353SBarry Smith   user->xvec = X;
444b3e6a353SBarry Smith 
445b3e6a353SBarry Smith   /* Initialize Hessian entries and work vector to zero */
446b3e6a353SBarry Smith   PetscCall(MatAssembled(H, &assembled));
447b3e6a353SBarry Smith   if (assembled) PetscCall(MatZeroEntries(H));
448b3e6a353SBarry Smith 
449b3e6a353SBarry Smith   PetscCall(VecSet(user->s, zero));
450b3e6a353SBarry Smith 
451b3e6a353SBarry Smith   /* Loop over matrix columns to compute entries of the Hessian */
452b3e6a353SBarry Smith   for (j = 0; j < ndim; j++) {
453b3e6a353SBarry Smith     PetscCall(VecSetValues(user->s, 1, &j, &one, INSERT_VALUES));
454b3e6a353SBarry Smith     PetscCall(VecAssemblyBegin(user->s));
455b3e6a353SBarry Smith     PetscCall(VecAssemblyEnd(user->s));
456b3e6a353SBarry Smith 
457b3e6a353SBarry Smith     PetscCall(HessianProduct(ptr, user->s, user->y));
458b3e6a353SBarry Smith 
459b3e6a353SBarry Smith     PetscCall(VecSetValues(user->s, 1, &j, &zero, INSERT_VALUES));
460b3e6a353SBarry Smith     PetscCall(VecAssemblyBegin(user->s));
461b3e6a353SBarry Smith     PetscCall(VecAssemblyEnd(user->s));
462b3e6a353SBarry Smith 
463b3e6a353SBarry Smith     PetscCall(VecGetArray(user->y, &y));
464b3e6a353SBarry Smith     for (i = 0; i < ndim; i++) {
465b3e6a353SBarry Smith       if (y[i] != zero) PetscCall(MatSetValues(H, 1, &i, 1, &j, &y[i], ADD_VALUES));
466b3e6a353SBarry Smith     }
467b3e6a353SBarry Smith     PetscCall(VecRestoreArray(user->y, &y));
468b3e6a353SBarry Smith   }
469b3e6a353SBarry Smith   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
470b3e6a353SBarry Smith   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
472b3e6a353SBarry Smith }
473b3e6a353SBarry Smith 
474b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
475b3e6a353SBarry Smith /*
476b3e6a353SBarry Smith    MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector
477b3e6a353SBarry Smith    products.
478b3e6a353SBarry Smith 
479b3e6a353SBarry Smith    Input Parameters:
480b3e6a353SBarry Smith .  tao - the Tao context
481b3e6a353SBarry Smith .  X    - the input vector
482b3e6a353SBarry Smith .  ptr  - optional user-defined context, as set by TaoSetHessian()
483b3e6a353SBarry Smith 
484b3e6a353SBarry Smith    Output Parameters:
485b3e6a353SBarry Smith .  H     - Hessian matrix
486b3e6a353SBarry Smith .  PrecH - optionally different preconditioning Hessian
487b3e6a353SBarry Smith */
MatrixFreeHessian(Tao tao,Vec X,Mat H,Mat PrecH,void * ptr)488b3e6a353SBarry Smith PetscErrorCode MatrixFreeHessian(Tao tao, Vec X, Mat H, Mat PrecH, void *ptr)
489b3e6a353SBarry Smith {
490b3e6a353SBarry Smith   AppCtx *user = (AppCtx *)ptr;
491b3e6a353SBarry Smith 
492b3e6a353SBarry Smith   /* Sets location of vector for use in computing matrix-vector products  of the form H(X)*y  */
493b3e6a353SBarry Smith   PetscFunctionBeginUser;
494b3e6a353SBarry Smith   user->xvec = X;
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
496b3e6a353SBarry Smith }
497b3e6a353SBarry Smith 
498b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
499b3e6a353SBarry Smith /*
500b3e6a353SBarry Smith    HessianProductMat - Computes the matrix-vector product
501b3e6a353SBarry Smith    y = mat*svec.
502b3e6a353SBarry Smith 
503b3e6a353SBarry Smith    Input Parameters:
504b3e6a353SBarry Smith .  mat  - input matrix
505b3e6a353SBarry Smith .  svec - input vector
506b3e6a353SBarry Smith 
507b3e6a353SBarry Smith    Output Parameters:
508b3e6a353SBarry Smith .  y    - solution vector
509b3e6a353SBarry Smith */
HessianProductMat(Mat mat,Vec svec,Vec y)510b3e6a353SBarry Smith PetscErrorCode HessianProductMat(Mat mat, Vec svec, Vec y)
511b3e6a353SBarry Smith {
512b3e6a353SBarry Smith   void *ptr;
513b3e6a353SBarry Smith 
514b3e6a353SBarry Smith   PetscFunctionBeginUser;
515b3e6a353SBarry Smith   PetscCall(MatShellGetContext(mat, &ptr));
516b3e6a353SBarry Smith   PetscCall(HessianProduct(ptr, svec, y));
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
518b3e6a353SBarry Smith }
519b3e6a353SBarry Smith 
520b3e6a353SBarry Smith /* ------------------------------------------------------------------- */
521b3e6a353SBarry Smith /*
522b3e6a353SBarry Smith    Hessian Product - Computes the matrix-vector product:
523b3e6a353SBarry Smith    y = f''(x)*svec.
524b3e6a353SBarry Smith 
525b3e6a353SBarry Smith    Input Parameters:
526b3e6a353SBarry Smith .  ptr  - pointer to the user-defined context
527b3e6a353SBarry Smith .  svec - input vector
528b3e6a353SBarry Smith 
529b3e6a353SBarry Smith    Output Parameters:
530b3e6a353SBarry Smith .  y    - product vector
531b3e6a353SBarry Smith */
HessianProduct(void * ptr,Vec svec,Vec y)532b3e6a353SBarry Smith PetscErrorCode HessianProduct(void *ptr, Vec svec, Vec y)
533b3e6a353SBarry Smith {
534b3e6a353SBarry Smith   AppCtx            *user = (AppCtx *)ptr;
535b3e6a353SBarry Smith   PetscReal          p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area;
536b3e6a353SBarry Smith   const PetscScalar *x, *s;
537b3e6a353SBarry Smith   PetscReal          v, vb, vl, vr, vt, hxhx, hyhy;
538b3e6a353SBarry Smith   PetscInt           nx, ny, i, j, k, ind;
539b3e6a353SBarry Smith 
540b3e6a353SBarry Smith   PetscFunctionBeginUser;
541b3e6a353SBarry Smith   nx   = user->mx;
542b3e6a353SBarry Smith   ny   = user->my;
543b3e6a353SBarry Smith   hx   = user->hx;
544b3e6a353SBarry Smith   hy   = user->hy;
545b3e6a353SBarry Smith   hxhx = one / (hx * hx);
546b3e6a353SBarry Smith   hyhy = one / (hy * hy);
547b3e6a353SBarry Smith 
548b3e6a353SBarry Smith   /* Get pointers to vector data */
549b3e6a353SBarry Smith   PetscCall(VecGetArrayRead(user->xvec, &x));
550b3e6a353SBarry Smith   PetscCall(VecGetArrayRead(svec, &s));
551b3e6a353SBarry Smith 
552b3e6a353SBarry Smith   /* Initialize product vector to zero */
553b3e6a353SBarry Smith   PetscCall(VecSet(y, zero));
554b3e6a353SBarry Smith 
555b3e6a353SBarry Smith   /* Compute f''(x)*s over the lower triangular elements */
556b3e6a353SBarry Smith   for (j = -1; j < ny; j++) {
557b3e6a353SBarry Smith     for (i = -1; i < nx; i++) {
558b3e6a353SBarry Smith       k  = nx * j + i;
559b3e6a353SBarry Smith       v  = zero;
560b3e6a353SBarry Smith       vr = zero;
561b3e6a353SBarry Smith       vt = zero;
562b3e6a353SBarry Smith       if (i != -1 && j != -1) v = s[k];
563b3e6a353SBarry Smith       if (i != nx - 1 && j != -1) {
564b3e6a353SBarry Smith         vr  = s[k + 1];
565b3e6a353SBarry Smith         ind = k + 1;
566b3e6a353SBarry Smith         val = hxhx * (vr - v);
567b3e6a353SBarry Smith         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
568b3e6a353SBarry Smith       }
569b3e6a353SBarry Smith       if (i != -1 && j != ny - 1) {
570b3e6a353SBarry Smith         vt  = s[k + nx];
571b3e6a353SBarry Smith         ind = k + nx;
572b3e6a353SBarry Smith         val = hyhy * (vt - v);
573b3e6a353SBarry Smith         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
574b3e6a353SBarry Smith       }
575b3e6a353SBarry Smith       if (i != -1 && j != -1) {
576b3e6a353SBarry Smith         ind = k;
577b3e6a353SBarry Smith         val = hxhx * (v - vr) + hyhy * (v - vt);
578b3e6a353SBarry Smith         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
579b3e6a353SBarry Smith       }
580b3e6a353SBarry Smith     }
581b3e6a353SBarry Smith   }
582b3e6a353SBarry Smith 
583b3e6a353SBarry Smith   /* Compute f''(x)*s over the upper triangular elements */
584b3e6a353SBarry Smith   for (j = 0; j <= ny; j++) {
585b3e6a353SBarry Smith     for (i = 0; i <= nx; i++) {
586b3e6a353SBarry Smith       k  = nx * j + i;
587b3e6a353SBarry Smith       v  = zero;
588b3e6a353SBarry Smith       vl = zero;
589b3e6a353SBarry Smith       vb = zero;
590b3e6a353SBarry Smith       if (i != nx && j != ny) v = s[k];
591b3e6a353SBarry Smith       if (i != nx && j != 0) {
592b3e6a353SBarry Smith         vb  = s[k - nx];
593b3e6a353SBarry Smith         ind = k - nx;
594b3e6a353SBarry Smith         val = hyhy * (vb - v);
595b3e6a353SBarry Smith         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
596b3e6a353SBarry Smith       }
597b3e6a353SBarry Smith       if (i != 0 && j != ny) {
598b3e6a353SBarry Smith         vl  = s[k - 1];
599b3e6a353SBarry Smith         ind = k - 1;
600b3e6a353SBarry Smith         val = hxhx * (vl - v);
601b3e6a353SBarry Smith         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
602b3e6a353SBarry Smith       }
603b3e6a353SBarry Smith       if (i != nx && j != ny) {
604b3e6a353SBarry Smith         ind = k;
605b3e6a353SBarry Smith         val = hxhx * (v - vl) + hyhy * (v - vb);
606b3e6a353SBarry Smith         PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES));
607b3e6a353SBarry Smith       }
608b3e6a353SBarry Smith     }
609b3e6a353SBarry Smith   }
610b3e6a353SBarry Smith   /* Restore vector data */
611b3e6a353SBarry Smith   PetscCall(VecRestoreArrayRead(svec, &s));
612b3e6a353SBarry Smith   PetscCall(VecRestoreArrayRead(user->xvec, &x));
613b3e6a353SBarry Smith 
614b3e6a353SBarry Smith   /* Assemble vector */
615b3e6a353SBarry Smith   PetscCall(VecAssemblyBegin(y));
616b3e6a353SBarry Smith   PetscCall(VecAssemblyEnd(y));
617b3e6a353SBarry Smith 
618b3e6a353SBarry Smith   /* Scale resulting vector by area */
619b3e6a353SBarry Smith   area = p5 * hx * hy;
620b3e6a353SBarry Smith   PetscCall(VecScale(y, area));
621b3e6a353SBarry Smith   PetscCall(PetscLogFlops(18.0 * nx * ny));
6223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
623b3e6a353SBarry Smith }
624b3e6a353SBarry Smith 
625b3e6a353SBarry Smith /*TEST
626b3e6a353SBarry Smith 
627b3e6a353SBarry Smith    build:
628b3e6a353SBarry Smith       requires: !complex
629b3e6a353SBarry Smith 
630b3e6a353SBarry Smith    test:
631b3e6a353SBarry Smith       args: -tao_monitor -tao_type bntr -tao_view -tao_bnk_ksp_monitor_short
632b3e6a353SBarry Smith 
633b3e6a353SBarry Smith TEST*/
634