xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision dadf0e6bf3fbf8619c99f30f0124ab4e08bd7289)
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){
21598bc09dSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2226be0446SHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
23598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2426be0446SHong Zhang   }
25598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
2601e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
27598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
285c66b693SKris Buschelman   PetscFunctionReturn(0);
295c66b693SKris Buschelman }
30ec01deb9SMatthew Knepley EXTERN_C_END
311c24bd37SHong Zhang 
3226be0446SHong Zhang #undef __FUNCT__
33b561aa9dSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ"
34b561aa9dSHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
3558c24d83SHong Zhang {
36dfbe8321SBarry Smith   PetscErrorCode     ierr;
37a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
38b561aa9dSHong Zhang   PetscInt           *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj;
39b561aa9dSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0;
40be0fcf8dSHong Zhang   PetscBT            lnkbt;
4158c24d83SHong Zhang 
4258c24d83SHong Zhang   PetscFunctionBegin;
4358c24d83SHong Zhang   /* Allocate ci array, arrays for fill computation and */
4458c24d83SHong Zhang   /* free space for accumulating nonzero column info */
4538baddfdSBarry Smith   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4658c24d83SHong Zhang   ci[0] = 0;
4758c24d83SHong Zhang 
48be0fcf8dSHong Zhang   /* create and initialize a linked list */
49be0fcf8dSHong Zhang   nlnk = bn+1;
50be0fcf8dSHong Zhang   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
5158c24d83SHong Zhang 
52c5db241fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
53a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5458c24d83SHong Zhang   current_space = free_space;
5558c24d83SHong Zhang 
5658c24d83SHong Zhang   /* Determine symbolic info for each row of the product: */
5758c24d83SHong Zhang   for (i=0; i<am; i++) {
5858c24d83SHong Zhang     anzi = ai[i+1] - ai[i];
5958c24d83SHong Zhang     cnzi = 0;
60b561aa9dSHong Zhang     aj   = Aj + ai[i];
6177584fe6SHong Zhang     for (j=0; j<anzi; j++){
62b561aa9dSHong Zhang       brow = aj[j];
6358c24d83SHong Zhang       bnzj = bi[brow+1] - bi[brow];
6458c24d83SHong Zhang       bjj  = bj + bi[brow];
651c239cc6SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
66dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
671c239cc6SHong Zhang       cnzi += nlnk;
6858c24d83SHong Zhang     }
6958c24d83SHong Zhang 
7058c24d83SHong Zhang     /* If free space is not available, make more free space */
7158c24d83SHong Zhang     /* Double the amount of total space in the list */
7258c24d83SHong Zhang     if (current_space->local_remaining<cnzi) {
734238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
74b561aa9dSHong Zhang       ndouble++;
7558c24d83SHong Zhang     }
7658c24d83SHong Zhang 
77c5db241fSHong Zhang     /* Copy data into free space, then initialize lnk */
78be0fcf8dSHong Zhang     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
79c5db241fSHong Zhang     current_space->array           += cnzi;
8058c24d83SHong Zhang     current_space->local_used      += cnzi;
8158c24d83SHong Zhang     current_space->local_remaining -= cnzi;
8258c24d83SHong Zhang     ci[i+1] = ci[i] + cnzi;
8358c24d83SHong Zhang   }
8458c24d83SHong Zhang 
8558c24d83SHong Zhang   /* Column indices are in the list of free space */
8658c24d83SHong Zhang   /* Allocate space for cj, initialize cj, and */
8758c24d83SHong Zhang   /* destroy list of free space and other temporary array(s) */
8838baddfdSBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
89a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
90be0fcf8dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
9158c24d83SHong Zhang 
92b561aa9dSHong Zhang   *Ci           = ci;
93b561aa9dSHong Zhang   *Cj           = cj;
94b561aa9dSHong Zhang   *nspacedouble = ndouble;
95b561aa9dSHong Zhang   PetscFunctionReturn(0);
96b561aa9dSHong Zhang }
97b561aa9dSHong Zhang 
98b561aa9dSHong Zhang #undef __FUNCT__
99b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
100b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
101b561aa9dSHong Zhang {
102b561aa9dSHong Zhang   PetscErrorCode ierr;
103b561aa9dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
104b561aa9dSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj;
105b561aa9dSHong Zhang   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
106b561aa9dSHong Zhang   MatScalar      *ca;
107b561aa9dSHong Zhang   PetscReal      afill;
10836ec6d2dSHong Zhang   PetscBool      dense_axpy; /* false: use sparse axpy; otherwise use dense axpy in MatMatMultNumeric_SeqAIJ_SeqAIJ() */
109b561aa9dSHong Zhang 
110b561aa9dSHong Zhang   PetscFunctionBegin;
111b561aa9dSHong Zhang   /* Get ci and cj */
112b561aa9dSHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
113b561aa9dSHong Zhang 
11458c24d83SHong Zhang   /* Allocate space for ca */
11558c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
11658c24d83SHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
11758c24d83SHong Zhang 
11826be0446SHong Zhang   /* put together the new symbolic matrix */
1197adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
12058c24d83SHong Zhang 
12158c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
12258c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
12358c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
124e6b907acSBarry Smith   c->free_a   = PETSC_TRUE;
125e6b907acSBarry Smith   c->free_ij  = PETSC_TRUE;
12658c24d83SHong Zhang   c->nonew    = 0;
1278cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqAIJ;
12858c24d83SHong Zhang 
12953565b12SHong Zhang   /* Determine which MatMatMultNumeric_SeqAIJ_SeqAIJ() to be used */
13036ec6d2dSHong Zhang   dense_axpy = PETSC_TRUE;
13136ec6d2dSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_denseaxpy",&dense_axpy,PETSC_NULL);CHKERRQ(ierr);
13236ec6d2dSHong Zhang   if (dense_axpy){
13336ec6d2dSHong Zhang     ierr = PetscMalloc(bn*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr);
134c58ca2e3SHong Zhang     ierr = PetscMemzero(c->matmult_abdense,dense_axpy*bn*sizeof(PetscScalar));CHKERRQ(ierr);
135c58ca2e3SHong Zhang     (*C)->ops->matmultnumeric =  MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, takes additional dense_axpy*bn*sizeof(PetscScalar) space */
136c58ca2e3SHong Zhang   } else { /* slower, but use less memory */
137c58ca2e3SHong Zhang     (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy; /* slower, less memory */
138c58ca2e3SHong Zhang   }
1390b7e3e3dSHong Zhang 
1407212da7cSHong Zhang   /* set MatInfo */
1417212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
142f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1437212da7cSHong Zhang   c->maxnz                     = ci[am];
1447212da7cSHong Zhang   c->nz                        = ci[am];
1457212da7cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
1467212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1477212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1487212da7cSHong Zhang 
1497212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1507212da7cSHong Zhang   if (ci[am]) {
151f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
152f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
153f2b054eeSHong Zhang   } else {
154f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
155be0fcf8dSHong Zhang   }
156f2b054eeSHong Zhang #endif
15758c24d83SHong Zhang   PetscFunctionReturn(0);
15858c24d83SHong Zhang }
159d50806bdSBarry Smith 
16026be0446SHong Zhang #undef __FUNCT__
16126be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
162dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
163d50806bdSBarry Smith {
164dfbe8321SBarry Smith   PetscErrorCode ierr;
165d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
166d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
167d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
168d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
16938baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
170c58ca2e3SHong Zhang   PetscInt       am=A->rmap->n,cm=C->rmap->n;
171519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
17236ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
1730b7e3e3dSHong Zhang   PetscScalar    *ab_dense=c->matmult_abdense;
174d50806bdSBarry Smith 
175d50806bdSBarry Smith   PetscFunctionBegin;
176c124e916SHong Zhang   /* clean old values in C */
177c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
178d50806bdSBarry Smith   /* Traverse A row-wise. */
179d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
180d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
181d50806bdSBarry Smith   for (i=0; i<am; i++) {
182d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
183d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
184519eb980SHong Zhang       brow = aj[j];
185d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
186d50806bdSBarry Smith       bjj  = bj + bi[brow];
187d50806bdSBarry Smith       baj  = ba + bi[brow];
188519eb980SHong Zhang       /* perform dense axpy */
18936ec6d2dSHong Zhang       valtmp = aa[j];
190519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
19136ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
192519eb980SHong Zhang       }
193519eb980SHong Zhang       flops += 2*bnzi;
194519eb980SHong Zhang     }
195c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
196c58ca2e3SHong Zhang 
197c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
198519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
199519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
200519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
201519eb980SHong Zhang     }
202519eb980SHong Zhang     flops += cnzi;
203519eb980SHong Zhang     cj += cnzi; ca += cnzi;
204519eb980SHong Zhang   }
205c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
206c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
207c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
208c58ca2e3SHong Zhang   C->ops->matmultnumeric =  MatMatMultNumeric_SeqAIJ_SeqAIJ;
209c58ca2e3SHong Zhang   PetscFunctionReturn(0);
210c58ca2e3SHong Zhang }
211c58ca2e3SHong Zhang 
212c58ca2e3SHong Zhang #undef __FUNCT__
213c58ca2e3SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
214c58ca2e3SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat B,Mat C)
215c58ca2e3SHong Zhang {
216c58ca2e3SHong Zhang   PetscErrorCode ierr;
217c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
218c58ca2e3SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
219c58ca2e3SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
220c58ca2e3SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
221c58ca2e3SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
222c58ca2e3SHong Zhang   PetscInt       am=A->rmap->N,cm=C->rmap->N;
223c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
22436ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
225c58ca2e3SHong Zhang   PetscInt       nextb;
226c58ca2e3SHong Zhang 
227c58ca2e3SHong Zhang   PetscFunctionBegin;
228c58ca2e3SHong Zhang   /* clean old values in C */
229c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
230c58ca2e3SHong Zhang   /* Traverse A row-wise. */
231c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
232c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
233519eb980SHong Zhang   for (i=0;i<am;i++) {
234519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
235519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
236519eb980SHong Zhang     for (j=0;j<anzi;j++) {
237519eb980SHong Zhang       brow = aj[j];
238519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
239519eb980SHong Zhang       bjj  = bj + bi[brow];
240519eb980SHong Zhang       baj  = ba + bi[brow];
241519eb980SHong Zhang       /* perform sparse axpy */
24236ec6d2dSHong Zhang       valtmp = aa[j];
243c124e916SHong Zhang       nextb  = 0;
244c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
245c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
24636ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
247c124e916SHong Zhang         }
248d50806bdSBarry Smith       }
249d50806bdSBarry Smith       flops += 2*bnzi;
250d50806bdSBarry Smith     }
251519eb980SHong Zhang     aj += anzi; aa += anzi;
252519eb980SHong Zhang     cj += cnzi; ca += cnzi;
253d50806bdSBarry Smith   }
254c58ca2e3SHong Zhang 
255716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
258d50806bdSBarry Smith   PetscFunctionReturn(0);
259d50806bdSBarry Smith }
260bc011b1eSHong Zhang 
261d2853540SHong Zhang /* This routine is not used. Should be removed! */
262bc011b1eSHong Zhang #undef __FUNCT__
2636fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
2646fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2655df89d91SHong Zhang {
266bc011b1eSHong Zhang   PetscErrorCode ierr;
267bc011b1eSHong Zhang 
268bc011b1eSHong Zhang   PetscFunctionBegin;
269bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
2706fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
271bc011b1eSHong Zhang   }
2726fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
273bc011b1eSHong Zhang   PetscFunctionReturn(0);
274bc011b1eSHong Zhang }
275bc011b1eSHong Zhang 
276bc011b1eSHong Zhang #undef __FUNCT__
277e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
278e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
2792128a86cSHong Zhang {
2802128a86cSHong Zhang   PetscErrorCode      ierr;
281e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
2822128a86cSHong Zhang 
2832128a86cSHong Zhang   PetscFunctionBegin;
284b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
2852128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
2862128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
2872128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
2882128a86cSHong Zhang   PetscFunctionReturn(0);
2892128a86cSHong Zhang }
2902128a86cSHong Zhang 
2912128a86cSHong Zhang #undef __FUNCT__
2922128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
2932128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
2942128a86cSHong Zhang {
2952128a86cSHong Zhang   PetscErrorCode      ierr;
2962128a86cSHong Zhang   PetscContainer      container;
297e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
2982128a86cSHong Zhang 
2992128a86cSHong Zhang   PetscFunctionBegin;
300e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
3012128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
3022128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
3032128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
3042128a86cSHong Zhang   if (A->ops->destroy) {
3052128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
3062128a86cSHong Zhang   }
307e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
3082128a86cSHong Zhang   PetscFunctionReturn(0);
3092128a86cSHong Zhang }
3102128a86cSHong Zhang 
3112128a86cSHong Zhang #undef __FUNCT__
3126fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
3136fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
314bc011b1eSHong Zhang {
3155df89d91SHong Zhang   PetscErrorCode      ierr;
31681d82fe4SHong Zhang   Mat                 Bt;
31781d82fe4SHong Zhang   PetscInt            *bti,*btj;
318e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
3192128a86cSHong Zhang   PetscContainer      container;
320d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
321d2853540SHong Zhang 
32281d82fe4SHong Zhang   PetscFunctionBegin;
323d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
32481d82fe4SHong Zhang    /* create symbolic Bt */
32581d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
32681d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
32781d82fe4SHong Zhang 
32881d82fe4SHong Zhang   /* get symbolic C=A*Bt */
32981d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
33081d82fe4SHong Zhang 
3312128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
332e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
3332128a86cSHong Zhang 
3342128a86cSHong Zhang   /* attach the supporting struct to C */
3352128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
3362128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
337e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
338e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
3392128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
3402128a86cSHong Zhang 
3412128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
3422128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
3432128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
3442128a86cSHong Zhang 
345d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
346d2853540SHong Zhang   etime2 += tf - t0;
347d2853540SHong Zhang 
348f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
3492128a86cSHong Zhang   if (multtrans->usecoloring){
350b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
351b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
3522128a86cSHong Zhang     ISColoring           iscoloring;
3532128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
354d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
355e8354b3eSHong Zhang 
356e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
357502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
358d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
359d2853540SHong Zhang     etime0 += tf - t0;
360d2853540SHong Zhang 
361d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
362b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
3632128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
3642128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
365e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
366d2853540SHong Zhang     etime01 += tf - t0;
3672128a86cSHong Zhang 
368e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
3692128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
3702128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
3712128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
3722128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
3732128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
3742128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
3752128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
3762128a86cSHong Zhang 
3772128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
3782128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
3792128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
3802128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
3812128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
3822128a86cSHong Zhang     multtrans->ABt_den = C_dense;
383e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
384e8354b3eSHong Zhang     etime1 += tf - t0;
385f94ccd6cSHong Zhang 
386f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
3871ce71dffSSatish Balay     {
388f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
389f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));
390f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
3911ce71dffSSatish Balay     }
392f94ccd6cSHong Zhang #endif
3932128a86cSHong Zhang   }
394*e99be685SHong Zhang   /* clean up */
395*e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
396*e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
3972128a86cSHong Zhang 
398f94ccd6cSHong Zhang 
399f94ccd6cSHong Zhang 
40081d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
40181d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
4021fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
4031fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
4041fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
4051fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
4061fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
4071fa1209cSHong Zhang   MatScalar          *ca;
4081fa1209cSHong Zhang   PetscBT            lnkbt;
4091fa1209cSHong Zhang   PetscReal          afill;
4105df89d91SHong Zhang 
4111fa1209cSHong Zhang   /* Allocate row pointer array ci  */
4121fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4131fa1209cSHong Zhang   ci[0] = 0;
4141fa1209cSHong Zhang 
4151fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
4161fa1209cSHong Zhang   nlnk = bm+1;
4171fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
4181fa1209cSHong Zhang 
4191fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
4201fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
4211fa1209cSHong Zhang   current_space = free_space;
4221fa1209cSHong Zhang 
4231fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
4241fa1209cSHong Zhang   for (i=0; i<am; i++) {
4251fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
4261fa1209cSHong Zhang     cnzi = 0;
4271fa1209cSHong Zhang     acol = aj + ai[i];
4281fa1209cSHong Zhang     for (j=0; j<bm; j++){
4291fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
4301fa1209cSHong Zhang       bcol= bj + bi[j];
43181d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
4321fa1209cSHong Zhang       ka = 0; kb = 0;
4331fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
43481d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
43581d82fe4SHong Zhang         if (ka == anzi) break;
43681d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
43781d82fe4SHong Zhang         if (kb == bnzj) break;
4381fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
4391fa1209cSHong Zhang           index[0] = j;
4401fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
4411fa1209cSHong Zhang           cnzi++;
4421fa1209cSHong Zhang           break;
4431fa1209cSHong Zhang         }
4441fa1209cSHong Zhang       }
4451fa1209cSHong Zhang     }
4461fa1209cSHong Zhang 
4471fa1209cSHong Zhang     /* If free space is not available, make more free space */
4481fa1209cSHong Zhang     /* Double the amount of total space in the list */
4491fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
4501fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
4511fa1209cSHong Zhang       nspacedouble++;
4521fa1209cSHong Zhang     }
4531fa1209cSHong Zhang 
4541fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
4551fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
4561fa1209cSHong Zhang     current_space->array           += cnzi;
4571fa1209cSHong Zhang     current_space->local_used      += cnzi;
4581fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
4591fa1209cSHong Zhang 
4601fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
4611fa1209cSHong Zhang   }
4621fa1209cSHong Zhang 
4631fa1209cSHong Zhang 
4641fa1209cSHong Zhang   /* Column indices are in the list of free space.
4651fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
4661fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
4671fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
4681fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
4691fa1209cSHong Zhang 
4701fa1209cSHong Zhang   /* Allocate space for ca */
4711fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4721fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4731fa1209cSHong Zhang 
4741fa1209cSHong Zhang   /* put together the new symbolic matrix */
4751fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
4761fa1209cSHong Zhang 
4771fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4781fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
4791fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
4801fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
4811fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
4821fa1209cSHong Zhang   c->nonew    = 0;
4831fa1209cSHong Zhang 
4841fa1209cSHong Zhang   /* set MatInfo */
4851fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
4861fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
4871fa1209cSHong Zhang   c->maxnz                     = ci[am];
4881fa1209cSHong Zhang   c->nz                        = ci[am];
4891fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
4901fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
4911fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
4921fa1209cSHong Zhang 
4931fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
4941fa1209cSHong Zhang   if (ci[am]) {
4951fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
4966fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
4971fa1209cSHong Zhang   } else {
4981fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
4991fa1209cSHong Zhang   }
5001fa1209cSHong Zhang #endif
50181d82fe4SHong Zhang #endif
5025df89d91SHong Zhang   PetscFunctionReturn(0);
5035df89d91SHong Zhang }
5045df89d91SHong Zhang 
5056973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
5065df89d91SHong Zhang #undef __FUNCT__
5076fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
5086fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
5095df89d91SHong Zhang {
5105df89d91SHong Zhang   PetscErrorCode ierr;
5115df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
512e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
51381d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
5145df89d91SHong Zhang   PetscLogDouble flops=0.0;
5156973769bSHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval;
516e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
5172128a86cSHong Zhang   PetscContainer      container;
5186973769bSHong Zhang #if defined(USE_ARRAY)
5196973769bSHong Zhang   MatScalar      *spdot;
5206973769bSHong Zhang #endif
5215df89d91SHong Zhang 
5225df89d91SHong Zhang   PetscFunctionBegin;
523e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
5242128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
5252128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
5262128a86cSHong Zhang   if (multtrans->usecoloring){
527b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
528c8db22f6SHong Zhang     Mat                   Bt_dense;
529c8db22f6SHong Zhang     PetscInt              m,n;
530b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
531a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
532a0b95ffbSSatish Balay 
5332128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
534c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
535c8db22f6SHong Zhang 
536b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
5372128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
538b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
5392128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
540b2d2b04fSHong Zhang     etime0 += tf - t0;
541c8db22f6SHong Zhang 
542c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
5432128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
5442128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
5452128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
5462128a86cSHong Zhang     etime2 += tf - t0;
547c8db22f6SHong Zhang 
5482128a86cSHong Zhang     /* Recover C from C_dense */
5492128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
550b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
5512128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
5522128a86cSHong Zhang     etime1 += tf - t0;
553f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
554f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
555f94ccd6cSHong Zhang #endif
556c8db22f6SHong Zhang     PetscFunctionReturn(0);
557c8db22f6SHong Zhang   }
558c8db22f6SHong Zhang 
5596973769bSHong Zhang #if defined(USE_ARRAY)
5606973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
561e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
562e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
5636973769bSHong Zhang #endif
5646973769bSHong Zhang 
56581d82fe4SHong Zhang   /* clear old values in C */
56681d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
5675df89d91SHong Zhang 
5681fa1209cSHong Zhang   for (i=0; i<cm; i++) {
56981d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
57081d82fe4SHong Zhang     acol = aj + ai[i];
5716973769bSHong Zhang     aval = aa + ai[i];
5721fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
5731fa1209cSHong Zhang     ccol = cj + ci[i];
5746973769bSHong Zhang     cval = ca + ci[i];
5751fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
57681d82fe4SHong Zhang       brow = ccol[j];
57781d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
57881d82fe4SHong Zhang       bcol = bj + bi[brow];
5796973769bSHong Zhang       bval = ba + bi[brow];
5806973769bSHong Zhang 
58181d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
5826973769bSHong Zhang #if defined(USE_ARRAY)
5836973769bSHong Zhang       /* put ba to spdot array */
5846973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
5856973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
5866973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
5876973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
5886973769bSHong Zhang       }
5896973769bSHong Zhang       /* zero spdot array */
5906973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
5916973769bSHong Zhang #else
59281d82fe4SHong Zhang       nexta = 0; nextb = 0;
59381d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
59481d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
59581d82fe4SHong Zhang         if (nexta == anzi) break;
59681d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
59781d82fe4SHong Zhang         if (nextb == bnzj) break;
59881d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
5996973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
60081d82fe4SHong Zhang           nexta++; nextb++;
60181d82fe4SHong Zhang           flops += 2;
6021fa1209cSHong Zhang         }
6031fa1209cSHong Zhang       }
6046973769bSHong Zhang #endif
60581d82fe4SHong Zhang     }
60681d82fe4SHong Zhang   }
60781d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60881d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60981d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
6106973769bSHong Zhang #if defined(USE_ARRAY)
6116973769bSHong Zhang   ierr = PetscFree(spdot);
6126973769bSHong Zhang #endif
6135df89d91SHong Zhang   PetscFunctionReturn(0);
6145df89d91SHong Zhang }
6155df89d91SHong Zhang 
6165df89d91SHong Zhang #undef __FUNCT__
61775648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
61875648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
6195df89d91SHong Zhang   PetscErrorCode ierr;
6205df89d91SHong Zhang 
6215df89d91SHong Zhang   PetscFunctionBegin;
6225df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
62375648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
6245df89d91SHong Zhang   }
62575648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
6265df89d91SHong Zhang   PetscFunctionReturn(0);
6275df89d91SHong Zhang }
6285df89d91SHong Zhang 
6295df89d91SHong Zhang #undef __FUNCT__
63075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
63175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
6325df89d91SHong Zhang {
633bc011b1eSHong Zhang   PetscErrorCode ierr;
634bc011b1eSHong Zhang   Mat            At;
63538baddfdSBarry Smith   PetscInt       *ati,*atj;
636bc011b1eSHong Zhang 
637bc011b1eSHong Zhang   PetscFunctionBegin;
638bc011b1eSHong Zhang   /* create symbolic At */
639bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
640d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
641bc011b1eSHong Zhang 
642bc011b1eSHong Zhang   /* get symbolic C=At*B */
643bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
644bc011b1eSHong Zhang 
645bc011b1eSHong Zhang   /* clean up */
6466bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
647bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
648bc011b1eSHong Zhang   PetscFunctionReturn(0);
649bc011b1eSHong Zhang }
650bc011b1eSHong Zhang 
651bc011b1eSHong Zhang #undef __FUNCT__
65275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
65375648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
654bc011b1eSHong Zhang {
655bc011b1eSHong Zhang   PetscErrorCode ierr;
6560fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
657d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
658d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
659d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
6600fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
661bc011b1eSHong Zhang 
662bc011b1eSHong Zhang   PetscFunctionBegin;
663bc011b1eSHong Zhang   /* clear old values in C */
664bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
665bc011b1eSHong Zhang 
666bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
667bc011b1eSHong Zhang   for (i=0;i<am;i++) {
668bc011b1eSHong Zhang     bj   = b->j + bi[i];
669bc011b1eSHong Zhang     ba   = b->a + bi[i];
670bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
671bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
672bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
673bc011b1eSHong Zhang       nextb = 0;
6740fbc74f4SHong Zhang       crow  = *aj++;
675bc011b1eSHong Zhang       cjj   = cj + ci[crow];
676bc011b1eSHong Zhang       caj   = ca + ci[crow];
677bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
678bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
6790fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
6800fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
681bc011b1eSHong Zhang           nextb++;
682bc011b1eSHong Zhang         }
683bc011b1eSHong Zhang       }
684bc011b1eSHong Zhang       flops += 2*bnzi;
6850fbc74f4SHong Zhang       aa++;
686bc011b1eSHong Zhang     }
687bc011b1eSHong Zhang   }
688bc011b1eSHong Zhang 
689bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
690bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
691bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
692bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
693bc011b1eSHong Zhang   PetscFunctionReturn(0);
694bc011b1eSHong Zhang }
6959123193aSHong Zhang 
696ec01deb9SMatthew Knepley EXTERN_C_BEGIN
6979123193aSHong Zhang #undef __FUNCT__
6989123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
6999123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
7009123193aSHong Zhang {
7019123193aSHong Zhang   PetscErrorCode ierr;
7029123193aSHong Zhang 
7039123193aSHong Zhang   PetscFunctionBegin;
7049123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
7059123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
7069123193aSHong Zhang   }
7079123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
7089123193aSHong Zhang   PetscFunctionReturn(0);
7099123193aSHong Zhang }
710ec01deb9SMatthew Knepley EXTERN_C_END
711edd81eecSMatthew Knepley 
7129123193aSHong Zhang #undef __FUNCT__
7139123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
7149123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
7159123193aSHong Zhang {
7169123193aSHong Zhang   PetscErrorCode ierr;
7179123193aSHong Zhang 
7189123193aSHong Zhang   PetscFunctionBegin;
7195a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
7208cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense;
7219123193aSHong Zhang   PetscFunctionReturn(0);
7229123193aSHong Zhang }
7239123193aSHong Zhang 
7249123193aSHong Zhang #undef __FUNCT__
7259123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
7269123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
7279123193aSHong Zhang {
728f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
7299123193aSHong Zhang   PetscErrorCode ierr;
730dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
731dd6ea824SBarry Smith   MatScalar      *aa;
732d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
733f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
7349123193aSHong Zhang 
7359123193aSHong Zhang   PetscFunctionBegin;
736f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
737e32f2f54SBarry 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);
738e32f2f54SBarry 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);
739e32f2f54SBarry 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);
740f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
741f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
742f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
743f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
744f32d5d43SBarry Smith     colam = col*am;
745f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
746f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
747f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
748f32d5d43SBarry Smith       aj  = a->j + a->i[i];
749f32d5d43SBarry Smith       aa  = a->a + a->i[i];
750f32d5d43SBarry Smith       for (j=0; j<n; j++) {
751f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
752f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
753f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
754f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
7559123193aSHong Zhang       }
756f32d5d43SBarry Smith       c[colam + i]       = r1;
757f32d5d43SBarry Smith       c[colam + am + i]  = r2;
758f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
759f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
760f32d5d43SBarry Smith     }
761f32d5d43SBarry Smith     b1 += bm4;
762f32d5d43SBarry Smith     b2 += bm4;
763f32d5d43SBarry Smith     b3 += bm4;
764f32d5d43SBarry Smith     b4 += bm4;
765f32d5d43SBarry Smith   }
766f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
767f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
768f32d5d43SBarry Smith       r1 = 0.0;
769f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
770f32d5d43SBarry Smith       aj  = a->j + a->i[i];
771f32d5d43SBarry Smith       aa  = a->a + a->i[i];
772f32d5d43SBarry Smith 
773f32d5d43SBarry Smith       for (j=0; j<n; j++) {
774f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
775f32d5d43SBarry Smith       }
776f32d5d43SBarry Smith       c[col*am + i]     = r1;
777f32d5d43SBarry Smith     }
778f32d5d43SBarry Smith     b1 += bm;
779f32d5d43SBarry Smith   }
780dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
781f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
782f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
783f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
784f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
785f32d5d43SBarry Smith   PetscFunctionReturn(0);
786f32d5d43SBarry Smith }
787f32d5d43SBarry Smith 
788f32d5d43SBarry Smith /*
7894324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
790f32d5d43SBarry Smith */
791f32d5d43SBarry Smith #undef __FUNCT__
792f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
793f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
794f32d5d43SBarry Smith {
795f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
796f32d5d43SBarry Smith   PetscErrorCode ierr;
797dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
798dd6ea824SBarry Smith   MatScalar      *aa;
799d0f46423SBarry 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;
8004324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
801f32d5d43SBarry Smith 
802f32d5d43SBarry Smith   PetscFunctionBegin;
803f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
804f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
805f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
806f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
8074324174eSBarry Smith 
8084324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
8094324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
8104324174eSBarry Smith       colam = col*am;
8114324174eSBarry Smith       arm   = a->compressedrow.nrows;
8124324174eSBarry Smith       ii    = a->compressedrow.i;
8134324174eSBarry Smith       ridx  = a->compressedrow.rindex;
8144324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
8154324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
8164324174eSBarry Smith 	n   = ii[i+1] - ii[i];
8174324174eSBarry Smith 	aj  = a->j + ii[i];
8184324174eSBarry Smith 	aa  = a->a + ii[i];
8194324174eSBarry Smith 	for (j=0; j<n; j++) {
8204324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
8214324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
8224324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
8234324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
8244324174eSBarry Smith 	}
8254324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
8264324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
8274324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
8284324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
8294324174eSBarry Smith       }
8304324174eSBarry Smith       b1 += bm4;
8314324174eSBarry Smith       b2 += bm4;
8324324174eSBarry Smith       b3 += bm4;
8334324174eSBarry Smith       b4 += bm4;
8344324174eSBarry Smith     }
8354324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
8364324174eSBarry Smith       colam = col*am;
8374324174eSBarry Smith       arm   = a->compressedrow.nrows;
8384324174eSBarry Smith       ii    = a->compressedrow.i;
8394324174eSBarry Smith       ridx  = a->compressedrow.rindex;
8404324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
8414324174eSBarry Smith 	r1 = 0.0;
8424324174eSBarry Smith 	n   = ii[i+1] - ii[i];
8434324174eSBarry Smith 	aj  = a->j + ii[i];
8444324174eSBarry Smith 	aa  = a->a + ii[i];
8454324174eSBarry Smith 
8464324174eSBarry Smith 	for (j=0; j<n; j++) {
8474324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
8484324174eSBarry Smith 	}
8494324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
8504324174eSBarry Smith       }
8514324174eSBarry Smith       b1 += bm;
8524324174eSBarry Smith     }
8534324174eSBarry Smith   } else {
854f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
855f32d5d43SBarry Smith       colam = col*am;
856f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
857f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
858f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
859f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
860f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
861f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
862f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
863f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
864f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
865f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
866f32d5d43SBarry Smith 	}
867f32d5d43SBarry Smith 	c[colam + i]       += r1;
868f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
869f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
870f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
871f32d5d43SBarry Smith       }
872f32d5d43SBarry Smith       b1 += bm4;
873f32d5d43SBarry Smith       b2 += bm4;
874f32d5d43SBarry Smith       b3 += bm4;
875f32d5d43SBarry Smith       b4 += bm4;
876f32d5d43SBarry Smith     }
877f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
878f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
879f32d5d43SBarry Smith 	r1 = 0.0;
880f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
881f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
882f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
883f32d5d43SBarry Smith 
884f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
885f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
886f32d5d43SBarry Smith 	}
887f32d5d43SBarry Smith 	c[col*am + i]     += r1;
888f32d5d43SBarry Smith       }
889f32d5d43SBarry Smith       b1 += bm;
890f32d5d43SBarry Smith     }
8914324174eSBarry Smith   }
892dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
893f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
894f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
8959123193aSHong Zhang   PetscFunctionReturn(0);
8969123193aSHong Zhang }
897b1683b59SHong Zhang 
898b1683b59SHong Zhang #undef __FUNCT__
899b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
900b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
901c8db22f6SHong Zhang {
902c8db22f6SHong Zhang   PetscErrorCode ierr;
9032f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
9042f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
9052f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
9062f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
9072f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
9082f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
909c8db22f6SHong Zhang 
910c8db22f6SHong Zhang   PetscFunctionBegin;
9112f5f1f90SHong Zhang   btval_den=btdense->v;
9122f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
9132f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
9142f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
9152f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
916d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
9172f5f1f90SHong Zhang       btcol = bj + bi[col];
9182f5f1f90SHong Zhang       btval = ba + bi[col];
9192f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
920d2853540SHong Zhang       for (j=0; j<anz; j++){
9212f5f1f90SHong Zhang         brow            = btcol[j];
9222f5f1f90SHong Zhang         btval_den[brow] = btval[j];
923c8db22f6SHong Zhang       }
924c8db22f6SHong Zhang     }
9252f5f1f90SHong Zhang     btval_den += m;
926c8db22f6SHong Zhang   }
927c8db22f6SHong Zhang   PetscFunctionReturn(0);
928c8db22f6SHong Zhang }
929c8db22f6SHong Zhang 
930c8db22f6SHong Zhang #undef __FUNCT__
931b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
932b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
9338972f759SHong Zhang {
9348972f759SHong Zhang   PetscErrorCode ierr;
935b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
9362f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
937b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
9382f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
9398972f759SHong Zhang 
9408972f759SHong Zhang   PetscFunctionBegin;
9418972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
942b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
943b2d2b04fSHong Zhang   cp_den = ca_den;
9442f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
9452f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
9462f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
9472f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
9482f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
9492f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
950b2d2b04fSHong Zhang     }
951b2d2b04fSHong Zhang     cp_den += m;
952b2d2b04fSHong Zhang   }
953b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
9548972f759SHong Zhang   PetscFunctionReturn(0);
9558972f759SHong Zhang }
9568972f759SHong Zhang 
957*e99be685SHong Zhang /*
958*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
959*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
960*e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
961*e99be685SHong Zhang  */
962*e99be685SHong Zhang #undef __FUNCT__
963*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
964*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
965*e99be685SHong Zhang {
966*e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
967*e99be685SHong Zhang   PetscErrorCode ierr;
968*e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
969*e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
970*e99be685SHong Zhang   PetscInt       *cspidx;
971*e99be685SHong Zhang 
972*e99be685SHong Zhang   PetscFunctionBegin;
973*e99be685SHong Zhang   *nn = n;
974*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
975*e99be685SHong Zhang   if (symmetric) {
976*e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
977*e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
978*e99be685SHong Zhang   } else {
979*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
980*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
981*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
982*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
983*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
984*e99be685SHong Zhang     jj = a->j;
985*e99be685SHong Zhang     for (i=0; i<nz; i++) {
986*e99be685SHong Zhang       collengths[jj[i]]++;
987*e99be685SHong Zhang     }
988*e99be685SHong Zhang     cia[0] = oshift;
989*e99be685SHong Zhang     for (i=0; i<n; i++) {
990*e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
991*e99be685SHong Zhang     }
992*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
993*e99be685SHong Zhang     jj   = a->j;
994*e99be685SHong Zhang     for (row=0; row<m; row++) {
995*e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
996*e99be685SHong Zhang       for (i=0; i<mr; i++) {
997*e99be685SHong Zhang         col = *jj++;
998*e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
999*e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1000*e99be685SHong Zhang       }
1001*e99be685SHong Zhang     }
1002*e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
1003*e99be685SHong Zhang     *ia = cia; *ja = cja;
1004*e99be685SHong Zhang     *spidx = cspidx;
1005*e99be685SHong Zhang   }
1006*e99be685SHong Zhang   PetscFunctionReturn(0);
1007*e99be685SHong Zhang }
1008*e99be685SHong Zhang 
1009*e99be685SHong Zhang #undef __FUNCT__
1010*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
1011*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1012*e99be685SHong Zhang {
1013*e99be685SHong Zhang   PetscErrorCode ierr;
1014*e99be685SHong Zhang 
1015*e99be685SHong Zhang   PetscFunctionBegin;
1016*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1017*e99be685SHong Zhang 
1018*e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1019*e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1020*e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1021*e99be685SHong Zhang   PetscFunctionReturn(0);
1022*e99be685SHong Zhang }
1023*e99be685SHong Zhang 
10248972f759SHong Zhang #undef __FUNCT__
1025b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1026b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1027b1683b59SHong Zhang {
1028b1683b59SHong Zhang   PetscErrorCode ierr;
1029d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
1030b1683b59SHong Zhang   const PetscInt *is;
1031b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1032b1683b59SHong Zhang   IS             *isa;
1033b1683b59SHong Zhang   PetscBool      done;
1034b1683b59SHong Zhang   PetscBool      flg1,flg2;
1035b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1036*e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1037d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1038b1683b59SHong Zhang 
1039b1683b59SHong Zhang   PetscFunctionBegin;
1040b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1041*e99be685SHong Zhang 
1042b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1043b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1044b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1045b1683b59SHong Zhang   if (flg1 || flg2) {
1046b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1047b1683b59SHong Zhang   }
1048b1683b59SHong Zhang 
1049b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1050b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1051b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1052b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1053b1683b59SHong Zhang   c->rstart = 0;
1054b1683b59SHong Zhang 
1055b1683b59SHong Zhang   c->ncolors = nis;
1056b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1057b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1058d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1059d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1060d2853540SHong Zhang   colorforrow[0]    = 0;
1061d2853540SHong Zhang   rows_i            = rows;
1062d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1063b1683b59SHong Zhang 
1064d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1065d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1066d2853540SHong Zhang   colorforcol[0] = 0;
1067d2853540SHong Zhang   columns_i      = columns;
1068d2853540SHong Zhang 
1069*e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1070b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1071b1683b59SHong Zhang 
1072b28d80bdSHong Zhang   cm = c->m;
1073b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1074*e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1075b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1076b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1077b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1078b1683b59SHong Zhang     c->ncolumns[i] = n;
1079b1683b59SHong Zhang     if (n) {
1080d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1081b1683b59SHong Zhang     }
1082d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1083d2853540SHong Zhang     columns_i       += n;
1084b1683b59SHong Zhang 
1085b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1086e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1087*e99be685SHong Zhang 
1088b1683b59SHong Zhang     /* loop over columns*/
1089b1683b59SHong Zhang     for (j=0; j<n; j++) {
1090b1683b59SHong Zhang       col     = is[j];
1091d2853540SHong Zhang       row_idx = cj + ci[col];
1092b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1093b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1094b1683b59SHong Zhang       for (k=0; k<m; k++) {
1095*e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1096d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1097b1683b59SHong Zhang       }
1098b1683b59SHong Zhang     }
1099b1683b59SHong Zhang     /* count the number of hits */
1100b1683b59SHong Zhang     nrows = 0;
1101e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1102b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1103b1683b59SHong Zhang     }
1104b1683b59SHong Zhang     c->nrows[i]      = nrows;
1105d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1106d2853540SHong Zhang 
1107b1683b59SHong Zhang     nrows       = 0;
1108e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1109b1683b59SHong Zhang       if (rowhit[j]) {
1110d2853540SHong Zhang         rows_i[nrows]            = j;
1111*e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1112b1683b59SHong Zhang         nrows++;
1113b1683b59SHong Zhang       }
1114b1683b59SHong Zhang     }
1115b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1116d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1117b1683b59SHong Zhang   }
1118*e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1119b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1120b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1121d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1122d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1123d2853540SHong Zhang #endif
1124b28d80bdSHong Zhang 
1125d2853540SHong Zhang   c->colorforrow     = colorforrow;
1126d2853540SHong Zhang   c->rows            = rows;
1127d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1128d2853540SHong Zhang   c->colorforcol     = colorforcol;
1129d2853540SHong Zhang   c->columns         = columns;
1130f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1131b1683b59SHong Zhang   PetscFunctionReturn(0);
1132b1683b59SHong Zhang }
1133