Actual source code: blmvm.c
petsc-3.9.4 2018-09-11
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: TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
12: PetscReal f, fold, gdx, gnorm;
13: PetscReal stepsize = 1.0,delta;
16: /* Project initial point onto bounds */
17: TaoComputeVariableBounds(tao);
18: VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
19: TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);
22: /* Check convergence criteria */
23: TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);
24: VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);
26: TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);
27: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
29: tao->reason = TAO_CONTINUE_ITERATING;
30: TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);
31: TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);
32: (*tao->ops->convergencetest)(tao,tao->cnvP);
33: if (tao->reason != TAO_CONTINUE_ITERATING) return(0);
35: /* Set initial scaling for the function */
36: if (f != 0.0) {
37: delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
38: } else {
39: delta = 2.0 / (gnorm*gnorm);
40: }
41: MatLMVMSetDelta(blmP->M,delta);
42: MatLMVMReset(blmP->M);
44: /* Set counter for gradient/reset steps */
45: blmP->grad = 0;
46: blmP->reset = 0;
48: /* Have not converged; continue with Newton method */
49: while (tao->reason == TAO_CONTINUE_ITERATING) {
50: /* Compute direction */
51: MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);
52: MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
53: VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);
55: /* Check for success (descent direction) */
56: VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);
57: if (gdx <= 0) {
58: /* Step is not descent or solve was not successful
59: Use steepest descent direction (scaled) */
60: ++blmP->grad;
62: if (f != 0.0) {
63: delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
64: } else {
65: delta = 2.0 / (gnorm*gnorm);
66: }
67: MatLMVMSetDelta(blmP->M,delta);
68: MatLMVMReset(blmP->M);
69: MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
70: MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);
71: }
72: VecScale(tao->stepdirection,-1.0);
74: /* Perform the linesearch */
75: fold = f;
76: VecCopy(tao->solution, blmP->Xold);
77: VecCopy(blmP->unprojected_gradient, blmP->Gold);
78: TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
79: TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);
80: TaoAddLineSearchCounts(tao);
82: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
83: /* Linesearch failed
84: Reset factors and use scaled (projected) gradient step */
85: ++blmP->reset;
87: f = fold;
88: VecCopy(blmP->Xold, tao->solution);
89: VecCopy(blmP->Gold, blmP->unprojected_gradient);
91: if (f != 0.0) {
92: delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
93: } else {
94: delta = 2.0/ (gnorm*gnorm);
95: }
96: MatLMVMSetDelta(blmP->M,delta);
97: MatLMVMReset(blmP->M);
98: MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);
99: MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);
100: VecScale(tao->stepdirection, -1.0);
102: /* This may be incorrect; linesearch has values for stepmax and stepmin
103: that should be reset. */
104: TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
105: TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);
106: TaoAddLineSearchCounts(tao);
108: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
109: tao->reason = TAO_DIVERGED_LS_FAILURE;
110: break;
111: }
112: }
114: /* Check for converged */
115: VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);
116: TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);
119: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
120: tao->niter++;
121: TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);
122: TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);
123: (*tao->ops->convergencetest)(tao,tao->cnvP);
124: }
125: return(0);
126: }
128: static PetscErrorCode TaoSetup_BLMVM(Tao tao)
129: {
130: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
131: PetscInt n,N;
133: KSP H0ksp;
136: /* Existence of tao->solution checked in TaoSetup() */
137: VecDuplicate(tao->solution,&blmP->Xold);
138: VecDuplicate(tao->solution,&blmP->Gold);
139: VecDuplicate(tao->solution, &blmP->unprojected_gradient);
141: if (!tao->stepdirection) {
142: VecDuplicate(tao->solution, &tao->stepdirection);
143: }
144: if (!tao->gradient) {
145: VecDuplicate(tao->solution,&tao->gradient);
146: }
147: if (!tao->XL) {
148: VecDuplicate(tao->solution,&tao->XL);
149: VecSet(tao->XL,PETSC_NINFINITY);
150: }
151: if (!tao->XU) {
152: VecDuplicate(tao->solution,&tao->XU);
153: VecSet(tao->XU,PETSC_INFINITY);
154: }
155: /* Create matrix for the limited memory approximation */
156: VecGetLocalSize(tao->solution,&n);
157: VecGetSize(tao->solution,&N);
158: MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);
159: MatLMVMAllocateVectors(blmP->M,tao->solution);
161: /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
162: if (blmP->H0) {
163: const char *prefix;
164: MatLMVMSetH0(blmP->M, blmP->H0);
165: MatLMVMGetH0KSP(blmP->M, &H0ksp);
167: TaoGetOptionsPrefix(tao, &prefix);
168: KSPSetOptionsPrefix(H0ksp, prefix);
169: KSPAppendOptionsPrefix(H0ksp, "tao_h0_");
170: KSPSetFromOptions(H0ksp);
171: KSPSetUp(H0ksp);
172: }
173: return(0);
174: }
176: /* ---------------------------------------------------------- */
177: static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
178: {
179: TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
183: if (tao->setupcalled) {
184: MatDestroy(&blmP->M);
185: VecDestroy(&blmP->unprojected_gradient);
186: VecDestroy(&blmP->Xold);
187: VecDestroy(&blmP->Gold);
188: }
190: if (blmP->H0) {
191: PetscObjectDereference((PetscObject)blmP->H0);
192: }
194: PetscFree(tao->data);
195: return(0);
196: }
198: /*------------------------------------------------------------*/
199: static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
200: {
204: PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");
205: TaoLineSearchSetFromOptions(tao->linesearch);
206: PetscOptionsTail();
207: return(0);
208: }
211: /*------------------------------------------------------------*/
212: static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
213: {
214: TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data;
215: PetscBool isascii;
219: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
220: if (isascii) {
221: PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);
222: }
223: return(0);
224: }
226: static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
227: {
228: TAO_BLMVM *blm = (TAO_BLMVM *) tao->data;
235: if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
237: VecCopy(tao->gradient,DXL);
238: VecAXPY(DXL,-1.0,blm->unprojected_gradient);
239: VecSet(DXU,0.0);
240: VecPointwiseMax(DXL,DXL,DXU);
242: VecCopy(blm->unprojected_gradient,DXU);
243: VecAXPY(DXU,-1.0,tao->gradient);
244: VecAXPY(DXU,1.0,DXL);
245: return(0);
246: }
248: /* ---------------------------------------------------------- */
249: /*MC
250: TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
251: for nonlinear minimization with bound constraints. It is an extension
252: of TAOLMVM
254: Options Database Keys:
255: + -tao_lmm_vectors - number of vectors to use for approximation
256: . -tao_lmm_scale_type - "none","scalar","broyden"
257: . -tao_lmm_limit_type - "none","average","relative","absolute"
258: . -tao_lmm_rescale_type - "none","scalar","gl"
259: . -tao_lmm_limit_mu - mu limiting factor
260: . -tao_lmm_limit_nu - nu limiting factor
261: . -tao_lmm_delta_min - minimum delta value
262: . -tao_lmm_delta_max - maximum delta value
263: . -tao_lmm_broyden_phi - phi factor for Broyden scaling
264: . -tao_lmm_scalar_alpha - alpha factor for scalar scaling
265: . -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
266: . -tao_lmm_rescale_beta - beta factor for rescaling diagonal
267: . -tao_lmm_scalar_history - amount of history for scalar scaling
268: . -tao_lmm_rescale_history - amount of history for rescaling diagonal
269: - -tao_lmm_eps - rejection tolerance
271: Level: beginner
272: M*/
273: PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
274: {
275: TAO_BLMVM *blmP;
276: const char *morethuente_type = TAOLINESEARCHMT;
280: tao->ops->setup = TaoSetup_BLMVM;
281: tao->ops->solve = TaoSolve_BLMVM;
282: tao->ops->view = TaoView_BLMVM;
283: tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
284: tao->ops->destroy = TaoDestroy_BLMVM;
285: tao->ops->computedual = TaoComputeDual_BLMVM;
287: PetscNewLog(tao,&blmP);
288: blmP->H0 = NULL;
289: tao->data = (void*)blmP;
291: /* Override default settings (unless already changed) */
292: if (!tao->max_it_changed) tao->max_it = 2000;
293: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
295: TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
296: PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
297: TaoLineSearchSetType(tao->linesearch, morethuente_type);
298: TaoLineSearchUseTaoRoutines(tao->linesearch,tao);
299: TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
300: return(0);
301: }
303: PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
304: {
305: TAO_LMVM *lmP;
306: TAO_BLMVM *blmP;
307: const TaoType type;
308: PetscBool is_lmvm, is_blmvm;
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: PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
328: {
329: TAO_LMVM *lmP;
330: TAO_BLMVM *blmP;
331: const TaoType type;
332: PetscBool is_lmvm, is_blmvm;
333: Mat M;
337: TaoGetType(tao, &type);
338: PetscStrcmp(type, TAOLMVM, &is_lmvm);
339: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
341: if (is_lmvm) {
342: lmP = (TAO_LMVM *)tao->data;
343: M = lmP->M;
344: } else if (is_blmvm) {
345: blmP = (TAO_BLMVM *)tao->data;
346: M = blmP->M;
347: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
348: MatLMVMGetH0(M, H0);
349: return(0);
350: }
352: PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
353: {
354: TAO_LMVM *lmP;
355: TAO_BLMVM *blmP;
356: const TaoType type;
357: PetscBool is_lmvm, is_blmvm;
358: Mat M;
361: TaoGetType(tao, &type);
362: PetscStrcmp(type, TAOLMVM, &is_lmvm);
363: PetscStrcmp(type, TAOBLMVM, &is_blmvm);
365: if (is_lmvm) {
366: lmP = (TAO_LMVM *)tao->data;
367: M = lmP->M;
368: } else if (is_blmvm) {
369: blmP = (TAO_BLMVM *)tao->data;
370: M = blmP->M;
371: } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
372: MatLMVMGetH0KSP(M, ksp);
373: return(0);
374: }