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