Actual source code: blmvm.c
petsc-3.10.5 2019-03-28
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(PETSC_COMM_SELF,1, "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: /* Compute direction */
44: gnorm2 = gnorm*gnorm;
45: if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
46: if (f == 0.0) {
47: delta = 2.0 / gnorm2;
48: } else {
49: delta = 2.0 * PetscAbsScalar(f) / gnorm2;
50: }
51: MatSymBrdnSetDelta(blmP->M, delta);
52: MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);
53: MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
54: VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);
56: /* Check for success (descent direction) */
57: VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);
58: if (gdx <= 0) {
59: /* Step is not descent or solve was not successful
60: Use steepest descent direction (scaled) */
61: ++blmP->grad;
63: MatLMVMReset(blmP->M, PETSC_FALSE);
64: MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
65: MatSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);
66: }
67: VecScale(tao->stepdirection,-1.0);
69: /* Perform the linesearch */
70: fold = f;
71: VecCopy(tao->solution, blmP->Xold);
72: VecCopy(blmP->unprojected_gradient, blmP->Gold);
73: TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
74: TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);
75: TaoAddLineSearchCounts(tao);
77: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
78: /* Linesearch failed
79: Reset factors and use scaled (projected) gradient step */
80: ++blmP->reset;
82: f = fold;
83: VecCopy(blmP->Xold, tao->solution);
84: VecCopy(blmP->Gold, blmP->unprojected_gradient);
86: MatLMVMReset(blmP->M, PETSC_FALSE);
87: MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
88: MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
89: VecScale(tao->stepdirection, -1.0);
91: /* This may be incorrect; linesearch has values for stepmax and stepmin
92: that should be reset. */
93: TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
94: TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);
95: TaoAddLineSearchCounts(tao);
97: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
98: tao->reason = TAO_DIVERGED_LS_FAILURE;
99: break;
100: }
101: }
103: /* Check for converged */
104: VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);
105: TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);
108: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
109: tao->niter++;
110: TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);
111: TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);
112: (*tao->ops->convergencetest)(tao,tao->cnvP);
113: }
114: return(0);
115: }
117: static PetscErrorCode TaoSetup_BLMVM(Tao tao)
118: {
119: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
123: /* Existence of tao->solution checked in TaoSetup() */
124: VecDuplicate(tao->solution,&blmP->Xold);
125: VecDuplicate(tao->solution,&blmP->Gold);
126: VecDuplicate(tao->solution, &blmP->unprojected_gradient);
128: if (!tao->stepdirection) {
129: VecDuplicate(tao->solution, &tao->stepdirection);
130: }
131: if (!tao->gradient) {
132: VecDuplicate(tao->solution,&tao->gradient);
133: }
134: if (!tao->XL) {
135: VecDuplicate(tao->solution,&tao->XL);
136: VecSet(tao->XL,PETSC_NINFINITY);
137: }
138: if (!tao->XU) {
139: VecDuplicate(tao->solution,&tao->XU);
140: VecSet(tao->XU,PETSC_INFINITY);
141: }
142: /* Allocate matrix for the limited memory approximation */
143: MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);
145: /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
146: if (blmP->H0) {
147: MatLMVMSetJ0(blmP->M, blmP->H0);
148: }
149: return(0);
150: }
152: /* ---------------------------------------------------------- */
153: static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
154: {
155: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
159: if (tao->setupcalled) {
160: VecDestroy(&blmP->unprojected_gradient);
161: VecDestroy(&blmP->Xold);
162: VecDestroy(&blmP->Gold);
163: }
164: MatDestroy(&blmP->M);
165: if (blmP->H0) {
166: PetscObjectDereference((PetscObject)blmP->H0);
167: }
168: PetscFree(tao->data);
169: return(0);
170: }
172: /*------------------------------------------------------------*/
173: static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
174: {
175: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
177: PetscBool is_spd;
180: PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");
181: PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);
182: PetscOptionsTail();
183: TaoLineSearchSetFromOptions(tao->linesearch);
184: MatSetFromOptions(blmP->M);
185: MatGetOption(blmP->M, MAT_SPD, &is_spd);
186: if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
187: return(0);
188: }
191: /*------------------------------------------------------------*/
192: static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
193: {
194: TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data;
195: PetscBool isascii;
199: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
200: if (isascii) {
201: PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);
202: PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);
203: MatView(lmP->M, viewer);
204: PetscViewerPopFormat(viewer);
205: }
206: return(0);
207: }
209: static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
210: {
211: TAO_BLMVM *blm = (TAO_BLMVM *) tao->data;
218: if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
220: VecCopy(tao->gradient,DXL);
221: VecAXPY(DXL,-1.0,blm->unprojected_gradient);
222: VecSet(DXU,0.0);
223: VecPointwiseMax(DXL,DXL,DXU);
225: VecCopy(blm->unprojected_gradient,DXU);
226: VecAXPY(DXU,-1.0,tao->gradient);
227: VecAXPY(DXU,1.0,DXL);
228: return(0);
229: }
231: /* ---------------------------------------------------------- */
232: /*MC
233: TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
234: for nonlinear minimization with bound constraints. It is an extension
235: of TAOLMVM
237: Options Database Keys:
238: . -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls
240: Level: beginner
241: M*/
242: PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
243: {
244: TAO_BLMVM *blmP;
245: const char *morethuente_type = TAOLINESEARCHMT;
249: tao->ops->setup = TaoSetup_BLMVM;
250: tao->ops->solve = TaoSolve_BLMVM;
251: tao->ops->view = TaoView_BLMVM;
252: tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
253: tao->ops->destroy = TaoDestroy_BLMVM;
254: tao->ops->computedual = TaoComputeDual_BLMVM;
256: PetscNewLog(tao,&blmP);
257: blmP->H0 = NULL;
258: blmP->recycle = PETSC_FALSE;
259: tao->data = (void*)blmP;
261: /* Override default settings (unless already changed) */
262: if (!tao->max_it_changed) tao->max_it = 2000;
263: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
265: TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
266: PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
267: TaoLineSearchSetType(tao->linesearch, morethuente_type);
268: TaoLineSearchUseTaoRoutines(tao->linesearch,tao);
269: TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
270:
271: KSPInitializePackage();
272: MatCreate(((PetscObject)tao)->comm, &blmP->M);
273: MatSetType(blmP->M, MATLMVMBFGS);
274: PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);
275: MatSetOptionsPrefix(blmP->M, "tao_blmvm_");
276: return(0);
277: }
279: PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg)
280: {
281: TAO_LMVM *lmP;
282: TAO_BLMVM *blmP;
283: TaoType type;
284: PetscBool is_lmvm, is_blmvm;
286:
288: TaoGetType(tao, &type);
289: PetscStrcmp(type, TAOLMVM, &is_lmvm);
290: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
291:
292: if (is_lmvm) {
293: lmP = (TAO_LMVM *)tao->data;
294: lmP->recycle = flg;
295: } else if (is_blmvm) {
296: blmP = (TAO_BLMVM *)tao->data;
297: blmP->recycle = flg;
298: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
299: return(0);
300: }
302: PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
303: {
304: TAO_LMVM *lmP;
305: TAO_BLMVM *blmP;
306: TaoType type;
307: PetscBool is_lmvm, is_blmvm;
309:
311: TaoGetType(tao, &type);
312: PetscStrcmp(type, TAOLMVM, &is_lmvm);
313: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
315: if (is_lmvm) {
316: lmP = (TAO_LMVM *)tao->data;
317: PetscObjectReference((PetscObject)H0);
318: lmP->H0 = H0;
319: } else if (is_blmvm) {
320: blmP = (TAO_BLMVM *)tao->data;
321: PetscObjectReference((PetscObject)H0);
322: blmP->H0 = H0;
323: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
324: return(0);
325: }
327: PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
328: {
329: TAO_LMVM *lmP;
330: TAO_BLMVM *blmP;
331: TaoType type;
332: PetscBool is_lmvm, is_blmvm;
333: Mat M;
336:
338: TaoGetType(tao, &type);
339: PetscStrcmp(type, TAOLMVM, &is_lmvm);
340: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
342: if (is_lmvm) {
343: lmP = (TAO_LMVM *)tao->data;
344: M = lmP->M;
345: } else if (is_blmvm) {
346: blmP = (TAO_BLMVM *)tao->data;
347: M = blmP->M;
348: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
349: MatLMVMGetJ0(M, H0);
350: return(0);
351: }
353: PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
354: {
355: TAO_LMVM *lmP;
356: TAO_BLMVM *blmP;
357: TaoType type;
358: PetscBool is_lmvm, is_blmvm;
359: Mat M;
362: TaoGetType(tao, &type);
363: PetscStrcmp(type, TAOLMVM, &is_lmvm);
364: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
366: if (is_lmvm) {
367: lmP = (TAO_LMVM *)tao->data;
368: M = lmP->M;
369: } else if (is_blmvm) {
370: blmP = (TAO_BLMVM *)tao->data;
371: M = blmP->M;
372: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
373: MatLMVMGetJ0KSP(M, ksp);
374: return(0);
375: }