Actual source code: crl.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

  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: {
 17:   Mat_AIJCRL     *aijcrl = (Mat_AIJCRL*) A->spptr;

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

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

 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;

 47:   aijcrl->nz   = a->nz;
 48:   aijcrl->m    = A->rmap->n;
 49:   aijcrl->rmax = rmax;

 51:   PetscFree2(aijcrl->acols,aijcrl->icols);
 52:   PetscMalloc2(rmax*m,&aijcrl->acols,rmax*m,&aijcrl->icols);
 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:   PetscInfo2(A,"Percentage of 0's introduced for vectorized multiply %g. Rmax= %D\n",1.0-((double)a->nz)/((double)(rmax*m)),rmax);
 66:   return(0);
 67: }

 69: PetscErrorCode MatAssemblyEnd_SeqAIJCRL(Mat A, MatAssemblyType mode)
 70: {
 72:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

 75:   a->inode.use = PETSC_FALSE;

 77:   MatAssemblyEnd_SeqAIJ(A,mode);
 78:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);

 80:   /* Now calculate the permutation and grouping information. */
 81:   MatSeqAIJCRL_create_aijcrl(A);
 82:   return(0);
 83: }

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

 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

 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

105: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
106: #pragma disjoint(*x,*y,*aa)
107: #endif

110:   if (aijcrl->xscat) {
111:     VecCopy(xx,aijcrl->xwork);
112:     /* get remote values needed for local part of multiply */
113:     VecScatterBegin(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);
114:     VecScatterEnd(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);
115:     xx   = aijcrl->xwork;
116:   }

118:   VecGetArrayRead(xx,&x);
119:   VecGetArray(yy,&y);

121: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
122:   fortranmultcrl_(&m,&rmax,x,y,icols,acols);
123: #else

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

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

143: #endif
144:   PetscLogFlops(2.0*aijcrl->nz - m);
145:   VecRestoreArrayRead(xx,&x);
146:   VecRestoreArray(yy,&y);
147:   return(0);
148: }


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: {
158:   Mat            B = *newmat;
159:   Mat_AIJCRL     *aijcrl;
160:   PetscBool      sametype;

163:   if (reuse == MAT_INITIAL_MATRIX) {
164:     MatDuplicate(A,MAT_COPY_VALUES,&B);
165:   }
166:   PetscObjectTypeCompare((PetscObject)A,type,&sametype);
167:   if (sametype) return(0);

169:   PetscNewLog(B,&aijcrl);
170:   B->spptr = (void*) aijcrl;

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;

178:   /* If A has already been assembled, compute the permutation. */
179:   if (A->assembled) {
180:     MatSeqAIJCRL_create_aijcrl(B);
181:   }
182:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCRL);
183:   *newmat = B;
184:   return(0);
185: }

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.

198:    Collective

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

208:    Output Parameter:
209: .  A - the matrix

211:    Notes:
212:    If nnz is given then nz is ignored

214:    Level: intermediate

216: .seealso: MatCreate(), MatCreateMPIAIJPERM(), MatSetValues()
217: @*/
218: PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
219: {

223:   MatCreate(comm,A);
224:   MatSetSizes(*A,m,n,m,n);
225:   MatSetType(*A,MATSEQAIJCRL);
226:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
227:   return(0);
228: }

230: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCRL(Mat A)
231: {

235:   MatSetType(A,MATSEQAIJ);
236:   MatConvert_SeqAIJ_SeqAIJCRL(A,MATSEQAIJCRL,MAT_INPLACE_MATRIX,&A);
237:   return(0);
238: }