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