xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision 65f8aed5f7eaa1e2ef2ddeffe666264e0669c876)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
3 
4 #define LMVM_STEP_BFGS     0
5 #define LMVM_STEP_GRAD     1
6 
7 static PetscErrorCode TaoSolve_LMVM(Tao tao)
8 {
9   TAO_LMVM                     *lmP = (TAO_LMVM *)tao->data;
10   PetscReal                    f, fold, gdx, gnorm;
11   PetscReal                    step = 1.0;
12   PetscErrorCode               ierr;
13   PetscInt                     stepType = LMVM_STEP_GRAD, nupdates;
14   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
15 
16   PetscFunctionBegin;
17 
18   if (tao->XL || tao->XU || tao->ops->computebounds) {
19     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr);
20   }
21 
22   /*  Check convergence criteria */
23   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
24   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
25 
26   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
27 
28   tao->reason = TAO_CONTINUE_ITERATING;
29   ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
30   ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr);
31   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
32   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
33 
34   /*  Set counter for gradient/reset steps */
35   if (!lmP->recycle) {
36     lmP->bfgs = 0;
37     lmP->grad = 0;
38     ierr = MatLMVMReset(lmP->M, PETSC_FALSE); CHKERRQ(ierr);
39   }
40 
41   /*  Have not converged; continue with Newton method */
42   while (tao->reason == TAO_CONTINUE_ITERATING) {
43     /*  Compute direction */
44     if (lmP->H0) {
45       ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr);
46       stepType = LMVM_STEP_BFGS;
47     }
48     ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr);
49     ierr = MatSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
50     ierr = MatLMVMGetUpdateCount(lmP->M, &nupdates); CHKERRQ(ierr);
51     if (nupdates > 0) stepType = LMVM_STEP_BFGS;
52 
53     /*  Check for success (descent direction) */
54     ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr);
55     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
56       /* Step is not descent or direction produced not a number
57          We can assert bfgsUpdates > 1 in this case because
58          the first solve produces the scaled gradient direction,
59          which is guaranteed to be descent
60 
61          Use steepest descent direction (scaled)
62       */
63 
64       ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr);
65       ierr = MatLMVMClearJ0(lmP->M);CHKERRQ(ierr);
66       ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
67       ierr = MatSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr);
68 
69       /* On a reset, the direction cannot be not a number; it is a
70          scaled gradient step.  No need to check for this condition. */
71       stepType = LMVM_STEP_GRAD;
72     }
73     ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
74 
75     /*  Perform the linesearch */
76     fold = f;
77     ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr);
78     ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr);
79 
80     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr);
81     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
82 
83     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) {
84       /*  Reset factors and use scaled gradient step */
85       f = fold;
86       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
87       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
88 
89       /*  Failed to obtain acceptable iterate with BFGS step */
90       /*  Attempt to use the scaled gradient direction */
91 
92       ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr);
93       ierr = MatLMVMClearJ0(lmP->M);CHKERRQ(ierr);
94       ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
95       ierr = MatSolve(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
96 
97       /* On a reset, the direction cannot be not a number; it is a
98           scaled gradient step.  No need to check for this condition. */
99       stepType = LMVM_STEP_GRAD;
100       ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
101 
102       /*  Perform the linesearch */
103       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr);
104       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
105     }
106 
107     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
108       /*  Failed to find an improving point */
109       f = fold;
110       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
111       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
112       step = 0.0;
113       tao->reason = TAO_DIVERGED_LS_FAILURE;
114     } else {
115       /* LS found valid step, so tally up step type */
116       switch (stepType) {
117       case LMVM_STEP_BFGS:
118         ++lmP->bfgs;
119         break;
120       case LMVM_STEP_GRAD:
121         ++lmP->grad;
122         break;
123       default:
124         break;
125       }
126       /*  Compute new gradient norm */
127       ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
128     }
129 
130     /* Check convergence */
131     tao->niter++;
132     ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
133     ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr);
134     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
135   }
136   PetscFunctionReturn(0);
137 }
138 
139 static PetscErrorCode TaoSetUp_LMVM(Tao tao)
140 {
141   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
142   PetscInt       n,N;
143   PetscErrorCode ierr;
144   PetscBool      is_spd;
145 
146   PetscFunctionBegin;
147   /* Existence of tao->solution checked in TaoSetUp() */
148   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);  }
149   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);  }
150   if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr);  }
151   if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr);  }
152   if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr);  }
153 
154   /*  Create matrix for the limited memory approximation */
155   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
156   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
157   ierr = MatSetSizes(lmP->M, n, n, N, N);CHKERRQ(ierr);
158   ierr = MatLMVMAllocate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr);
159   ierr = MatGetOption(lmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr);
160   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite.");
161 
162   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
163   if (lmP->H0) {
164     ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr);
165   }
166 
167   PetscFunctionReturn(0);
168 }
169 
170 /* ---------------------------------------------------------- */
171 static PetscErrorCode TaoDestroy_LMVM(Tao tao)
172 {
173   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
174   PetscErrorCode ierr;
175 
176   PetscFunctionBegin;
177   if (tao->setupcalled) {
178     ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr);
179     ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr);
180     ierr = VecDestroy(&lmP->D);CHKERRQ(ierr);
181   }
182   ierr = MatDestroy(&lmP->M);CHKERRQ(ierr);
183   if (lmP->H0) {
184     ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr);
185   }
186   ierr = PetscFree(tao->data);CHKERRQ(ierr);
187 
188   PetscFunctionReturn(0);
189 }
190 
191 /*------------------------------------------------------------*/
192 static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao)
193 {
194   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
195   PetscErrorCode ierr;
196 
197   PetscFunctionBegin;
198   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr);
199   ierr = PetscOptionsBool("-tao_lmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",lm->recycle,&lm->recycle,NULL);CHKERRQ(ierr);
200   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
201   ierr = MatSetFromOptions(lm->M);CHKERRQ(ierr);
202   ierr = PetscOptionsTail();CHKERRQ(ierr);
203   PetscFunctionReturn(0);
204 }
205 
206 /*------------------------------------------------------------*/
207 static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
208 {
209   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
210   PetscBool      isascii;
211   PetscInt       recycled_its;
212   PetscErrorCode ierr;
213 
214   PetscFunctionBegin;
215   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
216   if (isascii) {
217     ierr = PetscViewerASCIIPrintf(viewer, "  Gradient steps: %D\n", lm->grad);CHKERRQ(ierr);
218     if (lm->recycle) {
219       ierr = PetscViewerASCIIPrintf(viewer, "  Recycle: on\n");CHKERRQ(ierr);
220       recycled_its = lm->bfgs + lm->grad;
221       ierr = PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr);
222     }
223   }
224   PetscFunctionReturn(0);
225 }
226 
227 /* ---------------------------------------------------------- */
228 
229 /*MC
230   TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
231   optimization solver for unconstrained minimization. It solves
232   the Newton step
233           Hkdk = - gk
234 
235   using an approximation Bk in place of Hk, where Bk is composed using
236   the BFGS update formula. A More-Thuente line search is then used
237   to computed the steplength in the dk direction
238 
239   Options Database Keys:
240 .   -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls
241 .   -tao_lmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation
242 
243   Level: beginner
244 M*/
245 
246 PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
247 {
248   TAO_LMVM       *lmP;
249   const char     *morethuente_type = TAOLINESEARCHMT;
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   tao->ops->setup = TaoSetUp_LMVM;
254   tao->ops->solve = TaoSolve_LMVM;
255   tao->ops->view = TaoView_LMVM;
256   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
257   tao->ops->destroy = TaoDestroy_LMVM;
258 
259   ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr);
260   lmP->D = 0;
261   lmP->M = 0;
262   lmP->Xold = 0;
263   lmP->Gold = 0;
264   lmP->H0   = NULL;
265   lmP->recycle = PETSC_FALSE;
266 
267   tao->data = (void*)lmP;
268   /* Override default settings (unless already changed) */
269   if (!tao->max_it_changed) tao->max_it = 2000;
270   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
271 
272   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
273   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
274   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
275   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
276   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
277 
278   ierr = KSPInitializePackage();CHKERRQ(ierr);
279   ierr = MatCreate(((PetscObject)tao)->comm, &lmP->M);CHKERRQ(ierr);
280   ierr = PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1);CHKERRQ(ierr);
281   ierr = MatSetType(lmP->M, MATLMVMBFGS);CHKERRQ(ierr);
282   ierr = MatSetOptionsPrefix(lmP->M, "tao_lmvm_");CHKERRQ(ierr);
283   PetscFunctionReturn(0);
284 }
285