Actual source code: denseqn.h
1: #include <../src/ksp/ksp/utils/lmvm/lmvm.h>
3: /*
4: dense representation for the limited-memory BFGS/DFP method.
5: */
7: typedef struct {
8: Mat diag_qn; /* diagonalized Hessian init */
10: PetscInt num_updates;
11: PetscInt num_mult_updates;
12: Mat Sfull, Yfull, HY, BS; // Stored in recycled order
13: Vec StFprev;
14: Mat StY_triu; // triu(StY) is the R matrix
15: Mat StY_triu_strict; // strict_triu(YtS) is the R matrix
16: Mat YtS_triu_strict; // strict_triu(YtS) is the L^T matrix
17: Mat YtS_triu; // triu(YtS) is the L matrix
18: Mat YtHY;
19: Mat StBS;
20: Mat J;
21: Mat temp_mat;
22: Vec *PQ; /* P for BFGS, Q for DFP */
23: Vec diag_vec;
24: Vec diag_vec_recycle_order;
25: Vec inv_diag_vec;
26: Vec column_work, column_work2, rwork1, rwork2, rwork3;
27: Vec rwork2_local, rwork3_local;
28: Vec local_work_vec, local_work_vec_copy;
29: Vec cyclic_work_vec;
30: MatType dense_type;
31: MatLMVMDenseType strategy;
32: MatLMVMSymBroydenScaleType scale_type;
34: PetscReal *ytq, *stp, *yts;
35: PetscScalar *workscalar;
36: PetscInt S_count, St_count, Y_count, Yt_count;
37: PetscInt watchdog, max_seq_rejects; /* tracker to reset after a certain # of consecutive rejects */
38: PetscBool allocated, use_recursive, needPQ; /* P for BFGS, Q for DFP */
39: Vec Fprev_ref;
40: PetscObjectState Fprev_state;
41: } Mat_DQN;
43: PETSC_INTERN PetscErrorCode MatView_LMVMDDFP(Mat, PetscViewer);
44: PETSC_INTERN PetscErrorCode MatView_LMVMDBFGS(Mat, PetscViewer);
46: PETSC_INTERN PetscErrorCode MatUpperTriangularSolveInPlace_CUPM(PetscBool, PetscInt, const PetscScalar[], PetscInt, PetscScalar[], PetscInt);
47: PETSC_INTERN PetscErrorCode MatUpperTriangularSolveInPlaceCyclic_CUPM(PetscBool, PetscInt, PetscInt, const PetscScalar[], PetscInt, PetscScalar[], PetscInt);
48: PETSC_INTERN PetscErrorCode MatMultAddColumnRange(Mat, Vec, Vec, Vec, PetscInt, PetscInt);
49: PETSC_INTERN PetscErrorCode MatMultHermitianTransposeColumnRange(Mat, Vec, Vec, PetscInt, PetscInt);
50: PETSC_INTERN PetscErrorCode MatMultHermitianTransposeAddColumnRange(Mat, Vec, Vec, Vec, PetscInt, PetscInt);
51: PETSC_INTERN PetscErrorCode VecCyclicShift(Mat, Vec, PetscInt, Vec);
52: PETSC_INTERN PetscErrorCode VecRecycleOrderToHistoryOrder(Mat, Vec, PetscInt, Vec);
53: PETSC_INTERN PetscErrorCode VecHistoryOrderToRecycleOrder(Mat, Vec, PetscInt, Vec);
54: PETSC_INTERN PetscErrorCode MatUpperTriangularSolveInPlace(Mat, Mat, Vec, PetscBool, PetscInt, MatLMVMDenseType);
55: PETSC_INTERN PetscErrorCode MatMove_LR3(Mat, Mat, PetscInt);