xref: /petsc/src/tao/unconstrained/tutorials/eptorsion2.c (revision 70a7d78aacfbd24b2e31399a3d8e056944bb7de3)
1c4762a1bSJed Brown /* Program usage: mpiexec -n <proc> eptorsion2 [-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 uniprocessor version of this code is eptorsion1.c; the Fortran
16c4762a1bSJed Brown   version of this code is eptorsion2f.F.
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   This application solves the problem without calculating hessians
19c4762a1bSJed Brown ---------------------------------------------------------------------- */
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown   Include "petsctao.h" so that we can use TAO solvers.  Note that this
23c4762a1bSJed Brown   file automatically includes files for lower-level support, such as those
24c4762a1bSJed Brown   provided by the PETSc library:
25c4762a1bSJed Brown      petsc.h       - base PETSc routines   petscvec.h - vectors
26c4762a1bSJed Brown      petscsys.h    - sysem routines        petscmat.h - matrices
27c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
28c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
29c4762a1bSJed Brown   Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing
30c4762a1bSJed Brown   the parallel mesh.
31c4762a1bSJed Brown */
32c4762a1bSJed Brown 
33c4762a1bSJed Brown #include <petsctao.h>
34c4762a1bSJed Brown #include <petscdmda.h>
35c4762a1bSJed Brown 
36c4762a1bSJed Brown static  char help[] =
37c4762a1bSJed Brown "Demonstrates use of the TAO package to solve \n\
38c4762a1bSJed Brown unconstrained minimization problems in parallel.  This example is based on \n\
39c4762a1bSJed Brown the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 test suite.\n\
40c4762a1bSJed Brown The command line options are:\n\
41c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
42c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
43c4762a1bSJed Brown   -par <param>, where <param> = angle of twist per unit length\n\n";
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /*T
46c4762a1bSJed Brown    Concepts: TAO^Solving an unconstrained minimization problem
47c4762a1bSJed Brown    Routines: TaoCreate(); TaoSetType();
48c4762a1bSJed Brown    Routines: TaoSetInitialVector();
49c4762a1bSJed Brown    Routines: TaoSetObjectiveAndGradientRoutine();
50c4762a1bSJed Brown    Routines: TaoSetHessianRoutine(); TaoSetFromOptions();
51c4762a1bSJed Brown    Routines: TaoSolve();
52c4762a1bSJed Brown    Routines: TaoDestroy();
53c4762a1bSJed Brown    Processors: n
54c4762a1bSJed Brown T*/
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*
57c4762a1bSJed Brown    User-defined application context - contains data needed by the
58c4762a1bSJed Brown    application-provided call-back routines, FormFunction() and
59c4762a1bSJed Brown    FormGradient().
60c4762a1bSJed Brown */
61c4762a1bSJed Brown typedef struct {
62c4762a1bSJed Brown   /* parameters */
63c4762a1bSJed Brown    PetscInt      mx, my;         /* global discretization in x- and y-directions */
64c4762a1bSJed Brown    PetscReal     param;          /* nonlinearity parameter */
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* work space */
67c4762a1bSJed Brown    Vec           localX;         /* local vectors */
68c4762a1bSJed Brown    DM            dm;             /* distributed array data structure */
69c4762a1bSJed Brown } AppCtx;
70c4762a1bSJed Brown 
71c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx*, Vec);
72c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
73c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
74c4762a1bSJed Brown 
75c4762a1bSJed Brown int main(int argc, char **argv)
76c4762a1bSJed Brown {
77c4762a1bSJed Brown     PetscErrorCode     ierr;
78c4762a1bSJed Brown     Vec                x;
79c4762a1bSJed Brown     Mat                H;
80c4762a1bSJed Brown     PetscInt           Nx, Ny;
81c4762a1bSJed Brown     Tao                tao;
82c4762a1bSJed Brown     PetscBool          flg;
83c4762a1bSJed Brown     KSP                ksp;
84c4762a1bSJed Brown     PC                 pc;
85c4762a1bSJed Brown     AppCtx             user;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown     PetscInitialize(&argc, &argv, (char *)0, help);
88c4762a1bSJed Brown 
89c4762a1bSJed Brown     /* Specify default dimension of the problem */
90c4762a1bSJed Brown     user.param = 5.0; user.mx = 10; user.my = 10;
91c4762a1bSJed Brown     Nx = Ny = PETSC_DECIDE;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     /* Check for any command line arguments that override defaults */
94c4762a1bSJed Brown     ierr = PetscOptionsGetReal(NULL,NULL,"-par",&user.param,&flg);CHKERRQ(ierr);
95c4762a1bSJed Brown     ierr = PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg);CHKERRQ(ierr);
96c4762a1bSJed Brown     ierr = PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg);CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n---- Elastic-Plastic Torsion Problem -----\n");CHKERRQ(ierr);
99c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"mx: %D     my: %D   \n\n",user.mx,user.my);CHKERRQ(ierr);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown     /* Set up distributed array */
102c4762a1bSJed Brown     ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm);CHKERRQ(ierr);
103c4762a1bSJed Brown     ierr = DMSetFromOptions(user.dm);CHKERRQ(ierr);
104c4762a1bSJed Brown     ierr = DMSetUp(user.dm);CHKERRQ(ierr);
105c4762a1bSJed Brown 
106c4762a1bSJed Brown     /* Create vectors */
107c4762a1bSJed Brown     ierr = DMCreateGlobalVector(user.dm,&x);CHKERRQ(ierr);
108c4762a1bSJed Brown 
109c4762a1bSJed Brown     ierr = DMCreateLocalVector(user.dm,&user.localX);CHKERRQ(ierr);
110c4762a1bSJed Brown 
111c4762a1bSJed Brown     /* Create Hessian */
112c4762a1bSJed Brown     ierr = DMCreateMatrix(user.dm,&H);CHKERRQ(ierr);
113c4762a1bSJed Brown     ierr = MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
114c4762a1bSJed Brown 
115c4762a1bSJed Brown     /* The TAO code begins here */
116c4762a1bSJed Brown 
117c4762a1bSJed Brown     /* Create TAO solver and set desired solution method */
118c4762a1bSJed Brown     ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
119c4762a1bSJed Brown     ierr = TaoSetType(tao,TAOCG);CHKERRQ(ierr);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown     /* Set initial solution guess */
122c4762a1bSJed Brown     ierr = FormInitialGuess(&user,x);CHKERRQ(ierr);
123c4762a1bSJed Brown     ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
124c4762a1bSJed Brown 
125c4762a1bSJed Brown     /* Set routine for function and gradient evaluation */
126c4762a1bSJed Brown     ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);CHKERRQ(ierr);
127c4762a1bSJed Brown 
128c4762a1bSJed Brown     ierr = TaoSetHessianRoutine(tao,H,H,FormHessian,(void*)&user);CHKERRQ(ierr);
129c4762a1bSJed Brown 
130c4762a1bSJed Brown     /* Check for any TAO command line options */
131c4762a1bSJed Brown     ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown     ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
134c4762a1bSJed Brown     if (ksp) {
135c4762a1bSJed Brown       ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
136c4762a1bSJed Brown       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
137c4762a1bSJed Brown     }
138c4762a1bSJed Brown 
139c4762a1bSJed Brown     /* SOLVE THE APPLICATION */
140c4762a1bSJed Brown     ierr = TaoSolve(tao);CHKERRQ(ierr);
141c4762a1bSJed Brown 
142c4762a1bSJed Brown     /* Free TAO data structures */
143c4762a1bSJed Brown     ierr = TaoDestroy(&tao);CHKERRQ(ierr);
144c4762a1bSJed Brown 
145c4762a1bSJed Brown     /* Free PETSc data structures */
146c4762a1bSJed Brown     ierr = VecDestroy(&x);CHKERRQ(ierr);
147c4762a1bSJed Brown     ierr = MatDestroy(&H);CHKERRQ(ierr);
148c4762a1bSJed Brown 
149c4762a1bSJed Brown     ierr = VecDestroy(&user.localX);CHKERRQ(ierr);
150c4762a1bSJed Brown     ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
151c4762a1bSJed Brown 
152c4762a1bSJed Brown     PetscFinalize();
153c4762a1bSJed Brown     return 0;
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown /* ------------------------------------------------------------------- */
157c4762a1bSJed Brown /*
158c4762a1bSJed Brown     FormInitialGuess - Computes an initial approximation to the solution.
159c4762a1bSJed Brown 
160c4762a1bSJed Brown     Input Parameters:
161c4762a1bSJed Brown .   user - user-defined application context
162c4762a1bSJed Brown .   X    - vector
163c4762a1bSJed Brown 
164c4762a1bSJed Brown     Output Parameters:
165c4762a1bSJed Brown     X    - vector
166c4762a1bSJed Brown */
167c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
168c4762a1bSJed Brown {
169c4762a1bSJed Brown   PetscErrorCode ierr;
170c4762a1bSJed Brown   PetscInt       i, j, k, mx = user->mx, my = user->my;
171c4762a1bSJed Brown   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye;
172c4762a1bSJed Brown   PetscReal      hx = 1.0/(mx+1), hy = 1.0/(my+1), temp, val;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   PetscFunctionBegin;
175c4762a1bSJed Brown   /* Get local mesh boundaries */
176c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* Compute initial guess over locally owned part of mesh */
180c4762a1bSJed Brown   xe = xs+xm;
181c4762a1bSJed Brown   ye = ys+ym;
182c4762a1bSJed Brown   for (j=ys; j<ye; j++) {  /*  for (j=0; j<my; j++) */
183c4762a1bSJed Brown     temp = PetscMin(j+1,my-j)*hy;
184c4762a1bSJed Brown     for (i=xs; i<xe; i++) {  /*  for (i=0; i<mx; i++) */
185c4762a1bSJed Brown       k   = (j-gys)*gxm + i-gxs;
186c4762a1bSJed Brown       val = PetscMin((PetscMin(i+1,mx-i))*hx,temp);
187c4762a1bSJed Brown       ierr = VecSetValuesLocal(X,1,&k,&val,ADD_VALUES);CHKERRQ(ierr);
188c4762a1bSJed Brown     }
189c4762a1bSJed Brown   }
190c4762a1bSJed Brown   ierr = VecAssemblyBegin(X);CHKERRQ(ierr);
191c4762a1bSJed Brown   ierr = VecAssemblyEnd(X);CHKERRQ(ierr);
192c4762a1bSJed Brown   PetscFunctionReturn(0);
193c4762a1bSJed Brown }
194c4762a1bSJed Brown 
195c4762a1bSJed Brown /* ------------------------------------------------------------------ */
196c4762a1bSJed Brown /*
197c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    Input Parameters:
200c4762a1bSJed Brown    tao - the Tao context
201c4762a1bSJed Brown    X   - the input vector
202c4762a1bSJed Brown    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    Output Parameters:
205c4762a1bSJed Brown    f   - the newly evaluated function
206c4762a1bSJed Brown    G   - the newly evaluated gradient
207c4762a1bSJed Brown */
208*70a7d78aSStefano Zampini PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr)
209*70a7d78aSStefano Zampini {
210c4762a1bSJed Brown   AppCtx         *user = (AppCtx *)ptr;
211c4762a1bSJed Brown   PetscErrorCode ierr;
212c4762a1bSJed Brown   PetscInt       i,j,k,ind;
213c4762a1bSJed Brown   PetscInt       xe,ye,xsm,ysm,xep,yep;
214c4762a1bSJed Brown   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys;
215c4762a1bSJed Brown   PetscInt       mx = user->mx, my = user->my;
216c4762a1bSJed Brown   PetscReal      three = 3.0, zero = 0.0, *x, floc, cdiv3 = user->param/three;
217c4762a1bSJed Brown   PetscReal      p5 = 0.5, area, val, flin, fquad;
218c4762a1bSJed Brown   PetscReal      v,vb,vl,vr,vt,dvdx,dvdy;
219c4762a1bSJed Brown   PetscReal      hx = 1.0/(user->mx + 1);
220c4762a1bSJed Brown   PetscReal      hy = 1.0/(user->my + 1);
221c4762a1bSJed Brown   Vec            localX = user->localX;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   PetscFunctionBegin;
224c4762a1bSJed Brown   /* Initialize */
225c4762a1bSJed Brown   flin = fquad = zero;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   ierr = VecSet(G, zero);CHKERRQ(ierr);
228c4762a1bSJed Brown   /*
229c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
230c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
231c4762a1bSJed Brown      By placing code between these two statements, computations can be
232c4762a1bSJed Brown      done while messages are in transition.
233c4762a1bSJed Brown   */
234c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /* Get pointer to vector data */
238c4762a1bSJed Brown   ierr = VecGetArray(localX,&x);CHKERRQ(ierr);
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   /* Get local mesh boundaries */
241c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
242c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /* Set local loop dimensions */
245c4762a1bSJed Brown   xe = xs+xm;
246c4762a1bSJed Brown   ye = ys+ym;
247c4762a1bSJed Brown   if (xs == 0)  xsm = xs-1;
248c4762a1bSJed Brown   else          xsm = xs;
249c4762a1bSJed Brown   if (ys == 0)  ysm = ys-1;
250c4762a1bSJed Brown   else          ysm = ys;
251c4762a1bSJed Brown   if (xe == mx) xep = xe+1;
252c4762a1bSJed Brown   else          xep = xe;
253c4762a1bSJed Brown   if (ye == my) yep = ye+1;
254c4762a1bSJed Brown   else          yep = ye;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Compute local gradient contributions over the lower triangular elements */
257c4762a1bSJed Brown   for (j=ysm; j<ye; j++) {  /*  for (j=-1; j<my; j++) */
258c4762a1bSJed Brown     for (i=xsm; i<xe; i++) {  /*   for (i=-1; i<mx; i++) */
259c4762a1bSJed Brown       k = (j-gys)*gxm + i-gxs;
260c4762a1bSJed Brown       v = zero;
261c4762a1bSJed Brown       vr = zero;
262c4762a1bSJed Brown       vt = zero;
263c4762a1bSJed Brown       if (i >= 0 && j >= 0) v = x[k];
264c4762a1bSJed Brown       if (i < mx-1 && j > -1) vr = x[k+1];
265c4762a1bSJed Brown       if (i > -1 && j < my-1) vt = x[k+gxm];
266c4762a1bSJed Brown       dvdx = (vr-v)/hx;
267c4762a1bSJed Brown       dvdy = (vt-v)/hy;
268c4762a1bSJed Brown       if (i != -1 && j != -1) {
269c4762a1bSJed Brown         ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
270c4762a1bSJed Brown         ierr = VecSetValuesLocal(G,1,&k,&val,ADD_VALUES);CHKERRQ(ierr);
271c4762a1bSJed Brown       }
272c4762a1bSJed Brown       if (i != mx-1 && j != -1) {
273c4762a1bSJed Brown         ind = k+1; val =  dvdx/hx - cdiv3;
274c4762a1bSJed Brown         ierr = VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);CHKERRQ(ierr);
275c4762a1bSJed Brown       }
276c4762a1bSJed Brown       if (i != -1 && j != my-1) {
277c4762a1bSJed Brown         ind = k+gxm; val = dvdy/hy - cdiv3;
278c4762a1bSJed Brown         ierr = VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);CHKERRQ(ierr);
279c4762a1bSJed Brown       }
280c4762a1bSJed Brown       fquad += dvdx*dvdx + dvdy*dvdy;
281c4762a1bSJed Brown       flin -= cdiv3 * (v + vr + vt);
282c4762a1bSJed Brown     }
283c4762a1bSJed Brown   }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   /* Compute local gradient contributions over the upper triangular elements */
286c4762a1bSJed Brown   for (j=ys; j<yep; j++) { /*  for (j=0; j<=my; j++) */
287c4762a1bSJed Brown     for (i=xs; i<xep; i++) {  /*   for (i=0; i<=mx; i++) */
288c4762a1bSJed Brown       k = (j-gys)*gxm + i-gxs;
289c4762a1bSJed Brown       vb = zero;
290c4762a1bSJed Brown       vl = zero;
291c4762a1bSJed Brown       v  = zero;
292c4762a1bSJed Brown       if (i < mx && j > 0) vb = x[k-gxm];
293c4762a1bSJed Brown       if (i > 0 && j < my) vl = x[k-1];
294c4762a1bSJed Brown       if (i < mx && j < my) v = x[k];
295c4762a1bSJed Brown       dvdx = (v-vl)/hx;
296c4762a1bSJed Brown       dvdy = (v-vb)/hy;
297c4762a1bSJed Brown       if (i != mx && j != 0) {
298c4762a1bSJed Brown         ind = k-gxm; val = - dvdy/hy - cdiv3;
299c4762a1bSJed Brown         ierr = VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);CHKERRQ(ierr);
300c4762a1bSJed Brown       }
301c4762a1bSJed Brown       if (i != 0 && j != my) {
302c4762a1bSJed Brown         ind = k-1; val =  - dvdx/hx - cdiv3;
303c4762a1bSJed Brown         ierr = VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);CHKERRQ(ierr);
304c4762a1bSJed Brown       }
305c4762a1bSJed Brown       if (i != mx && j != my) {
306c4762a1bSJed Brown         ind = k; val =  dvdx/hx + dvdy/hy - cdiv3;
307c4762a1bSJed Brown         ierr = VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);CHKERRQ(ierr);
308c4762a1bSJed Brown       }
309c4762a1bSJed Brown       fquad += dvdx*dvdx + dvdy*dvdy;
310c4762a1bSJed Brown       flin -= cdiv3 * (vb + vl + v);
311c4762a1bSJed Brown     }
312c4762a1bSJed Brown   }
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   /* Restore vector */
315c4762a1bSJed Brown   ierr = VecRestoreArray(localX,&x);CHKERRQ(ierr);
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   /* Assemble gradient vector */
318c4762a1bSJed Brown   ierr = VecAssemblyBegin(G);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = VecAssemblyEnd(G);CHKERRQ(ierr);
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   /* Scale the gradient */
322c4762a1bSJed Brown   area = p5*hx*hy;
323c4762a1bSJed Brown   floc = area * (p5 * fquad + flin);
324c4762a1bSJed Brown   ierr = VecScale(G, area);CHKERRQ(ierr);
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /* Sum function contributions from all processes */
327c4762a1bSJed Brown   ierr = (PetscErrorCode)MPI_Allreduce((void*)&floc,(void*)f,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD);CHKERRQ(ierr);
328c4762a1bSJed Brown 
329c4762a1bSJed Brown   ierr=PetscLogFlops((ye-ysm)*(xe-xsm)*20+(xep-xs)*(yep-ys)*16);CHKERRQ(ierr);
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   PetscFunctionReturn(0);
332c4762a1bSJed Brown }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao, Vec X, Mat A, Mat Hpre, void*ctx)
335c4762a1bSJed Brown {
336c4762a1bSJed Brown   AppCtx         *user= (AppCtx*) ctx;
337c4762a1bSJed Brown   PetscErrorCode ierr;
338c4762a1bSJed Brown   PetscInt       i,j,k;
339c4762a1bSJed Brown   PetscInt       col[5],row;
340c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
341c4762a1bSJed Brown   PetscReal      v[5];
342c4762a1bSJed Brown   PetscReal      hx=1.0/(user->mx+1), hy=1.0/(user->my+1), hxhx=1.0/(hx*hx), hyhy=1.0/(hy*hy), area=0.5*hx*hy;
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   /* Compute the quadratic term in the objective function */
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   /*
347c4762a1bSJed Brown      Get local grid boundaries
348c4762a1bSJed Brown   */
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   PetscFunctionBegin;
351c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
352c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
353c4762a1bSJed Brown 
354c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
355c4762a1bSJed Brown 
356c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
357c4762a1bSJed Brown 
358c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
359c4762a1bSJed Brown 
360c4762a1bSJed Brown       k=0;
361c4762a1bSJed Brown       if (j>gys) {
362c4762a1bSJed Brown         v[k]=-2*hyhy; col[k]=row - gxm; k++;
363c4762a1bSJed Brown       }
364c4762a1bSJed Brown 
365c4762a1bSJed Brown       if (i>gxs) {
366c4762a1bSJed Brown         v[k]= -2*hxhx; col[k]=row - 1; k++;
367c4762a1bSJed Brown       }
368c4762a1bSJed Brown 
369c4762a1bSJed Brown       v[k]= 4.0*(hxhx+hyhy); col[k]=row; k++;
370c4762a1bSJed Brown 
371c4762a1bSJed Brown       if (i+1 < gxs+gxm) {
372c4762a1bSJed Brown         v[k]= -2.0*hxhx; col[k]=row+1; k++;
373c4762a1bSJed Brown       }
374c4762a1bSJed Brown 
375c4762a1bSJed Brown       if (j+1 <gys+gym) {
376c4762a1bSJed Brown         v[k]= -2*hyhy; col[k] = row+gxm; k++;
377c4762a1bSJed Brown       }
378c4762a1bSJed Brown 
379c4762a1bSJed Brown       ierr = MatSetValuesLocal(A,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr);
380c4762a1bSJed Brown 
381c4762a1bSJed Brown     }
382c4762a1bSJed Brown   }
383c4762a1bSJed Brown   /*
384c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
385c4762a1bSJed Brown      MatAssemblyBegin(), MatAssemblyEnd().
386c4762a1bSJed Brown      By placing code between these two statements, computations can be
387c4762a1bSJed Brown      done while messages are in transition.
388c4762a1bSJed Brown   */
389c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
390c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
391c4762a1bSJed Brown   /*
392c4762a1bSJed Brown     Tell the matrix we will never add a new nonzero location to the
393c4762a1bSJed Brown     matrix. If we do it will generate an error.
394c4762a1bSJed Brown   */
395c4762a1bSJed Brown   ierr = MatScale(A,area);CHKERRQ(ierr);
396c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
397c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
398c4762a1bSJed Brown   ierr = PetscLogFlops(9*xm*ym+49*xm);CHKERRQ(ierr);
399c4762a1bSJed Brown   PetscFunctionReturn(0);
400c4762a1bSJed Brown }
401c4762a1bSJed Brown 
402c4762a1bSJed Brown /*TEST
403c4762a1bSJed Brown 
404c4762a1bSJed Brown    build:
405c4762a1bSJed Brown       requires: !complex
406c4762a1bSJed Brown 
407c4762a1bSJed Brown    test:
408c4762a1bSJed Brown       suffix: 1
409c4762a1bSJed Brown       nsize: 2
410c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -mx 16 -my 16 -tao_gatol 1.e-4
411c4762a1bSJed Brown 
412c4762a1bSJed Brown TEST*/
413