Actual source code: owlqn.c
petsc-dev 2014-02-02
1: #include <petsctaolinesearch.h>
2: #include <../src/tao/matrix/lmvmmat.h>
3: #include <../src/tao/unconstrained/impls/owlqn/owlqn.h>
5: #define OWLQN_BFGS 0
6: #define OWLQN_SCALED_GRADIENT 1
7: #define OWLQN_GRADIENT 2
11: static PetscErrorCode ProjDirect_OWLQN(Vec d, Vec g)
12: {
14: PetscReal *gptr,*dptr;
15: PetscInt low,high,low1,high1,i;
18: ierr=VecGetOwnershipRange(d,&low,&high);
19: ierr=VecGetOwnershipRange(g,&low1,&high1);
21: VecGetArray(g,&gptr);
22: VecGetArray(d,&dptr);
23: for (i = 0; i < high-low; i++) {
24: if (dptr[i] * gptr[i] <= 0.0 ) {
25: dptr[i] = 0.0;
26: }
27: }
28: VecRestoreArray(d,&dptr);
29: VecRestoreArray(g,&gptr);
30: return(0);
31: }
35: static PetscErrorCode ComputePseudoGrad_OWLQN(Vec x, Vec gv, PetscReal lambda)
36: {
38: PetscReal *xptr,*gptr;
39: PetscInt low,high,low1,high1,i;
42: ierr=VecGetOwnershipRange(x,&low,&high);
43: ierr=VecGetOwnershipRange(gv,&low1,&high1);
45: VecGetArray(x,&xptr);
46: VecGetArray(gv,&gptr);
47: for (i = 0; i < high-low; i++) {
48: if (xptr[i] < 0.0) gptr[i] = gptr[i] - lambda;
49: else if (xptr[i] > 0.0) gptr[i] = gptr[i] + lambda;
50: else if (gptr[i] + lambda < 0.0) gptr[i] = gptr[i] + lambda;
51: else if (gptr[i] - lambda > 0.0) gptr[i] = gptr[i] - lambda;
52: else gptr[i] = 0.0;
53: }
54: VecRestoreArray(gv,&gptr);
55: VecRestoreArray(x,&xptr);
56: return(0);
57: }
61: static PetscErrorCode TaoSolve_OWLQN(Tao tao)
62: {
63: TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data;
64: PetscReal f, fold, gdx, gnorm;
65: PetscReal step = 1.0;
66: PetscReal delta;
67: PetscErrorCode ierr;
68: PetscInt stepType;
69: PetscInt iter = 0;
70: TaoTerminationReason reason = TAO_CONTINUE_ITERATING;
71: TaoLineSearchTerminationReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
75: if (tao->XL || tao->XU || tao->ops->computebounds) {
76: PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by owlqn algorithm\n");
77: }
79: /* Check convergence criteria */
80: TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);
82: VecCopy(tao->gradient, lmP->GV);
84: ComputePseudoGrad_OWLQN(tao->solution,lmP->GV,lmP->lambda);
86: VecNorm(lmP->GV,NORM_2,&gnorm);
88: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
90: TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);
91: if (reason != TAO_CONTINUE_ITERATING) return(0);
93: /* Set initial scaling for the function */
94: if (f != 0.0) {
95: delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
96: } else {
97: delta = 2.0 / (gnorm*gnorm);
98: }
99: MatLMVMSetDelta(lmP->M,delta);
101: /* Set counter for gradient/reset steps */
102: lmP->bfgs = 0;
103: lmP->sgrad = 0;
104: lmP->grad = 0;
106: /* Have not converged; continue with Newton method */
107: while (reason == TAO_CONTINUE_ITERATING) {
108: /* Compute direction */
109: MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);
110: MatLMVMSolve(lmP->M, lmP->GV, lmP->D);
112: ProjDirect_OWLQN(lmP->D,lmP->GV);
114: ++lmP->bfgs;
116: /* Check for success (descent direction) */
117: VecDot(lmP->D, lmP->GV , &gdx);
118: if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
120: /* Step is not descent or direction produced not a number
121: We can assert bfgsUpdates > 1 in this case because
122: the first solve produces the scaled gradient direction,
123: which is guaranteed to be descent
125: Use steepest descent direction (scaled) */
126: ++lmP->grad;
128: if (f != 0.0) {
129: delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
130: } else {
131: delta = 2.0 / (gnorm*gnorm);
132: }
133: MatLMVMSetDelta(lmP->M, delta);
134: MatLMVMReset(lmP->M);
135: MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);
136: MatLMVMSolve(lmP->M,lmP->GV, lmP->D);
138: ProjDirect_OWLQN(lmP->D,lmP->GV);
140: lmP->bfgs = 1;
141: ++lmP->sgrad;
142: stepType = OWLQN_SCALED_GRADIENT;
143: } else {
144: if (1 == lmP->bfgs) {
145: /* The first BFGS direction is always the scaled gradient */
146: ++lmP->sgrad;
147: stepType = OWLQN_SCALED_GRADIENT;
148: } else {
149: ++lmP->bfgs;
150: stepType = OWLQN_BFGS;
151: }
152: }
154: VecScale(lmP->D, -1.0);
156: /* Perform the linesearch */
157: fold = f;
158: VecCopy(tao->solution, lmP->Xold);
159: VecCopy(tao->gradient, lmP->Gold);
161: TaoLineSearchApply(tao->linesearch, tao->solution, &f, lmP->GV, lmP->D, &step,&ls_status);
162: TaoAddLineSearchCounts(tao);
164: while (((int)ls_status < 0) && (stepType != OWLQN_GRADIENT)) {
166: /* Reset factors and use scaled gradient step */
167: f = fold;
168: VecCopy(lmP->Xold, tao->solution);
169: VecCopy(lmP->Gold, tao->gradient);
170: VecCopy(tao->gradient, lmP->GV);
172: ComputePseudoGrad_OWLQN(tao->solution,lmP->GV,lmP->lambda);
174: switch(stepType) {
175: case OWLQN_BFGS:
176: /* Failed to obtain acceptable iterate with BFGS step
177: Attempt to use the scaled gradient direction */
179: if (f != 0.0) {
180: delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
181: } else {
182: delta = 2.0 / (gnorm*gnorm);
183: }
184: MatLMVMSetDelta(lmP->M, delta);
185: MatLMVMReset(lmP->M);
186: MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);
187: MatLMVMSolve(lmP->M, lmP->GV, lmP->D);
189: ProjDirect_OWLQN(lmP->D,lmP->GV);
191: lmP->bfgs = 1;
192: ++lmP->sgrad;
193: stepType = OWLQN_SCALED_GRADIENT;
194: break;
196: case OWLQN_SCALED_GRADIENT:
197: /* The scaled gradient step did not produce a new iterate;
198: attempt to use the gradient direction.
199: Need to make sure we are not using a different diagonal scaling */
200: MatLMVMSetDelta(lmP->M, 1.0);
201: MatLMVMReset(lmP->M);
202: MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);
203: MatLMVMSolve(lmP->M, lmP->GV, lmP->D);
205: ProjDirect_OWLQN(lmP->D,lmP->GV);
207: lmP->bfgs = 1;
208: ++lmP->grad;
209: stepType = OWLQN_GRADIENT;
210: break;
211: }
212: VecScale(lmP->D, -1.0);
215: /* Perform the linesearch */
216: TaoLineSearchApply(tao->linesearch, tao->solution, &f, lmP->GV, lmP->D, &step, &ls_status);
217: TaoAddLineSearchCounts(tao);
218: }
220: if ((int)ls_status < 0) {
221: /* Failed to find an improving point*/
222: f = fold;
223: VecCopy(lmP->Xold, tao->solution);
224: VecCopy(lmP->Gold, tao->gradient);
225: VecCopy(tao->gradient, lmP->GV);
226: step = 0.0;
227: } else {
228: /* a little hack here, because that gv is used to store g */
229: VecCopy(lmP->GV, tao->gradient);
230: }
232: ComputePseudoGrad_OWLQN(tao->solution,lmP->GV,lmP->lambda);
234: /* Check for termination */
236: VecNorm(lmP->GV,NORM_2,&gnorm);
238: iter++;
239: TaoMonitor(tao,iter,f,gnorm,0.0,step,&reason);
241: if ((int)ls_status < 0) break;
242: }
243: return(0);
244: }
248: static PetscErrorCode TaoSetUp_OWLQN(Tao tao)
249: {
250: TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data;
251: PetscInt n,N;
255: /* Existence of tao->solution checked in TaoSetUp() */
256: if (!tao->gradient) {VecDuplicate(tao->solution,&tao->gradient); }
257: if (!tao->stepdirection) {VecDuplicate(tao->solution,&tao->stepdirection); }
258: if (!lmP->D) {VecDuplicate(tao->solution,&lmP->D); }
259: if (!lmP->GV) {VecDuplicate(tao->solution,&lmP->GV); }
260: if (!lmP->Xold) {VecDuplicate(tao->solution,&lmP->Xold); }
261: if (!lmP->Gold) {VecDuplicate(tao->solution,&lmP->Gold); }
263: /* Create matrix for the limited memory approximation */
264: VecGetLocalSize(tao->solution,&n);
265: VecGetSize(tao->solution,&N);
266: MatCreateLMVM(((PetscObject)tao)->comm,n,N,&lmP->M);
267: MatLMVMAllocateVectors(lmP->M,tao->solution);
268: return(0);
269: }
271: /* ---------------------------------------------------------- */
274: static PetscErrorCode TaoDestroy_OWLQN(Tao tao)
275: {
276: TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data;
280: if (tao->setupcalled) {
281: VecDestroy(&lmP->Xold);
282: VecDestroy(&lmP->Gold);
283: VecDestroy(&lmP->D);
284: MatDestroy(&lmP->M);
285: VecDestroy(&lmP->GV);
286: }
287: PetscFree(tao->data);
288: return(0);
289: }
291: /*------------------------------------------------------------*/
294: static PetscErrorCode TaoSetFromOptions_OWLQN(Tao tao)
295: {
296: TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data;
297: PetscBool flg;
301: PetscOptionsHead("Orthant-Wise Limited-memory method for Quasi-Newton unconstrained optimization");
302: PetscOptionsReal("-tao_owlqn_lambda", "regulariser weight","", 100,&lmP->lambda,&flg);
303: PetscOptionsTail();
304: TaoLineSearchSetFromOptions(tao->linesearch);
305: return(0);
306: }
308: /*------------------------------------------------------------*/
311: static PetscErrorCode TaoView_OWLQN(Tao tao, PetscViewer viewer)
312: {
313: TAO_OWLQN *lm = (TAO_OWLQN *)tao->data;
314: PetscBool isascii;
318: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
319: if (isascii) {
320: PetscViewerASCIIPushTab(viewer);
321: PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);
322: PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);
323: PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);
324: PetscViewerASCIIPopTab(viewer);
325: }
326: return(0);
327: }
329: /* ---------------------------------------------------------- */
331: EXTERN_C_BEGIN
334: PetscErrorCode TaoCreate_OWLQN(Tao tao)
335: {
336: TAO_OWLQN *lmP;
337: const char *owarmijo_type = TAOLINESEARCH_OWARMIJO;
341: tao->ops->setup = TaoSetUp_OWLQN;
342: tao->ops->solve = TaoSolve_OWLQN;
343: tao->ops->view = TaoView_OWLQN;
344: tao->ops->setfromoptions = TaoSetFromOptions_OWLQN;
345: tao->ops->destroy = TaoDestroy_OWLQN;
347: PetscNewLog(tao,&lmP);
348: lmP->D = 0;
349: lmP->M = 0;
350: lmP->GV = 0;
351: lmP->Xold = 0;
352: lmP->Gold = 0;
353: lmP->lambda = 1.0;
355: tao->data = (void*)lmP;
356: tao->max_it = 2000;
357: tao->max_funcs = 4000;
358: tao->fatol = 1e-4;
359: tao->frtol = 1e-4;
361: TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);
362: TaoLineSearchSetType(tao->linesearch,owarmijo_type);
363: TaoLineSearchUseTaoRoutines(tao->linesearch,tao);
364: return(0);
365: }
366: EXTERN_C_END