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>
16: PetscErrorCode MatDestroy_SeqAIJCRL(Mat A) 17: {
19: Mat_AIJCRL *aijcrl = (Mat_AIJCRL *) A->spptr;
21: /* Free everything in the Mat_AIJCRL data structure. */
22: if (aijcrl) {
23: PetscFree2(aijcrl->acols,aijcrl->icols);
24: }
25: PetscFree(A->spptr);
26: PetscObjectChangeTypeName( (PetscObject)A, MATSEQAIJ);
27: MatDestroy_SeqAIJ(A);
28: return(0);
29: }
31: PetscErrorCode MatDuplicate_AIJCRL(Mat A, MatDuplicateOption op, Mat *M) 32: {
34: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot duplicate AIJCRL matrices yet");
35: return(0);
36: }
40: PetscErrorCode MatSeqAIJCRL_create_aijcrl(Mat A) 41: {
42: Mat_SeqAIJ *a = (Mat_SeqAIJ *)(A)->data;
43: Mat_AIJCRL *aijcrl = (Mat_AIJCRL*) A->spptr;
44: PetscInt m = A->rmap->n; /* Number of rows in the matrix. */
45: PetscInt *aj = a->j; /* From the CSR representation; points to the beginning of each row. */
46: PetscInt i, j,rmax = a->rmax,*icols, *ilen = a->ilen;
47: MatScalar *aa = a->a;
48: PetscScalar *acols;
52: aijcrl->nz = a->nz;
53: aijcrl->m = A->rmap->n;
54: aijcrl->rmax = rmax;
55: PetscFree2(aijcrl->acols,aijcrl->icols);
56: PetscMalloc2(rmax*m,PetscScalar,&aijcrl->acols,rmax*m,PetscInt,&aijcrl->icols);
57: acols = aijcrl->acols;
58: icols = aijcrl->icols;
59: for (i=0; i<m; i++) {
60: for (j=0; j<ilen[i]; j++) {
61: acols[j*m+i] = *aa++;
62: icols[j*m+i] = *aj++;
63: }
64: for (;j<rmax; j++) { /* empty column entries */
65: acols[j*m+i] = 0.0;
66: icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0; /* handle case where row is EMPTY */
67: }
68: }
69: PetscInfo2(A,"Percentage of 0's introduced for vectorized multiply %g. Rmax= %D\n",1.0-((double)a->nz)/((double)(rmax*m)),rmax);
70: return(0);
71: }
73: extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
77: PetscErrorCode MatAssemblyEnd_SeqAIJCRL(Mat A, MatAssemblyType mode) 78: {
80: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
83: a->inode.use = PETSC_FALSE;
84: MatAssemblyEnd_SeqAIJ(A,mode);
85: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
87: /* Now calculate the permutation and grouping information. */
88: MatSeqAIJCRL_create_aijcrl(A);
89: return(0);
90: }
92: #include <../src/mat/impls/aij/seq/crl/ftn-kernels/fmultcrl.h>
96: /*
97: Shared by both sequential and parallel versions of CRL matrix: MATMPIAIJCRL and MATSEQAIJCRL
98: - the scatter is used only in the parallel version
100: */
101: PetscErrorCode MatMult_AIJCRL(Mat A,Vec xx,Vec yy)102: {
103: Mat_AIJCRL *aijcrl = (Mat_AIJCRL*) A->spptr;
104: PetscInt m = aijcrl->m; /* Number of rows in the matrix. */
105: PetscInt rmax = aijcrl->rmax,*icols = aijcrl->icols;
106: PetscScalar *acols = aijcrl->acols;
108: PetscScalar *x,*y;
109: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
110: PetscInt i,j,ii;
111: #endif
114: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
115: #pragma disjoint(*x,*y,*aa)
116: #endif
119: if (aijcrl->xscat) {
120: VecCopy(xx,aijcrl->xwork);
121: /* get remote values needed for local part of multiply */
122: VecScatterBegin(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);
123: VecScatterEnd(aijcrl->xscat,xx,aijcrl->fwork,INSERT_VALUES,SCATTER_FORWARD);
124: xx = aijcrl->xwork;
125: };
127: VecGetArray(xx,&x);
128: VecGetArray(yy,&y);
130: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
131: fortranmultcrl_(&m,&rmax,x,y,icols,acols);
132: #else
134: /* first column */
135: for (j=0; j<m; j++) {
136: y[j] = acols[j]*x[icols[j]];
137: }
139: /* other columns */
140: #if defined(PETSC_HAVE_CRAY_VECTOR)
141: #pragma _CRI preferstream
142: #endif
143: for (i=1; i<rmax; i++) {
144: ii = i*m;
145: #if defined(PETSC_HAVE_CRAY_VECTOR)
146: #pragma _CRI prefervector
147: #endif
148: for (j=0; j<m; j++) {
149: y[j] = y[j] + acols[ii+j]*x[icols[ii+j]];
150: }
151: }
152: #if defined(PETSC_HAVE_CRAY_VECTOR)
153: #pragma _CRI ivdep
154: #endif
156: #endif
157: PetscLogFlops(2.0*aijcrl->nz - m);
158: VecRestoreArray(xx,&x);
159: VecRestoreArray(yy,&y);
160: return(0);
161: }
164: /* MatConvert_SeqAIJ_SeqAIJCRL converts a SeqAIJ matrix into a
165: * SeqAIJCRL matrix. This routine is called by the MatCreate_SeqAIJCRL()
166: * routine, but can also be used to convert an assembled SeqAIJ matrix
167: * into a SeqAIJCRL one. */
168: EXTERN_C_BEGIN
171: PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat A,const MatType type,MatReuse reuse,Mat *newmat)172: {
174: Mat B = *newmat;
175: Mat_AIJCRL *aijcrl;
178: if (reuse == MAT_INITIAL_MATRIX) {
179: MatDuplicate(A,MAT_COPY_VALUES,&B);
180: }
182: PetscNewLog(B,Mat_AIJCRL,&aijcrl);
183: B->spptr = (void *) aijcrl;
185: /* Set function pointers for methods that we inherit from AIJ but override. */
186: B->ops->duplicate = MatDuplicate_AIJCRL;
187: B->ops->assemblyend = MatAssemblyEnd_SeqAIJCRL;
188: B->ops->destroy = MatDestroy_SeqAIJCRL;
189: B->ops->mult = MatMult_AIJCRL;
191: /* If A has already been assembled, compute the permutation. */
192: if (A->assembled) {
193: MatSeqAIJCRL_create_aijcrl(B);
194: }
195: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCRL);
196: *newmat = B;
197: return(0);
198: }
199: EXTERN_C_END
204: /*@C
205: MatCreateSeqAIJCRL - Creates a sparse matrix of type SEQAIJCRL.
206: This type inherits from AIJ, but stores some additional
207: information that is used to allow better vectorization of
208: the matrix-vector product. At the cost of increased storage, the AIJ formatted
209: matrix can be copied to a format in which pieces of the matrix are
210: stored in ELLPACK format, allowing the vectorized matrix multiply
211: routine to use stride-1 memory accesses. As with the AIJ type, it is
212: important to preallocate matrix storage in order to get good assembly
213: performance.
214: 215: Collective on MPI_Comm217: Input Parameters:
218: + comm - MPI communicator, set to PETSC_COMM_SELF219: . m - number of rows
220: . n - number of columns
221: . nz - number of nonzeros per row (same for all rows)
222: - nnz - array containing the number of nonzeros in the various rows
223: (possibly different for each row) or PETSC_NULL225: Output Parameter:
226: . A - the matrix
228: Notes:
229: If nnz is given then nz is ignored
231: Level: intermediate
233: .keywords: matrix, cray, sparse, parallel
235: .seealso: MatCreate(), MatCreateMPIAIJPERM(), MatSetValues()
236: @*/
237: PetscErrorCodeMatCreateSeqAIJCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)238: {
242: MatCreate(comm,A);
243: MatSetSizes(*A,m,n,m,n);
244: MatSetType(*A,MATSEQAIJCRL);
245: MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
246: return(0);
247: }
250: EXTERN_C_BEGIN
253: PetscErrorCode MatCreate_SeqAIJCRL(Mat A)254: {
258: MatSetType(A,MATSEQAIJ);
259: MatConvert_SeqAIJ_SeqAIJCRL(A,MATSEQAIJCRL,MAT_REUSE_MATRIX,&A);
260: return(0);
261: }
262: EXTERN_C_END