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