Actual source code: rosenbrock1.c
1: /* Program usage: mpirun -np 1 rosenbrock1 [-help] [all TAO options] */
3: /* Include "tao.h" so we can use TAO solvers. */
4: #include "tao.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: TaoInitialize(); TaoFinalize();
14: Routines: TaoCreate();
15: Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
16: Routines: TaoSetHessianRoutine();
17: Routines: TaoSetInitialVector();
18: Routines: TaoSetFromOptions();
19: Routines: TaoSolve();
20: Routines: TaoGetTerminationReason(); TaoDestroy();
21: Processors: 1
22: 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: } AppCtx;
35: /* -------------- User-defined routines ---------- */
36: PetscErrorCode FormFunctionGradient(TaoSolver,Vec,PetscReal*,Vec,void*);
37: PetscErrorCode FormHessian(TaoSolver,Vec,Mat*,Mat*,MatStructure*,void*);
41: int main(int argc,char **argv)
42: {
43: PetscErrorCode ierr; /* used to check for functions returning nonzeros */
44: PetscReal zero=0.0;
45: Vec x; /* solution vector */
46: Mat H;
47: TaoSolver tao; /* TaoSolver solver context */
48: PetscBool flg;
49: int size,rank; /* number of processes running */
50: TaoSolverTerminationReason reason;
51: AppCtx user; /* user-defined application context */
53: /* Initialize TAO and PETSc */
54: TaoInitialize(&argc,&argv,(char*)0,help);
55: MPI_Comm_size(PETSC_COMM_WORLD,&size);
56: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
58: if (size >1) {
59: if (rank == 0) {
60: PetscPrintf(PETSC_COMM_SELF,"This example is intended for single processor use!\n");
61: SETERRQ(PETSC_COMM_SELF,1,"Incorrect number of processors");
62: }
63: }
66: /* Initialize problem parameters */
67: user.n = 2; user.alpha = 99.0;
68: /* Check for command line arguments to override defaults */
69: PetscOptionsGetInt(PETSC_NULL,"-n",&user.n,&flg);
70: PetscOptionsGetReal(PETSC_NULL,"-alpha",&user.alpha,&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,PETSC_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,"tao_lmvm");
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);
89:
90: /* Check for TAO command line options */
91: TaoSetFromOptions(tao);
93: /* SOLVE THE APPLICATION */
94: TaoSolve(tao);
96: /* Get termination information */
97: TaoGetTerminationReason(tao,&reason);
98: if (reason <= 0)
99: PetscPrintf(MPI_COMM_WORLD,"Try a different TAO type, adjust some parameters, or check the function evaluation routines\n");
102: /* Free TAO data structures */
103: TaoDestroy(&tao);
105: /* Free PETSc data structures */
106: VecDestroy(&x);
107: MatDestroy(&H);
109: TaoFinalize();
111: return 0;
112: }
114: /* -------------------------------------------------------------------- */
117: /*
118: FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
120: Input Parameters:
121: . tao - the TaoSolver context
122: . X - input vector
123: . ptr - optional user-defined context, as set by TaoSetFunctionGradient()
124:
125: Output Parameters:
126: . G - vector containing the newly evaluated gradient
127: . f - function value
129: Note:
130: Some optimization methods ask for the function and the gradient evaluation
131: at the same time. Evaluating both at once may be more efficient that
132: evaluating each separately.
133: */
134: PetscErrorCode FormFunctionGradient(TaoSolver tao,Vec X,PetscReal *f, Vec G,void *ptr)
135: {
136: AppCtx *user = (AppCtx *) ptr;
137: PetscInt i,nn=user->n/2;
139: PetscReal ff=0,t1,t2,alpha=user->alpha;
140: PetscReal *x,*g;
142: /* Get pointers to vector data */
143: VecGetArray(X,&x);
144: VecGetArray(G,&g);
146: /* Compute G(X) */
147: for (i=0; i<nn; i++){
148: t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
149: ff += alpha*t1*t1 + t2*t2;
150: g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
151: g[2*i+1] = 2*alpha*t1;
152: }
154: /* Restore vectors */
155: VecRestoreArray(X,&x);
156: VecRestoreArray(G,&g);
157: *f=ff;
159: PetscLogFlops(nn*15);
160: return 0;
161: }
163: /* ------------------------------------------------------------------- */
166: /*
167: FormHessian - Evaluates Hessian matrix.
169: Input Parameters:
170: . tao - the TaoSolver context
171: . x - input vector
172: . ptr - optional user-defined context, as set by TaoSetHessian()
174: Output Parameters:
175: . H - Hessian matrix
177: Note: Providing the Hessian may not be necessary. Only some solvers
178: require this matrix.
179: */
180: PetscErrorCode FormHessian(TaoSolver tao,Vec X,Mat *HH, Mat *Hpre, MatStructure *flag,void *ptr)
181: {
182: AppCtx *user = (AppCtx*)ptr;
184: PetscInt i, nn=user->n/2, ind[2];
185: PetscReal alpha=user->alpha;
186: PetscReal v[2][2],*x;
187: Mat H=*HH;
188: PetscBool assembled;
190: /* Zero existing matrix entries */
191: MatAssembled(H,&assembled);
192: if (assembled){MatZeroEntries(H); }
195: /* Get a pointer to vector data */
196: VecGetArray(X,&x);
198: /* Compute H(X) entries */
199: for (i=0; i<user->n/2; i++){
200: v[1][1] = 2*alpha;
201: v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
202: v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
203: ind[0]=2*i; ind[1]=2*i+1;
204: MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);
205: }
206: VecRestoreArray(X,&x);
208: /* Assemble matrix */
209: MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
210: MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
211: *flag=SAME_NONZERO_PATTERN;
212:
213: PetscLogFlops(nn*9);
214: return 0;
215: }