xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision bfb15611a2f62d5f3c65919701908022350ceac6)
1 #define PETSCMAT_DLL
2 
3 /*
4   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
5           C = A * B
6 */
7 
8 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
9 #include "src/mat/utils/freespace.h"
10 #include "petscbt.h"
11 
12 
13 #undef __FUNCT__
14 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
15 PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
16 {
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   if (scall == MAT_INITIAL_MATRIX){
21     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
22   }
23   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
24   PetscFunctionReturn(0);
25 }
26 
27 
28 #undef __FUNCT__
29 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
30 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
31 {
32   PetscErrorCode     ierr;
33   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
34   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
35   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci,*cj;
36   PetscInt           am=A->rmap.N,bn=B->cmap.N,bm=B->rmap.N;
37   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,nspacedouble=0;
38   MatScalar          *ca;
39   PetscBT            lnkbt;
40 
41   PetscFunctionBegin;
42   /* Set up */
43   /* Allocate ci array, arrays for fill computation and */
44   /* free space for accumulating nonzero column info */
45   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
46   ci[0] = 0;
47 
48   /* create and initialize a linked list */
49   nlnk = bn+1;
50   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
51 
52   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
53   if (fill < 1.0) fill = 1.0; /* In case user input a wrong fill, reset it to 1.0 */
54   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
55   current_space = free_space;
56 
57   /* Determine symbolic info for each row of the product: */
58   for (i=0;i<am;i++) {
59     anzi = ai[i+1] - ai[i];
60     cnzi = 0;
61     j    = anzi;
62     aj   = a->j + ai[i];
63     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
64       j--;
65       brow = *(aj + j);
66       bnzj = bi[brow+1] - bi[brow];
67       bjj  = bj + bi[brow];
68       /* add non-zero cols of B into the sorted linked list lnk */
69       ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
70       cnzi += nlnk;
71     }
72 
73     /* If free space is not available, make more free space */
74     /* Double the amount of total space in the list */
75     if (current_space->local_remaining<cnzi) {
76       ierr = PetscFreeSpaceGet(current_space->total_array_size,&current_space);CHKERRQ(ierr);
77       nspacedouble++;
78     }
79 
80     /* Copy data into free space, then initialize lnk */
81     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
82     current_space->array           += cnzi;
83     current_space->local_used      += cnzi;
84     current_space->local_remaining -= cnzi;
85 
86     ci[i+1] = ci[i] + cnzi;
87   }
88 
89   /* Column indices are in the list of free space */
90   /* Allocate space for cj, initialize cj, and */
91   /* destroy list of free space and other temporary array(s) */
92   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
93   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
94   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
95 
96   /* Allocate space for ca */
97   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
98   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
99 
100   /* put together the new symbolic matrix */
101   ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
102 
103   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
104   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
105   c = (Mat_SeqAIJ *)((*C)->data);
106   c->freedata = PETSC_TRUE;
107   c->nonew    = 0;
108 
109   if (nspacedouble){
110     ierr = PetscInfo5((*C),"nspacedouble:%D, nnz(A):%D, nnz(B):%D, fill:%G, nnz(C):%D\n",nspacedouble,ai[am],bi[bm],fill,ci[am]);CHKERRQ(ierr);
111   }
112   PetscFunctionReturn(0);
113 }
114 
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
118 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
119 {
120   PetscErrorCode ierr;
121   PetscInt       flops=0;
122   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
123   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
124   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
125   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
126   PetscInt       am=A->rmap.N,cm=C->rmap.N;
127   PetscInt       i,j,k,anzi,bnzi,cnzi,brow,nextb;
128   MatScalar      *aa=a->a,*ba=b->a,*baj,*ca=c->a;
129 
130   PetscFunctionBegin;
131   /* clean old values in C */
132   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
133   /* Traverse A row-wise. */
134   /* Build the ith row in C by summing over nonzero columns in A, */
135   /* the rows of B corresponding to nonzeros of A. */
136   for (i=0;i<am;i++) {
137     anzi = ai[i+1] - ai[i];
138     for (j=0;j<anzi;j++) {
139       brow = *aj++;
140       bnzi = bi[brow+1] - bi[brow];
141       bjj  = bj + bi[brow];
142       baj  = ba + bi[brow];
143       nextb = 0;
144       for (k=0; nextb<bnzi; k++) {
145         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
146           ca[k] += (*aa)*baj[nextb++];
147         }
148       }
149       flops += 2*bnzi;
150       aa++;
151     }
152     cnzi = ci[i+1] - ci[i];
153     ca += cnzi;
154     cj += cnzi;
155   }
156   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158 
159   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 
164 #undef __FUNCT__
165 #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ"
166 PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
167   PetscErrorCode ierr;
168 
169   PetscFunctionBegin;
170   if (scall == MAT_INITIAL_MATRIX){
171     ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
172   }
173   ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ"
179 PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
180 {
181   PetscErrorCode ierr;
182   Mat            At;
183   PetscInt       *ati,*atj;
184 
185   PetscFunctionBegin;
186   /* create symbolic At */
187   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
188   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap.n,A->rmap.n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
189 
190   /* get symbolic C=At*B */
191   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
192 
193   /* clean up */
194   ierr = MatDestroy(At);CHKERRQ(ierr);
195   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
196 
197   PetscFunctionReturn(0);
198 }
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ"
202 PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
203 {
204   PetscErrorCode ierr;
205   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
206   PetscInt       am=A->rmap.n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
207   PetscInt       cm=C->rmap.n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k,flops=0;
208   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
209 
210   PetscFunctionBegin;
211   /* clear old values in C */
212   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
213 
214   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
215   for (i=0;i<am;i++) {
216     bj   = b->j + bi[i];
217     ba   = b->a + bi[i];
218     bnzi = bi[i+1] - bi[i];
219     anzi = ai[i+1] - ai[i];
220     for (j=0; j<anzi; j++) {
221       nextb = 0;
222       crow  = *aj++;
223       cjj   = cj + ci[crow];
224       caj   = ca + ci[crow];
225       /* perform sparse axpy operation.  Note cjj includes bj. */
226       for (k=0; nextb<bnzi; k++) {
227         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
228           caj[k] += (*aa)*(*(ba+nextb));
229           nextb++;
230         }
231       }
232       flops += 2*bnzi;
233       aa++;
234     }
235   }
236 
237   /* Assemble the final matrix and clean up */
238   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
239   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
240   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243