MatCreateSELLCUDA#

Creates a sparse matrix in SELL format. This matrix will ultimately pushed down to NVIDIA GPUs.

Synopsis#

#include "petscmat.h" 
PetscErrorCode MatCreateSELLCUDA(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[], Mat *A)

Collective

Input Parameters#

  • comm - MPI communicator, set to PETSC_COMM_SELF

  • m - number of local rows (or PETSC_DECIDE to have calculated if M is given) This value should be the same as the local size used in creating the y vector for the matrix-vector product y = Ax.

  • n - This value should be the same as the local size used in creating the x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have calculated if N is given) For square matrices n is almost always m.

  • M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)

  • N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)

  • d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix (same value is used for all local rows)

  • d_nnz - array containing the number of nonzeros in the various rows of the DIAGONAL portion of the local submatrix (possibly different for each row) or NULL, if d_nz is used to specify the nonzero structure. The size of this array is equal to the number of local rows, i.e m. For matrices you plan to factor you must leave room for the diagonal entry and put in the entry even if it is zero.

  • o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local submatrix (same value is used for all local rows).

  • o_nnz - array containing the number of nonzeros in the various rows of the OFF-DIAGONAL portion of the local submatrix (possibly different for each row) or NULL, if o_nz is used to specify the nonzero structure. The size of this array is equal to the number of local rows, i.e m.

Output Parameter#

  • A - the matrix

Notes#

If nnz is given then nz is ignored

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#

Matrices, Mat, MatCreate(), MatCreateSELL(), MatSetValues(), MATMPISELLCUDA, MATSELLCUDA

Level#

intermediate

Location#

src/mat/impls/sell/mpi/mpicuda/mpisellcuda.cu


Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages