Actual source code: symbrdn.h
petsc-3.11.4 2019-09-28
1: #include <../src/ksp/ksp/utils/lmvm/lmvm.h>
3: /*
4: Limited-memory Symmetric Broyden method for approximating both
5: the forward product and inverse application of a Jacobian.
6: */
8: typedef struct {
9: Mat D; /* diagonal scaling term */
10: Vec *P, *Q; /* storage vectors for (B_i)*S[i] and (B_i)^{-1}*Y[i] */
11: Vec invDnew, invD, BFGS, DFP, U, V, W; /* work vectors for diagonal scaling */
12: Vec work;
13: PetscBool allocated, needP, needQ;
14: PetscReal *stp, *ytq, *yts, *yty, *sts; /* scalar arrays for recycling dot products */
15: PetscReal theta, phi, *psi; /* convex combination factors between DFP and BFGS */
16: PetscReal rho, alpha, beta; /* convex combination factors for the scalar or diagonal scaling */
17: PetscReal delta, delta_min, delta_max, sigma;
18: PetscInt sigma_hist; /* length of update history to be used for scaling */
19: PetscInt scale_type;
20: PetscInt watchdog, max_seq_rejects; /* tracker to reset after a certain # of consecutive rejects */
21: } Mat_SymBrdn;
23: #define SYMBRDN_SCALE_NONE 0
24: #define SYMBRDN_SCALE_SCALAR 1
25: #define SYMBRDN_SCALE_DIAG 2
26: #define SYMBRDN_SCALE_SIZE 3
28: static const char *Scale_Table[64] = {"none","scalar","diagonal"};
30: PETSC_INTERN PetscErrorCode MatSymBrdnApplyJ0Fwd(Mat, Vec, Vec);
31: PETSC_INTERN PetscErrorCode MatSymBrdnApplyJ0Inv(Mat, Vec, Vec);
32: PETSC_INTERN PetscErrorCode MatSymBrdnComputeJ0Diag(Mat);
33: PETSC_INTERN PetscErrorCode MatSymBrdnComputeJ0Scalar(Mat);
35: PETSC_INTERN PetscErrorCode MatView_LMVMSymBrdn(Mat, PetscViewer);
36: PETSC_INTERN PetscErrorCode MatSetFromOptions_LMVMSymBrdn(PetscOptionItems*, Mat);