Actual source code: rosenbrock1.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  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,rank;                  /* 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:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 63:   if (size >1) SETERRQ(PETSC_COMM_SELF,1,"Incorrect number of processors");

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

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

 77:   /* The TAO code begins here */

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

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

 87:   /* Set routines for function, gradient, hessian evaluation */
 88:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);
 89:   TaoSetHessianRoutine(tao,H,H,FormHessian,&user);
 90: 
 91:   /* Test the LMVM matrix */
 92:   if (test_lmvm) {
 93:     PetscOptionsSetValue(NULL, "-tao_type", "bqnktr");
 94:   }

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

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

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

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

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

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

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

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

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

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

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

190:   PetscLogFlops(nn*15);
191:   return(0);
192: }

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

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

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

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

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

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

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

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


258: /*TEST

260:    build:
261:       requires: !complex

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

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

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

276:    test:
277:       suffix: 4
278:       args: -tao_smonitor -tao_type ntr -tao_mf_hessian -tao_ntr_pc_type none -tao_gatol 1.e-4
279:       
280:    test:
281:       suffix: 5
282:       args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4
283:       
284:    test:
285:       suffix: 6
286:       args: -tao_smonitor -tao_type bntl -tao_gatol 1.e-4
287:    
288:    test:
289:       suffix: 7
290:       args: -tao_smonitor -tao_type bnls -tao_gatol 1.e-4
291:    
292:    test:
293:       suffix: 8
294:       args: -tao_smonitor -tao_type bntr -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
295:    
296:    test:
297:       suffix: 9
298:       args: -tao_smonitor -tao_type bntl -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
299:    
300:    test:
301:       suffix: 10
302:       args: -tao_smonitor -tao_type bnls -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
303:       
304:    test:
305:       suffix: 11
306:       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbrdn
307:       
308:    test:
309:       suffix: 12
310:       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbadbrdn
311:       
312:    test:
313:      suffix: 13
314:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbrdn

316:    test:
317:      suffix: 14
318:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbfgs
319:      
320:    test:
321:      suffix: 15
322:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmdfp
323:      
324:    test:
325:      suffix: 16
326:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsr1
327:      
328:    test:
329:      suffix: 17
330:      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnls
331:      
332:    test:
333:      suffix: 18
334:      args: -tao_smonitor -tao_gatol 1e-4 -tao_type blmvm
335:      
336:    test:
337:      suffix: 19
338:      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
339:      
340:    test:
341:      suffix: 20
342:      args: -tao_monitor -tao_gatol 1e-4 -tao_type blmvm -tao_ls_monitor
343:      
344:    test:
345:      suffix: 21
346:      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbadbrdn

348: TEST*/