Actual source code: aijviennacl.cxx

petsc-3.10.5 2019-03-28
Report Typos and Errors


  3: /*
  4:     Defines the basic matrix operations for the AIJ (compressed row)
  5:   matrix storage format.
  6: */

  8: #include <petscconf.h>
  9:  #include <../src/mat/impls/aij/seq/aij.h>
 10:  #include <petscbt.h>
 11:  #include <../src/vec/vec/impls/dvecimpl.h>
 12:  #include <petsc/private/vecimpl.h>

 14:  #include <../src/mat/impls/aij/seq/seqviennacl/viennaclmatimpl.h>


 17: #include <algorithm>
 18: #include <vector>
 19: #include <string>

 21: #include "viennacl/linalg/prod.hpp"

 23: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat A, MatType type, MatReuse reuse, Mat *newmat);
 24: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_petsc(Mat,MatFactorType,Mat*);


 27: PetscErrorCode MatViennaCLCopyToGPU(Mat A)
 28: {

 30:   Mat_SeqAIJViennaCL *viennaclstruct = (Mat_SeqAIJViennaCL*)A->spptr;
 31:   Mat_SeqAIJ         *a              = (Mat_SeqAIJ*)A->data;
 32:   PetscErrorCode     ierr;

 35:   if (A->rmap->n > 0 && A->cmap->n > 0 && a->nz) { //some OpenCL SDKs have issues with buffers of size 0
 36:     if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) {
 37:       PetscLogEventBegin(MAT_ViennaCLCopyToGPU,A,0,0,0);

 39:       try {
 40:         if (a->compressedrow.use) {
 41:           if (!viennaclstruct->compressed_mat) viennaclstruct->compressed_mat = new ViennaCLCompressedAIJMatrix();

 43:           // Since PetscInt is different from cl_uint, we have to convert:
 44:           viennacl::backend::mem_handle dummy;

 46:           viennacl::backend::typesafe_host_array<unsigned int> row_buffer; row_buffer.raw_resize(dummy, a->compressedrow.nrows+1);
 47:           for (PetscInt i=0; i<=a->compressedrow.nrows; ++i)
 48:             row_buffer.set(i, (a->compressedrow.i)[i]);

 50:           viennacl::backend::typesafe_host_array<unsigned int> row_indices; row_indices.raw_resize(dummy, a->compressedrow.nrows);
 51:           for (PetscInt i=0; i<a->compressedrow.nrows; ++i)
 52:             row_indices.set(i, (a->compressedrow.rindex)[i]);

 54:           viennacl::backend::typesafe_host_array<unsigned int> col_buffer; col_buffer.raw_resize(dummy, a->nz);
 55:           for (PetscInt i=0; i<a->nz; ++i)
 56:             col_buffer.set(i, (a->j)[i]);

 58:           viennaclstruct->compressed_mat->set(row_buffer.get(), row_indices.get(), col_buffer.get(), a->a, A->rmap->n, A->cmap->n, a->compressedrow.nrows, a->nz);
 59:         } else {
 60:           if (!viennaclstruct->mat) viennaclstruct->mat = new ViennaCLAIJMatrix();

 62:           // Since PetscInt is in general different from cl_uint, we have to convert:
 63:           viennacl::backend::mem_handle dummy;

 65:           viennacl::backend::typesafe_host_array<unsigned int> row_buffer; row_buffer.raw_resize(dummy, A->rmap->n+1);
 66:           for (PetscInt i=0; i<=A->rmap->n; ++i)
 67:             row_buffer.set(i, (a->i)[i]);

 69:           viennacl::backend::typesafe_host_array<unsigned int> col_buffer; col_buffer.raw_resize(dummy, a->nz);
 70:           for (PetscInt i=0; i<a->nz; ++i)
 71:             col_buffer.set(i, (a->j)[i]);

 73:           viennaclstruct->mat->set(row_buffer.get(), col_buffer.get(), a->a, A->rmap->n, A->cmap->n, a->nz);
 74:         }
 75:         ViennaCLWaitForGPU();
 76:       } catch(std::exception const & ex) {
 77:         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ViennaCL error: %s", ex.what());
 78:       }

 80:       // Create temporary vector for v += A*x:
 81:       if (viennaclstruct->tempvec) {
 82:         if (viennaclstruct->tempvec->size() != static_cast<std::size_t>(A->rmap->n)) {
 83:           delete (ViennaCLVector*)viennaclstruct->tempvec;
 84:           viennaclstruct->tempvec = new ViennaCLVector(A->rmap->n);
 85:         } else {
 86:           viennaclstruct->tempvec->clear();
 87:         }
 88:       } else {
 89:         viennaclstruct->tempvec = new ViennaCLVector(A->rmap->n);
 90:       }

 92:       A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH;

 94:       PetscLogEventEnd(MAT_ViennaCLCopyToGPU,A,0,0,0);
 95:     }
 96:   }
 97:   return(0);
 98: }

100: PetscErrorCode MatViennaCLCopyFromGPU(Mat A, const ViennaCLAIJMatrix *Agpu)
101: {
102:   Mat_SeqAIJ         *a              = (Mat_SeqAIJ*)A->data;
103:   PetscInt           m               = A->rmap->n;
104:   PetscErrorCode     ierr;


108:   if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED) {
109:     try {
110:       if (a->compressedrow.use) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "ViennaCL: Cannot handle row compression for GPU matrices");
111:       else {

113:         if ((PetscInt)Agpu->size1() != m) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "GPU matrix has %d rows, should be %d", Agpu->size1(), m);
114:         a->nz           = Agpu->nnz();
115:         a->maxnz        = a->nz; /* Since we allocate exactly the right amount */
116:         A->preallocated = PETSC_TRUE;
117:         if (a->singlemalloc) {
118:           if (a->a) {PetscFree3(a->a,a->j,a->i);}
119:         } else {
120:           if (a->i) {PetscFree(a->i);}
121:           if (a->j) {PetscFree(a->j);}
122:           if (a->a) {PetscFree(a->a);}
123:         }
124:         PetscMalloc3(a->nz,&a->a,a->nz,&a->j,m+1,&a->i);
125:         PetscLogObjectMemory((PetscObject)A, a->nz*(sizeof(PetscScalar)+sizeof(PetscInt))+(m+1)*sizeof(PetscInt));

127:         a->singlemalloc = PETSC_TRUE;

129:         /* Setup row lengths */
130:         if (a->imax) {PetscFree2(a->imax,a->ilen);}
131:         PetscMalloc2(m,&a->imax,m,&a->ilen);
132:         PetscLogObjectMemory((PetscObject)A, 2*m*sizeof(PetscInt));

134:         /* Copy data back from GPU */
135:         viennacl::backend::typesafe_host_array<unsigned int> row_buffer; row_buffer.raw_resize(Agpu->handle1(), Agpu->size1() + 1);

137:         // copy row array
138:         viennacl::backend::memory_read(Agpu->handle1(), 0, row_buffer.raw_size(), row_buffer.get());
139:         (a->i)[0] = row_buffer[0];
140:         for (PetscInt i = 0; i < (PetscInt)Agpu->size1(); ++i) {
141:           (a->i)[i+1] = row_buffer[i+1];
142:           a->imax[i]  = a->ilen[i] = a->i[i+1] - a->i[i];  //Set imax[] and ilen[] arrays at the same time as i[] for better cache reuse
143:         }

145:         // copy column indices
146:         viennacl::backend::typesafe_host_array<unsigned int> col_buffer; col_buffer.raw_resize(Agpu->handle2(), Agpu->nnz());
147:         viennacl::backend::memory_read(Agpu->handle2(), 0, col_buffer.raw_size(), col_buffer.get());
148:         for (PetscInt i=0; i < (PetscInt)Agpu->nnz(); ++i)
149:           (a->j)[i] = col_buffer[i];

151:         // copy nonzero entries directly to destination (no conversion required)
152:         viennacl::backend::memory_read(Agpu->handle(), 0, sizeof(PetscScalar)*Agpu->nnz(), a->a);

154:         ViennaCLWaitForGPU();
155:         /* TODO: Once a->diag is moved out of MatAssemblyEnd(), invalidate it here. */
156:       }
157:     } catch(std::exception const & ex) {
158:       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex.what());
159:     }

161:     /* This assembly prevents resetting the flag to PETSC_OFFLOAD_CPU and recopying */
162:     MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
163:     MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);

165:     A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH;
166:   } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "ViennaCL error: Only valid for unallocated GPU matrices");
167:   return(0);
168: }

170: PetscErrorCode MatMult_SeqAIJViennaCL(Mat A,Vec xx,Vec yy)
171: {
172:   Mat_SeqAIJ           *a = (Mat_SeqAIJ*)A->data;
173:   PetscErrorCode       ierr;
174:   Mat_SeqAIJViennaCL   *viennaclstruct = (Mat_SeqAIJViennaCL*)A->spptr;
175:   const ViennaCLVector *xgpu=NULL;
176:   ViennaCLVector       *ygpu=NULL;

179:   if (A->rmap->n > 0 && A->cmap->n > 0 && a->nz) {
180:     VecViennaCLGetArrayRead(xx,&xgpu);
181:     VecViennaCLGetArrayWrite(yy,&ygpu);
182:     try {
183:       if (a->compressedrow.use) {
184:         *ygpu = viennacl::linalg::prod(*viennaclstruct->compressed_mat, *xgpu);
185:       } else {
186:         *ygpu = viennacl::linalg::prod(*viennaclstruct->mat,*xgpu);
187:       }
188:       ViennaCLWaitForGPU();
189:     } catch (std::exception const & ex) {
190:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ViennaCL error: %s", ex.what());
191:     }
192:     VecViennaCLRestoreArrayRead(xx,&xgpu);
193:     VecViennaCLRestoreArrayWrite(yy,&ygpu);
194:     PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);
195:   } else {
196:     VecSet(yy,0);
197:   }
198:   return(0);
199: }

201: PetscErrorCode MatMultAdd_SeqAIJViennaCL(Mat A,Vec xx,Vec yy,Vec zz)
202: {
203:   Mat_SeqAIJ           *a = (Mat_SeqAIJ*)A->data;
204:   PetscErrorCode       ierr;
205:   Mat_SeqAIJViennaCL   *viennaclstruct = (Mat_SeqAIJViennaCL*)A->spptr;
206:   const ViennaCLVector *xgpu=NULL,*ygpu=NULL;
207:   ViennaCLVector       *zgpu=NULL;

210:   if (A->rmap->n > 0 && A->cmap->n > 0 && a->nz) {
211:     try {
212:       VecViennaCLGetArrayRead(xx,&xgpu);
213:       VecViennaCLGetArrayRead(yy,&ygpu);
214:       VecViennaCLGetArrayWrite(zz,&zgpu);

216:       if (a->compressedrow.use) {
217:         ViennaCLVector temp = viennacl::linalg::prod(*viennaclstruct->compressed_mat, *xgpu);
218:         *zgpu = *ygpu + temp;
219:         ViennaCLWaitForGPU();
220:       } else {
221:         if (zz == xx || zz == yy) { //temporary required
222:           ViennaCLVector temp = viennacl::linalg::prod(*viennaclstruct->mat, *xgpu);
223:           *zgpu = *ygpu;
224:           *zgpu += temp;
225:           ViennaCLWaitForGPU();
226:         } else {
227:           *viennaclstruct->tempvec = viennacl::linalg::prod(*viennaclstruct->mat, *xgpu);
228:           *zgpu = *ygpu + *viennaclstruct->tempvec;
229:           ViennaCLWaitForGPU();
230:         }
231:       }

233:       VecViennaCLRestoreArrayRead(xx,&xgpu);
234:       VecViennaCLRestoreArrayRead(yy,&ygpu);
235:       VecViennaCLRestoreArrayWrite(zz,&zgpu);

237:     } catch(std::exception const & ex) {
238:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ViennaCL error: %s", ex.what());
239:     }
240:     PetscLogFlops(2.0*a->nz);
241:   } else {
242:     VecCopy(yy,zz);
243:   }
244:   return(0);
245: }

247: PetscErrorCode MatAssemblyEnd_SeqAIJViennaCL(Mat A,MatAssemblyType mode)
248: {

252:   MatAssemblyEnd_SeqAIJ(A,mode);
253:   MatViennaCLCopyToGPU(A);
254:   return(0);
255: }

257: /* --------------------------------------------------------------------------------*/
258: /*@
259:    MatCreateSeqAIJViennaCL - Creates a sparse matrix in AIJ (compressed row) format
260:    (the default parallel PETSc format).  This matrix will ultimately be pushed down
261:    to GPUs and use the ViennaCL library for calculations. For good matrix
262:    assembly performance the user should preallocate the matrix storage by setting
263:    the parameter nz (or the array nnz).  By setting these parameters accurately,
264:    performance during matrix assembly can be increased substantially.


267:    Collective on MPI_Comm

269:    Input Parameters:
270: +  comm - MPI communicator, set to PETSC_COMM_SELF
271: .  m - number of rows
272: .  n - number of columns
273: .  nz - number of nonzeros per row (same for all rows)
274: -  nnz - array containing the number of nonzeros in the various rows
275:          (possibly different for each row) or NULL

277:    Output Parameter:
278: .  A - the matrix

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

284:    Notes:
285:    If nnz is given then nz is ignored

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

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

297:    Level: intermediate

299: .seealso: MatCreate(), MatCreateAIJ(), MatCreateAIJCUSPARSE(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ()

301: @*/
302: PetscErrorCode  MatCreateSeqAIJViennaCL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
303: {

307:   MatCreate(comm,A);
308:   MatSetSizes(*A,m,n,m,n);
309:   MatSetType(*A,MATSEQAIJVIENNACL);
310:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
311:   return(0);
312: }


315: PetscErrorCode MatDestroy_SeqAIJViennaCL(Mat A)
316: {
318:   Mat_SeqAIJViennaCL *viennaclcontainer = (Mat_SeqAIJViennaCL*)A->spptr;

321:   try {
322:     if (viennaclcontainer) {
323:       delete viennaclcontainer->tempvec;
324:       delete viennaclcontainer->mat;
325:       delete viennaclcontainer->compressed_mat;
326:       delete viennaclcontainer;
327:     }
328:     A->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED;
329:   } catch(std::exception const & ex) {
330:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ViennaCL error: %s", ex.what());
331:   }

333:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqaijviennacl_C",NULL);

335:   /* this next line is because MatDestroy tries to PetscFree spptr if it is not zero, and PetscFree only works if the memory was allocated with PetscNew or PetscMalloc, which don't call the constructor */
336:   A->spptr = 0;
337:   MatDestroy_SeqAIJ(A);
338:   return(0);
339: }


342: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJViennaCL(Mat B)
343: {

347:   MatCreate_SeqAIJ(B);
348:   MatConvert_SeqAIJ_SeqAIJViennaCL(B,MATSEQAIJVIENNACL,MAT_INPLACE_MATRIX,&B);
349:   return(0);
350: }

352: static PetscErrorCode MatDuplicate_SeqAIJViennaCL(Mat A,MatDuplicateOption cpvalues,Mat *B)
353: {
355:   Mat C;

358:   MatDuplicate_SeqAIJ(A,cpvalues,B);
359:   C = *B;

361:   C->ops->mult        = MatMult_SeqAIJViennaCL;
362:   C->ops->multadd     = MatMultAdd_SeqAIJViennaCL;
363:   C->ops->assemblyend = MatAssemblyEnd_SeqAIJViennaCL;
364:   C->ops->destroy     = MatDestroy_SeqAIJViennaCL;
365:   C->ops->duplicate   = MatDuplicate_SeqAIJViennaCL;

367:   C->spptr        = new Mat_SeqAIJViennaCL();
368:   ((Mat_SeqAIJViennaCL*)C->spptr)->tempvec        = NULL;
369:   ((Mat_SeqAIJViennaCL*)C->spptr)->mat            = NULL;
370:   ((Mat_SeqAIJViennaCL*)C->spptr)->compressed_mat = NULL;

372:   PetscObjectChangeTypeName((PetscObject)C,MATSEQAIJVIENNACL);

374:   C->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED;

376:   /* If the source matrix is already assembled, copy the destination matrix to the GPU */
377:   if (C->assembled) {
378:     MatViennaCLCopyToGPU(C);
379:   }

381:   return(0);
382: }

384: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
385: {
387:   Mat            B;
388:   Mat_SeqAIJ     *aij;


392:   if (reuse == MAT_REUSE_MATRIX) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MAT_REUSE_MATRIX is not supported. Consider using MAT_INPLACE_MATRIX instead");

394:   if (reuse == MAT_INITIAL_MATRIX) {
395:     MatDuplicate(A,MAT_COPY_VALUES,newmat);
396:   }

398:   B = *newmat;

400:   aij             = (Mat_SeqAIJ*)B->data;
401:   aij->inode.use  = PETSC_FALSE;

403:   B->ops->mult    = MatMult_SeqAIJViennaCL;
404:   B->ops->multadd = MatMultAdd_SeqAIJViennaCL;
405:   B->spptr        = new Mat_SeqAIJViennaCL();

407:   ((Mat_SeqAIJViennaCL*)B->spptr)->tempvec        = NULL;
408:   ((Mat_SeqAIJViennaCL*)B->spptr)->mat            = NULL;
409:   ((Mat_SeqAIJViennaCL*)B->spptr)->compressed_mat = NULL;

411:   B->ops->assemblyend    = MatAssemblyEnd_SeqAIJViennaCL;
412:   B->ops->destroy        = MatDestroy_SeqAIJViennaCL;
413:   B->ops->duplicate      = MatDuplicate_SeqAIJViennaCL;

415:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJVIENNACL);
416:   PetscFree(B->defaultvectype);
417:   PetscStrallocpy(VECVIENNACL,&B->defaultvectype);

419:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqaijviennacl_C",MatConvert_SeqAIJ_SeqAIJViennaCL);

421:   B->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED;

423:   /* If the source matrix is already assembled, copy the destination matrix to the GPU */
424:   if (B->assembled) {
425:     MatViennaCLCopyToGPU(B);
426:   }

428:   return(0);
429: }


432: /*MC
433:    MATSEQAIJVIENNACL - MATAIJVIENNACL = "aijviennacl" = "seqaijviennacl" - A matrix type to be used for sparse matrices.

435:    A matrix type type whose data resides on GPUs. These matrices are in CSR format by
436:    default. All matrix calculations are performed using the ViennaCL library.

438:    Options Database Keys:
439: +  -mat_type aijviennacl - sets the matrix type to "seqaijviennacl" during a call to MatSetFromOptions()
440: .  -mat_viennacl_storage_format csr - sets the storage format of matrices for MatMult during a call to MatSetFromOptions().
441: -  -mat_viennacl_mult_storage_format csr - sets the storage format of matrices for MatMult during a call to MatSetFromOptions().

443:   Level: beginner

445: .seealso: MatCreateSeqAIJViennaCL(), MATAIJVIENNACL, MatCreateAIJViennaCL()
446: M*/


449: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_ViennaCL(void)
450: {

454:   MatSolverTypeRegister(MATSOLVERPETSC, MATSEQAIJVIENNACL,    MAT_FACTOR_LU,MatGetFactor_seqaij_petsc);
455:   MatSolverTypeRegister(MATSOLVERPETSC, MATSEQAIJVIENNACL,    MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_petsc);
456:   MatSolverTypeRegister(MATSOLVERPETSC, MATSEQAIJVIENNACL,    MAT_FACTOR_ILU,MatGetFactor_seqaij_petsc);
457:   MatSolverTypeRegister(MATSOLVERPETSC, MATSEQAIJVIENNACL,    MAT_FACTOR_ICC,MatGetFactor_seqaij_petsc);
458:   return(0);
459: }