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 PetscReal zero=0.0; 43c4762a1bSJed Brown Vec x; /* solution vector */ 44c4762a1bSJed Brown Mat H; 45c4762a1bSJed Brown Tao tao; /* Tao solver context */ 46c4762a1bSJed Brown PetscBool flg, test_lmvm = PETSC_FALSE; 47c4762a1bSJed Brown PetscMPIInt size; /* number of processes running */ 48c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 49c4762a1bSJed Brown TaoConvergedReason reason; 50c4762a1bSJed Brown PetscInt its, recycled_its=0, oneshot_its=0; 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Initialize TAO and PETSc */ 53*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 545f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 553c859ba3SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Initialize problem parameters */ 58c4762a1bSJed Brown user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE; 59c4762a1bSJed Brown /* Check for command line arguments to override defaults */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Allocate vectors for the solution and gradient */ 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,user.n,&x)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* The TAO code begins here */ 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Create TAO solver with desired solution method */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOLMVM)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* Set solution vec and an initial guess */ 765f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x, zero)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao,x)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Set routines for function, gradient, hessian evaluation */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,H,H,FormHessian,&user)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Check for TAO command line options */ 845f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Solve the problem */ 875f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetMaximumIterations(tao, 5)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLMVMRecycle(tao, PETSC_TRUE)); 90c4762a1bSJed Brown reason = TAO_CONTINUE_ITERATING; 91c4762a1bSJed Brown while (reason != TAO_CONVERGED_GATOL) { 925f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetConvergedReason(tao, &reason)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetIterationNumber(tao, &its)); 95c4762a1bSJed Brown recycled_its += its; 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n")); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* Disable recycling and solve again! */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetMaximumIterations(tao, 100)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLMVMRecycle(tao, PETSC_FALSE)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x, zero)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetConvergedReason(tao, &reason)); 1053c859ba3SBarry Smith PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!"); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetIterationNumber(tao, &oneshot_its)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n")); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "recycled its: %D | oneshot its: %D\n", recycled_its, oneshot_its)); 1093c859ba3SBarry Smith PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "LMVM recycling does not work!"); 110c4762a1bSJed Brown 1115f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&H)); 114c4762a1bSJed Brown 115*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 116*b122ec5aSJacob Faibussowitsch return 0; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X). 122c4762a1bSJed Brown 123c4762a1bSJed Brown Input Parameters: 124c4762a1bSJed Brown . tao - the Tao context 125c4762a1bSJed Brown . X - input vector 126c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunctionGradient() 127c4762a1bSJed Brown 128c4762a1bSJed Brown Output Parameters: 129c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 130c4762a1bSJed Brown . f - function value 131c4762a1bSJed Brown 132c4762a1bSJed Brown Note: 133c4762a1bSJed Brown Some optimization methods ask for the function and the gradient evaluation 134c4762a1bSJed Brown at the same time. Evaluating both at once may be more efficient that 135c4762a1bSJed Brown evaluating each separately. 136c4762a1bSJed Brown */ 137c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr) 138c4762a1bSJed Brown { 139c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 140c4762a1bSJed Brown PetscInt i,nn=user->n/2; 141c4762a1bSJed Brown PetscReal ff=0,t1,t2,alpha=user->alpha; 142c4762a1bSJed Brown PetscScalar *g; 143c4762a1bSJed Brown const PetscScalar *x; 144c4762a1bSJed Brown 145c4762a1bSJed Brown PetscFunctionBeginUser; 146c4762a1bSJed Brown /* Get pointers to vector data */ 1475f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(G,&g)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* Compute G(X) */ 151c4762a1bSJed Brown if (user->chained) { 152c4762a1bSJed Brown g[0] = 0; 153c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 154c4762a1bSJed Brown t1 = x[i+1] - x[i]*x[i]; 155c4762a1bSJed Brown ff += PetscSqr(1 - x[i]) + alpha*t1*t1; 156c4762a1bSJed Brown g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]); 157c4762a1bSJed Brown g[i+1] = 2*alpha*t1; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } else { 160c4762a1bSJed Brown for (i=0; i<nn; i++) { 161c4762a1bSJed Brown t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i]; 162c4762a1bSJed Brown ff += alpha*t1*t1 + t2*t2; 163c4762a1bSJed Brown g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2; 164c4762a1bSJed Brown g[2*i+1] = 2*alpha*t1; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Restore vectors */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(G,&g)); 171c4762a1bSJed Brown *f = ff; 172c4762a1bSJed Brown 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(15.0*nn)); 174c4762a1bSJed Brown PetscFunctionReturn(0); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 180c4762a1bSJed Brown 181c4762a1bSJed Brown Input Parameters: 182c4762a1bSJed Brown . tao - the Tao context 183c4762a1bSJed Brown . x - input vector 184c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 185c4762a1bSJed Brown 186c4762a1bSJed Brown Output Parameters: 187c4762a1bSJed Brown . H - Hessian matrix 188c4762a1bSJed Brown 189c4762a1bSJed Brown Note: Providing the Hessian may not be necessary. Only some solvers 190c4762a1bSJed Brown require this matrix. 191c4762a1bSJed Brown */ 192c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 193c4762a1bSJed Brown { 194c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 195c4762a1bSJed Brown PetscInt i, ind[2]; 196c4762a1bSJed Brown PetscReal alpha=user->alpha; 197c4762a1bSJed Brown PetscReal v[2][2]; 198c4762a1bSJed Brown const PetscScalar *x; 199c4762a1bSJed Brown PetscBool assembled; 200c4762a1bSJed Brown 201c4762a1bSJed Brown PetscFunctionBeginUser; 202c4762a1bSJed Brown /* Zero existing matrix entries */ 2035f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssembled(H,&assembled)); 2045f80ce2aSJacob Faibussowitsch if (assembled) CHKERRQ(MatZeroEntries(H)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* Get a pointer to vector data */ 2075f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* Compute H(X) entries */ 210c4762a1bSJed Brown if (user->chained) { 2115f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(H)); 212c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 213c4762a1bSJed Brown PetscScalar t1 = x[i+1] - x[i]*x[i]; 214c4762a1bSJed Brown v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]); 215c4762a1bSJed Brown v[0][1] = 2*alpha*(-2*x[i]); 216c4762a1bSJed Brown v[1][0] = 2*alpha*(-2*x[i]); 217c4762a1bSJed Brown v[1][1] = 2*alpha*t1; 218c4762a1bSJed Brown ind[0] = i; ind[1] = i+1; 2195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES)); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown } else { 222c4762a1bSJed Brown for (i=0; i<user->n/2; i++) { 223c4762a1bSJed Brown v[1][1] = 2*alpha; 224c4762a1bSJed Brown v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2; 225c4762a1bSJed Brown v[1][0] = v[0][1] = -4.0*alpha*x[2*i]; 226c4762a1bSJed Brown ind[0]=2*i; ind[1]=2*i+1; 2275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES)); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown } 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* Assemble matrix */ 2335f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(9.0*user->n/2.0)); 236c4762a1bSJed Brown PetscFunctionReturn(0); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 239c4762a1bSJed Brown /*TEST 240c4762a1bSJed Brown 241c4762a1bSJed Brown build: 242c4762a1bSJed Brown requires: !complex 243c4762a1bSJed Brown 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown args: -tao_type lmvm -tao_monitor 246c4762a1bSJed Brown requires: !single 247c4762a1bSJed Brown 248c4762a1bSJed Brown TEST*/ 249