Actual source code: rosenbrock2.c

  1: /* Program usage: mpiexec -n 1 rosenbrock2 [-help] [all TAO options] */

  3: /*  Include "petsctao.h" so we can use TAO solvers.  */
  4: #include <petsctao.h>

  6: static  char help[] = "This example demonstrates use of the TAO package to \n\
  7: solve an unconstrained minimization problem on a single processor.  We \n\
  8: minimize the extended Rosenbrock function: \n\
  9:    sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) \n\
 10: or the chained Rosenbrock function:\n\
 11:    sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n";

 13: /*T
 14:    Concepts: TAO^Solving an unconstrained minimization problem
 15:    Routines: TaoCreate();
 16:    Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
 17:    Routines: TaoSetHessianRoutine();
 18:    Routines: TaoSetInitialVector();
 19:    Routines: TaoSetFromOptions();
 20:    Routines: TaoSolve();
 21:    Routines: TaoDestroy();
 22:    Processors: 1
 23: T*/

 25: /*
 26:    User-defined application context - contains data needed by the
 27:    application-provided call-back routines that evaluate the function,
 28:    gradient, and hessian.
 29: */
 30: typedef struct {
 31:   PetscInt  n;          /* dimension */
 32:   PetscReal alpha;   /* condition parameter */
 33:   PetscBool chained;
 34: } AppCtx;

 36: /* -------------- User-defined routines ---------- */
 37: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
 38: PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);

 40: int main(int argc,char **argv)
 41: {
 42:   PetscErrorCode     ierr;                  /* used to check for functions returning nonzeros */
 43:   PetscReal          zero=0.0;
 44:   Vec                x;                     /* solution vector */
 45:   Mat                H;
 46:   Tao                tao;                   /* Tao solver context */
 47:   PetscBool          flg, test_lmvm = PETSC_FALSE;
 48:   PetscMPIInt        size;                  /* number of processes running */
 49:   AppCtx             user;                  /* user-defined application context */
 50:   TaoConvergedReason reason;
 51:   PetscInt           its, recycled_its=0, oneshot_its=0;

 53:   /* Initialize TAO and PETSc */
 54:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 55:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 56:   if (size >1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors");

 58:   /* Initialize problem parameters */
 59:   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
 60:   /* Check for command line arguments to override defaults */
 61:   PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg);
 62:   PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg);
 63:   PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg);
 64:   PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg);

 66:   /* Allocate vectors for the solution and gradient */
 67:   VecCreateSeq(PETSC_COMM_SELF,user.n,&x);
 68:   MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H);

 70:   /* The TAO code begins here */

 72:   /* Create TAO solver with desired solution method */
 73:   TaoCreate(PETSC_COMM_SELF,&tao);
 74:   TaoSetType(tao,TAOLMVM);

 76:   /* Set solution vec and an initial guess */
 77:   VecSet(x, zero);
 78:   TaoSetInitialVector(tao,x);

 80:   /* Set routines for function, gradient, hessian evaluation */
 81:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);
 82:   TaoSetHessianRoutine(tao,H,H,FormHessian,&user);

 84:   /* Check for TAO command line options */
 85:   TaoSetFromOptions(tao);

 87:   /* Solve the problem */
 88:   TaoSetTolerances(tao, 1.e-5, 0.0, 0.0);
 89:   TaoSetMaximumIterations(tao, 5);
 90:   TaoLMVMRecycle(tao, PETSC_TRUE);
 91:   reason = TAO_CONTINUE_ITERATING;
 92:   while (reason != TAO_CONVERGED_GATOL) {
 93:     TaoSolve(tao);
 94:     TaoGetConvergedReason(tao, &reason);
 95:     TaoGetIterationNumber(tao, &its);
 96:     recycled_its += its;
 97:     PetscPrintf(PETSC_COMM_SELF, "-----------------------\n");
 98:   }

100:   /* Disable recycling and solve again! */
101:   TaoSetMaximumIterations(tao, 100);
102:   TaoLMVMRecycle(tao, PETSC_FALSE);
103:   VecSet(x, zero);
104:   TaoSolve(tao);
105:   TaoGetConvergedReason(tao, &reason);
106:   if (reason != TAO_CONVERGED_GATOL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
107:   TaoGetIterationNumber(tao, &oneshot_its);
108:   PetscPrintf(PETSC_COMM_SELF, "-----------------------\n");
109:   PetscPrintf(PETSC_COMM_SELF, "recycled its: %D | oneshot its: %D\n", recycled_its, oneshot_its);
110:   if (recycled_its != oneshot_its) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "LMVM recycling does not work!");

112:   TaoDestroy(&tao);
113:   VecDestroy(&x);
114:   MatDestroy(&H);

116:   PetscFinalize();
117:   return ierr;
118: }

120: /* -------------------------------------------------------------------- */
121: /*
122:     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).

124:     Input Parameters:
125: .   tao  - the Tao context
126: .   X    - input vector
127: .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()

129:     Output Parameters:
130: .   G - vector containing the newly evaluated gradient
131: .   f - function value

133:     Note:
134:     Some optimization methods ask for the function and the gradient evaluation
135:     at the same time.  Evaluating both at once may be more efficient that
136:     evaluating each separately.
137: */
138: PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
139: {
140:   AppCtx            *user = (AppCtx *) ptr;
141:   PetscInt          i,nn=user->n/2;
142:   PetscErrorCode    ierr;
143:   PetscReal         ff=0,t1,t2,alpha=user->alpha;
144:   PetscScalar       *g;
145:   const PetscScalar *x;

148:   /* Get pointers to vector data */
149:   VecGetArrayRead(X,&x);
150:   VecGetArray(G,&g);

152:   /* Compute G(X) */
153:   if (user->chained) {
154:     g[0] = 0;
155:     for (i=0; i<user->n-1; i++) {
156:       t1 = x[i+1] - x[i]*x[i];
157:       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
158:       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
159:       g[i+1] = 2*alpha*t1;
160:     }
161:   } else {
162:     for (i=0; i<nn; i++) {
163:       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
164:       ff += alpha*t1*t1 + t2*t2;
165:       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
166:       g[2*i+1] = 2*alpha*t1;
167:     }
168:   }

170:   /* Restore vectors */
171:   VecRestoreArrayRead(X,&x);
172:   VecRestoreArray(G,&g);
173:   *f   = ff;

175:   PetscLogFlops(15.0*nn);
176:   return(0);
177: }

179: /* ------------------------------------------------------------------- */
180: /*
181:    FormHessian - Evaluates Hessian matrix.

183:    Input Parameters:
184: .  tao   - the Tao context
185: .  x     - input vector
186: .  ptr   - optional user-defined context, as set by TaoSetHessian()

188:    Output Parameters:
189: .  H     - Hessian matrix

191:    Note:  Providing the Hessian may not be necessary.  Only some solvers
192:    require this matrix.
193: */
194: PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
195: {
196:   AppCtx            *user = (AppCtx*)ptr;
197:   PetscErrorCode    ierr;
198:   PetscInt          i, ind[2];
199:   PetscReal         alpha=user->alpha;
200:   PetscReal         v[2][2];
201:   const PetscScalar *x;
202:   PetscBool         assembled;

205:   /* Zero existing matrix entries */
206:   MatAssembled(H,&assembled);
207:   if (assembled) {MatZeroEntries(H);}

209:   /* Get a pointer to vector data */
210:   VecGetArrayRead(X,&x);

212:   /* Compute H(X) entries */
213:   if (user->chained) {
214:     MatZeroEntries(H);
215:     for (i=0; i<user->n-1; i++) {
216:       PetscScalar t1 = x[i+1] - x[i]*x[i];
217:       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
218:       v[0][1] = 2*alpha*(-2*x[i]);
219:       v[1][0] = 2*alpha*(-2*x[i]);
220:       v[1][1] = 2*alpha*t1;
221:       ind[0] = i; ind[1] = i+1;
222:       MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES);
223:     }
224:   } else {
225:     for (i=0; i<user->n/2; i++) {
226:       v[1][1] = 2*alpha;
227:       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
228:       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
229:       ind[0]=2*i; ind[1]=2*i+1;
230:       MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);
231:     }
232:   }
233:   VecRestoreArrayRead(X,&x);

235:   /* Assemble matrix */
236:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
237:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
238:   PetscLogFlops(9.0*user->n/2.0);
239:   return(0);
240: }

242: /*TEST

244:    build:
245:       requires: !complex

247:    test:
248:       args: -tao_type lmvm -tao_monitor
249:       requires: !single

251: TEST*/