Actual source code: mcrl.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: /*
  3:   Defines a matrix-vector product for the MATMPIAIJCRL matrix class.
  4:   This class is derived from the MATMPIAIJ 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.

 12:    See src/mat/impls/aij/seq/crl/crl.c for the sequential version
 13: */

 15: #include <../src/mat/impls/aij/mpi/mpiaij.h>
 16: #include <../src/mat/impls/aij/seq/crl/crl.h>

 18: extern PetscErrorCode MatDestroy_MPIAIJ(Mat);

 22: PetscErrorCode MatDestroy_MPIAIJCRL(Mat A)
 23: {
 25:   Mat_AIJCRL     *aijcrl = (Mat_AIJCRL*) A->spptr;

 27:   /* Free everything in the Mat_AIJCRL data structure. */
 28:   if (aijcrl) {
 29:     PetscFree2(aijcrl->acols,aijcrl->icols);
 30:     VecDestroy(&aijcrl->fwork);
 31:     VecDestroy(&aijcrl->xwork);
 32:     PetscFree(aijcrl->array);
 33:   }
 34:   PetscFree(A->spptr);

 36:   PetscObjectChangeTypeName((PetscObject)A, MATMPIAIJ);
 37:   MatDestroy_MPIAIJ(A);
 38:   return(0);
 39: }

 43: PetscErrorCode MatMPIAIJCRL_create_aijcrl(Mat A)
 44: {
 45:   Mat_MPIAIJ     *a      = (Mat_MPIAIJ*)(A)->data;
 46:   Mat_SeqAIJ     *Aij    = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->B->data);
 47:   Mat_AIJCRL     *aijcrl = (Mat_AIJCRL*) A->spptr;
 48:   PetscInt       m       = A->rmap->n; /* Number of rows in the matrix. */
 49:   PetscInt       nd      = a->A->cmap->n; /* number of columns in diagonal portion */
 50:   PetscInt       *aj     = Aij->j,*bj = Bij->j; /* From the CSR representation; points to the beginning  of each row. */
 51:   PetscInt       i, j,rmax = 0,*icols, *ailen = Aij->ilen, *bilen = Bij->ilen;
 52:   PetscScalar    *aa = Aij->a,*ba = Bij->a,*acols,*array;

 56:   /* determine the row with the most columns */
 57:   for (i=0; i<m; i++) {
 58:     rmax = PetscMax(rmax,ailen[i]+bilen[i]);
 59:   }
 60:   aijcrl->nz   = Aij->nz+Bij->nz;
 61:   aijcrl->m    = A->rmap->n;
 62:   aijcrl->rmax = rmax;

 64:   PetscFree2(aijcrl->acols,aijcrl->icols);
 65:   PetscMalloc2(rmax*m,&aijcrl->acols,rmax*m,&aijcrl->icols);
 66:   acols = aijcrl->acols;
 67:   icols = aijcrl->icols;
 68:   for (i=0; i<m; i++) {
 69:     for (j=0; j<ailen[i]; j++) {
 70:       acols[j*m+i] = *aa++;
 71:       icols[j*m+i] = *aj++;
 72:     }
 73:     for (; j<ailen[i]+bilen[i]; j++) {
 74:       acols[j*m+i] = *ba++;
 75:       icols[j*m+i] = nd + *bj++;
 76:     }
 77:     for (; j<rmax; j++) { /* empty column entries */
 78:       acols[j*m+i] = 0.0;
 79:       icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0;  /* handle case where row is EMPTY */
 80:     }
 81:   }
 82:   PetscInfo1(A,"Percentage of 0's introduced for vectorized multiply %g\n",1.0-((double)(aijcrl->nz))/((double)(rmax*m)));

 84:   PetscFree(aijcrl->array);
 85:   PetscMalloc1(a->B->cmap->n+nd,&array);
 86:   /* xwork array is actually B->n+nd long, but we define xwork this length so can copy into it */
 87:   VecDestroy(&aijcrl->xwork);
 88:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,nd,PETSC_DECIDE,array,&aijcrl->xwork);
 89:   VecDestroy(&aijcrl->fwork);
 90:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->B->cmap->n,array+nd,&aijcrl->fwork);

 92:   aijcrl->array = array;
 93:   aijcrl->xscat = a->Mvctx;
 94:   return(0);
 95: }

 97: extern PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat,MatAssemblyType);

101: PetscErrorCode MatAssemblyEnd_MPIAIJCRL(Mat A, MatAssemblyType mode)
102: {
104:   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data;
105:   Mat_SeqAIJ     *Aij = (Mat_SeqAIJ*)(a->A->data), *Bij = (Mat_SeqAIJ*)(a->A->data);

108:   Aij->inode.use = PETSC_FALSE;
109:   Bij->inode.use = PETSC_FALSE;

111:   MatAssemblyEnd_MPIAIJ(A,mode);
112:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);

114:   /* Now calculate the permutation and grouping information. */
115:   MatMPIAIJCRL_create_aijcrl(A);
116:   return(0);
117: }

119: extern PetscErrorCode MatMult_AIJCRL(Mat,Vec,Vec);
120: extern PetscErrorCode MatDuplicate_AIJCRL(Mat,MatDuplicateOption,Mat*);

122: /* MatConvert_MPIAIJ_MPIAIJCRL converts a MPIAIJ matrix into a
123:  * MPIAIJCRL matrix.  This routine is called by the MatCreate_MPIAIJCRL()
124:  * routine, but can also be used to convert an assembled MPIAIJ matrix
125:  * into a MPIAIJCRL one. */

129: PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPIAIJCRL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
130: {
132:   Mat            B = *newmat;
133:   Mat_AIJCRL     *aijcrl;

136:   if (reuse == MAT_INITIAL_MATRIX) {
137:     MatDuplicate(A,MAT_COPY_VALUES,&B);
138:   }

140:   PetscNewLog(B,&aijcrl);
141:   B->spptr = (void*) aijcrl;

143:   /* Set function pointers for methods that we inherit from AIJ but override. */
144:   B->ops->duplicate   = MatDuplicate_AIJCRL;
145:   B->ops->assemblyend = MatAssemblyEnd_MPIAIJCRL;
146:   B->ops->destroy     = MatDestroy_MPIAIJCRL;
147:   B->ops->mult        = MatMult_AIJCRL;

149:   /* If A has already been assembled, compute the permutation. */
150:   if (A->assembled) {
151:     MatMPIAIJCRL_create_aijcrl(B);
152:   }
153:   PetscObjectChangeTypeName((PetscObject)B,MATMPIAIJCRL);
154:   *newmat = B;
155:   return(0);
156: }

160: /*@C
161:    MatCreateMPIAIJCRL - Creates a sparse matrix of type MPIAIJCRL.
162:    This type inherits from AIJ, but stores some additional
163:    information that is used to allow better vectorization of
164:    the matrix-vector product. At the cost of increased storage, the AIJ formatted
165:    matrix can be copied to a format in which pieces of the matrix are
166:    stored in ELLPACK format, allowing the vectorized matrix multiply
167:    routine to use stride-1 memory accesses.  As with the AIJ type, it is
168:    important to preallocate matrix storage in order to get good assembly
169:    performance.

171:    Collective on MPI_Comm

173:    Input Parameters:
174: +  comm - MPI communicator, set to PETSC_COMM_SELF
175: .  m - number of rows
176: .  n - number of columns
177: .  nz - number of nonzeros per row (same for all rows)
178: -  nnz - array containing the number of nonzeros in the various rows
179:          (possibly different for each row) or NULL

181:    Output Parameter:
182: .  A - the matrix

184:    Notes:
185:    If nnz is given then nz is ignored

187:    Level: intermediate

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

191: .seealso: MatCreate(), MatCreateMPIAIJPERM(), MatSetValues()
192: @*/
193: PetscErrorCode  MatCreateMPIAIJCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[],Mat *A)
194: {

198:   MatCreate(comm,A);
199:   MatSetSizes(*A,m,n,m,n);
200:   MatSetType(*A,MATMPIAIJCRL);
201:   MatMPIAIJSetPreallocation_MPIAIJ(*A,nz,(PetscInt*)nnz,onz,(PetscInt*)onnz);
202:   return(0);
203: }

207: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJCRL(Mat A)
208: {

212:   MatSetType(A,MATMPIAIJ);
213:   MatConvert_MPIAIJ_MPIAIJCRL(A,MATMPIAIJCRL,MAT_INPLACE_MATRIX,&A);
214:   return(0);
215: }