xref: /petsc/src/tao/unconstrained/tutorials/eptorsion2.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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
26a5b23f4aSJose E. Roman      petscsys.h    - system 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 /*
46c4762a1bSJed Brown    User-defined application context - contains data needed by the
47c4762a1bSJed Brown    application-provided call-back routines, FormFunction() and
48c4762a1bSJed Brown    FormGradient().
49c4762a1bSJed Brown */
50c4762a1bSJed Brown typedef struct {
51c4762a1bSJed Brown   /* parameters */
52c4762a1bSJed Brown    PetscInt      mx, my;         /* global discretization in x- and y-directions */
53c4762a1bSJed Brown    PetscReal     param;          /* nonlinearity parameter */
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* work space */
56c4762a1bSJed Brown    Vec           localX;         /* local vectors */
57c4762a1bSJed Brown    DM            dm;             /* distributed array data structure */
58c4762a1bSJed Brown } AppCtx;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx*, Vec);
61c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
62c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
63c4762a1bSJed Brown 
64c4762a1bSJed Brown int main(int argc, char **argv)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown     Vec                x;
67c4762a1bSJed Brown     Mat                H;
68c4762a1bSJed Brown     PetscInt           Nx, Ny;
69c4762a1bSJed Brown     Tao                tao;
70c4762a1bSJed Brown     PetscBool          flg;
71c4762a1bSJed Brown     KSP                ksp;
72c4762a1bSJed Brown     PC                 pc;
73c4762a1bSJed Brown     AppCtx             user;
74c4762a1bSJed Brown 
75c4762a1bSJed Brown     PetscInitialize(&argc, &argv, (char *)0, help);
76c4762a1bSJed Brown 
77c4762a1bSJed Brown     /* Specify default dimension of the problem */
78c4762a1bSJed Brown     user.param = 5.0; user.mx = 10; user.my = 10;
79c4762a1bSJed Brown     Nx = Ny = PETSC_DECIDE;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     /* Check for any command line arguments that override defaults */
829566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,&flg));
839566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- Elastic-Plastic Torsion Problem -----\n"));
87*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n",user.mx,user.my));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown     /* Set up distributed array */
909566063dSJacob Faibussowitsch     PetscCall(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));
919566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(user.dm));
929566063dSJacob Faibussowitsch     PetscCall(DMSetUp(user.dm));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown     /* Create vectors */
959566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(user.dm,&x));
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(user.dm,&user.localX));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown     /* Create Hessian */
1009566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(user.dm,&H));
1019566063dSJacob Faibussowitsch     PetscCall(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown     /* The TAO code begins here */
104c4762a1bSJed Brown 
105c4762a1bSJed Brown     /* Create TAO solver and set desired solution method */
1069566063dSJacob Faibussowitsch     PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
1079566063dSJacob Faibussowitsch     PetscCall(TaoSetType(tao,TAOCG));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown     /* Set initial solution guess */
1109566063dSJacob Faibussowitsch     PetscCall(FormInitialGuess(&user,x));
1119566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(tao,x));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown     /* Set routine for function and gradient evaluation */
1149566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao,H,H,FormHessian,(void*)&user));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown     /* Check for any TAO command line options */
1199566063dSJacob Faibussowitsch     PetscCall(TaoSetFromOptions(tao));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch     PetscCall(TaoGetKSP(tao,&ksp));
122c4762a1bSJed Brown     if (ksp) {
1239566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp,&pc));
1249566063dSJacob Faibussowitsch       PetscCall(PCSetType(pc,PCNONE));
125c4762a1bSJed Brown     }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown     /* SOLVE THE APPLICATION */
1289566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown     /* Free TAO data structures */
1319566063dSJacob Faibussowitsch     PetscCall(TaoDestroy(&tao));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown     /* Free PETSc data structures */
1349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
1359566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&H));
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&user.localX));
1389566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&user.dm));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     PetscFinalize();
141c4762a1bSJed Brown     return 0;
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown /* ------------------------------------------------------------------- */
145c4762a1bSJed Brown /*
146c4762a1bSJed Brown     FormInitialGuess - Computes an initial approximation to the solution.
147c4762a1bSJed Brown 
148c4762a1bSJed Brown     Input Parameters:
149c4762a1bSJed Brown .   user - user-defined application context
150c4762a1bSJed Brown .   X    - vector
151c4762a1bSJed Brown 
152c4762a1bSJed Brown     Output Parameters:
153c4762a1bSJed Brown     X    - vector
154c4762a1bSJed Brown */
155c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
156c4762a1bSJed Brown {
157c4762a1bSJed Brown   PetscInt       i, j, k, mx = user->mx, my = user->my;
158c4762a1bSJed Brown   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye;
159c4762a1bSJed Brown   PetscReal      hx = 1.0/(mx+1), hy = 1.0/(my+1), temp, val;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   PetscFunctionBegin;
162c4762a1bSJed Brown   /* Get local mesh boundaries */
1639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
1649566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Compute initial guess over locally owned part of mesh */
167c4762a1bSJed Brown   xe = xs+xm;
168c4762a1bSJed Brown   ye = ys+ym;
169c4762a1bSJed Brown   for (j=ys; j<ye; j++) {  /*  for (j=0; j<my; j++) */
170c4762a1bSJed Brown     temp = PetscMin(j+1,my-j)*hy;
171c4762a1bSJed Brown     for (i=xs; i<xe; i++) {  /*  for (i=0; i<mx; i++) */
172c4762a1bSJed Brown       k   = (j-gys)*gxm + i-gxs;
173c4762a1bSJed Brown       val = PetscMin((PetscMin(i+1,mx-i))*hx,temp);
1749566063dSJacob Faibussowitsch       PetscCall(VecSetValuesLocal(X,1,&k,&val,ADD_VALUES));
175c4762a1bSJed Brown     }
176c4762a1bSJed Brown   }
1779566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
1789566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
179c4762a1bSJed Brown   PetscFunctionReturn(0);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown /* ------------------------------------------------------------------ */
183c4762a1bSJed Brown /*
184c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    Input Parameters:
187c4762a1bSJed Brown    tao - the Tao context
188c4762a1bSJed Brown    X   - the input vector
189a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
190c4762a1bSJed Brown 
191c4762a1bSJed Brown    Output Parameters:
192c4762a1bSJed Brown    f   - the newly evaluated function
193c4762a1bSJed Brown    G   - the newly evaluated gradient
194c4762a1bSJed Brown */
19570a7d78aSStefano Zampini PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr)
19670a7d78aSStefano Zampini {
197c4762a1bSJed Brown   AppCtx         *user = (AppCtx *)ptr;
198c4762a1bSJed Brown   PetscInt       i,j,k,ind;
199c4762a1bSJed Brown   PetscInt       xe,ye,xsm,ysm,xep,yep;
200c4762a1bSJed Brown   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys;
201c4762a1bSJed Brown   PetscInt       mx = user->mx, my = user->my;
202c4762a1bSJed Brown   PetscReal      three = 3.0, zero = 0.0, *x, floc, cdiv3 = user->param/three;
203c4762a1bSJed Brown   PetscReal      p5 = 0.5, area, val, flin, fquad;
204c4762a1bSJed Brown   PetscReal      v,vb,vl,vr,vt,dvdx,dvdy;
205c4762a1bSJed Brown   PetscReal      hx = 1.0/(user->mx + 1);
206c4762a1bSJed Brown   PetscReal      hy = 1.0/(user->my + 1);
207c4762a1bSJed Brown   Vec            localX = user->localX;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   PetscFunctionBegin;
210c4762a1bSJed Brown   /* Initialize */
211c4762a1bSJed Brown   flin = fquad = zero;
212c4762a1bSJed Brown 
2139566063dSJacob Faibussowitsch   PetscCall(VecSet(G, zero));
214c4762a1bSJed Brown   /*
215c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
216c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
217c4762a1bSJed Brown      By placing code between these two statements, computations can be
218c4762a1bSJed Brown      done while messages are in transition.
219c4762a1bSJed Brown   */
2209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
2219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* Get pointer to vector data */
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX,&x));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* Get local mesh boundaries */
2279566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
2289566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /* Set local loop dimensions */
231c4762a1bSJed Brown   xe = xs+xm;
232c4762a1bSJed Brown   ye = ys+ym;
233c4762a1bSJed Brown   if (xs == 0)  xsm = xs-1;
234c4762a1bSJed Brown   else          xsm = xs;
235c4762a1bSJed Brown   if (ys == 0)  ysm = ys-1;
236c4762a1bSJed Brown   else          ysm = ys;
237c4762a1bSJed Brown   if (xe == mx) xep = xe+1;
238c4762a1bSJed Brown   else          xep = xe;
239c4762a1bSJed Brown   if (ye == my) yep = ye+1;
240c4762a1bSJed Brown   else          yep = ye;
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Compute local gradient contributions over the lower triangular elements */
243c4762a1bSJed Brown   for (j=ysm; j<ye; j++) {  /*  for (j=-1; j<my; j++) */
244c4762a1bSJed Brown     for (i=xsm; i<xe; i++) {  /*   for (i=-1; i<mx; i++) */
245c4762a1bSJed Brown       k = (j-gys)*gxm + i-gxs;
246c4762a1bSJed Brown       v = zero;
247c4762a1bSJed Brown       vr = zero;
248c4762a1bSJed Brown       vt = zero;
249c4762a1bSJed Brown       if (i >= 0 && j >= 0) v = x[k];
250c4762a1bSJed Brown       if (i < mx-1 && j > -1) vr = x[k+1];
251c4762a1bSJed Brown       if (i > -1 && j < my-1) vt = x[k+gxm];
252c4762a1bSJed Brown       dvdx = (vr-v)/hx;
253c4762a1bSJed Brown       dvdy = (vt-v)/hy;
254c4762a1bSJed Brown       if (i != -1 && j != -1) {
255c4762a1bSJed Brown         ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
2569566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G,1,&k,&val,ADD_VALUES));
257c4762a1bSJed Brown       }
258c4762a1bSJed Brown       if (i != mx-1 && j != -1) {
259c4762a1bSJed Brown         ind = k+1; val =  dvdx/hx - cdiv3;
2609566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES));
261c4762a1bSJed Brown       }
262c4762a1bSJed Brown       if (i != -1 && j != my-1) {
263c4762a1bSJed Brown         ind = k+gxm; val = dvdy/hy - cdiv3;
2649566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES));
265c4762a1bSJed Brown       }
266c4762a1bSJed Brown       fquad += dvdx*dvdx + dvdy*dvdy;
267c4762a1bSJed Brown       flin -= cdiv3 * (v + vr + vt);
268c4762a1bSJed Brown     }
269c4762a1bSJed Brown   }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* Compute local gradient contributions over the upper triangular elements */
272c4762a1bSJed Brown   for (j=ys; j<yep; j++) { /*  for (j=0; j<=my; j++) */
273c4762a1bSJed Brown     for (i=xs; i<xep; i++) {  /*   for (i=0; i<=mx; i++) */
274c4762a1bSJed Brown       k = (j-gys)*gxm + i-gxs;
275c4762a1bSJed Brown       vb = zero;
276c4762a1bSJed Brown       vl = zero;
277c4762a1bSJed Brown       v  = zero;
278c4762a1bSJed Brown       if (i < mx && j > 0) vb = x[k-gxm];
279c4762a1bSJed Brown       if (i > 0 && j < my) vl = x[k-1];
280c4762a1bSJed Brown       if (i < mx && j < my) v = x[k];
281c4762a1bSJed Brown       dvdx = (v-vl)/hx;
282c4762a1bSJed Brown       dvdy = (v-vb)/hy;
283c4762a1bSJed Brown       if (i != mx && j != 0) {
284c4762a1bSJed Brown         ind = k-gxm; val = - dvdy/hy - cdiv3;
2859566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES));
286c4762a1bSJed Brown       }
287c4762a1bSJed Brown       if (i != 0 && j != my) {
288c4762a1bSJed Brown         ind = k-1; val =  - dvdx/hx - cdiv3;
2899566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES));
290c4762a1bSJed Brown       }
291c4762a1bSJed Brown       if (i != mx && j != my) {
292c4762a1bSJed Brown         ind = k; val =  dvdx/hx + dvdy/hy - cdiv3;
2939566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES));
294c4762a1bSJed Brown       }
295c4762a1bSJed Brown       fquad += dvdx*dvdx + dvdy*dvdy;
296c4762a1bSJed Brown       flin -= cdiv3 * (vb + vl + v);
297c4762a1bSJed Brown     }
298c4762a1bSJed Brown   }
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   /* Restore vector */
3019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX,&x));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   /* Assemble gradient vector */
3049566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(G));
3059566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(G));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /* Scale the gradient */
308c4762a1bSJed Brown   area = p5*hx*hy;
309c4762a1bSJed Brown   floc = area * (p5 * fquad + flin);
3109566063dSJacob Faibussowitsch   PetscCall(VecScale(G, area));
311c4762a1bSJed Brown 
3129566063dSJacob Faibussowitsch   /* Sum function contributions from all processes */  /* TODO: Change to PetscCallMPI() */
3139566063dSJacob Faibussowitsch   PetscCall((PetscErrorCode)MPI_Allreduce((void*)&floc,(void*)f,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD));
314c4762a1bSJed Brown 
3159566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((ye-ysm)*(xe-xsm)*20+(xep-xs)*(yep-ys)*16));
316c4762a1bSJed Brown   PetscFunctionReturn(0);
317c4762a1bSJed Brown }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao, Vec X, Mat A, Mat Hpre, void*ctx)
320c4762a1bSJed Brown {
321c4762a1bSJed Brown   AppCtx         *user= (AppCtx*) ctx;
322c4762a1bSJed Brown   PetscInt       i,j,k;
323c4762a1bSJed Brown   PetscInt       col[5],row;
324c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
325c4762a1bSJed Brown   PetscReal      v[5];
326c4762a1bSJed 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;
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   /* Compute the quadratic term in the objective function */
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   /*
331c4762a1bSJed Brown      Get local grid boundaries
332c4762a1bSJed Brown   */
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   PetscFunctionBegin;
3359566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
3369566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
339c4762a1bSJed Brown 
340c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
341c4762a1bSJed Brown 
342c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
343c4762a1bSJed Brown 
344c4762a1bSJed Brown       k=0;
345c4762a1bSJed Brown       if (j>gys) {
346c4762a1bSJed Brown         v[k]=-2*hyhy; col[k]=row - gxm; k++;
347c4762a1bSJed Brown       }
348c4762a1bSJed Brown 
349c4762a1bSJed Brown       if (i>gxs) {
350c4762a1bSJed Brown         v[k]= -2*hxhx; col[k]=row - 1; k++;
351c4762a1bSJed Brown       }
352c4762a1bSJed Brown 
353c4762a1bSJed Brown       v[k]= 4.0*(hxhx+hyhy); col[k]=row; k++;
354c4762a1bSJed Brown 
355c4762a1bSJed Brown       if (i+1 < gxs+gxm) {
356c4762a1bSJed Brown         v[k]= -2.0*hxhx; col[k]=row+1; k++;
357c4762a1bSJed Brown       }
358c4762a1bSJed Brown 
359c4762a1bSJed Brown       if (j+1 <gys+gym) {
360c4762a1bSJed Brown         v[k]= -2*hyhy; col[k] = row+gxm; k++;
361c4762a1bSJed Brown       }
362c4762a1bSJed Brown 
3639566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(A,1,&row,k,col,v,INSERT_VALUES));
364c4762a1bSJed Brown 
365c4762a1bSJed Brown     }
366c4762a1bSJed Brown   }
367c4762a1bSJed Brown   /*
368c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
369c4762a1bSJed Brown      MatAssemblyBegin(), MatAssemblyEnd().
370c4762a1bSJed Brown      By placing code between these two statements, computations can be
371c4762a1bSJed Brown      done while messages are in transition.
372c4762a1bSJed Brown   */
3739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
3749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
375c4762a1bSJed Brown   /*
376c4762a1bSJed Brown     Tell the matrix we will never add a new nonzero location to the
377c4762a1bSJed Brown     matrix. If we do it will generate an error.
378c4762a1bSJed Brown   */
3799566063dSJacob Faibussowitsch   PetscCall(MatScale(A,area));
3809566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
3819566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
3829566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(9*xm*ym+49*xm));
383c4762a1bSJed Brown   PetscFunctionReturn(0);
384c4762a1bSJed Brown }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown /*TEST
387c4762a1bSJed Brown 
388c4762a1bSJed Brown    build:
389c4762a1bSJed Brown       requires: !complex
390c4762a1bSJed Brown 
391c4762a1bSJed Brown    test:
392c4762a1bSJed Brown       suffix: 1
393c4762a1bSJed Brown       nsize: 2
394c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -mx 16 -my 16 -tao_gatol 1.e-4
395c4762a1bSJed Brown 
396c4762a1bSJed Brown TEST*/
397