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 */ 86*327415f7SBarry Smith PetscFunctionBeginUser; 879566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char *)0,help)); 889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size)); 893c859ba3SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Specify default parameters for the problem, check for command-line overrides */ 92c4762a1bSJed Brown user.param = 5.0; 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&my,&flg)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&mx,&flg)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,&flg)); 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg)); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n---- Elastic-Plastic Torsion Problem -----\n")); 9963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n",mx,my)); 100c4762a1bSJed Brown user.ndim = mx * my; user.mx = mx; user.my = my; 101c4762a1bSJed Brown user.hx = one/(mx+1); user.hy = one/(my+1); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Allocate vectors */ 1049566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.ndim,&user.y)); 1059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.y,&user.s)); 1069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.y,&x)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1099566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF,&tao)); 1109566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOLMVM)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1139566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user,x)); 1149566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 1179566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* From command line options, determine if using matrix-free hessian */ 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-my_tao_mf",&flg)); 121c4762a1bSJed Brown if (flg) { 1229566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_SELF,user.ndim,user.ndim,user.ndim,user.ndim,(void*)&user,&H)); 1239566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(H,MATOP_MULT,(void(*)(void))HessianProductMat)); 1249566063dSJacob Faibussowitsch PetscCall(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE)); 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,H,H,MatrixFreeHessian,(void *)&user)); 127c4762a1bSJed Brown } else { 1289566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,user.ndim,user.ndim,5,NULL,&H)); 1299566063dSJacob Faibussowitsch PetscCall(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE)); 1309566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,H,H,FormHessian,(void *)&user)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Test the LMVM matrix */ 134c4762a1bSJed Brown if (test_lmvm) { 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bntr")); 1369566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-tao_bnk_pc_type", "lmvm")); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Check for any TAO command line options */ 1409566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1439566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Test the LMVM matrix */ 146c4762a1bSJed Brown if (test_lmvm) { 1479566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 1489566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1499566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(pc, &M)); 1509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &in)); 1519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out)); 1529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &out2)); 1539566063dSJacob Faibussowitsch PetscCall(VecSet(in, 5.0)); 1549566063dSJacob Faibussowitsch PetscCall(MatMult(M, in, out)); 1559566063dSJacob Faibussowitsch PetscCall(MatSolve(M, out, out2)); 1569566063dSJacob Faibussowitsch PetscCall(VecAXPY(out2, -1.0, in)); 1579566063dSJacob Faibussowitsch PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist)); 15863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between MatMult and MatSolve: %e\n", (double)mult_solve_dist)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&in)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&out2)); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.s)); 1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.y)); 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 171b122ec5aSJacob Faibussowitsch return 0; 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 175c4762a1bSJed Brown /* 176c4762a1bSJed Brown FormInitialGuess - Computes an initial approximation to the solution. 177c4762a1bSJed Brown 178c4762a1bSJed Brown Input Parameters: 179c4762a1bSJed Brown . user - user-defined application context 180c4762a1bSJed Brown . X - vector 181c4762a1bSJed Brown 182c4762a1bSJed Brown Output Parameters: 183c4762a1bSJed Brown . X - vector 184c4762a1bSJed Brown */ 185c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 186c4762a1bSJed Brown { 187c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy, temp; 188c4762a1bSJed Brown PetscReal val; 189c4762a1bSJed Brown PetscInt i, j, k, nx = user->mx, ny = user->my; 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* Compute initial guess */ 192c4762a1bSJed Brown PetscFunctionBeginUser; 193c4762a1bSJed Brown for (j=0; j<ny; j++) { 194c4762a1bSJed Brown temp = PetscMin(j+1,ny-j)*hy; 195c4762a1bSJed Brown for (i=0; i<nx; i++) { 196c4762a1bSJed Brown k = nx*j + i; 197c4762a1bSJed Brown val = PetscMin((PetscMin(i+1,nx-i))*hx,temp); 1989566063dSJacob Faibussowitsch PetscCall(VecSetValues(X,1,&k,&val,ADD_VALUES)); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 2029566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 203c4762a1bSJed Brown PetscFunctionReturn(0); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 207c4762a1bSJed Brown /* 208c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 209c4762a1bSJed Brown 210c4762a1bSJed Brown Input Parameters: 211c4762a1bSJed Brown tao - the Tao context 212c4762a1bSJed Brown X - the input vector 213c4762a1bSJed Brown ptr - optional user-defined context, as set by TaoSetFunction() 214c4762a1bSJed Brown 215c4762a1bSJed Brown Output Parameters: 216c4762a1bSJed Brown f - the newly evaluated function 217c4762a1bSJed Brown G - the newly evaluated gradient 218c4762a1bSJed Brown */ 219c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown PetscFunctionBeginUser; 2229566063dSJacob Faibussowitsch PetscCall(FormFunction(tao,X,f,ptr)); 2239566063dSJacob Faibussowitsch PetscCall(FormGradient(tao,X,G,ptr)); 224c4762a1bSJed Brown PetscFunctionReturn(0); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 228c4762a1bSJed Brown /* 229c4762a1bSJed Brown FormFunction - Evaluates the function, f(X). 230c4762a1bSJed Brown 231c4762a1bSJed Brown Input Parameters: 232c4762a1bSJed Brown . tao - the Tao context 233c4762a1bSJed Brown . X - the input vector 234c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunction() 235c4762a1bSJed Brown 236c4762a1bSJed Brown Output Parameters: 237c4762a1bSJed Brown . f - the newly evaluated function 238c4762a1bSJed Brown */ 239c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr) 240c4762a1bSJed Brown { 241c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 242c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5; 243c4762a1bSJed Brown PetscReal zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0; 244c4762a1bSJed Brown PetscReal v, cdiv3 = user->param/three; 245c4762a1bSJed Brown const PetscScalar *x; 246c4762a1bSJed Brown PetscInt nx = user->mx, ny = user->my, i, j, k; 247c4762a1bSJed Brown 248c4762a1bSJed Brown PetscFunctionBeginUser; 249c4762a1bSJed Brown /* Get pointer to vector data */ 2509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* Compute function contributions over the lower triangular elements */ 253c4762a1bSJed Brown for (j=-1; j<ny; j++) { 254c4762a1bSJed Brown for (i=-1; i<nx; i++) { 255c4762a1bSJed Brown k = nx*j + i; 256c4762a1bSJed Brown v = zero; 257c4762a1bSJed Brown vr = zero; 258c4762a1bSJed Brown vt = zero; 259c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k]; 260c4762a1bSJed Brown if (i < nx-1 && j > -1) vr = x[k+1]; 261c4762a1bSJed Brown if (i > -1 && j < ny-1) vt = x[k+nx]; 262c4762a1bSJed Brown dvdx = (vr-v)/hx; 263c4762a1bSJed Brown dvdy = (vt-v)/hy; 264c4762a1bSJed Brown fquad += dvdx*dvdx + dvdy*dvdy; 265c4762a1bSJed Brown flin -= cdiv3*(v+vr+vt); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* Compute function contributions over the upper triangular elements */ 270c4762a1bSJed Brown for (j=0; j<=ny; j++) { 271c4762a1bSJed Brown for (i=0; i<=nx; i++) { 272c4762a1bSJed Brown k = nx*j + i; 273c4762a1bSJed Brown vb = zero; 274c4762a1bSJed Brown vl = zero; 275c4762a1bSJed Brown v = zero; 276c4762a1bSJed Brown if (i < nx && j > 0) vb = x[k-nx]; 277c4762a1bSJed Brown if (i > 0 && j < ny) vl = x[k-1]; 278c4762a1bSJed Brown if (i < nx && j < ny) v = x[k]; 279c4762a1bSJed Brown dvdx = (v-vl)/hx; 280c4762a1bSJed Brown dvdy = (v-vb)/hy; 281c4762a1bSJed Brown fquad = fquad + dvdx*dvdx + dvdy*dvdy; 282c4762a1bSJed Brown flin = flin - cdiv3*(vb+vl+v); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* Restore vector */ 2879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* Scale the function */ 290c4762a1bSJed Brown area = p5*hx*hy; 291c4762a1bSJed Brown *f = area*(p5*fquad+flin); 292c4762a1bSJed Brown 2939566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(24.0*nx*ny)); 294c4762a1bSJed Brown PetscFunctionReturn(0); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 298c4762a1bSJed Brown /* 299c4762a1bSJed Brown FormGradient - Evaluates the gradient, G(X). 300c4762a1bSJed Brown 301c4762a1bSJed Brown Input Parameters: 302c4762a1bSJed Brown . tao - the Tao context 303c4762a1bSJed Brown . X - input vector 304c4762a1bSJed Brown . ptr - optional user-defined context 305c4762a1bSJed Brown 306c4762a1bSJed Brown Output Parameters: 307c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 308c4762a1bSJed Brown */ 309c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr) 310c4762a1bSJed Brown { 311c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 312c4762a1bSJed Brown PetscReal zero=0.0, p5=0.5, three = 3.0, area, val; 313c4762a1bSJed Brown PetscInt nx = user->mx, ny = user->my, ind, i, j, k; 314c4762a1bSJed Brown PetscReal hx = user->hx, hy = user->hy; 315c4762a1bSJed Brown PetscReal vb, vl, vr, vt, dvdx, dvdy; 316c4762a1bSJed Brown PetscReal v, cdiv3 = user->param/three; 317c4762a1bSJed Brown const PetscScalar *x; 318c4762a1bSJed Brown 319c4762a1bSJed Brown PetscFunctionBeginUser; 320c4762a1bSJed Brown /* Initialize gradient to zero */ 3219566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* Get pointer to vector data */ 3249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* Compute gradient contributions over the lower triangular elements */ 327c4762a1bSJed Brown for (j=-1; j<ny; j++) { 328c4762a1bSJed Brown for (i=-1; i<nx; i++) { 329c4762a1bSJed Brown k = nx*j + i; 330c4762a1bSJed Brown v = zero; 331c4762a1bSJed Brown vr = zero; 332c4762a1bSJed Brown vt = zero; 333c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k]; 334c4762a1bSJed Brown if (i < nx-1 && j > -1) vr = x[k+1]; 335c4762a1bSJed Brown if (i > -1 && j < ny-1) vt = x[k+nx]; 336c4762a1bSJed Brown dvdx = (vr-v)/hx; 337c4762a1bSJed Brown dvdy = (vt-v)/hy; 338c4762a1bSJed Brown if (i != -1 && j != -1) { 339c4762a1bSJed Brown ind = k; val = - dvdx/hx - dvdy/hy - cdiv3; 3409566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown if (i != nx-1 && j != -1) { 343c4762a1bSJed Brown ind = k+1; val = dvdx/hx - cdiv3; 3449566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown if (i != -1 && j != ny-1) { 347c4762a1bSJed Brown ind = k+nx; val = dvdy/hy - cdiv3; 3489566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown } 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* Compute gradient contributions over the upper triangular elements */ 354c4762a1bSJed Brown for (j=0; j<=ny; j++) { 355c4762a1bSJed Brown for (i=0; i<=nx; i++) { 356c4762a1bSJed Brown k = nx*j + i; 357c4762a1bSJed Brown vb = zero; 358c4762a1bSJed Brown vl = zero; 359c4762a1bSJed Brown v = zero; 360c4762a1bSJed Brown if (i < nx && j > 0) vb = x[k-nx]; 361c4762a1bSJed Brown if (i > 0 && j < ny) vl = x[k-1]; 362c4762a1bSJed Brown if (i < nx && j < ny) v = x[k]; 363c4762a1bSJed Brown dvdx = (v-vl)/hx; 364c4762a1bSJed Brown dvdy = (v-vb)/hy; 365c4762a1bSJed Brown if (i != nx && j != 0) { 366c4762a1bSJed Brown ind = k-nx; val = - dvdy/hy - cdiv3; 3679566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 368c4762a1bSJed Brown } 369c4762a1bSJed Brown if (i != 0 && j != ny) { 370c4762a1bSJed Brown ind = k-1; val = - dvdx/hx - cdiv3; 3719566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown if (i != nx && j != ny) { 374c4762a1bSJed Brown ind = k; val = dvdx/hx + dvdy/hy - cdiv3; 3759566063dSJacob Faibussowitsch PetscCall(VecSetValues(G,1,&ind,&val,ADD_VALUES)); 376c4762a1bSJed Brown } 377c4762a1bSJed Brown } 378c4762a1bSJed Brown } 3799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 380c4762a1bSJed Brown 381c4762a1bSJed Brown /* Assemble gradient vector */ 3829566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(G)); 3839566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(G)); 384c4762a1bSJed Brown 385c4762a1bSJed Brown /* Scale the gradient */ 386c4762a1bSJed Brown area = p5*hx*hy; 3879566063dSJacob Faibussowitsch PetscCall(VecScale(G, area)); 3889566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(24.0*nx*ny)); 389c4762a1bSJed Brown PetscFunctionReturn(0); 390c4762a1bSJed Brown } 391c4762a1bSJed Brown 392c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 393c4762a1bSJed Brown /* 394c4762a1bSJed Brown FormHessian - Forms the Hessian matrix. 395c4762a1bSJed Brown 396c4762a1bSJed Brown Input Parameters: 397c4762a1bSJed Brown . tao - the Tao context 398c4762a1bSJed Brown . X - the input vector 399c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 400c4762a1bSJed Brown 401c4762a1bSJed Brown Output Parameters: 402c4762a1bSJed Brown . H - Hessian matrix 403c4762a1bSJed Brown . PrecH - optionally different preconditioning Hessian 404c4762a1bSJed Brown . flag - flag indicating matrix structure 405c4762a1bSJed Brown 406c4762a1bSJed Brown Notes: 407c4762a1bSJed Brown This routine is intended simply as an example of the interface 408c4762a1bSJed Brown to a Hessian evaluation routine. Since this example compute the 409c4762a1bSJed Brown Hessian a column at a time, it is not particularly efficient and 410c4762a1bSJed Brown is not recommended. 411c4762a1bSJed Brown */ 412c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H,Mat Hpre, void *ptr) 413c4762a1bSJed Brown { 414c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 415c4762a1bSJed Brown PetscInt i,j, ndim = user->ndim; 416c4762a1bSJed Brown PetscReal *y, zero = 0.0, one = 1.0; 417c4762a1bSJed Brown PetscBool assembled; 418c4762a1bSJed Brown 419c4762a1bSJed Brown PetscFunctionBeginUser; 420c4762a1bSJed Brown user->xvec = X; 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* Initialize Hessian entries and work vector to zero */ 4239566063dSJacob Faibussowitsch PetscCall(MatAssembled(H,&assembled)); 4249566063dSJacob Faibussowitsch if (assembled)PetscCall(MatZeroEntries(H)); 425c4762a1bSJed Brown 4269566063dSJacob Faibussowitsch PetscCall(VecSet(user->s, zero)); 427c4762a1bSJed Brown 428c4762a1bSJed Brown /* Loop over matrix columns to compute entries of the Hessian */ 429c4762a1bSJed Brown for (j=0; j<ndim; j++) { 4309566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->s,1,&j,&one,INSERT_VALUES)); 4319566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->s)); 4329566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->s)); 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(HessianProduct(ptr,user->s,user->y)); 435c4762a1bSJed Brown 4369566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->s,1,&j,&zero,INSERT_VALUES)); 4379566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->s)); 4389566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->s)); 439c4762a1bSJed Brown 4409566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->y,&y)); 441c4762a1bSJed Brown for (i=0; i<ndim; i++) { 442c4762a1bSJed Brown if (y[i] != zero) { 4439566063dSJacob Faibussowitsch PetscCall(MatSetValues(H,1,&i,1,&j,&y[i],ADD_VALUES)); 444c4762a1bSJed Brown } 445c4762a1bSJed Brown } 4469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->y,&y)); 447c4762a1bSJed Brown } 4489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 4499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 450c4762a1bSJed Brown PetscFunctionReturn(0); 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 454c4762a1bSJed Brown /* 455c4762a1bSJed Brown MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector 456c4762a1bSJed Brown products. 457c4762a1bSJed Brown 458c4762a1bSJed Brown Input Parameters: 459c4762a1bSJed Brown . tao - the Tao context 460c4762a1bSJed Brown . X - the input vector 461c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 462c4762a1bSJed Brown 463c4762a1bSJed Brown Output Parameters: 464c4762a1bSJed Brown . H - Hessian matrix 465c4762a1bSJed Brown . PrecH - optionally different preconditioning Hessian 466c4762a1bSJed Brown . flag - flag indicating matrix structure 467c4762a1bSJed Brown */ 468c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao tao,Vec X,Mat H,Mat PrecH, void *ptr) 469c4762a1bSJed Brown { 470c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 471c4762a1bSJed Brown 472c4762a1bSJed Brown /* Sets location of vector for use in computing matrix-vector products of the form H(X)*y */ 473362febeeSStefano Zampini PetscFunctionBeginUser; 474c4762a1bSJed Brown user->xvec = X; 475c4762a1bSJed Brown PetscFunctionReturn(0); 476c4762a1bSJed Brown } 477c4762a1bSJed Brown 478c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 479c4762a1bSJed Brown /* 480c4762a1bSJed Brown HessianProductMat - Computes the matrix-vector product 481c4762a1bSJed Brown y = mat*svec. 482c4762a1bSJed Brown 483c4762a1bSJed Brown Input Parameters: 484c4762a1bSJed Brown . mat - input matrix 485c4762a1bSJed Brown . svec - input vector 486c4762a1bSJed Brown 487c4762a1bSJed Brown Output Parameters: 488c4762a1bSJed Brown . y - solution vector 489c4762a1bSJed Brown */ 490c4762a1bSJed Brown PetscErrorCode HessianProductMat(Mat mat,Vec svec,Vec y) 491c4762a1bSJed Brown { 492c4762a1bSJed Brown void *ptr; 493c4762a1bSJed Brown 494c4762a1bSJed Brown PetscFunctionBeginUser; 4959566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ptr)); 4969566063dSJacob Faibussowitsch PetscCall(HessianProduct(ptr,svec,y)); 497c4762a1bSJed Brown PetscFunctionReturn(0); 498c4762a1bSJed Brown } 499c4762a1bSJed Brown 500c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 501c4762a1bSJed Brown /* 502c4762a1bSJed Brown Hessian Product - Computes the matrix-vector product: 503c4762a1bSJed Brown y = f''(x)*svec. 504c4762a1bSJed Brown 5057a7aea1fSJed Brown Input Parameters: 506c4762a1bSJed Brown . ptr - pointer to the user-defined context 507c4762a1bSJed Brown . svec - input vector 508c4762a1bSJed Brown 509c4762a1bSJed Brown Output Parameters: 510c4762a1bSJed Brown . y - product vector 511c4762a1bSJed Brown */ 512c4762a1bSJed Brown PetscErrorCode HessianProduct(void *ptr,Vec svec,Vec y) 513c4762a1bSJed Brown { 514c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 515c4762a1bSJed Brown PetscReal p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area; 516c4762a1bSJed Brown const PetscScalar *x, *s; 517c4762a1bSJed Brown PetscReal v, vb, vl, vr, vt, hxhx, hyhy; 518c4762a1bSJed Brown PetscInt nx, ny, i, j, k, ind; 519c4762a1bSJed Brown 520c4762a1bSJed Brown PetscFunctionBeginUser; 521c4762a1bSJed Brown nx = user->mx; 522c4762a1bSJed Brown ny = user->my; 523c4762a1bSJed Brown hx = user->hx; 524c4762a1bSJed Brown hy = user->hy; 525c4762a1bSJed Brown hxhx = one/(hx*hx); 526c4762a1bSJed Brown hyhy = one/(hy*hy); 527c4762a1bSJed Brown 528c4762a1bSJed Brown /* Get pointers to vector data */ 5299566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user->xvec,&x)); 5309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(svec,&s)); 531c4762a1bSJed Brown 532c4762a1bSJed Brown /* Initialize product vector to zero */ 5339566063dSJacob Faibussowitsch PetscCall(VecSet(y, zero)); 534c4762a1bSJed Brown 535c4762a1bSJed Brown /* Compute f''(x)*s over the lower triangular elements */ 536c4762a1bSJed Brown for (j=-1; j<ny; j++) { 537c4762a1bSJed Brown for (i=-1; i<nx; i++) { 538c4762a1bSJed Brown k = nx*j + i; 539c4762a1bSJed Brown v = zero; 540c4762a1bSJed Brown vr = zero; 541c4762a1bSJed Brown vt = zero; 542c4762a1bSJed Brown if (i != -1 && j != -1) v = s[k]; 543c4762a1bSJed Brown if (i != nx-1 && j != -1) { 544c4762a1bSJed Brown vr = s[k+1]; 545c4762a1bSJed Brown ind = k+1; val = hxhx*(vr-v); 5469566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 547c4762a1bSJed Brown } 548c4762a1bSJed Brown if (i != -1 && j != ny-1) { 549c4762a1bSJed Brown vt = s[k+nx]; 550c4762a1bSJed Brown ind = k+nx; val = hyhy*(vt-v); 5519566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 552c4762a1bSJed Brown } 553c4762a1bSJed Brown if (i != -1 && j != -1) { 554c4762a1bSJed Brown ind = k; val = hxhx*(v-vr) + hyhy*(v-vt); 5559566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 556c4762a1bSJed Brown } 557c4762a1bSJed Brown } 558c4762a1bSJed Brown } 559c4762a1bSJed Brown 560c4762a1bSJed Brown /* Compute f''(x)*s over the upper triangular elements */ 561c4762a1bSJed Brown for (j=0; j<=ny; j++) { 562c4762a1bSJed Brown for (i=0; i<=nx; i++) { 563c4762a1bSJed Brown k = nx*j + i; 564c4762a1bSJed Brown v = zero; 565c4762a1bSJed Brown vl = zero; 566c4762a1bSJed Brown vb = zero; 567c4762a1bSJed Brown if (i != nx && j != ny) v = s[k]; 568c4762a1bSJed Brown if (i != nx && j != 0) { 569c4762a1bSJed Brown vb = s[k-nx]; 570c4762a1bSJed Brown ind = k-nx; val = hyhy*(vb-v); 5719566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 572c4762a1bSJed Brown } 573c4762a1bSJed Brown if (i != 0 && j != ny) { 574c4762a1bSJed Brown vl = s[k-1]; 575c4762a1bSJed Brown ind = k-1; val = hxhx*(vl-v); 5769566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 577c4762a1bSJed Brown } 578c4762a1bSJed Brown if (i != nx && j != ny) { 579c4762a1bSJed Brown ind = k; val = hxhx*(v-vl) + hyhy*(v-vb); 5809566063dSJacob Faibussowitsch PetscCall(VecSetValues(y,1,&ind,&val,ADD_VALUES)); 581c4762a1bSJed Brown } 582c4762a1bSJed Brown } 583c4762a1bSJed Brown } 584c4762a1bSJed Brown /* Restore vector data */ 5859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(svec,&s)); 5869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user->xvec,&x)); 587c4762a1bSJed Brown 588c4762a1bSJed Brown /* Assemble vector */ 5899566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(y)); 5909566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(y)); 591c4762a1bSJed Brown 592c4762a1bSJed Brown /* Scale resulting vector by area */ 593c4762a1bSJed Brown area = p5*hx*hy; 5949566063dSJacob Faibussowitsch PetscCall(VecScale(y, area)); 5959566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(18.0*nx*ny)); 596c4762a1bSJed Brown PetscFunctionReturn(0); 597c4762a1bSJed Brown } 598c4762a1bSJed Brown 599c4762a1bSJed Brown /*TEST 600c4762a1bSJed Brown 601c4762a1bSJed Brown build: 602c4762a1bSJed Brown requires: !complex 603c4762a1bSJed Brown 604c4762a1bSJed Brown test: 605c4762a1bSJed Brown suffix: 1 606c4762a1bSJed Brown args: -tao_smonitor -tao_type ntl -tao_gatol 1.e-4 607c4762a1bSJed Brown 608c4762a1bSJed Brown test: 609c4762a1bSJed Brown suffix: 2 610c4762a1bSJed Brown args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4 611c4762a1bSJed Brown 612c4762a1bSJed Brown test: 613c4762a1bSJed Brown suffix: 3 614c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4 -my_tao_mf -tao_test_hessian 615c4762a1bSJed Brown 616c4762a1bSJed Brown test: 617c4762a1bSJed Brown suffix: 4 618c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnls 619c4762a1bSJed Brown 620c4762a1bSJed Brown test: 621c4762a1bSJed Brown suffix: 5 622c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-3 -tao_type blmvm 623c4762a1bSJed Brown 624c4762a1bSJed Brown test: 625c4762a1bSJed Brown suffix: 6 626c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 627c4762a1bSJed Brown 628c4762a1bSJed Brown TEST*/ 629