Actual source code: crl.c

  1: #define PETSCMAT_DLL

  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.

 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

 17: PetscErrorCode MatDestroy_SeqCRL(Mat A)
 18: {
 20:   Mat_CRL        *crl = (Mat_CRL *) A->spptr;

 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;

 31:   /* Free everything in the Mat_CRL data structure. */
 32:   PetscFree2(crl->acols,crl->icols);

 34:   /* Change the type of A back to SEQAIJ and use MatDestroy() 
 35:    * to destroy everything that remains. */
 36:   PetscObjectChangeTypeName( (PetscObject)A, MATSEQAIJ);
 37:   /* Note that I don't call MatSetType().  I believe this is because that 
 38:    * is only to be called when *building* a matrix. */
 39:   (*A->ops->destroy)(A);
 40:   return(0);
 41: }

 43: PetscErrorCode MatDuplicate_CRL(Mat A, MatDuplicateOption op, Mat *M)
 44: {
 46:   Mat_CRL        *crl = (Mat_CRL *) A->spptr;

 49:   (*crl->MatDuplicate)(A,op,M);
 50:   SETERRQ(PETSC_ERR_SUP,"Cannot duplicate CRL matrices yet");
 51:   return(0);
 52: }

 56: PetscErrorCode SeqCRL_create_crl(Mat A)
 57: {
 58:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)(A)->data;
 59:   Mat_CRL        *crl = (Mat_CRL*) A->spptr;
 60:   PetscInt       m = A->rmap->n;  /* Number of rows in the matrix. */
 61:   PetscInt       *aj = a->j;  /* From the CSR representation; points to the beginning  of each row. */
 62:   PetscInt       i, j,rmax = a->rmax,*icols, *ilen = a->ilen;
 63:   MatScalar      *aa = a->a;
 64:   PetscScalar    *acols;

 68:   crl->nz   = a->nz;
 69:   crl->m    = A->rmap->n;
 70:   crl->rmax = rmax;
 71:   PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);
 72:   acols = crl->acols;
 73:   icols = crl->icols;
 74:   for (i=0; i<m; i++) {
 75:     for (j=0; j<ilen[i]; j++) {
 76:       acols[j*m+i] = *aa++;
 77:       icols[j*m+i] = *aj++;
 78:     }
 79:     for (;j<rmax; j++) { /* empty column entries */
 80:       acols[j*m+i] = 0.0;
 81:       icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0;  /* handle case where row is EMPTY */
 82:     }
 83:   }
 84:   PetscInfo2(A,"Percentage of 0's introduced for vectorized multiply %G. Rmax= %D\n",1.0-((double)a->nz)/((double)(rmax*m)),rmax);
 85:   return(0);
 86: }

 90: PetscErrorCode MatAssemblyEnd_SeqCRL(Mat A, MatAssemblyType mode)
 91: {
 93:   Mat_CRL        *crl = (Mat_CRL*) A->spptr;
 94:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

 97:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);
 98: 
 99:   /* Since a MATSEQCRL matrix is really just a MATSEQAIJ with some 
100:    * extra information, call the AssemblyEnd routine for a MATSEQAIJ. 
101:    * I'm not sure if this is the best way to do this, but it avoids 
102:    * a lot of code duplication.
103:    * I also note that currently MATSEQCRL doesn't know anything about 
104:    * the Mat_CompressedRow data structure that SeqAIJ now uses when there 
105:    * are many zero rows.  If the SeqAIJ assembly end routine decides to use 
106:    * this, this may break things.  (Don't know... haven't looked at it.) */
107:   a->inode.use = PETSC_FALSE;
108:   (*crl->AssemblyEnd)(A, mode);

110:   /* Now calculate the permutation and grouping information. */
111:   SeqCRL_create_crl(A);
112:   return(0);
113: }

115:  #include ../src/inline/dot.h

119: PetscErrorCode MatMult_CRL(Mat A,Vec xx,Vec yy)
120: {
121:   Mat_CRL        *crl = (Mat_CRL*) A->spptr;
122:   PetscInt       m = crl->m;  /* Number of rows in the matrix. */
123:   PetscInt       rmax = crl->rmax,*icols = crl->icols;
124:   PetscScalar    *acols = crl->acols;
126:   PetscScalar    *x,*y;
127: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
128:   PetscInt       i,j,ii;
129: #endif


132: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
133: #pragma disjoint(*x,*y,*aa)
134: #endif

137:   if (crl->xscat) {
138:     VecCopy(xx,crl->xwork);
139:     /* get remote values needed for local part of multiply */
140:     VecScatterBegin(crl->xscat,xx,crl->fwork,INSERT_VALUES,SCATTER_FORWARD);
141:     VecScatterEnd(crl->xscat,xx,crl->fwork,INSERT_VALUES,SCATTER_FORWARD);
142:     xx = crl->xwork;
143:   };

145:   VecGetArray(xx,&x);
146:   VecGetArray(yy,&y);

148: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
149:   fortranmultcrl_(&m,&rmax,x,y,icols,acols);
150: #else

152:   /* first column */
153:   for (j=0; j<m; j++) {
154:     y[j] = acols[j]*x[icols[j]];
155:   }

157:   /* other columns */
158: #if defined(PETSC_HAVE_CRAYC)
159: #pragma _CRI preferstream
160: #endif
161:   for (i=1; i<rmax; i++) {
162:     ii = i*m;
163: #if defined(PETSC_HAVE_CRAYC)
164: #pragma _CRI prefervector
165: #endif
166:     for (j=0; j<m; j++) {
167:       y[j] = y[j] + acols[ii+j]*x[icols[ii+j]];
168:     }
169:   }
170: #if defined(PETSC_HAVE_CRAYC)
171: #pragma _CRI ivdep
172: #endif

174: #endif
175:   PetscLogFlops(2*crl->nz - m);
176:   VecRestoreArray(xx,&x);
177:   VecRestoreArray(yy,&y);
178:   return(0);
179: }


182: /* MatConvert_SeqAIJ_SeqCRL converts a SeqAIJ matrix into a 
183:  * SeqCRL matrix.  This routine is called by the MatCreate_SeqCRL() 
184:  * routine, but can also be used to convert an assembled SeqAIJ matrix 
185:  * into a SeqCRL one. */
189: PetscErrorCode  MatConvert_SeqAIJ_SeqCRL(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
190: {
192:   Mat            B = *newmat;
193:   Mat_CRL        *crl;

196:   if (reuse == MAT_INITIAL_MATRIX) {
197:     MatDuplicate(A,MAT_COPY_VALUES,&B);
198:   }

200:   PetscNewLog(B,Mat_CRL,&crl);
201:   B->spptr = (void *) crl;

203:   crl->AssemblyEnd  = A->ops->assemblyend;
204:   crl->MatDestroy   = A->ops->destroy;
205:   crl->MatDuplicate = A->ops->duplicate;

207:   /* Set function pointers for methods that we inherit from AIJ but override. */
208:   B->ops->duplicate   = MatDuplicate_CRL;
209:   B->ops->assemblyend = MatAssemblyEnd_SeqCRL;
210:   B->ops->destroy     = MatDestroy_SeqCRL;
211:   B->ops->mult        = MatMult_CRL;

213:   /* If A has already been assembled, compute the permutation. */
214:   if (A->assembled == PETSC_TRUE) {
215:     SeqCRL_create_crl(B);
216:   }
217:   PetscObjectChangeTypeName((PetscObject)B,MATSEQCRL);
218:   *newmat = B;
219:   return(0);
220: }


226: /*@C
227:    MatCreateSeqCRL - Creates a sparse matrix of type SEQCRL.
228:    This type inherits from AIJ, but stores some additional
229:    information that is used to allow better vectorization of 
230:    the matrix-vector product. At the cost of increased storage, the AIJ formatted 
231:    matrix can be copied to a format in which pieces of the matrix are 
232:    stored in ELLPACK format, allowing the vectorized matrix multiply 
233:    routine to use stride-1 memory accesses.  As with the AIJ type, it is 
234:    important to preallocate matrix storage in order to get good assembly 
235:    performance.
236:    
237:    Collective on MPI_Comm

239:    Input Parameters:
240: +  comm - MPI communicator, set to PETSC_COMM_SELF
241: .  m - number of rows
242: .  n - number of columns
243: .  nz - number of nonzeros per row (same for all rows)
244: -  nnz - array containing the number of nonzeros in the various rows 
245:          (possibly different for each row) or PETSC_NULL

247:    Output Parameter:
248: .  A - the matrix 

250:    Notes:
251:    If nnz is given then nz is ignored

253:    Level: intermediate

255: .keywords: matrix, cray, sparse, parallel

257: .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues()
258: @*/
259: PetscErrorCode  MatCreateSeqCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
260: {

264:   MatCreate(comm,A);
265:   MatSetSizes(*A,m,n,m,n);
266:   MatSetType(*A,MATSEQCRL);
267:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
268:   return(0);
269: }


275: PetscErrorCode  MatCreate_SeqCRL(Mat A)
276: {

280:   MatSetType(A,MATSEQAIJ);
281:   MatConvert_SeqAIJ_SeqCRL(A,MATSEQCRL,MAT_REUSE_MATRIX,&A);
282:   return(0);
283: }