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 489371c9d4SSatish Balay int main(int argc, char **argv) { 49c4762a1bSJed Brown PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 50c4762a1bSJed Brown PetscInt m; /* number of local elements in vectors */ 51c4762a1bSJed Brown Vec x; /* variables vector */ 52c4762a1bSJed Brown Vec xl, xu; /* bounds vectors */ 53c4762a1bSJed Brown PetscReal d1000 = 1000; 54c4762a1bSJed Brown PetscBool flg, testgetdiag; /* A return variable when checking for user options */ 55c4762a1bSJed Brown Tao tao; /* Tao solver context */ 56c4762a1bSJed Brown KSP ksp; 57c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 58c4762a1bSJed Brown PetscReal zero = 0.0; /* lower bound on all variables */ 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Initialize PETSC and TAO */ 61327415f7SBarry Smith PetscFunctionBeginUser; 629566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Set the default values for the problem parameters */ 659371c9d4SSatish Balay user.nx = 50; 669371c9d4SSatish Balay user.ny = 50; 679371c9d4SSatish Balay user.ecc = 0.1; 689371c9d4SSatish Balay user.b = 10.0; 69c4762a1bSJed Brown testgetdiag = PETSC_FALSE; 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.nx, &flg)); 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.ny, &flg)); 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-ecc", &user.ecc, &flg)); 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-b", &user.b, &flg)); 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_getdiagonal", &testgetdiag, NULL)); 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n---- Journal Bearing Problem SHB-----\n")); 7963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mx: %" PetscInt_FMT ", my: %" PetscInt_FMT ", ecc: %g \n\n", user.nx, user.ny, (double)user.ecc)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Let Petsc determine the grid division */ 829371c9d4SSatish Balay Nx = PETSC_DECIDE; 839371c9d4SSatish Balay Ny = PETSC_DECIDE; 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* 86c4762a1bSJed Brown A two dimensional distributed array will help define this problem, 87c4762a1bSJed Brown which derives from an elliptic PDE on two dimensional domain. From 88c4762a1bSJed Brown the distributed array, Create the vectors. 89c4762a1bSJed Brown */ 909566063dSJacob 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)); 919566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm)); 929566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* 95c4762a1bSJed Brown Extract global and local vectors from DM; the vector user.B is 96c4762a1bSJed Brown used solely as work space for the evaluation of the function, 97c4762a1bSJed Brown gradient, and Hessian. Duplicate for remaining vectors that are 98c4762a1bSJed Brown the same types. 99c4762a1bSJed Brown */ 1009566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm, &x)); /* Solution */ 1019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &user.B)); /* Linear objective */ 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Create matrix user.A to store quadratic, Create a local ordering scheme. */ 1049566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &m)); 1059566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm, &user.A)); 106c4762a1bSJed Brown 1071baa6e33SBarry Smith if (testgetdiag) PetscCall(MatSetOperation(user.A, MATOP_GET_DIAGONAL, NULL)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* User defined function -- compute linear term of quadratic */ 1109566063dSJacob Faibussowitsch PetscCall(ComputeB(&user)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* The TAO code begins here */ 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown Create the optimization solver 116c4762a1bSJed Brown Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM 117c4762a1bSJed Brown */ 1189566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 1199566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* Set the initial vector */ 1229566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero)); 1239566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Set the user function, gradient, hessian evaluation routines and data structures */ 1269566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.A, user.A, FormHessian, (void *)&user)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Set a routine that defines the bounds */ 1319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xl)); 1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xu)); 1339566063dSJacob Faibussowitsch PetscCall(VecSet(xl, zero)); 1349566063dSJacob Faibussowitsch PetscCall(VecSet(xu, d1000)); 1359566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, xl, xu)); 136c4762a1bSJed Brown 1379566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 1381baa6e33SBarry Smith if (ksp) PetscCall(KSPSetType(ksp, KSPCG)); 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-testmonitor", &flg)); 141*48a46eb9SPierre Jolivet if (flg) PetscCall(TaoSetMonitor(tao, Monitor, &user, NULL)); 1429566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-testconvergence", &flg)); 143*48a46eb9SPierre Jolivet if (flg) PetscCall(TaoSetConvergenceTest(tao, ConvergenceTest, &user)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Check for any tao command line options */ 1469566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* Solve the bound constrained problem */ 1499566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* Free PETSc data structures */ 1529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xl)); 1549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xu)); 1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A)); 1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.B)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Free TAO data structures */ 1599566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1609566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 1619566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 162b122ec5aSJacob Faibussowitsch return 0; 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 1659371c9d4SSatish Balay static PetscReal p(PetscReal xi, PetscReal ecc) { 166c4762a1bSJed Brown PetscReal t = 1.0 + ecc * PetscCosScalar(xi); 167c4762a1bSJed Brown return (t * t * t); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 1709371c9d4SSatish Balay PetscErrorCode ComputeB(AppCtx *user) { 171c4762a1bSJed Brown PetscInt i, j, k; 172c4762a1bSJed Brown PetscInt nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym; 173c4762a1bSJed Brown PetscReal two = 2.0, pi = 4.0 * atan(1.0); 174c4762a1bSJed Brown PetscReal hx, hy, ehxhy; 175c4762a1bSJed Brown PetscReal temp, *b; 176c4762a1bSJed Brown PetscReal ecc = user->ecc; 177c4762a1bSJed Brown 178780b99b1SStefano Zampini PetscFunctionBegin; 179c4762a1bSJed Brown nx = user->nx; 180c4762a1bSJed Brown ny = user->ny; 181c4762a1bSJed Brown hx = two * pi / (nx + 1.0); 182c4762a1bSJed Brown hy = two * user->b / (ny + 1.0); 183c4762a1bSJed Brown ehxhy = ecc * hx * hy; 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown Get local grid boundaries 187c4762a1bSJed Brown */ 1889566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 1899566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* Compute the linear term in the objective function */ 1929566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->B, &b)); 193c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 194c4762a1bSJed Brown temp = PetscSinScalar((i + 1) * hx); 195c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 196c4762a1bSJed Brown k = xm * (j - ys) + (i - xs); 197c4762a1bSJed Brown b[k] = -ehxhy * temp; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown } 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->B, &b)); 2019566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(5.0 * xm * ym + 3.0 * xm)); 202780b99b1SStefano Zampini PetscFunctionReturn(0); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 2059371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *ptr) { 206c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 207c4762a1bSJed Brown PetscInt i, j, k, kk; 208c4762a1bSJed Brown PetscInt col[5], row, nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym; 209c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, six = 6.0, pi = 4.0 * atan(1.0); 210c4762a1bSJed Brown PetscReal hx, hy, hxhy, hxhx, hyhy; 211c4762a1bSJed Brown PetscReal xi, v[5]; 212c4762a1bSJed Brown PetscReal ecc = user->ecc, trule1, trule2, trule3, trule4, trule5, trule6; 213c4762a1bSJed Brown PetscReal vmiddle, vup, vdown, vleft, vright; 214c4762a1bSJed Brown PetscReal tt, f1, f2; 215c4762a1bSJed Brown PetscReal *x, *g, zero = 0.0; 216c4762a1bSJed Brown Vec localX; 217c4762a1bSJed Brown 218780b99b1SStefano Zampini PetscFunctionBegin; 219c4762a1bSJed Brown nx = user->nx; 220c4762a1bSJed Brown ny = user->ny; 221c4762a1bSJed Brown hx = two * pi / (nx + 1.0); 222c4762a1bSJed Brown hy = two * user->b / (ny + 1.0); 223c4762a1bSJed Brown hxhy = hx * hy; 224c4762a1bSJed Brown hxhx = one / (hx * hx); 225c4762a1bSJed Brown hyhy = one / (hy * hy); 226c4762a1bSJed Brown 2279566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localX)); 228c4762a1bSJed Brown 2299566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 2309566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 231c4762a1bSJed Brown 2329566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 233c4762a1bSJed Brown /* 234c4762a1bSJed Brown Get local grid boundaries 235c4762a1bSJed Brown */ 2369566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 2379566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 238c4762a1bSJed Brown 2399566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x)); 2409566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 243c4762a1bSJed Brown xi = (i + 1) * hx; 244c4762a1bSJed Brown trule1 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi, ecc)) / six; /* L(i,j) */ 245c4762a1bSJed Brown trule2 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi, ecc)) / six; /* U(i,j) */ 246c4762a1bSJed Brown trule3 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi + hx, ecc)) / six; /* U(i+1,j) */ 247c4762a1bSJed Brown trule4 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi - hx, ecc)) / six; /* L(i-1,j) */ 248c4762a1bSJed Brown trule5 = trule1; /* L(i,j-1) */ 249c4762a1bSJed Brown trule6 = trule2; /* U(i,j+1) */ 250c4762a1bSJed Brown 251c4762a1bSJed Brown vdown = -(trule5 + trule2) * hyhy; 252c4762a1bSJed Brown vleft = -hxhx * (trule2 + trule4); 253c4762a1bSJed Brown vright = -hxhx * (trule1 + trule3); 254c4762a1bSJed Brown vup = -hyhy * (trule1 + trule6); 255c4762a1bSJed Brown vmiddle = (hxhx) * (trule1 + trule2 + trule3 + trule4) + hyhy * (trule1 + trule2 + trule5 + trule6); 256c4762a1bSJed Brown 257c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 258c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs); 2599371c9d4SSatish Balay v[0] = 0; 2609371c9d4SSatish Balay v[1] = 0; 2619371c9d4SSatish Balay v[2] = 0; 2629371c9d4SSatish Balay v[3] = 0; 2639371c9d4SSatish Balay v[4] = 0; 264c4762a1bSJed Brown 265c4762a1bSJed Brown k = 0; 266c4762a1bSJed Brown if (j > gys) { 2679371c9d4SSatish Balay v[k] = vdown; 2689371c9d4SSatish Balay col[k] = row - gxm; 2699371c9d4SSatish Balay k++; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown if (i > gxs) { 2739371c9d4SSatish Balay v[k] = vleft; 2749371c9d4SSatish Balay col[k] = row - 1; 2759371c9d4SSatish Balay k++; 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 2789371c9d4SSatish Balay v[k] = vmiddle; 2799371c9d4SSatish Balay col[k] = row; 2809371c9d4SSatish Balay k++; 281c4762a1bSJed Brown 282c4762a1bSJed Brown if (i + 1 < gxs + gxm) { 2839371c9d4SSatish Balay v[k] = vright; 2849371c9d4SSatish Balay col[k] = row + 1; 2859371c9d4SSatish Balay k++; 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown if (j + 1 < gys + gym) { 2899371c9d4SSatish Balay v[k] = vup; 2909371c9d4SSatish Balay col[k] = row + gxm; 2919371c9d4SSatish Balay k++; 292c4762a1bSJed Brown } 293c4762a1bSJed Brown tt = 0; 2949371c9d4SSatish Balay for (kk = 0; kk < k; kk++) { tt += v[kk] * x[col[kk]]; } 295c4762a1bSJed Brown row = (j - ys) * xm + (i - xs); 296c4762a1bSJed Brown g[row] = tt; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 3009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x)); 3019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g)); 302c4762a1bSJed Brown 3039566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localX)); 304c4762a1bSJed Brown 3059566063dSJacob Faibussowitsch PetscCall(VecDot(X, G, &f1)); 3069566063dSJacob Faibussowitsch PetscCall(VecDot(user->B, X, &f2)); 3079566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, one, user->B)); 308c4762a1bSJed Brown *fcn = f1 / 2.0 + f2; 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((91 + 10.0 * ym) * xm)); 311780b99b1SStefano Zampini PetscFunctionReturn(0); 312c4762a1bSJed Brown } 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* 315c4762a1bSJed Brown FormHessian computes the quadratic term in the quadratic objective function 316c4762a1bSJed Brown Notice that the objective function in this problem is quadratic (therefore a constant 317c4762a1bSJed Brown hessian). If using a nonquadratic solver, then you might want to reconsider this function 318c4762a1bSJed Brown */ 3199371c9d4SSatish Balay PetscErrorCode FormHessian(Tao tao, Vec X, Mat hes, Mat Hpre, void *ptr) { 320c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 321c4762a1bSJed Brown PetscInt i, j, k; 322c4762a1bSJed Brown PetscInt col[5], row, nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym; 323c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, six = 6.0, pi = 4.0 * atan(1.0); 324c4762a1bSJed Brown PetscReal hx, hy, hxhy, hxhx, hyhy; 325c4762a1bSJed Brown PetscReal xi, v[5]; 326c4762a1bSJed Brown PetscReal ecc = user->ecc, trule1, trule2, trule3, trule4, trule5, trule6; 327c4762a1bSJed Brown PetscReal vmiddle, vup, vdown, vleft, vright; 328c4762a1bSJed Brown PetscBool assembled; 329c4762a1bSJed Brown 330780b99b1SStefano Zampini PetscFunctionBegin; 331c4762a1bSJed Brown nx = user->nx; 332c4762a1bSJed Brown ny = user->ny; 333c4762a1bSJed Brown hx = two * pi / (nx + 1.0); 334c4762a1bSJed Brown hy = two * user->b / (ny + 1.0); 335c4762a1bSJed Brown hxhy = hx * hy; 336c4762a1bSJed Brown hxhx = one / (hx * hx); 337c4762a1bSJed Brown hyhy = one / (hy * hy); 338c4762a1bSJed Brown 339c4762a1bSJed Brown /* 340c4762a1bSJed Brown Get local grid boundaries 341c4762a1bSJed Brown */ 3429566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 3439566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 3449566063dSJacob Faibussowitsch PetscCall(MatAssembled(hes, &assembled)); 3459566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(hes)); 346c4762a1bSJed Brown 347c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 348c4762a1bSJed Brown xi = (i + 1) * hx; 349c4762a1bSJed Brown trule1 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi, ecc)) / six; /* L(i,j) */ 350c4762a1bSJed Brown trule2 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi, ecc)) / six; /* U(i,j) */ 351c4762a1bSJed Brown trule3 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi + hx, ecc)) / six; /* U(i+1,j) */ 352c4762a1bSJed Brown trule4 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi - hx, ecc)) / six; /* L(i-1,j) */ 353c4762a1bSJed Brown trule5 = trule1; /* L(i,j-1) */ 354c4762a1bSJed Brown trule6 = trule2; /* U(i,j+1) */ 355c4762a1bSJed Brown 356c4762a1bSJed Brown vdown = -(trule5 + trule2) * hyhy; 357c4762a1bSJed Brown vleft = -hxhx * (trule2 + trule4); 358c4762a1bSJed Brown vright = -hxhx * (trule1 + trule3); 359c4762a1bSJed Brown vup = -hyhy * (trule1 + trule6); 360c4762a1bSJed Brown vmiddle = (hxhx) * (trule1 + trule2 + trule3 + trule4) + hyhy * (trule1 + trule2 + trule5 + trule6); 3619371c9d4SSatish Balay v[0] = 0; 3629371c9d4SSatish Balay v[1] = 0; 3639371c9d4SSatish Balay v[2] = 0; 3649371c9d4SSatish Balay v[3] = 0; 3659371c9d4SSatish Balay v[4] = 0; 366c4762a1bSJed Brown 367c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 368c4762a1bSJed Brown row = (j - gys) * gxm + (i - gxs); 369c4762a1bSJed Brown 370c4762a1bSJed Brown k = 0; 371c4762a1bSJed Brown if (j > gys) { 3729371c9d4SSatish Balay v[k] = vdown; 3739371c9d4SSatish Balay col[k] = row - gxm; 3749371c9d4SSatish Balay k++; 375c4762a1bSJed Brown } 376c4762a1bSJed Brown 377c4762a1bSJed Brown if (i > gxs) { 3789371c9d4SSatish Balay v[k] = vleft; 3799371c9d4SSatish Balay col[k] = row - 1; 3809371c9d4SSatish Balay k++; 381c4762a1bSJed Brown } 382c4762a1bSJed Brown 3839371c9d4SSatish Balay v[k] = vmiddle; 3849371c9d4SSatish Balay col[k] = row; 3859371c9d4SSatish Balay k++; 386c4762a1bSJed Brown 387c4762a1bSJed Brown if (i + 1 < gxs + gxm) { 3889371c9d4SSatish Balay v[k] = vright; 3899371c9d4SSatish Balay col[k] = row + 1; 3909371c9d4SSatish Balay k++; 391c4762a1bSJed Brown } 392c4762a1bSJed Brown 393c4762a1bSJed Brown if (j + 1 < gys + gym) { 3949371c9d4SSatish Balay v[k] = vup; 3959371c9d4SSatish Balay col[k] = row + gxm; 3969371c9d4SSatish Balay k++; 397c4762a1bSJed Brown } 3989566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(hes, 1, &row, k, col, v, INSERT_VALUES)); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown } 401c4762a1bSJed Brown 402c4762a1bSJed Brown /* 403c4762a1bSJed Brown Assemble matrix, using the 2-step process: 404c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 405c4762a1bSJed Brown By placing code between these two statements, computations can be 406c4762a1bSJed Brown done while messages are in transition. 407c4762a1bSJed Brown */ 4089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(hes, MAT_FINAL_ASSEMBLY)); 4099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(hes, MAT_FINAL_ASSEMBLY)); 410c4762a1bSJed Brown 411c4762a1bSJed Brown /* 412c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 413c4762a1bSJed Brown matrix. If we do it will generate an error. 414c4762a1bSJed Brown */ 4159566063dSJacob Faibussowitsch PetscCall(MatSetOption(hes, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 4169566063dSJacob Faibussowitsch PetscCall(MatSetOption(hes, MAT_SYMMETRIC, PETSC_TRUE)); 417c4762a1bSJed Brown 4189566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0 * xm * ym + 49.0 * xm)); 419780b99b1SStefano Zampini PetscFunctionReturn(0); 420c4762a1bSJed Brown } 421c4762a1bSJed Brown 4229371c9d4SSatish Balay PetscErrorCode Monitor(Tao tao, void *ctx) { 423c4762a1bSJed Brown PetscInt its; 424c4762a1bSJed Brown PetscReal f, gnorm, cnorm, xdiff; 425c4762a1bSJed Brown TaoConvergedReason reason; 426c4762a1bSJed Brown 427c4762a1bSJed Brown PetscFunctionBegin; 4289566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason)); 429*48a46eb9SPierre Jolivet if (!(its % 5)) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration=%" PetscInt_FMT "\tf=%g\n", its, (double)f)); 430c4762a1bSJed Brown PetscFunctionReturn(0); 431c4762a1bSJed Brown } 432c4762a1bSJed Brown 4339371c9d4SSatish Balay PetscErrorCode ConvergenceTest(Tao tao, void *ctx) { 434c4762a1bSJed Brown PetscInt its; 435c4762a1bSJed Brown PetscReal f, gnorm, cnorm, xdiff; 436c4762a1bSJed Brown TaoConvergedReason reason; 437c4762a1bSJed Brown 438c4762a1bSJed Brown PetscFunctionBegin; 4399566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason)); 440*48a46eb9SPierre Jolivet if (its == 100) PetscCall(TaoSetConvergedReason(tao, TAO_DIVERGED_MAXITS)); 441c4762a1bSJed Brown PetscFunctionReturn(0); 442c4762a1bSJed Brown } 443c4762a1bSJed Brown 444c4762a1bSJed Brown /*TEST 445c4762a1bSJed Brown 446c4762a1bSJed Brown build: 447c4762a1bSJed Brown requires: !complex 448c4762a1bSJed Brown 449c4762a1bSJed Brown test: 450c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5 451c4762a1bSJed Brown requires: !single 452c4762a1bSJed Brown 453c4762a1bSJed Brown test: 454c4762a1bSJed Brown suffix: 2 455c4762a1bSJed Brown nsize: 2 456c4762a1bSJed Brown args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5 457c4762a1bSJed Brown requires: !single 458c4762a1bSJed Brown 459c4762a1bSJed Brown test: 460c4762a1bSJed Brown suffix: 3 461c4762a1bSJed Brown nsize: 2 462c4762a1bSJed Brown args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 463c4762a1bSJed Brown requires: !single 464c4762a1bSJed Brown 465c4762a1bSJed Brown test: 466c4762a1bSJed Brown suffix: 4 467c4762a1bSJed Brown nsize: 2 468c4762a1bSJed Brown args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal 469c4762a1bSJed Brown output_file: output/jbearing2_3.out 470c4762a1bSJed Brown requires: !single 471c4762a1bSJed Brown 472c4762a1bSJed Brown test: 473c4762a1bSJed Brown suffix: 5 474c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4 475c4762a1bSJed Brown requires: !single 476c4762a1bSJed Brown 477c4762a1bSJed Brown test: 478c4762a1bSJed Brown suffix: 6 479c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4 480c4762a1bSJed Brown requires: !single 481c4762a1bSJed Brown 482c4762a1bSJed Brown test: 483c4762a1bSJed Brown suffix: 7 484c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 485c4762a1bSJed Brown requires: !single 486c4762a1bSJed Brown 487c4762a1bSJed Brown test: 488c4762a1bSJed Brown suffix: 8 489c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 490c4762a1bSJed Brown requires: !single 491c4762a1bSJed Brown 492c4762a1bSJed Brown test: 493c4762a1bSJed Brown suffix: 9 494c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 495c4762a1bSJed Brown requires: !single 496c4762a1bSJed Brown 497c4762a1bSJed Brown test: 498c4762a1bSJed Brown suffix: 10 499c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 500c4762a1bSJed Brown requires: !single 501c4762a1bSJed Brown 502c4762a1bSJed Brown test: 503c4762a1bSJed Brown suffix: 11 504c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 505c4762a1bSJed Brown requires: !single 506c4762a1bSJed Brown 507c4762a1bSJed Brown test: 508c4762a1bSJed Brown suffix: 12 509c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 510c4762a1bSJed Brown requires: !single 511c4762a1bSJed Brown 512c4762a1bSJed Brown test: 513c4762a1bSJed Brown suffix: 13 514c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls 515c4762a1bSJed Brown requires: !single 516c4762a1bSJed Brown 517c4762a1bSJed Brown test: 518c4762a1bSJed Brown suffix: 14 519c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm 520c4762a1bSJed Brown requires: !single 521c4762a1bSJed Brown 522c4762a1bSJed Brown test: 523c4762a1bSJed Brown suffix: 15 524c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs 525c4762a1bSJed Brown requires: !single 526c4762a1bSJed Brown 527c4762a1bSJed Brown test: 528c4762a1bSJed Brown suffix: 16 529c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 530c4762a1bSJed Brown requires: !single 531c4762a1bSJed Brown 532c4762a1bSJed Brown test: 533c4762a1bSJed Brown suffix: 17 534864588a7SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view 535c4762a1bSJed Brown requires: !single 536c4762a1bSJed Brown 537c4762a1bSJed Brown test: 538c4762a1bSJed Brown suffix: 18 539864588a7SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view 540c4762a1bSJed Brown requires: !single 541c4762a1bSJed Brown 54234ad9904SAlp Dener test: 54334ad9904SAlp Dener suffix: 19 54434ad9904SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian 54534ad9904SAlp Dener requires: !single 54634ad9904SAlp Dener 54734ad9904SAlp Dener test: 54834ad9904SAlp Dener suffix: 20 54934ad9904SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian 55034ad9904SAlp Dener requires: !single 55134ad9904SAlp Dener 55234ad9904SAlp Dener test: 55334ad9904SAlp Dener suffix: 21 55434ad9904SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian 55534ad9904SAlp Dener requires: !single 556c4762a1bSJed Brown TEST*/ 557