xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscReal          zero=0.0;
43414d97d3SAlp Dener   Vec                x;                     /* solution vector */
44414d97d3SAlp Dener   Mat                H;
45414d97d3SAlp Dener   Tao                tao;                   /* Tao solver context */
46414d97d3SAlp Dener   PetscBool          flg, test_lmvm = PETSC_FALSE;
47414d97d3SAlp Dener   PetscMPIInt        size;                  /* number of processes running */
48414d97d3SAlp Dener   AppCtx             user;                  /* user-defined application context */
49414d97d3SAlp Dener   TaoConvergedReason reason;
50414d97d3SAlp Dener   PetscInt           its, recycled_its=0, oneshot_its=0;
51414d97d3SAlp Dener 
52414d97d3SAlp Dener   /* Initialize TAO and PETSc */
53*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
545f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
553c859ba3SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors");
56414d97d3SAlp Dener 
57414d97d3SAlp Dener   /* Initialize problem parameters */
58414d97d3SAlp Dener   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
59414d97d3SAlp Dener   /* Check for command line arguments to override defaults */
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg));
64414d97d3SAlp Dener 
65414d97d3SAlp Dener   /* Allocate vectors for the solution and gradient */
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,user.n,&x));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H));
68414d97d3SAlp Dener 
69414d97d3SAlp Dener   /* The TAO code begins here */
70414d97d3SAlp Dener 
71414d97d3SAlp Dener   /* Create TAO solver with desired solution method */
725f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBQNLS));
74414d97d3SAlp Dener 
75414d97d3SAlp Dener   /* Set solution vec and an initial guess */
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x, zero));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
78414d97d3SAlp Dener 
79414d97d3SAlp Dener   /* Set routines for function, gradient, hessian evaluation */
805f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetHessian(tao,H,H,FormHessian,&user));
82414d97d3SAlp Dener 
83414d97d3SAlp Dener   /* Check for TAO command line options */
845f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
85414d97d3SAlp Dener 
86414d97d3SAlp Dener   /* Solve the problem */
875f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetTolerances(tao, 1.e-5, 0.0, 0.0));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumIterations(tao, 5));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetRecycleHistory(tao, PETSC_TRUE));
90414d97d3SAlp Dener   reason = TAO_CONTINUE_ITERATING;
91414d97d3SAlp Dener   flg = PETSC_FALSE;
925f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetRecycleHistory(tao, &flg));
935f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Recycle: enabled\n"));
94414d97d3SAlp Dener   while (reason != TAO_CONVERGED_GATOL) {
955f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSolve(tao));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoGetConvergedReason(tao, &reason));
975f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoGetIterationNumber(tao, &its));
98414d97d3SAlp Dener     recycled_its += its;
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
100414d97d3SAlp Dener   }
101414d97d3SAlp Dener 
102414d97d3SAlp Dener   /* Disable recycling and solve again! */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumIterations(tao, 100));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetRecycleHistory(tao, PETSC_FALSE));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x, zero));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetRecycleHistory(tao, &flg));
1075f80ce2aSJacob Faibussowitsch   if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Recycle: disabled\n"));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetConvergedReason(tao, &reason));
1103c859ba3SBarry Smith   PetscCheck(reason == TAO_CONVERGED_GATOL,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Solution failed to converge!");
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetIterationNumber(tao, &oneshot_its));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-----------------------\n"));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "recycled its: %D | oneshot its: %D\n", recycled_its, oneshot_its));
1143c859ba3SBarry Smith   PetscCheck(recycled_its == oneshot_its,PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Recycled solution does not match oneshot solution!");
115414d97d3SAlp Dener 
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&H));
119414d97d3SAlp Dener 
120*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
121*b122ec5aSJacob Faibussowitsch   return 0;
122414d97d3SAlp Dener }
123414d97d3SAlp Dener 
124414d97d3SAlp Dener /* -------------------------------------------------------------------- */
125414d97d3SAlp Dener /*
126414d97d3SAlp Dener     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
127414d97d3SAlp Dener 
128414d97d3SAlp Dener     Input Parameters:
129414d97d3SAlp Dener .   tao  - the Tao context
130414d97d3SAlp Dener .   X    - input vector
131414d97d3SAlp Dener .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()
132414d97d3SAlp Dener 
133414d97d3SAlp Dener     Output Parameters:
134414d97d3SAlp Dener .   G - vector containing the newly evaluated gradient
135414d97d3SAlp Dener .   f - function value
136414d97d3SAlp Dener 
137414d97d3SAlp Dener     Note:
138414d97d3SAlp Dener     Some optimization methods ask for the function and the gradient evaluation
139414d97d3SAlp Dener     at the same time.  Evaluating both at once may be more efficient than
140414d97d3SAlp Dener     evaluating each separately.
141414d97d3SAlp Dener */
142414d97d3SAlp Dener PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
143414d97d3SAlp Dener {
144414d97d3SAlp Dener   AppCtx            *user = (AppCtx *) ptr;
145414d97d3SAlp Dener   PetscInt          i,nn=user->n/2;
146414d97d3SAlp Dener   PetscReal         ff=0,t1,t2,alpha=user->alpha;
147414d97d3SAlp Dener   PetscScalar       *g;
148414d97d3SAlp Dener   const PetscScalar *x;
149414d97d3SAlp Dener 
150414d97d3SAlp Dener   PetscFunctionBeginUser;
151414d97d3SAlp Dener   /* Get pointers to vector data */
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(G,&g));
154414d97d3SAlp Dener 
155414d97d3SAlp Dener   /* Compute G(X) */
156414d97d3SAlp Dener   if (user->chained) {
157414d97d3SAlp Dener     g[0] = 0;
158414d97d3SAlp Dener     for (i=0; i<user->n-1; i++) {
159414d97d3SAlp Dener       t1 = x[i+1] - x[i]*x[i];
160414d97d3SAlp Dener       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
161414d97d3SAlp Dener       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
162414d97d3SAlp Dener       g[i+1] = 2*alpha*t1;
163414d97d3SAlp Dener     }
164414d97d3SAlp Dener   } else {
165414d97d3SAlp Dener     for (i=0; i<nn; i++) {
166414d97d3SAlp Dener       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
167414d97d3SAlp Dener       ff += alpha*t1*t1 + t2*t2;
168414d97d3SAlp Dener       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
169414d97d3SAlp Dener       g[2*i+1] = 2*alpha*t1;
170414d97d3SAlp Dener     }
171414d97d3SAlp Dener   }
172414d97d3SAlp Dener 
173414d97d3SAlp Dener   /* Restore vectors */
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(G,&g));
176414d97d3SAlp Dener   *f   = ff;
177414d97d3SAlp Dener 
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(15.0*nn));
179414d97d3SAlp Dener   PetscFunctionReturn(0);
180414d97d3SAlp Dener }
181414d97d3SAlp Dener 
182414d97d3SAlp Dener /* ------------------------------------------------------------------- */
183414d97d3SAlp Dener /*
184414d97d3SAlp Dener    FormHessian - Evaluates Hessian matrix.
185414d97d3SAlp Dener 
186414d97d3SAlp Dener    Input Parameters:
187414d97d3SAlp Dener .  tao   - the Tao context
188414d97d3SAlp Dener .  x     - input vector
189414d97d3SAlp Dener .  ptr   - optional user-defined context, as set by TaoSetHessian()
190414d97d3SAlp Dener 
191414d97d3SAlp Dener    Output Parameters:
192414d97d3SAlp Dener .  H     - Hessian matrix
193414d97d3SAlp Dener 
194414d97d3SAlp Dener    Note:  Providing the Hessian may not be necessary.  Only some solvers
195414d97d3SAlp Dener    require this matrix.
196414d97d3SAlp Dener */
197414d97d3SAlp Dener PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
198414d97d3SAlp Dener {
199414d97d3SAlp Dener   AppCtx            *user = (AppCtx*)ptr;
200414d97d3SAlp Dener   PetscInt          i, ind[2];
201414d97d3SAlp Dener   PetscReal         alpha=user->alpha;
202414d97d3SAlp Dener   PetscReal         v[2][2];
203414d97d3SAlp Dener   const PetscScalar *x;
204414d97d3SAlp Dener   PetscBool         assembled;
205414d97d3SAlp Dener 
206414d97d3SAlp Dener   PetscFunctionBeginUser;
207414d97d3SAlp Dener   /* Zero existing matrix entries */
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssembled(H,&assembled));
2095f80ce2aSJacob Faibussowitsch   if (assembled || user->chained) CHKERRQ(MatZeroEntries(H));
210414d97d3SAlp Dener 
211414d97d3SAlp Dener   /* Get a pointer to vector data */
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
213414d97d3SAlp Dener 
214414d97d3SAlp Dener   /* Compute H(X) entries */
215414d97d3SAlp Dener   if (user->chained) {
216414d97d3SAlp Dener     for (i=0; i<user->n-1; i++) {
217414d97d3SAlp Dener       PetscScalar t1 = x[i+1] - x[i]*x[i];
218414d97d3SAlp Dener       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
219414d97d3SAlp Dener       v[0][1] = 2*alpha*(-2*x[i]);
220414d97d3SAlp Dener       v[1][0] = 2*alpha*(-2*x[i]);
221414d97d3SAlp Dener       v[1][1] = 2*alpha*t1;
222414d97d3SAlp Dener       ind[0] = i; ind[1] = i+1;
2235f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES));
224414d97d3SAlp Dener     }
225414d97d3SAlp Dener   } else {
226414d97d3SAlp Dener     for (i=0; i<user->n/2; i++) {
227414d97d3SAlp Dener       v[1][1] = 2*alpha;
228414d97d3SAlp Dener       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
229414d97d3SAlp Dener       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
230414d97d3SAlp Dener       ind[0]=2*i; ind[1]=2*i+1;
2315f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES));
232414d97d3SAlp Dener     }
233414d97d3SAlp Dener   }
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
235414d97d3SAlp Dener 
236414d97d3SAlp Dener   /* Assemble matrix */
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(9.0*user->n/2.0));
240414d97d3SAlp Dener   PetscFunctionReturn(0);
241414d97d3SAlp Dener }
242414d97d3SAlp Dener 
243414d97d3SAlp Dener /*TEST
244414d97d3SAlp Dener 
245414d97d3SAlp Dener    build:
246414d97d3SAlp Dener       requires: !complex
247414d97d3SAlp Dener 
248414d97d3SAlp Dener    test:
249414d97d3SAlp Dener       args: -tao_type bqnls -tao_monitor
250414d97d3SAlp Dener       requires: !single
251414d97d3SAlp Dener 
252414d97d3SAlp Dener TEST*/
253