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