1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 rosenbrock2 [-help] [all TAO options] */
2c4762a1bSJed Brown
3c4762a1bSJed Brown /* Include "petsctao.h" so we can use TAO solvers. */
4c4762a1bSJed Brown #include <petsctao.h>
5c4762a1bSJed Brown
6c4762a1bSJed Brown static char help[] = "This example demonstrates use of the TAO package to \n\
7c4762a1bSJed Brown solve an unconstrained minimization problem on a single processor. We \n\
8c4762a1bSJed Brown minimize the extended Rosenbrock function: \n\
9c4762a1bSJed Brown sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) \n\
10c4762a1bSJed Brown or the chained Rosenbrock function:\n\
11c4762a1bSJed Brown sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n";
12c4762a1bSJed Brown
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown User-defined application context - contains data needed by the
15c4762a1bSJed Brown application-provided call-back routines that evaluate the function,
16c4762a1bSJed Brown gradient, and hessian.
17c4762a1bSJed Brown */
18c4762a1bSJed Brown typedef struct {
19c4762a1bSJed Brown PetscInt n; /* dimension */
20c4762a1bSJed Brown PetscReal alpha; /* condition parameter */
21c4762a1bSJed Brown PetscBool chained;
22c4762a1bSJed Brown } AppCtx;
23c4762a1bSJed Brown
24c4762a1bSJed Brown /* -------------- User-defined routines ---------- */
25c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
26c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
27c4762a1bSJed Brown
main(int argc,char ** argv)28d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
29d71ae5a4SJacob Faibussowitsch {
30c4762a1bSJed Brown PetscReal zero = 0.0;
31c4762a1bSJed Brown Vec x; /* solution vector */
32c4762a1bSJed Brown Mat H;
33c4762a1bSJed Brown Tao tao; /* Tao solver context */
34c4762a1bSJed Brown PetscBool flg, test_lmvm = PETSC_FALSE;
35c4762a1bSJed Brown PetscMPIInt size; /* number of processes running */
36c4762a1bSJed Brown AppCtx user; /* user-defined application context */
37c4762a1bSJed Brown TaoConvergedReason reason;
38c4762a1bSJed Brown PetscInt its, recycled_its = 0, oneshot_its = 0;
39c4762a1bSJed Brown
40c4762a1bSJed Brown /* Initialize TAO and PETSc */
41327415f7SBarry Smith PetscFunctionBeginUser;
42*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
443c859ba3SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors");
45c4762a1bSJed Brown
46c4762a1bSJed Brown /* Initialize problem parameters */
479371c9d4SSatish Balay user.n = 2;
489371c9d4SSatish Balay user.alpha = 99.0;
499371c9d4SSatish Balay user.chained = PETSC_FALSE;
50c4762a1bSJed Brown /* Check for command line arguments to override defaults */
519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &user.n, &flg));
529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-alpha", &user.alpha, &flg));
539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-chained", &user.chained, &flg));
549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg));
55c4762a1bSJed Brown
56c4762a1bSJed Brown /* Allocate vectors for the solution and gradient */
579566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.n, &x));
589566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, 2, user.n, user.n, 1, NULL, &H));
59c4762a1bSJed Brown
60c4762a1bSJed Brown /* The TAO code begins here */
61c4762a1bSJed Brown
62c4762a1bSJed Brown /* Create TAO solver with desired solution method */
639566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
649566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLMVM));
65c4762a1bSJed Brown
66c4762a1bSJed Brown /* Set solution vec and an initial guess */
679566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero));
689566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
69c4762a1bSJed Brown
70c4762a1bSJed Brown /* Set routines for function, gradient, hessian evaluation */
719566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, &user));
729566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, FormHessian, &user));
73c4762a1bSJed Brown
74c4762a1bSJed Brown /* Check for TAO command line options */
759566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
76c4762a1bSJed Brown
77c4762a1bSJed Brown /* Solve the problem */
789566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0));
799566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(tao, 5));
809566063dSJacob Faibussowitsch PetscCall(TaoLMVMRecycle(tao, PETSC_TRUE));
81c4762a1bSJed Brown reason = TAO_CONTINUE_ITERATING;
82c4762a1bSJed Brown while (reason != TAO_CONVERGED_GATOL) {
839566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
849566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao, &reason));
859566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao, &its));
86c4762a1bSJed Brown recycled_its += its;
879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
88c4762a1bSJed Brown }
89c4762a1bSJed Brown
90c4762a1bSJed Brown /* Disable recycling and solve again! */
919566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(tao, 100));
929566063dSJacob Faibussowitsch PetscCall(TaoLMVMRecycle(tao, PETSC_FALSE));
939566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero));
949566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
959566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao, &reason));
963c859ba3SBarry Smith PetscCheck(reason == TAO_CONVERGED_GATOL, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
979566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao, &oneshot_its));
989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
9963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "recycled its: %" PetscInt_FMT " | oneshot its: %" PetscInt_FMT "\n", recycled_its, oneshot_its));
1003c859ba3SBarry Smith PetscCheck(recycled_its == oneshot_its, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "LMVM recycling does not work!");
101c4762a1bSJed Brown
1029566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
1039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H));
105c4762a1bSJed Brown
1069566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
107b122ec5aSJacob Faibussowitsch return 0;
108c4762a1bSJed Brown }
109c4762a1bSJed Brown
110c4762a1bSJed Brown /* -------------------------------------------------------------------- */
111c4762a1bSJed Brown /*
112c4762a1bSJed Brown FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
113c4762a1bSJed Brown
114c4762a1bSJed Brown Input Parameters:
115c4762a1bSJed Brown . tao - the Tao context
116c4762a1bSJed Brown . X - input vector
117c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunctionGradient()
118c4762a1bSJed Brown
119c4762a1bSJed Brown Output Parameters:
120c4762a1bSJed Brown . G - vector containing the newly evaluated gradient
121c4762a1bSJed Brown . f - function value
122c4762a1bSJed Brown
123c4762a1bSJed Brown Note:
124c4762a1bSJed Brown Some optimization methods ask for the function and the gradient evaluation
125c4762a1bSJed Brown at the same time. Evaluating both at once may be more efficient that
126c4762a1bSJed Brown evaluating each separately.
127c4762a1bSJed Brown */
FormFunctionGradient(Tao tao,Vec X,PetscReal * f,Vec G,void * ptr)128d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
129d71ae5a4SJacob Faibussowitsch {
130c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
131c4762a1bSJed Brown PetscInt i, nn = user->n / 2;
132c4762a1bSJed Brown PetscReal ff = 0, t1, t2, alpha = user->alpha;
133c4762a1bSJed Brown PetscScalar *g;
134c4762a1bSJed Brown const PetscScalar *x;
135c4762a1bSJed Brown
136c4762a1bSJed Brown PetscFunctionBeginUser;
137c4762a1bSJed Brown /* Get pointers to vector data */
1389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
1399566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g));
140c4762a1bSJed Brown
141c4762a1bSJed Brown /* Compute G(X) */
142c4762a1bSJed Brown if (user->chained) {
143c4762a1bSJed Brown g[0] = 0;
144c4762a1bSJed Brown for (i = 0; i < user->n - 1; i++) {
145c4762a1bSJed Brown t1 = x[i + 1] - x[i] * x[i];
146c4762a1bSJed Brown ff += PetscSqr(1 - x[i]) + alpha * t1 * t1;
147c4762a1bSJed Brown g[i] += -2 * (1 - x[i]) + 2 * alpha * t1 * (-2 * x[i]);
148c4762a1bSJed Brown g[i + 1] = 2 * alpha * t1;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown } else {
151c4762a1bSJed Brown for (i = 0; i < nn; i++) {
1529371c9d4SSatish Balay t1 = x[2 * i + 1] - x[2 * i] * x[2 * i];
1539371c9d4SSatish Balay t2 = 1 - x[2 * i];
154c4762a1bSJed Brown ff += alpha * t1 * t1 + t2 * t2;
155c4762a1bSJed Brown g[2 * i] = -4 * alpha * t1 * x[2 * i] - 2.0 * t2;
156c4762a1bSJed Brown g[2 * i + 1] = 2 * alpha * t1;
157c4762a1bSJed Brown }
158c4762a1bSJed Brown }
159c4762a1bSJed Brown
160c4762a1bSJed Brown /* Restore vectors */
1619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
1629566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g));
163c4762a1bSJed Brown *f = ff;
164c4762a1bSJed Brown
1659566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(15.0 * nn));
1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown
169c4762a1bSJed Brown /* ------------------------------------------------------------------- */
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix.
172c4762a1bSJed Brown
173c4762a1bSJed Brown Input Parameters:
174c4762a1bSJed Brown . tao - the Tao context
175c4762a1bSJed Brown . x - input vector
176c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian()
177c4762a1bSJed Brown
178c4762a1bSJed Brown Output Parameters:
179c4762a1bSJed Brown . H - Hessian matrix
180c4762a1bSJed Brown
181c4762a1bSJed Brown Note: Providing the Hessian may not be necessary. Only some solvers
182c4762a1bSJed Brown require this matrix.
183c4762a1bSJed Brown */
FormHessian(Tao tao,Vec X,Mat H,Mat Hpre,void * ptr)184d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
185d71ae5a4SJacob Faibussowitsch {
186c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
187c4762a1bSJed Brown PetscInt i, ind[2];
188c4762a1bSJed Brown PetscReal alpha = user->alpha;
189c4762a1bSJed Brown PetscReal v[2][2];
190c4762a1bSJed Brown const PetscScalar *x;
191c4762a1bSJed Brown PetscBool assembled;
192c4762a1bSJed Brown
193c4762a1bSJed Brown PetscFunctionBeginUser;
194c4762a1bSJed Brown /* Zero existing matrix entries */
1959566063dSJacob Faibussowitsch PetscCall(MatAssembled(H, &assembled));
1969566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(H));
197c4762a1bSJed Brown
198c4762a1bSJed Brown /* Get a pointer to vector data */
1999566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
200c4762a1bSJed Brown
201c4762a1bSJed Brown /* Compute H(X) entries */
202c4762a1bSJed Brown if (user->chained) {
2039566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H));
204c4762a1bSJed Brown for (i = 0; i < user->n - 1; i++) {
205c4762a1bSJed Brown PetscScalar t1 = x[i + 1] - x[i] * x[i];
206c4762a1bSJed Brown v[0][0] = 2 + 2 * alpha * (t1 * (-2) - 2 * x[i]);
207c4762a1bSJed Brown v[0][1] = 2 * alpha * (-2 * x[i]);
208c4762a1bSJed Brown v[1][0] = 2 * alpha * (-2 * x[i]);
209c4762a1bSJed Brown v[1][1] = 2 * alpha * t1;
2109371c9d4SSatish Balay ind[0] = i;
2119371c9d4SSatish Balay ind[1] = i + 1;
2129566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], ADD_VALUES));
213c4762a1bSJed Brown }
214c4762a1bSJed Brown } else {
215c4762a1bSJed Brown for (i = 0; i < user->n / 2; i++) {
216c4762a1bSJed Brown v[1][1] = 2 * alpha;
217c4762a1bSJed Brown v[0][0] = -4 * alpha * (x[2 * i + 1] - 3 * x[2 * i] * x[2 * i]) + 2;
218c4762a1bSJed Brown v[1][0] = v[0][1] = -4.0 * alpha * x[2 * i];
2199371c9d4SSatish Balay ind[0] = 2 * i;
2209371c9d4SSatish Balay ind[1] = 2 * i + 1;
2219566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], INSERT_VALUES));
222c4762a1bSJed Brown }
223c4762a1bSJed Brown }
2249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
225c4762a1bSJed Brown
226c4762a1bSJed Brown /* Assemble matrix */
2279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
2289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
2299566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0 * user->n / 2.0));
2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
231c4762a1bSJed Brown }
232c4762a1bSJed Brown
233c4762a1bSJed Brown /*TEST
234c4762a1bSJed Brown
235c4762a1bSJed Brown build:
236c4762a1bSJed Brown requires: !complex
237c4762a1bSJed Brown
238c4762a1bSJed Brown test:
239c4762a1bSJed Brown args: -tao_type lmvm -tao_monitor
240c4762a1bSJed Brown requires: !single
241c4762a1bSJed Brown
242c4762a1bSJed Brown TEST*/
243