:orphan: # MatCreateSeqAIJCUSPARSE Creates a sparse matrix in `MATAIJCUSPARSE` (compressed row) format (the default parallel PETSc format). This matrix will ultimately pushed down to NVIDIA GPUs and use the CuSPARSE library for calculations. For good matrix assembly performance the user should preallocate the matrix storage by setting the parameter `nz` (or the array `nnz`). ## Synopsis ``` #include "petscmat.h" PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A) ``` Collective ## Input Parameters - ***comm -*** MPI communicator, set to `PETSC_COMM_SELF` - ***m -*** number of rows - ***n -*** number of columns - ***nz -*** number of nonzeros per row (same for all rows), ignored if `nnz` is provide - ***nnz -*** array containing the number of nonzeros in the various rows (possibly different for each row) or `NULL` ## Output Parameter - ***A -*** the matrix ## Notes It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`, MatXXXXSetPreallocation() paradgm instead of this routine directly. [MatXXXXSetPreallocation() is, for example, `MatSeqAIJSetPreallocation()`] The AIJ format, also called compressed row storage, is fully compatible with standard Fortran storage. That is, the stored row and column indices can begin at either one (as in Fortran) or zero. Specify the preallocated storage with either nz or nnz (not both). Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory allocation. ## See Also [](ch_matrices), `Mat`, `MATSEQAIJCUSPARSE`, `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`, `MatCreateAIJ()`, `MATSEQAIJCUSPARSE`, `MATAIJCUSPARSE` ## Level intermediate ## Location src/mat/impls/aij/seq/seqcusparse/aijcusparse.cu --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/impls/aij/seq/seqcusparse/aijcusparse.cu) [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)