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