1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 eptorsion1 [-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 multiprocessor version of this code is eptorsion2.c. 16c4762a1bSJed Brown 17c4762a1bSJed Brown ---------------------------------------------------------------------- */ 18c4762a1bSJed Brown 19c4762a1bSJed Brown /* 20c4762a1bSJed Brown Include "petsctao.h" so that we can use TAO solvers. Note that this 21c4762a1bSJed Brown file automatically includes files for lower-level support, such as those 22c4762a1bSJed Brown provided by the PETSc library: 23c4762a1bSJed Brown petsc.h - base PETSc routines petscvec.h - vectors 24a5b23f4aSJose E. Roman petscsys.h - system routines petscmat.h - matrices 25c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 26c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 27c4762a1bSJed Brown */ 28c4762a1bSJed Brown 29c4762a1bSJed Brown #include <petsctao.h> 30c4762a1bSJed Brown 319371c9d4SSatish Balay static char help[] = "Demonstrates use of the TAO package to solve \n\ 32c4762a1bSJed Brown unconstrained minimization problems on a single processor. This example \n\ 33c4762a1bSJed Brown is based on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 \n\ 34c4762a1bSJed Brown test suite.\n\ 35c4762a1bSJed Brown The command line options are:\n\ 36c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 37c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 38c4762a1bSJed Brown -par <param>, where <param> = angle of twist per unit length\n\n"; 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* 41c4762a1bSJed Brown User-defined application context - contains data needed by the 42c4762a1bSJed Brown application-provided call-back routines, FormFunction(), 43c4762a1bSJed Brown FormGradient(), and FormHessian(). 44c4762a1bSJed Brown */ 45c4762a1bSJed Brown 46c4762a1bSJed Brown typedef struct { 47c4762a1bSJed Brown PetscReal param; /* nonlinearity parameter */ 48c4762a1bSJed Brown PetscInt mx, my; /* discretization in x- and y-directions */ 49c4762a1bSJed Brown PetscInt ndim; /* problem dimension */ 50c4762a1bSJed Brown Vec s, y, xvec; /* work space for computing Hessian */ 51c4762a1bSJed Brown PetscReal hx, hy; /* mesh spacing in x- and y-directions */ 52c4762a1bSJed Brown } AppCtx; 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *, Vec); 57c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 58c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 59c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 60c4762a1bSJed Brown PetscErrorCode HessianProductMat(Mat, Vec, Vec); 61c4762a1bSJed Brown PetscErrorCode HessianProduct(void *, Vec, Vec); 62c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *); 63c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 64c4762a1bSJed Brown 653ba16761SJacob Faibussowitsch int main(int argc, char **argv) 66d71ae5a4SJacob Faibussowitsch { 67c4762a1bSJed Brown PetscInt mx = 10; /* discretization in x-direction */ 68c4762a1bSJed Brown PetscInt my = 10; /* discretization in y-direction */ 69c4762a1bSJed Brown Vec x; /* solution, gradient vectors */ 70c4762a1bSJed Brown PetscBool flg; /* A return value when checking for use options */ 71c4762a1bSJed Brown Tao tao; /* Tao solver context */ 72c4762a1bSJed Brown Mat H; /* Hessian matrix */ 73c4762a1bSJed Brown AppCtx user; /* application context */ 74c4762a1bSJed Brown PetscMPIInt size; /* number of processes */ 75c4762a1bSJed Brown PetscReal one = 1.0; 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscBool test_lmvm = PETSC_FALSE; 78c4762a1bSJed Brown KSP ksp; 79c4762a1bSJed Brown PC pc; 80c4762a1bSJed Brown Mat M; 81c4762a1bSJed Brown Vec in, out, out2; 82c4762a1bSJed Brown PetscReal mult_solve_dist; 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Initialize TAO,PETSc */ 85327415f7SBarry Smith PetscFunctionBeginUser; 86c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size)); 883c859ba3SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors"); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Specify default parameters for the problem, check for command-line overrides */ 91c4762a1bSJed Brown user.param = 5.0; 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, &flg)); 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, &flg)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, &flg)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg)); 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Elastic-Plastic Torsion Problem -----\n")); 9863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n", mx, my)); 999371c9d4SSatish Balay user.ndim = mx * my; 1009371c9d4SSatish Balay user.mx = mx; 1019371c9d4SSatish Balay user.my = my; 1029371c9d4SSatish Balay user.hx = one / (mx + 1); 1039371c9d4SSatish Balay user.hy = one / (my + 1); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Allocate vectors */ 1069566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.ndim, &user.y)); 1079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.y, &user.s)); 1089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.y, &x)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1119566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 1129566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLMVM)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1159566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user, x)); 1169566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 1199566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* From command line options, determine if using matrix-free hessian */ 1229566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-my_tao_mf", &flg)); 123c4762a1bSJed Brown if (flg) { 1249566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_SELF, user.ndim, user.ndim, user.ndim, user.ndim, (void *)&user, &H)); 125*57d50842SBarry Smith PetscCall(MatShellSetOperation(H, MATOP_MULT, (PetscErrorCodeFn *)HessianProductMat)); 1269566063dSJacob Faibussowitsch PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE)); 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, MatrixFreeHessian, (void *)&user)); 129c4762a1bSJed Brown } else { 1309566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, user.ndim, user.ndim, 5, NULL, &H)); 1319566063dSJacob Faibussowitsch PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE)); 1329566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, FormHessian, (void *)&user)); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Test the LMVM matrix */ 136c4762a1bSJed Brown if (test_lmvm) { 1379566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bntr")); 1389566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-tao_bnk_pc_type", "lmvm")); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Check for any TAO command line options */ 1429566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1459566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Test the LMVM matrix */ 148c4762a1bSJed Brown if (test_lmvm) { 1499566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 1509566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1519566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(pc, &M)); 1529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &in)); 1539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out)); 1549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out2)); 1559566063dSJacob Faibussowitsch PetscCall(VecSet(in, 5.0)); 1569566063dSJacob Faibussowitsch PetscCall(MatMult(M, in, out)); 1579566063dSJacob Faibussowitsch PetscCall(MatSolve(M, out, out2)); 1589566063dSJacob Faibussowitsch PetscCall(VecAXPY(out2, -1.0, in)); 1599566063dSJacob Faibussowitsch PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist)); 16063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between MatMult and MatSolve: %e\n", (double)mult_solve_dist)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&in)); 1629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out)); 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out2)); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.s)); 1689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.y)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 171c4762a1bSJed Brown 1729566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 173b122ec5aSJacob Faibussowitsch return 0; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 177c4762a1bSJed Brown /* 178c4762a1bSJed Brown FormInitialGuess - Computes an initial approximation to the solution. 179c4762a1bSJed Brown 180c4762a1bSJed Brown Input Parameters: 181c4762a1bSJed Brown . user - user-defined application context 182c4762a1bSJed Brown . X - vector 183c4762a1bSJed Brown 184c4762a1bSJed Brown Output Parameters: 185c4762a1bSJed Brown . X - vector 186c4762a1bSJed Brown */ 187d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) 188d71ae5a4SJacob Faibussowitsch { 189c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy, temp; 190c4762a1bSJed Brown PetscReal val; 191c4762a1bSJed Brown PetscInt i, j, k, nx = user->mx, ny = user->my; 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* Compute initial guess */ 194c4762a1bSJed Brown PetscFunctionBeginUser; 195c4762a1bSJed Brown for (j = 0; j < ny; j++) { 196c4762a1bSJed Brown temp = PetscMin(j + 1, ny - j) * hy; 197c4762a1bSJed Brown for (i = 0; i < nx; i++) { 198c4762a1bSJed Brown k = nx * j + i; 199c4762a1bSJed Brown val = PetscMin((PetscMin(i + 1, nx - i)) * hx, temp); 2009566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 1, &k, &val, ADD_VALUES)); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown } 2039566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 2049566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 2053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 209c4762a1bSJed Brown /* 210c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 211c4762a1bSJed Brown 212c4762a1bSJed Brown Input Parameters: 213c4762a1bSJed Brown tao - the Tao context 214c4762a1bSJed Brown X - the input vector 215c4762a1bSJed Brown ptr - optional user-defined context, as set by TaoSetFunction() 216c4762a1bSJed Brown 217c4762a1bSJed Brown Output Parameters: 218c4762a1bSJed Brown f - the newly evaluated function 219c4762a1bSJed Brown G - the newly evaluated gradient 220c4762a1bSJed Brown */ 221d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) 222d71ae5a4SJacob Faibussowitsch { 223c4762a1bSJed Brown PetscFunctionBeginUser; 2249566063dSJacob Faibussowitsch PetscCall(FormFunction(tao, X, f, ptr)); 2259566063dSJacob Faibussowitsch PetscCall(FormGradient(tao, X, G, ptr)); 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 230c4762a1bSJed Brown /* 231c4762a1bSJed Brown FormFunction - Evaluates the function, f(X). 232c4762a1bSJed Brown 233c4762a1bSJed Brown Input Parameters: 234c4762a1bSJed Brown . tao - the Tao context 235c4762a1bSJed Brown . X - the input vector 236c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunction() 237c4762a1bSJed Brown 238c4762a1bSJed Brown Output Parameters: 239c4762a1bSJed Brown . f - the newly evaluated function 240c4762a1bSJed Brown */ 241d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) 242d71ae5a4SJacob Faibussowitsch { 243c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 244c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5; 245c4762a1bSJed Brown PetscReal zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0; 246c4762a1bSJed Brown PetscReal v, cdiv3 = user->param / three; 247c4762a1bSJed Brown const PetscScalar *x; 248c4762a1bSJed Brown PetscInt nx = user->mx, ny = user->my, i, j, k; 249c4762a1bSJed Brown 250c4762a1bSJed Brown PetscFunctionBeginUser; 251c4762a1bSJed Brown /* Get pointer to vector data */ 2529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* Compute function contributions over the lower triangular elements */ 255c4762a1bSJed Brown for (j = -1; j < ny; j++) { 256c4762a1bSJed Brown for (i = -1; i < nx; i++) { 257c4762a1bSJed Brown k = nx * j + i; 258c4762a1bSJed Brown v = zero; 259c4762a1bSJed Brown vr = zero; 260c4762a1bSJed Brown vt = zero; 261c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k]; 262c4762a1bSJed Brown if (i < nx - 1 && j > -1) vr = x[k + 1]; 263c4762a1bSJed Brown if (i > -1 && j < ny - 1) vt = x[k + nx]; 264c4762a1bSJed Brown dvdx = (vr - v) / hx; 265c4762a1bSJed Brown dvdy = (vt - v) / hy; 266c4762a1bSJed Brown fquad += dvdx * dvdx + dvdy * dvdy; 267c4762a1bSJed Brown flin -= cdiv3 * (v + vr + vt); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* Compute function contributions over the upper triangular elements */ 272c4762a1bSJed Brown for (j = 0; j <= ny; j++) { 273c4762a1bSJed Brown for (i = 0; i <= nx; i++) { 274c4762a1bSJed Brown k = nx * j + i; 275c4762a1bSJed Brown vb = zero; 276c4762a1bSJed Brown vl = zero; 277c4762a1bSJed Brown v = zero; 278c4762a1bSJed Brown if (i < nx && j > 0) vb = x[k - nx]; 279c4762a1bSJed Brown if (i > 0 && j < ny) vl = x[k - 1]; 280c4762a1bSJed Brown if (i < nx && j < ny) v = x[k]; 281c4762a1bSJed Brown dvdx = (v - vl) / hx; 282c4762a1bSJed Brown dvdy = (v - vb) / hy; 283c4762a1bSJed Brown fquad = fquad + dvdx * dvdx + dvdy * dvdy; 284c4762a1bSJed Brown flin = flin - cdiv3 * (vb + vl + v); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* Restore vector */ 2899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* Scale the function */ 292c4762a1bSJed Brown area = p5 * hx * hy; 293c4762a1bSJed Brown *f = area * (p5 * fquad + flin); 294c4762a1bSJed Brown 2959566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(24.0 * nx * ny)); 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 300c4762a1bSJed Brown /* 301c4762a1bSJed Brown FormGradient - Evaluates the gradient, G(X). 302c4762a1bSJed Brown 303c4762a1bSJed Brown Input Parameters: 304c4762a1bSJed Brown . tao - the Tao context 305c4762a1bSJed Brown . X - input vector 306c4762a1bSJed Brown . ptr - optional user-defined context 307c4762a1bSJed Brown 308c4762a1bSJed Brown Output Parameters: 309c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 310c4762a1bSJed Brown */ 311d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) 312d71ae5a4SJacob Faibussowitsch { 313c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 314c4762a1bSJed Brown PetscReal zero = 0.0, p5 = 0.5, three = 3.0, area, val; 315c4762a1bSJed Brown PetscInt nx = user->mx, ny = user->my, ind, i, j, k; 316c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy; 317c4762a1bSJed Brown PetscReal vb, vl, vr, vt, dvdx, dvdy; 318c4762a1bSJed Brown PetscReal v, cdiv3 = user->param / three; 319c4762a1bSJed Brown const PetscScalar *x; 320c4762a1bSJed Brown 321c4762a1bSJed Brown PetscFunctionBeginUser; 322c4762a1bSJed Brown /* Initialize gradient to zero */ 3239566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* Get pointer to vector data */ 3269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* Compute gradient contributions over the lower triangular elements */ 329c4762a1bSJed Brown for (j = -1; j < ny; j++) { 330c4762a1bSJed Brown for (i = -1; i < nx; i++) { 331c4762a1bSJed Brown k = nx * j + i; 332c4762a1bSJed Brown v = zero; 333c4762a1bSJed Brown vr = zero; 334c4762a1bSJed Brown vt = zero; 335c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k]; 336c4762a1bSJed Brown if (i < nx - 1 && j > -1) vr = x[k + 1]; 337c4762a1bSJed Brown if (i > -1 && j < ny - 1) vt = x[k + nx]; 338c4762a1bSJed Brown dvdx = (vr - v) / hx; 339c4762a1bSJed Brown dvdy = (vt - v) / hy; 340c4762a1bSJed Brown if (i != -1 && j != -1) { 3419371c9d4SSatish Balay ind = k; 3429371c9d4SSatish Balay val = -dvdx / hx - dvdy / hy - cdiv3; 3439566063dSJacob Faibussowitsch PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown if (i != nx - 1 && j != -1) { 3469371c9d4SSatish Balay ind = k + 1; 3479371c9d4SSatish Balay val = dvdx / hx - cdiv3; 3489566063dSJacob Faibussowitsch PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown if (i != -1 && j != ny - 1) { 3519371c9d4SSatish Balay ind = k + nx; 3529371c9d4SSatish Balay val = dvdy / hy - cdiv3; 3539566063dSJacob Faibussowitsch PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 354c4762a1bSJed Brown } 355c4762a1bSJed Brown } 356c4762a1bSJed Brown } 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* Compute gradient contributions over the upper triangular elements */ 359c4762a1bSJed Brown for (j = 0; j <= ny; j++) { 360c4762a1bSJed Brown for (i = 0; i <= nx; i++) { 361c4762a1bSJed Brown k = nx * j + i; 362c4762a1bSJed Brown vb = zero; 363c4762a1bSJed Brown vl = zero; 364c4762a1bSJed Brown v = zero; 365c4762a1bSJed Brown if (i < nx && j > 0) vb = x[k - nx]; 366c4762a1bSJed Brown if (i > 0 && j < ny) vl = x[k - 1]; 367c4762a1bSJed Brown if (i < nx && j < ny) v = x[k]; 368c4762a1bSJed Brown dvdx = (v - vl) / hx; 369c4762a1bSJed Brown dvdy = (v - vb) / hy; 370c4762a1bSJed Brown if (i != nx && j != 0) { 3719371c9d4SSatish Balay ind = k - nx; 3729371c9d4SSatish Balay val = -dvdy / hy - cdiv3; 3739566063dSJacob Faibussowitsch PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 374c4762a1bSJed Brown } 375c4762a1bSJed Brown if (i != 0 && j != ny) { 3769371c9d4SSatish Balay ind = k - 1; 3779371c9d4SSatish Balay val = -dvdx / hx - cdiv3; 3789566063dSJacob Faibussowitsch PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 379c4762a1bSJed Brown } 380c4762a1bSJed Brown if (i != nx && j != ny) { 3819371c9d4SSatish Balay ind = k; 3829371c9d4SSatish Balay val = dvdx / hx + dvdy / hy - cdiv3; 3839566063dSJacob Faibussowitsch PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 384c4762a1bSJed Brown } 385c4762a1bSJed Brown } 386c4762a1bSJed Brown } 3879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 388c4762a1bSJed Brown 389c4762a1bSJed Brown /* Assemble gradient vector */ 3909566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(G)); 3919566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(G)); 392c4762a1bSJed Brown 393c4762a1bSJed Brown /* Scale the gradient */ 394c4762a1bSJed Brown area = p5 * hx * hy; 3959566063dSJacob Faibussowitsch PetscCall(VecScale(G, area)); 3969566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(24.0 * nx * ny)); 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 401c4762a1bSJed Brown /* 402c4762a1bSJed Brown FormHessian - Forms the Hessian matrix. 403c4762a1bSJed Brown 404c4762a1bSJed Brown Input Parameters: 405c4762a1bSJed Brown . tao - the Tao context 406c4762a1bSJed Brown . X - the input vector 407c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 408c4762a1bSJed Brown 409c4762a1bSJed Brown Output Parameters: 410c4762a1bSJed Brown . H - Hessian matrix 411c4762a1bSJed Brown . PrecH - optionally different preconditioning Hessian 412c4762a1bSJed Brown 413c4762a1bSJed Brown Notes: 414c4762a1bSJed Brown This routine is intended simply as an example of the interface 415c4762a1bSJed Brown to a Hessian evaluation routine. Since this example compute the 416c4762a1bSJed Brown Hessian a column at a time, it is not particularly efficient and 417c4762a1bSJed Brown is not recommended. 418c4762a1bSJed Brown */ 419d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) 420d71ae5a4SJacob Faibussowitsch { 421c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 422c4762a1bSJed Brown PetscInt i, j, ndim = user->ndim; 423c4762a1bSJed Brown PetscReal *y, zero = 0.0, one = 1.0; 424c4762a1bSJed Brown PetscBool assembled; 425c4762a1bSJed Brown 426c4762a1bSJed Brown PetscFunctionBeginUser; 427c4762a1bSJed Brown user->xvec = X; 428c4762a1bSJed Brown 429c4762a1bSJed Brown /* Initialize Hessian entries and work vector to zero */ 4309566063dSJacob Faibussowitsch PetscCall(MatAssembled(H, &assembled)); 4319566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(H)); 432c4762a1bSJed Brown 4339566063dSJacob Faibussowitsch PetscCall(VecSet(user->s, zero)); 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* Loop over matrix columns to compute entries of the Hessian */ 436c4762a1bSJed Brown for (j = 0; j < ndim; j++) { 4379566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->s, 1, &j, &one, INSERT_VALUES)); 4389566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->s)); 4399566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->s)); 440c4762a1bSJed Brown 4419566063dSJacob Faibussowitsch PetscCall(HessianProduct(ptr, user->s, user->y)); 442c4762a1bSJed Brown 4439566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->s, 1, &j, &zero, INSERT_VALUES)); 4449566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->s)); 4459566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->s)); 446c4762a1bSJed Brown 4479566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->y, &y)); 448c4762a1bSJed Brown for (i = 0; i < ndim; i++) { 44948a46eb9SPierre Jolivet if (y[i] != zero) PetscCall(MatSetValues(H, 1, &i, 1, &j, &y[i], ADD_VALUES)); 450c4762a1bSJed Brown } 4519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->y, &y)); 452c4762a1bSJed Brown } 4539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 4549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 4553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 456c4762a1bSJed Brown } 457c4762a1bSJed Brown 458c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 459c4762a1bSJed Brown /* 460c4762a1bSJed Brown MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector 461c4762a1bSJed Brown products. 462c4762a1bSJed Brown 463c4762a1bSJed Brown Input Parameters: 464c4762a1bSJed Brown . tao - the Tao context 465c4762a1bSJed Brown . X - the input vector 466c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 467c4762a1bSJed Brown 468c4762a1bSJed Brown Output Parameters: 469c4762a1bSJed Brown . H - Hessian matrix 470c4762a1bSJed Brown . PrecH - optionally different preconditioning Hessian 471c4762a1bSJed Brown */ 472d71ae5a4SJacob Faibussowitsch PetscErrorCode MatrixFreeHessian(Tao tao, Vec X, Mat H, Mat PrecH, void *ptr) 473d71ae5a4SJacob Faibussowitsch { 474c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 475c4762a1bSJed Brown 476c4762a1bSJed Brown /* Sets location of vector for use in computing matrix-vector products of the form H(X)*y */ 477362febeeSStefano Zampini PetscFunctionBeginUser; 478c4762a1bSJed Brown user->xvec = X; 4793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 480c4762a1bSJed Brown } 481c4762a1bSJed Brown 482c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 483c4762a1bSJed Brown /* 484c4762a1bSJed Brown HessianProductMat - Computes the matrix-vector product 485c4762a1bSJed Brown y = mat*svec. 486c4762a1bSJed Brown 487c4762a1bSJed Brown Input Parameters: 488c4762a1bSJed Brown . mat - input matrix 489c4762a1bSJed Brown . svec - input vector 490c4762a1bSJed Brown 491c4762a1bSJed Brown Output Parameters: 492c4762a1bSJed Brown . y - solution vector 493c4762a1bSJed Brown */ 494d71ae5a4SJacob Faibussowitsch PetscErrorCode HessianProductMat(Mat mat, Vec svec, Vec y) 495d71ae5a4SJacob Faibussowitsch { 496c4762a1bSJed Brown void *ptr; 497c4762a1bSJed Brown 498c4762a1bSJed Brown PetscFunctionBeginUser; 4999566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ptr)); 5009566063dSJacob Faibussowitsch PetscCall(HessianProduct(ptr, svec, y)); 5013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 502c4762a1bSJed Brown } 503c4762a1bSJed Brown 504c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 505c4762a1bSJed Brown /* 506c4762a1bSJed Brown Hessian Product - Computes the matrix-vector product: 507c4762a1bSJed Brown y = f''(x)*svec. 508c4762a1bSJed Brown 5097a7aea1fSJed Brown Input Parameters: 510c4762a1bSJed Brown . ptr - pointer to the user-defined context 511c4762a1bSJed Brown . svec - input vector 512c4762a1bSJed Brown 513c4762a1bSJed Brown Output Parameters: 514c4762a1bSJed Brown . y - product vector 515c4762a1bSJed Brown */ 516d71ae5a4SJacob Faibussowitsch PetscErrorCode HessianProduct(void *ptr, Vec svec, Vec y) 517d71ae5a4SJacob Faibussowitsch { 518c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 519c4762a1bSJed Brown PetscReal p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area; 520c4762a1bSJed Brown const PetscScalar *x, *s; 521c4762a1bSJed Brown PetscReal v, vb, vl, vr, vt, hxhx, hyhy; 522c4762a1bSJed Brown PetscInt nx, ny, i, j, k, ind; 523c4762a1bSJed Brown 524c4762a1bSJed Brown PetscFunctionBeginUser; 525c4762a1bSJed Brown nx = user->mx; 526c4762a1bSJed Brown ny = user->my; 527c4762a1bSJed Brown hx = user->hx; 528c4762a1bSJed Brown hy = user->hy; 529c4762a1bSJed Brown hxhx = one / (hx * hx); 530c4762a1bSJed Brown hyhy = one / (hy * hy); 531c4762a1bSJed Brown 532c4762a1bSJed Brown /* Get pointers to vector data */ 5339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user->xvec, &x)); 5349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(svec, &s)); 535c4762a1bSJed Brown 536c4762a1bSJed Brown /* Initialize product vector to zero */ 5379566063dSJacob Faibussowitsch PetscCall(VecSet(y, zero)); 538c4762a1bSJed Brown 539c4762a1bSJed Brown /* Compute f''(x)*s over the lower triangular elements */ 540c4762a1bSJed Brown for (j = -1; j < ny; j++) { 541c4762a1bSJed Brown for (i = -1; i < nx; i++) { 542c4762a1bSJed Brown k = nx * j + i; 543c4762a1bSJed Brown v = zero; 544c4762a1bSJed Brown vr = zero; 545c4762a1bSJed Brown vt = zero; 546c4762a1bSJed Brown if (i != -1 && j != -1) v = s[k]; 547c4762a1bSJed Brown if (i != nx - 1 && j != -1) { 548c4762a1bSJed Brown vr = s[k + 1]; 5499371c9d4SSatish Balay ind = k + 1; 5509371c9d4SSatish Balay val = hxhx * (vr - v); 5519566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 552c4762a1bSJed Brown } 553c4762a1bSJed Brown if (i != -1 && j != ny - 1) { 554c4762a1bSJed Brown vt = s[k + nx]; 5559371c9d4SSatish Balay ind = k + nx; 5569371c9d4SSatish Balay val = hyhy * (vt - v); 5579566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 558c4762a1bSJed Brown } 559c4762a1bSJed Brown if (i != -1 && j != -1) { 5609371c9d4SSatish Balay ind = k; 5619371c9d4SSatish Balay val = hxhx * (v - vr) + hyhy * (v - vt); 5629566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 563c4762a1bSJed Brown } 564c4762a1bSJed Brown } 565c4762a1bSJed Brown } 566c4762a1bSJed Brown 567c4762a1bSJed Brown /* Compute f''(x)*s over the upper triangular elements */ 568c4762a1bSJed Brown for (j = 0; j <= ny; j++) { 569c4762a1bSJed Brown for (i = 0; i <= nx; i++) { 570c4762a1bSJed Brown k = nx * j + i; 571c4762a1bSJed Brown v = zero; 572c4762a1bSJed Brown vl = zero; 573c4762a1bSJed Brown vb = zero; 574c4762a1bSJed Brown if (i != nx && j != ny) v = s[k]; 575c4762a1bSJed Brown if (i != nx && j != 0) { 576c4762a1bSJed Brown vb = s[k - nx]; 5779371c9d4SSatish Balay ind = k - nx; 5789371c9d4SSatish Balay val = hyhy * (vb - v); 5799566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 580c4762a1bSJed Brown } 581c4762a1bSJed Brown if (i != 0 && j != ny) { 582c4762a1bSJed Brown vl = s[k - 1]; 5839371c9d4SSatish Balay ind = k - 1; 5849371c9d4SSatish Balay val = hxhx * (vl - v); 5859566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 586c4762a1bSJed Brown } 587c4762a1bSJed Brown if (i != nx && j != ny) { 5889371c9d4SSatish Balay ind = k; 5899371c9d4SSatish Balay val = hxhx * (v - vl) + hyhy * (v - vb); 5909566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 591c4762a1bSJed Brown } 592c4762a1bSJed Brown } 593c4762a1bSJed Brown } 594c4762a1bSJed Brown /* Restore vector data */ 5959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(svec, &s)); 5969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user->xvec, &x)); 597c4762a1bSJed Brown 598c4762a1bSJed Brown /* Assemble vector */ 5999566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(y)); 6009566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(y)); 601c4762a1bSJed Brown 602c4762a1bSJed Brown /* Scale resulting vector by area */ 603c4762a1bSJed Brown area = p5 * hx * hy; 6049566063dSJacob Faibussowitsch PetscCall(VecScale(y, area)); 6059566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(18.0 * nx * ny)); 6063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 607c4762a1bSJed Brown } 608c4762a1bSJed Brown 609c4762a1bSJed Brown /*TEST 610c4762a1bSJed Brown 611c4762a1bSJed Brown build: 612c4762a1bSJed Brown requires: !complex 613c4762a1bSJed Brown 614c4762a1bSJed Brown test: 615c4762a1bSJed Brown suffix: 1 61610978b7dSBarry Smith args: -tao_monitor_short -tao_type ntl -tao_gatol 1.e-4 617c4762a1bSJed Brown 618c4762a1bSJed Brown test: 619c4762a1bSJed Brown suffix: 2 62010978b7dSBarry Smith args: -tao_monitor_short -tao_type ntr -tao_gatol 1.e-4 621c4762a1bSJed Brown 622c4762a1bSJed Brown test: 623c4762a1bSJed Brown suffix: 3 62410978b7dSBarry Smith args: -tao_monitor_short -tao_type bntr -tao_gatol 1.e-4 -my_tao_mf -tao_test_hessian 625c4762a1bSJed Brown 626c4762a1bSJed Brown test: 627c4762a1bSJed Brown suffix: 4 62810978b7dSBarry Smith args: -tao_monitor_short -tao_gatol 1e-3 -tao_type bqnls 629c4762a1bSJed Brown 630c4762a1bSJed Brown test: 631c4762a1bSJed Brown suffix: 5 63210978b7dSBarry Smith args: -tao_monitor_short -tao_gatol 1e-3 -tao_type blmvm 633c4762a1bSJed Brown 634c4762a1bSJed Brown test: 635c4762a1bSJed Brown suffix: 6 63610978b7dSBarry Smith args: -tao_monitor_short -tao_gatol 1e-3 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 637c4762a1bSJed Brown 638f4f59681SStefano Zampini test: 639f4f59681SStefano Zampini suffix: snes 640f4f59681SStefano Zampini args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtontr -ksp_type cg -snes_atol 1.e-4 -tao_mf_hessian {{0 1}} -pc_type none 641f4f59681SStefano Zampini 642f4f59681SStefano Zampini test: 643f4f59681SStefano Zampini suffix: snes_2 644f4f59681SStefano Zampini args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtontr -snes_atol 5.e-4 -tao_mf_hessian -pc_type none -snes_tr_fallback_type cauchy 645f4f59681SStefano Zampini 646a0254a93SStefano Zampini test: 647a0254a93SStefano Zampini suffix: snes_3 648a0254a93SStefano Zampini args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtontr -snes_atol 5.e-4 -tao_mf_hessian -pc_type lmvm -snes_tr_fallback_type cauchy 649a0254a93SStefano Zampini 650c4762a1bSJed Brown TEST*/ 651