xref: /petsc/src/mat/impls/aij/seq/crl/crl.c (revision 009ddefdd03f050e7a6a2a9c71891d751851e81d)
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 
14 #include "src/mat/impls/aij/seq/aij.h"
15 
16 typedef struct {
17   PetscInt    ncols;    /* number of columns in each row */
18   PetscInt    *icols;   /* columns of nonzeros, stored one column at a time */
19   PetscScalar *acols;   /* values of nonzeros, stored as icols */
20 
21   /* We need to keep a pointer to MatAssemblyEnd_SeqAIJ because we
22    * actually want to call this function from within the
23    * MatAssemblyEnd_SeqCRL function.  Similarly, we also need
24    * MatDestroy_SeqAIJ and MatDuplicate_SeqAIJ. */
25   PetscErrorCode (*AssemblyEnd_SeqAIJ)(Mat,MatAssemblyType);
26   PetscErrorCode (*MatDestroy_SeqAIJ)(Mat);
27   PetscErrorCode (*MatDuplicate_SeqAIJ)(Mat,MatDuplicateOption,Mat*);
28 } Mat_SeqCRL;
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "MatDestroy_SeqCRL"
32 PetscErrorCode MatDestroy_SeqCRL(Mat A)
33 {
34   PetscErrorCode ierr;
35   Mat_SeqCRL     *crl = (Mat_SeqCRL *) A->spptr;
36 
37   /* We are going to convert A back into a SEQAIJ matrix, since we are
38    * eventually going to use MatDestroy_SeqAIJ() to destroy everything
39    * that is not specific to CRL.
40    * In preparation for this, reset the operations pointers in A to
41    * their SeqAIJ versions. */
42   A->ops->assemblyend = crl->AssemblyEnd_SeqAIJ;
43   A->ops->destroy     = crl->MatDestroy_SeqAIJ;
44   A->ops->duplicate   = crl->MatDuplicate_SeqAIJ;
45 
46   /* Free everything in the Mat_SeqCRL data structure. */
47   if (crl->icols) {
48     ierr = PetscFree2(crl->acols,crl->icols);CHKERRQ(ierr);
49   }
50   /* Free the Mat_SeqCRL struct itself. */
51   ierr = PetscFree(crl);CHKERRQ(ierr);
52 
53   /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
54    * to destroy everything that remains. */
55   ierr = PetscObjectChangeTypeName( (PetscObject)A, MATSEQAIJ);CHKERRQ(ierr);
56   /* Note that I don't call MatSetType().  I believe this is because that
57    * is only to be called when *building* a matrix. */
58   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
59   PetscFunctionReturn(0);
60 }
61 
62 PetscErrorCode MatDuplicate_SeqCRL(Mat A, MatDuplicateOption op, Mat *M)
63 {
64   PetscErrorCode ierr;
65   Mat_SeqCRL     *crl = (Mat_SeqCRL *) A->spptr;
66 
67   PetscFunctionBegin;
68   ierr = (*crl->MatDuplicate_SeqAIJ)(A,op,M);CHKERRQ(ierr);
69   SETERRQ(PETSC_ERR_SUP,"Cannot duplicate CRL matrices yet");
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "SeqCRL_create_crl"
75 PetscErrorCode SeqCRL_create_crl(Mat A)
76 {
77   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)(A)->data;
78   Mat_SeqCRL     *crl = (Mat_SeqCRL*) A->spptr;
79   PetscInt       m = A->m;  /* Number of rows in the matrix. */
80   PetscInt       *aj = a->j;  /* From the CSR representation; points to the beginning  of each row. */
81   PetscInt       i, j,rmax = a->rmax,*icols, *ilen = a->ilen;
82   PetscScalar    *aa = a->a,*acols;
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   ierr  = PetscMalloc2(rmax*m,PetscScalar,&crl->acols,rmax*m,PetscInt,&crl->icols);CHKERRQ(ierr);
87   acols = crl->acols;
88   icols = crl->icols;
89   for (i=0; i<m; i++) {
90     for (j=0; j<ilen[i]; j++) {
91       acols[j*m+i] = *aa++;
92       icols[j*m+i] = *aj++;
93     }
94     for (;j<rmax; j++) { /* empty column entries */
95       acols[j*m+i] = 0.0;
96       icols[j*m+i] = (j) ? icols[(j-1)*m+i] : 0;  /* handle case where row is EMPTY */
97     }
98   }
99   ierr = PetscLogInfo((A,"SeqCRL_create_crl: Percentage of 0's introduced for vectorized multiply %g\n",1.0-((double)a->nz)/((double)(rmax*m))));
100   PetscFunctionReturn(0);
101 }
102 
103 #undef __FUNCT__
104 #define __FUNCT__ "MatAssemblyEnd_SeqCRL"
105 PetscErrorCode MatAssemblyEnd_SeqCRL(Mat A, MatAssemblyType mode)
106 {
107   PetscErrorCode ierr;
108   Mat_SeqCRL     *crl = (Mat_SeqCRL*) A->spptr;
109   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
110 
111   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
112 
113   /* Since a MATSEQCRL matrix is really just a MATSEQAIJ with some
114    * extra information, call the AssemblyEnd routine for a MATSEQAIJ.
115    * I'm not sure if this is the best way to do this, but it avoids
116    * a lot of code duplication.
117    * I also note that currently MATSEQCRL doesn't know anything about
118    * the Mat_CompressedRow data structure that SeqAIJ now uses when there
119    * are many zero rows.  If the SeqAIJ assembly end routine decides to use
120    * this, this may break things.  (Don't know... haven't looked at it.) */
121   a->inode.use = PETSC_FALSE;
122   (*crl->AssemblyEnd_SeqAIJ)(A, mode);
123 
124   /* Now calculate the permutation and grouping information. */
125   ierr = SeqCRL_create_crl(A);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "MatMult_SeqCRL"
131 PetscErrorCode MatMult_SeqCRL(Mat A,Vec xx,Vec yy)
132 {
133   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)(A)->data;
134   Mat_SeqCRL     *crl = (Mat_SeqCRL*) A->spptr;
135   PetscInt       m = A->m;  /* Number of rows in the matrix. */
136   PetscInt       rmax = a->rmax,*icols = crl->icols;
137   PetscScalar    *acols = crl->acols;
138   PetscErrorCode ierr;
139   PetscScalar    *x,*y;
140 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
141   PetscInt       i,j,ii;
142 #endif
143 
144 
145 #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
146 #pragma disjoint(*x,*y,*aa)
147 #endif
148 
149   PetscFunctionBegin;
150   ierr = VecZeroEntries(yy);CHKERRQ(ierr);
151   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
152   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
153 
154 
155 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
156 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
157   fortranmultcrl_(&m,&rmax,x,y,icols,acols);
158 #else
159 
160 #if defined(PETSC_HAVE_CRAYC)
161 #pragma _CRI preferstream
162 #endif
163   for (i=0; i<rmax; i++) {
164     ii = i*m;
165 #if defined(PETSC_HAVE_CRAYC)
166 #pragma _CRI prefervector
167 #endif
168     for (j=0; j<m; j++) {
169       y[j] = y[j] + acols[ii+j]*x[icols[ii+j]];
170     }
171   }
172 #if defined(PETSC_HAVE_CRAYC)
173 #pragma _CRI ivdep
174 #endif
175 
176 #endif
177   ierr = PetscLogFlops(2*a->nz - m);CHKERRQ(ierr);
178   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
179   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
180   PetscFunctionReturn(0);
181 }
182 
183 
184 /* MatConvert_SeqAIJ_SeqCRL converts a SeqAIJ matrix into a
185  * SeqCRL matrix.  This routine is called by the MatCreate_SeqCRL()
186  * routine, but can also be used to convert an assembled SeqAIJ matrix
187  * into a SeqCRL one. */
188 EXTERN_C_BEGIN
189 #undef __FUNCT__
190 #define __FUNCT__ "MatConvert_SeqAIJ_SeqCRL"
191 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCRL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
192 {
193   /* This routine is only called to convert to MATSEQCRL
194    * from MATSEQAIJ, so we can ignore 'MatType Type'. */
195   PetscErrorCode ierr;
196   Mat            B = *newmat;
197   Mat_SeqCRL     *crl;
198 
199   PetscFunctionBegin;
200   if (reuse == MAT_INITIAL_MATRIX) {
201     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
202   }
203 
204   ierr = PetscNew(Mat_SeqCRL,&crl);CHKERRQ(ierr);
205   B->spptr = (void *) crl;
206 
207   /* Save a pointer to the original SeqAIJ assembly end routine, because we
208    * will want to use it later in the CRL assembly end routine.
209    * Also, save a pointer to the original SeqAIJ Destroy routine, because we
210    * will want to use it in the CRL destroy routine. */
211   crl->AssemblyEnd_SeqAIJ  = A->ops->assemblyend;
212   crl->MatDestroy_SeqAIJ   = A->ops->destroy;
213   crl->MatDuplicate_SeqAIJ = A->ops->duplicate;
214 
215   /* Set function pointers for methods that we inherit from AIJ but
216    * override. */
217   B->ops->duplicate   = MatDuplicate_SeqCRL;
218   B->ops->assemblyend = MatAssemblyEnd_SeqCRL;
219   B->ops->destroy     = MatDestroy_SeqCRL;
220   B->ops->mult        = MatMult_SeqCRL;
221 
222   /* If A has already been assembled, compute the permutation. */
223   if (A->assembled == PETSC_TRUE) {
224     ierr = SeqCRL_create_crl(B);CHKERRQ(ierr);
225   }
226   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQCRL);CHKERRQ(ierr);
227   *newmat = B;
228   PetscFunctionReturn(0);
229 }
230 EXTERN_C_END
231 
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "MatCreateSeqCRL"
235 /*@C
236    MatCreateSeqCRL - Creates a sparse matrix of type SEQCRL.
237    This type inherits from AIJ, but stores some additional
238    information that is used to allow better vectorization of
239    the matrix-vector product. At the cost of increased storage, the AIJ formatted
240    matrix can be copied to a format in which pieces of the matrix are
241    stored in ELLPACK format, allowing the vectorized matrix multiply
242    routine to use stride-1 memory accesses.  As with the AIJ type, it is
243    important to preallocate matrix storage in order to get good assembly
244    performance.
245 
246    Collective on MPI_Comm
247 
248    Input Parameters:
249 +  comm - MPI communicator, set to PETSC_COMM_SELF
250 .  m - number of rows
251 .  n - number of columns
252 .  nz - number of nonzeros per row (same for all rows)
253 -  nnz - array containing the number of nonzeros in the various rows
254          (possibly different for each row) or PETSC_NULL
255 
256    Output Parameter:
257 .  A - the matrix
258 
259    Notes:
260    If nnz is given then nz is ignored
261 
262    Level: intermediate
263 
264 .keywords: matrix, cray, sparse, parallel
265 
266 .seealso: MatCreate(), MatCreateMPICSRPERM(), MatSetValues()
267 @*/
268 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
269 {
270   PetscErrorCode ierr;
271 
272   PetscFunctionBegin;
273   ierr = MatCreate(comm,A);CHKERRQ(ierr);
274   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
275   ierr = MatSetType(*A,MATSEQCRL);CHKERRQ(ierr);
276   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 
281 EXTERN_C_BEGIN
282 #undef __FUNCT__
283 #define __FUNCT__ "MatCreate_SeqCRL"
284 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqCRL(Mat A)
285 {
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   /* Change the type name before calling MatSetType() to force proper construction of SeqAIJ
290      and MATSEQCRL types. */
291   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSEQCRL);CHKERRQ(ierr);
292   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
293   ierr = MatConvert_SeqAIJ_SeqCRL(A,MATSEQCRL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
294   PetscFunctionReturn(0);
295 }
296 EXTERN_C_END
297 
298