# MatCreateDense
Creates a matrix in `MATDENSE` format.
## Synopsis
```
#include "petscmat.h"
PetscErrorCode MatCreateDense(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar *data, Mat *A)
```
Collective
## Input Parameters
- ***comm -*** MPI communicator
- ***m -*** number of local rows (or `PETSC_DECIDE` to have calculated if M is given)
- ***n -*** number of local columns (or `PETSC_DECIDE` to have calculated if N is given)
- ***M -*** number of global rows (or `PETSC_DECIDE` to have calculated if m is given)
- ***N -*** number of global columns (or `PETSC_DECIDE` to have calculated if n is given)
- ***data -*** optional location of matrix data. Set data to NULL (`PETSC_NULL_SCALAR` for Fortran users) for PETSc
to control all matrix memory allocation.
## Output Parameter
- ***A -*** the matrix
## Notes
The dense format is fully compatible with standard Fortran 77
storage by columns.
Although local portions of the matrix are stored in column-major
order, the matrix is partitioned across MPI ranks by row.
The data input variable is intended primarily for Fortran programmers
who wish to allocate their own matrix memory space. Most users should
set data=NULL (`PETSC_NULL_SCALAR` for Fortran users).
The user MUST specify either the local or global matrix dimensions
(possibly both).
## See Also
`MATDENSE`, `MatCreate()`, `MatCreateSeqDense()`, `MatSetValues()`
## Level
intermediate
## Location
src/mat/impls/dense/mpi/mpidense.c
## Examples
src/ksp/ksp/tutorials/ex21.c.html
src/ksp/ksp/tutorials/ex76.c.html
src/ksp/ksp/tutorials/ex76f.F90.html
src/ksp/ksp/tutorials/ex77.c.html
src/ksp/ksp/tutorials/ex77f.F90.html
src/ksp/ksp/tutorials/ex79.c.html
src/ksp/ksp/tutorials/ex81.c.html
src/ksp/ksp/tutorials/ex82.c.html
src/ts/tutorials/ex16fwd.c.html
src/ts/tutorials/ex20fwd.c.html
src/ts/tutorials/ex23fwdadj.c.html
---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/impls/dense/mpi/mpidense.c)
[Index of all Mat routines](index.md)
[Table of Contents for all manual pages](/docs/manualpages/index.md)
[Index of all manual pages](/docs/manualpages/singleindex.md)