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 #undef __FUNCT__ 20 #define __FUNCT__ "MatDestroy_MPICRL" 21 PetscErrorCode MatDestroy_MPICRL(Mat A) 22 { 23 PetscErrorCode ierr; 24 Mat_CRL *crl = (Mat_CRL *) A->spptr; 25 26 /* We are going to convert A back into a MPIAIJ matrix, since we are 27 * eventually going to use MatDestroy_MPIAIJ() to destroy everything 28 * that is not specific to CRL. 29 * In preparation for this, reset the operations pointers in A to 30 * their MPIAIJ versions. */ 31 A->ops->assemblyend = crl->AssemblyEnd; 32 A->ops->destroy = crl->MatDestroy; 33 A->ops->duplicate = crl->MatDuplicate; 34 35 /* Free everything in the Mat_CRL data structure. */ 36 if (crl->icols) { 37 ierr = PetscFree2(crl->acols,crl->icols);CHKERRQ(ierr); 38 } 39 if (crl->fwork) { 40 ierr = VecDestroy(crl->fwork);CHKERRQ(ierr); 41 } 42 if (crl->xwork) { 43 ierr = VecDestroy(crl->xwork);CHKERRQ(ierr); 44 } 45 if (crl->array) { 46 ierr = PetscFree(crl->array);CHKERRQ(ierr); 47 } 48 /* Free the Mat_CRL struct itself. */ 49 ierr = PetscFree(crl);CHKERRQ(ierr); 50 51 /* Change the type of A back to MPIAIJ and use MatDestroy_MPIAIJ() 52 * to destroy everything that remains. */ 53 ierr = PetscObjectChangeTypeName( (PetscObject)A, MATMPIAIJ);CHKERRQ(ierr); 54 /* Note that I don't call MatSetType(). I believe this is because that 55 * is only to be called when *building* a matrix. */ 56 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 57 PetscFunctionReturn(0); 58 } 59 60 #undef __FUNCT__ 61 #define __FUNCT__ "MPICRL_create_crl" 62 PetscErrorCode MPICRL_create_crl(Mat A) 63 { 64 Mat_MPIAIJ *a = (Mat_MPIAIJ *)(A)->data; 65 Mat_SeqAIJ *Aij = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->B->data); 66 Mat_CRL *crl = (Mat_CRL*) A->spptr; 67 PetscInt m = A->m; /* Number of rows in the matrix. */ 68 PetscInt nd = a->A->n; /* number of columns in diagonal portion */ 69 PetscInt *aj = Aij->j,*bj = Bij->j; /* From the CSR representation; points to the beginning of each row. */ 70 PetscInt i, j,rmax = 0,*icols, *ailen = Aij->ilen, *bilen = Bij->ilen; 71 PetscScalar *aa = Aij->a,*ba = Bij->a,*acols,*array; 72 PetscErrorCode ierr; 73 74 PetscFunctionBegin; 75 /* determine the row with the most columns */ 76 for (i=0; i<m; i++) { 77 rmax = PetscMax(rmax,ailen[i]+bilen[i]); 78 } 79 crl->nz = Aij->nz+Bij->nz; 80 crl->m = A->m; 81 crl->rmax = rmax; 82 ierr = PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);CHKERRQ(ierr); 83 acols = crl->acols; 84 icols = crl->icols; 85 for (i=0; i<m; i++) { 86 for (j=0; j<ailen[i]; j++) { 87 acols[j*m+i] = *aa++; 88 icols[j*m+i] = *aj++; 89 } 90 for (;j<ailen[i]+bilen[i]; j++) { 91 acols[j*m+i] = *ba++; 92 icols[j*m+i] = nd + *bj++; 93 } 94 for (;j<rmax; j++) { /* empty column entries */ 95 acols[j*m+i] = 0.0; 96 icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0; /* handle case where row is EMPTY */ 97 } 98 } 99 ierr = PetscVerboseInfo((A,"MPICRL_create_crl: Percentage of 0's introduced for vectorized multiply %g\n",1.0-((double)(crl->nz))/((double)(rmax*m)))); 100 101 ierr = PetscMalloc((a->B->n+nd)*sizeof(PetscScalar),&array);CHKERRQ(ierr); 102 /* xwork array is actually B->n+nd long, but we define xwork this length so can copy into it */ 103 ierr = VecCreateMPIWithArray(A->comm,nd,PETSC_DECIDE,array,&crl->xwork);CHKERRQ(ierr); 104 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,array+nd,&crl->fwork);CHKERRQ(ierr); 105 crl->array = array; 106 crl->xscat = a->Mvctx; 107 PetscFunctionReturn(0); 108 } 109 110 #undef __FUNCT__ 111 #define __FUNCT__ "MatAssemblyEnd_MPICRL" 112 PetscErrorCode MatAssemblyEnd_MPICRL(Mat A, MatAssemblyType mode) 113 { 114 PetscErrorCode ierr; 115 Mat_CRL *crl = (Mat_CRL*) A->spptr; 116 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 117 Mat_SeqAIJ *Aij = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->A->data); 118 119 PetscFunctionBegin; 120 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 121 122 /* Since a MATMPICRL matrix is really just a MATMPIAIJ with some 123 * extra information, call the AssemblyEnd routine for a MATMPIAIJ. 124 * I'm not sure if this is the best way to do this, but it avoids 125 * a lot of code duplication. 126 * I also note that currently MATMPICRL doesn't know anything about 127 * the Mat_CompressedRow data structure that MPIAIJ now uses when there 128 * are many zero rows. If the MPIAIJ assembly end routine decides to use 129 * this, this may break things. (Don't know... haven't looked at it.) */ 130 Aij->inode.use = PETSC_FALSE; 131 Bij->inode.use = PETSC_FALSE; 132 (*crl->AssemblyEnd)(A, mode); 133 134 /* Now calculate the permutation and grouping information. */ 135 ierr = MPICRL_create_crl(A);CHKERRQ(ierr); 136 PetscFunctionReturn(0); 137 } 138 139 extern PetscErrorCode MatMult_CRL(Mat,Vec,Vec); 140 extern PetscErrorCode MatDuplicate_CRL(Mat,MatDuplicateOption,Mat*); 141 142 /* MatConvert_MPIAIJ_MPICRL converts a MPIAIJ matrix into a 143 * MPICRL matrix. This routine is called by the MatCreate_MPICRL() 144 * routine, but can also be used to convert an assembled MPIAIJ matrix 145 * into a MPICRL one. */ 146 EXTERN_C_BEGIN 147 #undef __FUNCT__ 148 #define __FUNCT__ "MatConvert_MPIAIJ_MPICRL" 149 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MPIAIJ_MPICRL(Mat A,MatType type,MatReuse reuse,Mat *newmat) 150 { 151 /* This routine is only called to convert to MATMPICRL 152 * from MATMPIAIJ, so we can ignore 'MatType Type'. */ 153 PetscErrorCode ierr; 154 Mat B = *newmat; 155 Mat_CRL *crl; 156 157 PetscFunctionBegin; 158 if (reuse == MAT_INITIAL_MATRIX) { 159 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 160 } 161 162 ierr = PetscNew(Mat_CRL,&crl);CHKERRQ(ierr); 163 B->spptr = (void *) crl; 164 165 /* Save a pointer to the original MPIAIJ assembly end routine, because we 166 * will want to use it later in the CRL assembly end routine. 167 * Also, save a pointer to the original MPIAIJ Destroy routine, because we 168 * will want to use it in the CRL destroy routine. */ 169 crl->AssemblyEnd = A->ops->assemblyend; 170 crl->MatDestroy = A->ops->destroy; 171 crl->MatDuplicate = A->ops->duplicate; 172 173 /* Set function pointers for methods that we inherit from AIJ but 174 * override. */ 175 B->ops->duplicate = MatDuplicate_CRL; 176 B->ops->assemblyend = MatAssemblyEnd_MPICRL; 177 B->ops->destroy = MatDestroy_MPICRL; 178 B->ops->mult = MatMult_CRL; 179 180 /* If A has already been assembled, compute the permutation. */ 181 if (A->assembled == PETSC_TRUE) { 182 ierr = MPICRL_create_crl(B);CHKERRQ(ierr); 183 } 184 ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPICRL);CHKERRQ(ierr); 185 *newmat = B; 186 PetscFunctionReturn(0); 187 } 188 EXTERN_C_END 189 190 191 #undef __FUNCT__ 192 #define __FUNCT__ "MatCreateMPICRL" 193 /*@C 194 MatCreateMPICRL - Creates a sparse matrix of type MPICRL. 195 This type inherits from AIJ, but stores some additional 196 information that is used to allow better vectorization of 197 the matrix-vector product. At the cost of increased storage, the AIJ formatted 198 matrix can be copied to a format in which pieces of the matrix are 199 stored in ELLPACK format, allowing the vectorized matrix multiply 200 routine to use stride-1 memory accesses. As with the AIJ type, it is 201 important to preallocate matrix storage in order to get good assembly 202 performance. 203 204 Collective on MPI_Comm 205 206 Input Parameters: 207 + comm - MPI communicator, set to PETSC_COMM_SELF 208 . m - number of rows 209 . n - number of columns 210 . nz - number of nonzeros per row (same for all rows) 211 - nnz - array containing the number of nonzeros in the various rows 212 (possibly different for each row) or PETSC_NULL 213 214 Output Parameter: 215 . A - the matrix 216 217 Notes: 218 If nnz is given then nz is ignored 219 220 Level: intermediate 221 222 .keywords: matrix, cray, sparse, parallel 223 224 .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues() 225 @*/ 226 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[],Mat *A) 227 { 228 PetscErrorCode ierr; 229 230 PetscFunctionBegin; 231 ierr = MatCreate(comm,A);CHKERRQ(ierr); 232 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 233 ierr = MatSetType(*A,MATMPICRL);CHKERRQ(ierr); 234 ierr = MatMPIAIJSetPreallocation_MPIAIJ(*A,nz,(PetscInt*)nnz,onz,(PetscInt*)onnz);CHKERRQ(ierr); 235 PetscFunctionReturn(0); 236 } 237 238 239 EXTERN_C_BEGIN 240 #undef __FUNCT__ 241 #define __FUNCT__ "MatCreate_MPICRL" 242 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPICRL(Mat A) 243 { 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 /* Change the type name before calling MatSetType() to force proper construction of MPIAIJ 248 and MATMPICRL types. */ 249 ierr = PetscObjectChangeTypeName((PetscObject)A,MATMPICRL);CHKERRQ(ierr); 250 ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 251 ierr = MatConvert_MPIAIJ_MPICRL(A,MATMPICRL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 252 PetscFunctionReturn(0); 253 } 254 EXTERN_C_END 255 256