Actual source code: blmvm.c

  1: #include <petsctaolinesearch.h>
  2: #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
  3: #include <../src/tao/bound/impls/blmvm/blmvm.h>

  5: /*------------------------------------------------------------*/
  6: static PetscErrorCode TaoSolve_BLMVM(Tao tao)
  7: {
  8:   PetscErrorCode               ierr;
  9:   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
 10:   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
 11:   PetscReal                    f, fold, gdx, gnorm, gnorm2;
 12:   PetscReal                    stepsize = 1.0,delta;

 15:   /*  Project initial point onto bounds */
 16:   TaoComputeVariableBounds(tao);
 17:   VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
 18:   TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);


 21:   /* Check convergence criteria */
 22:   TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);
 23:   VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);

 25:   TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);
 26:   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");

 28:   tao->reason = TAO_CONTINUE_ITERATING;
 29:   TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);
 30:   TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);
 31:   (*tao->ops->convergencetest)(tao,tao->cnvP);
 32:   if (tao->reason != TAO_CONTINUE_ITERATING) return(0);

 34:   /* Set counter for gradient/reset steps */
 35:   if (!blmP->recycle) {
 36:     blmP->grad = 0;
 37:     blmP->reset = 0;
 38:     MatLMVMReset(blmP->M, PETSC_FALSE);
 39:   }

 41:   /* Have not converged; continue with Newton method */
 42:   while (tao->reason == TAO_CONTINUE_ITERATING) {
 43:     /* Call general purpose update function */
 44:     if (tao->ops->update) {
 45:       (*tao->ops->update)(tao, tao->niter, tao->user_update);
 46:     }
 47:     /* Compute direction */
 48:     gnorm2 = gnorm*gnorm;
 49:     if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
 50:     if (f == 0.0) {
 51:       delta = 2.0 / gnorm2;
 52:     } else {
 53:       delta = 2.0 * PetscAbsScalar(f) / gnorm2;
 54:     }
 55:     MatLMVMSymBroydenSetDelta(blmP->M, delta);
 56:     MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);
 57:     MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
 58:     VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);

 60:     /* Check for success (descent direction) */
 61:     VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);
 62:     if (gdx <= 0) {
 63:       /* Step is not descent or solve was not successful
 64:          Use steepest descent direction (scaled) */
 65:       ++blmP->grad;

 67:       MatLMVMReset(blmP->M, PETSC_FALSE);
 68:       MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
 69:       MatSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);
 70:     }
 71:     VecScale(tao->stepdirection,-1.0);

 73:     /* Perform the linesearch */
 74:     fold = f;
 75:     VecCopy(tao->solution, blmP->Xold);
 76:     VecCopy(blmP->unprojected_gradient, blmP->Gold);
 77:     TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
 78:     TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);
 79:     TaoAddLineSearchCounts(tao);

 81:     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
 82:       /* Linesearch failed
 83:          Reset factors and use scaled (projected) gradient step */
 84:       ++blmP->reset;

 86:       f = fold;
 87:       VecCopy(blmP->Xold, tao->solution);
 88:       VecCopy(blmP->Gold, blmP->unprojected_gradient);

 90:       MatLMVMReset(blmP->M, PETSC_FALSE);
 91:       MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
 92:       MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
 93:       VecScale(tao->stepdirection, -1.0);

 95:       /* This may be incorrect; linesearch has values for stepmax and stepmin
 96:          that should be reset. */
 97:       TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
 98:       TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);
 99:       TaoAddLineSearchCounts(tao);

101:       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
102:         tao->reason = TAO_DIVERGED_LS_FAILURE;
103:         break;
104:       }
105:     }

107:     /* Check for converged */
108:     VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);
109:     TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);
110:     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Not-a-Number");
111:     tao->niter++;
112:     TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);
113:     TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);
114:     (*tao->ops->convergencetest)(tao,tao->cnvP);
115:   }
116:   return(0);
117: }

119: static PetscErrorCode TaoSetup_BLMVM(Tao tao)
120: {
121:   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;

125:   /* Existence of tao->solution checked in TaoSetup() */
126:   VecDuplicate(tao->solution,&blmP->Xold);
127:   VecDuplicate(tao->solution,&blmP->Gold);
128:   VecDuplicate(tao->solution, &blmP->unprojected_gradient);

130:   if (!tao->stepdirection) {
131:     VecDuplicate(tao->solution, &tao->stepdirection);
132:   }
133:   if (!tao->gradient) {
134:     VecDuplicate(tao->solution,&tao->gradient);
135:   }
136:   if (!tao->XL) {
137:     VecDuplicate(tao->solution,&tao->XL);
138:     VecSet(tao->XL,PETSC_NINFINITY);
139:   }
140:   if (!tao->XU) {
141:     VecDuplicate(tao->solution,&tao->XU);
142:     VecSet(tao->XU,PETSC_INFINITY);
143:   }
144:   /* Allocate matrix for the limited memory approximation */
145:   MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);

147:   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
148:   if (blmP->H0) {
149:     MatLMVMSetJ0(blmP->M, blmP->H0);
150:   }
151:   return(0);
152: }

154: /* ---------------------------------------------------------- */
155: static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
156: {
157:   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;

161:   if (tao->setupcalled) {
162:     VecDestroy(&blmP->unprojected_gradient);
163:     VecDestroy(&blmP->Xold);
164:     VecDestroy(&blmP->Gold);
165:   }
166:   MatDestroy(&blmP->M);
167:   if (blmP->H0) {
168:     PetscObjectDereference((PetscObject)blmP->H0);
169:   }
170:   PetscFree(tao->data);
171:   return(0);
172: }

174: /*------------------------------------------------------------*/
175: static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
176: {
177:   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
179:   PetscBool      is_spd;

182:   PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");
183:   PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);
184:   PetscOptionsTail();
185:   TaoLineSearchSetFromOptions(tao->linesearch);
186:   MatSetFromOptions(blmP->M);
187:   MatGetOption(blmP->M, MAT_SPD, &is_spd);
188:   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
189:   return(0);
190: }


193: /*------------------------------------------------------------*/
194: static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
195: {
196:   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
197:   PetscBool      isascii;

201:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
202:   if (isascii) {
203:     PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);
204:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);
205:     MatView(lmP->M, viewer);
206:     PetscViewerPopFormat(viewer);
207:   }
208:   return(0);
209: }

211: static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
212: {
213:   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;

220:   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");

222:   VecCopy(tao->gradient,DXL);
223:   VecAXPY(DXL,-1.0,blm->unprojected_gradient);
224:   VecSet(DXU,0.0);
225:   VecPointwiseMax(DXL,DXL,DXU);

227:   VecCopy(blm->unprojected_gradient,DXU);
228:   VecAXPY(DXU,-1.0,tao->gradient);
229:   VecAXPY(DXU,1.0,DXL);
230:   return(0);
231: }

233: /* ---------------------------------------------------------- */
234: /*MC
235:   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
236:          for nonlinear minimization with bound constraints. It is an extension
237:          of TAOLMVM

239:   Options Database Keys:
240: .     -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls

242:   Level: beginner
243: M*/
244: PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
245: {
246:   TAO_BLMVM      *blmP;
247:   const char     *morethuente_type = TAOLINESEARCHMT;

251:   tao->ops->setup = TaoSetup_BLMVM;
252:   tao->ops->solve = TaoSolve_BLMVM;
253:   tao->ops->view = TaoView_BLMVM;
254:   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
255:   tao->ops->destroy = TaoDestroy_BLMVM;
256:   tao->ops->computedual = TaoComputeDual_BLMVM;

258:   PetscNewLog(tao,&blmP);
259:   blmP->H0 = NULL;
260:   blmP->recycle = PETSC_FALSE;
261:   tao->data = (void*)blmP;

263:   /* Override default settings (unless already changed) */
264:   if (!tao->max_it_changed) tao->max_it = 2000;
265:   if (!tao->max_funcs_changed) tao->max_funcs = 4000;

267:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
268:   PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
269:   TaoLineSearchSetType(tao->linesearch, morethuente_type);
270:   TaoLineSearchUseTaoRoutines(tao->linesearch,tao);
271:   TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);

273:   KSPInitializePackage();
274:   MatCreate(((PetscObject)tao)->comm, &blmP->M);
275:   MatSetType(blmP->M, MATLMVMBFGS);
276:   PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);
277:   MatSetOptionsPrefix(blmP->M, "tao_blmvm_");
278:   return(0);
279: }

281: /*@
282:   TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls.

284:   Input Parameters:
285: +  tao  - the Tao solver context
286: -  flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE)

288:   Level: intermediate
289: @*/
290: PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg)
291: {
292:   TAO_LMVM       *lmP;
293:   TAO_BLMVM      *blmP;
294:   PetscBool      is_lmvm, is_blmvm;

298:   PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);
299:   PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);
300:   if (is_lmvm) {
301:     lmP = (TAO_LMVM *)tao->data;
302:     lmP->recycle = flg;
303:   } else if (is_blmvm) {
304:     blmP = (TAO_BLMVM *)tao->data;
305:     blmP->recycle = flg;
306:   }
307:   return(0);
308: }

310: /*@
311:   TaoLMVMSetH0 - Set the initial Hessian for the QN approximation

313:   Input Parameters:
314: +  tao  - the Tao solver context
315: -  H0 - Mat object for the initial Hessian

317:   Level: advanced

319: .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP()
320: @*/
321: PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
322: {
323:   TAO_LMVM       *lmP;
324:   TAO_BLMVM      *blmP;
325:   PetscBool      is_lmvm, is_blmvm;

329:   PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);
330:   PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);
331:   if (is_lmvm) {
332:     lmP = (TAO_LMVM *)tao->data;
333:     PetscObjectReference((PetscObject)H0);
334:     lmP->H0 = H0;
335:   } else if (is_blmvm) {
336:     blmP = (TAO_BLMVM *)tao->data;
337:     PetscObjectReference((PetscObject)H0);
338:     blmP->H0 = H0;
339:   }
340:   return(0);
341: }

343: /*@
344:   TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian

346:   Input Parameters:
347: .  tao  - the Tao solver context

349:   Output Parameters:
350: .  H0 - Mat object for the initial Hessian

352:   Level: advanced

354: .seealso: TaoLMVMSetH0(), TaoLMVMGetH0KSP()
355: @*/
356: PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
357: {
358:   TAO_LMVM       *lmP;
359:   TAO_BLMVM      *blmP;
360:   PetscBool      is_lmvm, is_blmvm;
361:   Mat            M;

365:   PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);
366:   PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);
367:   if (is_lmvm) {
368:     lmP = (TAO_LMVM *)tao->data;
369:     M = lmP->M;
370:   } else if (is_blmvm) {
371:     blmP = (TAO_BLMVM *)tao->data;
372:     M = blmP->M;
373:   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM.");
374:   MatLMVMGetJ0(M, H0);
375:   return(0);
376: }

378: /*@
379:   TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian

381:   Input Parameters:
382: .  tao  - the Tao solver context

384:   Output Parameters:
385: .  ksp - KSP solver context for the initial Hessian

387:   Level: advanced

389: .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP()
390: @*/
391: PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
392: {
393:   TAO_LMVM       *lmP;
394:   TAO_BLMVM      *blmP;
395:   PetscBool      is_lmvm, is_blmvm;
396:   Mat            M;

400:   PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm);
401:   PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm);
402:   if (is_lmvm) {
403:     lmP = (TAO_LMVM *)tao->data;
404:     M = lmP->M;
405:   } else if (is_blmvm) {
406:     blmP = (TAO_BLMVM *)tao->data;
407:     M = blmP->M;
408:   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM.");
409:   MatLMVMGetJ0KSP(M, ksp);
410:   return(0);
411: }