xref: /petsc/src/mat/impls/aij/seq/crl/crl.c (revision d32f9abdbc052d6e1fd06679b17a55415c3aae30)
1 #define PETSCMAT_DLL
2 
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.
9 
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"
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "MatDestroy_SeqCRL"
17 PetscErrorCode MatDestroy_SeqCRL(Mat A)
18 {
19   PetscErrorCode ierr;
20   Mat_CRL        *crl = (Mat_CRL *) A->spptr;
21 
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;
30 
31   /* Free everything in the Mat_CRL data structure. */
32   ierr = PetscFree2(crl->acols,crl->icols);CHKERRQ(ierr);
33 
34   /* Change the type of A back to SEQAIJ and use MatDestroy()
35    * to destroy everything that remains. */
36   ierr = PetscObjectChangeTypeName( (PetscObject)A, MATSEQAIJ);CHKERRQ(ierr);
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   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
40   PetscFunctionReturn(0);
41 }
42 
43 PetscErrorCode MatDuplicate_CRL(Mat A, MatDuplicateOption op, Mat *M)
44 {
45   PetscErrorCode ierr;
46   Mat_CRL        *crl = (Mat_CRL *) A->spptr;
47 
48   PetscFunctionBegin;
49   ierr = (*crl->MatDuplicate)(A,op,M);CHKERRQ(ierr);
50   SETERRQ(PETSC_ERR_SUP,"Cannot duplicate CRL matrices yet");
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "SeqCRL_create_crl"
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;
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68   crl->nz   = a->nz;
69   crl->m    = A->rmap.n;
70   crl->rmax = rmax;
71   ierr  = PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);CHKERRQ(ierr);
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   ierr = PetscInfo2(A,"Percentage of 0's introduced for vectorized multiply %G. Rmax= %D\n",1.0-((double)a->nz)/((double)(rmax*m)),rmax);
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "MatAssemblyEnd_SeqCRL"
90 PetscErrorCode MatAssemblyEnd_SeqCRL(Mat A, MatAssemblyType mode)
91 {
92   PetscErrorCode ierr;
93   Mat_CRL        *crl = (Mat_CRL*) A->spptr;
94   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
95 
96   PetscFunctionBegin;
97   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(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);
109 
110   /* Now calculate the permutation and grouping information. */
111   ierr = SeqCRL_create_crl(A);CHKERRQ(ierr);
112   PetscFunctionReturn(0);
113 }
114 
115 #include "src/inline/dot.h"
116 
117 #undef __FUNCT__
118 #define __FUNCT__ "MatMult_CRL"
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;
125   PetscErrorCode ierr;
126   PetscScalar    *x,*y;
127 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
128   PetscInt       i,j,ii;
129 #endif
130 
131 
132 #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
133 #pragma disjoint(*x,*y,*aa)
134 #endif
135 
136   PetscFunctionBegin;
137   if (crl->xscat) {
138     ierr = VecCopy(xx,crl->xwork);CHKERRQ(ierr);
139     /* get remote values needed for local part of multiply */
140     ierr = VecScatterBegin(crl->xscat,xx,crl->fwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
141     ierr = VecScatterEnd(crl->xscat,xx,crl->fwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
142     xx = crl->xwork;
143   };
144 
145   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
146   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
147 
148 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
149   fortranmultcrl_(&m,&rmax,x,y,icols,acols);
150 #else
151 
152   /* first column */
153   for (j=0; j<m; j++) {
154     y[j] = acols[j]*x[icols[j]];
155   }
156 
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
173 
174 #endif
175   ierr = PetscLogFlops(2*crl->nz - m);CHKERRQ(ierr);
176   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
177   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 
181 
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. */
186 EXTERN_C_BEGIN
187 #undef __FUNCT__
188 #define __FUNCT__ "MatConvert_SeqAIJ_SeqCRL"
189 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCRL(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
190 {
191   PetscErrorCode ierr;
192   Mat            B = *newmat;
193   Mat_CRL        *crl;
194 
195   PetscFunctionBegin;
196   if (reuse == MAT_INITIAL_MATRIX) {
197     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
198   }
199 
200   ierr = PetscNewLog(B,Mat_CRL,&crl);CHKERRQ(ierr);
201   B->spptr = (void *) crl;
202 
203   crl->AssemblyEnd  = A->ops->assemblyend;
204   crl->MatDestroy   = A->ops->destroy;
205   crl->MatDuplicate = A->ops->duplicate;
206 
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;
212 
213   /* If A has already been assembled, compute the permutation. */
214   if (A->assembled == PETSC_TRUE) {
215     ierr = SeqCRL_create_crl(B);CHKERRQ(ierr);
216   }
217   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQCRL);CHKERRQ(ierr);
218   *newmat = B;
219   PetscFunctionReturn(0);
220 }
221 EXTERN_C_END
222 
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "MatCreateSeqCRL"
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
238 
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
246 
247    Output Parameter:
248 .  A - the matrix
249 
250    Notes:
251    If nnz is given then nz is ignored
252 
253    Level: intermediate
254 
255 .keywords: matrix, cray, sparse, parallel
256 
257 .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues()
258 @*/
259 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
260 {
261   PetscErrorCode ierr;
262 
263   PetscFunctionBegin;
264   ierr = MatCreate(comm,A);CHKERRQ(ierr);
265   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
266   ierr = MatSetType(*A,MATSEQCRL);CHKERRQ(ierr);
267   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 
272 EXTERN_C_BEGIN
273 #undef __FUNCT__
274 #define __FUNCT__ "MatCreate_SeqCRL"
275 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqCRL(Mat A)
276 {
277   PetscErrorCode ierr;
278 
279   PetscFunctionBegin;
280   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
281   ierr = MatConvert_SeqAIJ_SeqCRL(A,MATSEQCRL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
282   PetscFunctionReturn(0);
283 }
284 EXTERN_C_END
285 
286