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