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 /* 14c4762a1bSJed Brown User-defined application context - contains data needed by the 15c4762a1bSJed Brown application-provided call-back routines that evaluate the function, 16c4762a1bSJed Brown gradient, and hessian. 17c4762a1bSJed Brown */ 18c4762a1bSJed Brown typedef struct { 19c4762a1bSJed Brown PetscInt n; /* dimension */ 20c4762a1bSJed Brown PetscReal alpha; /* condition parameter */ 21c4762a1bSJed Brown PetscBool chained; 22c4762a1bSJed Brown } AppCtx; 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* -------------- User-defined routines ---------- */ 25c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 26c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 27c4762a1bSJed Brown 28*9371c9d4SSatish Balay int main(int argc, char **argv) { 29c4762a1bSJed Brown PetscReal zero = 0.0; 30c4762a1bSJed Brown Vec x; /* solution vector */ 31c4762a1bSJed Brown Mat H; 32c4762a1bSJed Brown Tao tao; /* Tao solver context */ 33c4762a1bSJed Brown PetscBool flg, test_lmvm = PETSC_FALSE; 34c4762a1bSJed Brown PetscMPIInt size; /* number of processes running */ 35c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 36c4762a1bSJed Brown TaoConvergedReason reason; 37c4762a1bSJed Brown PetscInt its, recycled_its = 0, oneshot_its = 0; 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Initialize TAO and PETSc */ 40327415f7SBarry Smith PetscFunctionBeginUser; 419566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 433c859ba3SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors"); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Initialize problem parameters */ 46*9371c9d4SSatish Balay user.n = 2; 47*9371c9d4SSatish Balay user.alpha = 99.0; 48*9371c9d4SSatish Balay user.chained = PETSC_FALSE; 49c4762a1bSJed Brown /* Check for command line arguments to override defaults */ 509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &user.n, &flg)); 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-alpha", &user.alpha, &flg)); 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-chained", &user.chained, &flg)); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Allocate vectors for the solution and gradient */ 569566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.n, &x)); 579566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, 2, user.n, user.n, 1, NULL, &H)); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* The TAO code begins here */ 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Create TAO solver with desired solution method */ 629566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 639566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLMVM)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Set solution vec and an initial guess */ 669566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero)); 679566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Set routines for function, gradient, hessian evaluation */ 709566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, &user)); 719566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, FormHessian, &user)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Check for TAO command line options */ 749566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Solve the problem */ 779566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0)); 789566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(tao, 5)); 799566063dSJacob Faibussowitsch PetscCall(TaoLMVMRecycle(tao, PETSC_TRUE)); 80c4762a1bSJed Brown reason = TAO_CONTINUE_ITERATING; 81c4762a1bSJed Brown while (reason != TAO_CONVERGED_GATOL) { 829566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 839566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao, &reason)); 849566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao, &its)); 85c4762a1bSJed Brown recycled_its += its; 869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n")); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* Disable recycling and solve again! */ 909566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(tao, 100)); 919566063dSJacob Faibussowitsch PetscCall(TaoLMVMRecycle(tao, PETSC_FALSE)); 929566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero)); 939566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 949566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao, &reason)); 953c859ba3SBarry Smith PetscCheck(reason == TAO_CONVERGED_GATOL, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!"); 969566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao, &oneshot_its)); 979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n")); 9863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "recycled its: %" PetscInt_FMT " | oneshot its: %" PetscInt_FMT "\n", recycled_its, oneshot_its)); 993c859ba3SBarry Smith PetscCheck(recycled_its == oneshot_its, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "LMVM recycling does not work!"); 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 104c4762a1bSJed Brown 1059566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 106b122ec5aSJacob Faibussowitsch return 0; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 110c4762a1bSJed Brown /* 111c4762a1bSJed Brown FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X). 112c4762a1bSJed Brown 113c4762a1bSJed Brown Input Parameters: 114c4762a1bSJed Brown . tao - the Tao context 115c4762a1bSJed Brown . X - input vector 116c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunctionGradient() 117c4762a1bSJed Brown 118c4762a1bSJed Brown Output Parameters: 119c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 120c4762a1bSJed Brown . f - function value 121c4762a1bSJed Brown 122c4762a1bSJed Brown Note: 123c4762a1bSJed Brown Some optimization methods ask for the function and the gradient evaluation 124c4762a1bSJed Brown at the same time. Evaluating both at once may be more efficient that 125c4762a1bSJed Brown evaluating each separately. 126c4762a1bSJed Brown */ 127*9371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) { 128c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 129c4762a1bSJed Brown PetscInt i, nn = user->n / 2; 130c4762a1bSJed Brown PetscReal ff = 0, t1, t2, alpha = user->alpha; 131c4762a1bSJed Brown PetscScalar *g; 132c4762a1bSJed Brown const PetscScalar *x; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBeginUser; 135c4762a1bSJed Brown /* Get pointers to vector data */ 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Compute G(X) */ 140c4762a1bSJed Brown if (user->chained) { 141c4762a1bSJed Brown g[0] = 0; 142c4762a1bSJed Brown for (i = 0; i < user->n - 1; i++) { 143c4762a1bSJed Brown t1 = x[i + 1] - x[i] * x[i]; 144c4762a1bSJed Brown ff += PetscSqr(1 - x[i]) + alpha * t1 * t1; 145c4762a1bSJed Brown g[i] += -2 * (1 - x[i]) + 2 * alpha * t1 * (-2 * x[i]); 146c4762a1bSJed Brown g[i + 1] = 2 * alpha * t1; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } else { 149c4762a1bSJed Brown for (i = 0; i < nn; i++) { 150*9371c9d4SSatish Balay t1 = x[2 * i + 1] - x[2 * i] * x[2 * i]; 151*9371c9d4SSatish Balay t2 = 1 - x[2 * i]; 152c4762a1bSJed Brown ff += alpha * t1 * t1 + t2 * t2; 153c4762a1bSJed Brown g[2 * i] = -4 * alpha * t1 * x[2 * i] - 2.0 * t2; 154c4762a1bSJed Brown g[2 * i + 1] = 2 * alpha * t1; 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Restore vectors */ 1599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g)); 161c4762a1bSJed Brown *f = ff; 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(15.0 * nn)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 170c4762a1bSJed Brown 171c4762a1bSJed Brown Input Parameters: 172c4762a1bSJed Brown . tao - the Tao context 173c4762a1bSJed Brown . x - input vector 174c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 175c4762a1bSJed Brown 176c4762a1bSJed Brown Output Parameters: 177c4762a1bSJed Brown . H - Hessian matrix 178c4762a1bSJed Brown 179c4762a1bSJed Brown Note: Providing the Hessian may not be necessary. Only some solvers 180c4762a1bSJed Brown require this matrix. 181c4762a1bSJed Brown */ 182*9371c9d4SSatish Balay PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) { 183c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 184c4762a1bSJed Brown PetscInt i, ind[2]; 185c4762a1bSJed Brown PetscReal alpha = user->alpha; 186c4762a1bSJed Brown PetscReal v[2][2]; 187c4762a1bSJed Brown const PetscScalar *x; 188c4762a1bSJed Brown PetscBool assembled; 189c4762a1bSJed Brown 190c4762a1bSJed Brown PetscFunctionBeginUser; 191c4762a1bSJed Brown /* Zero existing matrix entries */ 1929566063dSJacob Faibussowitsch PetscCall(MatAssembled(H, &assembled)); 1939566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(H)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* Get a pointer to vector data */ 1969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* Compute H(X) entries */ 199c4762a1bSJed Brown if (user->chained) { 2009566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 201c4762a1bSJed Brown for (i = 0; i < user->n - 1; i++) { 202c4762a1bSJed Brown PetscScalar t1 = x[i + 1] - x[i] * x[i]; 203c4762a1bSJed Brown v[0][0] = 2 + 2 * alpha * (t1 * (-2) - 2 * x[i]); 204c4762a1bSJed Brown v[0][1] = 2 * alpha * (-2 * x[i]); 205c4762a1bSJed Brown v[1][0] = 2 * alpha * (-2 * x[i]); 206c4762a1bSJed Brown v[1][1] = 2 * alpha * t1; 207*9371c9d4SSatish Balay ind[0] = i; 208*9371c9d4SSatish Balay ind[1] = i + 1; 2099566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], ADD_VALUES)); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown } else { 212c4762a1bSJed Brown for (i = 0; i < user->n / 2; i++) { 213c4762a1bSJed Brown v[1][1] = 2 * alpha; 214c4762a1bSJed Brown v[0][0] = -4 * alpha * (x[2 * i + 1] - 3 * x[2 * i] * x[2 * i]) + 2; 215c4762a1bSJed Brown v[1][0] = v[0][1] = -4.0 * alpha * x[2 * i]; 216*9371c9d4SSatish Balay ind[0] = 2 * i; 217*9371c9d4SSatish Balay ind[1] = 2 * i + 1; 2189566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], INSERT_VALUES)); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 2219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* Assemble matrix */ 2249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 2259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 2269566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0 * user->n / 2.0)); 227c4762a1bSJed Brown PetscFunctionReturn(0); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown /*TEST 231c4762a1bSJed Brown 232c4762a1bSJed Brown build: 233c4762a1bSJed Brown requires: !complex 234c4762a1bSJed Brown 235c4762a1bSJed Brown test: 236c4762a1bSJed Brown args: -tao_type lmvm -tao_monitor 237c4762a1bSJed Brown requires: !single 238c4762a1bSJed Brown 239c4762a1bSJed Brown TEST*/ 240