Actual source code: crl.c


  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.

  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>

 14: PetscErrorCode MatDestroy_SeqAIJCRL(Mat A)
 15: {
 16:   Mat_AIJCRL     *aijcrl = (Mat_AIJCRL*) A->spptr;

 18:   /* Free everything in the Mat_AIJCRL data structure. */
 19:   if (aijcrl) {
 20:     PetscFree2(aijcrl->acols,aijcrl->icols);
 21:   }
 22:   PetscFree(A->spptr);
 23:   PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);
 24:   MatDestroy_SeqAIJ(A);
 25:   return 0;
 26: }

 28: PetscErrorCode MatDuplicate_AIJCRL(Mat A, MatDuplicateOption op, Mat *M)
 29: {
 30:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot duplicate AIJCRL matrices yet");
 31: }

 33: PetscErrorCode MatSeqAIJCRL_create_aijcrl(Mat A)
 34: {
 35:   Mat_SeqAIJ     *a      = (Mat_SeqAIJ*)(A)->data;
 36:   Mat_AIJCRL     *aijcrl = (Mat_AIJCRL*) A->spptr;
 37:   PetscInt       m       = A->rmap->n; /* Number of rows in the matrix. */
 38:   PetscInt       *aj     = a->j; /* From the CSR representation; points to the beginning  of each row. */
 39:   PetscInt       i, j,rmax = a->rmax,*icols, *ilen = a->ilen;
 40:   MatScalar      *aa = a->a;
 41:   PetscScalar    *acols;

 43:   aijcrl->nz   = a->nz;
 44:   aijcrl->m    = A->rmap->n;
 45:   aijcrl->rmax = rmax;

 47:   PetscFree2(aijcrl->acols,aijcrl->icols);
 48:   PetscMalloc2(rmax*m,&aijcrl->acols,rmax*m,&aijcrl->icols);
 49:   acols = aijcrl->acols;
 50:   icols = aijcrl->icols;
 51:   for (i=0; i<m; i++) {
 52:     for (j=0; j<ilen[i]; j++) {
 53:       acols[j*m+i] = *aa++;
 54:       icols[j*m+i] = *aj++;
 55:     }
 56:     for (; j<rmax; j++) { /* empty column entries */
 57:       acols[j*m+i] = 0.0;
 58:       icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0;  /* handle case where row is EMPTY */
 59:     }
 60:   }
 61:   PetscInfo(A,"Percentage of 0's introduced for vectorized multiply %g. Rmax= %" PetscInt_FMT "\n",1.0-((double)a->nz)/((double)(rmax*m)),rmax);
 62:   return 0;
 63: }

 65: PetscErrorCode MatAssemblyEnd_SeqAIJCRL(Mat A, MatAssemblyType mode)
 66: {
 67:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

 69:   a->inode.use = PETSC_FALSE;

 71:   MatAssemblyEnd_SeqAIJ(A,mode);
 72:   if (mode == MAT_FLUSH_ASSEMBLY) return 0;

 74:   /* Now calculate the permutation and grouping information. */
 75:   MatSeqAIJCRL_create_aijcrl(A);
 76:   return 0;
 77: }

 79: #include <../src/mat/impls/aij/seq/crl/ftn-kernels/fmultcrl.h>

 81: /*
 82:     Shared by both sequential and parallel versions of CRL matrix: MATMPIAIJCRL and MATSEQAIJCRL
 83:     - the scatter is used only in the parallel version

 85: */
 86: PetscErrorCode MatMult_AIJCRL(Mat A,Vec xx,Vec yy)
 87: {
 88:   Mat_AIJCRL        *aijcrl = (Mat_AIJCRL*) A->spptr;
 89:   PetscInt          m       = aijcrl->m; /* Number of rows in the matrix. */
 90:   PetscInt          rmax    = aijcrl->rmax,*icols = aijcrl->icols;
 91:   PetscScalar       *acols  = aijcrl->acols;
 92:   PetscScalar       *y;
 93:   const PetscScalar *x;
 94: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
 95:   PetscInt          i,j,ii;
 96: #endif

 98: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
 99: #pragma disjoint(*x,*y,*aa)
100: #endif

102:   if (aijcrl->xscat) {
103:     VecCopy(xx,aijcrl->xwork);
104:     /* get remote values needed for local part of multiply */
105:     VecScatterBegin(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);
106:     VecScatterEnd(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);
107:     xx   = aijcrl->xwork;
108:   }

110:   VecGetArrayRead(xx,&x);
111:   VecGetArray(yy,&y);

113: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
114:   fortranmultcrl_(&m,&rmax,x,y,icols,acols);
115: #else

117:   /* first column */
118:   for (j=0; j<m; j++) y[j] = acols[j]*x[icols[j]];

120:   /* other columns */
121: #if defined(PETSC_HAVE_CRAY_VECTOR)
122: #pragma _CRI preferstream
123: #endif
124:   for (i=1; i<rmax; i++) {
125:     ii = i*m;
126: #if defined(PETSC_HAVE_CRAY_VECTOR)
127: #pragma _CRI prefervector
128: #endif
129:     for (j=0; j<m; j++) y[j] = y[j] + acols[ii+j]*x[icols[ii+j]];
130:   }
131: #if defined(PETSC_HAVE_CRAY_VECTOR)
132: #pragma _CRI ivdep
133: #endif

135: #endif
136:   PetscLogFlops(2.0*aijcrl->nz - m);
137:   VecRestoreArrayRead(xx,&x);
138:   VecRestoreArray(yy,&y);
139:   return 0;
140: }

142: /* MatConvert_SeqAIJ_SeqAIJCRL converts a SeqAIJ matrix into a
143:  * SeqAIJCRL matrix.  This routine is called by the MatCreate_SeqAIJCRL()
144:  * routine, but can also be used to convert an assembled SeqAIJ matrix
145:  * into a SeqAIJCRL one. */
146: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
147: {
148:   Mat            B = *newmat;
149:   Mat_AIJCRL     *aijcrl;
150:   PetscBool      sametype;

152:   if (reuse == MAT_INITIAL_MATRIX) {
153:     MatDuplicate(A,MAT_COPY_VALUES,&B);
154:   }
155:   PetscObjectTypeCompare((PetscObject)A,type,&sametype);
156:   if (sametype) return 0;

158:   PetscNewLog(B,&aijcrl);
159:   B->spptr = (void*) aijcrl;

161:   /* Set function pointers for methods that we inherit from AIJ but override. */
162:   B->ops->duplicate   = MatDuplicate_AIJCRL;
163:   B->ops->assemblyend = MatAssemblyEnd_SeqAIJCRL;
164:   B->ops->destroy     = MatDestroy_SeqAIJCRL;
165:   B->ops->mult        = MatMult_AIJCRL;

167:   /* If A has already been assembled, compute the permutation. */
168:   if (A->assembled) {
169:     MatSeqAIJCRL_create_aijcrl(B);
170:   }
171:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCRL);
172:   *newmat = B;
173:   return 0;
174: }

176: /*@C
177:    MatCreateSeqAIJCRL - Creates a sparse matrix of type SEQAIJCRL.
178:    This type inherits from AIJ, but stores some additional
179:    information that is used to allow better vectorization of
180:    the matrix-vector product. At the cost of increased storage, the AIJ formatted
181:    matrix can be copied to a format in which pieces of the matrix are
182:    stored in ELLPACK format, allowing the vectorized matrix multiply
183:    routine to use stride-1 memory accesses.  As with the AIJ type, it is
184:    important to preallocate matrix storage in order to get good assembly
185:    performance.

187:    Collective

189:    Input Parameters:
190: +  comm - MPI communicator, set to PETSC_COMM_SELF
191: .  m - number of rows
192: .  n - number of columns
193: .  nz - number of nonzeros per row (same for all rows)
194: -  nnz - array containing the number of nonzeros in the various rows
195:          (possibly different for each row) or NULL

197:    Output Parameter:
198: .  A - the matrix

200:    Notes:
201:    If nnz is given then nz is ignored

203:    Level: intermediate

205: .seealso: MatCreate(), MatCreateMPIAIJPERM(), MatSetValues()
206: @*/
207: PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
208: {
209:   MatCreate(comm,A);
210:   MatSetSizes(*A,m,n,m,n);
211:   MatSetType(*A,MATSEQAIJCRL);
212:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
213:   return 0;
214: }

216: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCRL(Mat A)
217: {
218:   MatSetType(A,MATSEQAIJ);
219:   MatConvert_SeqAIJ_SeqAIJCRL(A,MATSEQAIJCRL,MAT_INPLACE_MATRIX,&A);
220:   return 0;
221: }