Actual source code: rosenbrock1.c

  1: /* Program usage: mpiexec -n 1 rosenbrock1 [-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*/




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

 39: /* -------------- User-defined routines ---------- */
 40: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
 41: PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);

 43: int main(int argc,char **argv)
 44: {
 45:   PetscErrorCode     ierr;                  /* used to check for functions returning nonzeros */
 46:   PetscReal          zero=0.0;
 47:   Vec                x;                     /* solution vector */
 48:   Mat                H;
 49:   Tao                tao;                   /* Tao solver context */
 50:   PetscBool          flg, test_lmvm = PETSC_FALSE;
 51:   PetscMPIInt        size;                  /* number of processes running */
 52:   AppCtx             user;                  /* user-defined application context */
 53:   KSP                ksp;
 54:   PC                 pc;
 55:   Mat                M;
 56:   Vec                in, out, out2;
 57:   PetscReal          mult_solve_dist;

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

 64:   /* Initialize problem parameters */
 65:   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
 66:   /* Check for command line arguments to override defaults */
 67:   PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg);
 68:   PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg);
 69:   PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg);
 70:   PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg);

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

 76:   /* The TAO code begins here */

 78:   /* Create TAO solver with desired solution method */
 79:   TaoCreate(PETSC_COMM_SELF,&tao);
 80:   TaoSetType(tao,TAOLMVM);

 82:   /* Set solution vec and an initial guess */
 83:   VecSet(x, zero);
 84:   TaoSetInitialVector(tao,x);

 86:   /* Set routines for function, gradient, hessian evaluation */
 87:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);
 88:   TaoSetHessianRoutine(tao,H,H,FormHessian,&user);

 90:   /* Test the LMVM matrix */
 91:   if (test_lmvm) {
 92:     PetscOptionsSetValue(NULL, "-tao_type", "bqnktr");
 93:   }

 95:   /* Check for TAO command line options */
 96:   TaoSetFromOptions(tao);

 98:   /* SOLVE THE APPLICATION */
 99:   TaoSolve(tao);

101:   /* Test the LMVM matrix */
102:   if (test_lmvm) {
103:     TaoGetKSP(tao, &ksp);
104:     KSPGetPC(ksp, &pc);
105:     PCLMVMGetMatLMVM(pc, &M);
106:     VecDuplicate(x, &in);
107:     VecDuplicate(x, &out);
108:     VecDuplicate(x, &out2);
109:     VecSet(in, 1.0);
110:     MatMult(M, in, out);
111:     MatSolve(M, out, out2);
112:     VecAXPY(out2, -1.0, in);
113:     VecNorm(out2, NORM_2, &mult_solve_dist);
114:     if (mult_solve_dist < 1.e-11) {
115:       PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-11\n");
116:     } else if (mult_solve_dist < 1.e-6) {
117:       PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-6\n");
118:     } else {
119:       PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: %e\n", (double)mult_solve_dist);
120:     }
121:     VecDestroy(&in);
122:     VecDestroy(&out);
123:     VecDestroy(&out2);
124:   }

126:   TaoDestroy(&tao);
127:   VecDestroy(&x);
128:   MatDestroy(&H);

130:   PetscFinalize();
131:   return ierr;
132: }

134: /* -------------------------------------------------------------------- */
135: /*
136:     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).

138:     Input Parameters:
139: .   tao  - the Tao context
140: .   X    - input vector
141: .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()

143:     Output Parameters:
144: .   G - vector containing the newly evaluated gradient
145: .   f - function value

147:     Note:
148:     Some optimization methods ask for the function and the gradient evaluation
149:     at the same time.  Evaluating both at once may be more efficient that
150:     evaluating each separately.
151: */
152: PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
153: {
154:   AppCtx            *user = (AppCtx *) ptr;
155:   PetscInt          i,nn=user->n/2;
156:   PetscErrorCode    ierr;
157:   PetscReal         ff=0,t1,t2,alpha=user->alpha;
158:   PetscScalar       *g;
159:   const PetscScalar *x;

162:   /* Get pointers to vector data */
163:   VecGetArrayRead(X,&x);
164:   VecGetArray(G,&g);

166:   /* Compute G(X) */
167:   if (user->chained) {
168:     g[0] = 0;
169:     for (i=0; i<user->n-1; i++) {
170:       t1 = x[i+1] - x[i]*x[i];
171:       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
172:       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
173:       g[i+1] = 2*alpha*t1;
174:     }
175:   } else {
176:     for (i=0; i<nn; i++){
177:       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
178:       ff += alpha*t1*t1 + t2*t2;
179:       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
180:       g[2*i+1] = 2*alpha*t1;
181:     }
182:   }

184:   /* Restore vectors */
185:   VecRestoreArrayRead(X,&x);
186:   VecRestoreArray(G,&g);
187:   *f   = ff;

189:   PetscLogFlops(15.0*nn);
190:   return(0);
191: }

193: /* ------------------------------------------------------------------- */
194: /*
195:    FormHessian - Evaluates Hessian matrix.

197:    Input Parameters:
198: .  tao   - the Tao context
199: .  x     - input vector
200: .  ptr   - optional user-defined context, as set by TaoSetHessian()

202:    Output Parameters:
203: .  H     - Hessian matrix

205:    Note:  Providing the Hessian may not be necessary.  Only some solvers
206:    require this matrix.
207: */
208: PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
209: {
210:   AppCtx            *user = (AppCtx*)ptr;
211:   PetscErrorCode    ierr;
212:   PetscInt          i, ind[2];
213:   PetscReal         alpha=user->alpha;
214:   PetscReal         v[2][2];
215:   const PetscScalar *x;
216:   PetscBool         assembled;

219:   /* Zero existing matrix entries */
220:   MatAssembled(H,&assembled);
221:   if (assembled){MatZeroEntries(H);}

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

226:   /* Compute H(X) entries */
227:   if (user->chained) {
228:     MatZeroEntries(H);
229:     for (i=0; i<user->n-1; i++) {
230:       PetscScalar t1 = x[i+1] - x[i]*x[i];
231:       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
232:       v[0][1] = 2*alpha*(-2*x[i]);
233:       v[1][0] = 2*alpha*(-2*x[i]);
234:       v[1][1] = 2*alpha*t1;
235:       ind[0] = i; ind[1] = i+1;
236:       MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES);
237:     }
238:   } else {
239:     for (i=0; i<user->n/2; i++){
240:       v[1][1] = 2*alpha;
241:       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
242:       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
243:       ind[0]=2*i; ind[1]=2*i+1;
244:       MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);
245:     }
246:   }
247:   VecRestoreArrayRead(X,&x);

249:   /* Assemble matrix */
250:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
251:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
252:   PetscLogFlops(9.0*user->n/2.0);
253:   return(0);
254: }


257: /*TEST

259:    build:
260:       requires: !complex

262:    test:
263:       args: -tao_smonitor -tao_type nls -tao_gatol 1.e-4
264:       requires: !single

266:    test:
267:       suffix: 2
268:       args: -tao_smonitor -tao_type lmvm -tao_gatol 1.e-3

270:    test:
271:       suffix: 3
272:       args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4
273:       requires: !single

275:    test:
276:       suffix: 4
277:       args: -tao_smonitor -tao_type ntr -tao_mf_hessian -tao_ntr_pc_type none -tao_gatol 1.e-4

279:    test:
280:       suffix: 5
281:       args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4

283:    test:
284:       suffix: 6
285:       args: -tao_smonitor -tao_type bntl -tao_gatol 1.e-4

287:    test:
288:       suffix: 7
289:       args: -tao_smonitor -tao_type bnls -tao_gatol 1.e-4

291:    test:
292:       suffix: 8
293:       args: -tao_smonitor -tao_type bntr -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4

295:    test:
296:       suffix: 9
297:       args: -tao_smonitor -tao_type bntl -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4

299:    test:
300:       suffix: 10
301:       args: -tao_smonitor -tao_type bnls -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4

303:    test:
304:       suffix: 11
305:       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbroyden

307:    test:
308:       suffix: 12
309:       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbadbroyden

311:    test:
312:      suffix: 13
313:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbroyden

315:    test:
316:      suffix: 14
317:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbfgs

319:    test:
320:      suffix: 15
321:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmdfp

323:    test:
324:      suffix: 16
325:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsr1

327:    test:
328:      suffix: 17
329:      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnls

331:    test:
332:      suffix: 18
333:      args: -tao_smonitor -tao_gatol 1e-4 -tao_type blmvm

335:    test:
336:      suffix: 19
337:      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1

339:    test:
340:      suffix: 20
341:      args: -tao_monitor -tao_gatol 1e-4 -tao_type blmvm -tao_ls_monitor

343:    test:
344:      suffix: 21
345:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbadbroyden

347:    test:
348:      suffix: 22
349:      args: -tao_max_it 1 -tao_converged_reason

351:    test:
352:      suffix: 23
353:      args: -tao_max_funcs 0 -tao_converged_reason

355:    test:
356:      suffix: 24
357:      args: -tao_gatol 10 -tao_converged_reason

359:    test:
360:      suffix: 25
361:      args: -tao_grtol 10 -tao_converged_reason

363:    test:
364:      suffix: 26
365:      args: -tao_gttol 10 -tao_converged_reason

367:    test:
368:      suffix: 27
369:      args: -tao_steptol 10 -tao_converged_reason

371:    test:
372:      suffix: 28
373:      args: -tao_fmin 10 -tao_converged_reason

375: TEST*/