xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock2.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 /*T
14c4762a1bSJed Brown    Concepts: TAO^Solving an unconstrained minimization problem
15c4762a1bSJed Brown    Routines: TaoCreate();
16a82e8c82SStefano Zampini    Routines: TaoSetType(); TaoSetObjectiveAndGradient();
17a82e8c82SStefano Zampini    Routines: TaoSetHessian();
18a82e8c82SStefano Zampini    Routines: TaoSetSolution();
19c4762a1bSJed Brown    Routines: TaoSetFromOptions();
20c4762a1bSJed Brown    Routines: TaoSolve();
21c4762a1bSJed Brown    Routines: TaoDestroy();
22c4762a1bSJed Brown    Processors: 1
23c4762a1bSJed Brown T*/
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    User-defined application context - contains data needed by the
27c4762a1bSJed Brown    application-provided call-back routines that evaluate the function,
28c4762a1bSJed Brown    gradient, and hessian.
29c4762a1bSJed Brown */
30c4762a1bSJed Brown typedef struct {
31c4762a1bSJed Brown   PetscInt  n;          /* dimension */
32c4762a1bSJed Brown   PetscReal alpha;   /* condition parameter */
33c4762a1bSJed Brown   PetscBool chained;
34c4762a1bSJed Brown } AppCtx;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown /* -------------- User-defined routines ---------- */
37c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
38c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
39c4762a1bSJed Brown 
40c4762a1bSJed Brown int main(int argc,char **argv)
41c4762a1bSJed Brown {
42c4762a1bSJed Brown   PetscErrorCode     ierr;                  /* used to check for functions returning nonzeros */
43c4762a1bSJed Brown   PetscReal          zero=0.0;
44c4762a1bSJed Brown   Vec                x;                     /* solution vector */
45c4762a1bSJed Brown   Mat                H;
46c4762a1bSJed Brown   Tao                tao;                   /* Tao solver context */
47c4762a1bSJed Brown   PetscBool          flg, test_lmvm = PETSC_FALSE;
48c4762a1bSJed Brown   PetscMPIInt        size;                  /* number of processes running */
49c4762a1bSJed Brown   AppCtx             user;                  /* user-defined application context */
50c4762a1bSJed Brown   TaoConvergedReason reason;
51c4762a1bSJed Brown   PetscInt           its, recycled_its=0, oneshot_its=0;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Initialize TAO and PETSc */
54c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
55*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
563c859ba3SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors");
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Initialize problem parameters */
59c4762a1bSJed Brown   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
60c4762a1bSJed Brown   /* Check for command line arguments to override defaults */
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* Allocate vectors for the solution and gradient */
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,user.n,&x));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* The TAO code begins here */
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Create TAO solver with desired solution method */
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOLMVM));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Set solution vec and an initial guess */
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x, zero));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Set routines for function, gradient, hessian evaluation */
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetHessian(tao,H,H,FormHessian,&user));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* Check for TAO command line options */
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Solve the problem */
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumIterations(tao, 5));
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoLMVMRecycle(tao, PETSC_TRUE));
91c4762a1bSJed Brown   reason = TAO_CONTINUE_ITERATING;
92c4762a1bSJed Brown   while (reason != TAO_CONVERGED_GATOL) {
93*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSolve(tao));
94*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoGetConvergedReason(tao, &reason));
95*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoGetIterationNumber(tao, &its));
96c4762a1bSJed Brown     recycled_its += its;
97*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
98c4762a1bSJed Brown   }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* Disable recycling and solve again! */
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumIterations(tao, 100));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoLMVMRecycle(tao, PETSC_FALSE));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x, zero));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetConvergedReason(tao, &reason));
1063c859ba3SBarry Smith   PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetIterationNumber(tao, &oneshot_its));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "recycled its: %D | oneshot its: %D\n", recycled_its, oneshot_its));
1103c859ba3SBarry Smith   PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "LMVM recycling does not work!");
111c4762a1bSJed Brown 
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&H));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   ierr = PetscFinalize();
117c4762a1bSJed Brown   return ierr;
118c4762a1bSJed Brown }
119c4762a1bSJed Brown 
120c4762a1bSJed Brown /* -------------------------------------------------------------------- */
121c4762a1bSJed Brown /*
122c4762a1bSJed Brown     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
123c4762a1bSJed Brown 
124c4762a1bSJed Brown     Input Parameters:
125c4762a1bSJed Brown .   tao  - the Tao context
126c4762a1bSJed Brown .   X    - input vector
127c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()
128c4762a1bSJed Brown 
129c4762a1bSJed Brown     Output Parameters:
130c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
131c4762a1bSJed Brown .   f - function value
132c4762a1bSJed Brown 
133c4762a1bSJed Brown     Note:
134c4762a1bSJed Brown     Some optimization methods ask for the function and the gradient evaluation
135c4762a1bSJed Brown     at the same time.  Evaluating both at once may be more efficient that
136c4762a1bSJed Brown     evaluating each separately.
137c4762a1bSJed Brown */
138c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
139c4762a1bSJed Brown {
140c4762a1bSJed Brown   AppCtx            *user = (AppCtx *) ptr;
141c4762a1bSJed Brown   PetscInt          i,nn=user->n/2;
142c4762a1bSJed Brown   PetscReal         ff=0,t1,t2,alpha=user->alpha;
143c4762a1bSJed Brown   PetscScalar       *g;
144c4762a1bSJed Brown   const PetscScalar *x;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   PetscFunctionBeginUser;
147c4762a1bSJed Brown   /* Get pointers to vector data */
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(G,&g));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Compute G(X) */
152c4762a1bSJed Brown   if (user->chained) {
153c4762a1bSJed Brown     g[0] = 0;
154c4762a1bSJed Brown     for (i=0; i<user->n-1; i++) {
155c4762a1bSJed Brown       t1 = x[i+1] - x[i]*x[i];
156c4762a1bSJed Brown       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
157c4762a1bSJed Brown       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
158c4762a1bSJed Brown       g[i+1] = 2*alpha*t1;
159c4762a1bSJed Brown     }
160c4762a1bSJed Brown   } else {
161c4762a1bSJed Brown     for (i=0; i<nn; i++) {
162c4762a1bSJed Brown       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
163c4762a1bSJed Brown       ff += alpha*t1*t1 + t2*t2;
164c4762a1bSJed Brown       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
165c4762a1bSJed Brown       g[2*i+1] = 2*alpha*t1;
166c4762a1bSJed Brown     }
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* Restore vectors */
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(G,&g));
172c4762a1bSJed Brown   *f   = ff;
173c4762a1bSJed Brown 
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(15.0*nn));
175c4762a1bSJed Brown   PetscFunctionReturn(0);
176c4762a1bSJed Brown }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown /* ------------------------------------------------------------------- */
179c4762a1bSJed Brown /*
180c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
181c4762a1bSJed Brown 
182c4762a1bSJed Brown    Input Parameters:
183c4762a1bSJed Brown .  tao   - the Tao context
184c4762a1bSJed Brown .  x     - input vector
185c4762a1bSJed Brown .  ptr   - optional user-defined context, as set by TaoSetHessian()
186c4762a1bSJed Brown 
187c4762a1bSJed Brown    Output Parameters:
188c4762a1bSJed Brown .  H     - Hessian matrix
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    Note:  Providing the Hessian may not be necessary.  Only some solvers
191c4762a1bSJed Brown    require this matrix.
192c4762a1bSJed Brown */
193c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
194c4762a1bSJed Brown {
195c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
196c4762a1bSJed Brown   PetscInt          i, ind[2];
197c4762a1bSJed Brown   PetscReal         alpha=user->alpha;
198c4762a1bSJed Brown   PetscReal         v[2][2];
199c4762a1bSJed Brown   const PetscScalar *x;
200c4762a1bSJed Brown   PetscBool         assembled;
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   PetscFunctionBeginUser;
203c4762a1bSJed Brown   /* Zero existing matrix entries */
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssembled(H,&assembled));
205*5f80ce2aSJacob Faibussowitsch   if (assembled) CHKERRQ(MatZeroEntries(H));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /* Get a pointer to vector data */
208*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* Compute H(X) entries */
211c4762a1bSJed Brown   if (user->chained) {
212*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(H));
213c4762a1bSJed Brown     for (i=0; i<user->n-1; i++) {
214c4762a1bSJed Brown       PetscScalar t1 = x[i+1] - x[i]*x[i];
215c4762a1bSJed Brown       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
216c4762a1bSJed Brown       v[0][1] = 2*alpha*(-2*x[i]);
217c4762a1bSJed Brown       v[1][0] = 2*alpha*(-2*x[i]);
218c4762a1bSJed Brown       v[1][1] = 2*alpha*t1;
219c4762a1bSJed Brown       ind[0] = i; ind[1] = i+1;
220*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES));
221c4762a1bSJed Brown     }
222c4762a1bSJed Brown   } else {
223c4762a1bSJed Brown     for (i=0; i<user->n/2; i++) {
224c4762a1bSJed Brown       v[1][1] = 2*alpha;
225c4762a1bSJed Brown       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
226c4762a1bSJed Brown       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
227c4762a1bSJed Brown       ind[0]=2*i; ind[1]=2*i+1;
228*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES));
229c4762a1bSJed Brown     }
230c4762a1bSJed Brown   }
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* Assemble matrix */
234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
235*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
236*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(9.0*user->n/2.0));
237c4762a1bSJed Brown   PetscFunctionReturn(0);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown /*TEST
241c4762a1bSJed Brown 
242c4762a1bSJed Brown    build:
243c4762a1bSJed Brown       requires: !complex
244c4762a1bSJed Brown 
245c4762a1bSJed Brown    test:
246c4762a1bSJed Brown       args: -tao_type lmvm -tao_monitor
247c4762a1bSJed Brown       requires: !single
248c4762a1bSJed Brown 
249c4762a1bSJed Brown TEST*/
250