:orphan: # MatCreateLMVMSymBadBroyden Creates a limited-memory Symmetric "Bad" Broyden-type matrix used for approximating Jacobians. L-SymBadBrdn is a convex combination of L-DFP and L-BFGS such that SymBadBrdn^{-1} = (1 - phi)*BFGS^{-1} + phi*DFP^{-1}. The combination factor phi is restricted to the range [0, 1], where the L-SymBadBrdn matrix is guaranteed to be symmetric positive-definite. Note that this combination is on the inverses and not on the forwards. For forward convex combinations, use the L-SymBrdn matrix. ## Synopsis ``` #include "petscksp.h" PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm comm, PetscInt n, PetscInt N, Mat *B) ``` To use the L-SymBrdn matrix with other vector types, the matrix must be created using MatCreate() and MatSetType(), followed by `MatLMVMAllocate()`. This ensures that the internal storage and work vectors are duplicated from the correct type of vector. Collective ## Input Parameters - ***comm -*** MPI communicator - ***n -*** number of local rows for storage vectors - ***N -*** global size of the storage vectors ## Output Parameter - ***B -*** the matrix ## Options Database Keys - ***-mat_lmvm_phi -*** (developer) convex ratio between BFGS and DFP components of the update - ***-mat_lmvm_scale_type -*** (developer) type of scaling applied to J0 (none, scalar, diagonal) - ***-mat_lmvm_theta -*** (developer) convex ratio between BFGS and DFP components of the diagonal J0 scaling - ***-mat_lmvm_rho -*** (developer) update limiter for the J0 scaling - ***-mat_lmvm_alpha -*** (developer) coefficient factor for the quadratic subproblem in J0 scaling - ***-mat_lmvm_beta -*** (developer) exponential factor for the diagonal J0 scaling - ***-mat_lmvm_sigma_hist -*** (developer) number of past updates to use in J0 scaling ## Note It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()` paradigm instead of this routine directly. ## See Also [](ch_ksp), [LMVM Matrices](sec_matlmvm), `MatCreate()`, `MATLMVM`, `MATLMVMSYMBROYDEN`, `MatCreateLMVMDFP()`, `MatCreateLMVMSR1()`, `MatCreateLMVMBFGS()`, `MatCreateLMVMBrdn()`, `MatCreateLMVMBadBrdn()` ## Level intermediate ## Location src/ksp/ksp/utils/lmvm/symbrdn/symbadbrdn.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/utils/lmvm/symbrdn/symbadbrdn.c) [Index of all KSP routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)