xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 6bf464f92cc51e6fd6163850774a6badb2f63b6b)
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 
179bc011b1eSHong Zhang #undef __FUNCT__
180bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ"
181bc011b1eSHong Zhang PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
182bc011b1eSHong Zhang   PetscErrorCode ierr;
183bc011b1eSHong Zhang 
184bc011b1eSHong Zhang   PetscFunctionBegin;
185bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
186bc011b1eSHong Zhang     ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
187bc011b1eSHong Zhang   }
188bc011b1eSHong Zhang   ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
189bc011b1eSHong Zhang   PetscFunctionReturn(0);
190bc011b1eSHong Zhang }
191bc011b1eSHong Zhang 
192bc011b1eSHong Zhang #undef __FUNCT__
193bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ"
194bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
195bc011b1eSHong Zhang {
196bc011b1eSHong Zhang   PetscErrorCode ierr;
197bc011b1eSHong Zhang   Mat            At;
19838baddfdSBarry Smith   PetscInt       *ati,*atj;
199bc011b1eSHong Zhang 
200bc011b1eSHong Zhang   PetscFunctionBegin;
201bc011b1eSHong Zhang   /* create symbolic At */
202bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
203d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
204bc011b1eSHong Zhang 
205bc011b1eSHong Zhang   /* get symbolic C=At*B */
206bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
207bc011b1eSHong Zhang 
208bc011b1eSHong Zhang   /* clean up */
209*6bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
210bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
211bc011b1eSHong Zhang   PetscFunctionReturn(0);
212bc011b1eSHong Zhang }
213bc011b1eSHong Zhang 
214bc011b1eSHong Zhang #undef __FUNCT__
215bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ"
216bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
217bc011b1eSHong Zhang {
218bc011b1eSHong Zhang   PetscErrorCode ierr;
2190fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
220d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
221d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
222d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
2230fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
224bc011b1eSHong Zhang 
225bc011b1eSHong Zhang   PetscFunctionBegin;
226bc011b1eSHong Zhang   /* clear old values in C */
227bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
228bc011b1eSHong Zhang 
229bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
230bc011b1eSHong Zhang   for (i=0;i<am;i++) {
231bc011b1eSHong Zhang     bj   = b->j + bi[i];
232bc011b1eSHong Zhang     ba   = b->a + bi[i];
233bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
234bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
235bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
236bc011b1eSHong Zhang       nextb = 0;
2370fbc74f4SHong Zhang       crow  = *aj++;
238bc011b1eSHong Zhang       cjj   = cj + ci[crow];
239bc011b1eSHong Zhang       caj   = ca + ci[crow];
240bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
241bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
2420fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
2430fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
244bc011b1eSHong Zhang           nextb++;
245bc011b1eSHong Zhang         }
246bc011b1eSHong Zhang       }
247bc011b1eSHong Zhang       flops += 2*bnzi;
2480fbc74f4SHong Zhang       aa++;
249bc011b1eSHong Zhang     }
250bc011b1eSHong Zhang   }
251bc011b1eSHong Zhang 
252bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
253bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
254bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
256bc011b1eSHong Zhang   PetscFunctionReturn(0);
257bc011b1eSHong Zhang }
2589123193aSHong Zhang 
259ec01deb9SMatthew Knepley EXTERN_C_BEGIN
2609123193aSHong Zhang #undef __FUNCT__
2619123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
2629123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2639123193aSHong Zhang {
2649123193aSHong Zhang   PetscErrorCode ierr;
2659123193aSHong Zhang 
2669123193aSHong Zhang   PetscFunctionBegin;
2679123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
2689123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
2699123193aSHong Zhang   }
2709123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
2719123193aSHong Zhang   PetscFunctionReturn(0);
2729123193aSHong Zhang }
273ec01deb9SMatthew Knepley EXTERN_C_END
274edd81eecSMatthew Knepley 
2759123193aSHong Zhang #undef __FUNCT__
2769123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
2779123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2789123193aSHong Zhang {
2799123193aSHong Zhang   PetscErrorCode ierr;
2809123193aSHong Zhang 
2819123193aSHong Zhang   PetscFunctionBegin;
2825a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
2839123193aSHong Zhang   PetscFunctionReturn(0);
2849123193aSHong Zhang }
2859123193aSHong Zhang 
2869123193aSHong Zhang #undef __FUNCT__
2879123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
2889123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
2899123193aSHong Zhang {
290f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
2919123193aSHong Zhang   PetscErrorCode ierr;
292dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
293dd6ea824SBarry Smith   MatScalar      *aa;
294d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
295f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
2969123193aSHong Zhang 
2979123193aSHong Zhang   PetscFunctionBegin;
298f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
299e32f2f54SBarry 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);
300e32f2f54SBarry 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);
301e32f2f54SBarry 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);
302f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
303f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
304f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
305f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
306f32d5d43SBarry Smith     colam = col*am;
307f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
308f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
309f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
310f32d5d43SBarry Smith       aj  = a->j + a->i[i];
311f32d5d43SBarry Smith       aa  = a->a + a->i[i];
312f32d5d43SBarry Smith       for (j=0; j<n; j++) {
313f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
314f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
315f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
316f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
3179123193aSHong Zhang       }
318f32d5d43SBarry Smith       c[colam + i]       = r1;
319f32d5d43SBarry Smith       c[colam + am + i]  = r2;
320f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
321f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
322f32d5d43SBarry Smith     }
323f32d5d43SBarry Smith     b1 += bm4;
324f32d5d43SBarry Smith     b2 += bm4;
325f32d5d43SBarry Smith     b3 += bm4;
326f32d5d43SBarry Smith     b4 += bm4;
327f32d5d43SBarry Smith   }
328f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
329f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
330f32d5d43SBarry Smith       r1 = 0.0;
331f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
332f32d5d43SBarry Smith       aj  = a->j + a->i[i];
333f32d5d43SBarry Smith       aa  = a->a + a->i[i];
334f32d5d43SBarry Smith 
335f32d5d43SBarry Smith       for (j=0; j<n; j++) {
336f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
337f32d5d43SBarry Smith       }
338f32d5d43SBarry Smith       c[col*am + i]     = r1;
339f32d5d43SBarry Smith     }
340f32d5d43SBarry Smith     b1 += bm;
341f32d5d43SBarry Smith   }
342dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
343f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
344f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
345f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
346f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
347f32d5d43SBarry Smith   PetscFunctionReturn(0);
348f32d5d43SBarry Smith }
349f32d5d43SBarry Smith 
350f32d5d43SBarry Smith /*
3514324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
352f32d5d43SBarry Smith */
353f32d5d43SBarry Smith #undef __FUNCT__
354f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
355f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
356f32d5d43SBarry Smith {
357f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
358f32d5d43SBarry Smith   PetscErrorCode ierr;
359dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
360dd6ea824SBarry Smith   MatScalar      *aa;
361d0f46423SBarry 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;
3624324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
363f32d5d43SBarry Smith 
364f32d5d43SBarry Smith   PetscFunctionBegin;
365f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
366f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
367f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
368f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
3694324174eSBarry Smith 
3704324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
3714324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
3724324174eSBarry Smith       colam = col*am;
3734324174eSBarry Smith       arm   = a->compressedrow.nrows;
3744324174eSBarry Smith       ii    = a->compressedrow.i;
3754324174eSBarry Smith       ridx  = a->compressedrow.rindex;
3764324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
3774324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
3784324174eSBarry Smith 	n   = ii[i+1] - ii[i];
3794324174eSBarry Smith 	aj  = a->j + ii[i];
3804324174eSBarry Smith 	aa  = a->a + ii[i];
3814324174eSBarry Smith 	for (j=0; j<n; j++) {
3824324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
3834324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
3844324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
3854324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
3864324174eSBarry Smith 	}
3874324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
3884324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
3894324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
3904324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
3914324174eSBarry Smith       }
3924324174eSBarry Smith       b1 += bm4;
3934324174eSBarry Smith       b2 += bm4;
3944324174eSBarry Smith       b3 += bm4;
3954324174eSBarry Smith       b4 += bm4;
3964324174eSBarry Smith     }
3974324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
3984324174eSBarry Smith       colam = col*am;
3994324174eSBarry Smith       arm   = a->compressedrow.nrows;
4004324174eSBarry Smith       ii    = a->compressedrow.i;
4014324174eSBarry Smith       ridx  = a->compressedrow.rindex;
4024324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
4034324174eSBarry Smith 	r1 = 0.0;
4044324174eSBarry Smith 	n   = ii[i+1] - ii[i];
4054324174eSBarry Smith 	aj  = a->j + ii[i];
4064324174eSBarry Smith 	aa  = a->a + ii[i];
4074324174eSBarry Smith 
4084324174eSBarry Smith 	for (j=0; j<n; j++) {
4094324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
4104324174eSBarry Smith 	}
4114324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
4124324174eSBarry Smith       }
4134324174eSBarry Smith       b1 += bm;
4144324174eSBarry Smith     }
4154324174eSBarry Smith   } else {
416f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
417f32d5d43SBarry Smith       colam = col*am;
418f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
419f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
420f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
421f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
422f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
423f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
424f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
425f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
426f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
427f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
428f32d5d43SBarry Smith 	}
429f32d5d43SBarry Smith 	c[colam + i]       += r1;
430f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
431f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
432f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
433f32d5d43SBarry Smith       }
434f32d5d43SBarry Smith       b1 += bm4;
435f32d5d43SBarry Smith       b2 += bm4;
436f32d5d43SBarry Smith       b3 += bm4;
437f32d5d43SBarry Smith       b4 += bm4;
438f32d5d43SBarry Smith     }
439f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
440f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
441f32d5d43SBarry Smith 	r1 = 0.0;
442f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
443f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
444f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
445f32d5d43SBarry Smith 
446f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
447f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
448f32d5d43SBarry Smith 	}
449f32d5d43SBarry Smith 	c[col*am + i]     += r1;
450f32d5d43SBarry Smith       }
451f32d5d43SBarry Smith       b1 += bm;
452f32d5d43SBarry Smith     }
4534324174eSBarry Smith   }
454dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
455f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
456f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
4579123193aSHong Zhang   PetscFunctionReturn(0);
4589123193aSHong Zhang }
459