xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock2.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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 
28c4762a1bSJed Brown int main(int argc,char **argv)
29c4762a1bSJed Brown {
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 */
419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
433c859ba3SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors");
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Initialize problem parameters */
46c4762a1bSJed Brown   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
47c4762a1bSJed Brown   /* Check for command line arguments to override defaults */
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg));
509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Allocate vectors for the solution and gradient */
549566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.n,&x));
559566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* The TAO code begins here */
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Create TAO solver with desired solution method */
609566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF,&tao));
619566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOLMVM));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* Set solution vec and an initial guess */
649566063dSJacob Faibussowitsch   PetscCall(VecSet(x, zero));
659566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,x));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Set routines for function, gradient, hessian evaluation */
689566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user));
699566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao,H,H,FormHessian,&user));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Check for TAO command line options */
729566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Solve the problem */
759566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0));
769566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(tao, 5));
779566063dSJacob Faibussowitsch   PetscCall(TaoLMVMRecycle(tao, PETSC_TRUE));
78c4762a1bSJed Brown   reason = TAO_CONTINUE_ITERATING;
79c4762a1bSJed Brown   while (reason != TAO_CONVERGED_GATOL) {
809566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
819566063dSJacob Faibussowitsch     PetscCall(TaoGetConvergedReason(tao, &reason));
829566063dSJacob Faibussowitsch     PetscCall(TaoGetIterationNumber(tao, &its));
83c4762a1bSJed Brown     recycled_its += its;
849566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Disable recycling and solve again! */
889566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(tao, 100));
899566063dSJacob Faibussowitsch   PetscCall(TaoLMVMRecycle(tao, PETSC_FALSE));
909566063dSJacob Faibussowitsch   PetscCall(VecSet(x, zero));
919566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
929566063dSJacob Faibussowitsch   PetscCall(TaoGetConvergedReason(tao, &reason));
933c859ba3SBarry Smith   PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
949566063dSJacob Faibussowitsch   PetscCall(TaoGetIterationNumber(tao, &oneshot_its));
959566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
96*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "recycled its: %" PetscInt_FMT " | oneshot its: %" PetscInt_FMT "\n", recycled_its, oneshot_its));
973c859ba3SBarry Smith   PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "LMVM recycling does not work!");
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&H));
102c4762a1bSJed Brown 
1039566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
104b122ec5aSJacob Faibussowitsch   return 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /* -------------------------------------------------------------------- */
108c4762a1bSJed Brown /*
109c4762a1bSJed Brown     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
110c4762a1bSJed Brown 
111c4762a1bSJed Brown     Input Parameters:
112c4762a1bSJed Brown .   tao  - the Tao context
113c4762a1bSJed Brown .   X    - input vector
114c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     Output Parameters:
117c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
118c4762a1bSJed Brown .   f - function value
119c4762a1bSJed Brown 
120c4762a1bSJed Brown     Note:
121c4762a1bSJed Brown     Some optimization methods ask for the function and the gradient evaluation
122c4762a1bSJed Brown     at the same time.  Evaluating both at once may be more efficient that
123c4762a1bSJed Brown     evaluating each separately.
124c4762a1bSJed Brown */
125c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   AppCtx            *user = (AppCtx *) ptr;
128c4762a1bSJed Brown   PetscInt          i,nn=user->n/2;
129c4762a1bSJed Brown   PetscReal         ff=0,t1,t2,alpha=user->alpha;
130c4762a1bSJed Brown   PetscScalar       *g;
131c4762a1bSJed Brown   const PetscScalar *x;
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   PetscFunctionBeginUser;
134c4762a1bSJed Brown   /* Get pointers to vector data */
1359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
1369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G,&g));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* Compute G(X) */
139c4762a1bSJed Brown   if (user->chained) {
140c4762a1bSJed Brown     g[0] = 0;
141c4762a1bSJed Brown     for (i=0; i<user->n-1; i++) {
142c4762a1bSJed Brown       t1 = x[i+1] - x[i]*x[i];
143c4762a1bSJed Brown       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
144c4762a1bSJed Brown       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
145c4762a1bSJed Brown       g[i+1] = 2*alpha*t1;
146c4762a1bSJed Brown     }
147c4762a1bSJed Brown   } else {
148c4762a1bSJed Brown     for (i=0; i<nn; i++) {
149c4762a1bSJed Brown       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
150c4762a1bSJed Brown       ff += alpha*t1*t1 + t2*t2;
151c4762a1bSJed Brown       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
152c4762a1bSJed Brown       g[2*i+1] = 2*alpha*t1;
153c4762a1bSJed Brown     }
154c4762a1bSJed Brown   }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Restore vectors */
1579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
1589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G,&g));
159c4762a1bSJed Brown   *f   = ff;
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(15.0*nn));
162c4762a1bSJed Brown   PetscFunctionReturn(0);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown /* ------------------------------------------------------------------- */
166c4762a1bSJed Brown /*
167c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
168c4762a1bSJed Brown 
169c4762a1bSJed Brown    Input Parameters:
170c4762a1bSJed Brown .  tao   - the Tao context
171c4762a1bSJed Brown .  x     - input vector
172c4762a1bSJed Brown .  ptr   - optional user-defined context, as set by TaoSetHessian()
173c4762a1bSJed Brown 
174c4762a1bSJed Brown    Output Parameters:
175c4762a1bSJed Brown .  H     - Hessian matrix
176c4762a1bSJed Brown 
177c4762a1bSJed Brown    Note:  Providing the Hessian may not be necessary.  Only some solvers
178c4762a1bSJed Brown    require this matrix.
179c4762a1bSJed Brown */
180c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
181c4762a1bSJed Brown {
182c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
183c4762a1bSJed Brown   PetscInt          i, ind[2];
184c4762a1bSJed Brown   PetscReal         alpha=user->alpha;
185c4762a1bSJed Brown   PetscReal         v[2][2];
186c4762a1bSJed Brown   const PetscScalar *x;
187c4762a1bSJed Brown   PetscBool         assembled;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   PetscFunctionBeginUser;
190c4762a1bSJed Brown   /* Zero existing matrix entries */
1919566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H,&assembled));
1929566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(H));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* Get a pointer to vector data */
1959566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* Compute H(X) entries */
198c4762a1bSJed Brown   if (user->chained) {
1999566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(H));
200c4762a1bSJed Brown     for (i=0; i<user->n-1; i++) {
201c4762a1bSJed Brown       PetscScalar t1 = x[i+1] - x[i]*x[i];
202c4762a1bSJed Brown       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
203c4762a1bSJed Brown       v[0][1] = 2*alpha*(-2*x[i]);
204c4762a1bSJed Brown       v[1][0] = 2*alpha*(-2*x[i]);
205c4762a1bSJed Brown       v[1][1] = 2*alpha*t1;
206c4762a1bSJed Brown       ind[0] = i; ind[1] = i+1;
2079566063dSJacob Faibussowitsch       PetscCall(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES));
208c4762a1bSJed Brown     }
209c4762a1bSJed Brown   } else {
210c4762a1bSJed Brown     for (i=0; i<user->n/2; i++) {
211c4762a1bSJed Brown       v[1][1] = 2*alpha;
212c4762a1bSJed Brown       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
213c4762a1bSJed Brown       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
214c4762a1bSJed Brown       ind[0]=2*i; ind[1]=2*i+1;
2159566063dSJacob Faibussowitsch       PetscCall(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES));
216c4762a1bSJed Brown     }
217c4762a1bSJed Brown   }
2189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* Assemble matrix */
2219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
2229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
2239566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(9.0*user->n/2.0));
224c4762a1bSJed Brown   PetscFunctionReturn(0);
225c4762a1bSJed Brown }
226c4762a1bSJed Brown 
227c4762a1bSJed Brown /*TEST
228c4762a1bSJed Brown 
229c4762a1bSJed Brown    build:
230c4762a1bSJed Brown       requires: !complex
231c4762a1bSJed Brown 
232c4762a1bSJed Brown    test:
233c4762a1bSJed Brown       args: -tao_type lmvm -tao_monitor
234c4762a1bSJed Brown       requires: !single
235c4762a1bSJed Brown 
236c4762a1bSJed Brown TEST*/
237