1 #define PETSCMAT_DLL 2 3 /* 4 Defines a matrix-vector product for the MATMPIAIJCRL matrix class. 5 This class is derived from the MATMPIAIJ class and retains the 6 compressed row storage (aka Yale sparse matrix format) but augments 7 it with a column oriented storage that is more efficient for 8 matrix vector products on Vector machines. 9 10 CRL stands for constant row length (that is the same number of columns 11 is kept (padded with zeros) for each row of the sparse matrix. 12 13 See src/mat/impls/aij/seq/crl/crl.c for the sequential version 14 */ 15 16 #include "../src/mat/impls/aij/mpi/mpiaij.h" 17 #include "../src/mat/impls/aij/seq/crl/crl.h" 18 19 extern PetscErrorCode MatDestroy_MPIAIJ(Mat); 20 21 #undef __FUNCT__ 22 #define __FUNCT__ "MatDestroy_MPICRL" 23 PetscErrorCode MatDestroy_MPICRL(Mat A) 24 { 25 PetscErrorCode ierr; 26 Mat_CRL *crl = (Mat_CRL *) A->spptr; 27 28 /* Free everything in the Mat_CRL data structure. */ 29 ierr = PetscFree2(crl->acols,crl->icols);CHKERRQ(ierr); 30 if (crl->fwork) { 31 ierr = VecDestroy(crl->fwork);CHKERRQ(ierr); 32 } 33 if (crl->xwork) { 34 ierr = VecDestroy(crl->xwork);CHKERRQ(ierr); 35 } 36 ierr = PetscFree(crl->array);CHKERRQ(ierr); 37 ierr = PetscFree(crl);CHKERRQ(ierr); 38 A->spptr = 0; 39 40 ierr = PetscObjectChangeTypeName( (PetscObject)A, MATMPIAIJ);CHKERRQ(ierr); 41 ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 42 PetscFunctionReturn(0); 43 } 44 45 #undef __FUNCT__ 46 #define __FUNCT__ "MPICRL_create_crl" 47 PetscErrorCode MPICRL_create_crl(Mat A) 48 { 49 Mat_MPIAIJ *a = (Mat_MPIAIJ *)(A)->data; 50 Mat_SeqAIJ *Aij = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->B->data); 51 Mat_CRL *crl = (Mat_CRL*) A->spptr; 52 PetscInt m = A->rmap->n; /* Number of rows in the matrix. */ 53 PetscInt nd = a->A->cmap->n; /* number of columns in diagonal portion */ 54 PetscInt *aj = Aij->j,*bj = Bij->j; /* From the CSR representation; points to the beginning of each row. */ 55 PetscInt i, j,rmax = 0,*icols, *ailen = Aij->ilen, *bilen = Bij->ilen; 56 PetscScalar *aa = Aij->a,*ba = Bij->a,*acols,*array; 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 /* determine the row with the most columns */ 61 for (i=0; i<m; i++) { 62 rmax = PetscMax(rmax,ailen[i]+bilen[i]); 63 } 64 crl->nz = Aij->nz+Bij->nz; 65 crl->m = A->rmap->n; 66 crl->rmax = rmax; 67 ierr = PetscFree2(crl->acols,crl->icols);CHKERRQ(ierr); 68 ierr = PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);CHKERRQ(ierr); 69 acols = crl->acols; 70 icols = crl->icols; 71 for (i=0; i<m; i++) { 72 for (j=0; j<ailen[i]; j++) { 73 acols[j*m+i] = *aa++; 74 icols[j*m+i] = *aj++; 75 } 76 for (;j<ailen[i]+bilen[i]; j++) { 77 acols[j*m+i] = *ba++; 78 icols[j*m+i] = nd + *bj++; 79 } 80 for (;j<rmax; j++) { /* empty column entries */ 81 acols[j*m+i] = 0.0; 82 icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0; /* handle case where row is EMPTY */ 83 } 84 } 85 ierr = PetscInfo1(A,"Percentage of 0's introduced for vectorized multiply %g\n",1.0-((double)(crl->nz))/((double)(rmax*m))); 86 87 ierr = PetscFree(crl->array);CHKERRQ(ierr); 88 ierr = PetscMalloc((a->B->cmap->n+nd)*sizeof(PetscScalar),&array);CHKERRQ(ierr); 89 /* xwork array is actually B->n+nd long, but we define xwork this length so can copy into it */ 90 if (crl->xwork) {ierr = VecDestroy(crl->xwork);CHKERRQ(ierr);} 91 ierr = VecCreateMPIWithArray(((PetscObject)A)->comm,nd,PETSC_DECIDE,array,&crl->xwork);CHKERRQ(ierr); 92 if (crl->fwork) {ierr = VecDestroy(crl->fwork);CHKERRQ(ierr);} 93 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->cmap->n,array+nd,&crl->fwork);CHKERRQ(ierr); 94 crl->array = array; 95 crl->xscat = a->Mvctx; 96 PetscFunctionReturn(0); 97 } 98 99 extern PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat,MatAssemblyType); 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "MatAssemblyEnd_MPICRL" 103 PetscErrorCode MatAssemblyEnd_MPICRL(Mat A, MatAssemblyType mode) 104 { 105 PetscErrorCode ierr; 106 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 107 Mat_SeqAIJ *Aij = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->A->data); 108 109 PetscFunctionBegin; 110 Aij->inode.use = PETSC_FALSE; 111 Bij->inode.use = PETSC_FALSE; 112 ierr = MatAssemblyEnd_MPIAIJ(A,mode);CHKERRQ(ierr); 113 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 114 115 /* Now calculate the permutation and grouping information. */ 116 ierr = MPICRL_create_crl(A);CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119 120 extern PetscErrorCode MatMult_CRL(Mat,Vec,Vec); 121 extern PetscErrorCode MatDuplicate_CRL(Mat,MatDuplicateOption,Mat*); 122 123 /* MatConvert_MPIAIJ_MPICRL converts a MPIAIJ matrix into a 124 * MPICRL matrix. This routine is called by the MatCreate_MPICRL() 125 * routine, but can also be used to convert an assembled MPIAIJ matrix 126 * into a MPICRL one. */ 127 EXTERN_C_BEGIN 128 #undef __FUNCT__ 129 #define __FUNCT__ "MatConvert_MPIAIJ_MPICRL" 130 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MPIAIJ_MPICRL(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 131 { 132 PetscErrorCode ierr; 133 Mat B = *newmat; 134 Mat_CRL *crl; 135 136 PetscFunctionBegin; 137 if (reuse == MAT_INITIAL_MATRIX) { 138 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 139 } 140 141 ierr = PetscNewLog(B,Mat_CRL,&crl);CHKERRQ(ierr); 142 B->spptr = (void *) crl; 143 144 /* Set function pointers for methods that we inherit from AIJ but override. */ 145 B->ops->duplicate = MatDuplicate_CRL; 146 B->ops->assemblyend = MatAssemblyEnd_MPICRL; 147 B->ops->destroy = MatDestroy_MPICRL; 148 B->ops->mult = MatMult_CRL; 149 150 /* If A has already been assembled, compute the permutation. */ 151 if (A->assembled) { 152 ierr = MPICRL_create_crl(B);CHKERRQ(ierr); 153 } 154 ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPICRL);CHKERRQ(ierr); 155 *newmat = B; 156 PetscFunctionReturn(0); 157 } 158 EXTERN_C_END 159 160 161 #undef __FUNCT__ 162 #define __FUNCT__ "MatCreateMPICRL" 163 /*@C 164 MatCreateMPICRL - Creates a sparse matrix of type MPICRL. 165 This type inherits from AIJ, but stores some additional 166 information that is used to allow better vectorization of 167 the matrix-vector product. At the cost of increased storage, the AIJ formatted 168 matrix can be copied to a format in which pieces of the matrix are 169 stored in ELLPACK format, allowing the vectorized matrix multiply 170 routine to use stride-1 memory accesses. As with the AIJ type, it is 171 important to preallocate matrix storage in order to get good assembly 172 performance. 173 174 Collective on MPI_Comm 175 176 Input Parameters: 177 + comm - MPI communicator, set to PETSC_COMM_SELF 178 . m - number of rows 179 . n - number of columns 180 . nz - number of nonzeros per row (same for all rows) 181 - nnz - array containing the number of nonzeros in the various rows 182 (possibly different for each row) or PETSC_NULL 183 184 Output Parameter: 185 . A - the matrix 186 187 Notes: 188 If nnz is given then nz is ignored 189 190 Level: intermediate 191 192 .keywords: matrix, cray, sparse, parallel 193 194 .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues() 195 @*/ 196 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[],Mat *A) 197 { 198 PetscErrorCode ierr; 199 200 PetscFunctionBegin; 201 ierr = MatCreate(comm,A);CHKERRQ(ierr); 202 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 203 ierr = MatSetType(*A,MATMPICRL);CHKERRQ(ierr); 204 ierr = MatMPIAIJSetPreallocation_MPIAIJ(*A,nz,(PetscInt*)nnz,onz,(PetscInt*)onnz);CHKERRQ(ierr); 205 PetscFunctionReturn(0); 206 } 207 208 209 EXTERN_C_BEGIN 210 #undef __FUNCT__ 211 #define __FUNCT__ "MatCreate_MPICRL" 212 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPICRL(Mat A) 213 { 214 PetscErrorCode ierr; 215 216 PetscFunctionBegin; 217 ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 218 ierr = MatConvert_MPIAIJ_MPICRL(A,MATMPICRL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 219 PetscFunctionReturn(0); 220 } 221 EXTERN_C_END 222 223