xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 5df89d910620d6cb6d1f55622404bc9966ce1990)
1be1d678aSKris Buschelman 
2d50806bdSBarry Smith /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4d50806bdSBarry Smith           C = A * B
5d50806bdSBarry Smith */
6d50806bdSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
116284ec50SHong Zhang 
12ec01deb9SMatthew Knepley EXTERN_C_BEGIN
136284ec50SHong Zhang #undef __FUNCT__
145c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1538baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1638baddfdSBarry Smith {
17dfbe8321SBarry Smith   PetscErrorCode ierr;
185c66b693SKris Buschelman 
195c66b693SKris Buschelman   PetscFunctionBegin;
2026be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
2126be0446SHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
2226be0446SHong Zhang   }
2326be0446SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
245c66b693SKris Buschelman   PetscFunctionReturn(0);
255c66b693SKris Buschelman }
26ec01deb9SMatthew Knepley EXTERN_C_END
271c24bd37SHong Zhang 
2826be0446SHong Zhang #undef __FUNCT__
2926be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
30dfbe8321SBarry Smith PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
3158c24d83SHong Zhang {
32dfbe8321SBarry Smith   PetscErrorCode     ierr;
33a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
3458c24d83SHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
3538baddfdSBarry Smith   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci,*cj;
36d0f46423SBarry Smith   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
3738baddfdSBarry Smith   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,nspacedouble=0;
3858c24d83SHong Zhang   MatScalar          *ca;
39be0fcf8dSHong Zhang   PetscBT            lnkbt;
407212da7cSHong Zhang   PetscReal          afill;
4158c24d83SHong Zhang 
4258c24d83SHong Zhang   PetscFunctionBegin;
4358c24d83SHong Zhang   /* Set up */
4458c24d83SHong Zhang   /* Allocate ci array, arrays for fill computation and */
4558c24d83SHong Zhang   /* free space for accumulating nonzero column info */
4638baddfdSBarry Smith   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4758c24d83SHong Zhang   ci[0] = 0;
4858c24d83SHong Zhang 
49be0fcf8dSHong Zhang   /* create and initialize a linked list */
50be0fcf8dSHong Zhang   nlnk = bn+1;
51be0fcf8dSHong Zhang   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
5258c24d83SHong Zhang 
53c5db241fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
54a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5558c24d83SHong Zhang   current_space = free_space;
5658c24d83SHong Zhang 
5758c24d83SHong Zhang   /* Determine symbolic info for each row of the product: */
5858c24d83SHong Zhang   for (i=0;i<am;i++) {
5958c24d83SHong Zhang     anzi = ai[i+1] - ai[i];
6058c24d83SHong Zhang     cnzi = 0;
612d09714cSHong Zhang     j    = anzi;
622d09714cSHong Zhang     aj   = a->j + ai[i];
632d09714cSHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
642d09714cSHong Zhang       j--;
652d09714cSHong Zhang       brow = *(aj + j);
6658c24d83SHong Zhang       bnzj = bi[brow+1] - bi[brow];
6758c24d83SHong Zhang       bjj  = bj + bi[brow];
681c239cc6SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
69be0fcf8dSHong Zhang       ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
701c239cc6SHong Zhang       cnzi += nlnk;
7158c24d83SHong Zhang     }
7258c24d83SHong Zhang 
7358c24d83SHong Zhang     /* If free space is not available, make more free space */
7458c24d83SHong Zhang     /* Double the amount of total space in the list */
7558c24d83SHong Zhang     if (current_space->local_remaining<cnzi) {
764238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
77c5db241fSHong Zhang       nspacedouble++;
7858c24d83SHong Zhang     }
7958c24d83SHong Zhang 
80c5db241fSHong Zhang     /* Copy data into free space, then initialize lnk */
81be0fcf8dSHong Zhang     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
82c5db241fSHong Zhang     current_space->array           += cnzi;
8358c24d83SHong Zhang     current_space->local_used      += cnzi;
8458c24d83SHong Zhang     current_space->local_remaining -= cnzi;
8558c24d83SHong Zhang 
8658c24d83SHong Zhang     ci[i+1] = ci[i] + cnzi;
8758c24d83SHong Zhang   }
8858c24d83SHong Zhang 
8958c24d83SHong Zhang   /* Column indices are in the list of free space */
9058c24d83SHong Zhang   /* Allocate space for cj, initialize cj, and */
9158c24d83SHong Zhang   /* destroy list of free space and other temporary array(s) */
9238baddfdSBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
93a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
94be0fcf8dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
9558c24d83SHong Zhang 
9658c24d83SHong Zhang   /* Allocate space for ca */
9758c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
9858c24d83SHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
9958c24d83SHong Zhang 
10026be0446SHong Zhang   /* put together the new symbolic matrix */
1017adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
10258c24d83SHong Zhang 
10358c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
10458c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
10558c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
106e6b907acSBarry Smith   c->free_a   = PETSC_TRUE;
107e6b907acSBarry Smith   c->free_ij  = PETSC_TRUE;
10858c24d83SHong Zhang   c->nonew    = 0;
10958c24d83SHong Zhang 
1107212da7cSHong Zhang   /* set MatInfo */
1117212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
112f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1137212da7cSHong Zhang   c->maxnz                     = ci[am];
1147212da7cSHong Zhang   c->nz                        = ci[am];
1157212da7cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
1167212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1177212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1187212da7cSHong Zhang 
1197212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1207212da7cSHong Zhang   if (ci[am]) {
121f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
122f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
123f2b054eeSHong Zhang   } else {
124f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
125be0fcf8dSHong Zhang   }
126f2b054eeSHong Zhang #endif
12758c24d83SHong Zhang   PetscFunctionReturn(0);
12858c24d83SHong Zhang }
129d50806bdSBarry Smith 
1305c66b693SKris Buschelman 
13126be0446SHong Zhang #undef __FUNCT__
13226be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
133dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
134d50806bdSBarry Smith {
135dfbe8321SBarry Smith   PetscErrorCode ierr;
136d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
137d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
138d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
139d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
14038baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
141d0f46423SBarry Smith   PetscInt       am=A->rmap->N,cm=C->rmap->N;
142c124e916SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow,nextb;
143c124e916SHong Zhang   MatScalar      *aa=a->a,*ba=b->a,*baj,*ca=c->a;
144d50806bdSBarry Smith 
145d50806bdSBarry Smith   PetscFunctionBegin;
146c124e916SHong Zhang   /* clean old values in C */
147c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
148d50806bdSBarry Smith   /* Traverse A row-wise. */
149d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
150d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
151d50806bdSBarry Smith   for (i=0;i<am;i++) {
152d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
153d50806bdSBarry Smith     for (j=0;j<anzi;j++) {
154d50806bdSBarry Smith       brow = *aj++;
155d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
156d50806bdSBarry Smith       bjj  = bj + bi[brow];
157d50806bdSBarry Smith       baj  = ba + bi[brow];
158c124e916SHong Zhang       nextb = 0;
159c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
160c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
161c124e916SHong Zhang           ca[k] += (*aa)*baj[nextb++];
162c124e916SHong Zhang         }
163d50806bdSBarry Smith       }
164d50806bdSBarry Smith       flops += 2*bnzi;
165d50806bdSBarry Smith       aa++;
166d50806bdSBarry Smith     }
167d50806bdSBarry Smith     cnzi = ci[i+1] - ci[i];
168d50806bdSBarry Smith     ca += cnzi;
169d50806bdSBarry Smith     cj += cnzi;
170d50806bdSBarry Smith   }
171716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
172716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
173716bacf3SKris Buschelman 
174d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
175d50806bdSBarry Smith   PetscFunctionReturn(0);
176d50806bdSBarry Smith }
177bc011b1eSHong Zhang 
178bc011b1eSHong Zhang #undef __FUNCT__
179bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ"
180*5df89d91SHong Zhang PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
181*5df89d91SHong Zhang {
182bc011b1eSHong Zhang   PetscErrorCode ierr;
183bc011b1eSHong Zhang 
184bc011b1eSHong Zhang   PetscFunctionBegin;
185*5df89d91SHong Zhang   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Not done yet");
186bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
187bc011b1eSHong Zhang     ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
188bc011b1eSHong Zhang   }
189bc011b1eSHong Zhang   ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
190bc011b1eSHong Zhang   PetscFunctionReturn(0);
191bc011b1eSHong Zhang }
192bc011b1eSHong Zhang 
193bc011b1eSHong Zhang #undef __FUNCT__
194bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ"
195bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
196bc011b1eSHong Zhang {
197*5df89d91SHong Zhang   /*
198*5df89d91SHong Zhang   PetscErrorCode ierr;
199*5df89d91SHong Zhang   Mat            At;
200*5df89d91SHong Zhang   PetscInt       *ati,*atj;
201*5df89d91SHong Zhang    */
202*5df89d91SHong Zhang   PetscFunctionBegin;
203*5df89d91SHong Zhang 
204*5df89d91SHong Zhang   PetscFunctionReturn(0);
205*5df89d91SHong Zhang }
206*5df89d91SHong Zhang 
207*5df89d91SHong Zhang #undef __FUNCT__
208*5df89d91SHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ"
209*5df89d91SHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
210*5df89d91SHong Zhang {
211*5df89d91SHong Zhang   /*
212*5df89d91SHong Zhang   PetscErrorCode ierr;
213*5df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
214*5df89d91SHong Zhang   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
215*5df89d91SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
216*5df89d91SHong Zhang   PetscLogDouble flops=0.0;
217*5df89d91SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
218*5df89d91SHong Zhang    */
219*5df89d91SHong Zhang 
220*5df89d91SHong Zhang   PetscFunctionBegin;
221*5df89d91SHong Zhang 
222*5df89d91SHong Zhang   PetscFunctionReturn(0);
223*5df89d91SHong Zhang }
224*5df89d91SHong Zhang 
225*5df89d91SHong Zhang #undef __FUNCT__
226*5df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
227*5df89d91SHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
228*5df89d91SHong Zhang   PetscErrorCode ierr;
229*5df89d91SHong Zhang 
230*5df89d91SHong Zhang   PetscFunctionBegin;
231*5df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
232*5df89d91SHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
233*5df89d91SHong Zhang   }
234*5df89d91SHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
235*5df89d91SHong Zhang   PetscFunctionReturn(0);
236*5df89d91SHong Zhang }
237*5df89d91SHong Zhang 
238*5df89d91SHong Zhang #undef __FUNCT__
239*5df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
240*5df89d91SHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
241*5df89d91SHong Zhang {
242bc011b1eSHong Zhang   PetscErrorCode ierr;
243bc011b1eSHong Zhang   Mat            At;
24438baddfdSBarry Smith   PetscInt       *ati,*atj;
245bc011b1eSHong Zhang 
246bc011b1eSHong Zhang   PetscFunctionBegin;
247bc011b1eSHong Zhang   /* create symbolic At */
248bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
249d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
250bc011b1eSHong Zhang 
251bc011b1eSHong Zhang   /* get symbolic C=At*B */
252bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
253bc011b1eSHong Zhang 
254bc011b1eSHong Zhang   /* clean up */
2556bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
256bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
257bc011b1eSHong Zhang   PetscFunctionReturn(0);
258bc011b1eSHong Zhang }
259bc011b1eSHong Zhang 
260bc011b1eSHong Zhang #undef __FUNCT__
261*5df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
262*5df89d91SHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
263bc011b1eSHong Zhang {
264bc011b1eSHong Zhang   PetscErrorCode ierr;
2650fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
266d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
267d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
268d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
2690fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
270bc011b1eSHong Zhang 
271bc011b1eSHong Zhang   PetscFunctionBegin;
272bc011b1eSHong Zhang   /* clear old values in C */
273bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
274bc011b1eSHong Zhang 
275bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
276bc011b1eSHong Zhang   for (i=0;i<am;i++) {
277bc011b1eSHong Zhang     bj   = b->j + bi[i];
278bc011b1eSHong Zhang     ba   = b->a + bi[i];
279bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
280bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
281bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
282bc011b1eSHong Zhang       nextb = 0;
2830fbc74f4SHong Zhang       crow  = *aj++;
284bc011b1eSHong Zhang       cjj   = cj + ci[crow];
285bc011b1eSHong Zhang       caj   = ca + ci[crow];
286bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
287bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
2880fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
2890fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
290bc011b1eSHong Zhang           nextb++;
291bc011b1eSHong Zhang         }
292bc011b1eSHong Zhang       }
293bc011b1eSHong Zhang       flops += 2*bnzi;
2940fbc74f4SHong Zhang       aa++;
295bc011b1eSHong Zhang     }
296bc011b1eSHong Zhang   }
297bc011b1eSHong Zhang 
298bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
299bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
300bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
301bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
302bc011b1eSHong Zhang   PetscFunctionReturn(0);
303bc011b1eSHong Zhang }
3049123193aSHong Zhang 
305ec01deb9SMatthew Knepley EXTERN_C_BEGIN
3069123193aSHong Zhang #undef __FUNCT__
3079123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
3089123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3099123193aSHong Zhang {
3109123193aSHong Zhang   PetscErrorCode ierr;
3119123193aSHong Zhang 
3129123193aSHong Zhang   PetscFunctionBegin;
3139123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
3149123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
3159123193aSHong Zhang   }
3169123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
3179123193aSHong Zhang   PetscFunctionReturn(0);
3189123193aSHong Zhang }
319ec01deb9SMatthew Knepley EXTERN_C_END
320edd81eecSMatthew Knepley 
3219123193aSHong Zhang #undef __FUNCT__
3229123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
3239123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
3249123193aSHong Zhang {
3259123193aSHong Zhang   PetscErrorCode ierr;
3269123193aSHong Zhang 
3279123193aSHong Zhang   PetscFunctionBegin;
3285a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
3299123193aSHong Zhang   PetscFunctionReturn(0);
3309123193aSHong Zhang }
3319123193aSHong Zhang 
3329123193aSHong Zhang #undef __FUNCT__
3339123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
3349123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
3359123193aSHong Zhang {
336f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
3379123193aSHong Zhang   PetscErrorCode ierr;
338dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
339dd6ea824SBarry Smith   MatScalar      *aa;
340d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
341f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
3429123193aSHong Zhang 
3439123193aSHong Zhang   PetscFunctionBegin;
344f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
345e32f2f54SBarry Smith   if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm);
346e32f2f54SBarry Smith   if (A->rmap->n != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap->n,A->rmap->n);
347e32f2f54SBarry Smith   if (B->cmap->n != C->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap->n,C->cmap->n);
348f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
349f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
350f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
351f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
352f32d5d43SBarry Smith     colam = col*am;
353f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
354f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
355f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
356f32d5d43SBarry Smith       aj  = a->j + a->i[i];
357f32d5d43SBarry Smith       aa  = a->a + a->i[i];
358f32d5d43SBarry Smith       for (j=0; j<n; j++) {
359f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
360f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
361f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
362f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
3639123193aSHong Zhang       }
364f32d5d43SBarry Smith       c[colam + i]       = r1;
365f32d5d43SBarry Smith       c[colam + am + i]  = r2;
366f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
367f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
368f32d5d43SBarry Smith     }
369f32d5d43SBarry Smith     b1 += bm4;
370f32d5d43SBarry Smith     b2 += bm4;
371f32d5d43SBarry Smith     b3 += bm4;
372f32d5d43SBarry Smith     b4 += bm4;
373f32d5d43SBarry Smith   }
374f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
375f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
376f32d5d43SBarry Smith       r1 = 0.0;
377f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
378f32d5d43SBarry Smith       aj  = a->j + a->i[i];
379f32d5d43SBarry Smith       aa  = a->a + a->i[i];
380f32d5d43SBarry Smith 
381f32d5d43SBarry Smith       for (j=0; j<n; j++) {
382f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
383f32d5d43SBarry Smith       }
384f32d5d43SBarry Smith       c[col*am + i]     = r1;
385f32d5d43SBarry Smith     }
386f32d5d43SBarry Smith     b1 += bm;
387f32d5d43SBarry Smith   }
388dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
389f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
390f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
391f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
392f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
393f32d5d43SBarry Smith   PetscFunctionReturn(0);
394f32d5d43SBarry Smith }
395f32d5d43SBarry Smith 
396f32d5d43SBarry Smith /*
3974324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
398f32d5d43SBarry Smith */
399f32d5d43SBarry Smith #undef __FUNCT__
400f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
401f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
402f32d5d43SBarry Smith {
403f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
404f32d5d43SBarry Smith   PetscErrorCode ierr;
405dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
406dd6ea824SBarry Smith   MatScalar      *aa;
407d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
4084324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
409f32d5d43SBarry Smith 
410f32d5d43SBarry Smith   PetscFunctionBegin;
411f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
412f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
413f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
414f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
4154324174eSBarry Smith 
4164324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
4174324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
4184324174eSBarry Smith       colam = col*am;
4194324174eSBarry Smith       arm   = a->compressedrow.nrows;
4204324174eSBarry Smith       ii    = a->compressedrow.i;
4214324174eSBarry Smith       ridx  = a->compressedrow.rindex;
4224324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
4234324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
4244324174eSBarry Smith 	n   = ii[i+1] - ii[i];
4254324174eSBarry Smith 	aj  = a->j + ii[i];
4264324174eSBarry Smith 	aa  = a->a + ii[i];
4274324174eSBarry Smith 	for (j=0; j<n; j++) {
4284324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
4294324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
4304324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
4314324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
4324324174eSBarry Smith 	}
4334324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
4344324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
4354324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
4364324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
4374324174eSBarry Smith       }
4384324174eSBarry Smith       b1 += bm4;
4394324174eSBarry Smith       b2 += bm4;
4404324174eSBarry Smith       b3 += bm4;
4414324174eSBarry Smith       b4 += bm4;
4424324174eSBarry Smith     }
4434324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
4444324174eSBarry Smith       colam = col*am;
4454324174eSBarry Smith       arm   = a->compressedrow.nrows;
4464324174eSBarry Smith       ii    = a->compressedrow.i;
4474324174eSBarry Smith       ridx  = a->compressedrow.rindex;
4484324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
4494324174eSBarry Smith 	r1 = 0.0;
4504324174eSBarry Smith 	n   = ii[i+1] - ii[i];
4514324174eSBarry Smith 	aj  = a->j + ii[i];
4524324174eSBarry Smith 	aa  = a->a + ii[i];
4534324174eSBarry Smith 
4544324174eSBarry Smith 	for (j=0; j<n; j++) {
4554324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
4564324174eSBarry Smith 	}
4574324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
4584324174eSBarry Smith       }
4594324174eSBarry Smith       b1 += bm;
4604324174eSBarry Smith     }
4614324174eSBarry Smith   } else {
462f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
463f32d5d43SBarry Smith       colam = col*am;
464f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
465f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
466f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
467f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
468f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
469f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
470f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
471f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
472f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
473f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
474f32d5d43SBarry Smith 	}
475f32d5d43SBarry Smith 	c[colam + i]       += r1;
476f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
477f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
478f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
479f32d5d43SBarry Smith       }
480f32d5d43SBarry Smith       b1 += bm4;
481f32d5d43SBarry Smith       b2 += bm4;
482f32d5d43SBarry Smith       b3 += bm4;
483f32d5d43SBarry Smith       b4 += bm4;
484f32d5d43SBarry Smith     }
485f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
486f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
487f32d5d43SBarry Smith 	r1 = 0.0;
488f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
489f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
490f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
491f32d5d43SBarry Smith 
492f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
493f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
494f32d5d43SBarry Smith 	}
495f32d5d43SBarry Smith 	c[col*am + i]     += r1;
496f32d5d43SBarry Smith       }
497f32d5d43SBarry Smith       b1 += bm;
498f32d5d43SBarry Smith     }
4994324174eSBarry Smith   }
500dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
501f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
502f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
5039123193aSHong Zhang   PetscFunctionReturn(0);
5049123193aSHong Zhang }
505