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