Actual source code: symbrdn.h
1: #pragma once
3: #include <../src/ksp/ksp/utils/lmvm/lmvm.h>
4: #include <../src/ksp/ksp/utils/lmvm/rescale/symbrdnrescale.h>
6: /*
7: Limited-memory Symmetric Broyden method for approximating both
8: the forward product and inverse application of a Jacobian.
9: */
11: // bases needed by symmetric [bad] Broyden algorithms beyond those in Mat_LMVM
12: enum {
13: SYMBROYDEN_BASIS_BKS = 0, // B_k S_k for recursive algorithms
14: SYMBROYDEN_BASIS_HKY = 1, // dual to the above, H_k Y_k
15: SYMBROYDEN_BASIS_COUNT
16: };
18: typedef PetscInt SymBroydenBasisType;
20: // products needed by symmetric [bad] Broyden algorithms beyond those in Mat_LMVM
21: enum {
22: SYMBROYDEN_PRODUCTS_PHI = 0, // diagonal: either phi_k = phi_scalar (symm. Broyden), or phi_k is different for every k (symm. bad Broyden)
23: SYMBROYDEN_PRODUCTS_PSI = 1, // diagonal: either psi_k = psi_scalar (symm. bad Broyden), or psi_k is different for every k (symm. Broyden)
24: SYMBROYDEN_PRODUCTS_STBKS = 2, // diagonal S_k^T B_k S_k values for recursive algorithms
25: SYMBROYDEN_PRODUCTS_YTHKY = 3, // dual to the above: diagonal Y_k^T H_k Y_k values
26: SYMBROYDEN_PRODUCTS_M00 = 4, // matrix that appears in (B_* S) M_00 (B_* S)^T rank-m updates, either diagonal (recursive) or full (compact)
27: SYMBROYDEN_PRODUCTS_N00 = 5, // dual to the above, appears in (H_* Y) N_00 (H_* Y)^T rank-m updates
28: SYMBROYDEN_PRODUCTS_M01 = 6, // matrix that appears in (B_* S) M_01 Y^T rank-m updates, either diagonal (recursive) or full (compact)
29: SYMBROYDEN_PRODUCTS_N01 = 7, // dual to the above, appears in (H_* Y) N_01 S^T rank-m updates
30: SYMBROYDEN_PRODUCTS_M11 = 8, // matrix that appers in Y M_11 Y^T rank-m updates, either diagonal (recursive) or full (compact)
31: SYMBROYDEN_PRODUCTS_N11 = 9, // dual to the above, appers in S N_11 S^T rank-m updates
32: SYMBROYDEN_PRODUCTS_COUNT
33: };
35: typedef PetscInt SymBroydenProductsType;
37: typedef struct {
38: PetscReal phi_scalar, psi_scalar;
39: PetscInt watchdog, max_seq_rejects; /* tracker to reset after a certain # of consecutive rejects */
40: SymBroydenRescale rescale; /* context for diagonal or scalar rescaling */
41: LMBasis basis[SYMBROYDEN_BASIS_COUNT];
42: LMProducts products[SYMBROYDEN_PRODUCTS_COUNT];
43: Vec StFprev, YtH0Fprev;
44: } Mat_SymBrdn;
46: PETSC_INTERN PetscErrorCode SymBroydenKernel_Recursive(Mat, MatLMVMMode, Vec, Vec, PetscBool);
47: PETSC_INTERN PetscErrorCode SymBroydenKernel_CompactDense(Mat, MatLMVMMode, Vec, Vec, PetscBool);
49: PETSC_INTERN PetscErrorCode DFPKernel_Recursive(Mat, MatLMVMMode, Vec, Vec);
50: PETSC_INTERN PetscErrorCode DFPKernel_CompactDense(Mat, MatLMVMMode, Vec, Vec);
51: PETSC_INTERN PetscErrorCode DFPKernel_Dense(Mat, MatLMVMMode, Vec, Vec);
53: PETSC_INTERN PetscErrorCode BFGSKernel_Recursive(Mat, MatLMVMMode, Vec, Vec);
54: PETSC_INTERN PetscErrorCode BFGSKernel_CompactDense(Mat, MatLMVMMode, Vec, Vec);
56: PETSC_INTERN PetscErrorCode SymBroydenCompactDenseKernelUseB0S(Mat, MatLMVMMode, Vec, PetscBool *);