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: PetscScalar *aa = a->a,*acols;
67: crl->nz = a->nz;
68: crl->m = A->rmap.n;
69: crl->rmax = rmax;
70: PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);
71: acols = crl->acols;
72: icols = crl->icols;
73: for (i=0; i<m; i++) {
74: for (j=0; j<ilen[i]; j++) {
75: acols[j*m+i] = *aa++;
76: icols[j*m+i] = *aj++;
77: }
78: for (;j<rmax; j++) { /* empty column entries */
79: acols[j*m+i] = 0.0;
80: icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0; /* handle case where row is EMPTY */
81: }
82: }
83: PetscInfo2(A,"Percentage of 0's introduced for vectorized multiply %G. Rmax= %D\n",1.0-((double)a->nz)/((double)(rmax*m)),rmax);
84: return(0);
85: }
89: PetscErrorCode MatAssemblyEnd_SeqCRL(Mat A, MatAssemblyType mode)
90: {
92: Mat_CRL *crl = (Mat_CRL*) A->spptr;
93: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
96: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
97:
98: /* Since a MATSEQCRL matrix is really just a MATSEQAIJ with some
99: * extra information, call the AssemblyEnd routine for a MATSEQAIJ.
100: * I'm not sure if this is the best way to do this, but it avoids
101: * a lot of code duplication.
102: * I also note that currently MATSEQCRL doesn't know anything about
103: * the Mat_CompressedRow data structure that SeqAIJ now uses when there
104: * are many zero rows. If the SeqAIJ assembly end routine decides to use
105: * this, this may break things. (Don't know... haven't looked at it.) */
106: a->inode.use = PETSC_FALSE;
107: (*crl->AssemblyEnd)(A, mode);
109: /* Now calculate the permutation and grouping information. */
110: SeqCRL_create_crl(A);
111: return(0);
112: }
114: #include src/inline/dot.h
118: PetscErrorCode MatMult_CRL(Mat A,Vec xx,Vec yy)
119: {
120: Mat_CRL *crl = (Mat_CRL*) A->spptr;
121: PetscInt m = crl->m; /* Number of rows in the matrix. */
122: PetscInt rmax = crl->rmax,*icols = crl->icols;
123: PetscScalar *acols = crl->acols;
125: PetscScalar *x,*y;
126: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
127: PetscInt i,j,ii;
128: #endif
131: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
132: #pragma disjoint(*x,*y,*aa)
133: #endif
136: if (crl->xscat) {
137: VecCopy(xx,crl->xwork);
138: /* get remote values needed for local part of multiply */
139: VecScatterBegin(xx,crl->fwork,INSERT_VALUES,SCATTER_FORWARD,crl->xscat);
140: VecScatterEnd(xx,crl->fwork,INSERT_VALUES,SCATTER_FORWARD,crl->xscat);
141: xx = crl->xwork;
142: };
144: VecGetArray(xx,&x);
145: VecGetArray(yy,&y);
147: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
148: fortranmultcrl_(&m,&rmax,x,y,icols,acols);
149: #else
151: /* first column */
152: for (j=0; j<m; j++) {
153: y[j] = acols[j]*x[icols[j]];
154: }
156: /* other columns */
157: #if defined(PETSC_HAVE_CRAYC)
158: #pragma _CRI preferstream
159: #endif
160: for (i=1; i<rmax; i++) {
161: ii = i*m;
162: #if defined(PETSC_HAVE_CRAYC)
163: #pragma _CRI prefervector
164: #endif
165: for (j=0; j<m; j++) {
166: y[j] = y[j] + acols[ii+j]*x[icols[ii+j]];
167: }
168: }
169: #if defined(PETSC_HAVE_CRAYC)
170: #pragma _CRI ivdep
171: #endif
173: #endif
174: PetscLogFlops(2*crl->nz - m);
175: VecRestoreArray(xx,&x);
176: VecRestoreArray(yy,&y);
177: return(0);
178: }
181: /* MatConvert_SeqAIJ_SeqCRL converts a SeqAIJ matrix into a
182: * SeqCRL matrix. This routine is called by the MatCreate_SeqCRL()
183: * routine, but can also be used to convert an assembled SeqAIJ matrix
184: * into a SeqCRL one. */
188: PetscErrorCode MatConvert_SeqAIJ_SeqCRL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
189: {
191: Mat B = *newmat;
192: Mat_CRL *crl;
195: if (reuse == MAT_INITIAL_MATRIX) {
196: MatDuplicate(A,MAT_COPY_VALUES,&B);
197: }
199: PetscNew(Mat_CRL,&crl);
200: B->spptr = (void *) crl;
202: crl->AssemblyEnd = A->ops->assemblyend;
203: crl->MatDestroy = A->ops->destroy;
204: crl->MatDuplicate = A->ops->duplicate;
206: /* Set function pointers for methods that we inherit from AIJ but override. */
207: B->ops->duplicate = MatDuplicate_CRL;
208: B->ops->assemblyend = MatAssemblyEnd_SeqCRL;
209: B->ops->destroy = MatDestroy_SeqCRL;
210: B->ops->mult = MatMult_CRL;
212: /* If A has already been assembled, compute the permutation. */
213: if (A->assembled == PETSC_TRUE) {
214: SeqCRL_create_crl(B);
215: }
216: PetscObjectChangeTypeName((PetscObject)B,MATSEQCRL);
217: *newmat = B;
218: return(0);
219: }
225: /*@C
226: MatCreateSeqCRL - Creates a sparse matrix of type SEQCRL.
227: This type inherits from AIJ, but stores some additional
228: information that is used to allow better vectorization of
229: the matrix-vector product. At the cost of increased storage, the AIJ formatted
230: matrix can be copied to a format in which pieces of the matrix are
231: stored in ELLPACK format, allowing the vectorized matrix multiply
232: routine to use stride-1 memory accesses. As with the AIJ type, it is
233: important to preallocate matrix storage in order to get good assembly
234: performance.
235:
236: Collective on MPI_Comm
238: Input Parameters:
239: + comm - MPI communicator, set to PETSC_COMM_SELF
240: . m - number of rows
241: . n - number of columns
242: . nz - number of nonzeros per row (same for all rows)
243: - nnz - array containing the number of nonzeros in the various rows
244: (possibly different for each row) or PETSC_NULL
246: Output Parameter:
247: . A - the matrix
249: Notes:
250: If nnz is given then nz is ignored
252: Level: intermediate
254: .keywords: matrix, cray, sparse, parallel
256: .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues()
257: @*/
258: PetscErrorCode MatCreateSeqCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
259: {
263: MatCreate(comm,A);
264: MatSetSizes(*A,m,n,m,n);
265: MatSetType(*A,MATSEQCRL);
266: MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
267: return(0);
268: }
274: PetscErrorCode MatCreate_SeqCRL(Mat A)
275: {
279: MatSetType(A,MATSEQAIJ);
280: MatConvert_SeqAIJ_SeqCRL(A,MATSEQCRL,MAT_REUSE_MATRIX,&A);
281: return(0);
282: }