Actual source code: blmvm.c
petsc-3.8.4 2018-03-24
1: #include <petsctaolinesearch.h>
2: #include <../src/tao/matrix/lmvmmat.h>
3: #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
4: #include <../src/tao/bound/impls/blmvm/blmvm.h>
6: /*------------------------------------------------------------*/
7: static PetscErrorCode TaoSolve_BLMVM(Tao tao)
8: {
9: PetscErrorCode ierr;
10: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
11: TaoConvergedReason reason = TAO_CONTINUE_ITERATING;
12: TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
13: PetscReal f, fold, gdx, gnorm;
14: PetscReal stepsize = 1.0,delta;
17: /* Project initial point onto bounds */
18: TaoComputeVariableBounds(tao);
19: VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
20: TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);
23: /* Check convergence criteria */
24: TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);
25: VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);
27: TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);
28: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
30: TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);
31: if (reason != TAO_CONTINUE_ITERATING) return(0);
33: /* Set initial scaling for the function */
34: if (f != 0.0) {
35: delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
36: } else {
37: delta = 2.0 / (gnorm*gnorm);
38: }
39: MatLMVMSetDelta(blmP->M,delta);
40: MatLMVMReset(blmP->M);
42: /* Set counter for gradient/reset steps */
43: blmP->grad = 0;
44: blmP->reset = 0;
46: /* Have not converged; continue with Newton method */
47: while (reason == TAO_CONTINUE_ITERATING) {
48: /* Compute direction */
49: MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);
50: MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
51: VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);
53: /* Check for success (descent direction) */
54: VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);
55: if (gdx <= 0) {
56: /* Step is not descent or solve was not successful
57: Use steepest descent direction (scaled) */
58: ++blmP->grad;
60: if (f != 0.0) {
61: delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
62: } else {
63: delta = 2.0 / (gnorm*gnorm);
64: }
65: MatLMVMSetDelta(blmP->M,delta);
66: MatLMVMReset(blmP->M);
67: MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
68: MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);
69: }
70: VecScale(tao->stepdirection,-1.0);
72: /* Perform the linesearch */
73: fold = f;
74: VecCopy(tao->solution, blmP->Xold);
75: VecCopy(blmP->unprojected_gradient, blmP->Gold);
76: TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
77: TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);
78: TaoAddLineSearchCounts(tao);
80: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
81: /* Linesearch failed
82: Reset factors and use scaled (projected) gradient step */
83: ++blmP->reset;
85: f = fold;
86: VecCopy(blmP->Xold, tao->solution);
87: VecCopy(blmP->Gold, blmP->unprojected_gradient);
89: if (f != 0.0) {
90: delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
91: } else {
92: delta = 2.0/ (gnorm*gnorm);
93: }
94: MatLMVMSetDelta(blmP->M,delta);
95: MatLMVMReset(blmP->M);
96: MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
97: MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
98: VecScale(tao->stepdirection, -1.0);
100: /* This may be incorrect; linesearch has values for stepmax and stepmin
101: that should be reset. */
102: TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
103: TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);
104: TaoAddLineSearchCounts(tao);
106: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
107: tao->reason = TAO_DIVERGED_LS_FAILURE;
108: break;
109: }
110: }
112: /* Check for converged */
113: VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);
114: TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);
117: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
118: tao->niter++;
119: TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);
120: }
121: return(0);
122: }
124: static PetscErrorCode TaoSetup_BLMVM(Tao tao)
125: {
126: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
127: PetscInt n,N;
129: KSP H0ksp;
132: /* Existence of tao->solution checked in TaoSetup() */
133: VecDuplicate(tao->solution,&blmP->Xold);
134: VecDuplicate(tao->solution,&blmP->Gold);
135: VecDuplicate(tao->solution, &blmP->unprojected_gradient);
137: if (!tao->stepdirection) {
138: VecDuplicate(tao->solution, &tao->stepdirection);
139: }
140: if (!tao->gradient) {
141: VecDuplicate(tao->solution,&tao->gradient);
142: }
143: if (!tao->XL) {
144: VecDuplicate(tao->solution,&tao->XL);
145: VecSet(tao->XL,PETSC_NINFINITY);
146: }
147: if (!tao->XU) {
148: VecDuplicate(tao->solution,&tao->XU);
149: VecSet(tao->XU,PETSC_INFINITY);
150: }
151: /* Create matrix for the limited memory approximation */
152: VecGetLocalSize(tao->solution,&n);
153: VecGetSize(tao->solution,&N);
154: MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);
155: MatLMVMAllocateVectors(blmP->M,tao->solution);
157: /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
158: if (blmP->H0) {
159: const char *prefix;
160: MatLMVMSetH0(blmP->M, blmP->H0);
161: MatLMVMGetH0KSP(blmP->M, &H0ksp);
163: TaoGetOptionsPrefix(tao, &prefix);
164: KSPSetOptionsPrefix(H0ksp, prefix);
165: KSPAppendOptionsPrefix(H0ksp, "tao_h0_");
166: KSPSetFromOptions(H0ksp);
167: KSPSetUp(H0ksp);
168: }
169: return(0);
170: }
172: /* ---------------------------------------------------------- */
173: static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
174: {
175: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
179: if (tao->setupcalled) {
180: MatDestroy(&blmP->M);
181: VecDestroy(&blmP->unprojected_gradient);
182: VecDestroy(&blmP->Xold);
183: VecDestroy(&blmP->Gold);
184: }
186: if (blmP->H0) {
187: PetscObjectDereference((PetscObject)blmP->H0);
188: }
190: PetscFree(tao->data);
191: return(0);
192: }
194: /*------------------------------------------------------------*/
195: static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
196: {
200: PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");
201: TaoLineSearchSetFromOptions(tao->linesearch);
202: PetscOptionsTail();
203: return(0);
204: }
207: /*------------------------------------------------------------*/
208: static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
209: {
210: TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data;
211: PetscBool isascii;
215: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
216: if (isascii) {
217: PetscViewerASCIIPushTab(viewer);
218: PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);
219: PetscViewerASCIIPopTab(viewer);
220: }
221: return(0);
222: }
224: static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
225: {
226: TAO_BLMVM *blm = (TAO_BLMVM *) tao->data;
233: if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
235: VecCopy(tao->gradient,DXL);
236: VecAXPY(DXL,-1.0,blm->unprojected_gradient);
237: VecSet(DXU,0.0);
238: VecPointwiseMax(DXL,DXL,DXU);
240: VecCopy(blm->unprojected_gradient,DXU);
241: VecAXPY(DXU,-1.0,tao->gradient);
242: VecAXPY(DXU,1.0,DXL);
243: return(0);
244: }
246: /* ---------------------------------------------------------- */
247: /*MC
248: TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
249: for nonlinear minimization with bound constraints. It is an extension
250: of TAOLMVM
252: Options Database Keys:
253: + -tao_lmm_vectors - number of vectors to use for approximation
254: . -tao_lmm_scale_type - "none","scalar","broyden"
255: . -tao_lmm_limit_type - "none","average","relative","absolute"
256: . -tao_lmm_rescale_type - "none","scalar","gl"
257: . -tao_lmm_limit_mu - mu limiting factor
258: . -tao_lmm_limit_nu - nu limiting factor
259: . -tao_lmm_delta_min - minimum delta value
260: . -tao_lmm_delta_max - maximum delta value
261: . -tao_lmm_broyden_phi - phi factor for Broyden scaling
262: . -tao_lmm_scalar_alpha - alpha factor for scalar scaling
263: . -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
264: . -tao_lmm_rescale_beta - beta factor for rescaling diagonal
265: . -tao_lmm_scalar_history - amount of history for scalar scaling
266: . -tao_lmm_rescale_history - amount of history for rescaling diagonal
267: - -tao_lmm_eps - rejection tolerance
269: Level: beginner
270: M*/
271: PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
272: {
273: TAO_BLMVM *blmP;
274: const char *morethuente_type = TAOLINESEARCHMT;
278: tao->ops->setup = TaoSetup_BLMVM;
279: tao->ops->solve = TaoSolve_BLMVM;
280: tao->ops->view = TaoView_BLMVM;
281: tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
282: tao->ops->destroy = TaoDestroy_BLMVM;
283: tao->ops->computedual = TaoComputeDual_BLMVM;
285: PetscNewLog(tao,&blmP);
286: blmP->H0 = NULL;
287: tao->data = (void*)blmP;
289: /* Override default settings (unless already changed) */
290: if (!tao->max_it_changed) tao->max_it = 2000;
291: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
293: TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
294: TaoLineSearchSetType(tao->linesearch, morethuente_type);
295: TaoLineSearchUseTaoRoutines(tao->linesearch,tao);
296: TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
297: return(0);
298: }
300: PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
301: {
302: TAO_LMVM *lmP;
303: TAO_BLMVM *blmP;
304: const TaoType type;
305: PetscBool is_lmvm, is_blmvm;
308: TaoGetType(tao, &type);
309: PetscStrcmp(type, TAOLMVM, &is_lmvm);
310: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
312: if (is_lmvm) {
313: lmP = (TAO_LMVM *)tao->data;
314: PetscObjectReference((PetscObject)H0);
315: lmP->H0 = H0;
316: } else if (is_blmvm) {
317: blmP = (TAO_BLMVM *)tao->data;
318: PetscObjectReference((PetscObject)H0);
319: blmP->H0 = H0;
320: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
321: return(0);
322: }
324: PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
325: {
326: TAO_LMVM *lmP;
327: TAO_BLMVM *blmP;
328: const TaoType type;
329: PetscBool is_lmvm, is_blmvm;
330: Mat M;
334: TaoGetType(tao, &type);
335: PetscStrcmp(type, TAOLMVM, &is_lmvm);
336: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
338: if (is_lmvm) {
339: lmP = (TAO_LMVM *)tao->data;
340: M = lmP->M;
341: } else if (is_blmvm) {
342: blmP = (TAO_BLMVM *)tao->data;
343: M = blmP->M;
344: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
345: MatLMVMGetH0(M, H0);
346: return(0);
347: }
349: PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
350: {
351: TAO_LMVM *lmP;
352: TAO_BLMVM *blmP;
353: const TaoType type;
354: PetscBool is_lmvm, is_blmvm;
355: Mat M;
358: TaoGetType(tao, &type);
359: PetscStrcmp(type, TAOLMVM, &is_lmvm);
360: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
362: if (is_lmvm) {
363: lmP = (TAO_LMVM *)tao->data;
364: M = lmP->M;
365: } else if (is_blmvm) {
366: blmP = (TAO_BLMVM *)tao->data;
367: M = blmP->M;
368: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
369: MatLMVMGetH0KSP(M, ksp);
370: return(0);
371: }