# MatCreateLMVMBroyden
Creates a limited-memory "good" Broyden-type approximation matrix used for a Jacobian. L-Brdn is not guaranteed to be symmetric or positive-definite. 
## Synopsis
```
#include "petscksp.h" 
PetscErrorCode MatCreateLMVMBroyden(MPI_Comm comm, PetscInt n, PetscInt N, Mat *B)
```
The provided local and global sizes must match the solution and function vectors
used with `MatLMVMUpdate()` and `MatSolve()`. The resulting L-Brdn matrix will have
storage vectors allocated with `VecCreateSeq()` in serial and `VecCreateMPI()` in
parallel. To use the L-Brdn 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





## Note
It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`
paradigm instead of this routine directly.


## See Also
 [](chapter_ksp), `MatCreate()`, `MATLMVM`, `MATLMVMBRDN`, `MatCreateLMVMDFP()`, `MatCreateLMVMSR1()`,
`MatCreateLMVMBFGS()`, `MatCreateLMVMBadBrdn()`, `MatCreateLMVMSymBrdn()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/utils/lmvm/brdn/brdn.c.html#MatCreateLMVMBroyden">src/ksp/ksp/utils/lmvm/brdn/brdn.c</A>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/utils/lmvm/brdn/brdn.c)


[Index of all KSP routines](index.md)  
[Table of Contents for all manual pages](/docs/manualpages/index.md)  
[Index of all manual pages](/docs/manualpages/singleindex.md)  
