xref: /petsc/src/tao/unconstrained/tutorials/eptorsion2.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
36*9371c9d4SSatish Balay static char help[] = "Demonstrates use of the TAO package to solve \n\
37c4762a1bSJed Brown unconstrained minimization problems in parallel.  This example is based on \n\
38c4762a1bSJed Brown the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 test suite.\n\
39c4762a1bSJed Brown The command line options are:\n\
40c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
41c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
42c4762a1bSJed Brown   -par <param>, where <param> = angle of twist per unit length\n\n";
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*
45c4762a1bSJed Brown    User-defined application context - contains data needed by the
46c4762a1bSJed Brown    application-provided call-back routines, FormFunction() and
47c4762a1bSJed Brown    FormGradient().
48c4762a1bSJed Brown */
49c4762a1bSJed Brown typedef struct {
50c4762a1bSJed Brown   /* parameters */
51c4762a1bSJed Brown   PetscInt  mx, my; /* global discretization in x- and y-directions */
52c4762a1bSJed Brown   PetscReal param;  /* nonlinearity parameter */
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* work space */
55c4762a1bSJed Brown   Vec localX; /* local vectors */
56c4762a1bSJed Brown   DM  dm;     /* distributed array data structure */
57c4762a1bSJed Brown } AppCtx;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *, Vec);
60c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
61c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
62c4762a1bSJed Brown 
63*9371c9d4SSatish Balay int main(int argc, char **argv) {
64c4762a1bSJed Brown   Vec       x;
65c4762a1bSJed Brown   Mat       H;
66c4762a1bSJed Brown   PetscInt  Nx, Ny;
67c4762a1bSJed Brown   Tao       tao;
68c4762a1bSJed Brown   PetscBool flg;
69c4762a1bSJed Brown   KSP       ksp;
70c4762a1bSJed Brown   PC        pc;
71c4762a1bSJed Brown   AppCtx    user;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscInitialize(&argc, &argv, (char *)0, help);
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* Specify default dimension of the problem */
76*9371c9d4SSatish Balay   user.param = 5.0;
77*9371c9d4SSatish Balay   user.mx    = 10;
78*9371c9d4SSatish Balay   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"));
8763a3b9bcSJacob 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 */
155*9371c9d4SSatish Balay PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) {
156c4762a1bSJed Brown   PetscInt  i, j, k, mx = user->mx, my = user->my;
157c4762a1bSJed Brown   PetscInt  xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye;
158c4762a1bSJed Brown   PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), temp, val;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBegin;
161c4762a1bSJed Brown   /* Get local mesh boundaries */
1629566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
1639566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* Compute initial guess over locally owned part of mesh */
166c4762a1bSJed Brown   xe = xs + xm;
167c4762a1bSJed Brown   ye = ys + ym;
168c4762a1bSJed Brown   for (j = ys; j < ye; j++) { /*  for (j=0; j<my; j++) */
169c4762a1bSJed Brown     temp = PetscMin(j + 1, my - j) * hy;
170c4762a1bSJed Brown     for (i = xs; i < xe; i++) { /*  for (i=0; i<mx; i++) */
171c4762a1bSJed Brown       k   = (j - gys) * gxm + i - gxs;
172c4762a1bSJed Brown       val = PetscMin((PetscMin(i + 1, mx - i)) * hx, temp);
1739566063dSJacob Faibussowitsch       PetscCall(VecSetValuesLocal(X, 1, &k, &val, ADD_VALUES));
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown   }
1769566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
1779566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
178c4762a1bSJed Brown   PetscFunctionReturn(0);
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown /* ------------------------------------------------------------------ */
182c4762a1bSJed Brown /*
183c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    Input Parameters:
186c4762a1bSJed Brown    tao - the Tao context
187c4762a1bSJed Brown    X   - the input vector
188a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    Output Parameters:
191c4762a1bSJed Brown    f   - the newly evaluated function
192c4762a1bSJed Brown    G   - the newly evaluated gradient
193c4762a1bSJed Brown */
194*9371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) {
195c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ptr;
196c4762a1bSJed Brown   PetscInt  i, j, k, ind;
197c4762a1bSJed Brown   PetscInt  xe, ye, xsm, ysm, xep, yep;
198c4762a1bSJed Brown   PetscInt  xs, ys, xm, ym, gxm, gym, gxs, gys;
199c4762a1bSJed Brown   PetscInt  mx = user->mx, my = user->my;
200c4762a1bSJed Brown   PetscReal three = 3.0, zero = 0.0, *x, floc, cdiv3 = user->param / three;
201c4762a1bSJed Brown   PetscReal p5 = 0.5, area, val, flin, fquad;
202c4762a1bSJed Brown   PetscReal v, vb, vl, vr, vt, dvdx, dvdy;
203c4762a1bSJed Brown   PetscReal hx     = 1.0 / (user->mx + 1);
204c4762a1bSJed Brown   PetscReal hy     = 1.0 / (user->my + 1);
205c4762a1bSJed Brown   Vec       localX = user->localX;
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   PetscFunctionBegin;
208c4762a1bSJed Brown   /* Initialize */
209c4762a1bSJed Brown   flin = fquad = zero;
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch   PetscCall(VecSet(G, zero));
212c4762a1bSJed Brown   /*
213c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
214c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
215c4762a1bSJed Brown      By placing code between these two statements, computations can be
216c4762a1bSJed Brown      done while messages are in transition.
217c4762a1bSJed Brown   */
2189566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
2199566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Get pointer to vector data */
2229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &x));
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   /* Get local mesh boundaries */
2259566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
2269566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* Set local loop dimensions */
229c4762a1bSJed Brown   xe = xs + xm;
230c4762a1bSJed Brown   ye = ys + ym;
231c4762a1bSJed Brown   if (xs == 0) xsm = xs - 1;
232c4762a1bSJed Brown   else xsm = xs;
233c4762a1bSJed Brown   if (ys == 0) ysm = ys - 1;
234c4762a1bSJed Brown   else ysm = ys;
235c4762a1bSJed Brown   if (xe == mx) xep = xe + 1;
236c4762a1bSJed Brown   else xep = xe;
237c4762a1bSJed Brown   if (ye == my) yep = ye + 1;
238c4762a1bSJed Brown   else yep = ye;
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   /* Compute local gradient contributions over the lower triangular elements */
241c4762a1bSJed Brown   for (j = ysm; j < ye; j++) {   /*  for (j=-1; j<my; j++) */
242c4762a1bSJed Brown     for (i = xsm; i < xe; i++) { /*   for (i=-1; i<mx; i++) */
243c4762a1bSJed Brown       k  = (j - gys) * gxm + i - gxs;
244c4762a1bSJed Brown       v  = zero;
245c4762a1bSJed Brown       vr = zero;
246c4762a1bSJed Brown       vt = zero;
247c4762a1bSJed Brown       if (i >= 0 && j >= 0) v = x[k];
248c4762a1bSJed Brown       if (i < mx - 1 && j > -1) vr = x[k + 1];
249c4762a1bSJed Brown       if (i > -1 && j < my - 1) vt = x[k + gxm];
250c4762a1bSJed Brown       dvdx = (vr - v) / hx;
251c4762a1bSJed Brown       dvdy = (vt - v) / hy;
252c4762a1bSJed Brown       if (i != -1 && j != -1) {
253*9371c9d4SSatish Balay         ind = k;
254*9371c9d4SSatish Balay         val = -dvdx / hx - dvdy / hy - cdiv3;
2559566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G, 1, &k, &val, ADD_VALUES));
256c4762a1bSJed Brown       }
257c4762a1bSJed Brown       if (i != mx - 1 && j != -1) {
258*9371c9d4SSatish Balay         ind = k + 1;
259*9371c9d4SSatish Balay         val = dvdx / hx - cdiv3;
2609566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
261c4762a1bSJed Brown       }
262c4762a1bSJed Brown       if (i != -1 && j != my - 1) {
263*9371c9d4SSatish Balay         ind = k + gxm;
264*9371c9d4SSatish Balay         val = dvdy / hy - cdiv3;
2659566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
266c4762a1bSJed Brown       }
267c4762a1bSJed Brown       fquad += dvdx * dvdx + dvdy * dvdy;
268c4762a1bSJed Brown       flin -= cdiv3 * (v + vr + vt);
269c4762a1bSJed Brown     }
270c4762a1bSJed Brown   }
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   /* Compute local gradient contributions over the upper triangular elements */
273c4762a1bSJed Brown   for (j = ys; j < yep; j++) {   /*  for (j=0; j<=my; j++) */
274c4762a1bSJed Brown     for (i = xs; i < xep; i++) { /*   for (i=0; i<=mx; i++) */
275c4762a1bSJed Brown       k  = (j - gys) * gxm + i - gxs;
276c4762a1bSJed Brown       vb = zero;
277c4762a1bSJed Brown       vl = zero;
278c4762a1bSJed Brown       v  = zero;
279c4762a1bSJed Brown       if (i < mx && j > 0) vb = x[k - gxm];
280c4762a1bSJed Brown       if (i > 0 && j < my) vl = x[k - 1];
281c4762a1bSJed Brown       if (i < mx && j < my) v = x[k];
282c4762a1bSJed Brown       dvdx = (v - vl) / hx;
283c4762a1bSJed Brown       dvdy = (v - vb) / hy;
284c4762a1bSJed Brown       if (i != mx && j != 0) {
285*9371c9d4SSatish Balay         ind = k - gxm;
286*9371c9d4SSatish Balay         val = -dvdy / hy - cdiv3;
2879566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
288c4762a1bSJed Brown       }
289c4762a1bSJed Brown       if (i != 0 && j != my) {
290*9371c9d4SSatish Balay         ind = k - 1;
291*9371c9d4SSatish Balay         val = -dvdx / hx - cdiv3;
2929566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
293c4762a1bSJed Brown       }
294c4762a1bSJed Brown       if (i != mx && j != my) {
295*9371c9d4SSatish Balay         ind = k;
296*9371c9d4SSatish Balay         val = dvdx / hx + dvdy / hy - cdiv3;
2979566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
298c4762a1bSJed Brown       }
299c4762a1bSJed Brown       fquad += dvdx * dvdx + dvdy * dvdy;
300c4762a1bSJed Brown       flin -= cdiv3 * (vb + vl + v);
301c4762a1bSJed Brown     }
302c4762a1bSJed Brown   }
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /* Restore vector */
3059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &x));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /* Assemble gradient vector */
3089566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(G));
3099566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(G));
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   /* Scale the gradient */
312c4762a1bSJed Brown   area = p5 * hx * hy;
313c4762a1bSJed Brown   floc = area * (p5 * fquad + flin);
3149566063dSJacob Faibussowitsch   PetscCall(VecScale(G, area));
315c4762a1bSJed Brown 
3169566063dSJacob Faibussowitsch   /* Sum function contributions from all processes */ /* TODO: Change to PetscCallMPI() */
3179566063dSJacob Faibussowitsch   PetscCall((PetscErrorCode)MPI_Allreduce((void *)&floc, (void *)f, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
318c4762a1bSJed Brown 
3199566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((ye - ysm) * (xe - xsm) * 20 + (xep - xs) * (yep - ys) * 16));
320c4762a1bSJed Brown   PetscFunctionReturn(0);
321c4762a1bSJed Brown }
322c4762a1bSJed Brown 
323*9371c9d4SSatish Balay PetscErrorCode FormHessian(Tao tao, Vec X, Mat A, Mat Hpre, void *ctx) {
324c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ctx;
325c4762a1bSJed Brown   PetscInt  i, j, k;
326c4762a1bSJed Brown   PetscInt  col[5], row;
327c4762a1bSJed Brown   PetscInt  xs, xm, gxs, gxm, ys, ym, gys, gym;
328c4762a1bSJed Brown   PetscReal v[5];
329c4762a1bSJed 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;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   /* Compute the quadratic term in the objective function */
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   /*
334c4762a1bSJed Brown      Get local grid boundaries
335c4762a1bSJed Brown   */
336c4762a1bSJed Brown 
337c4762a1bSJed Brown   PetscFunctionBegin;
3389566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
3399566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
342c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
343c4762a1bSJed Brown       row = (j - gys) * gxm + (i - gxs);
344c4762a1bSJed Brown 
345c4762a1bSJed Brown       k = 0;
346c4762a1bSJed Brown       if (j > gys) {
347*9371c9d4SSatish Balay         v[k]   = -2 * hyhy;
348*9371c9d4SSatish Balay         col[k] = row - gxm;
349*9371c9d4SSatish Balay         k++;
350c4762a1bSJed Brown       }
351c4762a1bSJed Brown 
352c4762a1bSJed Brown       if (i > gxs) {
353*9371c9d4SSatish Balay         v[k]   = -2 * hxhx;
354*9371c9d4SSatish Balay         col[k] = row - 1;
355*9371c9d4SSatish Balay         k++;
356c4762a1bSJed Brown       }
357c4762a1bSJed Brown 
358*9371c9d4SSatish Balay       v[k]   = 4.0 * (hxhx + hyhy);
359*9371c9d4SSatish Balay       col[k] = row;
360*9371c9d4SSatish Balay       k++;
361c4762a1bSJed Brown 
362c4762a1bSJed Brown       if (i + 1 < gxs + gxm) {
363*9371c9d4SSatish Balay         v[k]   = -2.0 * hxhx;
364*9371c9d4SSatish Balay         col[k] = row + 1;
365*9371c9d4SSatish Balay         k++;
366c4762a1bSJed Brown       }
367c4762a1bSJed Brown 
368c4762a1bSJed Brown       if (j + 1 < gys + gym) {
369*9371c9d4SSatish Balay         v[k]   = -2 * hyhy;
370*9371c9d4SSatish Balay         col[k] = row + gxm;
371*9371c9d4SSatish Balay         k++;
372c4762a1bSJed Brown       }
373c4762a1bSJed Brown 
3749566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(A, 1, &row, k, col, v, INSERT_VALUES));
375c4762a1bSJed Brown     }
376c4762a1bSJed Brown   }
377c4762a1bSJed Brown   /*
378c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
379c4762a1bSJed Brown      MatAssemblyBegin(), MatAssemblyEnd().
380c4762a1bSJed Brown      By placing code between these two statements, computations can be
381c4762a1bSJed Brown      done while messages are in transition.
382c4762a1bSJed Brown   */
3839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
385c4762a1bSJed Brown   /*
386c4762a1bSJed Brown     Tell the matrix we will never add a new nonzero location to the
387c4762a1bSJed Brown     matrix. If we do it will generate an error.
388c4762a1bSJed Brown   */
3899566063dSJacob Faibussowitsch   PetscCall(MatScale(A, area));
3909566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
3919566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
3929566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(9 * xm * ym + 49 * xm));
393c4762a1bSJed Brown   PetscFunctionReturn(0);
394c4762a1bSJed Brown }
395c4762a1bSJed Brown 
396c4762a1bSJed Brown /*TEST
397c4762a1bSJed Brown 
398c4762a1bSJed Brown    build:
399c4762a1bSJed Brown       requires: !complex
400c4762a1bSJed Brown 
401c4762a1bSJed Brown    test:
402c4762a1bSJed Brown       suffix: 1
403c4762a1bSJed Brown       nsize: 2
404c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -mx 16 -my 16 -tao_gatol 1.e-4
405c4762a1bSJed Brown 
406c4762a1bSJed Brown TEST*/
407