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