Actual source code: mpiaijcusp.cu

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <petscconf.h>
  2: PETSC_CUDA_EXTERN_C_BEGIN
  3: #include <../src/mat/impls/aij/mpi/mpiaij.h>   /*I "petscmat.h" I*/
  4: PETSC_CUDA_EXTERN_C_END
  5: #include <../src/mat/impls/aij/mpi/mpicusp/mpicuspmatimpl.h>

  9: PetscErrorCode  MatMPIAIJSetPreallocation_MPIAIJCUSP(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
 10: {
 11:   Mat_MPIAIJ *b = (Mat_MPIAIJ*)B->data;
 12:   Mat_MPIAIJCUSP * cuspStruct = (Mat_MPIAIJCUSP*)b->spptr;
 14:   PetscInt       i;

 17:   PetscLayoutSetUp(B->rmap);
 18:   PetscLayoutSetUp(B->cmap);
 19:   if (d_nnz) {
 20:     for (i=0; i<B->rmap->n; i++) {
 21:       if (d_nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"d_nnz cannot be less than 0: local row %D value %D",i,d_nnz[i]);
 22:     }
 23:   }
 24:   if (o_nnz) {
 25:     for (i=0; i<B->rmap->n; i++) {
 26:       if (o_nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"o_nnz cannot be less than 0: local row %D value %D",i,o_nnz[i]);
 27:     }
 28:   }
 29:   if (!B->preallocated) {
 30:     /* Explicitly create 2 MATSEQAIJCUSP matrices. */
 31:     MatCreate(PETSC_COMM_SELF,&b->A);
 32:     MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
 33:     MatSetType(b->A,MATSEQAIJCUSP);
 34:     PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
 35:     MatCreate(PETSC_COMM_SELF,&b->B);
 36:     MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
 37:     MatSetType(b->B,MATSEQAIJCUSP);
 38:     PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
 39:   }
 40:   MatSeqAIJSetPreallocation(b->A,d_nz,d_nnz);
 41:   MatSeqAIJSetPreallocation(b->B,o_nz,o_nnz);
 42:   MatCUSPSetFormat(b->A,MAT_CUSP_MULT,cuspStruct->diagGPUMatFormat);
 43:   MatCUSPSetFormat(b->B,MAT_CUSP_MULT,cuspStruct->offdiagGPUMatFormat);
 44:   MatCUSPSetStream(b->A,cuspStruct->stream);
 45:   MatCUSPSetStream(b->B,cuspStruct->stream);
 46:   B->preallocated = PETSC_TRUE;
 47:   return(0);
 48: }

 52: PetscErrorCode  MatGetVecs_MPIAIJCUSP(Mat mat,Vec *right,Vec *left)
 53: {
 55:   PetscInt rbs,cbs;

 58:   MatGetBlockSizes(mat,&rbs,&cbs);
 59:   if (right) {
 60:     VecCreate(PetscObjectComm((PetscObject)mat),right);
 61:     VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);
 62:     VecSetBlockSize(*right,cbs);
 63:     VecSetType(*right,VECCUSP);
 64:     VecSetLayout(*right,mat->cmap);
 65:   }
 66:   if (left) {
 67:     VecCreate(PetscObjectComm((PetscObject)mat),left);
 68:     VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);
 69:     VecSetBlockSize(*left,rbs);
 70:     VecSetType(*left,VECCUSP);
 71:     VecSetLayout(*left,mat->rmap);
 72:   }
 73:   return(0);
 74: }


 79: PetscErrorCode MatMult_MPIAIJCUSP(Mat A,Vec xx,Vec yy)
 80: {
 81:   /* This multiplication sequence is different sequence
 82:      than the CPU version. In particular, the diagonal block
 83:      multiplication kernel is launched in one stream. Then,
 84:      in a separate stream, the data transfers from DeviceToHost
 85:      (with MPI messaging in between), then HostToDevice are
 86:      launched. Once the data transfer stream is synchronized,
 87:      to ensure messaging is complete, the MatMultAdd kernel
 88:      is launched in the original (MatMult) stream to protect
 89:      against race conditions.

 91:      This sequence should only be called for GPU computation. */
 92:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
 94:   PetscInt       nt;

 97:   VecGetLocalSize(xx,&nt);
 98:   if (nt != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A (%D) and xx (%D)",A->cmap->n,nt);
 99:   VecScatterInitializeForGPU(a->Mvctx,xx,SCATTER_FORWARD);
100:   (*a->A->ops->mult)(a->A,xx,yy);
101:   VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
102:   VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
103:   (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
104:   VecScatterFinalizeForGPU(a->Mvctx);
105:   return(0);
106: }

108: PetscErrorCode MatSetValuesBatch_MPIAIJCUSP(Mat J, PetscInt Ne, PetscInt Nl, PetscInt *elemRows, const PetscScalar *elemMats);

112: PetscErrorCode MatCUSPSetFormat_MPIAIJCUSP(Mat A,MatCUSPFormatOperation op,MatCUSPStorageFormat format)
113: {
114:   Mat_MPIAIJ     *a           = (Mat_MPIAIJ*)A->data;
115:   Mat_MPIAIJCUSP * cuspStruct = (Mat_MPIAIJCUSP*)a->spptr;

118:   switch (op) {
119:   case MAT_CUSP_MULT_DIAG:
120:     cuspStruct->diagGPUMatFormat = format;
121:     break;
122:   case MAT_CUSP_MULT_OFFDIAG:
123:     cuspStruct->offdiagGPUMatFormat = format;
124:     break;
125:   case MAT_CUSP_ALL:
126:     cuspStruct->diagGPUMatFormat    = format;
127:     cuspStruct->offdiagGPUMatFormat = format;
128:     break;
129:   default:
130:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPFormatOperation. Only MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_DIAG, and MAT_CUSP_MULT_ALL are currently supported.",op);
131:   }
132:   return(0);
133: }

137: PetscErrorCode MatSetFromOptions_MPIAIJCUSP(Mat A)
138: {
139:   MatCUSPStorageFormat format;
140:   PetscErrorCode       ierr;
141:   PetscBool            flg;

144:   PetscOptionsHead("MPIAIJCUSP options");
145:   PetscObjectOptionsBegin((PetscObject)A);
146:   if (A->factortype==MAT_FACTOR_NONE) {
147:     PetscOptionsEnum("-mat_cusp_mult_diag_storage_format","sets storage format of the diagonal blocks of (mpi)aijcusp gpu matrices for SpMV",
148:                             "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)MAT_CUSP_CSR,(PetscEnum*)&format,&flg);
149:     if (flg) {
150:       MatCUSPSetFormat(A,MAT_CUSP_MULT_DIAG,format);
151:     }
152:     PetscOptionsEnum("-mat_cusp_mult_offdiag_storage_format","sets storage format of the off-diagonal blocks (mpi)aijcusp gpu matrices for SpMV",
153:                             "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)MAT_CUSP_CSR,(PetscEnum*)&format,&flg);
154:     if (flg) {
155:       MatCUSPSetFormat(A,MAT_CUSP_MULT_OFFDIAG,format);
156:     }
157:     PetscOptionsEnum("-mat_cusp_storage_format","sets storage format of the diagonal and off-diagonal blocks (mpi)aijcusp gpu matrices for SpMV",
158:                             "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)MAT_CUSP_CSR,(PetscEnum*)&format,&flg);
159:     if (flg) {
160:       MatCUSPSetFormat(A,MAT_CUSP_ALL,format);
161:     }
162:   }
163:   PetscOptionsEnd();
164:   return(0);
165: }

169: PetscErrorCode MatDestroy_MPIAIJCUSP(Mat A)
170: {
172:   Mat_MPIAIJ     *a           = (Mat_MPIAIJ*)A->data;
173:   Mat_MPIAIJCUSP *cuspStruct = (Mat_MPIAIJCUSP*)a->spptr;
174:   cudaError_t    err=cudaSuccess;

177:   try {
178:     err = cudaStreamDestroy(cuspStruct->stream);
179:     if (err!=cudaSuccess)
180:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Mat_MPIAIJCUSP error: %s", cudaGetErrorString(err));
181:     delete cuspStruct;
182:   } catch(char *ex) {
183:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Mat_MPIAIJCUSP error: %s", ex);
184:   }
185:   cuspStruct = 0;
186:   MatDestroy_MPIAIJ(A);
187:   return(0);
188: }

192: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJCUSP(Mat A)
193: {
195:   Mat_MPIAIJ     *a;
196:   Mat_MPIAIJCUSP * cuspStruct;
197:   cudaError_t    err=cudaSuccess;

200:   MatCreate_MPIAIJ(A);
201:   PetscObjectComposeFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",MatMPIAIJSetPreallocation_MPIAIJCUSP);
202:   A->ops->getvecs        = MatGetVecs_MPIAIJCUSP;
203:   A->ops->setvaluesbatch = MatSetValuesBatch_MPIAIJCUSP;

205:   a          = (Mat_MPIAIJ*)A->data;
206:   a->spptr   = new Mat_MPIAIJCUSP;
207:   cuspStruct = (Mat_MPIAIJCUSP*)a->spptr;

209:   cuspStruct->diagGPUMatFormat    = MAT_CUSP_CSR;
210:   cuspStruct->offdiagGPUMatFormat = MAT_CUSP_CSR;
211:   err = cudaStreamCreate(&(cuspStruct->stream));
212:   if (err!=cudaSuccess)
213:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Mat_MPIAIJCUSP error: %s", cudaGetErrorString(err));

215:   A->ops->mult           = MatMult_MPIAIJCUSP;
216:   A->ops->setfromoptions = MatSetFromOptions_MPIAIJCUSP;
217:   A->ops->destroy        = MatDestroy_MPIAIJCUSP;

219:   PetscObjectComposeFunction((PetscObject)A,"MatCUSPSetFormat_C", MatCUSPSetFormat_MPIAIJCUSP);
220:   PetscObjectChangeTypeName((PetscObject)A,MATMPIAIJCUSP);
221:   return(0);
222: }


225: /*@
226:    MatCreateAIJCUSP - Creates a sparse matrix in AIJ (compressed row) format
227:    (the default parallel PETSc format).  This matrix will ultimately pushed down
228:    to NVidia GPUs and use the CUSP library for calculations. For good matrix
229:    assembly performance the user should preallocate the matrix storage by setting
230:    the parameter nz (or the array nnz).  By setting these parameters accurately,
231:    performance during matrix assembly can be increased by more than a factor of 50.


234:    Collective on MPI_Comm

236:    Input Parameters:
237: +  comm - MPI communicator, set to PETSC_COMM_SELF
238: .  m - number of rows
239: .  n - number of columns
240: .  nz - number of nonzeros per row (same for all rows)
241: -  nnz - array containing the number of nonzeros in the various rows
242:          (possibly different for each row) or NULL

244:    Output Parameter:
245: .  A - the matrix

247:    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
248:    MatXXXXSetPreallocation() paradigm instead of this routine directly.
249:    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]

251:    Notes:
252:    If nnz is given then nz is ignored

254:    The AIJ format (also called the Yale sparse matrix format or
255:    compressed row storage), is fully compatible with standard Fortran 77
256:    storage.  That is, the stored row and column indices can begin at
257:    either one (as in Fortran) or zero.  See the users' manual for details.

259:    Specify the preallocated storage with either nz or nnz (not both).
260:    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
261:    allocation.  For large problems you MUST preallocate memory or you
262:    will get TERRIBLE performance, see the users' manual chapter on matrices.

264:    By default, this format uses inodes (identical nodes) when possible, to
265:    improve numerical efficiency of matrix-vector products and solves. We
266:    search for consecutive rows with the same nonzero structure, thereby
267:    reusing matrix information to achieve increased efficiency.

269:    Level: intermediate

271: .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATMPIAIJCUSP, MATAIJCUSP
272: @*/
275: PetscErrorCode  MatCreateAIJCUSP(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)
276: {
278:   PetscMPIInt    size;

281:   MatCreate(comm,A);
282:   MatSetSizes(*A,m,n,M,N);
283:   MPI_Comm_size(comm,&size);
284:   if (size > 1) {
285:     MatSetType(*A,MATMPIAIJCUSP);
286:     MatMPIAIJSetPreallocation(*A,d_nz,d_nnz,o_nz,o_nnz);
287:   } else {
288:     MatSetType(*A,MATSEQAIJCUSP);
289:     MatSeqAIJSetPreallocation(*A,d_nz,d_nnz);
290:   }
291:   return(0);
292: }

294: /*M
295:    MATAIJCUSP - MATMPIAIJCUSP= "aijcusp" = "mpiaijcusp" - A matrix type to be used for sparse matrices.

297:    A matrix type type whose data resides on Nvidia GPUs. These matrices can be CSR format.
298:    All matrix calculations are performed using the CUSP library. DIA and ELL
299:    formats are also available

301:    This matrix type is identical to MATSEQAIJCUSP when constructed with a single process communicator,
302:    and MATMPIAIJCUSP otherwise.  As a result, for single process communicators,
303:    MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
304:    for communicators controlling multiple processes.  It is recommended that you call both of
305:    the above preallocation routines for simplicity.

307:    Options Database Keys:
308: +  -mat_type mpiaijcusp - sets the matrix type to "mpiaijcusp" during a call to MatSetFromOptions()
309: .  -mat_cusp_storage_format csr - sets the storage format of diagonal and off-diagonal matrices during a call to MatSetFromOptions(). Other storage formats include dia (diagonal) or ell (ellpack).
310: .  -mat_cusp_mult_diag_storage_format csr - sets the storage format of diagonal matrix during a call to MatSetFromOptions(). Other storage formats include dia (diagonal) or ell (ellpack).
311: -  -mat_cusp_mult_offdiag_storage_format csr - sets the storage format of off-diagonal matrix during a call to MatSetFromOptions(). Other storage formats include dia (diagonal) or ell (ellpack).

313:   Level: beginner

315:  .seealso: MatCreateAIJCUSP(), MATSEQAIJCUSP, MatCreateSeqAIJCUSP(), MatCUSPSetFormat(), MatCUSPStorageFormat, MatCUSPFormatOperation
316: M*/