1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 rosenbrock2 [-help] [all TAO options] */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* Include "petsctao.h" so we can use TAO solvers. */ 4c4762a1bSJed Brown #include <petsctao.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown static char help[] = "This example demonstrates use of the TAO package to \n\ 7c4762a1bSJed Brown solve an unconstrained minimization problem on a single processor. We \n\ 8c4762a1bSJed Brown minimize the extended Rosenbrock function: \n\ 9c4762a1bSJed Brown sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) \n\ 10c4762a1bSJed Brown or the chained Rosenbrock function:\n\ 11c4762a1bSJed Brown sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n"; 12c4762a1bSJed Brown 13c4762a1bSJed Brown /*T 14c4762a1bSJed Brown Concepts: TAO^Solving an unconstrained minimization problem 15c4762a1bSJed Brown Routines: TaoCreate(); 16a82e8c82SStefano Zampini Routines: TaoSetType(); TaoSetObjectiveAndGradient(); 17a82e8c82SStefano Zampini Routines: TaoSetHessian(); 18a82e8c82SStefano Zampini Routines: TaoSetSolution(); 19c4762a1bSJed Brown Routines: TaoSetFromOptions(); 20c4762a1bSJed Brown Routines: TaoSolve(); 21c4762a1bSJed Brown Routines: TaoDestroy(); 22c4762a1bSJed Brown Processors: 1 23c4762a1bSJed Brown T*/ 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* 26c4762a1bSJed Brown User-defined application context - contains data needed by the 27c4762a1bSJed Brown application-provided call-back routines that evaluate the function, 28c4762a1bSJed Brown gradient, and hessian. 29c4762a1bSJed Brown */ 30c4762a1bSJed Brown typedef struct { 31c4762a1bSJed Brown PetscInt n; /* dimension */ 32c4762a1bSJed Brown PetscReal alpha; /* condition parameter */ 33c4762a1bSJed Brown PetscBool chained; 34c4762a1bSJed Brown } AppCtx; 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* -------------- User-defined routines ---------- */ 37c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 38c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 39c4762a1bSJed Brown 40c4762a1bSJed Brown int main(int argc,char **argv) 41c4762a1bSJed Brown { 42c4762a1bSJed Brown PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 43c4762a1bSJed Brown PetscReal zero=0.0; 44c4762a1bSJed Brown Vec x; /* solution vector */ 45c4762a1bSJed Brown Mat H; 46c4762a1bSJed Brown Tao tao; /* Tao solver context */ 47c4762a1bSJed Brown PetscBool flg, test_lmvm = PETSC_FALSE; 48c4762a1bSJed Brown PetscMPIInt size; /* number of processes running */ 49c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 50c4762a1bSJed Brown TaoConvergedReason reason; 51c4762a1bSJed Brown PetscInt its, recycled_its=0, oneshot_its=0; 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Initialize TAO and PETSc */ 54c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 55ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 56*3c859ba3SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* Initialize problem parameters */ 59c4762a1bSJed Brown user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE; 60c4762a1bSJed Brown /* Check for command line arguments to override defaults */ 61c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg);CHKERRQ(ierr); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Allocate vectors for the solution and gradient */ 67c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,user.n,&x);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H);CHKERRQ(ierr); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* The TAO code begins here */ 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Create TAO solver with desired solution method */ 73c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_SELF,&tao);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = TaoSetType(tao,TAOLMVM);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Set solution vec and an initial guess */ 77c4762a1bSJed Brown ierr = VecSet(x, zero);CHKERRQ(ierr); 78a82e8c82SStefano Zampini ierr = TaoSetSolution(tao,x);CHKERRQ(ierr); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Set routines for function, gradient, hessian evaluation */ 81a82e8c82SStefano Zampini ierr = TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user);CHKERRQ(ierr); 82a82e8c82SStefano Zampini ierr = TaoSetHessian(tao,H,H,FormHessian,&user);CHKERRQ(ierr); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Check for TAO command line options */ 85c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Solve the problem */ 88c4762a1bSJed Brown ierr = TaoSetTolerances(tao, 1.e-5, 0.0, 0.0);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = TaoSetMaximumIterations(tao, 5);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = TaoLMVMRecycle(tao, PETSC_TRUE);CHKERRQ(ierr); 91c4762a1bSJed Brown reason = TAO_CONTINUE_ITERATING; 92c4762a1bSJed Brown while (reason != TAO_CONVERGED_GATOL) { 93c4762a1bSJed Brown ierr = TaoSolve(tao);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = TaoGetConvergedReason(tao, &reason);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = TaoGetIterationNumber(tao, &its);CHKERRQ(ierr); 96c4762a1bSJed Brown recycled_its += its; 97c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "-----------------------\n");CHKERRQ(ierr); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* Disable recycling and solve again! */ 101c4762a1bSJed Brown ierr = TaoSetMaximumIterations(tao, 100);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = TaoLMVMRecycle(tao, PETSC_FALSE);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = VecSet(x, zero);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = TaoSolve(tao);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = TaoGetConvergedReason(tao, &reason);CHKERRQ(ierr); 106*3c859ba3SBarry Smith PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!"); 107c4762a1bSJed Brown ierr = TaoGetIterationNumber(tao, &oneshot_its);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "-----------------------\n");CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "recycled its: %D | oneshot its: %D\n", recycled_its, oneshot_its);CHKERRQ(ierr); 110*3c859ba3SBarry Smith PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "LMVM recycling does not work!"); 111c4762a1bSJed Brown 112c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = MatDestroy(&H);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown ierr = PetscFinalize(); 117c4762a1bSJed Brown return ierr; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X). 123c4762a1bSJed Brown 124c4762a1bSJed Brown Input Parameters: 125c4762a1bSJed Brown . tao - the Tao context 126c4762a1bSJed Brown . X - input vector 127c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunctionGradient() 128c4762a1bSJed Brown 129c4762a1bSJed Brown Output Parameters: 130c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 131c4762a1bSJed Brown . f - function value 132c4762a1bSJed Brown 133c4762a1bSJed Brown Note: 134c4762a1bSJed Brown Some optimization methods ask for the function and the gradient evaluation 135c4762a1bSJed Brown at the same time. Evaluating both at once may be more efficient that 136c4762a1bSJed Brown evaluating each separately. 137c4762a1bSJed Brown */ 138c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr) 139c4762a1bSJed Brown { 140c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 141c4762a1bSJed Brown PetscInt i,nn=user->n/2; 142c4762a1bSJed Brown PetscErrorCode ierr; 143c4762a1bSJed Brown PetscReal ff=0,t1,t2,alpha=user->alpha; 144c4762a1bSJed Brown PetscScalar *g; 145c4762a1bSJed Brown const PetscScalar *x; 146c4762a1bSJed Brown 147c4762a1bSJed Brown PetscFunctionBeginUser; 148c4762a1bSJed Brown /* Get pointers to vector data */ 149c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = VecGetArray(G,&g);CHKERRQ(ierr); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Compute G(X) */ 153c4762a1bSJed Brown if (user->chained) { 154c4762a1bSJed Brown g[0] = 0; 155c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 156c4762a1bSJed Brown t1 = x[i+1] - x[i]*x[i]; 157c4762a1bSJed Brown ff += PetscSqr(1 - x[i]) + alpha*t1*t1; 158c4762a1bSJed Brown g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]); 159c4762a1bSJed Brown g[i+1] = 2*alpha*t1; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } else { 162c4762a1bSJed Brown for (i=0; i<nn; i++) { 163c4762a1bSJed Brown t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i]; 164c4762a1bSJed Brown ff += alpha*t1*t1 + t2*t2; 165c4762a1bSJed Brown g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2; 166c4762a1bSJed Brown g[2*i+1] = 2*alpha*t1; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* Restore vectors */ 171c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = VecRestoreArray(G,&g);CHKERRQ(ierr); 173c4762a1bSJed Brown *f = ff; 174c4762a1bSJed Brown 175ca0c957dSBarry Smith ierr = PetscLogFlops(15.0*nn);CHKERRQ(ierr); 176c4762a1bSJed Brown PetscFunctionReturn(0); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 180c4762a1bSJed Brown /* 181c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 182c4762a1bSJed Brown 183c4762a1bSJed Brown Input Parameters: 184c4762a1bSJed Brown . tao - the Tao context 185c4762a1bSJed Brown . x - input vector 186c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 187c4762a1bSJed Brown 188c4762a1bSJed Brown Output Parameters: 189c4762a1bSJed Brown . H - Hessian matrix 190c4762a1bSJed Brown 191c4762a1bSJed Brown Note: Providing the Hessian may not be necessary. Only some solvers 192c4762a1bSJed Brown require this matrix. 193c4762a1bSJed Brown */ 194c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 195c4762a1bSJed Brown { 196c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 197c4762a1bSJed Brown PetscErrorCode ierr; 198c4762a1bSJed Brown PetscInt i, ind[2]; 199c4762a1bSJed Brown PetscReal alpha=user->alpha; 200c4762a1bSJed Brown PetscReal v[2][2]; 201c4762a1bSJed Brown const PetscScalar *x; 202c4762a1bSJed Brown PetscBool assembled; 203c4762a1bSJed Brown 204c4762a1bSJed Brown PetscFunctionBeginUser; 205c4762a1bSJed Brown /* Zero existing matrix entries */ 206c4762a1bSJed Brown ierr = MatAssembled(H,&assembled);CHKERRQ(ierr); 207c4762a1bSJed Brown if (assembled) {ierr = MatZeroEntries(H);CHKERRQ(ierr);} 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* Get a pointer to vector data */ 210c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Compute H(X) entries */ 213c4762a1bSJed Brown if (user->chained) { 214c4762a1bSJed Brown ierr = MatZeroEntries(H);CHKERRQ(ierr); 215c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 216c4762a1bSJed Brown PetscScalar t1 = x[i+1] - x[i]*x[i]; 217c4762a1bSJed Brown v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]); 218c4762a1bSJed Brown v[0][1] = 2*alpha*(-2*x[i]); 219c4762a1bSJed Brown v[1][0] = 2*alpha*(-2*x[i]); 220c4762a1bSJed Brown v[1][1] = 2*alpha*t1; 221c4762a1bSJed Brown ind[0] = i; ind[1] = i+1; 222c4762a1bSJed Brown ierr = MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES);CHKERRQ(ierr); 223c4762a1bSJed Brown } 224c4762a1bSJed Brown } else { 225c4762a1bSJed Brown for (i=0; i<user->n/2; i++) { 226c4762a1bSJed Brown v[1][1] = 2*alpha; 227c4762a1bSJed Brown v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2; 228c4762a1bSJed Brown v[1][0] = v[0][1] = -4.0*alpha*x[2*i]; 229c4762a1bSJed Brown ind[0]=2*i; ind[1]=2*i+1; 230c4762a1bSJed Brown ierr = MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);CHKERRQ(ierr); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown } 233c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* Assemble matrix */ 236c4762a1bSJed Brown ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = PetscLogFlops(9.0*user->n/2.0);CHKERRQ(ierr); 239c4762a1bSJed Brown PetscFunctionReturn(0); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown /*TEST 243c4762a1bSJed Brown 244c4762a1bSJed Brown build: 245c4762a1bSJed Brown requires: !complex 246c4762a1bSJed Brown 247c4762a1bSJed Brown test: 248c4762a1bSJed Brown args: -tao_type lmvm -tao_monitor 249c4762a1bSJed Brown requires: !single 250c4762a1bSJed Brown 251c4762a1bSJed Brown TEST*/ 252