xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock3.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 /*T
14414d97d3SAlp Dener    Concepts: TAO^Solving an unconstrained minimization problem
15414d97d3SAlp Dener    Routines: TaoCreate();
16a82e8c82SStefano Zampini    Routines: TaoSetType(); TaoSetObjectiveAndGradient();
17a82e8c82SStefano Zampini    Routines: TaoSetHessian();
18a82e8c82SStefano Zampini    Routines: TaoSetSolution();
19414d97d3SAlp Dener    Routines: TaoSetFromOptions();
20414d97d3SAlp Dener    Routines: TaoSolve();
21414d97d3SAlp Dener    Routines: TaoDestroy();
22414d97d3SAlp Dener    Processors: 1
23414d97d3SAlp Dener T*/
24414d97d3SAlp Dener 
25414d97d3SAlp Dener /*
26414d97d3SAlp Dener    User-defined application context - contains data needed by the
27414d97d3SAlp Dener    application-provided call-back routines that evaluate the function,
28414d97d3SAlp Dener    gradient, and hessian.
29414d97d3SAlp Dener */
30414d97d3SAlp Dener typedef struct {
31414d97d3SAlp Dener   PetscInt  n;          /* dimension */
32414d97d3SAlp Dener   PetscReal alpha;   /* condition parameter */
33414d97d3SAlp Dener   PetscBool chained;
34414d97d3SAlp Dener } AppCtx;
35414d97d3SAlp Dener 
36414d97d3SAlp Dener /* -------------- User-defined routines ---------- */
37414d97d3SAlp Dener PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
38414d97d3SAlp Dener PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
39414d97d3SAlp Dener 
40414d97d3SAlp Dener int main(int argc,char **argv)
41414d97d3SAlp Dener {
42414d97d3SAlp Dener   PetscErrorCode     ierr;                  /* used to check for functions returning nonzeros */
43414d97d3SAlp Dener   PetscReal          zero=0.0;
44414d97d3SAlp Dener   Vec                x;                     /* solution vector */
45414d97d3SAlp Dener   Mat                H;
46414d97d3SAlp Dener   Tao                tao;                   /* Tao solver context */
47414d97d3SAlp Dener   PetscBool          flg, test_lmvm = PETSC_FALSE;
48414d97d3SAlp Dener   PetscMPIInt        size;                  /* number of processes running */
49414d97d3SAlp Dener   AppCtx             user;                  /* user-defined application context */
50414d97d3SAlp Dener   TaoConvergedReason reason;
51414d97d3SAlp Dener   PetscInt           its, recycled_its=0, oneshot_its=0;
52414d97d3SAlp Dener 
53414d97d3SAlp Dener   /* Initialize TAO and PETSc */
54414d97d3SAlp Dener   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");
57414d97d3SAlp Dener 
58414d97d3SAlp Dener   /* Initialize problem parameters */
59414d97d3SAlp Dener   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
60414d97d3SAlp Dener   /* 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));
65414d97d3SAlp Dener 
66414d97d3SAlp Dener   /* 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));
69414d97d3SAlp Dener 
70414d97d3SAlp Dener   /* The TAO code begins here */
71414d97d3SAlp Dener 
72414d97d3SAlp Dener   /* Create TAO solver with desired solution method */
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBQNLS));
75414d97d3SAlp Dener 
76414d97d3SAlp Dener   /* Set solution vec and an initial guess */
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x, zero));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
79414d97d3SAlp Dener 
80414d97d3SAlp Dener   /* 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));
83414d97d3SAlp Dener 
84414d97d3SAlp Dener   /* Check for TAO command line options */
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
86414d97d3SAlp Dener 
87414d97d3SAlp Dener   /* 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(TaoSetRecycleHistory(tao, PETSC_TRUE));
91414d97d3SAlp Dener   reason = TAO_CONTINUE_ITERATING;
92414d97d3SAlp Dener   flg = PETSC_FALSE;
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetRecycleHistory(tao, &flg));
94*5f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Recycle: enabled\n"));
95414d97d3SAlp Dener   while (reason != TAO_CONVERGED_GATOL) {
96*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSolve(tao));
97*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoGetConvergedReason(tao, &reason));
98*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoGetIterationNumber(tao, &its));
99414d97d3SAlp Dener     recycled_its += its;
100*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
101414d97d3SAlp Dener   }
102414d97d3SAlp Dener 
103414d97d3SAlp Dener   /* Disable recycling and solve again! */
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumIterations(tao, 100));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetRecycleHistory(tao, PETSC_FALSE));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x, zero));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetRecycleHistory(tao, &flg));
108*5f80ce2aSJacob Faibussowitsch   if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Recycle: disabled\n"));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetConvergedReason(tao, &reason));
1113c859ba3SBarry Smith   PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetIterationNumber(tao, &oneshot_its));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "recycled its: %D | oneshot its: %D\n", recycled_its, oneshot_its));
1153c859ba3SBarry Smith   PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Recycled solution does not match oneshot solution!");
116414d97d3SAlp Dener 
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&H));
120414d97d3SAlp Dener 
121414d97d3SAlp Dener   ierr = PetscFinalize();
122414d97d3SAlp Dener   return ierr;
123414d97d3SAlp Dener }
124414d97d3SAlp Dener 
125414d97d3SAlp Dener /* -------------------------------------------------------------------- */
126414d97d3SAlp Dener /*
127414d97d3SAlp Dener     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
128414d97d3SAlp Dener 
129414d97d3SAlp Dener     Input Parameters:
130414d97d3SAlp Dener .   tao  - the Tao context
131414d97d3SAlp Dener .   X    - input vector
132414d97d3SAlp Dener .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()
133414d97d3SAlp Dener 
134414d97d3SAlp Dener     Output Parameters:
135414d97d3SAlp Dener .   G - vector containing the newly evaluated gradient
136414d97d3SAlp Dener .   f - function value
137414d97d3SAlp Dener 
138414d97d3SAlp Dener     Note:
139414d97d3SAlp Dener     Some optimization methods ask for the function and the gradient evaluation
140414d97d3SAlp Dener     at the same time.  Evaluating both at once may be more efficient than
141414d97d3SAlp Dener     evaluating each separately.
142414d97d3SAlp Dener */
143414d97d3SAlp Dener PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
144414d97d3SAlp Dener {
145414d97d3SAlp Dener   AppCtx            *user = (AppCtx *) ptr;
146414d97d3SAlp Dener   PetscInt          i,nn=user->n/2;
147414d97d3SAlp Dener   PetscReal         ff=0,t1,t2,alpha=user->alpha;
148414d97d3SAlp Dener   PetscScalar       *g;
149414d97d3SAlp Dener   const PetscScalar *x;
150414d97d3SAlp Dener 
151414d97d3SAlp Dener   PetscFunctionBeginUser;
152414d97d3SAlp Dener   /* Get pointers to vector data */
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(G,&g));
155414d97d3SAlp Dener 
156414d97d3SAlp Dener   /* Compute G(X) */
157414d97d3SAlp Dener   if (user->chained) {
158414d97d3SAlp Dener     g[0] = 0;
159414d97d3SAlp Dener     for (i=0; i<user->n-1; i++) {
160414d97d3SAlp Dener       t1 = x[i+1] - x[i]*x[i];
161414d97d3SAlp Dener       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
162414d97d3SAlp Dener       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
163414d97d3SAlp Dener       g[i+1] = 2*alpha*t1;
164414d97d3SAlp Dener     }
165414d97d3SAlp Dener   } else {
166414d97d3SAlp Dener     for (i=0; i<nn; i++) {
167414d97d3SAlp Dener       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
168414d97d3SAlp Dener       ff += alpha*t1*t1 + t2*t2;
169414d97d3SAlp Dener       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
170414d97d3SAlp Dener       g[2*i+1] = 2*alpha*t1;
171414d97d3SAlp Dener     }
172414d97d3SAlp Dener   }
173414d97d3SAlp Dener 
174414d97d3SAlp Dener   /* Restore vectors */
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(G,&g));
177414d97d3SAlp Dener   *f   = ff;
178414d97d3SAlp Dener 
179*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(15.0*nn));
180414d97d3SAlp Dener   PetscFunctionReturn(0);
181414d97d3SAlp Dener }
182414d97d3SAlp Dener 
183414d97d3SAlp Dener /* ------------------------------------------------------------------- */
184414d97d3SAlp Dener /*
185414d97d3SAlp Dener    FormHessian - Evaluates Hessian matrix.
186414d97d3SAlp Dener 
187414d97d3SAlp Dener    Input Parameters:
188414d97d3SAlp Dener .  tao   - the Tao context
189414d97d3SAlp Dener .  x     - input vector
190414d97d3SAlp Dener .  ptr   - optional user-defined context, as set by TaoSetHessian()
191414d97d3SAlp Dener 
192414d97d3SAlp Dener    Output Parameters:
193414d97d3SAlp Dener .  H     - Hessian matrix
194414d97d3SAlp Dener 
195414d97d3SAlp Dener    Note:  Providing the Hessian may not be necessary.  Only some solvers
196414d97d3SAlp Dener    require this matrix.
197414d97d3SAlp Dener */
198414d97d3SAlp Dener PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
199414d97d3SAlp Dener {
200414d97d3SAlp Dener   AppCtx            *user = (AppCtx*)ptr;
201414d97d3SAlp Dener   PetscInt          i, ind[2];
202414d97d3SAlp Dener   PetscReal         alpha=user->alpha;
203414d97d3SAlp Dener   PetscReal         v[2][2];
204414d97d3SAlp Dener   const PetscScalar *x;
205414d97d3SAlp Dener   PetscBool         assembled;
206414d97d3SAlp Dener 
207414d97d3SAlp Dener   PetscFunctionBeginUser;
208414d97d3SAlp Dener   /* Zero existing matrix entries */
209*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssembled(H,&assembled));
210*5f80ce2aSJacob Faibussowitsch   if (assembled || user->chained) CHKERRQ(MatZeroEntries(H));
211414d97d3SAlp Dener 
212414d97d3SAlp Dener   /* Get a pointer to vector data */
213*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
214414d97d3SAlp Dener 
215414d97d3SAlp Dener   /* Compute H(X) entries */
216414d97d3SAlp Dener   if (user->chained) {
217414d97d3SAlp Dener     for (i=0; i<user->n-1; i++) {
218414d97d3SAlp Dener       PetscScalar t1 = x[i+1] - x[i]*x[i];
219414d97d3SAlp Dener       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
220414d97d3SAlp Dener       v[0][1] = 2*alpha*(-2*x[i]);
221414d97d3SAlp Dener       v[1][0] = 2*alpha*(-2*x[i]);
222414d97d3SAlp Dener       v[1][1] = 2*alpha*t1;
223414d97d3SAlp Dener       ind[0] = i; ind[1] = i+1;
224*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES));
225414d97d3SAlp Dener     }
226414d97d3SAlp Dener   } else {
227414d97d3SAlp Dener     for (i=0; i<user->n/2; i++) {
228414d97d3SAlp Dener       v[1][1] = 2*alpha;
229414d97d3SAlp Dener       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
230414d97d3SAlp Dener       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
231414d97d3SAlp Dener       ind[0]=2*i; ind[1]=2*i+1;
232*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES));
233414d97d3SAlp Dener     }
234414d97d3SAlp Dener   }
235*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
236414d97d3SAlp Dener 
237414d97d3SAlp Dener   /* Assemble matrix */
238*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(9.0*user->n/2.0));
241414d97d3SAlp Dener   PetscFunctionReturn(0);
242414d97d3SAlp Dener }
243414d97d3SAlp Dener 
244414d97d3SAlp Dener /*TEST
245414d97d3SAlp Dener 
246414d97d3SAlp Dener    build:
247414d97d3SAlp Dener       requires: !complex
248414d97d3SAlp Dener 
249414d97d3SAlp Dener    test:
250414d97d3SAlp Dener       args: -tao_type bqnls -tao_monitor
251414d97d3SAlp Dener       requires: !single
252414d97d3SAlp Dener 
253414d97d3SAlp Dener TEST*/
254