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