1 #define PETSCMAT_DLL 2 3 /* 4 Defines a matrix-vector product for the MATSEQAIJCRL matrix class. 5 This class is derived from the MATSEQAIJ 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 14 #include "src/mat/impls/aij/seq/aij.h" 15 16 typedef struct { 17 PetscInt ncols; /* number of columns in each row */ 18 PetscInt *icols; /* columns of nonzeros, stored one column at a time */ 19 PetscScalar *acols; /* values of nonzeros, stored as icols */ 20 21 /* We need to keep a pointer to MatAssemblyEnd_SeqAIJ because we 22 * actually want to call this function from within the 23 * MatAssemblyEnd_SeqCRL function. Similarly, we also need 24 * MatDestroy_SeqAIJ and MatDuplicate_SeqAIJ. */ 25 PetscErrorCode (*AssemblyEnd_SeqAIJ)(Mat,MatAssemblyType); 26 PetscErrorCode (*MatDestroy_SeqAIJ)(Mat); 27 PetscErrorCode (*MatDuplicate_SeqAIJ)(Mat,MatDuplicateOption,Mat*); 28 } Mat_SeqCRL; 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "MatDestroy_SeqCRL" 32 PetscErrorCode MatDestroy_SeqCRL(Mat A) 33 { 34 PetscErrorCode ierr; 35 Mat_SeqCRL *crl = (Mat_SeqCRL *) A->spptr; 36 37 /* We are going to convert A back into a SEQAIJ matrix, since we are 38 * eventually going to use MatDestroy_SeqAIJ() to destroy everything 39 * that is not specific to CRL. 40 * In preparation for this, reset the operations pointers in A to 41 * their SeqAIJ versions. */ 42 A->ops->assemblyend = crl->AssemblyEnd_SeqAIJ; 43 A->ops->destroy = crl->MatDestroy_SeqAIJ; 44 A->ops->duplicate = crl->MatDuplicate_SeqAIJ; 45 46 /* Free everything in the Mat_SeqCRL data structure. */ 47 if (crl->icols) { 48 ierr = PetscFree2(crl->acols,crl->icols);CHKERRQ(ierr); 49 } 50 /* Free the Mat_SeqCRL struct itself. */ 51 ierr = PetscFree(crl);CHKERRQ(ierr); 52 53 /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ() 54 * to destroy everything that remains. */ 55 ierr = PetscObjectChangeTypeName( (PetscObject)A, MATSEQAIJ);CHKERRQ(ierr); 56 /* Note that I don't call MatSetType(). I believe this is because that 57 * is only to be called when *building* a matrix. */ 58 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 59 PetscFunctionReturn(0); 60 } 61 62 PetscErrorCode MatDuplicate_SeqCRL(Mat A, MatDuplicateOption op, Mat *M) 63 { 64 PetscErrorCode ierr; 65 Mat_SeqCRL *crl = (Mat_SeqCRL *) A->spptr; 66 67 PetscFunctionBegin; 68 ierr = (*crl->MatDuplicate_SeqAIJ)(A,op,M);CHKERRQ(ierr); 69 SETERRQ(PETSC_ERR_SUP,"Cannot duplicate CRL matrices yet"); 70 PetscFunctionReturn(0); 71 } 72 73 #undef __FUNCT__ 74 #define __FUNCT__ "SeqCRL_create_crl" 75 PetscErrorCode SeqCRL_create_crl(Mat A) 76 { 77 Mat_SeqAIJ *a = (Mat_SeqAIJ *)(A)->data; 78 Mat_SeqCRL *crl = (Mat_SeqCRL*) A->spptr; 79 PetscInt m = A->m; /* Number of rows in the matrix. */ 80 PetscInt *aj = a->j; /* From the CSR representation; points to the beginning of each row. */ 81 PetscInt i, j,rmax = a->rmax,*icols, *ilen = a->ilen; 82 PetscScalar *aa = a->a,*acols; 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 ierr = PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);CHKERRQ(ierr); 87 acols = crl->acols; 88 icols = crl->icols; 89 for (i=0; i<m; i++) { 90 for (j=0; j<ilen[i]; j++) { 91 acols[j*m+i] = *aa++; 92 icols[j*m+i] = *aj++; 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 = PetscLogInfo((A,"SeqCRL_create_crl: Percentage of 0's introduced for vectorized multiply %g\n",1.0-((double)a->nz)/((double)(rmax*m)))); 100 PetscFunctionReturn(0); 101 } 102 103 #undef __FUNCT__ 104 #define __FUNCT__ "MatAssemblyEnd_SeqCRL" 105 PetscErrorCode MatAssemblyEnd_SeqCRL(Mat A, MatAssemblyType mode) 106 { 107 PetscErrorCode ierr; 108 Mat_SeqCRL *crl = (Mat_SeqCRL*) A->spptr; 109 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 110 111 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 112 113 /* Since a MATSEQCRL matrix is really just a MATSEQAIJ with some 114 * extra information, call the AssemblyEnd routine for a MATSEQAIJ. 115 * I'm not sure if this is the best way to do this, but it avoids 116 * a lot of code duplication. 117 * I also note that currently MATSEQCRL doesn't know anything about 118 * the Mat_CompressedRow data structure that SeqAIJ now uses when there 119 * are many zero rows. If the SeqAIJ assembly end routine decides to use 120 * this, this may break things. (Don't know... haven't looked at it.) */ 121 a->inode.use = PETSC_FALSE; 122 (*crl->AssemblyEnd_SeqAIJ)(A, mode); 123 124 /* Now calculate the permutation and grouping information. */ 125 ierr = SeqCRL_create_crl(A);CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 #include "src/inline/dot.h" 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "MatMult_SeqCRL" 133 PetscErrorCode MatMult_SeqCRL(Mat A,Vec xx,Vec yy) 134 { 135 Mat_SeqAIJ *a = (Mat_SeqAIJ *)(A)->data; 136 Mat_SeqCRL *crl = (Mat_SeqCRL*) A->spptr; 137 PetscInt m = A->m; /* Number of rows in the matrix. */ 138 PetscInt rmax = a->rmax,*icols = crl->icols; 139 PetscScalar *acols = crl->acols; 140 PetscErrorCode ierr; 141 PetscScalar *x,*y; 142 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 143 PetscInt i,j,ii; 144 #endif 145 146 147 #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 148 #pragma disjoint(*x,*y,*aa) 149 #endif 150 151 PetscFunctionBegin; 152 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 153 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 154 155 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 156 fortranmultcrl_(&m,&rmax,x,y,icols,acols); 157 #else 158 159 /* first column */ 160 for (j=0; j<m; j++) { 161 y[j] = acols[j]*x[icols[j]]; 162 } 163 164 /* other columns */ 165 #if defined(PETSC_HAVE_CRAYC) 166 #pragma _CRI preferstream 167 #endif 168 for (i=1; i<rmax; i++) { 169 ii = i*m; 170 #if defined(PETSC_HAVE_CRAYC) 171 #pragma _CRI prefervector 172 #endif 173 for (j=0; j<m; j++) { 174 y[j] = y[j] + acols[ii+j]*x[icols[ii+j]]; 175 } 176 } 177 #if defined(PETSC_HAVE_CRAYC) 178 #pragma _CRI ivdep 179 #endif 180 181 #endif 182 ierr = PetscLogFlops(2*a->nz - m);CHKERRQ(ierr); 183 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 184 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 185 PetscFunctionReturn(0); 186 } 187 188 189 /* MatConvert_SeqAIJ_SeqCRL converts a SeqAIJ matrix into a 190 * SeqCRL matrix. This routine is called by the MatCreate_SeqCRL() 191 * routine, but can also be used to convert an assembled SeqAIJ matrix 192 * into a SeqCRL one. */ 193 EXTERN_C_BEGIN 194 #undef __FUNCT__ 195 #define __FUNCT__ "MatConvert_SeqAIJ_SeqCRL" 196 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCRL(Mat A,MatType type,MatReuse reuse,Mat *newmat) 197 { 198 /* This routine is only called to convert to MATSEQCRL 199 * from MATSEQAIJ, so we can ignore 'MatType Type'. */ 200 PetscErrorCode ierr; 201 Mat B = *newmat; 202 Mat_SeqCRL *crl; 203 204 PetscFunctionBegin; 205 if (reuse == MAT_INITIAL_MATRIX) { 206 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 207 } 208 209 ierr = PetscNew(Mat_SeqCRL,&crl);CHKERRQ(ierr); 210 B->spptr = (void *) crl; 211 212 /* Save a pointer to the original SeqAIJ assembly end routine, because we 213 * will want to use it later in the CRL assembly end routine. 214 * Also, save a pointer to the original SeqAIJ Destroy routine, because we 215 * will want to use it in the CRL destroy routine. */ 216 crl->AssemblyEnd_SeqAIJ = A->ops->assemblyend; 217 crl->MatDestroy_SeqAIJ = A->ops->destroy; 218 crl->MatDuplicate_SeqAIJ = A->ops->duplicate; 219 220 /* Set function pointers for methods that we inherit from AIJ but 221 * override. */ 222 B->ops->duplicate = MatDuplicate_SeqCRL; 223 B->ops->assemblyend = MatAssemblyEnd_SeqCRL; 224 B->ops->destroy = MatDestroy_SeqCRL; 225 B->ops->mult = MatMult_SeqCRL; 226 227 /* If A has already been assembled, compute the permutation. */ 228 if (A->assembled == PETSC_TRUE) { 229 ierr = SeqCRL_create_crl(B);CHKERRQ(ierr); 230 } 231 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQCRL);CHKERRQ(ierr); 232 *newmat = B; 233 PetscFunctionReturn(0); 234 } 235 EXTERN_C_END 236 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "MatCreateSeqCRL" 240 /*@C 241 MatCreateSeqCRL - Creates a sparse matrix of type SEQCRL. 242 This type inherits from AIJ, but stores some additional 243 information that is used to allow better vectorization of 244 the matrix-vector product. At the cost of increased storage, the AIJ formatted 245 matrix can be copied to a format in which pieces of the matrix are 246 stored in ELLPACK format, allowing the vectorized matrix multiply 247 routine to use stride-1 memory accesses. As with the AIJ type, it is 248 important to preallocate matrix storage in order to get good assembly 249 performance. 250 251 Collective on MPI_Comm 252 253 Input Parameters: 254 + comm - MPI communicator, set to PETSC_COMM_SELF 255 . m - number of rows 256 . n - number of columns 257 . nz - number of nonzeros per row (same for all rows) 258 - nnz - array containing the number of nonzeros in the various rows 259 (possibly different for each row) or PETSC_NULL 260 261 Output Parameter: 262 . A - the matrix 263 264 Notes: 265 If nnz is given then nz is ignored 266 267 Level: intermediate 268 269 .keywords: matrix, cray, sparse, parallel 270 271 .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues() 272 @*/ 273 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 274 { 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 ierr = MatCreate(comm,A);CHKERRQ(ierr); 279 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 280 ierr = MatSetType(*A,MATSEQCRL);CHKERRQ(ierr); 281 ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 286 EXTERN_C_BEGIN 287 #undef __FUNCT__ 288 #define __FUNCT__ "MatCreate_SeqCRL" 289 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqCRL(Mat A) 290 { 291 PetscErrorCode ierr; 292 293 PetscFunctionBegin; 294 /* Change the type name before calling MatSetType() to force proper construction of SeqAIJ 295 and MATSEQCRL types. */ 296 ierr = PetscObjectChangeTypeName((PetscObject)A,MATSEQCRL);CHKERRQ(ierr); 297 ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 298 ierr = MatConvert_SeqAIJ_SeqCRL(A,MATSEQCRL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 299 PetscFunctionReturn(0); 300 } 301 EXTERN_C_END 302 303