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 } 34 35 PetscErrorCode MatSeqAIJCRL_create_aijcrl(Mat A) 36 { 37 Mat_SeqAIJ *a = (Mat_SeqAIJ*)(A)->data; 38 Mat_AIJCRL *aijcrl = (Mat_AIJCRL*) A->spptr; 39 PetscInt m = A->rmap->n; /* Number of rows in the matrix. */ 40 PetscInt *aj = a->j; /* From the CSR representation; points to the beginning of each row. */ 41 PetscInt i, j,rmax = a->rmax,*icols, *ilen = a->ilen; 42 MatScalar *aa = a->a; 43 PetscScalar *acols; 44 PetscErrorCode ierr; 45 46 PetscFunctionBegin; 47 aijcrl->nz = a->nz; 48 aijcrl->m = A->rmap->n; 49 aijcrl->rmax = rmax; 50 51 ierr = PetscFree2(aijcrl->acols,aijcrl->icols);CHKERRQ(ierr); 52 ierr = PetscMalloc2(rmax*m,&aijcrl->acols,rmax*m,&aijcrl->icols);CHKERRQ(ierr); 53 acols = aijcrl->acols; 54 icols = aijcrl->icols; 55 for (i=0; i<m; i++) { 56 for (j=0; j<ilen[i]; j++) { 57 acols[j*m+i] = *aa++; 58 icols[j*m+i] = *aj++; 59 } 60 for (; j<rmax; j++) { /* empty column entries */ 61 acols[j*m+i] = 0.0; 62 icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0; /* handle case where row is EMPTY */ 63 } 64 } 65 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); 66 PetscFunctionReturn(0); 67 } 68 69 PetscErrorCode MatAssemblyEnd_SeqAIJCRL(Mat A, MatAssemblyType mode) 70 { 71 PetscErrorCode ierr; 72 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 73 74 PetscFunctionBegin; 75 a->inode.use = PETSC_FALSE; 76 77 ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr); 78 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 79 80 /* Now calculate the permutation and grouping information. */ 81 ierr = MatSeqAIJCRL_create_aijcrl(A);CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 #include <../src/mat/impls/aij/seq/crl/ftn-kernels/fmultcrl.h> 86 87 /* 88 Shared by both sequential and parallel versions of CRL matrix: MATMPIAIJCRL and MATSEQAIJCRL 89 - the scatter is used only in the parallel version 90 91 */ 92 PetscErrorCode MatMult_AIJCRL(Mat A,Vec xx,Vec yy) 93 { 94 Mat_AIJCRL *aijcrl = (Mat_AIJCRL*) A->spptr; 95 PetscInt m = aijcrl->m; /* Number of rows in the matrix. */ 96 PetscInt rmax = aijcrl->rmax,*icols = aijcrl->icols; 97 PetscScalar *acols = aijcrl->acols; 98 PetscErrorCode ierr; 99 PetscScalar *y; 100 const PetscScalar *x; 101 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 102 PetscInt i,j,ii; 103 #endif 104 105 #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 106 #pragma disjoint(*x,*y,*aa) 107 #endif 108 109 PetscFunctionBegin; 110 if (aijcrl->xscat) { 111 ierr = VecCopy(xx,aijcrl->xwork);CHKERRQ(ierr); 112 /* get remote values needed for local part of multiply */ 113 ierr = VecScatterBegin(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 114 ierr = VecScatterEnd(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 115 xx = aijcrl->xwork; 116 } 117 118 ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 119 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 120 121 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 122 fortranmultcrl_(&m,&rmax,x,y,icols,acols); 123 #else 124 125 /* first column */ 126 for (j=0; j<m; j++) y[j] = acols[j]*x[icols[j]]; 127 128 /* other columns */ 129 #if defined(PETSC_HAVE_CRAY_VECTOR) 130 #pragma _CRI preferstream 131 #endif 132 for (i=1; i<rmax; i++) { 133 ii = i*m; 134 #if defined(PETSC_HAVE_CRAY_VECTOR) 135 #pragma _CRI prefervector 136 #endif 137 for (j=0; j<m; j++) y[j] = y[j] + acols[ii+j]*x[icols[ii+j]]; 138 } 139 #if defined(PETSC_HAVE_CRAY_VECTOR) 140 #pragma _CRI ivdep 141 #endif 142 143 #endif 144 ierr = PetscLogFlops(2.0*aijcrl->nz - m);CHKERRQ(ierr); 145 ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 146 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 147 PetscFunctionReturn(0); 148 } 149 150 151 /* MatConvert_SeqAIJ_SeqAIJCRL converts a SeqAIJ matrix into a 152 * SeqAIJCRL matrix. This routine is called by the MatCreate_SeqAIJCRL() 153 * routine, but can also be used to convert an assembled SeqAIJ matrix 154 * into a SeqAIJCRL one. */ 155 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat A,MatType type,MatReuse reuse,Mat *newmat) 156 { 157 PetscErrorCode ierr; 158 Mat B = *newmat; 159 Mat_AIJCRL *aijcrl; 160 PetscBool sametype; 161 162 PetscFunctionBegin; 163 if (reuse == MAT_INITIAL_MATRIX) { 164 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 165 } 166 ierr = PetscObjectTypeCompare((PetscObject)A,type,&sametype);CHKERRQ(ierr); 167 if (sametype) PetscFunctionReturn(0); 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 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 .seealso: MatCreate(), MatCreateMPIAIJPERM(), MatSetValues() 217 @*/ 218 PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 219 { 220 PetscErrorCode ierr; 221 222 PetscFunctionBegin; 223 ierr = MatCreate(comm,A);CHKERRQ(ierr); 224 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 225 ierr = MatSetType(*A,MATSEQAIJCRL);CHKERRQ(ierr); 226 ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);CHKERRQ(ierr); 227 PetscFunctionReturn(0); 228 } 229 230 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCRL(Mat A) 231 { 232 PetscErrorCode ierr; 233 234 PetscFunctionBegin; 235 ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 236 ierr = MatConvert_SeqAIJ_SeqAIJCRL(A,MATSEQAIJCRL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 237 PetscFunctionReturn(0); 238 } 239 240 241