1 2 /* 3 Defines a matrix-vector product for the MATSEQAIJCRL matrix class. 4 This class is derived from the MATSEQAIJ class and retains the 5 compressed row storage (aka Yale sparse matrix format) but augments 6 it with a column oriented storage that is more efficient for 7 matrix vector products on Vector machines. 8 9 CRL stands for constant row length (that is the same number of columns 10 is kept (padded with zeros) for each row of the sparse matrix. 11 */ 12 #include <../src/mat/impls/aij/seq/crl/crl.h> 13 14 PetscErrorCode MatDestroy_SeqAIJCRL(Mat A) 15 { 16 PetscErrorCode ierr; 17 Mat_AIJCRL *aijcrl = (Mat_AIJCRL*) A->spptr; 18 19 /* Free everything in the Mat_AIJCRL data structure. */ 20 if (aijcrl) { 21 ierr = PetscFree2(aijcrl->acols,aijcrl->icols);CHKERRQ(ierr); 22 } 23 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 24 ierr = PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);CHKERRQ(ierr); 25 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 26 PetscFunctionReturn(0); 27 } 28 29 PetscErrorCode MatDuplicate_AIJCRL(Mat A, MatDuplicateOption op, Mat *M) 30 { 31 PetscFunctionBegin; 32 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot duplicate AIJCRL matrices yet"); 33 PetscFunctionReturn(0); 34 } 35 36 PetscErrorCode MatSeqAIJCRL_create_aijcrl(Mat A) 37 { 38 Mat_SeqAIJ *a = (Mat_SeqAIJ*)(A)->data; 39 Mat_AIJCRL *aijcrl = (Mat_AIJCRL*) A->spptr; 40 PetscInt m = A->rmap->n; /* Number of rows in the matrix. */ 41 PetscInt *aj = a->j; /* From the CSR representation; points to the beginning of each row. */ 42 PetscInt i, j,rmax = a->rmax,*icols, *ilen = a->ilen; 43 MatScalar *aa = a->a; 44 PetscScalar *acols; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 aijcrl->nz = a->nz; 49 aijcrl->m = A->rmap->n; 50 aijcrl->rmax = rmax; 51 52 ierr = PetscFree2(aijcrl->acols,aijcrl->icols);CHKERRQ(ierr); 53 ierr = PetscMalloc2(rmax*m,&aijcrl->acols,rmax*m,&aijcrl->icols);CHKERRQ(ierr); 54 acols = aijcrl->acols; 55 icols = aijcrl->icols; 56 for (i=0; i<m; i++) { 57 for (j=0; j<ilen[i]; j++) { 58 acols[j*m+i] = *aa++; 59 icols[j*m+i] = *aj++; 60 } 61 for (; j<rmax; j++) { /* empty column entries */ 62 acols[j*m+i] = 0.0; 63 icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0; /* handle case where row is EMPTY */ 64 } 65 } 66 ierr = PetscInfo2(A,"Percentage of 0's introduced for vectorized multiply %g. Rmax= %D\n",1.0-((double)a->nz)/((double)(rmax*m)),rmax);CHKERRQ(ierr); 67 PetscFunctionReturn(0); 68 } 69 70 extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 71 72 PetscErrorCode MatAssemblyEnd_SeqAIJCRL(Mat A, MatAssemblyType mode) 73 { 74 PetscErrorCode ierr; 75 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 76 77 PetscFunctionBegin; 78 a->inode.use = PETSC_FALSE; 79 80 ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr); 81 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 82 83 /* Now calculate the permutation and grouping information. */ 84 ierr = MatSeqAIJCRL_create_aijcrl(A);CHKERRQ(ierr); 85 PetscFunctionReturn(0); 86 } 87 88 #include <../src/mat/impls/aij/seq/crl/ftn-kernels/fmultcrl.h> 89 90 /* 91 Shared by both sequential and parallel versions of CRL matrix: MATMPIAIJCRL and MATSEQAIJCRL 92 - the scatter is used only in the parallel version 93 94 */ 95 PetscErrorCode MatMult_AIJCRL(Mat A,Vec xx,Vec yy) 96 { 97 Mat_AIJCRL *aijcrl = (Mat_AIJCRL*) A->spptr; 98 PetscInt m = aijcrl->m; /* Number of rows in the matrix. */ 99 PetscInt rmax = aijcrl->rmax,*icols = aijcrl->icols; 100 PetscScalar *acols = aijcrl->acols; 101 PetscErrorCode ierr; 102 PetscScalar *y; 103 const PetscScalar *x; 104 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 105 PetscInt i,j,ii; 106 #endif 107 108 #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 109 #pragma disjoint(*x,*y,*aa) 110 #endif 111 112 PetscFunctionBegin; 113 if (aijcrl->xscat) { 114 ierr = VecCopy(xx,aijcrl->xwork);CHKERRQ(ierr); 115 /* get remote values needed for local part of multiply */ 116 ierr = VecScatterBegin(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 117 ierr = VecScatterEnd(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 118 xx = aijcrl->xwork; 119 } 120 121 ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 122 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 123 124 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 125 fortranmultcrl_(&m,&rmax,x,y,icols,acols); 126 #else 127 128 /* first column */ 129 for (j=0; j<m; j++) y[j] = acols[j]*x[icols[j]]; 130 131 /* other columns */ 132 #if defined(PETSC_HAVE_CRAY_VECTOR) 133 #pragma _CRI preferstream 134 #endif 135 for (i=1; i<rmax; i++) { 136 ii = i*m; 137 #if defined(PETSC_HAVE_CRAY_VECTOR) 138 #pragma _CRI prefervector 139 #endif 140 for (j=0; j<m; j++) y[j] = y[j] + acols[ii+j]*x[icols[ii+j]]; 141 } 142 #if defined(PETSC_HAVE_CRAY_VECTOR) 143 #pragma _CRI ivdep 144 #endif 145 146 #endif 147 ierr = PetscLogFlops(2.0*aijcrl->nz - m);CHKERRQ(ierr); 148 ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 149 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 150 PetscFunctionReturn(0); 151 } 152 153 154 /* MatConvert_SeqAIJ_SeqAIJCRL converts a SeqAIJ matrix into a 155 * SeqAIJCRL matrix. This routine is called by the MatCreate_SeqAIJCRL() 156 * routine, but can also be used to convert an assembled SeqAIJ matrix 157 * into a SeqAIJCRL one. */ 158 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat A,MatType type,MatReuse reuse,Mat *newmat) 159 { 160 PetscErrorCode ierr; 161 Mat B = *newmat; 162 Mat_AIJCRL *aijcrl; 163 164 PetscFunctionBegin; 165 if (reuse == MAT_INITIAL_MATRIX) { 166 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 167 } 168 169 ierr = PetscNewLog(B,&aijcrl);CHKERRQ(ierr); 170 B->spptr = (void*) aijcrl; 171 172 /* Set function pointers for methods that we inherit from AIJ but override. */ 173 B->ops->duplicate = MatDuplicate_AIJCRL; 174 B->ops->assemblyend = MatAssemblyEnd_SeqAIJCRL; 175 B->ops->destroy = MatDestroy_SeqAIJCRL; 176 B->ops->mult = MatMult_AIJCRL; 177 178 /* If A has already been assembled, compute the permutation. */ 179 if (A->assembled) { 180 ierr = MatSeqAIJCRL_create_aijcrl(B);CHKERRQ(ierr); 181 } 182 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCRL);CHKERRQ(ierr); 183 *newmat = B; 184 PetscFunctionReturn(0); 185 } 186 187 /*@C 188 MatCreateSeqAIJCRL - Creates a sparse matrix of type SEQAIJCRL. 189 This type inherits from AIJ, but stores some additional 190 information that is used to allow better vectorization of 191 the matrix-vector product. At the cost of increased storage, the AIJ formatted 192 matrix can be copied to a format in which pieces of the matrix are 193 stored in ELLPACK format, allowing the vectorized matrix multiply 194 routine to use stride-1 memory accesses. As with the AIJ type, it is 195 important to preallocate matrix storage in order to get good assembly 196 performance. 197 198 Collective on MPI_Comm 199 200 Input Parameters: 201 + comm - MPI communicator, set to PETSC_COMM_SELF 202 . m - number of rows 203 . n - number of columns 204 . nz - number of nonzeros per row (same for all rows) 205 - nnz - array containing the number of nonzeros in the various rows 206 (possibly different for each row) or NULL 207 208 Output Parameter: 209 . A - the matrix 210 211 Notes: 212 If nnz is given then nz is ignored 213 214 Level: intermediate 215 216 .keywords: matrix, cray, sparse, parallel 217 218 .seealso: MatCreate(), MatCreateMPIAIJPERM(), MatSetValues() 219 @*/ 220 PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 221 { 222 PetscErrorCode ierr; 223 224 PetscFunctionBegin; 225 ierr = MatCreate(comm,A);CHKERRQ(ierr); 226 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 227 ierr = MatSetType(*A,MATSEQAIJCRL);CHKERRQ(ierr); 228 ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);CHKERRQ(ierr); 229 PetscFunctionReturn(0); 230 } 231 232 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCRL(Mat A) 233 { 234 PetscErrorCode ierr; 235 236 PetscFunctionBegin; 237 ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 238 ierr = MatConvert_SeqAIJ_SeqAIJCRL(A,MATSEQAIJCRL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 239 PetscFunctionReturn(0); 240 } 241 242 243