1*b3e6a353SBarry Smith /* Program usage: mpiexec -n 1 eptorsion1 [-help] [all TAO options] */ 2*b3e6a353SBarry Smith 3*b3e6a353SBarry Smith /* ---------------------------------------------------------------------- 4*b3e6a353SBarry Smith 5*b3e6a353SBarry Smith Elastic-plastic torsion problem. 6*b3e6a353SBarry Smith 7*b3e6a353SBarry Smith The elastic plastic torsion problem arises from the determination 8*b3e6a353SBarry Smith of the stress field on an infinitely long cylindrical bar, which is 9*b3e6a353SBarry Smith equivalent to the solution of the following problem: 10*b3e6a353SBarry Smith 11*b3e6a353SBarry Smith min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)} 12*b3e6a353SBarry Smith 13*b3e6a353SBarry Smith where C is the torsion angle per unit length. 14*b3e6a353SBarry Smith 15*b3e6a353SBarry Smith The multiprocessor version of this code is eptorsion2.c. 16*b3e6a353SBarry Smith 17*b3e6a353SBarry Smith ---------------------------------------------------------------------- */ 18*b3e6a353SBarry Smith 19*b3e6a353SBarry Smith /* 20*b3e6a353SBarry Smith Include "petsctao.h" so that we can use TAO solvers. Note that this 21*b3e6a353SBarry Smith file automatically includes files for lower-level support, such as those 22*b3e6a353SBarry Smith provided by the PETSc library: 23*b3e6a353SBarry Smith petsc.h - base PETSc routines petscvec.h - vectors 24*b3e6a353SBarry Smith petscsys.h - system routines petscmat.h - matrices 25*b3e6a353SBarry Smith petscis.h - index sets petscksp.h - Krylov subspace methods 26*b3e6a353SBarry Smith petscviewer.h - viewers petscpc.h - preconditioners 27*b3e6a353SBarry Smith */ 28*b3e6a353SBarry Smith 29*b3e6a353SBarry Smith #include <petsctao.h> 30*b3e6a353SBarry Smith 31*b3e6a353SBarry Smith static char help[] = "Demonstrates use of the TAO package to solve \n\ 32*b3e6a353SBarry Smith unconstrained minimization problems on a single processor. This example \n\ 33*b3e6a353SBarry Smith is based on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 \n\ 34*b3e6a353SBarry Smith test suite.\n\ 35*b3e6a353SBarry Smith The command line options are:\n\ 36*b3e6a353SBarry Smith -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 37*b3e6a353SBarry Smith -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 38*b3e6a353SBarry Smith -par <param>, where <param> = angle of twist per unit length\n\n"; 39*b3e6a353SBarry Smith 40*b3e6a353SBarry Smith /* 41*b3e6a353SBarry Smith User-defined application context - contains data needed by the 42*b3e6a353SBarry Smith application-provided call-back routines, FormFunction(), 43*b3e6a353SBarry Smith FormGradient(), and FormHessian(). 44*b3e6a353SBarry Smith */ 45*b3e6a353SBarry Smith 46*b3e6a353SBarry Smith typedef struct { 47*b3e6a353SBarry Smith PetscReal param; /* nonlinearity parameter */ 48*b3e6a353SBarry Smith PetscInt mx, my; /* discretization in x- and y-directions */ 49*b3e6a353SBarry Smith PetscInt ndim; /* problem dimension */ 50*b3e6a353SBarry Smith Vec s, y, xvec; /* work space for computing Hessian */ 51*b3e6a353SBarry Smith PetscReal hx, hy; /* mesh spacing in x- and y-directions */ 52*b3e6a353SBarry Smith } AppCtx; 53*b3e6a353SBarry Smith 54*b3e6a353SBarry Smith /* -------- User-defined Routines --------- */ 55*b3e6a353SBarry Smith 56*b3e6a353SBarry Smith PetscErrorCode FormInitialGuess(AppCtx *, Vec); 57*b3e6a353SBarry Smith PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 58*b3e6a353SBarry Smith PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 59*b3e6a353SBarry Smith PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 60*b3e6a353SBarry Smith PetscErrorCode HessianProductMat(Mat, Vec, Vec); 61*b3e6a353SBarry Smith PetscErrorCode HessianProduct(void *, Vec, Vec); 62*b3e6a353SBarry Smith PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *); 63*b3e6a353SBarry Smith PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 64*b3e6a353SBarry Smith 65*b3e6a353SBarry Smith PetscErrorCode main(int argc, char **argv) 66*b3e6a353SBarry Smith { 67*b3e6a353SBarry Smith PetscInt mx = 10; /* discretization in x-direction */ 68*b3e6a353SBarry Smith PetscInt my = 10; /* discretization in y-direction */ 69*b3e6a353SBarry Smith Vec x; /* solution, gradient vectors */ 70*b3e6a353SBarry Smith PetscBool flg; /* A return value when checking for use options */ 71*b3e6a353SBarry Smith Tao tao; /* Tao solver context */ 72*b3e6a353SBarry Smith Mat H; /* Hessian matrix */ 73*b3e6a353SBarry Smith AppCtx user; /* application context */ 74*b3e6a353SBarry Smith PetscMPIInt size; /* number of processes */ 75*b3e6a353SBarry Smith PetscReal one = 1.0; 76*b3e6a353SBarry Smith 77*b3e6a353SBarry Smith PetscBool test_lmvm = PETSC_FALSE; 78*b3e6a353SBarry Smith KSP ksp; 79*b3e6a353SBarry Smith PC pc; 80*b3e6a353SBarry Smith Mat M; 81*b3e6a353SBarry Smith Vec in, out, out2; 82*b3e6a353SBarry Smith PetscReal mult_solve_dist; 83*b3e6a353SBarry Smith Vec ub, lb; 84*b3e6a353SBarry Smith PetscInt i = 3; 85*b3e6a353SBarry Smith 86*b3e6a353SBarry Smith /* Initialize TAO,PETSc */ 87*b3e6a353SBarry Smith PetscFunctionBeginUser; 88*b3e6a353SBarry Smith PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 89*b3e6a353SBarry Smith PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size)); 90*b3e6a353SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors"); 91*b3e6a353SBarry Smith 92*b3e6a353SBarry Smith /* Specify default parameters for the problem, check for command-line overrides */ 93*b3e6a353SBarry Smith user.param = 5.0; 94*b3e6a353SBarry Smith PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, &flg)); 95*b3e6a353SBarry Smith PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, &flg)); 96*b3e6a353SBarry Smith PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, &flg)); 97*b3e6a353SBarry Smith PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg)); 98*b3e6a353SBarry Smith 99*b3e6a353SBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Elastic-Plastic Torsion Problem -----\n")); 100*b3e6a353SBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n", mx, my)); 101*b3e6a353SBarry Smith user.ndim = mx * my; 102*b3e6a353SBarry Smith user.mx = mx; 103*b3e6a353SBarry Smith user.my = my; 104*b3e6a353SBarry Smith user.hx = one / (mx + 1); 105*b3e6a353SBarry Smith user.hy = one / (my + 1); 106*b3e6a353SBarry Smith 107*b3e6a353SBarry Smith /* Allocate vectors */ 108*b3e6a353SBarry Smith PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.ndim, &user.y)); 109*b3e6a353SBarry Smith PetscCall(VecDuplicate(user.y, &user.s)); 110*b3e6a353SBarry Smith PetscCall(VecDuplicate(user.y, &x)); 111*b3e6a353SBarry Smith PetscCall(VecDuplicate(user.y, &lb)); 112*b3e6a353SBarry Smith PetscCall(VecDuplicate(user.y, &ub)); 113*b3e6a353SBarry Smith 114*b3e6a353SBarry Smith PetscCall(VecSet(ub, 0.1)); 115*b3e6a353SBarry Smith PetscCall(VecSet(lb, -0.1)); 116*b3e6a353SBarry Smith PetscCall(VecSetValue(ub, i, 0., INSERT_VALUES)); 117*b3e6a353SBarry Smith PetscCall(VecSetValue(lb, i, 0., INSERT_VALUES)); 118*b3e6a353SBarry Smith PetscCall(VecAssemblyBegin(ub)); 119*b3e6a353SBarry Smith PetscCall(VecAssemblyBegin(lb)); 120*b3e6a353SBarry Smith PetscCall(VecAssemblyEnd(ub)); 121*b3e6a353SBarry Smith PetscCall(VecAssemblyEnd(lb)); 122*b3e6a353SBarry Smith 123*b3e6a353SBarry Smith /* Create TAO solver and set desired solution method */ 124*b3e6a353SBarry Smith PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 125*b3e6a353SBarry Smith PetscCall(TaoSetType(tao, TAOLMVM)); 126*b3e6a353SBarry Smith 127*b3e6a353SBarry Smith /* Set solution vector with an initial guess */ 128*b3e6a353SBarry Smith PetscCall(FormInitialGuess(&user, x)); 129*b3e6a353SBarry Smith PetscCall(TaoSetSolution(tao, x)); 130*b3e6a353SBarry Smith PetscCall(TaoSetVariableBounds(tao, lb, ub)); 131*b3e6a353SBarry Smith 132*b3e6a353SBarry Smith /* Set routine for function and gradient evaluation */ 133*b3e6a353SBarry Smith PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 134*b3e6a353SBarry Smith 135*b3e6a353SBarry Smith /* From command line options, determine if using matrix-free hessian */ 136*b3e6a353SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-my_tao_mf", &flg)); 137*b3e6a353SBarry Smith if (flg) { 138*b3e6a353SBarry Smith PetscCall(MatCreateShell(PETSC_COMM_SELF, user.ndim, user.ndim, user.ndim, user.ndim, (void *)&user, &H)); 139*b3e6a353SBarry Smith PetscCall(MatShellSetOperation(H, MATOP_MULT, (void (*)(void))HessianProductMat)); 140*b3e6a353SBarry Smith PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE)); 141*b3e6a353SBarry Smith 142*b3e6a353SBarry Smith PetscCall(TaoSetHessian(tao, H, H, MatrixFreeHessian, (void *)&user)); 143*b3e6a353SBarry Smith } else { 144*b3e6a353SBarry Smith PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, user.ndim, user.ndim, 5, NULL, &H)); 145*b3e6a353SBarry Smith PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE)); 146*b3e6a353SBarry Smith PetscCall(TaoSetHessian(tao, H, H, FormHessian, (void *)&user)); 147*b3e6a353SBarry Smith } 148*b3e6a353SBarry Smith 149*b3e6a353SBarry Smith /* Test the LMVM matrix */ 150*b3e6a353SBarry Smith if (test_lmvm) { 151*b3e6a353SBarry Smith PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bntr")); 152*b3e6a353SBarry Smith PetscCall(PetscOptionsSetValue(NULL, "-tao_bnk_pc_type", "lmvm")); 153*b3e6a353SBarry Smith } 154*b3e6a353SBarry Smith 155*b3e6a353SBarry Smith /* Check for any TAO command line options */ 156*b3e6a353SBarry Smith PetscCall(TaoSetFromOptions(tao)); 157*b3e6a353SBarry Smith 158*b3e6a353SBarry Smith /* SOLVE THE APPLICATION */ 159*b3e6a353SBarry Smith PetscCall(TaoSolve(tao)); 160*b3e6a353SBarry Smith 161*b3e6a353SBarry Smith /* Test the LMVM matrix */ 162*b3e6a353SBarry Smith if (test_lmvm) { 163*b3e6a353SBarry Smith PetscCall(TaoGetKSP(tao, &ksp)); 164*b3e6a353SBarry Smith PetscCall(KSPGetPC(ksp, &pc)); 165*b3e6a353SBarry Smith PetscCall(PCLMVMGetMatLMVM(pc, &M)); 166*b3e6a353SBarry Smith PetscCall(VecDuplicate(x, &in)); 167*b3e6a353SBarry Smith PetscCall(VecDuplicate(x, &out)); 168*b3e6a353SBarry Smith PetscCall(VecDuplicate(x, &out2)); 169*b3e6a353SBarry Smith PetscCall(VecSet(in, 5.0)); 170*b3e6a353SBarry Smith PetscCall(MatMult(M, in, out)); 171*b3e6a353SBarry Smith PetscCall(MatSolve(M, out, out2)); 172*b3e6a353SBarry Smith PetscCall(VecAXPY(out2, -1.0, in)); 173*b3e6a353SBarry Smith PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist)); 174*b3e6a353SBarry Smith PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between MatMult and MatSolve: %e\n", (double)mult_solve_dist)); 175*b3e6a353SBarry Smith PetscCall(VecDestroy(&in)); 176*b3e6a353SBarry Smith PetscCall(VecDestroy(&out)); 177*b3e6a353SBarry Smith PetscCall(VecDestroy(&out2)); 178*b3e6a353SBarry Smith } 179*b3e6a353SBarry Smith 180*b3e6a353SBarry Smith PetscCall(TaoDestroy(&tao)); 181*b3e6a353SBarry Smith PetscCall(VecDestroy(&user.s)); 182*b3e6a353SBarry Smith PetscCall(VecDestroy(&user.y)); 183*b3e6a353SBarry Smith PetscCall(VecDestroy(&x)); 184*b3e6a353SBarry Smith PetscCall(MatDestroy(&H)); 185*b3e6a353SBarry Smith PetscCall(VecDestroy(&lb)); 186*b3e6a353SBarry Smith PetscCall(VecDestroy(&ub)); 187*b3e6a353SBarry Smith 188*b3e6a353SBarry Smith PetscCall(PetscFinalize()); 189*b3e6a353SBarry Smith return 0; 190*b3e6a353SBarry Smith } 191*b3e6a353SBarry Smith 192*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 193*b3e6a353SBarry Smith /* 194*b3e6a353SBarry Smith FormInitialGuess - Computes an initial approximation to the solution. 195*b3e6a353SBarry Smith 196*b3e6a353SBarry Smith Input Parameters: 197*b3e6a353SBarry Smith . user - user-defined application context 198*b3e6a353SBarry Smith . X - vector 199*b3e6a353SBarry Smith 200*b3e6a353SBarry Smith Output Parameters: 201*b3e6a353SBarry Smith . X - vector 202*b3e6a353SBarry Smith */ 203*b3e6a353SBarry Smith PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) 204*b3e6a353SBarry Smith { 205*b3e6a353SBarry Smith PetscReal hx = user->hx, hy = user->hy, temp; 206*b3e6a353SBarry Smith PetscReal val; 207*b3e6a353SBarry Smith PetscInt i, j, k, nx = user->mx, ny = user->my; 208*b3e6a353SBarry Smith 209*b3e6a353SBarry Smith /* Compute initial guess */ 210*b3e6a353SBarry Smith PetscFunctionBeginUser; 211*b3e6a353SBarry Smith for (j = 0; j < ny; j++) { 212*b3e6a353SBarry Smith temp = PetscMin(j + 1, ny - j) * hy; 213*b3e6a353SBarry Smith for (i = 0; i < nx; i++) { 214*b3e6a353SBarry Smith k = nx * j + i; 215*b3e6a353SBarry Smith val = PetscMin((PetscMin(i + 1, nx - i)) * hx, temp); 216*b3e6a353SBarry Smith PetscCall(VecSetValues(X, 1, &k, &val, ADD_VALUES)); 217*b3e6a353SBarry Smith } 218*b3e6a353SBarry Smith } 219*b3e6a353SBarry Smith PetscCall(VecAssemblyBegin(X)); 220*b3e6a353SBarry Smith PetscCall(VecAssemblyEnd(X)); 221*b3e6a353SBarry Smith PetscFunctionReturn(0); 222*b3e6a353SBarry Smith } 223*b3e6a353SBarry Smith 224*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 225*b3e6a353SBarry Smith /* 226*b3e6a353SBarry Smith FormFunctionGradient - Evaluates the function and corresponding gradient. 227*b3e6a353SBarry Smith 228*b3e6a353SBarry Smith Input Parameters: 229*b3e6a353SBarry Smith tao - the Tao context 230*b3e6a353SBarry Smith X - the input vector 231*b3e6a353SBarry Smith ptr - optional user-defined context, as set by TaoSetFunction() 232*b3e6a353SBarry Smith 233*b3e6a353SBarry Smith Output Parameters: 234*b3e6a353SBarry Smith f - the newly evaluated function 235*b3e6a353SBarry Smith G - the newly evaluated gradient 236*b3e6a353SBarry Smith */ 237*b3e6a353SBarry Smith PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) 238*b3e6a353SBarry Smith { 239*b3e6a353SBarry Smith PetscFunctionBeginUser; 240*b3e6a353SBarry Smith PetscCall(FormFunction(tao, X, f, ptr)); 241*b3e6a353SBarry Smith PetscCall(FormGradient(tao, X, G, ptr)); 242*b3e6a353SBarry Smith PetscFunctionReturn(0); 243*b3e6a353SBarry Smith } 244*b3e6a353SBarry Smith 245*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 246*b3e6a353SBarry Smith /* 247*b3e6a353SBarry Smith FormFunction - Evaluates the function, f(X). 248*b3e6a353SBarry Smith 249*b3e6a353SBarry Smith Input Parameters: 250*b3e6a353SBarry Smith . tao - the Tao context 251*b3e6a353SBarry Smith . X - the input vector 252*b3e6a353SBarry Smith . ptr - optional user-defined context, as set by TaoSetFunction() 253*b3e6a353SBarry Smith 254*b3e6a353SBarry Smith Output Parameters: 255*b3e6a353SBarry Smith . f - the newly evaluated function 256*b3e6a353SBarry Smith */ 257*b3e6a353SBarry Smith PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) 258*b3e6a353SBarry Smith { 259*b3e6a353SBarry Smith AppCtx *user = (AppCtx *)ptr; 260*b3e6a353SBarry Smith PetscReal hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5; 261*b3e6a353SBarry Smith PetscReal zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0; 262*b3e6a353SBarry Smith PetscReal v, cdiv3 = user->param / three; 263*b3e6a353SBarry Smith const PetscScalar *x; 264*b3e6a353SBarry Smith PetscInt nx = user->mx, ny = user->my, i, j, k; 265*b3e6a353SBarry Smith 266*b3e6a353SBarry Smith PetscFunctionBeginUser; 267*b3e6a353SBarry Smith /* Get pointer to vector data */ 268*b3e6a353SBarry Smith PetscCall(VecGetArrayRead(X, &x)); 269*b3e6a353SBarry Smith 270*b3e6a353SBarry Smith /* Compute function contributions over the lower triangular elements */ 271*b3e6a353SBarry Smith for (j = -1; j < ny; j++) { 272*b3e6a353SBarry Smith for (i = -1; i < nx; i++) { 273*b3e6a353SBarry Smith k = nx * j + i; 274*b3e6a353SBarry Smith v = zero; 275*b3e6a353SBarry Smith vr = zero; 276*b3e6a353SBarry Smith vt = zero; 277*b3e6a353SBarry Smith if (i >= 0 && j >= 0) v = x[k]; 278*b3e6a353SBarry Smith if (i < nx - 1 && j > -1) vr = x[k + 1]; 279*b3e6a353SBarry Smith if (i > -1 && j < ny - 1) vt = x[k + nx]; 280*b3e6a353SBarry Smith dvdx = (vr - v) / hx; 281*b3e6a353SBarry Smith dvdy = (vt - v) / hy; 282*b3e6a353SBarry Smith fquad += dvdx * dvdx + dvdy * dvdy; 283*b3e6a353SBarry Smith flin -= cdiv3 * (v + vr + vt); 284*b3e6a353SBarry Smith } 285*b3e6a353SBarry Smith } 286*b3e6a353SBarry Smith 287*b3e6a353SBarry Smith /* Compute function contributions over the upper triangular elements */ 288*b3e6a353SBarry Smith for (j = 0; j <= ny; j++) { 289*b3e6a353SBarry Smith for (i = 0; i <= nx; i++) { 290*b3e6a353SBarry Smith k = nx * j + i; 291*b3e6a353SBarry Smith vb = zero; 292*b3e6a353SBarry Smith vl = zero; 293*b3e6a353SBarry Smith v = zero; 294*b3e6a353SBarry Smith if (i < nx && j > 0) vb = x[k - nx]; 295*b3e6a353SBarry Smith if (i > 0 && j < ny) vl = x[k - 1]; 296*b3e6a353SBarry Smith if (i < nx && j < ny) v = x[k]; 297*b3e6a353SBarry Smith dvdx = (v - vl) / hx; 298*b3e6a353SBarry Smith dvdy = (v - vb) / hy; 299*b3e6a353SBarry Smith fquad = fquad + dvdx * dvdx + dvdy * dvdy; 300*b3e6a353SBarry Smith flin = flin - cdiv3 * (vb + vl + v); 301*b3e6a353SBarry Smith } 302*b3e6a353SBarry Smith } 303*b3e6a353SBarry Smith 304*b3e6a353SBarry Smith /* Restore vector */ 305*b3e6a353SBarry Smith PetscCall(VecRestoreArrayRead(X, &x)); 306*b3e6a353SBarry Smith 307*b3e6a353SBarry Smith /* Scale the function */ 308*b3e6a353SBarry Smith area = p5 * hx * hy; 309*b3e6a353SBarry Smith *f = area * (p5 * fquad + flin); 310*b3e6a353SBarry Smith 311*b3e6a353SBarry Smith PetscCall(PetscLogFlops(24.0 * nx * ny)); 312*b3e6a353SBarry Smith PetscFunctionReturn(0); 313*b3e6a353SBarry Smith } 314*b3e6a353SBarry Smith 315*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 316*b3e6a353SBarry Smith /* 317*b3e6a353SBarry Smith FormGradient - Evaluates the gradient, G(X). 318*b3e6a353SBarry Smith 319*b3e6a353SBarry Smith Input Parameters: 320*b3e6a353SBarry Smith . tao - the Tao context 321*b3e6a353SBarry Smith . X - input vector 322*b3e6a353SBarry Smith . ptr - optional user-defined context 323*b3e6a353SBarry Smith 324*b3e6a353SBarry Smith Output Parameters: 325*b3e6a353SBarry Smith . G - vector containing the newly evaluated gradient 326*b3e6a353SBarry Smith */ 327*b3e6a353SBarry Smith PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) 328*b3e6a353SBarry Smith { 329*b3e6a353SBarry Smith AppCtx *user = (AppCtx *)ptr; 330*b3e6a353SBarry Smith PetscReal zero = 0.0, p5 = 0.5, three = 3.0, area, val; 331*b3e6a353SBarry Smith PetscInt nx = user->mx, ny = user->my, ind, i, j, k; 332*b3e6a353SBarry Smith PetscReal hx = user->hx, hy = user->hy; 333*b3e6a353SBarry Smith PetscReal vb, vl, vr, vt, dvdx, dvdy; 334*b3e6a353SBarry Smith PetscReal v, cdiv3 = user->param / three; 335*b3e6a353SBarry Smith const PetscScalar *x; 336*b3e6a353SBarry Smith 337*b3e6a353SBarry Smith PetscFunctionBeginUser; 338*b3e6a353SBarry Smith /* Initialize gradient to zero */ 339*b3e6a353SBarry Smith PetscCall(VecSet(G, zero)); 340*b3e6a353SBarry Smith 341*b3e6a353SBarry Smith /* Get pointer to vector data */ 342*b3e6a353SBarry Smith PetscCall(VecGetArrayRead(X, &x)); 343*b3e6a353SBarry Smith 344*b3e6a353SBarry Smith /* Compute gradient contributions over the lower triangular elements */ 345*b3e6a353SBarry Smith for (j = -1; j < ny; j++) { 346*b3e6a353SBarry Smith for (i = -1; i < nx; i++) { 347*b3e6a353SBarry Smith k = nx * j + i; 348*b3e6a353SBarry Smith v = zero; 349*b3e6a353SBarry Smith vr = zero; 350*b3e6a353SBarry Smith vt = zero; 351*b3e6a353SBarry Smith if (i >= 0 && j >= 0) v = x[k]; 352*b3e6a353SBarry Smith if (i < nx - 1 && j > -1) vr = x[k + 1]; 353*b3e6a353SBarry Smith if (i > -1 && j < ny - 1) vt = x[k + nx]; 354*b3e6a353SBarry Smith dvdx = (vr - v) / hx; 355*b3e6a353SBarry Smith dvdy = (vt - v) / hy; 356*b3e6a353SBarry Smith if (i != -1 && j != -1) { 357*b3e6a353SBarry Smith ind = k; 358*b3e6a353SBarry Smith val = -dvdx / hx - dvdy / hy - cdiv3; 359*b3e6a353SBarry Smith PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 360*b3e6a353SBarry Smith } 361*b3e6a353SBarry Smith if (i != nx - 1 && j != -1) { 362*b3e6a353SBarry Smith ind = k + 1; 363*b3e6a353SBarry Smith val = dvdx / hx - cdiv3; 364*b3e6a353SBarry Smith PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 365*b3e6a353SBarry Smith } 366*b3e6a353SBarry Smith if (i != -1 && j != ny - 1) { 367*b3e6a353SBarry Smith ind = k + nx; 368*b3e6a353SBarry Smith val = dvdy / hy - cdiv3; 369*b3e6a353SBarry Smith PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 370*b3e6a353SBarry Smith } 371*b3e6a353SBarry Smith } 372*b3e6a353SBarry Smith } 373*b3e6a353SBarry Smith 374*b3e6a353SBarry Smith /* Compute gradient contributions over the upper triangular elements */ 375*b3e6a353SBarry Smith for (j = 0; j <= ny; j++) { 376*b3e6a353SBarry Smith for (i = 0; i <= nx; i++) { 377*b3e6a353SBarry Smith k = nx * j + i; 378*b3e6a353SBarry Smith vb = zero; 379*b3e6a353SBarry Smith vl = zero; 380*b3e6a353SBarry Smith v = zero; 381*b3e6a353SBarry Smith if (i < nx && j > 0) vb = x[k - nx]; 382*b3e6a353SBarry Smith if (i > 0 && j < ny) vl = x[k - 1]; 383*b3e6a353SBarry Smith if (i < nx && j < ny) v = x[k]; 384*b3e6a353SBarry Smith dvdx = (v - vl) / hx; 385*b3e6a353SBarry Smith dvdy = (v - vb) / hy; 386*b3e6a353SBarry Smith if (i != nx && j != 0) { 387*b3e6a353SBarry Smith ind = k - nx; 388*b3e6a353SBarry Smith val = -dvdy / hy - cdiv3; 389*b3e6a353SBarry Smith PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 390*b3e6a353SBarry Smith } 391*b3e6a353SBarry Smith if (i != 0 && j != ny) { 392*b3e6a353SBarry Smith ind = k - 1; 393*b3e6a353SBarry Smith val = -dvdx / hx - cdiv3; 394*b3e6a353SBarry Smith PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 395*b3e6a353SBarry Smith } 396*b3e6a353SBarry Smith if (i != nx && j != ny) { 397*b3e6a353SBarry Smith ind = k; 398*b3e6a353SBarry Smith val = dvdx / hx + dvdy / hy - cdiv3; 399*b3e6a353SBarry Smith PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 400*b3e6a353SBarry Smith } 401*b3e6a353SBarry Smith } 402*b3e6a353SBarry Smith } 403*b3e6a353SBarry Smith PetscCall(VecRestoreArrayRead(X, &x)); 404*b3e6a353SBarry Smith 405*b3e6a353SBarry Smith /* Assemble gradient vector */ 406*b3e6a353SBarry Smith PetscCall(VecAssemblyBegin(G)); 407*b3e6a353SBarry Smith PetscCall(VecAssemblyEnd(G)); 408*b3e6a353SBarry Smith 409*b3e6a353SBarry Smith /* Scale the gradient */ 410*b3e6a353SBarry Smith area = p5 * hx * hy; 411*b3e6a353SBarry Smith PetscCall(VecScale(G, area)); 412*b3e6a353SBarry Smith PetscCall(PetscLogFlops(24.0 * nx * ny)); 413*b3e6a353SBarry Smith PetscFunctionReturn(0); 414*b3e6a353SBarry Smith } 415*b3e6a353SBarry Smith 416*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 417*b3e6a353SBarry Smith /* 418*b3e6a353SBarry Smith FormHessian - Forms the Hessian matrix. 419*b3e6a353SBarry Smith 420*b3e6a353SBarry Smith Input Parameters: 421*b3e6a353SBarry Smith . tao - the Tao context 422*b3e6a353SBarry Smith . X - the input vector 423*b3e6a353SBarry Smith . ptr - optional user-defined context, as set by TaoSetHessian() 424*b3e6a353SBarry Smith 425*b3e6a353SBarry Smith Output Parameters: 426*b3e6a353SBarry Smith . H - Hessian matrix 427*b3e6a353SBarry Smith . PrecH - optionally different preconditioning Hessian 428*b3e6a353SBarry Smith . flag - flag indicating matrix structure 429*b3e6a353SBarry Smith 430*b3e6a353SBarry Smith Notes: 431*b3e6a353SBarry Smith This routine is intended simply as an example of the interface 432*b3e6a353SBarry Smith to a Hessian evaluation routine. Since this example compute the 433*b3e6a353SBarry Smith Hessian a column at a time, it is not particularly efficient and 434*b3e6a353SBarry Smith is not recommended. 435*b3e6a353SBarry Smith */ 436*b3e6a353SBarry Smith PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) 437*b3e6a353SBarry Smith { 438*b3e6a353SBarry Smith AppCtx *user = (AppCtx *)ptr; 439*b3e6a353SBarry Smith PetscInt i, j, ndim = user->ndim; 440*b3e6a353SBarry Smith PetscReal *y, zero = 0.0, one = 1.0; 441*b3e6a353SBarry Smith PetscBool assembled; 442*b3e6a353SBarry Smith 443*b3e6a353SBarry Smith PetscFunctionBeginUser; 444*b3e6a353SBarry Smith user->xvec = X; 445*b3e6a353SBarry Smith 446*b3e6a353SBarry Smith /* Initialize Hessian entries and work vector to zero */ 447*b3e6a353SBarry Smith PetscCall(MatAssembled(H, &assembled)); 448*b3e6a353SBarry Smith if (assembled) PetscCall(MatZeroEntries(H)); 449*b3e6a353SBarry Smith 450*b3e6a353SBarry Smith PetscCall(VecSet(user->s, zero)); 451*b3e6a353SBarry Smith 452*b3e6a353SBarry Smith /* Loop over matrix columns to compute entries of the Hessian */ 453*b3e6a353SBarry Smith for (j = 0; j < ndim; j++) { 454*b3e6a353SBarry Smith PetscCall(VecSetValues(user->s, 1, &j, &one, INSERT_VALUES)); 455*b3e6a353SBarry Smith PetscCall(VecAssemblyBegin(user->s)); 456*b3e6a353SBarry Smith PetscCall(VecAssemblyEnd(user->s)); 457*b3e6a353SBarry Smith 458*b3e6a353SBarry Smith PetscCall(HessianProduct(ptr, user->s, user->y)); 459*b3e6a353SBarry Smith 460*b3e6a353SBarry Smith PetscCall(VecSetValues(user->s, 1, &j, &zero, INSERT_VALUES)); 461*b3e6a353SBarry Smith PetscCall(VecAssemblyBegin(user->s)); 462*b3e6a353SBarry Smith PetscCall(VecAssemblyEnd(user->s)); 463*b3e6a353SBarry Smith 464*b3e6a353SBarry Smith PetscCall(VecGetArray(user->y, &y)); 465*b3e6a353SBarry Smith for (i = 0; i < ndim; i++) { 466*b3e6a353SBarry Smith if (y[i] != zero) PetscCall(MatSetValues(H, 1, &i, 1, &j, &y[i], ADD_VALUES)); 467*b3e6a353SBarry Smith } 468*b3e6a353SBarry Smith PetscCall(VecRestoreArray(user->y, &y)); 469*b3e6a353SBarry Smith } 470*b3e6a353SBarry Smith PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 471*b3e6a353SBarry Smith PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 472*b3e6a353SBarry Smith PetscFunctionReturn(0); 473*b3e6a353SBarry Smith } 474*b3e6a353SBarry Smith 475*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 476*b3e6a353SBarry Smith /* 477*b3e6a353SBarry Smith MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector 478*b3e6a353SBarry Smith products. 479*b3e6a353SBarry Smith 480*b3e6a353SBarry Smith Input Parameters: 481*b3e6a353SBarry Smith . tao - the Tao context 482*b3e6a353SBarry Smith . X - the input vector 483*b3e6a353SBarry Smith . ptr - optional user-defined context, as set by TaoSetHessian() 484*b3e6a353SBarry Smith 485*b3e6a353SBarry Smith Output Parameters: 486*b3e6a353SBarry Smith . H - Hessian matrix 487*b3e6a353SBarry Smith . PrecH - optionally different preconditioning Hessian 488*b3e6a353SBarry Smith . flag - flag indicating matrix structure 489*b3e6a353SBarry Smith */ 490*b3e6a353SBarry Smith PetscErrorCode MatrixFreeHessian(Tao tao, Vec X, Mat H, Mat PrecH, void *ptr) 491*b3e6a353SBarry Smith { 492*b3e6a353SBarry Smith AppCtx *user = (AppCtx *)ptr; 493*b3e6a353SBarry Smith 494*b3e6a353SBarry Smith /* Sets location of vector for use in computing matrix-vector products of the form H(X)*y */ 495*b3e6a353SBarry Smith PetscFunctionBeginUser; 496*b3e6a353SBarry Smith user->xvec = X; 497*b3e6a353SBarry Smith PetscFunctionReturn(0); 498*b3e6a353SBarry Smith } 499*b3e6a353SBarry Smith 500*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 501*b3e6a353SBarry Smith /* 502*b3e6a353SBarry Smith HessianProductMat - Computes the matrix-vector product 503*b3e6a353SBarry Smith y = mat*svec. 504*b3e6a353SBarry Smith 505*b3e6a353SBarry Smith Input Parameters: 506*b3e6a353SBarry Smith . mat - input matrix 507*b3e6a353SBarry Smith . svec - input vector 508*b3e6a353SBarry Smith 509*b3e6a353SBarry Smith Output Parameters: 510*b3e6a353SBarry Smith . y - solution vector 511*b3e6a353SBarry Smith */ 512*b3e6a353SBarry Smith PetscErrorCode HessianProductMat(Mat mat, Vec svec, Vec y) 513*b3e6a353SBarry Smith { 514*b3e6a353SBarry Smith void *ptr; 515*b3e6a353SBarry Smith 516*b3e6a353SBarry Smith PetscFunctionBeginUser; 517*b3e6a353SBarry Smith PetscCall(MatShellGetContext(mat, &ptr)); 518*b3e6a353SBarry Smith PetscCall(HessianProduct(ptr, svec, y)); 519*b3e6a353SBarry Smith PetscFunctionReturn(0); 520*b3e6a353SBarry Smith } 521*b3e6a353SBarry Smith 522*b3e6a353SBarry Smith /* ------------------------------------------------------------------- */ 523*b3e6a353SBarry Smith /* 524*b3e6a353SBarry Smith Hessian Product - Computes the matrix-vector product: 525*b3e6a353SBarry Smith y = f''(x)*svec. 526*b3e6a353SBarry Smith 527*b3e6a353SBarry Smith Input Parameters: 528*b3e6a353SBarry Smith . ptr - pointer to the user-defined context 529*b3e6a353SBarry Smith . svec - input vector 530*b3e6a353SBarry Smith 531*b3e6a353SBarry Smith Output Parameters: 532*b3e6a353SBarry Smith . y - product vector 533*b3e6a353SBarry Smith */ 534*b3e6a353SBarry Smith PetscErrorCode HessianProduct(void *ptr, Vec svec, Vec y) 535*b3e6a353SBarry Smith { 536*b3e6a353SBarry Smith AppCtx *user = (AppCtx *)ptr; 537*b3e6a353SBarry Smith PetscReal p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area; 538*b3e6a353SBarry Smith const PetscScalar *x, *s; 539*b3e6a353SBarry Smith PetscReal v, vb, vl, vr, vt, hxhx, hyhy; 540*b3e6a353SBarry Smith PetscInt nx, ny, i, j, k, ind; 541*b3e6a353SBarry Smith 542*b3e6a353SBarry Smith PetscFunctionBeginUser; 543*b3e6a353SBarry Smith nx = user->mx; 544*b3e6a353SBarry Smith ny = user->my; 545*b3e6a353SBarry Smith hx = user->hx; 546*b3e6a353SBarry Smith hy = user->hy; 547*b3e6a353SBarry Smith hxhx = one / (hx * hx); 548*b3e6a353SBarry Smith hyhy = one / (hy * hy); 549*b3e6a353SBarry Smith 550*b3e6a353SBarry Smith /* Get pointers to vector data */ 551*b3e6a353SBarry Smith PetscCall(VecGetArrayRead(user->xvec, &x)); 552*b3e6a353SBarry Smith PetscCall(VecGetArrayRead(svec, &s)); 553*b3e6a353SBarry Smith 554*b3e6a353SBarry Smith /* Initialize product vector to zero */ 555*b3e6a353SBarry Smith PetscCall(VecSet(y, zero)); 556*b3e6a353SBarry Smith 557*b3e6a353SBarry Smith /* Compute f''(x)*s over the lower triangular elements */ 558*b3e6a353SBarry Smith for (j = -1; j < ny; j++) { 559*b3e6a353SBarry Smith for (i = -1; i < nx; i++) { 560*b3e6a353SBarry Smith k = nx * j + i; 561*b3e6a353SBarry Smith v = zero; 562*b3e6a353SBarry Smith vr = zero; 563*b3e6a353SBarry Smith vt = zero; 564*b3e6a353SBarry Smith if (i != -1 && j != -1) v = s[k]; 565*b3e6a353SBarry Smith if (i != nx - 1 && j != -1) { 566*b3e6a353SBarry Smith vr = s[k + 1]; 567*b3e6a353SBarry Smith ind = k + 1; 568*b3e6a353SBarry Smith val = hxhx * (vr - v); 569*b3e6a353SBarry Smith PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 570*b3e6a353SBarry Smith } 571*b3e6a353SBarry Smith if (i != -1 && j != ny - 1) { 572*b3e6a353SBarry Smith vt = s[k + nx]; 573*b3e6a353SBarry Smith ind = k + nx; 574*b3e6a353SBarry Smith val = hyhy * (vt - v); 575*b3e6a353SBarry Smith PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 576*b3e6a353SBarry Smith } 577*b3e6a353SBarry Smith if (i != -1 && j != -1) { 578*b3e6a353SBarry Smith ind = k; 579*b3e6a353SBarry Smith val = hxhx * (v - vr) + hyhy * (v - vt); 580*b3e6a353SBarry Smith PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 581*b3e6a353SBarry Smith } 582*b3e6a353SBarry Smith } 583*b3e6a353SBarry Smith } 584*b3e6a353SBarry Smith 585*b3e6a353SBarry Smith /* Compute f''(x)*s over the upper triangular elements */ 586*b3e6a353SBarry Smith for (j = 0; j <= ny; j++) { 587*b3e6a353SBarry Smith for (i = 0; i <= nx; i++) { 588*b3e6a353SBarry Smith k = nx * j + i; 589*b3e6a353SBarry Smith v = zero; 590*b3e6a353SBarry Smith vl = zero; 591*b3e6a353SBarry Smith vb = zero; 592*b3e6a353SBarry Smith if (i != nx && j != ny) v = s[k]; 593*b3e6a353SBarry Smith if (i != nx && j != 0) { 594*b3e6a353SBarry Smith vb = s[k - nx]; 595*b3e6a353SBarry Smith ind = k - nx; 596*b3e6a353SBarry Smith val = hyhy * (vb - v); 597*b3e6a353SBarry Smith PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 598*b3e6a353SBarry Smith } 599*b3e6a353SBarry Smith if (i != 0 && j != ny) { 600*b3e6a353SBarry Smith vl = s[k - 1]; 601*b3e6a353SBarry Smith ind = k - 1; 602*b3e6a353SBarry Smith val = hxhx * (vl - v); 603*b3e6a353SBarry Smith PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 604*b3e6a353SBarry Smith } 605*b3e6a353SBarry Smith if (i != nx && j != ny) { 606*b3e6a353SBarry Smith ind = k; 607*b3e6a353SBarry Smith val = hxhx * (v - vl) + hyhy * (v - vb); 608*b3e6a353SBarry Smith PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 609*b3e6a353SBarry Smith } 610*b3e6a353SBarry Smith } 611*b3e6a353SBarry Smith } 612*b3e6a353SBarry Smith /* Restore vector data */ 613*b3e6a353SBarry Smith PetscCall(VecRestoreArrayRead(svec, &s)); 614*b3e6a353SBarry Smith PetscCall(VecRestoreArrayRead(user->xvec, &x)); 615*b3e6a353SBarry Smith 616*b3e6a353SBarry Smith /* Assemble vector */ 617*b3e6a353SBarry Smith PetscCall(VecAssemblyBegin(y)); 618*b3e6a353SBarry Smith PetscCall(VecAssemblyEnd(y)); 619*b3e6a353SBarry Smith 620*b3e6a353SBarry Smith /* Scale resulting vector by area */ 621*b3e6a353SBarry Smith area = p5 * hx * hy; 622*b3e6a353SBarry Smith PetscCall(VecScale(y, area)); 623*b3e6a353SBarry Smith PetscCall(PetscLogFlops(18.0 * nx * ny)); 624*b3e6a353SBarry Smith PetscFunctionReturn(0); 625*b3e6a353SBarry Smith } 626*b3e6a353SBarry Smith 627*b3e6a353SBarry Smith /*TEST 628*b3e6a353SBarry Smith 629*b3e6a353SBarry Smith build: 630*b3e6a353SBarry Smith requires: !complex 631*b3e6a353SBarry Smith 632*b3e6a353SBarry Smith test: 633*b3e6a353SBarry Smith args: -tao_monitor -tao_type bntr -tao_view -tao_bnk_ksp_monitor_short 634*b3e6a353SBarry Smith 635*b3e6a353SBarry Smith TEST*/ 636