1c4762a1bSJed Brown /*
2c4762a1bSJed Brown Include "petsctao.h" so we can use TAO solvers
3c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing
4c4762a1bSJed Brown Include "petscksp.h" so we can set KSP type
5c4762a1bSJed Brown the parallel mesh.
6c4762a1bSJed Brown */
7c4762a1bSJed Brown
8c4762a1bSJed Brown #include <petsctao.h>
9c4762a1bSJed Brown #include <petscdmda.h>
10c4762a1bSJed Brown
119371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to \n\
12c4762a1bSJed Brown solve a bound constrained minimization problem. This example is based on \n\
13c4762a1bSJed Brown the problem DPJB from the MINPACK-2 test suite. This pressure journal \n\
14c4762a1bSJed Brown bearing problem is an example of elliptic variational problem defined over \n\
15c4762a1bSJed Brown a two dimensional rectangle. By discretizing the domain into triangular \n\
16c4762a1bSJed Brown elements, the pressure surrounding the journal bearing is defined as the \n\
17c4762a1bSJed Brown minimum of a quadratic function whose variables are bounded below by zero.\n\
18c4762a1bSJed Brown The command line options are:\n\
19c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
20c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
21c4762a1bSJed Brown \n";
22c4762a1bSJed Brown
23c4762a1bSJed Brown /*
24c4762a1bSJed Brown User-defined application context - contains data needed by the
25c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient(),
26c4762a1bSJed Brown FormHessian().
27c4762a1bSJed Brown */
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown /* problem parameters */
30c4762a1bSJed Brown PetscReal ecc; /* test problem parameter */
31c4762a1bSJed Brown PetscReal b; /* A dimension of journal bearing */
32c4762a1bSJed Brown PetscInt nx, ny; /* discretization in x, y directions */
33c4762a1bSJed Brown
34c4762a1bSJed Brown /* Working space */
35c4762a1bSJed Brown DM dm; /* distributed array data structure */
36c4762a1bSJed Brown Mat A; /* Quadratic Objective term */
37c4762a1bSJed Brown Vec B; /* Linear Objective term */
38c4762a1bSJed Brown } AppCtx;
39c4762a1bSJed Brown
40c4762a1bSJed Brown /* User-defined routines */
41c4762a1bSJed Brown static PetscReal p(PetscReal xi, PetscReal ecc);
42c4762a1bSJed Brown static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
43c4762a1bSJed Brown static PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
44c4762a1bSJed Brown static PetscErrorCode ComputeB(AppCtx *);
45c4762a1bSJed Brown static PetscErrorCode Monitor(Tao, void *);
46c4762a1bSJed Brown static PetscErrorCode ConvergenceTest(Tao, void *);
47c4762a1bSJed Brown
main(int argc,char ** argv)48d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
49d71ae5a4SJacob Faibussowitsch {
50c4762a1bSJed Brown PetscInt Nx, Ny; /* number of processors in x- and y- directions */
51c4762a1bSJed Brown PetscInt m; /* number of local elements in vectors */
52c4762a1bSJed Brown Vec x; /* variables vector */
53c4762a1bSJed Brown Vec xl, xu; /* bounds vectors */
54c4762a1bSJed Brown PetscReal d1000 = 1000;
55c4762a1bSJed Brown PetscBool flg, testgetdiag; /* A return variable when checking for user options */
56c4762a1bSJed Brown Tao tao; /* Tao solver context */
57c4762a1bSJed Brown KSP ksp;
58c4762a1bSJed Brown AppCtx user; /* user-defined work context */
59c4762a1bSJed Brown PetscReal zero = 0.0; /* lower bound on all variables */
60c4762a1bSJed Brown
61f605775fSPierre Jolivet /* Initialize PETSc and TAO */
62327415f7SBarry Smith PetscFunctionBeginUser;
63c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
64c4762a1bSJed Brown
65c4762a1bSJed Brown /* Set the default values for the problem parameters */
669371c9d4SSatish Balay user.nx = 50;
679371c9d4SSatish Balay user.ny = 50;
689371c9d4SSatish Balay user.ecc = 0.1;
699371c9d4SSatish Balay user.b = 10.0;
70c4762a1bSJed Brown testgetdiag = PETSC_FALSE;
71c4762a1bSJed Brown
72c4762a1bSJed Brown /* Check for any command line arguments that override defaults */
739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.nx, &flg));
749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.ny, &flg));
759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-ecc", &user.ecc, &flg));
769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-b", &user.b, &flg));
779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_getdiagonal", &testgetdiag, NULL));
78c4762a1bSJed Brown
799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n---- Journal Bearing Problem SHB-----\n"));
8063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mx: %" PetscInt_FMT ", my: %" PetscInt_FMT ", ecc: %g \n\n", user.nx, user.ny, (double)user.ecc));
81c4762a1bSJed Brown
82f0b74427SPierre Jolivet /* Let PETSc determine the grid division */
839371c9d4SSatish Balay Nx = PETSC_DECIDE;
849371c9d4SSatish Balay Ny = PETSC_DECIDE;
85c4762a1bSJed Brown
86c4762a1bSJed Brown /*
87c4762a1bSJed Brown A two dimensional distributed array will help define this problem,
88c4762a1bSJed Brown which derives from an elliptic PDE on two dimensional domain. From
89c4762a1bSJed Brown the distributed array, Create the vectors.
90c4762a1bSJed Brown */
919566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.nx, user.ny, Nx, Ny, 1, 1, NULL, NULL, &user.dm));
929566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm));
939566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm));
94c4762a1bSJed Brown
95c4762a1bSJed Brown /*
96c4762a1bSJed Brown Extract global and local vectors from DM; the vector user.B is
97c4762a1bSJed Brown used solely as work space for the evaluation of the function,
98c4762a1bSJed Brown gradient, and Hessian. Duplicate for remaining vectors that are
99c4762a1bSJed Brown the same types.
100c4762a1bSJed Brown */
1019566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm, &x)); /* Solution */
1029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &user.B)); /* Linear objective */
103c4762a1bSJed Brown
104c4762a1bSJed Brown /* Create matrix user.A to store quadratic, Create a local ordering scheme. */
1059566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &m));
1069566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm, &user.A));
107c4762a1bSJed Brown
1081baa6e33SBarry Smith if (testgetdiag) PetscCall(MatSetOperation(user.A, MATOP_GET_DIAGONAL, NULL));
109c4762a1bSJed Brown
110c4762a1bSJed Brown /* User defined function -- compute linear term of quadratic */
1119566063dSJacob Faibussowitsch PetscCall(ComputeB(&user));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /* The TAO code begins here */
114c4762a1bSJed Brown
115c4762a1bSJed Brown /*
116c4762a1bSJed Brown Create the optimization solver
117c4762a1bSJed Brown Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
118c4762a1bSJed Brown */
1199566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1209566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM));
121c4762a1bSJed Brown
122c4762a1bSJed Brown /* Set the initial vector */
1239566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero));
1249566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
125c4762a1bSJed Brown
126c4762a1bSJed Brown /* Set the user function, gradient, hessian evaluation routines and data structures */
1279566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
128c4762a1bSJed Brown
1299566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.A, user.A, FormHessian, (void *)&user));
130c4762a1bSJed Brown
131c4762a1bSJed Brown /* Set a routine that defines the bounds */
1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xl));
1339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xu));
1349566063dSJacob Faibussowitsch PetscCall(VecSet(xl, zero));
1359566063dSJacob Faibussowitsch PetscCall(VecSet(xu, d1000));
1369566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, xl, xu));
137c4762a1bSJed Brown
1389566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp));
1391baa6e33SBarry Smith if (ksp) PetscCall(KSPSetType(ksp, KSPCG));
140c4762a1bSJed Brown
1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-testmonitor", &flg));
14210978b7dSBarry Smith if (flg) PetscCall(TaoMonitorSet(tao, Monitor, &user, NULL));
1439566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-testconvergence", &flg));
14448a46eb9SPierre Jolivet if (flg) PetscCall(TaoSetConvergenceTest(tao, ConvergenceTest, &user));
145c4762a1bSJed Brown
146c4762a1bSJed Brown /* Check for any tao command line options */
1479566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
148c4762a1bSJed Brown
149c4762a1bSJed Brown /* Solve the bound constrained problem */
1509566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
151c4762a1bSJed Brown
152c4762a1bSJed Brown /* Free PETSc data structures */
1539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xl));
1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xu));
1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A));
1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.B));
158c4762a1bSJed Brown
159c4762a1bSJed Brown /* Free TAO data structures */
1609566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
1619566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm));
1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
163b122ec5aSJacob Faibussowitsch return 0;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown
p(PetscReal xi,PetscReal ecc)166d71ae5a4SJacob Faibussowitsch static PetscReal p(PetscReal xi, PetscReal ecc)
167d71ae5a4SJacob Faibussowitsch {
168c4762a1bSJed Brown PetscReal t = 1.0 + ecc * PetscCosScalar(xi);
1694ad8454bSPierre Jolivet return t * t * t;
170c4762a1bSJed Brown }
171c4762a1bSJed Brown
ComputeB(AppCtx * user)172d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeB(AppCtx *user)
173d71ae5a4SJacob Faibussowitsch {
174c4762a1bSJed Brown PetscInt i, j, k;
175c4762a1bSJed Brown PetscInt nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym;
176c4762a1bSJed Brown PetscReal two = 2.0, pi = 4.0 * atan(1.0);
177c4762a1bSJed Brown PetscReal hx, hy, ehxhy;
178c4762a1bSJed Brown PetscReal temp, *b;
179c4762a1bSJed Brown PetscReal ecc = user->ecc;
180c4762a1bSJed Brown
181780b99b1SStefano Zampini PetscFunctionBegin;
182c4762a1bSJed Brown nx = user->nx;
183c4762a1bSJed Brown ny = user->ny;
184c4762a1bSJed Brown hx = two * pi / (nx + 1.0);
185c4762a1bSJed Brown hy = two * user->b / (ny + 1.0);
186c4762a1bSJed Brown ehxhy = ecc * hx * hy;
187c4762a1bSJed Brown
188c4762a1bSJed Brown /*
189c4762a1bSJed Brown Get local grid boundaries
190c4762a1bSJed Brown */
1919566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
1929566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
193c4762a1bSJed Brown
194c4762a1bSJed Brown /* Compute the linear term in the objective function */
1959566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->B, &b));
196c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
197c4762a1bSJed Brown temp = PetscSinScalar((i + 1) * hx);
198c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
199c4762a1bSJed Brown k = xm * (j - ys) + (i - xs);
200c4762a1bSJed Brown b[k] = -ehxhy * temp;
201c4762a1bSJed Brown }
202c4762a1bSJed Brown }
2039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->B, &b));
2049566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(5.0 * xm * ym + 3.0 * xm));
2053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown
FormFunctionGradient(Tao tao,Vec X,PetscReal * fcn,Vec G,void * ptr)208d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *ptr)
209d71ae5a4SJacob Faibussowitsch {
210c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
211c4762a1bSJed Brown PetscInt i, j, k, kk;
212c4762a1bSJed Brown PetscInt col[5], row, nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym;
213c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, six = 6.0, pi = 4.0 * atan(1.0);
214c4762a1bSJed Brown PetscReal hx, hy, hxhy, hxhx, hyhy;
215c4762a1bSJed Brown PetscReal xi, v[5];
216c4762a1bSJed Brown PetscReal ecc = user->ecc, trule1, trule2, trule3, trule4, trule5, trule6;
217c4762a1bSJed Brown PetscReal vmiddle, vup, vdown, vleft, vright;
218c4762a1bSJed Brown PetscReal tt, f1, f2;
219c4762a1bSJed Brown PetscReal *x, *g, zero = 0.0;
220c4762a1bSJed Brown Vec localX;
221c4762a1bSJed Brown
222780b99b1SStefano Zampini PetscFunctionBegin;
223c4762a1bSJed Brown nx = user->nx;
224c4762a1bSJed Brown ny = user->ny;
225c4762a1bSJed Brown hx = two * pi / (nx + 1.0);
226c4762a1bSJed Brown hy = two * user->b / (ny + 1.0);
227c4762a1bSJed Brown hxhy = hx * hy;
228c4762a1bSJed Brown hxhx = one / (hx * hx);
229c4762a1bSJed Brown hyhy = one / (hy * hy);
230c4762a1bSJed Brown
2319566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localX));
232c4762a1bSJed Brown
2339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
2349566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
235c4762a1bSJed Brown
2369566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero));
237c4762a1bSJed Brown /*
238c4762a1bSJed Brown Get local grid boundaries
239c4762a1bSJed Brown */
2409566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
2419566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
242c4762a1bSJed Brown
2439566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x));
2449566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g));
245c4762a1bSJed Brown
246c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
247c4762a1bSJed Brown xi = (i + 1) * hx;
248c4762a1bSJed Brown trule1 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi, ecc)) / six; /* L(i,j) */
249c4762a1bSJed Brown trule2 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi, ecc)) / six; /* U(i,j) */
250c4762a1bSJed Brown trule3 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi + hx, ecc)) / six; /* U(i+1,j) */
251c4762a1bSJed Brown trule4 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi - hx, ecc)) / six; /* L(i-1,j) */
252c4762a1bSJed Brown trule5 = trule1; /* L(i,j-1) */
253c4762a1bSJed Brown trule6 = trule2; /* U(i,j+1) */
254c4762a1bSJed Brown
255c4762a1bSJed Brown vdown = -(trule5 + trule2) * hyhy;
256c4762a1bSJed Brown vleft = -hxhx * (trule2 + trule4);
257c4762a1bSJed Brown vright = -hxhx * (trule1 + trule3);
258c4762a1bSJed Brown vup = -hyhy * (trule1 + trule6);
259c4762a1bSJed Brown vmiddle = (hxhx) * (trule1 + trule2 + trule3 + trule4) + hyhy * (trule1 + trule2 + trule5 + trule6);
260c4762a1bSJed Brown
261c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
262c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs);
2639371c9d4SSatish Balay v[0] = 0;
2649371c9d4SSatish Balay v[1] = 0;
2659371c9d4SSatish Balay v[2] = 0;
2669371c9d4SSatish Balay v[3] = 0;
2679371c9d4SSatish Balay v[4] = 0;
268c4762a1bSJed Brown
269c4762a1bSJed Brown k = 0;
270c4762a1bSJed Brown if (j > gys) {
2719371c9d4SSatish Balay v[k] = vdown;
2729371c9d4SSatish Balay col[k] = row - gxm;
2739371c9d4SSatish Balay k++;
274c4762a1bSJed Brown }
275c4762a1bSJed Brown
276c4762a1bSJed Brown if (i > gxs) {
2779371c9d4SSatish Balay v[k] = vleft;
2789371c9d4SSatish Balay col[k] = row - 1;
2799371c9d4SSatish Balay k++;
280c4762a1bSJed Brown }
281c4762a1bSJed Brown
2829371c9d4SSatish Balay v[k] = vmiddle;
2839371c9d4SSatish Balay col[k] = row;
2849371c9d4SSatish Balay k++;
285c4762a1bSJed Brown
286c4762a1bSJed Brown if (i + 1 < gxs + gxm) {
2879371c9d4SSatish Balay v[k] = vright;
2889371c9d4SSatish Balay col[k] = row + 1;
2899371c9d4SSatish Balay k++;
290c4762a1bSJed Brown }
291c4762a1bSJed Brown
292c4762a1bSJed Brown if (j + 1 < gys + gym) {
2939371c9d4SSatish Balay v[k] = vup;
2949371c9d4SSatish Balay col[k] = row + gxm;
2959371c9d4SSatish Balay k++;
296c4762a1bSJed Brown }
297c4762a1bSJed Brown tt = 0;
298ad540459SPierre Jolivet for (kk = 0; kk < k; kk++) tt += v[kk] * x[col[kk]];
299c4762a1bSJed Brown row = (j - ys) * xm + (i - xs);
300c4762a1bSJed Brown g[row] = tt;
301c4762a1bSJed Brown }
302c4762a1bSJed Brown }
303c4762a1bSJed Brown
3049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x));
3059566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g));
306c4762a1bSJed Brown
3079566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localX));
308c4762a1bSJed Brown
3099566063dSJacob Faibussowitsch PetscCall(VecDot(X, G, &f1));
3109566063dSJacob Faibussowitsch PetscCall(VecDot(user->B, X, &f2));
3119566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, one, user->B));
312c4762a1bSJed Brown *fcn = f1 / 2.0 + f2;
313c4762a1bSJed Brown
3149566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((91 + 10.0 * ym) * xm));
3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
316c4762a1bSJed Brown }
317c4762a1bSJed Brown
318c4762a1bSJed Brown /*
319c4762a1bSJed Brown FormHessian computes the quadratic term in the quadratic objective function
320c4762a1bSJed Brown Notice that the objective function in this problem is quadratic (therefore a constant
321c4762a1bSJed Brown hessian). If using a nonquadratic solver, then you might want to reconsider this function
322c4762a1bSJed Brown */
FormHessian(Tao tao,Vec X,Mat hes,Mat Hpre,void * ptr)323d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat hes, Mat Hpre, void *ptr)
324d71ae5a4SJacob Faibussowitsch {
325c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
326c4762a1bSJed Brown PetscInt i, j, k;
327c4762a1bSJed Brown PetscInt col[5], row, nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym;
328c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, six = 6.0, pi = 4.0 * atan(1.0);
329c4762a1bSJed Brown PetscReal hx, hy, hxhy, hxhx, hyhy;
330c4762a1bSJed Brown PetscReal xi, v[5];
331c4762a1bSJed Brown PetscReal ecc = user->ecc, trule1, trule2, trule3, trule4, trule5, trule6;
332c4762a1bSJed Brown PetscReal vmiddle, vup, vdown, vleft, vright;
333c4762a1bSJed Brown PetscBool assembled;
334c4762a1bSJed Brown
335780b99b1SStefano Zampini PetscFunctionBegin;
336c4762a1bSJed Brown nx = user->nx;
337c4762a1bSJed Brown ny = user->ny;
338c4762a1bSJed Brown hx = two * pi / (nx + 1.0);
339c4762a1bSJed Brown hy = two * user->b / (ny + 1.0);
340c4762a1bSJed Brown hxhy = hx * hy;
341c4762a1bSJed Brown hxhx = one / (hx * hx);
342c4762a1bSJed Brown hyhy = one / (hy * hy);
343c4762a1bSJed Brown
344c4762a1bSJed Brown /*
345c4762a1bSJed Brown Get local grid boundaries
346c4762a1bSJed Brown */
3479566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
3489566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
3499566063dSJacob Faibussowitsch PetscCall(MatAssembled(hes, &assembled));
3509566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(hes));
351c4762a1bSJed Brown
352c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
353c4762a1bSJed Brown xi = (i + 1) * hx;
354c4762a1bSJed Brown trule1 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi, ecc)) / six; /* L(i,j) */
355c4762a1bSJed Brown trule2 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi, ecc)) / six; /* U(i,j) */
356c4762a1bSJed Brown trule3 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi + hx, ecc)) / six; /* U(i+1,j) */
357c4762a1bSJed Brown trule4 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi - hx, ecc)) / six; /* L(i-1,j) */
358c4762a1bSJed Brown trule5 = trule1; /* L(i,j-1) */
359c4762a1bSJed Brown trule6 = trule2; /* U(i,j+1) */
360c4762a1bSJed Brown
361c4762a1bSJed Brown vdown = -(trule5 + trule2) * hyhy;
362c4762a1bSJed Brown vleft = -hxhx * (trule2 + trule4);
363c4762a1bSJed Brown vright = -hxhx * (trule1 + trule3);
364c4762a1bSJed Brown vup = -hyhy * (trule1 + trule6);
365c4762a1bSJed Brown vmiddle = (hxhx) * (trule1 + trule2 + trule3 + trule4) + hyhy * (trule1 + trule2 + trule5 + trule6);
3669371c9d4SSatish Balay v[0] = 0;
3679371c9d4SSatish Balay v[1] = 0;
3689371c9d4SSatish Balay v[2] = 0;
3699371c9d4SSatish Balay v[3] = 0;
3709371c9d4SSatish Balay v[4] = 0;
371c4762a1bSJed Brown
372c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
373c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs);
374c4762a1bSJed Brown
375c4762a1bSJed Brown k = 0;
376c4762a1bSJed Brown if (j > gys) {
3779371c9d4SSatish Balay v[k] = vdown;
3789371c9d4SSatish Balay col[k] = row - gxm;
3799371c9d4SSatish Balay k++;
380c4762a1bSJed Brown }
381c4762a1bSJed Brown
382c4762a1bSJed Brown if (i > gxs) {
3839371c9d4SSatish Balay v[k] = vleft;
3849371c9d4SSatish Balay col[k] = row - 1;
3859371c9d4SSatish Balay k++;
386c4762a1bSJed Brown }
387c4762a1bSJed Brown
3889371c9d4SSatish Balay v[k] = vmiddle;
3899371c9d4SSatish Balay col[k] = row;
3909371c9d4SSatish Balay k++;
391c4762a1bSJed Brown
392c4762a1bSJed Brown if (i + 1 < gxs + gxm) {
3939371c9d4SSatish Balay v[k] = vright;
3949371c9d4SSatish Balay col[k] = row + 1;
3959371c9d4SSatish Balay k++;
396c4762a1bSJed Brown }
397c4762a1bSJed Brown
398c4762a1bSJed Brown if (j + 1 < gys + gym) {
3999371c9d4SSatish Balay v[k] = vup;
4009371c9d4SSatish Balay col[k] = row + gxm;
4019371c9d4SSatish Balay k++;
402c4762a1bSJed Brown }
4039566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(hes, 1, &row, k, col, v, INSERT_VALUES));
404c4762a1bSJed Brown }
405c4762a1bSJed Brown }
406c4762a1bSJed Brown
407c4762a1bSJed Brown /*
408c4762a1bSJed Brown Assemble matrix, using the 2-step process:
409c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd().
410c4762a1bSJed Brown By placing code between these two statements, computations can be
411c4762a1bSJed Brown done while messages are in transition.
412c4762a1bSJed Brown */
4139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(hes, MAT_FINAL_ASSEMBLY));
4149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(hes, MAT_FINAL_ASSEMBLY));
415c4762a1bSJed Brown
416c4762a1bSJed Brown /*
417c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the
418c4762a1bSJed Brown matrix. If we do it will generate an error.
419c4762a1bSJed Brown */
4209566063dSJacob Faibussowitsch PetscCall(MatSetOption(hes, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
4219566063dSJacob Faibussowitsch PetscCall(MatSetOption(hes, MAT_SYMMETRIC, PETSC_TRUE));
422c4762a1bSJed Brown
4239566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0 * xm * ym + 49.0 * xm));
4243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
425c4762a1bSJed Brown }
426c4762a1bSJed Brown
Monitor(Tao tao,PetscCtx ctx)4272a8381b2SBarry Smith PetscErrorCode Monitor(Tao tao, PetscCtx ctx)
428d71ae5a4SJacob Faibussowitsch {
429c4762a1bSJed Brown PetscInt its;
430c4762a1bSJed Brown PetscReal f, gnorm, cnorm, xdiff;
431c4762a1bSJed Brown TaoConvergedReason reason;
432c4762a1bSJed Brown
433c4762a1bSJed Brown PetscFunctionBegin;
4349566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
43548a46eb9SPierre Jolivet if (!(its % 5)) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration=%" PetscInt_FMT "\tf=%g\n", its, (double)f));
4363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
437c4762a1bSJed Brown }
438c4762a1bSJed Brown
ConvergenceTest(Tao tao,PetscCtx ctx)4392a8381b2SBarry Smith PetscErrorCode ConvergenceTest(Tao tao, PetscCtx ctx)
440d71ae5a4SJacob Faibussowitsch {
441c4762a1bSJed Brown PetscInt its;
442c4762a1bSJed Brown PetscReal f, gnorm, cnorm, xdiff;
443c4762a1bSJed Brown TaoConvergedReason reason;
444c4762a1bSJed Brown
445c4762a1bSJed Brown PetscFunctionBegin;
4469566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
44748a46eb9SPierre Jolivet if (its == 100) PetscCall(TaoSetConvergedReason(tao, TAO_DIVERGED_MAXITS));
4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
449c4762a1bSJed Brown }
450c4762a1bSJed Brown
451c4762a1bSJed Brown /*TEST
452c4762a1bSJed Brown
453c4762a1bSJed Brown build:
454c4762a1bSJed Brown requires: !complex
455c4762a1bSJed Brown
456c4762a1bSJed Brown test:
45710978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5
458c4762a1bSJed Brown requires: !single
459c4762a1bSJed Brown
460c4762a1bSJed Brown test:
461c4762a1bSJed Brown suffix: 2
462c4762a1bSJed Brown nsize: 2
46310978b7dSBarry Smith args: -tao_monitor_short -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5
464c4762a1bSJed Brown requires: !single
465c4762a1bSJed Brown
466c4762a1bSJed Brown test:
467c4762a1bSJed Brown suffix: 3
468c4762a1bSJed Brown nsize: 2
46910978b7dSBarry Smith args: -tao_monitor_short -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4
470c4762a1bSJed Brown requires: !single
471c4762a1bSJed Brown
472c4762a1bSJed Brown test:
473c4762a1bSJed Brown suffix: 4
474c4762a1bSJed Brown nsize: 2
47510978b7dSBarry Smith args: -tao_monitor_short -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal
476c4762a1bSJed Brown output_file: output/jbearing2_3.out
477c4762a1bSJed Brown requires: !single
478c4762a1bSJed Brown
479c4762a1bSJed Brown test:
480c4762a1bSJed Brown suffix: 5
48110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
482c4762a1bSJed Brown requires: !single
483c4762a1bSJed Brown
484c4762a1bSJed Brown test:
485c4762a1bSJed Brown suffix: 6
48610978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4
487c4762a1bSJed Brown requires: !single
488c4762a1bSJed Brown
489c4762a1bSJed Brown test:
490c4762a1bSJed Brown suffix: 7
49110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5
492c4762a1bSJed Brown requires: !single
493c4762a1bSJed Brown
494c4762a1bSJed Brown test:
495c4762a1bSJed Brown suffix: 8
49610978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5
497c4762a1bSJed Brown requires: !single
498c4762a1bSJed Brown
499c4762a1bSJed Brown test:
500c4762a1bSJed Brown suffix: 9
50110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5
502c4762a1bSJed Brown requires: !single
503c4762a1bSJed Brown
504c4762a1bSJed Brown test:
505c4762a1bSJed Brown suffix: 10
506*a336c150SZach Atkins args: -tao_monitor_short -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 -tao_bnk_cg_tao_monitor_short
507c4762a1bSJed Brown requires: !single
508c4762a1bSJed Brown
509c4762a1bSJed Brown test:
510c4762a1bSJed Brown suffix: 11
511*a336c150SZach Atkins args: -tao_monitor_short -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 -tao_bnk_cg_tao_monitor_short
512c4762a1bSJed Brown requires: !single
513c4762a1bSJed Brown
514c4762a1bSJed Brown test:
515c4762a1bSJed Brown suffix: 12
516*a336c150SZach Atkins args: -tao_monitor_short -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 -tao_bnk_cg_tao_monitor_short
517c4762a1bSJed Brown requires: !single
518c4762a1bSJed Brown
519c4762a1bSJed Brown test:
520c4762a1bSJed Brown suffix: 13
52110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls
522c4762a1bSJed Brown requires: !single
523c4762a1bSJed Brown
524c4762a1bSJed Brown test:
525c4762a1bSJed Brown suffix: 14
52610978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm
527c4762a1bSJed Brown requires: !single
528c4762a1bSJed Brown
529c4762a1bSJed Brown test:
530c4762a1bSJed Brown suffix: 15
53110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
532c4762a1bSJed Brown requires: !single
533c4762a1bSJed Brown
534c4762a1bSJed Brown test:
535c4762a1bSJed Brown suffix: 16
53610978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
537c4762a1bSJed Brown requires: !single
538c4762a1bSJed Brown
539c4762a1bSJed Brown test:
540c4762a1bSJed Brown suffix: 17
54110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view
542c4762a1bSJed Brown requires: !single
543c4762a1bSJed Brown
544c4762a1bSJed Brown test:
545c4762a1bSJed Brown suffix: 18
54610978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view
547c4762a1bSJed Brown requires: !single
548c4762a1bSJed Brown
54934ad9904SAlp Dener test:
55034ad9904SAlp Dener suffix: 19
55110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
55234ad9904SAlp Dener requires: !single
55334ad9904SAlp Dener
55434ad9904SAlp Dener test:
55534ad9904SAlp Dener suffix: 20
55610978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
55734ad9904SAlp Dener requires: !single
55834ad9904SAlp Dener
55934ad9904SAlp Dener test:
56034ad9904SAlp Dener suffix: 21
56110978b7dSBarry Smith args: -tao_monitor_short -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
56234ad9904SAlp Dener requires: !single
563c4762a1bSJed Brown TEST*/
564