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