Actual source code: bfgs.c
1: #include <../src/ksp/ksp/utils/lmvm/symbrdn/symbrdn.h>
2: #include <../src/ksp/ksp/utils/lmvm/diagbrdn/diagbrdn.h>
4: /*
5: Limited-memory Broyden-Fletcher-Goldfarb-Shano method for approximating both
6: the forward product and inverse application of a Jacobian.
7: */
9: /*------------------------------------------------------------*/
11: /*
12: The solution method (approximate inverse Jacobian application) is adapted
13: from Algorithm 7.4 on page 178 of Nocedal and Wright "Numerical Optimization"
14: 2nd edition (https://doi.org/10.1007/978-0-387-40065-5). The initial inverse
15: Jacobian application falls back onto the gamma scaling recommended in equation
16: (7.20) if the user has not provided any estimation of the initial Jacobian or
17: its inverse.
19: work <- F
21: for i = k,k-1,k-2,...,0
22: rho[i] = 1 / (Y[i]^T S[i])
23: alpha[i] = rho[i] * (S[i]^T work)
24: Fwork <- work - (alpha[i] * Y[i])
25: end
27: dX <- J0^{-1} * work
29: for i = 0,1,2,...,k
30: beta = rho[i] * (Y[i]^T dX)
31: dX <- dX + ((alpha[i] - beta) * S[i])
32: end
33: */
34: PetscErrorCode MatSolve_LMVMBFGS(Mat B, Vec F, Vec dX)
35: {
36: Mat_LMVM *lmvm = (Mat_LMVM *)B->data;
37: Mat_SymBrdn *lbfgs = (Mat_SymBrdn *)lmvm->ctx;
38: PetscInt i;
39: PetscReal *alpha, beta;
40: PetscScalar stf, ytx;
42: VecCheckSameSize(F, 2, dX, 3);
43: VecCheckMatCompatible(B, dX, 3, F, 2);
45: /* Copy the function into the work vector for the first loop */
46: VecCopy(F, lbfgs->work);
48: /* Start the first loop */
49: PetscMalloc1(lmvm->k + 1, &alpha);
50: for (i = lmvm->k; i >= 0; --i) {
51: VecDot(lmvm->S[i], lbfgs->work, &stf);
52: alpha[i] = PetscRealPart(stf) / lbfgs->yts[i];
53: VecAXPY(lbfgs->work, -alpha[i], lmvm->Y[i]);
54: }
56: /* Invert the initial Jacobian onto the work vector (or apply scaling) */
57: MatSymBrdnApplyJ0Inv(B, lbfgs->work, dX);
59: /* Start the second loop */
60: for (i = 0; i <= lmvm->k; ++i) {
61: VecDot(lmvm->Y[i], dX, &ytx);
62: beta = PetscRealPart(ytx) / lbfgs->yts[i];
63: VecAXPY(dX, alpha[i] - beta, lmvm->S[i]);
64: }
65: PetscFree(alpha);
66: return 0;
67: }
69: /*------------------------------------------------------------*/
71: /*
72: The forward product for the approximate Jacobian is the matrix-free
73: implementation of Equation (6.19) in Nocedal and Wright "Numerical
74: Optimization" 2nd Edition, pg 140.
76: This forward product has the same structure as the inverse Jacobian
77: application in the DFP formulation, except with S and Y exchanging
78: roles.
80: Note: P[i] = (B_i)*S[i] terms are computed ahead of time whenever
81: the matrix is updated with a new (S[i], Y[i]) pair. This allows
82: repeated calls of MatMult inside KSP solvers without unnecessarily
83: recomputing P[i] terms in expensive nested-loops.
85: Z <- J0 * X
87: for i = 0,1,2,...,k
88: P[i] <- J0 * S[i]
89: for j = 0,1,2,...,(i-1)
90: gamma = (Y[j]^T S[i]) / (Y[j]^T S[j])
91: zeta = (S[j]^ P[i]) / (S[j]^T P[j])
92: P[i] <- P[i] - (zeta * P[j]) + (gamma * Y[j])
93: end
94: gamma = (Y[i]^T X) / (Y[i]^T S[i])
95: zeta = (S[i]^T Z) / (S[i]^T P[i])
96: Z <- Z - (zeta * P[i]) + (gamma * Y[i])
97: end
98: */
99: PetscErrorCode MatMult_LMVMBFGS(Mat B, Vec X, Vec Z)
100: {
101: Mat_LMVM *lmvm = (Mat_LMVM *)B->data;
102: Mat_SymBrdn *lbfgs = (Mat_SymBrdn *)lmvm->ctx;
103: PetscInt i, j;
104: PetscScalar sjtpi, yjtsi, ytx, stz, stp;
106: VecCheckSameSize(X, 2, Z, 3);
107: VecCheckMatCompatible(B, X, 2, Z, 3);
109: if (lbfgs->needP) {
110: /* Pre-compute (P[i] = B_i * S[i]) */
111: for (i = 0; i <= lmvm->k; ++i) {
112: MatSymBrdnApplyJ0Fwd(B, lmvm->S[i], lbfgs->P[i]);
113: for (j = 0; j <= i - 1; ++j) {
114: /* Compute the necessary dot products */
115: VecDotBegin(lmvm->S[j], lbfgs->P[i], &sjtpi);
116: VecDotBegin(lmvm->Y[j], lmvm->S[i], &yjtsi);
117: VecDotEnd(lmvm->S[j], lbfgs->P[i], &sjtpi);
118: VecDotEnd(lmvm->Y[j], lmvm->S[i], &yjtsi);
119: /* Compute the pure BFGS component of the forward product */
120: VecAXPBYPCZ(lbfgs->P[i], -PetscRealPart(sjtpi) / lbfgs->stp[j], PetscRealPart(yjtsi) / lbfgs->yts[j], 1.0, lbfgs->P[j], lmvm->Y[j]);
121: }
122: VecDot(lmvm->S[i], lbfgs->P[i], &stp);
123: lbfgs->stp[i] = PetscRealPart(stp);
124: }
125: lbfgs->needP = PETSC_FALSE;
126: }
128: /* Start the outer loop (i) for the recursive formula */
129: MatSymBrdnApplyJ0Fwd(B, X, Z);
130: for (i = 0; i <= lmvm->k; ++i) {
131: /* Get all the dot products we need */
132: VecDotBegin(lmvm->S[i], Z, &stz);
133: VecDotBegin(lmvm->Y[i], X, &ytx);
134: VecDotEnd(lmvm->S[i], Z, &stz);
135: VecDotEnd(lmvm->Y[i], X, &ytx);
136: /* Update Z_{i+1} = B_{i+1} * X */
137: VecAXPBYPCZ(Z, -PetscRealPart(stz) / lbfgs->stp[i], PetscRealPart(ytx) / lbfgs->yts[i], 1.0, lbfgs->P[i], lmvm->Y[i]);
138: }
139: return 0;
140: }
142: /*------------------------------------------------------------*/
144: static PetscErrorCode MatUpdate_LMVMBFGS(Mat B, Vec X, Vec F)
145: {
146: Mat_LMVM *lmvm = (Mat_LMVM *)B->data;
147: Mat_SymBrdn *lbfgs = (Mat_SymBrdn *)lmvm->ctx;
148: Mat_LMVM *dbase;
149: Mat_DiagBrdn *dctx;
150: PetscInt old_k, i;
151: PetscReal curvtol;
152: PetscScalar curvature, ytytmp, ststmp;
154: if (!lmvm->m) return 0;
155: if (lmvm->prev_set) {
156: /* Compute the new (S = X - Xprev) and (Y = F - Fprev) vectors */
157: VecAYPX(lmvm->Xprev, -1.0, X);
158: VecAYPX(lmvm->Fprev, -1.0, F);
159: /* Test if the updates can be accepted */
160: VecDotBegin(lmvm->Xprev, lmvm->Fprev, &curvature);
161: VecDotBegin(lmvm->Xprev, lmvm->Xprev, &ststmp);
162: VecDotEnd(lmvm->Xprev, lmvm->Fprev, &curvature);
163: VecDotEnd(lmvm->Xprev, lmvm->Xprev, &ststmp);
164: if (PetscRealPart(ststmp) < lmvm->eps) {
165: curvtol = 0.0;
166: } else {
167: curvtol = lmvm->eps * PetscRealPart(ststmp);
168: }
169: if (PetscRealPart(curvature) > curvtol) {
170: /* Update is good, accept it */
171: lbfgs->watchdog = 0;
172: lbfgs->needP = PETSC_TRUE;
173: old_k = lmvm->k;
174: MatUpdateKernel_LMVM(B, lmvm->Xprev, lmvm->Fprev);
175: /* If we hit the memory limit, shift the yts, yty and sts arrays */
176: if (old_k == lmvm->k) {
177: for (i = 0; i <= lmvm->k - 1; ++i) {
178: lbfgs->yts[i] = lbfgs->yts[i + 1];
179: lbfgs->yty[i] = lbfgs->yty[i + 1];
180: lbfgs->sts[i] = lbfgs->sts[i + 1];
181: }
182: }
183: /* Update history of useful scalars */
184: VecDot(lmvm->Y[lmvm->k], lmvm->Y[lmvm->k], &ytytmp);
185: lbfgs->yts[lmvm->k] = PetscRealPart(curvature);
186: lbfgs->yty[lmvm->k] = PetscRealPart(ytytmp);
187: lbfgs->sts[lmvm->k] = PetscRealPart(ststmp);
188: /* Compute the scalar scale if necessary */
189: if (lbfgs->scale_type == MAT_LMVM_SYMBROYDEN_SCALE_SCALAR) MatSymBrdnComputeJ0Scalar(B);
190: } else {
191: /* Update is bad, skip it */
192: ++lmvm->nrejects;
193: ++lbfgs->watchdog;
194: }
195: } else {
196: switch (lbfgs->scale_type) {
197: case MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL:
198: dbase = (Mat_LMVM *)lbfgs->D->data;
199: dctx = (Mat_DiagBrdn *)dbase->ctx;
200: VecSet(dctx->invD, lbfgs->delta);
201: break;
202: case MAT_LMVM_SYMBROYDEN_SCALE_SCALAR:
203: lbfgs->sigma = lbfgs->delta;
204: break;
205: case MAT_LMVM_SYMBROYDEN_SCALE_NONE:
206: lbfgs->sigma = 1.0;
207: break;
208: default:
209: break;
210: }
211: }
213: /* Update the scaling */
214: if (lbfgs->scale_type == MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL) MatLMVMUpdate(lbfgs->D, X, F);
216: if (lbfgs->watchdog > lbfgs->max_seq_rejects) {
217: MatLMVMReset(B, PETSC_FALSE);
218: if (lbfgs->scale_type == MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL) MatLMVMReset(lbfgs->D, PETSC_FALSE);
219: }
221: /* Save the solution and function to be used in the next update */
222: VecCopy(X, lmvm->Xprev);
223: VecCopy(F, lmvm->Fprev);
224: lmvm->prev_set = PETSC_TRUE;
225: return 0;
226: }
228: /*------------------------------------------------------------*/
230: static PetscErrorCode MatCopy_LMVMBFGS(Mat B, Mat M, MatStructure str)
231: {
232: Mat_LMVM *bdata = (Mat_LMVM *)B->data;
233: Mat_SymBrdn *bctx = (Mat_SymBrdn *)bdata->ctx;
234: Mat_LMVM *mdata = (Mat_LMVM *)M->data;
235: Mat_SymBrdn *mctx = (Mat_SymBrdn *)mdata->ctx;
236: PetscInt i;
238: mctx->needP = bctx->needP;
239: for (i = 0; i <= bdata->k; ++i) {
240: mctx->stp[i] = bctx->stp[i];
241: mctx->yts[i] = bctx->yts[i];
242: VecCopy(bctx->P[i], mctx->P[i]);
243: }
244: mctx->scale_type = bctx->scale_type;
245: mctx->alpha = bctx->alpha;
246: mctx->beta = bctx->beta;
247: mctx->rho = bctx->rho;
248: mctx->delta = bctx->delta;
249: mctx->sigma_hist = bctx->sigma_hist;
250: mctx->watchdog = bctx->watchdog;
251: mctx->max_seq_rejects = bctx->max_seq_rejects;
252: switch (bctx->scale_type) {
253: case MAT_LMVM_SYMBROYDEN_SCALE_SCALAR:
254: mctx->sigma = bctx->sigma;
255: break;
256: case MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL:
257: MatCopy(bctx->D, mctx->D, SAME_NONZERO_PATTERN);
258: break;
259: case MAT_LMVM_SYMBROYDEN_SCALE_NONE:
260: mctx->sigma = 1.0;
261: break;
262: default:
263: break;
264: }
265: return 0;
266: }
268: /*------------------------------------------------------------*/
270: static PetscErrorCode MatReset_LMVMBFGS(Mat B, PetscBool destructive)
271: {
272: Mat_LMVM *lmvm = (Mat_LMVM *)B->data;
273: Mat_SymBrdn *lbfgs = (Mat_SymBrdn *)lmvm->ctx;
274: Mat_LMVM *dbase;
275: Mat_DiagBrdn *dctx;
277: lbfgs->watchdog = 0;
278: lbfgs->needP = PETSC_TRUE;
279: if (lbfgs->allocated) {
280: if (destructive) {
281: VecDestroy(&lbfgs->work);
282: PetscFree4(lbfgs->stp, lbfgs->yts, lbfgs->yty, lbfgs->sts);
283: VecDestroyVecs(lmvm->m, &lbfgs->P);
284: switch (lbfgs->scale_type) {
285: case MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL:
286: MatLMVMReset(lbfgs->D, PETSC_TRUE);
287: break;
288: default:
289: break;
290: }
291: lbfgs->allocated = PETSC_FALSE;
292: } else {
293: switch (lbfgs->scale_type) {
294: case MAT_LMVM_SYMBROYDEN_SCALE_SCALAR:
295: lbfgs->sigma = lbfgs->delta;
296: break;
297: case MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL:
298: MatLMVMReset(lbfgs->D, PETSC_FALSE);
299: dbase = (Mat_LMVM *)lbfgs->D->data;
300: dctx = (Mat_DiagBrdn *)dbase->ctx;
301: VecSet(dctx->invD, lbfgs->delta);
302: break;
303: case MAT_LMVM_SYMBROYDEN_SCALE_NONE:
304: lbfgs->sigma = 1.0;
305: break;
306: default:
307: break;
308: }
309: }
310: }
311: MatReset_LMVM(B, destructive);
312: return 0;
313: }
315: /*------------------------------------------------------------*/
317: static PetscErrorCode MatAllocate_LMVMBFGS(Mat B, Vec X, Vec F)
318: {
319: Mat_LMVM *lmvm = (Mat_LMVM *)B->data;
320: Mat_SymBrdn *lbfgs = (Mat_SymBrdn *)lmvm->ctx;
322: MatAllocate_LMVM(B, X, F);
323: if (!lbfgs->allocated) {
324: VecDuplicate(X, &lbfgs->work);
325: PetscMalloc4(lmvm->m, &lbfgs->stp, lmvm->m, &lbfgs->yts, lmvm->m, &lbfgs->yty, lmvm->m, &lbfgs->sts);
326: if (lmvm->m > 0) VecDuplicateVecs(X, lmvm->m, &lbfgs->P);
327: switch (lbfgs->scale_type) {
328: case MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL:
329: MatLMVMAllocate(lbfgs->D, X, F);
330: break;
331: default:
332: break;
333: }
334: lbfgs->allocated = PETSC_TRUE;
335: }
336: return 0;
337: }
339: /*------------------------------------------------------------*/
341: static PetscErrorCode MatDestroy_LMVMBFGS(Mat B)
342: {
343: Mat_LMVM *lmvm = (Mat_LMVM *)B->data;
344: Mat_SymBrdn *lbfgs = (Mat_SymBrdn *)lmvm->ctx;
346: if (lbfgs->allocated) {
347: VecDestroy(&lbfgs->work);
348: PetscFree4(lbfgs->stp, lbfgs->yts, lbfgs->yty, lbfgs->sts);
349: VecDestroyVecs(lmvm->m, &lbfgs->P);
350: lbfgs->allocated = PETSC_FALSE;
351: }
352: MatDestroy(&lbfgs->D);
353: PetscFree(lmvm->ctx);
354: MatDestroy_LMVM(B);
355: return 0;
356: }
358: /*------------------------------------------------------------*/
360: static PetscErrorCode MatSetUp_LMVMBFGS(Mat B)
361: {
362: Mat_LMVM *lmvm = (Mat_LMVM *)B->data;
363: Mat_SymBrdn *lbfgs = (Mat_SymBrdn *)lmvm->ctx;
364: PetscInt n, N;
366: MatSetUp_LMVM(B);
367: lbfgs->max_seq_rejects = lmvm->m / 2;
368: if (!lbfgs->allocated) {
369: VecDuplicate(lmvm->Xprev, &lbfgs->work);
370: PetscMalloc4(lmvm->m, &lbfgs->stp, lmvm->m, &lbfgs->yts, lmvm->m, &lbfgs->yty, lmvm->m, &lbfgs->sts);
371: if (lmvm->m > 0) VecDuplicateVecs(lmvm->Xprev, lmvm->m, &lbfgs->P);
372: switch (lbfgs->scale_type) {
373: case MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL:
374: MatGetLocalSize(B, &n, &n);
375: MatGetSize(B, &N, &N);
376: MatSetSizes(lbfgs->D, n, n, N, N);
377: MatSetUp(lbfgs->D);
378: break;
379: default:
380: break;
381: }
382: lbfgs->allocated = PETSC_TRUE;
383: }
384: return 0;
385: }
387: /*------------------------------------------------------------*/
389: static PetscErrorCode MatSetFromOptions_LMVMBFGS(Mat B, PetscOptionItems *PetscOptionsObject)
390: {
391: MatSetFromOptions_LMVM(B, PetscOptionsObject);
392: PetscOptionsHeadBegin(PetscOptionsObject, "L-BFGS method for approximating SPD Jacobian actions (MATLMVMBFGS)");
393: MatSetFromOptions_LMVMSymBrdn_Private(B, PetscOptionsObject);
394: PetscOptionsHeadEnd();
395: return 0;
396: }
398: /*------------------------------------------------------------*/
400: PetscErrorCode MatCreate_LMVMBFGS(Mat B)
401: {
402: Mat_LMVM *lmvm;
403: Mat_SymBrdn *lbfgs;
405: MatCreate_LMVMSymBrdn(B);
406: PetscObjectChangeTypeName((PetscObject)B, MATLMVMBFGS);
407: B->ops->setup = MatSetUp_LMVMBFGS;
408: B->ops->destroy = MatDestroy_LMVMBFGS;
409: B->ops->setfromoptions = MatSetFromOptions_LMVMBFGS;
410: B->ops->solve = MatSolve_LMVMBFGS;
412: lmvm = (Mat_LMVM *)B->data;
413: lmvm->ops->allocate = MatAllocate_LMVMBFGS;
414: lmvm->ops->reset = MatReset_LMVMBFGS;
415: lmvm->ops->update = MatUpdate_LMVMBFGS;
416: lmvm->ops->mult = MatMult_LMVMBFGS;
417: lmvm->ops->copy = MatCopy_LMVMBFGS;
419: lbfgs = (Mat_SymBrdn *)lmvm->ctx;
420: lbfgs->needQ = PETSC_FALSE;
421: lbfgs->phi = 0.0;
422: return 0;
423: }
425: /*------------------------------------------------------------*/
427: /*@
428: MatCreateLMVMBFGS - Creates a limited-memory Broyden-Fletcher-Goldfarb-Shano (BFGS)
429: matrix used for approximating Jacobians. L-BFGS is symmetric positive-definite by
430: construction, and is commonly used to approximate Hessians in optimization
431: problems.
433: The provided local and global sizes must match the solution and function vectors
434: used with `MatLMVMUpdate()` and `MatSolve()`. The resulting L-BFGS matrix will have
435: storage vectors allocated with `VecCreateSeq()` in serial and `VecCreateMPI()` in
436: parallel. To use the L-BFGS matrix with other vector types, the matrix must be
437: created using `MatCreate()` and `MatSetType()`, followed by `MatLMVMAllocate()`.
438: This ensures that the internal storage and work vectors are duplicated from the
439: correct type of vector.
441: Collective
443: Input Parameters:
444: + comm - MPI communicator
445: . n - number of local rows for storage vectors
446: - N - global size of the storage vectors
448: Output Parameter:
449: . B - the matrix
451: Options Database Keys:
452: + -mat_lmvm_scale_type - (developer) type of scaling applied to J0 (none, scalar, diagonal)
453: . -mat_lmvm_theta - (developer) convex ratio between BFGS and DFP components of the diagonal J0 scaling
454: . -mat_lmvm_rho - (developer) update limiter for the J0 scaling
455: . -mat_lmvm_alpha - (developer) coefficient factor for the quadratic subproblem in J0 scaling
456: . -mat_lmvm_beta - (developer) exponential factor for the diagonal J0 scaling
457: - -mat_lmvm_sigma_hist - (developer) number of past updates to use in J0 scaling
459: Level: intermediate
461: Note:
462: It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`
463: paradigm instead of this routine directly.
465: .seealso: [](chapter_ksp), `MatCreate()`, `MATLMVM`, `MATLMVMBFGS`, `MatCreateLMVMDFP()`, `MatCreateLMVMSR1()`,
466: `MatCreateLMVMBrdn()`, `MatCreateLMVMBadBrdn()`, `MatCreateLMVMSymBrdn()`
467: @*/
468: PetscErrorCode MatCreateLMVMBFGS(MPI_Comm comm, PetscInt n, PetscInt N, Mat *B)
469: {
470: KSPInitializePackage();
471: MatCreate(comm, B);
472: MatSetSizes(*B, n, n, N, N);
473: MatSetType(*B, MATLMVMBFGS);
474: MatSetUp(*B);
475: return 0;
476: }