Actual source code: mpiaijcusp.cu
petsc-3.7.3 2016-08-01
1: #define PETSC_SKIP_COMPLEX
2: #define PETSC_SKIP_SPINLOCK
4: #include <petscconf.h>
5: #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I "petscmat.h" I*/
6: #include <../src/mat/impls/aij/mpi/mpicusp/mpicuspmatimpl.h>
10: PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJCUSP(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
11: {
12: Mat_MPIAIJ *b = (Mat_MPIAIJ*)B->data;
13: Mat_MPIAIJCUSP * cuspStruct = (Mat_MPIAIJCUSP*)b->spptr;
15: PetscInt i;
18: PetscLayoutSetUp(B->rmap);
19: PetscLayoutSetUp(B->cmap);
20: if (d_nnz) {
21: for (i=0; i<B->rmap->n; i++) {
22: 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]);
23: }
24: }
25: if (o_nnz) {
26: for (i=0; i<B->rmap->n; i++) {
27: 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]);
28: }
29: }
30: if (!B->preallocated) {
31: /* Explicitly create 2 MATSEQAIJCUSP matrices. */
32: MatCreate(PETSC_COMM_SELF,&b->A);
33: MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
34: MatSetType(b->A,MATSEQAIJCUSP);
35: PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
36: MatCreate(PETSC_COMM_SELF,&b->B);
37: MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
38: MatSetType(b->B,MATSEQAIJCUSP);
39: PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
40: }
41: MatSeqAIJSetPreallocation(b->A,d_nz,d_nnz);
42: MatSeqAIJSetPreallocation(b->B,o_nz,o_nnz);
43: MatCUSPSetFormat(b->A,MAT_CUSP_MULT,cuspStruct->diagGPUMatFormat);
44: MatCUSPSetFormat(b->B,MAT_CUSP_MULT,cuspStruct->offdiagGPUMatFormat);
45: MatCUSPSetStream(b->A,cuspStruct->stream);
46: MatCUSPSetStream(b->B,cuspStruct->stream);
47: B->preallocated = PETSC_TRUE;
48: return(0);
49: }
53: PetscErrorCode MatCreateVecs_MPIAIJCUSP(Mat mat,Vec *right,Vec *left)
54: {
56: PetscInt rbs,cbs;
59: MatGetBlockSizes(mat,&rbs,&cbs);
60: if (right) {
61: VecCreate(PetscObjectComm((PetscObject)mat),right);
62: VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);
63: VecSetBlockSize(*right,cbs);
64: VecSetType(*right,VECCUSP);
65: VecSetLayout(*right,mat->cmap);
66: }
67: if (left) {
68: VecCreate(PetscObjectComm((PetscObject)mat),left);
69: VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);
70: VecSetBlockSize(*left,rbs);
71: VecSetType(*left,VECCUSP);
72: VecSetLayout(*left,mat->rmap);
73: }
74: return(0);
75: }
80: PetscErrorCode MatMult_MPIAIJCUSP(Mat A,Vec xx,Vec yy)
81: {
82: /* This multiplication sequence is different sequence
83: than the CPU version. In particular, the diagonal block
84: multiplication kernel is launched in one stream. Then,
85: in a separate stream, the data transfers from DeviceToHost
86: (with MPI messaging in between), then HostToDevice are
87: launched. Once the data transfer stream is synchronized,
88: to ensure messaging is complete, the MatMultAdd kernel
89: is launched in the original (MatMult) stream to protect
90: against race conditions.
92: This sequence should only be called for GPU computation. */
93: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
95: PetscInt nt;
98: VecGetLocalSize(xx,&nt);
99: 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);
100: VecScatterInitializeForGPU(a->Mvctx,xx,SCATTER_FORWARD);
101: (*a->A->ops->mult)(a->A,xx,yy);
102: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
103: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
104: (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
105: VecScatterFinalizeForGPU(a->Mvctx);
106: return(0);
107: }
109: PetscErrorCode MatSetValuesBatch_MPIAIJCUSP(Mat J, PetscInt Ne, PetscInt Nl, PetscInt *elemRows, const PetscScalar *elemMats);
113: PetscErrorCode MatCUSPSetFormat_MPIAIJCUSP(Mat A,MatCUSPFormatOperation op,MatCUSPStorageFormat format)
114: {
115: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
116: Mat_MPIAIJCUSP * cuspStruct = (Mat_MPIAIJCUSP*)a->spptr;
119: switch (op) {
120: case MAT_CUSP_MULT_DIAG:
121: cuspStruct->diagGPUMatFormat = format;
122: break;
123: case MAT_CUSP_MULT_OFFDIAG:
124: cuspStruct->offdiagGPUMatFormat = format;
125: break;
126: case MAT_CUSP_ALL:
127: cuspStruct->diagGPUMatFormat = format;
128: cuspStruct->offdiagGPUMatFormat = format;
129: break;
130: default:
131: 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);
132: }
133: return(0);
134: }
138: PetscErrorCode MatSetFromOptions_MPIAIJCUSP(PetscOptionItems *PetscOptionsObject,Mat A)
139: {
140: MatCUSPStorageFormat format;
141: PetscErrorCode ierr;
142: PetscBool flg;
143: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
144: Mat_MPIAIJCUSP *cuspStruct = (Mat_MPIAIJCUSP*)a->spptr;
147: MatSetFromOptions_MPIAIJ(PetscOptionsObject,A);
149: PetscOptionsHead(PetscOptionsObject,"MPIAIJCUSP options");
150: PetscObjectOptionsBegin((PetscObject)A);
151: if (A->factortype==MAT_FACTOR_NONE) {
152: PetscOptionsEnum("-mat_cusp_mult_diag_storage_format","sets storage format of the diagonal blocks of (mpi)aijcusp gpu matrices for SpMV",
153: "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)cuspStruct->diagGPUMatFormat,(PetscEnum*)&format,&flg);
154: if (flg) {
155: MatCUSPSetFormat(A,MAT_CUSP_MULT_DIAG,format);
156: }
157: PetscOptionsEnum("-mat_cusp_mult_offdiag_storage_format","sets storage format of the off-diagonal blocks (mpi)aijcusp gpu matrices for SpMV",
158: "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)cuspStruct->offdiagGPUMatFormat,(PetscEnum*)&format,&flg);
159: if (flg) {
160: MatCUSPSetFormat(A,MAT_CUSP_MULT_OFFDIAG,format);
161: }
162: PetscOptionsEnum("-mat_cusp_storage_format","sets storage format of the diagonal and off-diagonal blocks (mpi)aijcusp gpu matrices for SpMV",
163: "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)cuspStruct->diagGPUMatFormat,(PetscEnum*)&format,&flg);
164: if (flg) {
165: MatCUSPSetFormat(A,MAT_CUSP_ALL,format);
166: }
167: }
168: PetscOptionsEnd();
169: return(0);
170: }
174: PetscErrorCode MatDestroy_MPIAIJCUSP(Mat A)
175: {
177: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
178: Mat_MPIAIJCUSP *cuspStruct = (Mat_MPIAIJCUSP*)a->spptr;
179: cudaError_t err=cudaSuccess;
182: try {
183: err = cudaStreamDestroy(cuspStruct->stream);
184: if (err!=cudaSuccess) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Mat_MPIAIJCUSP error: %s", cudaGetErrorString(err));
185: delete cuspStruct;
186: } catch(char *ex) {
187: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Mat_MPIAIJCUSP error: %s", ex);
188: }
189: cuspStruct = 0;
190: MatDestroy_MPIAIJ(A);
191: return(0);
192: }
196: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJCUSP(Mat A)
197: {
199: Mat_MPIAIJ *a;
200: Mat_MPIAIJCUSP * cuspStruct;
201: cudaError_t err=cudaSuccess;
204: MatCreate_MPIAIJ(A);
205: PetscObjectComposeFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",MatMPIAIJSetPreallocation_MPIAIJCUSP);
206: A->ops->getvecs = MatCreateVecs_MPIAIJCUSP;
207: A->ops->setvaluesbatch = MatSetValuesBatch_MPIAIJCUSP;
209: a = (Mat_MPIAIJ*)A->data;
210: a->spptr = new Mat_MPIAIJCUSP;
211: cuspStruct = (Mat_MPIAIJCUSP*)a->spptr;
213: cuspStruct->diagGPUMatFormat = MAT_CUSP_CSR;
214: cuspStruct->offdiagGPUMatFormat = MAT_CUSP_CSR;
215: err = cudaStreamCreate(&(cuspStruct->stream));
216: if (err!=cudaSuccess) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Mat_MPIAIJCUSP error: %s", cudaGetErrorString(err));
218: A->ops->mult = MatMult_MPIAIJCUSP;
219: A->ops->setfromoptions = MatSetFromOptions_MPIAIJCUSP;
220: A->ops->destroy = MatDestroy_MPIAIJCUSP;
222: PetscObjectComposeFunction((PetscObject)A,"MatCUSPSetFormat_C", MatCUSPSetFormat_MPIAIJCUSP);
223: PetscObjectChangeTypeName((PetscObject)A,MATMPIAIJCUSP);
224: return(0);
225: }
228: /*@
229: MatCreateAIJCUSP - Creates a sparse matrix in AIJ (compressed row) format
230: (the default parallel PETSc format). This matrix will ultimately pushed down
231: to NVidia GPUs and use the CUSP library for calculations. For good matrix
232: assembly performance the user should preallocate the matrix storage by setting
233: the parameter nz (or the array nnz). By setting these parameters accurately,
234: performance during matrix assembly can be increased by more than a factor of 50.
237: Collective on MPI_Comm
239: Input Parameters:
240: + comm - MPI communicator, set to PETSC_COMM_SELF
241: . m - number of rows
242: . n - number of columns
243: . nz - number of nonzeros per row (same for all rows)
244: - nnz - array containing the number of nonzeros in the various rows
245: (possibly different for each row) or NULL
247: Output Parameter:
248: . A - the matrix
250: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
251: MatXXXXSetPreallocation() paradigm instead of this routine directly.
252: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
254: Notes:
255: If nnz is given then nz is ignored
257: The AIJ format (also called the Yale sparse matrix format or
258: compressed row storage), is fully compatible with standard Fortran 77
259: storage. That is, the stored row and column indices can begin at
260: either one (as in Fortran) or zero. See the users' manual for details.
262: Specify the preallocated storage with either nz or nnz (not both).
263: Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
264: allocation. For large problems you MUST preallocate memory or you
265: will get TERRIBLE performance, see the users' manual chapter on matrices.
267: By default, this format uses inodes (identical nodes) when possible, to
268: improve numerical efficiency of matrix-vector products and solves. We
269: search for consecutive rows with the same nonzero structure, thereby
270: reusing matrix information to achieve increased efficiency.
272: Level: intermediate
274: .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATMPIAIJCUSP, MATAIJCUSP
275: @*/
278: 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)
279: {
281: PetscMPIInt size;
284: MatCreate(comm,A);
285: MatSetSizes(*A,m,n,M,N);
286: MPI_Comm_size(comm,&size);
287: if (size > 1) {
288: MatSetType(*A,MATMPIAIJCUSP);
289: MatMPIAIJSetPreallocation(*A,d_nz,d_nnz,o_nz,o_nnz);
290: } else {
291: MatSetType(*A,MATSEQAIJCUSP);
292: MatSeqAIJSetPreallocation(*A,d_nz,d_nnz);
293: }
294: return(0);
295: }
297: /*M
298: MATAIJCUSP - MATMPIAIJCUSP= "aijcusp" = "mpiaijcusp" - A matrix type to be used for sparse matrices.
300: A matrix type type whose data resides on Nvidia GPUs. These matrices can be CSR format.
301: All matrix calculations are performed using the CUSP library. DIA and ELL
302: formats are also available
304: This matrix type is identical to MATSEQAIJCUSP when constructed with a single process communicator,
305: and MATMPIAIJCUSP otherwise. As a result, for single process communicators,
306: MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
307: for communicators controlling multiple processes. It is recommended that you call both of
308: the above preallocation routines for simplicity.
310: Options Database Keys:
311: + -mat_type mpiaijcusp - sets the matrix type to "mpiaijcusp" during a call to MatSetFromOptions()
312: . -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).
313: . -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).
314: - -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).
316: Level: beginner
318: .seealso: MatCreateAIJCUSP(), MATSEQAIJCUSP, MatCreateSeqAIJCUSP(), MatCUSPSetFormat(), MatCUSPStorageFormat, MatCUSPFormatOperation
319: M*/