:orphan:
# MatCreateMPIAdj
Creates a sparse matrix representing an adjacency list. The matrix does not have numerical values associated with it, but is intended for ordering (to reduce bandwidth etc) and partitioning. 
## Synopsis
```
#include "petscmat.h" 
PetscErrorCode MatCreateMPIAdj(MPI_Comm comm, PetscInt m, PetscInt N, PetscInt *i, PetscInt *j, PetscInt *values, Mat *A)
```
Collective


## Input Parameters

- ***comm -*** MPI communicator
- ***m -*** number of local rows
- ***N -*** number of global columns
- ***i -*** the indices into `j` for the start of each row
- ***j -*** the column indices for each row (sorted for each row).
The indices in `i` and `j` start with zero (NOT with one).
- ***values -*** [use `NULL` if not provided] edge weights



## Output Parameter

- ***A -*** the matrix





## Notes
You must NOT free the `i`, `values` and `j` arrays yourself. PETSc will free them
when the matrix is destroyed; you must allocate them with `PetscMalloc()`. If you
call from Fortran you need not create the arrays with `PetscMalloc()`.

You should not include the matrix diagonals.

If you already have a matrix, you can create its adjacency matrix by a call
to `MatConvert()`, specifying a type of `MATMPIADJ`.

Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC`


## See Also
 [](ch_matrices), `Mat`, `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/adj/mpi/mpiadj.c.html#MatCreateMPIAdj">src/mat/impls/adj/mpi/mpiadj.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex11.c.html">src/mat/tutorials/ex11.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/impls/adj/mpi/mpiadj.c)


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