Actual source code: rosenbrock1.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: /* Program usage: mpirun -np 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";

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


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

 34: /* -------------- User-defined routines ---------- */
 35: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
 36: 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;
 48:   PetscMPIInt        size,rank;                  /* number of processes running */
 49:   TaoConvergedReason reason;
 50:   AppCtx             user;                  /* user-defined application context */

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

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

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

 68:   /* The TAO code begins here */

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

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

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

 82:   /* Check for TAO command line options */
 83:   TaoSetFromOptions(tao);

 85:   /* SOLVE THE APPLICATION */
 86:   TaoSolve(tao);

 88:   /* Get termination information */
 89:   TaoGetConvergedReason(tao,&reason);
 90:   if (reason <= 0) {
 91:     PetscPrintf(MPI_COMM_WORLD,"Try a different TAO type, adjust some parameters, or check the function evaluation routines\n");
 92:   }

 94:   TaoDestroy(&tao);
 95:   VecDestroy(&x);
 96:   MatDestroy(&H);

 98:   return 0;
 99: }

101: /* -------------------------------------------------------------------- */
104: /*
105:     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).

107:     Input Parameters:
108: .   tao  - the Tao context
109: .   X    - input vector
110: .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()

112:     Output Parameters:
113: .   G - vector containing the newly evaluated gradient
114: .   f - function value

116:     Note:
117:     Some optimization methods ask for the function and the gradient evaluation
118:     at the same time.  Evaluating both at once may be more efficient that
119:     evaluating each separately.
120: */
121: PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
122: {
123:   AppCtx         *user = (AppCtx *) ptr;
124:   PetscInt       i,nn=user->n/2;
126:   PetscReal      ff=0,t1,t2,alpha=user->alpha;
127:   PetscReal      *x,*g;

129:   /* Get pointers to vector data */
130:   VecGetArray(X,&x);
131:   VecGetArray(G,&g);

133:   /* Compute G(X) */
134:   for (i=0; i<nn; i++){
135:     t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
136:     ff += alpha*t1*t1 + t2*t2;
137:     g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
138:     g[2*i+1] = 2*alpha*t1;
139:   }

141:   /* Restore vectors */
142:   VecRestoreArray(X,&x);
143:   VecRestoreArray(G,&g);
144:   *f=ff;

146:   PetscLogFlops(nn*15);
147:   return 0;
148: }

150: /* ------------------------------------------------------------------- */
153: /*
154:    FormHessian - Evaluates Hessian matrix.

156:    Input Parameters:
157: .  tao   - the Tao context
158: .  x     - input vector
159: .  ptr   - optional user-defined context, as set by TaoSetHessian()

161:    Output Parameters:
162: .  H     - Hessian matrix

164:    Note:  Providing the Hessian may not be necessary.  Only some solvers
165:    require this matrix.
166: */
167: PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
168: {
169:   AppCtx         *user = (AppCtx*)ptr;
171:   PetscInt       i, nn=user->n/2, ind[2];
172:   PetscReal      alpha=user->alpha;
173:   PetscReal      v[2][2],*x;
174:   PetscBool      assembled;

176:   /* Zero existing matrix entries */
177:   MatAssembled(H,&assembled);
178:   if (assembled){MatZeroEntries(H); }

180:   /* Get a pointer to vector data */
181:   VecGetArray(X,&x);

183:   /* Compute H(X) entries */
184:   for (i=0; i<user->n/2; i++){
185:     v[1][1] = 2*alpha;
186:     v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
187:     v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
188:     ind[0]=2*i; ind[1]=2*i+1;
189:     MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);
190:   }
191:   VecRestoreArray(X,&x);

193:   /* Assemble matrix */
194:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
195:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
196:   PetscLogFlops(nn*9);
197:   return 0;
198: }