xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 6284ec503f7d48a443257752306f1a524cdcb213)
1 /*
2   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
3           C = A * B
4 */
5 
6 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
7 #include "src/mat/utils/freespace.h"
8 
9 /*
10   Add a index set into a sorted linked list
11   input:
12     nidx      - number of input indices
13     indices   - interger array
14     idx_head  - the header of the list
15     idx_unset - the value indicating the entry in the list is not set yet
16     lnk       - linked list(an integer array) that is created
17   output:
18     nlnk      - number of newly added indices
19     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
20  */
21 #undef __FUNCT__
22 #define __FUNCT__ "LnklistAdd"
23 int LnklistAdd(int nidx,int *indices,int idx_head,int idx_unset,int *nlnk,int *lnk)
24 {
25   int i,idx,lidx,entry,n;
26 
27   PetscFunctionBegin;
28   n    = 0;
29   lidx = idx_head;
30   i    = nidx;
31   while (i){
32     /* assume indices is almost in increasing order, starting from its end saves computation */
33     entry = indices[--i];
34     if (lnk[entry] == idx_unset) { /* new entry */
35       do {
36         idx = lidx;
37         lidx  = lnk[idx];
38       } while (entry > lidx);
39       lnk[idx] = entry;
40       lnk[entry] = lidx;
41       n++;
42     }
43   }
44   *nlnk = n;
45   PetscFunctionReturn(0);
46 }
47 
48 
49 static int logkey_matmatmult          = 0;
50 static int logkey_matmatmult_symbolic = 0;
51 static int logkey_matmatmult_numeric  = 0;
52 
53 #undef __FUNCT__
54 #define __FUNCT__ "MatMatMult"
55 /*@
56    MatMatMult - Performs Matrix-Matrix Multiplication C=A*B.
57 
58    Collective on Mat
59 
60    Input Parameters:
61 +  A - the left matrix
62 -  B - the right matrix
63 
64    Output Parameters:
65 .  C - the product matrix
66 
67    Notes:
68    C will be created and must be destroyed by the user with MatDestroy().
69 
70    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
71    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
72 
73    Level: intermediate
74 
75 .seealso: MatMatMultSymbolic(),MatMatMultNumeric()
76 @*/
77 int MatMatMult(Mat A,Mat B, Mat *C) {
78   /* Perhaps this "interface" routine should be moved into the interface directory.*/
79   /* To facilitate implementations with varying types, QueryFunction is used.*/
80   /* It is assumed that implementations will be composed as "MatMatMult_<type of A><type of B>". */
81   int  ierr;
82 #ifdef OLD
83   char funct[80];
84   int  (*mult)(Mat,Mat,Mat*);
85 #endif
86   PetscFunctionBegin;
87   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
88   PetscValidType(A,1);
89   MatPreallocated(A);
90   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
91   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
92 
93   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
94   PetscValidType(B,2);
95   MatPreallocated(B);
96   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
97   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
98 
99   PetscValidPointer(C,3);
100 
101   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
102 
103   ierr = PetscLogEventBegin(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
104   ierr = (*A->ops->matmult)(A,B,C);CHKERRQ(ierr);
105   ierr = PetscLogEventEnd(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
106 
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
112 int MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B, Mat *C) {
113   int ierr;
114 
115   PetscFunctionBegin;
116   SETERRQ(PETSC_ERR_SUP,"Not written yet");
117 
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
123 int MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B, Mat *C) {
124   int ierr;
125   char symfunct[80],numfunct[80];
126   int (*symbolic)(Mat,Mat,Mat*),(*numeric)(Mat,Mat,Mat);
127 
128   PetscFunctionBegin;
129   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
130   /* When other implementations exist, attack the multiple dispatch problem. */
131   ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr);
132   ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
133   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
134   ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
135   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
136 
137   ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr);
138   ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
139   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
140   ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
141   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
142 
143   ierr = PetscLogEventBegin(logkey_matmatmult,A,B,0,0);CHKERRQ(ierr);
144   ierr = (*symbolic)(A,B,C);CHKERRQ(ierr);
145   ierr = (*numeric)(A,B,*C);CHKERRQ(ierr);
146   ierr = PetscLogEventEnd(logkey_matmatmult,A,B,0,0);CHKERRQ(ierr);
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "MatMatMultSymbolic"
152 /*@
153    MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure
154    of the matrix-matrix product C=A*B.  Call this routine before calling MatMatMultNumeric().
155 
156    Collective on Mat
157 
158    Input Parameters:
159 +  A - the left matrix
160 -  B - the right matrix
161 
162    Output Parameters:
163 .  C - the matrix containing the ij structure of product matrix
164 
165    Notes:
166    C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy().
167 
168    This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ.
169 
170    Level: intermediate
171 
172 .seealso: MatMatMult(),MatMatMultNumeric()
173 @*/
174 int MatMatMultSymbolic(Mat A,Mat B,Mat *C) {
175   /* Perhaps this "interface" routine should be moved into the interface directory.*/
176   /* To facilitate implementations with varying types, QueryFunction is used.*/
177   /* It is assumed that implementations will be composed as "MatMatMultSymbolic_<type of A><type of B>". */
178   int  ierr;
179   char symfunct[80];
180   int  (*symbolic)(Mat,Mat,Mat *);
181 
182   PetscFunctionBegin;
183   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
184   PetscValidType(A,1);
185   MatPreallocated(A);
186   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
187   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
188 
189   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
190   PetscValidType(B,2);
191   MatPreallocated(B);
192   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
193   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
194   PetscValidPointer(C,3);
195 
196 
197   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
198 
199   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
200   /* When other implementations exist, attack the multiple dispatch problem. */
201   ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr);
202   ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
203   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
204   ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
205   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
206   ierr = (*symbolic)(A,B,C);CHKERRQ(ierr);
207 
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "MatMatMult_Symbolic_SeqAIJ_SeqAIJ"
213 int MatMatMult_Symbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat *C)
214 {
215   int            ierr;
216   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
217   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
218   int            *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj;
219   int            *ci,*cj,*lnk,idx0,idx;
220   int            am=A->M,bn=B->N,bm=B->M;
221   int            i,j,anzi,brow,bnzj,cnzi,nlnk;
222   MatScalar      *ca;
223 
224   PetscFunctionBegin;
225   /* Start timers */
226   ierr = PetscLogEventBegin(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
227   /* Set up */
228   /* Allocate ci array, arrays for fill computation and */
229   /* free space for accumulating nonzero column info */
230   ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr);
231   ci[0] = 0;
232 
233   ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr);
234   for (i=0; i<bn; i++) lnk[i] = -1;
235 
236   /* Initial FreeSpace size is nnz(B)=4*bi[bm] */
237   ierr = GetMoreSpace(4*bi[bm],&free_space);CHKERRQ(ierr);
238   current_space = free_space;
239 
240   /* Determine symbolic info for each row of the product: */
241   for (i=0;i<am;i++) {
242     anzi = ai[i+1] - ai[i];
243     cnzi = 0;
244     lnk[bn] = bn;
245     for (j=0;j<anzi;j++) {
246       brow = *aj++;
247       bnzj = bi[brow+1] - bi[brow];
248       bjj  = bj + bi[brow];
249       /* add non-zero cols of B into the sorted linked list lnk */
250       ierr = LnklistAdd(bnzj,bjj,bn,-1,&nlnk,lnk);CHKERRQ(ierr);
251       cnzi += nlnk;
252     }
253 
254     /* If free space is not available, make more free space */
255     /* Double the amount of total space in the list */
256     if (current_space->local_remaining<cnzi) {
257       printf("...%d -th row, double space ...\n",i);
258       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
259     }
260 
261     /* Copy data into free space, and zero out denserow and lnk */
262     idx = bn;
263     for (j=0; j<cnzi; j++){
264       idx0 = idx;
265       idx  = lnk[idx0];
266       *current_space->array++ = idx;
267       lnk[idx0] = -1;
268     }
269     lnk[idx] = -1;
270 
271     current_space->local_used      += cnzi;
272     current_space->local_remaining -= cnzi;
273 
274     ci[i+1] = ci[i] + cnzi;
275   }
276 
277   /* Column indices are in the list of free space */
278   /* Allocate space for cj, initialize cj, and */
279   /* destroy list of free space and other temporary array(s) */
280   ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr);
281   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
282   ierr = PetscFree(lnk);CHKERRQ(ierr);
283 
284   /* Allocate space for ca */
285   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
286   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
287 
288   /* put together the new matrix */
289   ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
290 
291   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
292   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
293   c = (Mat_SeqAIJ *)((*C)->data);
294   c->freedata = PETSC_TRUE;
295   c->nonew    = 0;
296 
297   ierr = PetscLogEventEnd(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
298   PetscFunctionReturn(0);
299 }
300 
301 #undef __FUNCT__
302 #define __FUNCT__ "MatMatMultNumeric"
303 /*@
304    MatMatMultNumeric - Performs the numeric matrix-matrix product.
305    Call this routine after first calling MatMatMultSymbolic().
306 
307    Collective on Mat
308 
309    Input Parameters:
310 +  A - the left matrix
311 -  B - the right matrix
312 
313    Output Parameters:
314 .  C - the product matrix, whose ij structure was defined from MatMatMultSymbolic().
315 
316    Notes:
317    C must have been created with MatMatMultSymbolic.
318 
319    This routine is currently only implemented for SeqAIJ type matrices.
320 
321    Level: intermediate
322 
323 .seealso: MatMatMult(),MatMatMultSymbolic()
324 @*/
325 int MatMatMultNumeric(Mat A,Mat B,Mat C){
326   /* Perhaps this "interface" routine should be moved into the interface directory.*/
327   /* To facilitate implementations with varying types, QueryFunction is used.*/
328   /* It is assumed that implementations will be composed as "MatMatMultNumeric_<type of A><type of B>". */
329   int ierr;
330   char numfunct[80];
331   int (*numeric)(Mat,Mat,Mat);
332 
333   PetscFunctionBegin;
334 
335   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
336   PetscValidType(A,1);
337   MatPreallocated(A);
338   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
339   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
340 
341   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
342   PetscValidType(B,2);
343   MatPreallocated(B);
344   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
345   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
346 
347   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
348   PetscValidType(C,3);
349   MatPreallocated(C);
350   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
351   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
352 
353   if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N);
354   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
355   if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M);
356 
357   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
358   /* When other implementations exist, attack the multiple dispatch problem. */
359   ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr);
360   ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
361   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
362   ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
363   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
364 
365   ierr = (*numeric)(A,B,C);CHKERRQ(ierr);
366 
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "MatMatMult_Numeric_SeqAIJ_SeqAIJ"
372 int MatMatMult_Numeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
373 {
374   int        ierr,flops=0;
375   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
376   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
377   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
378   int        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
379   int        am=A->M,cn=C->N;
380   int        i,j,k,anzi,bnzi,cnzi,brow;
381   MatScalar  *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp;
382 
383   PetscFunctionBegin;
384 
385   /* Start timers */
386   ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
387 
388   /* Allocate temp accumulation space to avoid searching for nonzero columns in C */
389   ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr);
390   ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr);
391   /* Traverse A row-wise. */
392   /* Build the ith row in C by summing over nonzero columns in A, */
393   /* the rows of B corresponding to nonzeros of A. */
394   for (i=0;i<am;i++) {
395     anzi = ai[i+1] - ai[i];
396     for (j=0;j<anzi;j++) {
397       brow = *aj++;
398       bnzi = bi[brow+1] - bi[brow];
399       bjj  = bj + bi[brow];
400       baj  = ba + bi[brow];
401       for (k=0;k<bnzi;k++) {
402         temp[bjj[k]] += (*aa)*baj[k];
403       }
404       flops += 2*bnzi;
405       aa++;
406     }
407     /* Store row back into C, and re-zero temp */
408     cnzi = ci[i+1] - ci[i];
409     for (j=0;j<cnzi;j++) {
410       ca[j] = temp[cj[j]];
411       temp[cj[j]] = 0.0;
412     }
413     ca += cnzi;
414     cj += cnzi;
415   }
416   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
417   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
418 
419   /* Free temp */
420   ierr = PetscFree(temp);CHKERRQ(ierr);
421   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
422   ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
423   PetscFunctionReturn(0);
424 }
425 
426 #undef __FUNCT__
427 #define __FUNCT__ "RegisterMatMatMultRoutines_Private"
428 int RegisterMatMatMultRoutines_Private(Mat A) {
429   int ierr;
430 
431   PetscFunctionBegin;
432 #ifdef OLD
433   if (!logkey_matmatmult) {
434     ierr = PetscLogEventRegister(&logkey_matmatmult,"MatMatMult",MAT_COOKIE);CHKERRQ(ierr);
435   }
436   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMult_seqaijseqaij",
437                                            "MatMatMult_SeqAIJ_SeqAIJ",
438                                            MatMatMult_SeqAIJ_SeqAIJ);CHKERRQ(ierr);
439 #endif
440   if (!logkey_matmatmult_symbolic) {
441     ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMult_Symbolic",MAT_COOKIE);CHKERRQ(ierr);
442   }
443   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultSymbolic_seqaijseqaij",
444                                            "MatMatMult_Symbolic_SeqAIJ_SeqAIJ",
445                                            MatMatMult_Symbolic_SeqAIJ_SeqAIJ);CHKERRQ(ierr);
446   if (!logkey_matmatmult_numeric) {
447     ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMult_Numeric",MAT_COOKIE);CHKERRQ(ierr);
448   }
449   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultNumeric_seqaijseqaij",
450                                            "MatMatMult_Numeric_SeqAIJ_SeqAIJ",
451                                            MatMatMult_Numeric_SeqAIJ_SeqAIJ);CHKERRQ(ierr);
452   PetscFunctionReturn(0);
453 }
454