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 31c4762a1bSJed Brown static char help[]= 32c4762a1bSJed Brown "Demonstrates use of the TAO package to solve \n\ 33c4762a1bSJed Brown unconstrained minimization problems on a single processor. This example \n\ 34c4762a1bSJed Brown is based on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 \n\ 35c4762a1bSJed Brown test suite.\n\ 36c4762a1bSJed Brown The command line options are:\n\ 37c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 38c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 39c4762a1bSJed Brown -par <param>, where <param> = angle of twist per unit length\n\n"; 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown User-defined application context - contains data needed by the 43c4762a1bSJed Brown application-provided call-back routines, FormFunction(), 44c4762a1bSJed Brown FormGradient(), and FormHessian(). 45c4762a1bSJed Brown */ 46c4762a1bSJed Brown 47c4762a1bSJed Brown typedef struct { 48c4762a1bSJed Brown PetscReal param; /* nonlinearity parameter */ 49c4762a1bSJed Brown PetscInt mx, my; /* discretization in x- and y-directions */ 50c4762a1bSJed Brown PetscInt ndim; /* problem dimension */ 51c4762a1bSJed Brown Vec s, y, xvec; /* work space for computing Hessian */ 52c4762a1bSJed Brown PetscReal hx, hy; /* mesh spacing in x- and y-directions */ 53c4762a1bSJed Brown } AppCtx; 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 56c4762a1bSJed Brown 57c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx*,Vec); 58c4762a1bSJed Brown PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*); 59c4762a1bSJed Brown PetscErrorCode FormGradient(Tao,Vec,Vec,void*); 60c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*); 61c4762a1bSJed Brown PetscErrorCode HessianProductMat(Mat,Vec,Vec); 62c4762a1bSJed Brown PetscErrorCode HessianProduct(void*,Vec,Vec); 63c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao,Vec,Mat,Mat,void*); 64c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *); 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscErrorCode main(int argc,char **argv) 67c4762a1bSJed Brown { 68c4762a1bSJed Brown PetscInt mx=10; /* discretization in x-direction */ 69c4762a1bSJed Brown PetscInt my=10; /* discretization in y-direction */ 70c4762a1bSJed Brown Vec x; /* solution, gradient vectors */ 71c4762a1bSJed Brown PetscBool flg; /* A return value when checking for use options */ 72c4762a1bSJed Brown Tao tao; /* Tao solver context */ 73c4762a1bSJed Brown Mat H; /* Hessian matrix */ 74c4762a1bSJed Brown AppCtx user; /* application context */ 75c4762a1bSJed Brown PetscMPIInt size; /* number of processes */ 76c4762a1bSJed Brown PetscReal one=1.0; 77c4762a1bSJed Brown 78c4762a1bSJed Brown PetscBool test_lmvm = PETSC_FALSE; 79c4762a1bSJed Brown KSP ksp; 80c4762a1bSJed Brown PC pc; 81c4762a1bSJed Brown Mat M; 82c4762a1bSJed Brown Vec in, out, out2; 83c4762a1bSJed Brown PetscReal mult_solve_dist; 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Initialize TAO,PETSc */ 869566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char *)0,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")); 98*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n",mx,my)); 99c4762a1bSJed Brown user.ndim = mx * my; user.mx = mx; user.my = my; 100c4762a1bSJed Brown user.hx = one/(mx+1); user.hy = one/(my+1); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Allocate vectors */ 1039566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.ndim,&user.y)); 1049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.y,&user.s)); 1059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.y,&x)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1089566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF,&tao)); 1099566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOLMVM)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1129566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user,x)); 1139566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 1169566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* From command line options, determine if using matrix-free hessian */ 1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-my_tao_mf",&flg)); 120c4762a1bSJed Brown if (flg) { 1219566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_SELF,user.ndim,user.ndim,user.ndim,user.ndim,(void*)&user,&H)); 1229566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(H,MATOP_MULT,(void(*)(void))HessianProductMat)); 1239566063dSJacob Faibussowitsch PetscCall(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE)); 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,H,H,MatrixFreeHessian,(void *)&user)); 126c4762a1bSJed Brown } else { 1279566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,user.ndim,user.ndim,5,NULL,&H)); 1289566063dSJacob Faibussowitsch PetscCall(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE)); 1299566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,H,H,FormHessian,(void *)&user)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Test the LMVM matrix */ 133c4762a1bSJed Brown if (test_lmvm) { 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bntr")); 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-tao_bnk_pc_type", "lmvm")); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Check for any TAO command line options */ 1399566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1429566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Test the LMVM matrix */ 145c4762a1bSJed Brown if (test_lmvm) { 1469566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 1479566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1489566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(pc, &M)); 1499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &in)); 1509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out)); 1519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out2)); 1529566063dSJacob Faibussowitsch PetscCall(VecSet(in, 5.0)); 1539566063dSJacob Faibussowitsch PetscCall(MatMult(M, in, out)); 1549566063dSJacob Faibussowitsch PetscCall(MatSolve(M, out, out2)); 1559566063dSJacob Faibussowitsch PetscCall(VecAXPY(out2, -1.0, in)); 1569566063dSJacob Faibussowitsch PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist)); 157*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between MatMult and MatSolve: %e\n", (double)mult_solve_dist)); 1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&in)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out2)); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.s)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.y)); 1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 168c4762a1bSJed Brown 1699566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 170b122ec5aSJacob Faibussowitsch return 0; 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 174c4762a1bSJed Brown /* 175c4762a1bSJed Brown FormInitialGuess - Computes an initial approximation to the solution. 176c4762a1bSJed Brown 177c4762a1bSJed Brown Input Parameters: 178c4762a1bSJed Brown . user - user-defined application context 179c4762a1bSJed Brown . X - vector 180c4762a1bSJed Brown 181c4762a1bSJed Brown Output Parameters: 182c4762a1bSJed Brown . X - vector 183c4762a1bSJed Brown */ 184c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy, temp; 187c4762a1bSJed Brown PetscReal val; 188c4762a1bSJed Brown PetscInt i, j, k, nx = user->mx, ny = user->my; 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* Compute initial guess */ 191c4762a1bSJed Brown PetscFunctionBeginUser; 192c4762a1bSJed Brown for (j=0; j<ny; j++) { 193c4762a1bSJed Brown temp = PetscMin(j+1,ny-j)*hy; 194c4762a1bSJed Brown for (i=0; i<nx; i++) { 195c4762a1bSJed Brown k = nx*j + i; 196c4762a1bSJed Brown val = PetscMin((PetscMin(i+1,nx-i))*hx,temp); 1979566063dSJacob Faibussowitsch PetscCall(VecSetValues(X,1,&k,&val,ADD_VALUES)); 198c4762a1bSJed Brown } 199c4762a1bSJed Brown } 2009566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 2019566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 202c4762a1bSJed Brown PetscFunctionReturn(0); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 208c4762a1bSJed Brown 209c4762a1bSJed Brown Input Parameters: 210c4762a1bSJed Brown tao - the Tao context 211c4762a1bSJed Brown X - the input vector 212c4762a1bSJed Brown ptr - optional user-defined context, as set by TaoSetFunction() 213c4762a1bSJed Brown 214c4762a1bSJed Brown Output Parameters: 215c4762a1bSJed Brown f - the newly evaluated function 216c4762a1bSJed Brown G - the newly evaluated gradient 217c4762a1bSJed Brown */ 218c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr) 219c4762a1bSJed Brown { 220c4762a1bSJed Brown PetscFunctionBeginUser; 2219566063dSJacob Faibussowitsch PetscCall(FormFunction(tao,X,f,ptr)); 2229566063dSJacob Faibussowitsch PetscCall(FormGradient(tao,X,G,ptr)); 223c4762a1bSJed Brown PetscFunctionReturn(0); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 227c4762a1bSJed Brown /* 228c4762a1bSJed Brown FormFunction - Evaluates the function, f(X). 229c4762a1bSJed Brown 230c4762a1bSJed Brown Input Parameters: 231c4762a1bSJed Brown . tao - the Tao context 232c4762a1bSJed Brown . X - the input vector 233c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunction() 234c4762a1bSJed Brown 235c4762a1bSJed Brown Output Parameters: 236c4762a1bSJed Brown . f - the newly evaluated function 237c4762a1bSJed Brown */ 238c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr) 239c4762a1bSJed Brown { 240c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 241c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5; 242c4762a1bSJed Brown PetscReal zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0; 243c4762a1bSJed Brown PetscReal v, cdiv3 = user->param/three; 244c4762a1bSJed Brown const PetscScalar *x; 245c4762a1bSJed Brown PetscInt nx = user->mx, ny = user->my, i, j, k; 246c4762a1bSJed Brown 247c4762a1bSJed Brown PetscFunctionBeginUser; 248c4762a1bSJed Brown /* Get pointer to vector data */ 2499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* Compute function contributions over the lower triangular elements */ 252c4762a1bSJed Brown for (j=-1; j<ny; j++) { 253c4762a1bSJed Brown for (i=-1; i<nx; i++) { 254c4762a1bSJed Brown k = nx*j + i; 255c4762a1bSJed Brown v = zero; 256c4762a1bSJed Brown vr = zero; 257c4762a1bSJed Brown vt = zero; 258c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k]; 259c4762a1bSJed Brown if (i < nx-1 && j > -1) vr = x[k+1]; 260c4762a1bSJed Brown if (i > -1 && j < ny-1) vt = x[k+nx]; 261c4762a1bSJed Brown dvdx = (vr-v)/hx; 262c4762a1bSJed Brown dvdy = (vt-v)/hy; 263c4762a1bSJed Brown fquad += dvdx*dvdx + dvdy*dvdy; 264c4762a1bSJed Brown flin -= cdiv3*(v+vr+vt); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* Compute function contributions over the upper triangular elements */ 269c4762a1bSJed Brown for (j=0; j<=ny; j++) { 270c4762a1bSJed Brown for (i=0; i<=nx; i++) { 271c4762a1bSJed Brown k = nx*j + i; 272c4762a1bSJed Brown vb = zero; 273c4762a1bSJed Brown vl = zero; 274c4762a1bSJed Brown v = zero; 275c4762a1bSJed Brown if (i < nx && j > 0) vb = x[k-nx]; 276c4762a1bSJed Brown if (i > 0 && j < ny) vl = x[k-1]; 277c4762a1bSJed Brown if (i < nx && j < ny) v = x[k]; 278c4762a1bSJed Brown dvdx = (v-vl)/hx; 279c4762a1bSJed Brown dvdy = (v-vb)/hy; 280c4762a1bSJed Brown fquad = fquad + dvdx*dvdx + dvdy*dvdy; 281c4762a1bSJed Brown flin = flin - cdiv3*(vb+vl+v); 282c4762a1bSJed Brown } 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* Restore vector */ 2869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* Scale the function */ 289c4762a1bSJed Brown area = p5*hx*hy; 290c4762a1bSJed Brown *f = area*(p5*fquad+flin); 291c4762a1bSJed Brown 2929566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(24.0*nx*ny)); 293c4762a1bSJed Brown PetscFunctionReturn(0); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 297c4762a1bSJed Brown /* 298c4762a1bSJed Brown FormGradient - Evaluates the gradient, G(X). 299c4762a1bSJed Brown 300c4762a1bSJed Brown Input Parameters: 301c4762a1bSJed Brown . tao - the Tao context 302c4762a1bSJed Brown . X - input vector 303c4762a1bSJed Brown . ptr - optional user-defined context 304c4762a1bSJed Brown 305c4762a1bSJed Brown Output Parameters: 306c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 307c4762a1bSJed Brown */ 308c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr) 309c4762a1bSJed Brown { 310c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 311c4762a1bSJed Brown PetscReal zero=0.0, p5=0.5, three = 3.0, area, val; 312c4762a1bSJed Brown PetscInt nx = user->mx, ny = user->my, ind, i, j, k; 313c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy; 314c4762a1bSJed Brown PetscReal vb, vl, vr, vt, dvdx, dvdy; 315c4762a1bSJed Brown PetscReal v, cdiv3 = user->param/three; 316c4762a1bSJed Brown const PetscScalar *x; 317c4762a1bSJed Brown 318c4762a1bSJed Brown PetscFunctionBeginUser; 319c4762a1bSJed Brown /* Initialize gradient to zero */ 3209566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 321c4762a1bSJed Brown 322c4762a1bSJed Brown /* Get pointer to vector data */ 3239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* Compute gradient contributions over the lower triangular elements */ 326c4762a1bSJed Brown for (j=-1; j<ny; j++) { 327c4762a1bSJed Brown for (i=-1; i<nx; i++) { 328c4762a1bSJed Brown k = nx*j + i; 329c4762a1bSJed Brown v = zero; 330c4762a1bSJed Brown vr = zero; 331c4762a1bSJed Brown vt = zero; 332c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k]; 333c4762a1bSJed Brown if (i < nx-1 && j > -1) vr = x[k+1]; 334c4762a1bSJed Brown if (i > -1 && j < ny-1) vt = x[k+nx]; 335c4762a1bSJed Brown dvdx = (vr-v)/hx; 336c4762a1bSJed Brown dvdy = (vt-v)/hy; 337c4762a1bSJed Brown if (i != -1 && j != -1) { 338c4762a1bSJed Brown ind = k; val = - dvdx/hx - dvdy/hy - cdiv3; 3399566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 340c4762a1bSJed Brown } 341c4762a1bSJed Brown if (i != nx-1 && j != -1) { 342c4762a1bSJed Brown ind = k+1; val = dvdx/hx - cdiv3; 3439566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown if (i != -1 && j != ny-1) { 346c4762a1bSJed Brown ind = k+nx; val = dvdy/hy - cdiv3; 3479566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 350c4762a1bSJed Brown } 351c4762a1bSJed Brown 352c4762a1bSJed Brown /* Compute gradient contributions over the upper triangular elements */ 353c4762a1bSJed Brown for (j=0; j<=ny; j++) { 354c4762a1bSJed Brown for (i=0; i<=nx; i++) { 355c4762a1bSJed Brown k = nx*j + i; 356c4762a1bSJed Brown vb = zero; 357c4762a1bSJed Brown vl = zero; 358c4762a1bSJed Brown v = zero; 359c4762a1bSJed Brown if (i < nx && j > 0) vb = x[k-nx]; 360c4762a1bSJed Brown if (i > 0 && j < ny) vl = x[k-1]; 361c4762a1bSJed Brown if (i < nx && j < ny) v = x[k]; 362c4762a1bSJed Brown dvdx = (v-vl)/hx; 363c4762a1bSJed Brown dvdy = (v-vb)/hy; 364c4762a1bSJed Brown if (i != nx && j != 0) { 365c4762a1bSJed Brown ind = k-nx; val = - dvdy/hy - cdiv3; 3669566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 367c4762a1bSJed Brown } 368c4762a1bSJed Brown if (i != 0 && j != ny) { 369c4762a1bSJed Brown ind = k-1; val = - dvdx/hx - cdiv3; 3709566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 371c4762a1bSJed Brown } 372c4762a1bSJed Brown if (i != nx && j != ny) { 373c4762a1bSJed Brown ind = k; val = dvdx/hx + dvdy/hy - cdiv3; 3749566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 375c4762a1bSJed Brown } 376c4762a1bSJed Brown } 377c4762a1bSJed Brown } 3789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 379c4762a1bSJed Brown 380c4762a1bSJed Brown /* Assemble gradient vector */ 3819566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(G)); 3829566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(G)); 383c4762a1bSJed Brown 384c4762a1bSJed Brown /* Scale the gradient */ 385c4762a1bSJed Brown area = p5*hx*hy; 3869566063dSJacob Faibussowitsch PetscCall(VecScale(G, area)); 3879566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(24.0*nx*ny)); 388c4762a1bSJed Brown PetscFunctionReturn(0); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 391c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 392c4762a1bSJed Brown /* 393c4762a1bSJed Brown FormHessian - Forms the Hessian matrix. 394c4762a1bSJed Brown 395c4762a1bSJed Brown Input Parameters: 396c4762a1bSJed Brown . tao - the Tao context 397c4762a1bSJed Brown . X - the input vector 398c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 399c4762a1bSJed Brown 400c4762a1bSJed Brown Output Parameters: 401c4762a1bSJed Brown . H - Hessian matrix 402c4762a1bSJed Brown . PrecH - optionally different preconditioning Hessian 403c4762a1bSJed Brown . flag - flag indicating matrix structure 404c4762a1bSJed Brown 405c4762a1bSJed Brown Notes: 406c4762a1bSJed Brown This routine is intended simply as an example of the interface 407c4762a1bSJed Brown to a Hessian evaluation routine. Since this example compute the 408c4762a1bSJed Brown Hessian a column at a time, it is not particularly efficient and 409c4762a1bSJed Brown is not recommended. 410c4762a1bSJed Brown */ 411c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H,Mat Hpre, void *ptr) 412c4762a1bSJed Brown { 413c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 414c4762a1bSJed Brown PetscInt i,j, ndim = user->ndim; 415c4762a1bSJed Brown PetscReal *y, zero = 0.0, one = 1.0; 416c4762a1bSJed Brown PetscBool assembled; 417c4762a1bSJed Brown 418c4762a1bSJed Brown PetscFunctionBeginUser; 419c4762a1bSJed Brown user->xvec = X; 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* Initialize Hessian entries and work vector to zero */ 4229566063dSJacob Faibussowitsch PetscCall(MatAssembled(H,&assembled)); 4239566063dSJacob Faibussowitsch if (assembled)PetscCall(MatZeroEntries(H)); 424c4762a1bSJed Brown 4259566063dSJacob Faibussowitsch PetscCall(VecSet(user->s, zero)); 426c4762a1bSJed Brown 427c4762a1bSJed Brown /* Loop over matrix columns to compute entries of the Hessian */ 428c4762a1bSJed Brown for (j=0; j<ndim; j++) { 4299566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->s,1,&j,&one,INSERT_VALUES)); 4309566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->s)); 4319566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->s)); 432c4762a1bSJed Brown 4339566063dSJacob Faibussowitsch PetscCall(HessianProduct(ptr,user->s,user->y)); 434c4762a1bSJed Brown 4359566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->s,1,&j,&zero,INSERT_VALUES)); 4369566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->s)); 4379566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->s)); 438c4762a1bSJed Brown 4399566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->y,&y)); 440c4762a1bSJed Brown for (i=0; i<ndim; i++) { 441c4762a1bSJed Brown if (y[i] != zero) { 4429566063dSJacob Faibussowitsch PetscCall(MatSetValues(H,1,&i,1,&j,&y[i],ADD_VALUES)); 443c4762a1bSJed Brown } 444c4762a1bSJed Brown } 4459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->y,&y)); 446c4762a1bSJed Brown } 4479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 4489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 449c4762a1bSJed Brown PetscFunctionReturn(0); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 453c4762a1bSJed Brown /* 454c4762a1bSJed Brown MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector 455c4762a1bSJed Brown products. 456c4762a1bSJed Brown 457c4762a1bSJed Brown Input Parameters: 458c4762a1bSJed Brown . tao - the Tao context 459c4762a1bSJed Brown . X - the input vector 460c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 461c4762a1bSJed Brown 462c4762a1bSJed Brown Output Parameters: 463c4762a1bSJed Brown . H - Hessian matrix 464c4762a1bSJed Brown . PrecH - optionally different preconditioning Hessian 465c4762a1bSJed Brown . flag - flag indicating matrix structure 466c4762a1bSJed Brown */ 467c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao tao,Vec X,Mat H,Mat PrecH, void *ptr) 468c4762a1bSJed Brown { 469c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 470c4762a1bSJed Brown 471c4762a1bSJed Brown /* Sets location of vector for use in computing matrix-vector products of the form H(X)*y */ 472362febeeSStefano Zampini PetscFunctionBeginUser; 473c4762a1bSJed Brown user->xvec = X; 474c4762a1bSJed Brown PetscFunctionReturn(0); 475c4762a1bSJed Brown } 476c4762a1bSJed Brown 477c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 478c4762a1bSJed Brown /* 479c4762a1bSJed Brown HessianProductMat - Computes the matrix-vector product 480c4762a1bSJed Brown y = mat*svec. 481c4762a1bSJed Brown 482c4762a1bSJed Brown Input Parameters: 483c4762a1bSJed Brown . mat - input matrix 484c4762a1bSJed Brown . svec - input vector 485c4762a1bSJed Brown 486c4762a1bSJed Brown Output Parameters: 487c4762a1bSJed Brown . y - solution vector 488c4762a1bSJed Brown */ 489c4762a1bSJed Brown PetscErrorCode HessianProductMat(Mat mat,Vec svec,Vec y) 490c4762a1bSJed Brown { 491c4762a1bSJed Brown void *ptr; 492c4762a1bSJed Brown 493c4762a1bSJed Brown PetscFunctionBeginUser; 4949566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ptr)); 4959566063dSJacob Faibussowitsch PetscCall(HessianProduct(ptr,svec,y)); 496c4762a1bSJed Brown PetscFunctionReturn(0); 497c4762a1bSJed Brown } 498c4762a1bSJed Brown 499c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 500c4762a1bSJed Brown /* 501c4762a1bSJed Brown Hessian Product - Computes the matrix-vector product: 502c4762a1bSJed Brown y = f''(x)*svec. 503c4762a1bSJed Brown 5047a7aea1fSJed Brown Input Parameters: 505c4762a1bSJed Brown . ptr - pointer to the user-defined context 506c4762a1bSJed Brown . svec - input vector 507c4762a1bSJed Brown 508c4762a1bSJed Brown Output Parameters: 509c4762a1bSJed Brown . y - product vector 510c4762a1bSJed Brown */ 511c4762a1bSJed Brown PetscErrorCode HessianProduct(void *ptr,Vec svec,Vec y) 512c4762a1bSJed Brown { 513c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 514c4762a1bSJed Brown PetscReal p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area; 515c4762a1bSJed Brown const PetscScalar *x, *s; 516c4762a1bSJed Brown PetscReal v, vb, vl, vr, vt, hxhx, hyhy; 517c4762a1bSJed Brown PetscInt nx, ny, i, j, k, ind; 518c4762a1bSJed Brown 519c4762a1bSJed Brown PetscFunctionBeginUser; 520c4762a1bSJed Brown nx = user->mx; 521c4762a1bSJed Brown ny = user->my; 522c4762a1bSJed Brown hx = user->hx; 523c4762a1bSJed Brown hy = user->hy; 524c4762a1bSJed Brown hxhx = one/(hx*hx); 525c4762a1bSJed Brown hyhy = one/(hy*hy); 526c4762a1bSJed Brown 527c4762a1bSJed Brown /* Get pointers to vector data */ 5289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user->xvec,&x)); 5299566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(svec,&s)); 530c4762a1bSJed Brown 531c4762a1bSJed Brown /* Initialize product vector to zero */ 5329566063dSJacob Faibussowitsch PetscCall(VecSet(y, zero)); 533c4762a1bSJed Brown 534c4762a1bSJed Brown /* Compute f''(x)*s over the lower triangular elements */ 535c4762a1bSJed Brown for (j=-1; j<ny; j++) { 536c4762a1bSJed Brown for (i=-1; i<nx; i++) { 537c4762a1bSJed Brown k = nx*j + i; 538c4762a1bSJed Brown v = zero; 539c4762a1bSJed Brown vr = zero; 540c4762a1bSJed Brown vt = zero; 541c4762a1bSJed Brown if (i != -1 && j != -1) v = s[k]; 542c4762a1bSJed Brown if (i != nx-1 && j != -1) { 543c4762a1bSJed Brown vr = s[k+1]; 544c4762a1bSJed Brown ind = k+1; val = hxhx*(vr-v); 5459566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 546c4762a1bSJed Brown } 547c4762a1bSJed Brown if (i != -1 && j != ny-1) { 548c4762a1bSJed Brown vt = s[k+nx]; 549c4762a1bSJed Brown ind = k+nx; val = hyhy*(vt-v); 5509566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 551c4762a1bSJed Brown } 552c4762a1bSJed Brown if (i != -1 && j != -1) { 553c4762a1bSJed Brown ind = k; val = hxhx*(v-vr) + hyhy*(v-vt); 5549566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 555c4762a1bSJed Brown } 556c4762a1bSJed Brown } 557c4762a1bSJed Brown } 558c4762a1bSJed Brown 559c4762a1bSJed Brown /* Compute f''(x)*s over the upper triangular elements */ 560c4762a1bSJed Brown for (j=0; j<=ny; j++) { 561c4762a1bSJed Brown for (i=0; i<=nx; i++) { 562c4762a1bSJed Brown k = nx*j + i; 563c4762a1bSJed Brown v = zero; 564c4762a1bSJed Brown vl = zero; 565c4762a1bSJed Brown vb = zero; 566c4762a1bSJed Brown if (i != nx && j != ny) v = s[k]; 567c4762a1bSJed Brown if (i != nx && j != 0) { 568c4762a1bSJed Brown vb = s[k-nx]; 569c4762a1bSJed Brown ind = k-nx; val = hyhy*(vb-v); 5709566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 571c4762a1bSJed Brown } 572c4762a1bSJed Brown if (i != 0 && j != ny) { 573c4762a1bSJed Brown vl = s[k-1]; 574c4762a1bSJed Brown ind = k-1; val = hxhx*(vl-v); 5759566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 576c4762a1bSJed Brown } 577c4762a1bSJed Brown if (i != nx && j != ny) { 578c4762a1bSJed Brown ind = k; val = hxhx*(v-vl) + hyhy*(v-vb); 5799566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 580c4762a1bSJed Brown } 581c4762a1bSJed Brown } 582c4762a1bSJed Brown } 583c4762a1bSJed Brown /* Restore vector data */ 5849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(svec,&s)); 5859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user->xvec,&x)); 586c4762a1bSJed Brown 587c4762a1bSJed Brown /* Assemble vector */ 5889566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(y)); 5899566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(y)); 590c4762a1bSJed Brown 591c4762a1bSJed Brown /* Scale resulting vector by area */ 592c4762a1bSJed Brown area = p5*hx*hy; 5939566063dSJacob Faibussowitsch PetscCall(VecScale(y, area)); 5949566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(18.0*nx*ny)); 595c4762a1bSJed Brown PetscFunctionReturn(0); 596c4762a1bSJed Brown } 597c4762a1bSJed Brown 598c4762a1bSJed Brown /*TEST 599c4762a1bSJed Brown 600c4762a1bSJed Brown build: 601c4762a1bSJed Brown requires: !complex 602c4762a1bSJed Brown 603c4762a1bSJed Brown test: 604c4762a1bSJed Brown suffix: 1 605c4762a1bSJed Brown args: -tao_smonitor -tao_type ntl -tao_gatol 1.e-4 606c4762a1bSJed Brown 607c4762a1bSJed Brown test: 608c4762a1bSJed Brown suffix: 2 609c4762a1bSJed Brown args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4 610c4762a1bSJed Brown 611c4762a1bSJed Brown test: 612c4762a1bSJed Brown suffix: 3 613c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4 -my_tao_mf -tao_test_hessian 614c4762a1bSJed Brown 615c4762a1bSJed Brown test: 616c4762a1bSJed Brown suffix: 4 617c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnls 618c4762a1bSJed Brown 619c4762a1bSJed Brown test: 620c4762a1bSJed Brown suffix: 5 621c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-3 -tao_type blmvm 622c4762a1bSJed Brown 623c4762a1bSJed Brown test: 624c4762a1bSJed Brown suffix: 6 625c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 626c4762a1bSJed Brown 627c4762a1bSJed Brown TEST*/ 628