xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock3.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
28414d97d3SAlp Dener int main(int argc,char **argv)
29414d97d3SAlp Dener {
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 */
41*327415f7SBarry Smith   PetscFunctionBeginUser;
429566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,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 */
47414d97d3SAlp Dener   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
48414d97d3SAlp Dener   /* Check for command line arguments to override defaults */
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg));
509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg));
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg));
53414d97d3SAlp Dener 
54414d97d3SAlp Dener   /* Allocate vectors for the solution and gradient */
559566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.n,&x));
569566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H));
57414d97d3SAlp Dener 
58414d97d3SAlp Dener   /* The TAO code begins here */
59414d97d3SAlp Dener 
60414d97d3SAlp Dener   /* Create TAO solver with desired solution method */
619566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF,&tao));
629566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBQNLS));
63414d97d3SAlp Dener 
64414d97d3SAlp Dener   /* Set solution vec and an initial guess */
659566063dSJacob Faibussowitsch   PetscCall(VecSet(x, zero));
669566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,x));
67414d97d3SAlp Dener 
68414d97d3SAlp Dener   /* Set routines for function, gradient, hessian evaluation */
699566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user));
709566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao,H,H,FormHessian,&user));
71414d97d3SAlp Dener 
72414d97d3SAlp Dener   /* Check for TAO command line options */
739566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
74414d97d3SAlp Dener 
75414d97d3SAlp Dener   /* Solve the problem */
769566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0));
779566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(tao, 5));
789566063dSJacob Faibussowitsch   PetscCall(TaoSetRecycleHistory(tao, PETSC_TRUE));
79414d97d3SAlp Dener   reason = TAO_CONTINUE_ITERATING;
80414d97d3SAlp Dener   flg = PETSC_FALSE;
819566063dSJacob Faibussowitsch   PetscCall(TaoGetRecycleHistory(tao, &flg));
829566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Recycle: enabled\n"));
83414d97d3SAlp Dener   while (reason != TAO_CONVERGED_GATOL) {
849566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
859566063dSJacob Faibussowitsch     PetscCall(TaoGetConvergedReason(tao, &reason));
869566063dSJacob Faibussowitsch     PetscCall(TaoGetIterationNumber(tao, &its));
87414d97d3SAlp Dener     recycled_its += its;
889566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
89414d97d3SAlp Dener   }
90414d97d3SAlp Dener 
91414d97d3SAlp Dener   /* Disable recycling and solve again! */
929566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(tao, 100));
939566063dSJacob Faibussowitsch   PetscCall(TaoSetRecycleHistory(tao, PETSC_FALSE));
949566063dSJacob Faibussowitsch   PetscCall(VecSet(x, zero));
959566063dSJacob Faibussowitsch   PetscCall(TaoGetRecycleHistory(tao, &flg));
969566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Recycle: disabled\n"));
979566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
989566063dSJacob Faibussowitsch   PetscCall(TaoGetConvergedReason(tao, &reason));
993c859ba3SBarry Smith   PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
1009566063dSJacob Faibussowitsch   PetscCall(TaoGetIterationNumber(tao, &oneshot_its));
1019566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
10263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "recycled its: %" PetscInt_FMT " | oneshot its: %" PetscInt_FMT "\n", recycled_its, oneshot_its));
1033c859ba3SBarry Smith   PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Recycled solution does not match oneshot solution!");
104414d97d3SAlp Dener 
1059566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&H));
108414d97d3SAlp Dener 
1099566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
110b122ec5aSJacob Faibussowitsch   return 0;
111414d97d3SAlp Dener }
112414d97d3SAlp Dener 
113414d97d3SAlp Dener /* -------------------------------------------------------------------- */
114414d97d3SAlp Dener /*
115414d97d3SAlp Dener     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
116414d97d3SAlp Dener 
117414d97d3SAlp Dener     Input Parameters:
118414d97d3SAlp Dener .   tao  - the Tao context
119414d97d3SAlp Dener .   X    - input vector
120414d97d3SAlp Dener .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()
121414d97d3SAlp Dener 
122414d97d3SAlp Dener     Output Parameters:
123414d97d3SAlp Dener .   G - vector containing the newly evaluated gradient
124414d97d3SAlp Dener .   f - function value
125414d97d3SAlp Dener 
126414d97d3SAlp Dener     Note:
127414d97d3SAlp Dener     Some optimization methods ask for the function and the gradient evaluation
128414d97d3SAlp Dener     at the same time.  Evaluating both at once may be more efficient than
129414d97d3SAlp Dener     evaluating each separately.
130414d97d3SAlp Dener */
131414d97d3SAlp Dener PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
132414d97d3SAlp Dener {
133414d97d3SAlp Dener   AppCtx            *user = (AppCtx *) ptr;
134414d97d3SAlp Dener   PetscInt          i,nn=user->n/2;
135414d97d3SAlp Dener   PetscReal         ff=0,t1,t2,alpha=user->alpha;
136414d97d3SAlp Dener   PetscScalar       *g;
137414d97d3SAlp Dener   const PetscScalar *x;
138414d97d3SAlp Dener 
139414d97d3SAlp Dener   PetscFunctionBeginUser;
140414d97d3SAlp Dener   /* Get pointers to vector data */
1419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
1429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(G,&g));
143414d97d3SAlp Dener 
144414d97d3SAlp Dener   /* Compute G(X) */
145414d97d3SAlp Dener   if (user->chained) {
146414d97d3SAlp Dener     g[0] = 0;
147414d97d3SAlp Dener     for (i=0; i<user->n-1; i++) {
148414d97d3SAlp Dener       t1 = x[i+1] - x[i]*x[i];
149414d97d3SAlp Dener       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
150414d97d3SAlp Dener       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
151414d97d3SAlp Dener       g[i+1] = 2*alpha*t1;
152414d97d3SAlp Dener     }
153414d97d3SAlp Dener   } else {
154414d97d3SAlp Dener     for (i=0; i<nn; i++) {
155414d97d3SAlp Dener       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
156414d97d3SAlp Dener       ff += alpha*t1*t1 + t2*t2;
157414d97d3SAlp Dener       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
158414d97d3SAlp Dener       g[2*i+1] = 2*alpha*t1;
159414d97d3SAlp Dener     }
160414d97d3SAlp Dener   }
161414d97d3SAlp Dener 
162414d97d3SAlp Dener   /* Restore vectors */
1639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(G,&g));
165414d97d3SAlp Dener   *f   = ff;
166414d97d3SAlp Dener 
1679566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(15.0*nn));
168414d97d3SAlp Dener   PetscFunctionReturn(0);
169414d97d3SAlp Dener }
170414d97d3SAlp Dener 
171414d97d3SAlp Dener /* ------------------------------------------------------------------- */
172414d97d3SAlp Dener /*
173414d97d3SAlp Dener    FormHessian - Evaluates Hessian matrix.
174414d97d3SAlp Dener 
175414d97d3SAlp Dener    Input Parameters:
176414d97d3SAlp Dener .  tao   - the Tao context
177414d97d3SAlp Dener .  x     - input vector
178414d97d3SAlp Dener .  ptr   - optional user-defined context, as set by TaoSetHessian()
179414d97d3SAlp Dener 
180414d97d3SAlp Dener    Output Parameters:
181414d97d3SAlp Dener .  H     - Hessian matrix
182414d97d3SAlp Dener 
183414d97d3SAlp Dener    Note:  Providing the Hessian may not be necessary.  Only some solvers
184414d97d3SAlp Dener    require this matrix.
185414d97d3SAlp Dener */
186414d97d3SAlp Dener PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
187414d97d3SAlp Dener {
188414d97d3SAlp Dener   AppCtx            *user = (AppCtx*)ptr;
189414d97d3SAlp Dener   PetscInt          i, ind[2];
190414d97d3SAlp Dener   PetscReal         alpha=user->alpha;
191414d97d3SAlp Dener   PetscReal         v[2][2];
192414d97d3SAlp Dener   const PetscScalar *x;
193414d97d3SAlp Dener   PetscBool         assembled;
194414d97d3SAlp Dener 
195414d97d3SAlp Dener   PetscFunctionBeginUser;
196414d97d3SAlp Dener   /* Zero existing matrix entries */
1979566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H,&assembled));
1989566063dSJacob Faibussowitsch   if (assembled || user->chained) PetscCall(MatZeroEntries(H));
199414d97d3SAlp Dener 
200414d97d3SAlp Dener   /* Get a pointer to vector data */
2019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
202414d97d3SAlp Dener 
203414d97d3SAlp Dener   /* Compute H(X) entries */
204414d97d3SAlp Dener   if (user->chained) {
205414d97d3SAlp Dener     for (i=0; i<user->n-1; i++) {
206414d97d3SAlp Dener       PetscScalar t1 = x[i+1] - x[i]*x[i];
207414d97d3SAlp Dener       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
208414d97d3SAlp Dener       v[0][1] = 2*alpha*(-2*x[i]);
209414d97d3SAlp Dener       v[1][0] = 2*alpha*(-2*x[i]);
210414d97d3SAlp Dener       v[1][1] = 2*alpha*t1;
211414d97d3SAlp Dener       ind[0] = i; ind[1] = i+1;
2129566063dSJacob Faibussowitsch       PetscCall(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES));
213414d97d3SAlp Dener     }
214414d97d3SAlp Dener   } else {
215414d97d3SAlp Dener     for (i=0; i<user->n/2; i++) {
216414d97d3SAlp Dener       v[1][1] = 2*alpha;
217414d97d3SAlp Dener       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
218414d97d3SAlp Dener       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
219414d97d3SAlp Dener       ind[0]=2*i; ind[1]=2*i+1;
2209566063dSJacob Faibussowitsch       PetscCall(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES));
221414d97d3SAlp Dener     }
222414d97d3SAlp Dener   }
2239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
224414d97d3SAlp Dener 
225414d97d3SAlp Dener   /* Assemble matrix */
2269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
2279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
2289566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(9.0*user->n/2.0));
229414d97d3SAlp Dener   PetscFunctionReturn(0);
230414d97d3SAlp Dener }
231414d97d3SAlp Dener 
232414d97d3SAlp Dener /*TEST
233414d97d3SAlp Dener 
234414d97d3SAlp Dener    build:
235414d97d3SAlp Dener       requires: !complex
236414d97d3SAlp Dener 
237414d97d3SAlp Dener    test:
238414d97d3SAlp Dener       args: -tao_type bqnls -tao_monitor
239414d97d3SAlp Dener       requires: !single
240414d97d3SAlp Dener 
241414d97d3SAlp Dener TEST*/
242