Actual source code: rosenbrock1.c

petsc-3.7.5 2017-01-01
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";

 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: 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:   AppCtx             user;                  /* user-defined application context */

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

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

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

 67:   /* The TAO code begins here */

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

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

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

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

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

 87:   TaoDestroy(&tao);
 88:   VecDestroy(&x);
 89:   MatDestroy(&H);

 91:   PetscFinalize();
 92:   return 0;
 93: }

 95: /* -------------------------------------------------------------------- */
 98: /*
 99:     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).

101:     Input Parameters:
102: .   tao  - the Tao context
103: .   X    - input vector
104: .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()

106:     Output Parameters:
107: .   G - vector containing the newly evaluated gradient
108: .   f - function value

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

123:   /* Get pointers to vector data */
124:   VecGetArray(X,&x);
125:   VecGetArray(G,&g);

127:   /* Compute G(X) */
128:   for (i=0; i<nn; i++){
129:     t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
130:     ff += alpha*t1*t1 + t2*t2;
131:     g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
132:     g[2*i+1] = 2*alpha*t1;
133:   }

135:   /* Restore vectors */
136:   VecRestoreArray(X,&x);
137:   VecRestoreArray(G,&g);
138:   *f=ff;

140:   PetscLogFlops(nn*15);
141:   return 0;
142: }

144: /* ------------------------------------------------------------------- */
147: /*
148:    FormHessian - Evaluates Hessian matrix.

150:    Input Parameters:
151: .  tao   - the Tao context
152: .  x     - input vector
153: .  ptr   - optional user-defined context, as set by TaoSetHessian()

155:    Output Parameters:
156: .  H     - Hessian matrix

158:    Note:  Providing the Hessian may not be necessary.  Only some solvers
159:    require this matrix.
160: */
161: PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
162: {
163:   AppCtx         *user = (AppCtx*)ptr;
165:   PetscInt       i, ind[2];
166:   PetscReal      alpha=user->alpha;
167:   PetscReal      v[2][2],*x;
168:   PetscBool      assembled;

170:   /* Zero existing matrix entries */
171:   MatAssembled(H,&assembled);
172:   if (assembled){MatZeroEntries(H); }

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

177:   /* Compute H(X) entries */
178:   for (i=0; i<user->n/2; i++){
179:     v[1][1] = 2*alpha;
180:     v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
181:     v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
182:     ind[0]=2*i; ind[1]=2*i+1;
183:     MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);
184:   }
185:   VecRestoreArray(X,&x);

187:   /* Assemble matrix */
188:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
189:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
190:   PetscLogFlops(9.0*user->n/2.0);
191:   return 0;
192: }