1414d97d3SAlp Dener /* Program usage: mpiexec -n 1 rosenbrock2 [-help] [all TAO options] */ 2414d97d3SAlp Dener 3414d97d3SAlp Dener /* Include "petsctao.h" so we can use TAO solvers. */ 4414d97d3SAlp Dener #include <petsctao.h> 5414d97d3SAlp Dener 6414d97d3SAlp Dener static char help[] = "This example demonstrates use of the TAO package to \n\ 7414d97d3SAlp Dener solve an unconstrained minimization problem on a single processor. We \n\ 8414d97d3SAlp Dener minimize the extended Rosenbrock function: \n\ 9414d97d3SAlp Dener sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) \n\ 10414d97d3SAlp Dener or the chained Rosenbrock function:\n\ 11414d97d3SAlp Dener sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n"; 12414d97d3SAlp Dener 13414d97d3SAlp Dener /*T 14414d97d3SAlp Dener Concepts: TAO^Solving an unconstrained minimization problem 15414d97d3SAlp Dener Routines: TaoCreate(); 16a82e8c82SStefano Zampini Routines: TaoSetType(); TaoSetObjectiveAndGradient(); 17a82e8c82SStefano Zampini Routines: TaoSetHessian(); 18a82e8c82SStefano Zampini Routines: TaoSetSolution(); 19414d97d3SAlp Dener Routines: TaoSetFromOptions(); 20414d97d3SAlp Dener Routines: TaoSolve(); 21414d97d3SAlp Dener Routines: TaoDestroy(); 22414d97d3SAlp Dener Processors: 1 23414d97d3SAlp Dener T*/ 24414d97d3SAlp Dener 25414d97d3SAlp Dener /* 26414d97d3SAlp Dener User-defined application context - contains data needed by the 27414d97d3SAlp Dener application-provided call-back routines that evaluate the function, 28414d97d3SAlp Dener gradient, and hessian. 29414d97d3SAlp Dener */ 30414d97d3SAlp Dener typedef struct { 31414d97d3SAlp Dener PetscInt n; /* dimension */ 32414d97d3SAlp Dener PetscReal alpha; /* condition parameter */ 33414d97d3SAlp Dener PetscBool chained; 34414d97d3SAlp Dener } AppCtx; 35414d97d3SAlp Dener 36414d97d3SAlp Dener /* -------------- User-defined routines ---------- */ 37414d97d3SAlp Dener PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 38414d97d3SAlp Dener PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 39414d97d3SAlp Dener 40414d97d3SAlp Dener int main(int argc,char **argv) 41414d97d3SAlp Dener { 42414d97d3SAlp Dener PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 43414d97d3SAlp Dener PetscReal zero=0.0; 44414d97d3SAlp Dener Vec x; /* solution vector */ 45414d97d3SAlp Dener Mat H; 46414d97d3SAlp Dener Tao tao; /* Tao solver context */ 47414d97d3SAlp Dener PetscBool flg, test_lmvm = PETSC_FALSE; 48414d97d3SAlp Dener PetscMPIInt size; /* number of processes running */ 49414d97d3SAlp Dener AppCtx user; /* user-defined application context */ 50414d97d3SAlp Dener TaoConvergedReason reason; 51414d97d3SAlp Dener PetscInt its, recycled_its=0, oneshot_its=0; 52414d97d3SAlp Dener 53414d97d3SAlp Dener /* Initialize TAO and PETSc */ 54414d97d3SAlp Dener ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 55414d97d3SAlp Dener 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"); 57414d97d3SAlp Dener 58414d97d3SAlp Dener /* Initialize problem parameters */ 59414d97d3SAlp Dener user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE; 60414d97d3SAlp Dener /* Check for command line arguments to override defaults */ 61414d97d3SAlp Dener ierr = PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg);CHKERRQ(ierr); 62414d97d3SAlp Dener ierr = PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg);CHKERRQ(ierr); 63414d97d3SAlp Dener ierr = PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg);CHKERRQ(ierr); 64414d97d3SAlp Dener ierr = PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg);CHKERRQ(ierr); 65414d97d3SAlp Dener 66414d97d3SAlp Dener /* Allocate vectors for the solution and gradient */ 67414d97d3SAlp Dener ierr = VecCreateSeq(PETSC_COMM_SELF,user.n,&x);CHKERRQ(ierr); 68414d97d3SAlp Dener ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H);CHKERRQ(ierr); 69414d97d3SAlp Dener 70414d97d3SAlp Dener /* The TAO code begins here */ 71414d97d3SAlp Dener 72414d97d3SAlp Dener /* Create TAO solver with desired solution method */ 73414d97d3SAlp Dener ierr = TaoCreate(PETSC_COMM_SELF,&tao);CHKERRQ(ierr); 74414d97d3SAlp Dener ierr = TaoSetType(tao,TAOBQNLS);CHKERRQ(ierr); 75414d97d3SAlp Dener 76414d97d3SAlp Dener /* Set solution vec and an initial guess */ 77414d97d3SAlp Dener ierr = VecSet(x, zero);CHKERRQ(ierr); 78a82e8c82SStefano Zampini ierr = TaoSetSolution(tao,x);CHKERRQ(ierr); 79414d97d3SAlp Dener 80414d97d3SAlp Dener /* 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); 83414d97d3SAlp Dener 84414d97d3SAlp Dener /* Check for TAO command line options */ 85414d97d3SAlp Dener ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 86414d97d3SAlp Dener 87414d97d3SAlp Dener /* Solve the problem */ 88414d97d3SAlp Dener ierr = TaoSetTolerances(tao, 1.e-5, 0.0, 0.0);CHKERRQ(ierr); 89414d97d3SAlp Dener ierr = TaoSetMaximumIterations(tao, 5);CHKERRQ(ierr); 90414d97d3SAlp Dener ierr = TaoSetRecycleHistory(tao, PETSC_TRUE);CHKERRQ(ierr); 91414d97d3SAlp Dener reason = TAO_CONTINUE_ITERATING; 92414d97d3SAlp Dener flg = PETSC_FALSE; 93414d97d3SAlp Dener ierr = TaoGetRecycleHistory(tao, &flg);CHKERRQ(ierr); 9460d4fc61SSatish Balay if (flg) {ierr = PetscPrintf(PETSC_COMM_SELF, "Recycle: enabled\n");CHKERRQ(ierr);} 95414d97d3SAlp Dener while (reason != TAO_CONVERGED_GATOL) { 96414d97d3SAlp Dener ierr = TaoSolve(tao);CHKERRQ(ierr); 97414d97d3SAlp Dener ierr = TaoGetConvergedReason(tao, &reason);CHKERRQ(ierr); 98414d97d3SAlp Dener ierr = TaoGetIterationNumber(tao, &its);CHKERRQ(ierr); 99414d97d3SAlp Dener recycled_its += its; 100414d97d3SAlp Dener ierr = PetscPrintf(PETSC_COMM_SELF, "-----------------------\n");CHKERRQ(ierr); 101414d97d3SAlp Dener } 102414d97d3SAlp Dener 103414d97d3SAlp Dener /* Disable recycling and solve again! */ 104414d97d3SAlp Dener ierr = TaoSetMaximumIterations(tao, 100);CHKERRQ(ierr); 105414d97d3SAlp Dener ierr = TaoSetRecycleHistory(tao, PETSC_FALSE);CHKERRQ(ierr); 106414d97d3SAlp Dener ierr = VecSet(x, zero);CHKERRQ(ierr); 107414d97d3SAlp Dener ierr = TaoGetRecycleHistory(tao, &flg);CHKERRQ(ierr); 10860d4fc61SSatish Balay if (!flg) {ierr = PetscPrintf(PETSC_COMM_SELF, "Recycle: disabled\n");CHKERRQ(ierr);} 109414d97d3SAlp Dener ierr = TaoSolve(tao);CHKERRQ(ierr); 110414d97d3SAlp Dener ierr = TaoGetConvergedReason(tao, &reason);CHKERRQ(ierr); 111*3c859ba3SBarry Smith PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!"); 112414d97d3SAlp Dener ierr = TaoGetIterationNumber(tao, &oneshot_its);CHKERRQ(ierr); 113414d97d3SAlp Dener ierr = PetscPrintf(PETSC_COMM_SELF, "-----------------------\n");CHKERRQ(ierr); 114414d97d3SAlp Dener ierr = PetscPrintf(PETSC_COMM_SELF, "recycled its: %D | oneshot its: %D\n", recycled_its, oneshot_its);CHKERRQ(ierr); 115*3c859ba3SBarry Smith PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Recycled solution does not match oneshot solution!"); 116414d97d3SAlp Dener 117414d97d3SAlp Dener ierr = TaoDestroy(&tao);CHKERRQ(ierr); 118414d97d3SAlp Dener ierr = VecDestroy(&x);CHKERRQ(ierr); 119414d97d3SAlp Dener ierr = MatDestroy(&H);CHKERRQ(ierr); 120414d97d3SAlp Dener 121414d97d3SAlp Dener ierr = PetscFinalize(); 122414d97d3SAlp Dener return ierr; 123414d97d3SAlp Dener } 124414d97d3SAlp Dener 125414d97d3SAlp Dener /* -------------------------------------------------------------------- */ 126414d97d3SAlp Dener /* 127414d97d3SAlp Dener FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X). 128414d97d3SAlp Dener 129414d97d3SAlp Dener Input Parameters: 130414d97d3SAlp Dener . tao - the Tao context 131414d97d3SAlp Dener . X - input vector 132414d97d3SAlp Dener . ptr - optional user-defined context, as set by TaoSetFunctionGradient() 133414d97d3SAlp Dener 134414d97d3SAlp Dener Output Parameters: 135414d97d3SAlp Dener . G - vector containing the newly evaluated gradient 136414d97d3SAlp Dener . f - function value 137414d97d3SAlp Dener 138414d97d3SAlp Dener Note: 139414d97d3SAlp Dener Some optimization methods ask for the function and the gradient evaluation 140414d97d3SAlp Dener at the same time. Evaluating both at once may be more efficient than 141414d97d3SAlp Dener evaluating each separately. 142414d97d3SAlp Dener */ 143414d97d3SAlp Dener PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr) 144414d97d3SAlp Dener { 145414d97d3SAlp Dener AppCtx *user = (AppCtx *) ptr; 146414d97d3SAlp Dener PetscInt i,nn=user->n/2; 147414d97d3SAlp Dener PetscErrorCode ierr; 148414d97d3SAlp Dener PetscReal ff=0,t1,t2,alpha=user->alpha; 149414d97d3SAlp Dener PetscScalar *g; 150414d97d3SAlp Dener const PetscScalar *x; 151414d97d3SAlp Dener 152414d97d3SAlp Dener PetscFunctionBeginUser; 153414d97d3SAlp Dener /* Get pointers to vector data */ 154414d97d3SAlp Dener ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 155414d97d3SAlp Dener ierr = VecGetArrayWrite(G,&g);CHKERRQ(ierr); 156414d97d3SAlp Dener 157414d97d3SAlp Dener /* Compute G(X) */ 158414d97d3SAlp Dener if (user->chained) { 159414d97d3SAlp Dener g[0] = 0; 160414d97d3SAlp Dener for (i=0; i<user->n-1; i++) { 161414d97d3SAlp Dener t1 = x[i+1] - x[i]*x[i]; 162414d97d3SAlp Dener ff += PetscSqr(1 - x[i]) + alpha*t1*t1; 163414d97d3SAlp Dener g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]); 164414d97d3SAlp Dener g[i+1] = 2*alpha*t1; 165414d97d3SAlp Dener } 166414d97d3SAlp Dener } else { 167414d97d3SAlp Dener for (i=0; i<nn; i++) { 168414d97d3SAlp Dener t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i]; 169414d97d3SAlp Dener ff += alpha*t1*t1 + t2*t2; 170414d97d3SAlp Dener g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2; 171414d97d3SAlp Dener g[2*i+1] = 2*alpha*t1; 172414d97d3SAlp Dener } 173414d97d3SAlp Dener } 174414d97d3SAlp Dener 175414d97d3SAlp Dener /* Restore vectors */ 176414d97d3SAlp Dener ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 177414d97d3SAlp Dener ierr = VecRestoreArrayWrite(G,&g);CHKERRQ(ierr); 178414d97d3SAlp Dener *f = ff; 179414d97d3SAlp Dener 180414d97d3SAlp Dener ierr = PetscLogFlops(15.0*nn);CHKERRQ(ierr); 181414d97d3SAlp Dener PetscFunctionReturn(0); 182414d97d3SAlp Dener } 183414d97d3SAlp Dener 184414d97d3SAlp Dener /* ------------------------------------------------------------------- */ 185414d97d3SAlp Dener /* 186414d97d3SAlp Dener FormHessian - Evaluates Hessian matrix. 187414d97d3SAlp Dener 188414d97d3SAlp Dener Input Parameters: 189414d97d3SAlp Dener . tao - the Tao context 190414d97d3SAlp Dener . x - input vector 191414d97d3SAlp Dener . ptr - optional user-defined context, as set by TaoSetHessian() 192414d97d3SAlp Dener 193414d97d3SAlp Dener Output Parameters: 194414d97d3SAlp Dener . H - Hessian matrix 195414d97d3SAlp Dener 196414d97d3SAlp Dener Note: Providing the Hessian may not be necessary. Only some solvers 197414d97d3SAlp Dener require this matrix. 198414d97d3SAlp Dener */ 199414d97d3SAlp Dener PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 200414d97d3SAlp Dener { 201414d97d3SAlp Dener AppCtx *user = (AppCtx*)ptr; 202414d97d3SAlp Dener PetscErrorCode ierr; 203414d97d3SAlp Dener PetscInt i, ind[2]; 204414d97d3SAlp Dener PetscReal alpha=user->alpha; 205414d97d3SAlp Dener PetscReal v[2][2]; 206414d97d3SAlp Dener const PetscScalar *x; 207414d97d3SAlp Dener PetscBool assembled; 208414d97d3SAlp Dener 209414d97d3SAlp Dener PetscFunctionBeginUser; 210414d97d3SAlp Dener /* Zero existing matrix entries */ 211414d97d3SAlp Dener ierr = MatAssembled(H,&assembled);CHKERRQ(ierr); 212414d97d3SAlp Dener if (assembled || user->chained) {ierr = MatZeroEntries(H);CHKERRQ(ierr);} 213414d97d3SAlp Dener 214414d97d3SAlp Dener /* Get a pointer to vector data */ 215414d97d3SAlp Dener ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 216414d97d3SAlp Dener 217414d97d3SAlp Dener /* Compute H(X) entries */ 218414d97d3SAlp Dener if (user->chained) { 219414d97d3SAlp Dener for (i=0; i<user->n-1; i++) { 220414d97d3SAlp Dener PetscScalar t1 = x[i+1] - x[i]*x[i]; 221414d97d3SAlp Dener v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]); 222414d97d3SAlp Dener v[0][1] = 2*alpha*(-2*x[i]); 223414d97d3SAlp Dener v[1][0] = 2*alpha*(-2*x[i]); 224414d97d3SAlp Dener v[1][1] = 2*alpha*t1; 225414d97d3SAlp Dener ind[0] = i; ind[1] = i+1; 226414d97d3SAlp Dener ierr = MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES);CHKERRQ(ierr); 227414d97d3SAlp Dener } 228414d97d3SAlp Dener } else { 229414d97d3SAlp Dener for (i=0; i<user->n/2; i++) { 230414d97d3SAlp Dener v[1][1] = 2*alpha; 231414d97d3SAlp Dener v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2; 232414d97d3SAlp Dener v[1][0] = v[0][1] = -4.0*alpha*x[2*i]; 233414d97d3SAlp Dener ind[0]=2*i; ind[1]=2*i+1; 234414d97d3SAlp Dener ierr = MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);CHKERRQ(ierr); 235414d97d3SAlp Dener } 236414d97d3SAlp Dener } 237414d97d3SAlp Dener ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 238414d97d3SAlp Dener 239414d97d3SAlp Dener /* Assemble matrix */ 240414d97d3SAlp Dener ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 241414d97d3SAlp Dener ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 242414d97d3SAlp Dener ierr = PetscLogFlops(9.0*user->n/2.0);CHKERRQ(ierr); 243414d97d3SAlp Dener PetscFunctionReturn(0); 244414d97d3SAlp Dener } 245414d97d3SAlp Dener 246414d97d3SAlp Dener /*TEST 247414d97d3SAlp Dener 248414d97d3SAlp Dener build: 249414d97d3SAlp Dener requires: !complex 250414d97d3SAlp Dener 251414d97d3SAlp Dener test: 252414d97d3SAlp Dener args: -tao_type bqnls -tao_monitor 253414d97d3SAlp Dener requires: !single 254414d97d3SAlp Dener 255414d97d3SAlp Dener TEST*/ 256