1414d97d3SAlp Dener /* Program usage: mpiexec -n 1 rosenbrock2 [-help] [all TAO options] */
2414d97d3SAlp Dener
3414d97d3SAlp Dener /* Include "petsctao.h" so we can use TAO solvers. */
4414d97d3SAlp Dener #include <petsctao.h>
5414d97d3SAlp Dener
6414d97d3SAlp Dener static char help[] = "This example demonstrates use of the TAO package to \n\
7414d97d3SAlp Dener solve an unconstrained minimization problem on a single processor. We \n\
8414d97d3SAlp Dener minimize the extended Rosenbrock function: \n\
9414d97d3SAlp Dener sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) \n\
10414d97d3SAlp Dener or the chained Rosenbrock function:\n\
11414d97d3SAlp Dener sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n";
12414d97d3SAlp Dener
13414d97d3SAlp Dener /*
14414d97d3SAlp Dener User-defined application context - contains data needed by the
15414d97d3SAlp Dener application-provided call-back routines that evaluate the function,
16414d97d3SAlp Dener gradient, and hessian.
17414d97d3SAlp Dener */
18414d97d3SAlp Dener typedef struct {
19414d97d3SAlp Dener PetscInt n; /* dimension */
20414d97d3SAlp Dener PetscReal alpha; /* condition parameter */
21414d97d3SAlp Dener PetscBool chained;
22414d97d3SAlp Dener } AppCtx;
23414d97d3SAlp Dener
24414d97d3SAlp Dener /* -------------- User-defined routines ---------- */
25414d97d3SAlp Dener PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
26414d97d3SAlp Dener PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
27414d97d3SAlp Dener
main(int argc,char ** argv)28d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
29d71ae5a4SJacob Faibussowitsch {
30414d97d3SAlp Dener PetscReal zero = 0.0;
31414d97d3SAlp Dener Vec x; /* solution vector */
32414d97d3SAlp Dener Mat H;
33414d97d3SAlp Dener Tao tao; /* Tao solver context */
34414d97d3SAlp Dener PetscBool flg, test_lmvm = PETSC_FALSE;
35414d97d3SAlp Dener PetscMPIInt size; /* number of processes running */
36414d97d3SAlp Dener AppCtx user; /* user-defined application context */
37414d97d3SAlp Dener TaoConvergedReason reason;
38414d97d3SAlp Dener PetscInt its, recycled_its = 0, oneshot_its = 0;
39414d97d3SAlp Dener
40414d97d3SAlp Dener /* 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");
45414d97d3SAlp Dener
46414d97d3SAlp Dener /* Initialize problem parameters */
479371c9d4SSatish Balay user.n = 2;
489371c9d4SSatish Balay user.alpha = 99.0;
499371c9d4SSatish Balay user.chained = PETSC_FALSE;
50414d97d3SAlp Dener /* 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));
55414d97d3SAlp Dener
56414d97d3SAlp Dener /* 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));
59414d97d3SAlp Dener
60414d97d3SAlp Dener /* The TAO code begins here */
61414d97d3SAlp Dener
62414d97d3SAlp Dener /* Create TAO solver with desired solution method */
639566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
649566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBQNLS));
65414d97d3SAlp Dener
66414d97d3SAlp Dener /* Set solution vec and an initial guess */
679566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero));
689566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
69414d97d3SAlp Dener
70414d97d3SAlp Dener /* Set routines for function, gradient, hessian evaluation */
719566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, &user));
729566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, FormHessian, &user));
73414d97d3SAlp Dener
74414d97d3SAlp Dener /* Check for TAO command line options */
759566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
76414d97d3SAlp Dener
77414d97d3SAlp Dener /* Solve the problem */
789566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0));
799566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(tao, 5));
809566063dSJacob Faibussowitsch PetscCall(TaoSetRecycleHistory(tao, PETSC_TRUE));
81414d97d3SAlp Dener reason = TAO_CONTINUE_ITERATING;
82414d97d3SAlp Dener flg = PETSC_FALSE;
839566063dSJacob Faibussowitsch PetscCall(TaoGetRecycleHistory(tao, &flg));
849566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Recycle: enabled\n"));
85414d97d3SAlp Dener while (reason != TAO_CONVERGED_GATOL) {
869566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
879566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao, &reason));
889566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao, &its));
89414d97d3SAlp Dener recycled_its += its;
909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
91414d97d3SAlp Dener }
92414d97d3SAlp Dener
93414d97d3SAlp Dener /* Disable recycling and solve again! */
949566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(tao, 100));
959566063dSJacob Faibussowitsch PetscCall(TaoSetRecycleHistory(tao, PETSC_FALSE));
969566063dSJacob Faibussowitsch PetscCall(VecSet(x, zero));
979566063dSJacob Faibussowitsch PetscCall(TaoGetRecycleHistory(tao, &flg));
989566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Recycle: disabled\n"));
999566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
1009566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao, &reason));
1013c859ba3SBarry Smith PetscCheck(reason == TAO_CONVERGED_GATOL, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
1029566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao, &oneshot_its));
1039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
10463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "recycled its: %" PetscInt_FMT " | oneshot its: %" PetscInt_FMT "\n", recycled_its, oneshot_its));
1053c859ba3SBarry Smith PetscCheck(recycled_its == oneshot_its, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Recycled solution does not match oneshot solution!");
106414d97d3SAlp Dener
1079566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H));
110414d97d3SAlp Dener
1119566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
112b122ec5aSJacob Faibussowitsch return 0;
113414d97d3SAlp Dener }
114414d97d3SAlp Dener
115414d97d3SAlp Dener /* -------------------------------------------------------------------- */
116414d97d3SAlp Dener /*
117414d97d3SAlp Dener FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
118414d97d3SAlp Dener
119414d97d3SAlp Dener Input Parameters:
120414d97d3SAlp Dener . tao - the Tao context
121414d97d3SAlp Dener . X - input vector
122414d97d3SAlp Dener . ptr - optional user-defined context, as set by TaoSetFunctionGradient()
123414d97d3SAlp Dener
124414d97d3SAlp Dener Output Parameters:
125414d97d3SAlp Dener . G - vector containing the newly evaluated gradient
126414d97d3SAlp Dener . f - function value
127414d97d3SAlp Dener
128414d97d3SAlp Dener Note:
129414d97d3SAlp Dener Some optimization methods ask for the function and the gradient evaluation
130414d97d3SAlp Dener at the same time. Evaluating both at once may be more efficient than
131414d97d3SAlp Dener evaluating each separately.
132414d97d3SAlp Dener */
FormFunctionGradient(Tao tao,Vec X,PetscReal * f,Vec G,void * ptr)133d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
134d71ae5a4SJacob Faibussowitsch {
135414d97d3SAlp Dener AppCtx *user = (AppCtx *)ptr;
136414d97d3SAlp Dener PetscInt i, nn = user->n / 2;
137414d97d3SAlp Dener PetscReal ff = 0, t1, t2, alpha = user->alpha;
138414d97d3SAlp Dener PetscScalar *g;
139414d97d3SAlp Dener const PetscScalar *x;
140414d97d3SAlp Dener
141414d97d3SAlp Dener PetscFunctionBeginUser;
142414d97d3SAlp Dener /* Get pointers to vector data */
1439566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
1449566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(G, &g));
145414d97d3SAlp Dener
146414d97d3SAlp Dener /* Compute G(X) */
147414d97d3SAlp Dener if (user->chained) {
148414d97d3SAlp Dener g[0] = 0;
149414d97d3SAlp Dener for (i = 0; i < user->n - 1; i++) {
150414d97d3SAlp Dener t1 = x[i + 1] - x[i] * x[i];
151414d97d3SAlp Dener ff += PetscSqr(1 - x[i]) + alpha * t1 * t1;
152414d97d3SAlp Dener g[i] += -2 * (1 - x[i]) + 2 * alpha * t1 * (-2 * x[i]);
153414d97d3SAlp Dener g[i + 1] = 2 * alpha * t1;
154414d97d3SAlp Dener }
155414d97d3SAlp Dener } else {
156414d97d3SAlp Dener for (i = 0; i < nn; i++) {
1579371c9d4SSatish Balay t1 = x[2 * i + 1] - x[2 * i] * x[2 * i];
1589371c9d4SSatish Balay t2 = 1 - x[2 * i];
159414d97d3SAlp Dener ff += alpha * t1 * t1 + t2 * t2;
160414d97d3SAlp Dener g[2 * i] = -4 * alpha * t1 * x[2 * i] - 2.0 * t2;
161414d97d3SAlp Dener g[2 * i + 1] = 2 * alpha * t1;
162414d97d3SAlp Dener }
163414d97d3SAlp Dener }
164414d97d3SAlp Dener
165414d97d3SAlp Dener /* Restore vectors */
1669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
1679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(G, &g));
168414d97d3SAlp Dener *f = ff;
169414d97d3SAlp Dener
1709566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(15.0 * nn));
1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
172414d97d3SAlp Dener }
173414d97d3SAlp Dener
174414d97d3SAlp Dener /* ------------------------------------------------------------------- */
175414d97d3SAlp Dener /*
176414d97d3SAlp Dener FormHessian - Evaluates Hessian matrix.
177414d97d3SAlp Dener
178414d97d3SAlp Dener Input Parameters:
179414d97d3SAlp Dener . tao - the Tao context
180414d97d3SAlp Dener . x - input vector
181414d97d3SAlp Dener . ptr - optional user-defined context, as set by TaoSetHessian()
182414d97d3SAlp Dener
183414d97d3SAlp Dener Output Parameters:
184414d97d3SAlp Dener . H - Hessian matrix
185414d97d3SAlp Dener
186414d97d3SAlp Dener Note: Providing the Hessian may not be necessary. Only some solvers
187414d97d3SAlp Dener require this matrix.
188414d97d3SAlp Dener */
FormHessian(Tao tao,Vec X,Mat H,Mat Hpre,void * ptr)189d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
190d71ae5a4SJacob Faibussowitsch {
191414d97d3SAlp Dener AppCtx *user = (AppCtx *)ptr;
192414d97d3SAlp Dener PetscInt i, ind[2];
193414d97d3SAlp Dener PetscReal alpha = user->alpha;
194414d97d3SAlp Dener PetscReal v[2][2];
195414d97d3SAlp Dener const PetscScalar *x;
196414d97d3SAlp Dener PetscBool assembled;
197414d97d3SAlp Dener
198414d97d3SAlp Dener PetscFunctionBeginUser;
199414d97d3SAlp Dener /* Zero existing matrix entries */
2009566063dSJacob Faibussowitsch PetscCall(MatAssembled(H, &assembled));
2019566063dSJacob Faibussowitsch if (assembled || user->chained) PetscCall(MatZeroEntries(H));
202414d97d3SAlp Dener
203414d97d3SAlp Dener /* Get a pointer to vector data */
2049566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
205414d97d3SAlp Dener
206414d97d3SAlp Dener /* Compute H(X) entries */
207414d97d3SAlp Dener if (user->chained) {
208414d97d3SAlp Dener for (i = 0; i < user->n - 1; i++) {
209414d97d3SAlp Dener PetscScalar t1 = x[i + 1] - x[i] * x[i];
210414d97d3SAlp Dener v[0][0] = 2 + 2 * alpha * (t1 * (-2) - 2 * x[i]);
211414d97d3SAlp Dener v[0][1] = 2 * alpha * (-2 * x[i]);
212414d97d3SAlp Dener v[1][0] = 2 * alpha * (-2 * x[i]);
213414d97d3SAlp Dener v[1][1] = 2 * alpha * t1;
2149371c9d4SSatish Balay ind[0] = i;
2159371c9d4SSatish Balay ind[1] = i + 1;
2169566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], ADD_VALUES));
217414d97d3SAlp Dener }
218414d97d3SAlp Dener } else {
219414d97d3SAlp Dener for (i = 0; i < user->n / 2; i++) {
220414d97d3SAlp Dener v[1][1] = 2 * alpha;
221414d97d3SAlp Dener v[0][0] = -4 * alpha * (x[2 * i + 1] - 3 * x[2 * i] * x[2 * i]) + 2;
222414d97d3SAlp Dener v[1][0] = v[0][1] = -4.0 * alpha * x[2 * i];
2239371c9d4SSatish Balay ind[0] = 2 * i;
2249371c9d4SSatish Balay ind[1] = 2 * i + 1;
2259566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], INSERT_VALUES));
226414d97d3SAlp Dener }
227414d97d3SAlp Dener }
2289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
229414d97d3SAlp Dener
230414d97d3SAlp Dener /* Assemble matrix */
2319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
2329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
2339566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0 * user->n / 2.0));
2343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
235414d97d3SAlp Dener }
236414d97d3SAlp Dener
237414d97d3SAlp Dener /*TEST
238414d97d3SAlp Dener
239414d97d3SAlp Dener build:
240414d97d3SAlp Dener requires: !complex
241414d97d3SAlp Dener
242414d97d3SAlp Dener test:
243414d97d3SAlp Dener args: -tao_type bqnls -tao_monitor
244414d97d3SAlp Dener requires: !single
245414d97d3SAlp Dener
246414d97d3SAlp Dener TEST*/
247