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
369371c9d4SSatish 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
main(int argc,char ** argv)63d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
64d71ae5a4SJacob Faibussowitsch {
65c4762a1bSJed Brown Vec x;
66c4762a1bSJed Brown Mat H;
67c4762a1bSJed Brown PetscInt Nx, Ny;
68c4762a1bSJed Brown Tao tao;
69c4762a1bSJed Brown PetscBool flg;
70c4762a1bSJed Brown KSP ksp;
71c4762a1bSJed Brown PC pc;
72c4762a1bSJed Brown AppCtx user;
73c4762a1bSJed Brown
74c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
75c4762a1bSJed Brown
76c4762a1bSJed Brown /* Specify default dimension of the problem */
779371c9d4SSatish Balay user.param = 5.0;
789371c9d4SSatish Balay user.mx = 10;
799371c9d4SSatish Balay user.my = 10;
80c4762a1bSJed Brown Nx = Ny = PETSC_DECIDE;
81c4762a1bSJed Brown
82c4762a1bSJed Brown /* Check for any command line arguments that override defaults */
839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, &flg));
849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg));
859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg));
86c4762a1bSJed Brown
879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n---- Elastic-Plastic Torsion Problem -----\n"));
8863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n", user.mx, user.my));
89c4762a1bSJed Brown
90c4762a1bSJed Brown /* Set up distributed array */
919566063dSJacob 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));
929566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm));
939566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm));
94c4762a1bSJed Brown
95c4762a1bSJed Brown /* Create vectors */
969566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm, &x));
97c4762a1bSJed Brown
989566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(user.dm, &user.localX));
99c4762a1bSJed Brown
100c4762a1bSJed Brown /* Create Hessian */
1019566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm, &H));
1029566063dSJacob Faibussowitsch PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE));
103c4762a1bSJed Brown
104c4762a1bSJed Brown /* The TAO code begins here */
105c4762a1bSJed Brown
106c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
1079566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1089566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOCG));
109c4762a1bSJed Brown
110c4762a1bSJed Brown /* Set initial solution guess */
1119566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user, x));
1129566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
113c4762a1bSJed Brown
114c4762a1bSJed Brown /* Set routine for function and gradient evaluation */
1159566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
116c4762a1bSJed Brown
1179566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, FormHessian, (void *)&user));
118c4762a1bSJed Brown
119c4762a1bSJed Brown /* Check for any TAO command line options */
1209566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
121c4762a1bSJed Brown
1229566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp));
123c4762a1bSJed Brown if (ksp) {
1249566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc));
1259566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE));
126c4762a1bSJed Brown }
127c4762a1bSJed Brown
128c4762a1bSJed Brown /* SOLVE THE APPLICATION */
1299566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
130c4762a1bSJed Brown
131c4762a1bSJed Brown /* Free TAO data structures */
1329566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
133c4762a1bSJed Brown
134c4762a1bSJed Brown /* Free PETSc data structures */
1359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H));
137c4762a1bSJed Brown
1389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.localX));
1399566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm));
140c4762a1bSJed Brown
1413ba16761SJacob Faibussowitsch PetscCall(PetscFinalize());
142c4762a1bSJed Brown return 0;
143c4762a1bSJed Brown }
144c4762a1bSJed Brown
145c4762a1bSJed Brown /* ------------------------------------------------------------------- */
146c4762a1bSJed Brown /*
147c4762a1bSJed Brown FormInitialGuess - Computes an initial approximation to the solution.
148c4762a1bSJed Brown
149c4762a1bSJed Brown Input Parameters:
150c4762a1bSJed Brown . user - user-defined application context
151c4762a1bSJed Brown . X - vector
152c4762a1bSJed Brown
153c4762a1bSJed Brown Output Parameters:
154c4762a1bSJed Brown X - vector
155c4762a1bSJed Brown */
FormInitialGuess(AppCtx * user,Vec X)156d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X)
157d71ae5a4SJacob Faibussowitsch {
158c4762a1bSJed Brown PetscInt i, j, k, mx = user->mx, my = user->my;
159c4762a1bSJed Brown PetscInt xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye;
160c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), temp, val;
161c4762a1bSJed Brown
162c4762a1bSJed Brown PetscFunctionBegin;
163c4762a1bSJed Brown /* Get local mesh boundaries */
1649566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
1659566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
166c4762a1bSJed Brown
167c4762a1bSJed Brown /* Compute initial guess over locally owned part of mesh */
168c4762a1bSJed Brown xe = xs + xm;
169c4762a1bSJed Brown ye = ys + ym;
170c4762a1bSJed Brown for (j = ys; j < ye; j++) { /* for (j=0; j<my; j++) */
171c4762a1bSJed Brown temp = PetscMin(j + 1, my - j) * hy;
172c4762a1bSJed Brown for (i = xs; i < xe; i++) { /* for (i=0; i<mx; i++) */
173c4762a1bSJed Brown k = (j - gys) * gxm + i - gxs;
174c4762a1bSJed Brown val = PetscMin((PetscMin(i + 1, mx - i)) * hx, temp);
1759566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(X, 1, &k, &val, ADD_VALUES));
176c4762a1bSJed Brown }
177c4762a1bSJed Brown }
1789566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X));
1799566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X));
1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
181c4762a1bSJed Brown }
182c4762a1bSJed Brown
183c4762a1bSJed Brown /* ------------------------------------------------------------------ */
184c4762a1bSJed Brown /*
185c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient.
186c4762a1bSJed Brown
187c4762a1bSJed Brown Input Parameters:
188c4762a1bSJed Brown tao - the Tao context
189c4762a1bSJed Brown X - the input vector
190a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
191c4762a1bSJed Brown
192c4762a1bSJed Brown Output Parameters:
193c4762a1bSJed Brown f - the newly evaluated function
194c4762a1bSJed Brown G - the newly evaluated gradient
195c4762a1bSJed Brown */
FormFunctionGradient(Tao tao,Vec X,PetscReal * f,Vec G,void * ptr)196d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
197d71ae5a4SJacob Faibussowitsch {
198c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
199c4762a1bSJed Brown PetscInt i, j, k, ind;
200c4762a1bSJed Brown PetscInt xe, ye, xsm, ysm, xep, yep;
201c4762a1bSJed Brown PetscInt xs, ys, xm, ym, gxm, gym, gxs, gys;
202c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
203c4762a1bSJed Brown PetscReal three = 3.0, zero = 0.0, *x, floc, cdiv3 = user->param / three;
204c4762a1bSJed Brown PetscReal p5 = 0.5, area, val, flin, fquad;
205c4762a1bSJed Brown PetscReal v, vb, vl, vr, vt, dvdx, dvdy;
206c4762a1bSJed Brown PetscReal hx = 1.0 / (user->mx + 1);
207c4762a1bSJed Brown PetscReal hy = 1.0 / (user->my + 1);
208c4762a1bSJed Brown Vec localX = user->localX;
209c4762a1bSJed Brown
210c4762a1bSJed Brown PetscFunctionBegin;
211c4762a1bSJed Brown /* Initialize */
212c4762a1bSJed Brown flin = fquad = zero;
213c4762a1bSJed Brown
2149566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero));
215c4762a1bSJed Brown /*
216c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
217c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
218c4762a1bSJed Brown By placing code between these two statements, computations can be
219c4762a1bSJed Brown done while messages are in transition.
220c4762a1bSJed Brown */
2219566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
2229566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
223c4762a1bSJed Brown
224c4762a1bSJed Brown /* Get pointer to vector data */
2259566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x));
226c4762a1bSJed Brown
227c4762a1bSJed Brown /* Get local mesh boundaries */
2289566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
2299566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
230c4762a1bSJed Brown
231c4762a1bSJed Brown /* Set local loop dimensions */
232c4762a1bSJed Brown xe = xs + xm;
233c4762a1bSJed Brown ye = ys + ym;
234c4762a1bSJed Brown if (xs == 0) xsm = xs - 1;
235c4762a1bSJed Brown else xsm = xs;
236c4762a1bSJed Brown if (ys == 0) ysm = ys - 1;
237c4762a1bSJed Brown else ysm = ys;
238c4762a1bSJed Brown if (xe == mx) xep = xe + 1;
239c4762a1bSJed Brown else xep = xe;
240c4762a1bSJed Brown if (ye == my) yep = ye + 1;
241c4762a1bSJed Brown else yep = ye;
242c4762a1bSJed Brown
243c4762a1bSJed Brown /* Compute local gradient contributions over the lower triangular elements */
244c4762a1bSJed Brown for (j = ysm; j < ye; j++) { /* for (j=-1; j<my; j++) */
245c4762a1bSJed Brown for (i = xsm; i < xe; i++) { /* for (i=-1; i<mx; i++) */
246c4762a1bSJed Brown k = (j - gys) * gxm + i - gxs;
247c4762a1bSJed Brown v = zero;
248c4762a1bSJed Brown vr = zero;
249c4762a1bSJed Brown vt = zero;
250c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k];
251c4762a1bSJed Brown if (i < mx - 1 && j > -1) vr = x[k + 1];
252c4762a1bSJed Brown if (i > -1 && j < my - 1) vt = x[k + gxm];
253c4762a1bSJed Brown dvdx = (vr - v) / hx;
254c4762a1bSJed Brown dvdy = (vt - v) / hy;
255c4762a1bSJed Brown if (i != -1 && j != -1) {
2569371c9d4SSatish Balay ind = k;
2579371c9d4SSatish Balay val = -dvdx / hx - dvdy / hy - cdiv3;
2589566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G, 1, &k, &val, ADD_VALUES));
259c4762a1bSJed Brown }
260c4762a1bSJed Brown if (i != mx - 1 && j != -1) {
2619371c9d4SSatish Balay ind = k + 1;
2629371c9d4SSatish Balay val = dvdx / hx - cdiv3;
2639566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
264c4762a1bSJed Brown }
265c4762a1bSJed Brown if (i != -1 && j != my - 1) {
2669371c9d4SSatish Balay ind = k + gxm;
2679371c9d4SSatish Balay val = dvdy / hy - cdiv3;
2689566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
269c4762a1bSJed Brown }
270c4762a1bSJed Brown fquad += dvdx * dvdx + dvdy * dvdy;
271c4762a1bSJed Brown flin -= cdiv3 * (v + vr + vt);
272c4762a1bSJed Brown }
273c4762a1bSJed Brown }
274c4762a1bSJed Brown
275c4762a1bSJed Brown /* Compute local gradient contributions over the upper triangular elements */
276c4762a1bSJed Brown for (j = ys; j < yep; j++) { /* for (j=0; j<=my; j++) */
277c4762a1bSJed Brown for (i = xs; i < xep; i++) { /* for (i=0; i<=mx; i++) */
278c4762a1bSJed Brown k = (j - gys) * gxm + i - gxs;
279c4762a1bSJed Brown vb = zero;
280c4762a1bSJed Brown vl = zero;
281c4762a1bSJed Brown v = zero;
282c4762a1bSJed Brown if (i < mx && j > 0) vb = x[k - gxm];
283c4762a1bSJed Brown if (i > 0 && j < my) vl = x[k - 1];
284c4762a1bSJed Brown if (i < mx && j < my) v = x[k];
285c4762a1bSJed Brown dvdx = (v - vl) / hx;
286c4762a1bSJed Brown dvdy = (v - vb) / hy;
287c4762a1bSJed Brown if (i != mx && j != 0) {
2889371c9d4SSatish Balay ind = k - gxm;
2899371c9d4SSatish Balay val = -dvdy / hy - cdiv3;
2909566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
291c4762a1bSJed Brown }
292c4762a1bSJed Brown if (i != 0 && j != my) {
2939371c9d4SSatish Balay ind = k - 1;
2949371c9d4SSatish Balay val = -dvdx / hx - cdiv3;
2959566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
296c4762a1bSJed Brown }
297c4762a1bSJed Brown if (i != mx && j != my) {
2989371c9d4SSatish Balay ind = k;
2999371c9d4SSatish Balay val = dvdx / hx + dvdy / hy - cdiv3;
3009566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G, 1, &ind, &val, ADD_VALUES));
301c4762a1bSJed Brown }
302c4762a1bSJed Brown fquad += dvdx * dvdx + dvdy * dvdy;
303c4762a1bSJed Brown flin -= cdiv3 * (vb + vl + v);
304c4762a1bSJed Brown }
305c4762a1bSJed Brown }
306c4762a1bSJed Brown
307c4762a1bSJed Brown /* Restore vector */
3089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x));
309c4762a1bSJed Brown
310c4762a1bSJed Brown /* Assemble gradient vector */
3119566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(G));
3129566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(G));
313c4762a1bSJed Brown
314c4762a1bSJed Brown /* Scale the gradient */
315c4762a1bSJed Brown area = p5 * hx * hy;
316c4762a1bSJed Brown floc = area * (p5 * fquad + flin);
3179566063dSJacob Faibussowitsch PetscCall(VecScale(G, area));
318c4762a1bSJed Brown
3199566063dSJacob Faibussowitsch /* Sum function contributions from all processes */ /* TODO: Change to PetscCallMPI() */
3209566063dSJacob Faibussowitsch PetscCall((PetscErrorCode)MPI_Allreduce((void *)&floc, (void *)f, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
321c4762a1bSJed Brown
3229566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((ye - ysm) * (xe - xsm) * 20 + (xep - xs) * (yep - ys) * 16));
3233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
324c4762a1bSJed Brown }
325c4762a1bSJed Brown
FormHessian(Tao tao,Vec X,Mat A,Mat Hpre,PetscCtx ctx)326*2a8381b2SBarry Smith PetscErrorCode FormHessian(Tao tao, Vec X, Mat A, Mat Hpre, PetscCtx ctx)
327d71ae5a4SJacob Faibussowitsch {
328c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
329c4762a1bSJed Brown PetscInt i, j, k;
330c4762a1bSJed Brown PetscInt col[5], row;
331c4762a1bSJed Brown PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym;
332c4762a1bSJed Brown PetscReal v[5];
333c4762a1bSJed 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;
334c4762a1bSJed Brown
335c4762a1bSJed Brown /* Compute the quadratic term in the objective function */
336c4762a1bSJed Brown
337c4762a1bSJed Brown /*
338c4762a1bSJed Brown Get local grid boundaries
339c4762a1bSJed Brown */
340c4762a1bSJed Brown
341c4762a1bSJed Brown PetscFunctionBegin;
3429566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
3439566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
344c4762a1bSJed Brown
345c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
346c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
347c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs);
348c4762a1bSJed Brown
349c4762a1bSJed Brown k = 0;
350c4762a1bSJed Brown if (j > gys) {
3519371c9d4SSatish Balay v[k] = -2 * hyhy;
3529371c9d4SSatish Balay col[k] = row - gxm;
3539371c9d4SSatish Balay k++;
354c4762a1bSJed Brown }
355c4762a1bSJed Brown
356c4762a1bSJed Brown if (i > gxs) {
3579371c9d4SSatish Balay v[k] = -2 * hxhx;
3589371c9d4SSatish Balay col[k] = row - 1;
3599371c9d4SSatish Balay k++;
360c4762a1bSJed Brown }
361c4762a1bSJed Brown
3629371c9d4SSatish Balay v[k] = 4.0 * (hxhx + hyhy);
3639371c9d4SSatish Balay col[k] = row;
3649371c9d4SSatish Balay k++;
365c4762a1bSJed Brown
366c4762a1bSJed Brown if (i + 1 < gxs + gxm) {
3679371c9d4SSatish Balay v[k] = -2.0 * hxhx;
3689371c9d4SSatish Balay col[k] = row + 1;
3699371c9d4SSatish Balay k++;
370c4762a1bSJed Brown }
371c4762a1bSJed Brown
372c4762a1bSJed Brown if (j + 1 < gys + gym) {
3739371c9d4SSatish Balay v[k] = -2 * hyhy;
3749371c9d4SSatish Balay col[k] = row + gxm;
3759371c9d4SSatish Balay k++;
376c4762a1bSJed Brown }
377c4762a1bSJed Brown
3789566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, 1, &row, k, col, v, INSERT_VALUES));
379c4762a1bSJed Brown }
380c4762a1bSJed Brown }
381c4762a1bSJed Brown /*
382c4762a1bSJed Brown Assemble matrix, using the 2-step process:
383c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd().
384c4762a1bSJed Brown By placing code between these two statements, computations can be
385c4762a1bSJed Brown done while messages are in transition.
386c4762a1bSJed Brown */
3879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
389c4762a1bSJed Brown /*
390c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the
391c4762a1bSJed Brown matrix. If we do it will generate an error.
392c4762a1bSJed Brown */
3939566063dSJacob Faibussowitsch PetscCall(MatScale(A, area));
3949566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
3959566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
3969566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9 * xm * ym + 49 * xm));
3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
398c4762a1bSJed Brown }
399c4762a1bSJed Brown
400c4762a1bSJed Brown /*TEST
401c4762a1bSJed Brown
402c4762a1bSJed Brown build:
403c4762a1bSJed Brown requires: !complex
404c4762a1bSJed Brown
405c4762a1bSJed Brown test:
406c4762a1bSJed Brown suffix: 1
407c4762a1bSJed Brown nsize: 2
40810978b7dSBarry Smith args: -tao_monitor_short -tao_type nls -mx 16 -my 16 -tao_gatol 1.e-4
409c4762a1bSJed Brown
410c4762a1bSJed Brown TEST*/
411