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